FE-Project
Loading...
Searching...
No Matches
scale_atm_dyn_dgm_nonhydro3d_rhot_hevi Module Reference

module FElib / Fluid dyn solver / Atmosphere / Regional nonhydrostatic model / HEVI More...

Functions/Subroutines

subroutine, public atm_dyn_dgm_nonhydro3d_rhot_hevi_init (mesh)
subroutine, public atm_dyn_dgm_nonhydro3d_rhot_hevi_final ()
subroutine, public atm_dyn_dgm_nonhydro3d_rhot_hevi_cal_tend_asis (dens_dt, momx_dt, momy_dt, momz_dt, rhot_dt, ddens_, momx_, momy_, momz_, drhot_, dpres_, dens_hyd, pres_hyd, pres_hyd_ref, therm_hyd, coriolis, rtot, cvtot, cptot, dphyddx, dphyddy, element3d_operation, dx, dy, dz, sx, sy, sz, lift, lmesh, elem, lmesh2d, elem2d)
subroutine, public atm_dyn_dgm_nonhydro3d_rhot_hevi_cal_tend (dens_dt, momx_dt, momy_dt, momz_dt, rhot_dt, ddens_, momx_, momy_, momz_, drhot_, dpres_, dens_hyd, pres_hyd, pres_hyd_ref, therm_hyd, coriolis, rtot, cvtot, cptot, dphyddx, dphyddy, element3d_operation, dx, dy, dz, sx, sy, sz, lift, lmesh, elem, lmesh2d, elem2d)
subroutine, public atm_dyn_dgm_nonhydro3d_rhot_hevi_cal_vi_asis (dens_dt, momx_dt, momy_dt, momz_dt, rhot_dt, ddens_, momx_, momy_, momz_, drhot_, dens_hyd, pres_hyd, ddens0_, momx0_, momy0_, momz0_, drhot0_, rtot, cvtot, cptot, element3d_operation, dz, lift, impl_fac, dt, lmesh, elem, lmesh2d, elem2d)
subroutine, public atm_dyn_dgm_nonhydro3d_rhot_hevi_cal_vi (dens_dt, momx_dt, momy_dt, momz_dt, rhot_dt, ddens_, momx_, momy_, momz_, drhot_, dens_hyd, pres_hyd, ddens0_, momx0_, momy0_, momz0_, drhot0_, rtot, cvtot, cptot, element3d_operation, dz, lift, impl_fac, dt, lmesh, elem, lmesh2d, elem2d)

Detailed Description

module FElib / Fluid dyn solver / Atmosphere / Regional nonhydrostatic model / HEVI

Description
HEVI DGM scheme for Atmospheric dynamical process The governing equations is a fully compressible nonhydrostatic equations, which consist of mass, momentum, and thermodynamics (density * potential temperature conservation) equations.
Author
Yuta Kawai, Team SCALE

Function/Subroutine Documentation

◆ atm_dyn_dgm_nonhydro3d_rhot_hevi_init()

subroutine, public scale_atm_dyn_dgm_nonhydro3d_rhot_hevi::atm_dyn_dgm_nonhydro3d_rhot_hevi_init ( class(meshbase3d), intent(in) mesh)

Definition at line 75 of file scale_atm_dyn_dgm_nonhydro3d_rhot_hevi.F90.

78 implicit none
79 class(MeshBase3D), intent(in) :: mesh
80 !--------------------------------------------
81
82 call atm_dyn_dgm_nonhydro3d_common_init( mesh )
84
85 return
module FElib / Fluid dyn solver / Atmosphere / Nonhydrostatic model / HEVI / Common

References scale_atm_dyn_dgm_nonhydro3d_common::atm_dyn_dgm_nonhydro3d_common_init(), and scale_atm_dyn_dgm_nonhydro3d_rhot_hevi_common::atm_dyn_dgm_nonhydro3d_rhot_hevi_common_init().

◆ atm_dyn_dgm_nonhydro3d_rhot_hevi_final()

subroutine, public scale_atm_dyn_dgm_nonhydro3d_rhot_hevi::atm_dyn_dgm_nonhydro3d_rhot_hevi_final

◆ atm_dyn_dgm_nonhydro3d_rhot_hevi_cal_tend_asis()

subroutine, public scale_atm_dyn_dgm_nonhydro3d_rhot_hevi::atm_dyn_dgm_nonhydro3d_rhot_hevi_cal_tend_asis ( real(rp), dimension(elem%np,lmesh%nea), intent(out) dens_dt,
real(rp), dimension(elem%np,lmesh%nea), intent(out) momx_dt,
real(rp), dimension(elem%np,lmesh%nea), intent(out) momy_dt,
real(rp), dimension(elem%np,lmesh%nea), intent(out) momz_dt,
real(rp), dimension(elem%np,lmesh%nea), intent(out) rhot_dt,
real(rp), dimension(elem%np,lmesh%nea), intent(in) ddens_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) momx_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) momy_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) momz_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) drhot_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) dpres_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) dens_hyd,
real(rp), dimension(elem%np,lmesh%nea), intent(in) pres_hyd,
real(rp), dimension(elem%np,lmesh%nea), intent(in) pres_hyd_ref,
real(rp), dimension(elem%np,lmesh%nea), intent(in) therm_hyd,
real(rp), dimension(elem2d%np,lmesh2d%nea), intent(in) coriolis,
real(rp), dimension(elem%np,lmesh%nea), intent(in) rtot,
real(rp), dimension(elem%np,lmesh%nea), intent(in) cvtot,
real(rp), dimension(elem%np,lmesh%nea), intent(in) cptot,
real(rp), dimension(elem%np,lmesh%nea), intent(in) dphyddx,
real(rp), dimension(elem%np,lmesh%nea), intent(in) dphyddy,
class(elementoperationbase3d), intent(in) element3d_operation,
type(sparsemat), intent(in) dx,
type(sparsemat), intent(in) dy,
type(sparsemat), intent(in) dz,
type(sparsemat), intent(in) sx,
type(sparsemat), intent(in) sy,
type(sparsemat), intent(in) sz,
type(sparsemat), intent(in) lift,
class(localmesh3d), intent(in) lmesh,
class(elementbase3d), intent(in) elem,
class(localmesh2d), intent(in) lmesh2d,
class(elementbase2d), intent(in) elem2d )

Definition at line 104 of file scale_atm_dyn_dgm_nonhydro3d_rhot_hevi.F90.

111
114
115 implicit none
116
117 class(LocalMesh3D), intent(in) :: lmesh
118 class(ElementBase3D), intent(in) :: elem
119 class(LocalMesh2D), intent(in) :: lmesh2D
120 class(ElementBase2D), intent(in) :: elem2D
121 class(ElementOperationBase3D), intent(in) :: element3D_operation
122 type(SparseMat), intent(in) :: Dx, Dy, Dz, Sx, Sy, Sz, Lift
123 real(RP), intent(out) :: DENS_dt(elem%Np,lmesh%NeA)
124 real(RP), intent(out) :: MOMX_dt(elem%Np,lmesh%NeA)
125 real(RP), intent(out) :: MOMY_dt(elem%Np,lmesh%NeA)
126 real(RP), intent(out) :: MOMZ_dt(elem%Np,lmesh%NeA)
127 real(RP), intent(out) :: RHOT_dt(elem%Np,lmesh%NeA)
128 real(RP), intent(in) :: DDENS_(elem%Np,lmesh%NeA)
129 real(RP), intent(in) :: MOMX_(elem%Np,lmesh%NeA)
130 real(RP), intent(in) :: MOMY_(elem%Np,lmesh%NeA)
131 real(RP), intent(in) :: MOMZ_(elem%Np,lmesh%NeA)
132 real(RP), intent(in) :: DRHOT_(elem%Np,lmesh%NeA)
133 real(RP), intent(in) :: DPRES_(elem%Np,lmesh%NeA)
134 real(RP), intent(in) :: DENS_hyd(elem%Np,lmesh%NeA)
135 real(RP), intent(in) :: PRES_hyd(elem%Np,lmesh%NeA)
136 real(RP), intent(in) :: PRES_hyd_ref(elem%Np,lmesh%NeA)
137 real(RP), intent(in) :: THERM_hyd(elem%Np,lmesh%NeA)
138 real(RP), intent(in) :: CORIOLIS(elem2D%Np,lmesh2D%NeA)
139 real(RP), intent(in) :: Rtot(elem%Np,lmesh%NeA)
140 real(RP), intent(in) :: CVtot(elem%Np,lmesh%NeA)
141 real(RP), intent(in) :: CPtot(elem%Np,lmesh%NeA)
142 real(RP), intent(in) :: DPhydDx(elem%Np,lmesh%NeA)
143 real(RP), intent(in) :: DPhydDy(elem%Np,lmesh%NeA)
144
145 real(RP) :: Fx(elem%Np), Fy(elem%Np), Fz(elem%Np), LiftDelFlx(elem%Np)
146 real(RP) :: DPRES_hyd(elem%Np), GradPhyd_x(elem%Np), GradPhyd_y(elem%Np)
147 real(RP) :: del_flux(elem%NfpTot,lmesh%Ne,PRGVAR_NUM)
148 real(RP) :: del_flux_hyd(elem%NfpTot,lmesh%Ne,2)
149 real(RP) :: RHOT_(elem%Np)
150 real(RP) :: rdens_(elem%Np), u_(elem%Np), v_(elem%Np), w_(elem%Np), wt_(elem%Np)
151 real(RP) :: GsqrtV(elem%Np), RGsqrtV(elem%Np)
152 real(RP) :: Cori(elem%Np)
153
154 integer :: ke, ke2d
155
156 real(RP) :: gamm, rgamm
157 real(RP) :: rP0
158 real(RP) :: RovP0, P0ovR
159 !------------------------------------------------------------------------
160
161 call prof_rapstart( 'cal_dyn_tend_bndflux', 3)
162 call get_ebnd_flux( &
163 del_flux, del_flux_hyd, & ! (out)
164 ddens_, momx_, momy_, momz_, drhot_, dpres_, dens_hyd, pres_hyd, & ! (in)
165 rtot, cvtot, cptot, & ! (in)
166 lmesh%Gsqrt, lmesh%GI3(:,:,1), lmesh%GI3(:,:,2), & ! (in)
167 lmesh%normal_fn(:,:,1), lmesh%normal_fn(:,:,2), lmesh%normal_fn(:,:,3), & ! (in)
168 lmesh%vmapM, lmesh%vmapP, & ! (in)
169 lmesh, elem, lmesh2d, elem2d ) ! (in)
170 call prof_rapend( 'cal_dyn_tend_bndflux', 3)
171
172 !-----
173 call prof_rapstart( 'cal_dyn_tend_interior', 3)
174 gamm = cpdry / cvdry
175 rgamm = cvdry / cpdry
176 rp0 = 1.0_rp / pres00
177 rovp0 = rdry * rp0
178 p0ovr = pres00 / rdry
179
180 !$omp parallel do private( ke2d, Cori, &
181 !$omp RHOT_, rdens_, u_, v_, w_, wt_, &
182 !$omp DPRES_hyd, GradPhyd_x, GradPhyd_y, &
183 !$omp GsqrtV, RGsqrtV, &
184 !$omp Fx, Fy, Fz, LiftDelFlx )
185 do ke = lmesh%NeS, lmesh%NeE
186 !--
187 ke2d = lmesh%EMap3Dto2D(ke)
188 cori(:) = coriolis(elem%IndexH2Dto3D(:),ke2d)
189
190 gsqrtv(:) = lmesh%Gsqrt(:,ke) / lmesh%GsqrtH(elem%IndexH2Dto3D,ke2d)
191 rgsqrtv(:) = 1.0_rp / gsqrtv(:)
192
193 !--
194 rhot_(:) = p0ovr * (pres_hyd(:,ke) * rp0)**rgamm + drhot_(:,ke)
195 ! DPRES_(:) = PRES00 * ( Rtot(:,ke) * rP0 * RHOT_(:) )**( CPtot(:,ke) / CVtot(:,ke) ) &
196 ! - PRES_hyd(:,ke)
197
198 rdens_(:) = 1.0_rp / (ddens_(:,ke) + dens_hyd(:,ke))
199 u_(:) = momx_(:,ke) * rdens_(:)
200 v_(:) = momy_(:,ke) * rdens_(:)
201 w_(:) = momz_(:,ke) * rdens_(:)
202 wt_(:) = w_(:) * rgsqrtv(:) + lmesh%GI3(:,ke,1) * u_(:) + lmesh%GI3(:,ke,2) * v_(:)
203
204 !-- Gradient hydrostatic pressure
205
206 dpres_hyd(:) = pres_hyd(:,ke) - pres_hyd_ref(:,ke)
207
208 call sparsemat_matmul(dx, gsqrtv(:) * dpres_hyd(:), fx)
209 call sparsemat_matmul(dz, gsqrtv(:) * lmesh%GI3(:,ke,1) * dpres_hyd(:), fz)
210 call sparsemat_matmul(lift, lmesh%Fscale(:,ke) * del_flux_hyd(:,ke,1), liftdelflx)
211 gradphyd_x(:) = lmesh%Escale(:,ke,1,1) * fx(:) &
212 + lmesh%Escale(:,ke,3,3) * fz(:) &
213 + liftdelflx(:)
214
215 call sparsemat_matmul(dy, gsqrtv(:) * dpres_hyd(:), fy)
216 call sparsemat_matmul(dz, gsqrtv(:) * lmesh%GI3(:,ke,2) * dpres_hyd(:), fz)
217 call sparsemat_matmul(lift, lmesh%Fscale(:,ke) * del_flux_hyd(:,ke,2), liftdelflx)
218 gradphyd_y(:) = lmesh%Escale(:,ke,2,2) * fy(:) &
219 + lmesh%Escale(:,ke,3,3) * fz(:) &
220 + liftdelflx(:)
221 !-- DENS
222 call sparsemat_matmul(dx, lmesh%Gsqrt(:,ke) * momx_(:,ke), fx)
223 call sparsemat_matmul(dy, lmesh%Gsqrt(:,ke) * momy_(:,ke), fy)
224 call sparsemat_matmul(lift, lmesh%Fscale(:,ke) * del_flux(:,ke,dens_vid), liftdelflx)
225
226 dens_dt(:,ke) = - ( &
227 lmesh%Escale(:,ke,1,1) * fx(:) &
228 + lmesh%Escale(:,ke,2,2) * fy(:) &
229 + liftdelflx(:) ) / lmesh%Gsqrt(:,ke)
230
231 !-- MOMX
232 call sparsemat_matmul(dx, lmesh%Gsqrt(:,ke) * ( u_(:) * momx_(:,ke) + dpres_(:,ke) ), fx)
233 call sparsemat_matmul(dy, lmesh%Gsqrt(:,ke) * v_(:) * momx_(:,ke) , fy)
234 call sparsemat_matmul(dz, lmesh%Gsqrt(:,ke) * ( wt_(:) * momx_(:,ke) &
235 + lmesh%GI3(:,ke,1) * dpres_(:,ke) ), fz)
236 call sparsemat_matmul(lift, lmesh%Fscale(:,ke) * del_flux(:,ke,momx_vid), liftdelflx)
237
238 momx_dt(:,ke) = &
239 - ( lmesh%Escale(:,ke,1,1) * fx(:) &
240 + lmesh%Escale(:,ke,2,2) * fy(:) &
241 + lmesh%Escale(:,ke,3,3) * fz(:) &
242 + liftdelflx(:) ) / lmesh%Gsqrt(:,ke) &
243 - gradphyd_x(:) * rgsqrtv(:) &
244 + cori(:) * momy_(:,ke)
245
246 !-- MOMY
247 call sparsemat_matmul(dx, lmesh%Gsqrt(:,ke) * u_(:) * momy_(:,ke) , fx)
248 call sparsemat_matmul(dy, lmesh%Gsqrt(:,ke) * ( v_(:) * momy_(:,ke) + dpres_(:,ke) ), fy)
249 call sparsemat_matmul(dz, lmesh%Gsqrt(:,ke) * ( wt_(:) * momy_(:,ke) &
250 + lmesh%GI3(:,ke,2) * dpres_(:,ke) ), fz)
251 call sparsemat_matmul(lift, lmesh%Fscale(:,ke) * del_flux(:,ke,momy_vid), liftdelflx)
252
253 momy_dt(:,ke) = &
254 - ( lmesh%Escale(:,ke,1,1) * fx(:) &
255 + lmesh%Escale(:,ke,2,2) * fy(:) &
256 + lmesh%Escale(:,ke,3,3) * fz(:) &
257 + liftdelflx(:) ) / lmesh%Gsqrt(:,ke) &
258 - gradphyd_y(:) * rgsqrtv(:) &
259 - cori(:) * momx_(:,ke)
260
261 !-- MOMZ
262 call sparsemat_matmul(dx, lmesh%Gsqrt(:,ke) * u_(:) * momz_(:,ke), fx)
263 call sparsemat_matmul(dy, lmesh%Gsqrt(:,ke) * v_(:) * momz_(:,ke), fy)
264 call sparsemat_matmul(dz, lmesh%Gsqrt(:,ke) * wt_(:) * momz_(:,ke), fz)
265 call sparsemat_matmul(lift, lmesh%Fscale(:,ke) * del_flux(:,ke,momz_vid), liftdelflx)
266
267 momz_dt(:,ke) = &
268 - ( lmesh%Escale(:,ke,1,1) * fx(:) &
269 + lmesh%Escale(:,ke,2,2) * fy(:) &
270 + lmesh%Escale(:,ke,3,3) * fz(:) &
271 + liftdelflx(:) ) / lmesh%Gsqrt(:,ke)
272
273 !-- RHOT
274 call sparsemat_matmul(dx, lmesh%Gsqrt(:,ke) * u_(:) * rhot_(:), fx)
275 call sparsemat_matmul(dy, lmesh%Gsqrt(:,ke) * v_(:) * rhot_(:), fy)
276 call sparsemat_matmul(lift, lmesh%Fscale(:,ke) * del_flux(:,ke,rhot_vid), liftdelflx)
277
278 rhot_dt(:,ke) = &
279 - ( lmesh%Escale(:,ke,1,1) * fx(:) &
280 + lmesh%Escale(:,ke,2,2) * fy(:) &
281 + liftdelflx(:) ) / lmesh%Gsqrt(:,ke)
282 end do
283 call prof_rapend( 'cal_dyn_tend_interior', 3)
284
285 return
module FElib / Fluid dyn solver / Atmosphere / Nonhydrostatic model / HEVI / Numflux
subroutine, public atm_dyn_dgm_nonhydro3d_rhot_hevi_numflux_get_generalvc_asis(del_flux, del_flux_hyd, ddens_, momx_, momy_, momz_, drhot_, dpres_, dens_hyd, pres_hyd, rtot, cvtot, cptot, gsqrt, g13, g23, nx, ny, nz, vmapm, vmapp, lmesh, elem, lmesh2d, elem2d)

References scale_atm_dyn_dgm_nonhydro3d_rhot_hevi_numflux::atm_dyn_dgm_nonhydro3d_rhot_hevi_numflux_get_generalvc_asis().

◆ atm_dyn_dgm_nonhydro3d_rhot_hevi_cal_tend()

subroutine, public scale_atm_dyn_dgm_nonhydro3d_rhot_hevi::atm_dyn_dgm_nonhydro3d_rhot_hevi_cal_tend ( real(rp), dimension(elem%np,lmesh%nea), intent(out) dens_dt,
real(rp), dimension(elem%np,lmesh%nea), intent(out) momx_dt,
real(rp), dimension(elem%np,lmesh%nea), intent(out) momy_dt,
real(rp), dimension(elem%np,lmesh%nea), intent(out) momz_dt,
real(rp), dimension(elem%np,lmesh%nea), intent(out) rhot_dt,
real(rp), dimension(elem%np,lmesh%nea), intent(in) ddens_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) momx_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) momy_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) momz_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) drhot_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) dpres_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) dens_hyd,
real(rp), dimension(elem%np,lmesh%nea), intent(in) pres_hyd,
real(rp), dimension(elem%np,lmesh%nea), intent(in) pres_hyd_ref,
real(rp), dimension(elem%np,lmesh%nea), intent(in) therm_hyd,
real(rp), dimension(elem2d%np,lmesh2d%nea), intent(in) coriolis,
real(rp), dimension (elem%np,lmesh%nea), intent(in) rtot,
real(rp), dimension(elem%np,lmesh%nea), intent(in) cvtot,
real(rp), dimension(elem%np,lmesh%nea), intent(in) cptot,
real(rp), dimension(elem%np,lmesh%nea), intent(in) dphyddx,
real(rp), dimension(elem%np,lmesh%nea), intent(in) dphyddy,
class(elementoperationbase3d), intent(in) element3d_operation,
type(sparsemat), intent(in) dx,
type(sparsemat), intent(in) dy,
type(sparsemat), intent(in) dz,
type(sparsemat), intent(in) sx,
type(sparsemat), intent(in) sy,
type(sparsemat), intent(in) sz,
type(sparsemat), intent(in) lift,
class(localmesh3d), intent(in) lmesh,
class(elementbase3d), intent(in) elem,
class(localmesh2d), intent(in) lmesh2d,
class(elementbase2d), intent(in) elem2d )

Definition at line 289 of file scale_atm_dyn_dgm_nonhydro3d_rhot_hevi.F90.

296
299 use scale_const, only: &
300 ohm => const_ohm
301 implicit none
302
303 class(LocalMesh3D), intent(in) :: lmesh
304 class(ElementBase3D), intent(in) :: elem
305 class(LocalMesh2D), intent(in) :: lmesh2D
306 class(ElementBase2D), intent(in) :: elem2D
307 class(ElementOperationBase3D), intent(in) :: element3D_operation
308 type(SparseMat), intent(in) :: Dx, Dy, Dz, Sx, Sy, Sz, Lift
309 real(RP), intent(out) :: DENS_dt(elem%Np,lmesh%NeA)
310 real(RP), intent(out) :: MOMX_dt(elem%Np,lmesh%NeA)
311 real(RP), intent(out) :: MOMY_dt(elem%Np,lmesh%NeA)
312 real(RP), intent(out) :: MOMZ_dt(elem%Np,lmesh%NeA)
313 real(RP), intent(out) :: RHOT_dt(elem%Np,lmesh%NeA)
314 real(RP), intent(in) :: DDENS_(elem%Np,lmesh%NeA)
315 real(RP), intent(in) :: MOMX_(elem%Np,lmesh%NeA)
316 real(RP), intent(in) :: MOMY_(elem%Np,lmesh%NeA)
317 real(RP), intent(in) :: MOMZ_(elem%Np,lmesh%NeA)
318 real(RP), intent(in) :: DRHOT_(elem%Np,lmesh%NeA)
319 real(RP), intent(in) :: DPRES_(elem%Np,lmesh%NeA)
320 real(RP), intent(in) :: DENS_hyd(elem%Np,lmesh%NeA)
321 real(RP), intent(in) :: PRES_hyd(elem%Np,lmesh%NeA)
322 real(RP), intent(in) :: PRES_hyd_ref(elem%Np,lmesh%NeA)
323 real(RP), intent(in) :: THERM_hyd(elem%Np,lmesh%NeA)
324 real(RP), intent(in) :: CORIOLIS(elem2D%Np,lmesh2D%NeA)
325 real(RP), intent(in) :: Rtot (elem%Np,lmesh%NeA)
326 real(RP), intent(in) :: CVtot(elem%Np,lmesh%NeA)
327 real(RP), intent(in) :: CPtot(elem%Np,lmesh%NeA)
328 real(RP), intent(in) :: DPhydDx(elem%Np,lmesh%NeA)
329 real(RP), intent(in) :: DPhydDy(elem%Np,lmesh%NeA)
330
331 real(RP) :: Flux(elem%Np,3,5), DFlux(elem%Np,4,5)
332 real(RP) :: del_flux(elem%NfpTot,PRGVAR_NUM,lmesh%Ne)
333 real(RP) :: u_, v_, w_, pt_
334 real(RP) :: RHOT_(elem%Np)
335 real(RP) :: RDENS_(elem%Np), GsqrtV(elem%Np), RGsqrtV(elem%Np), RGsqrt(elem%Np)
336 real(RP) :: Gsqrt_, GsqrtDPRES_, E11, E22, E33
337 real(RP) :: cor
338
339 integer :: ke, ke2d
340 integer :: p
341
342 real(RP) :: gamm, rgamm
343 real(RP) :: rP0
344 real(RP) :: RovP0, P0ovR
345 !------------------------------------------------------------------------
346
347 call prof_rapstart('cal_dyn_tend_bndflux', 3)
348 call get_ebnd_flux( &
349 del_flux, & ! (out)
350 ddens_, momx_, momy_, momz_, drhot_, dpres_, & ! (in)
351 dens_hyd, pres_hyd, therm_hyd, rtot, cvtot, cptot, & ! (in)
352 lmesh%Gsqrt, lmesh%GI3(:,:,1), lmesh%GI3(:,:,2), & ! (in)
353 lmesh%normal_fn(:,:,1), lmesh%normal_fn(:,:,2), lmesh%normal_fn(:,:,3), & ! (in)
354 lmesh%vmapM, lmesh%vmapP, & ! (in)
355 lmesh, elem, lmesh2d, elem2d ) ! (in)
356 call prof_rapend('cal_dyn_tend_bndflux', 3)
357
358 !-----
359 call prof_rapstart('cal_dyn_tend_interior', 3)
360 gamm = cpdry / cvdry
361 rgamm = cvdry / cpdry
362 rp0 = 1.0_rp / pres00
363 rovp0 = rdry * rp0
364 p0ovr = pres00 / rdry
365
366 !$omp parallel private( ke, ke2d, p, &
367 !$omp u_, v_, w_, pt_, &
368 !$omp cor, GsqrtDPRES_, &
369 !$omp RDENS_, RGsqrt, GsqrtV, RGsqrtV, Gsqrt_, &
370 !$omp E11, E22, E33, DFlux, Flux )
371
372 !$omp do
373 do ke = lmesh%NeS, lmesh%NeE
374 !--
375 ke2d = lmesh%EMap3Dto2D(ke)
376
377 do p=1, elem%Np
378 gsqrtv(p) = lmesh%Gsqrt(p,ke) / lmesh%GsqrtH(elem%IndexH2Dto3D(p),ke2d)
379 rgsqrtv(p) = 1.0_rp / gsqrtv(p)
380 rgsqrt(p) = 1.0_rp / lmesh%Gsqrt(p,ke)
381 rdens_(p) = 1.0_rp / ( ddens_(p,ke) + dens_hyd(p,ke) )
382 end do
383
384 !--
385 do p=1, elem%Np
386 gsqrt_ = lmesh%Gsqrt(p,ke)
387 flux(p,1,dens_vid) = gsqrt_ * momx_(p,ke)
388 flux(p,2,dens_vid) = gsqrt_ * momy_(p,ke)
389 flux(p,3,dens_vid) = gsqrt_ * ( &
390 momz_(p,ke) * rgsqrtv(p) &
391 + lmesh%GI3(p,ke,1) * momx_(p,ke) &
392 + lmesh%GI3(p,ke,2) * momy_(p,ke) )
393 end do
394 do p=1, elem%Np
395 pt_ = ( therm_hyd(p,ke) + drhot_(p,ke) ) * rdens_(p)
396 flux(p,1,rhot_vid) = flux(p,1,dens_vid) * pt_
397 flux(p,2,rhot_vid) = flux(p,2,dens_vid) * pt_
398! Flux(p,3,RHOT_VID) = 0.0_RP ! for HEVI
399
400 w_ = momz_(p,ke) * rdens_(p)
401 flux(p,1,momz_vid) = flux(p,1,dens_vid) * w_
402 flux(p,2,momz_vid) = flux(p,2,dens_vid) * w_
403 flux(p,3,momz_vid) = flux(p,3,dens_vid) * w_
404 end do
405
406 do p=1, elem%Np
407 gsqrtdpres_ = lmesh%Gsqrt(p,ke) * dpres_(p,ke)
408 u_ = momx_(p,ke) * rdens_(p)
409 v_ = momy_(p,ke) * rdens_(p)
410
411 flux(p,1,momx_vid) = flux(p,1,dens_vid) * u_ + gsqrtdpres_
412 flux(p,2,momx_vid) = flux(p,2,dens_vid) * u_
413 flux(p,3,momx_vid) = flux(p,3,dens_vid) * u_ + gsqrtdpres_ * lmesh%GI3(p,ke,1)
414
415 flux(p,1,momy_vid) = flux(p,1,dens_vid) * v_
416 flux(p,2,momy_vid) = flux(p,2,dens_vid) * v_ + gsqrtdpres_
417 flux(p,3,momy_vid) = flux(p,3,dens_vid) * v_ + gsqrtdpres_ * lmesh%GI3(p,ke,2)
418 end do
419
420 call element3d_operation%Div_var5( &
421 flux, del_flux(:,:,ke), & ! (in)
422 dflux ) ! (out)
423
424 do p=1, elem%Np
425 e11 = lmesh%Escale(p,ke,1,1)
426 e22 = lmesh%Escale(p,ke,2,2)
427 e33 = lmesh%Escale(p,ke,3,3)
428
429 dens_dt(p,ke) = - ( &
430 e11 * dflux(p,1,dens_vid) &
431 + e22 * dflux(p,2,dens_vid) &
432 + dflux(p,4,dens_vid) ) * rgsqrt(p)
433
434 rhot_dt(p,ke) = - ( &
435 e11 * dflux(p,1,rhot_vid) &
436 + e22 * dflux(p,2,rhot_vid) &
437 + dflux(p,4,rhot_vid) ) * rgsqrt(p)
438 end do
439
440 do p=1, elem%Np
441 e11 = lmesh%Escale(p,ke,1,1)
442 e22 = lmesh%Escale(p,ke,2,2)
443 e33 = lmesh%Escale(p,ke,3,3)
444
445 momz_dt(p,ke) = - ( &
446 e11 * dflux(p,1,momz_vid) &
447 + e22 * dflux(p,2,momz_vid) &
448 + e33 * dflux(p,3,momz_vid) &
449 + dflux(p,4,momz_vid) ) * rgsqrt(p)
450 end do
451
452 !--
453 do p=1, elem%Np
454 cor = coriolis(elem%IndexH2Dto3D(p),ke2d)
455 momx_dt(p,ke) = - dphyddx(p,ke) &
456 + cor * momy_(p,ke)
457 momy_dt(p,ke) = - dphyddy(p,ke) &
458 - cor * momx_(p,ke)
459 end do
460
461 do p=1, elem%Np
462 e11 = lmesh%Escale(p,ke,1,1)
463 e22 = lmesh%Escale(p,ke,2,2)
464 e33 = lmesh%Escale(p,ke,3,3)
465
466 momx_dt(p,ke) = momx_dt(p,ke) - ( &
467 e11 * dflux(p,1,momx_vid) &
468 + e22 * dflux(p,2,momx_vid) &
469 + e33 * dflux(p,3,momx_vid) &
470 + dflux(p,4,momx_vid) ) * rgsqrt(p)
471
472 momy_dt(p,ke) = momy_dt(p,ke) - ( &
473 e11 * dflux(p,1,momy_vid) &
474 + e22 * dflux(p,2,momy_vid) &
475 + e33 * dflux(p,3,momy_vid) &
476 + dflux(p,4,momy_vid) ) * rgsqrt(p)
477 end do
478 end do
479 !$omp end parallel
480 call prof_rapend('cal_dyn_tend_interior', 3)
481 return
subroutine, public atm_dyn_dgm_nonhydro3d_rhot_hevi_numflux_get_generalvc(del_flux, ddens_, momx_, momy_, momz_, drhot_, dpres, dens_hyd, pres_hyd, therm_hyd, rtot, cvtot, cptot, gsqrt, g13, g23, nx, ny, nz, vmapm, vmapp, lmesh, elem, lmesh2d, elem2d)

References scale_atm_dyn_dgm_nonhydro3d_rhot_hevi_numflux::atm_dyn_dgm_nonhydro3d_rhot_hevi_numflux_get_generalvc().

◆ atm_dyn_dgm_nonhydro3d_rhot_hevi_cal_vi_asis()

subroutine, public scale_atm_dyn_dgm_nonhydro3d_rhot_hevi::atm_dyn_dgm_nonhydro3d_rhot_hevi_cal_vi_asis ( real(rp), dimension(elem%np,lmesh%nea), intent(out) dens_dt,
real(rp), dimension(elem%np,lmesh%nea), intent(out) momx_dt,
real(rp), dimension(elem%np,lmesh%nea), intent(out) momy_dt,
real(rp), dimension(elem%np,lmesh%nea), intent(out) momz_dt,
real(rp), dimension(elem%np,lmesh%nea), intent(out) rhot_dt,
real(rp), dimension(elem%np,lmesh%nea), intent(in) ddens_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) momx_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) momy_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) momz_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) drhot_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) dens_hyd,
real(rp), dimension(elem%np,lmesh%nea), intent(in) pres_hyd,
real(rp), dimension(elem%np,lmesh%nea), intent(in) ddens0_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) momx0_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) momy0_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) momz0_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) drhot0_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) rtot,
real(rp), dimension(elem%np,lmesh%nea), intent(in) cvtot,
real(rp), dimension(elem%np,lmesh%nea), intent(in) cptot,
class(elementoperationbase3d), intent(in) element3d_operation,
class(sparsemat), intent(in) dz,
class(sparsemat), intent(in) lift,
real(rp), intent(in) impl_fac,
real(rp), intent(in) dt,
class(localmesh3d), intent(in) lmesh,
class(elementbase3d), intent(in) elem,
class(localmesh2d), intent(in) lmesh2d,
class(elementbase2d), intent(in) elem2d )

Definition at line 487 of file scale_atm_dyn_dgm_nonhydro3d_rhot_hevi.F90.

495
503 implicit none
504
505 class(LocalMesh3D), intent(in) :: lmesh
506 class(ElementBase3D), intent(in) :: elem
507 class(LocalMesh2D), intent(in) :: lmesh2D
508 class(ElementBase2D), intent(in) :: elem2D
509 real(RP), intent(out) :: DENS_dt(elem%Np,lmesh%NeA)
510 real(RP), intent(out) :: MOMX_dt(elem%Np,lmesh%NeA)
511 real(RP), intent(out) :: MOMY_dt(elem%Np,lmesh%NeA)
512 real(RP), intent(out) :: MOMZ_dt(elem%Np,lmesh%NeA)
513 real(RP), intent(out) :: RHOT_dt(elem%Np,lmesh%NeA)
514 real(RP), intent(in) :: DDENS_(elem%Np,lmesh%NeA)
515 real(RP), intent(in) :: MOMX_(elem%Np,lmesh%NeA)
516 real(RP), intent(in) :: MOMY_(elem%Np,lmesh%NeA)
517 real(RP), intent(in) :: MOMZ_(elem%Np,lmesh%NeA)
518 real(RP), intent(in) :: DRHOT_(elem%Np,lmesh%NeA)
519 real(RP), intent(in) :: DENS_hyd(elem%Np,lmesh%NeA)
520 real(RP), intent(in) :: PRES_hyd(elem%Np,lmesh%NeA)
521 real(RP), intent(in) :: DDENS0_(elem%Np,lmesh%NeA)
522 real(RP), intent(in) :: MOMX0_(elem%Np,lmesh%NeA)
523 real(RP), intent(in) :: MOMY0_(elem%Np,lmesh%NeA)
524 real(RP), intent(in) :: MOMZ0_(elem%Np,lmesh%NeA)
525 real(RP), intent(in) :: DRHOT0_(elem%Np,lmesh%NeA)
526 real(RP), intent(in) :: Rtot(elem%Np,lmesh%NeA)
527 real(RP), intent(in) :: CVtot(elem%Np,lmesh%NeA)
528 real(RP), intent(in) :: CPtot(elem%Np,lmesh%NeA)
529 class(ElementOperationBase3D), intent(in) :: element3D_operation
530 class(SparseMat), intent(in) :: Dz, Lift
531 real(RP), intent(in) :: impl_fac
532 real(RP), intent(in) :: dt
533
534 real(RP) :: PROG_VARS (elem%Np,lmesh%NeZ,PRGVAR_NUM,lmesh%NeX*lmesh%NeY)
535 real(RP) :: PROG_VARS0(elem%Np,lmesh%NeZ,PRGVAR_NUM,lmesh%NeX*lmesh%NeY)
536 real(RP) :: b1D(3,elem%Nnode_v,lmesh%NeZ,elem%Nnode_h1D**2,lmesh%NeX*lmesh%NeY)
537 integer :: ipiv(elem%Nnode_v*3*lmesh%NeZ,elem%Nnode_h1D**2)
538 real(RP) :: b1D_uv(elem%Nnode_v,lmesh%NeZ,2,elem%Nnode_h1D**2,lmesh%NeX*lmesh%NeY)
539 integer :: ipiv_uv(elem%Nnode_v*1*lmesh%NeZ,elem%Nnode_h1D**2)
540 real(RP) :: alph(elem%NfpTot,lmesh%NeZ,lmesh%NeX*lmesh%NeY)
541 real(RP) :: Rtot_z(elem%Np,lmesh%NeZ,lmesh%NeX*lmesh%NeY)
542 real(RP) :: CPtot_ov_CVtot(elem%Np,lmesh%NeZ,lmesh%NeX*lmesh%NeY)
543 real(RP) :: DENS_hyd_z(elem%Np,lmesh%NeZ,lmesh%NeX*lmesh%NeY)
544 real(RP) :: PRES_hyd_z(elem%Np,lmesh%NeZ,lmesh%NeX*lmesh%NeY)
545 real(RP) :: GnnM_z(elem%Np,lmesh%NeZ,lmesh%NeX*lmesh%NeY)
546 real(RP) :: G13_z(elem%Np,lmesh%NeZ,lmesh%NeX*lmesh%NeY)
547 real(RP) :: G23_z(elem%Np,lmesh%NeZ,lmesh%NeX*lmesh%NeY)
548 real(RP) :: GsqrtV_z(elem%Np,lmesh%NeZ,lmesh%NeX*lmesh%NeY)
549 real(RP) :: nz(elem%NfpTot,lmesh%NeZ,lmesh%NeX*lmesh%NeY)
550 integer :: vmapM(elem%NfpTot,lmesh%NeZ)
551 integer :: vmapP(elem%NfpTot,lmesh%NeZ)
552 integer :: ColMask(elem%Nnode_v)
553 integer :: ke_xy, ke_z, ke, ke2D, v
554 integer :: itr_nlin
555 integer :: kl, ku, nz_1D
556 integer :: kl_uv, ku_uv, nz_1D_uv
557 integer :: ij
558 logical :: is_converged
559
560 real(RP), allocatable :: PmatBnd(:,:,:)
561 real(RP), allocatable :: PmatBnd_uv(:,:,:)
562 integer :: info, info_uv
563 !------------------------------------------------------------------------
564
565 call prof_rapstart( 'hevi_cal_vi_prep', 3)
566
567 nz_1d = elem%Nnode_v * 3 * lmesh%NeZ
568 kl = ( elem%Nnode_v + 1 ) * 3 - 1
569 ku = kl
570 nz_1d_uv = elem%Nnode_v * 1 * lmesh%NeZ
571 kl_uv = elem%Nnode_v
572 ku_uv = kl_uv
573 allocate( pmatbnd(2*kl+ku+1,nz_1d,elem%Nnode_h1D**2) )
574 allocate( pmatbnd_uv(2*kl_uv+ku_uv+1,nz_1d_uv,elem%Nnode_h1D**2) )
575
576 call lmesh%GetVmapZ1D( vmapm, vmapp ) ! (out)
577
578 !-
579
580 !$omp parallel private( ke_xy, ke_z, ke, ke2D )
581 !$omp do collapse(2)
582 do ke_xy=1, lmesh%NeX*lmesh%NeY
583 do ke_z=1, lmesh%NeZ
584 ke = ke_xy + (ke_z-1)*lmesh%NeX*lmesh%NeY
585 ke2d = lmesh%EMap3Dto2D(ke)
586
587 prog_vars(:,ke_z,dens_vid,ke_xy) = ddens0_(:,ke)
588 prog_vars(:,ke_z,momx_vid,ke_xy) = momx0_(:,ke)
589 prog_vars(:,ke_z,momy_vid,ke_xy) = momy0_(:,ke)
590 prog_vars(:,ke_z,momz_vid,ke_xy) = momz0_(:,ke)
591 prog_vars(:,ke_z,rhot_vid,ke_xy) = drhot0_(:,ke)
592
593 dens_hyd_z(:,ke_z,ke_xy) = dens_hyd(:,ke)
594 pres_hyd_z(:,ke_z,ke_xy) = pres_hyd(:,ke)
595
596 rtot_z(:,ke_z,ke_xy) = rtot(:,ke)
597 cptot_ov_cvtot(:,ke_z,ke_xy) = cptot(:,ke) / cvtot(:,ke)
598
599 nz(:,ke_z,ke_xy) = lmesh%normal_fn(:,ke,3)
600 g13_z(:,ke_z,ke_xy) = lmesh%GI3(:,ke,1)
601 g23_z(:,ke_z,ke_xy) = lmesh%GI3(:,ke,2)
602 gsqrtv_z(:,ke_z,ke_xy) = lmesh%Gsqrt(:,ke) / lmesh%GsqrtH(elem%IndexH2Dto3D,ke2d)
603
604 gnnm_z(:,ke_z,ke_xy) = ( 1.0_rp / gsqrtv_z(:,ke_z,ke_xy)**2 &
605 + g13_z(:,ke_z,ke_xy)**2 + g23_z(:,ke_z,ke_xy) )
606 end do
607 end do
608 !$omp end do
609 !$omp workshare
610 prog_vars0(:,:,:,:) = prog_vars(:,:,:,:)
611 !$omp end workshare
612 !$omp end parallel
613
614 call prof_rapend( 'hevi_cal_vi_prep', 3)
615
616 !--
617
618 if ( abs(impl_fac) > 0.0_rp ) then
619 call prof_rapstart( 'hevi_cal_vi_itr', 3)
620
621 ! G = (q^n+1 - q^n*) + impl_fac * A(q^n+1) = 0
622 ! dG/dq^n+1 del[q] = - G(q^n*)
623 do itr_nlin = 1, 1
624 call prof_rapstart( 'hevi_cal_vi_ax', 3)
625
626 call vi_eval_ax_uv( &
627 momx_dt(:,:), momy_dt(:,:), alph(:,:,:), & ! (out)
628 prog_vars, prog_vars0, & ! (in)
629 ddens_, momx_, momy_, momz_, drhot_, & ! (in)
630 dens_hyd_z, pres_hyd_z, & ! (in)
631 rtot_z, cptot_ov_cvtot, & ! (in)
632 dz, lift, intrpmat_vpordm1, & ! (in)
633 gnnm_z, g13_z, g23_z, gsqrtv_z, & ! (in)
634 impl_fac, dt, & ! (in)
635 lmesh, elem, nz, vmapm, vmapp, & ! (in)
636 b1d_uv(:,:,:,:,:) ) ! (out)
637
638 call prof_rapend( 'hevi_cal_vi_ax', 3)
639
640 do ke_xy=1, lmesh%NeX * lmesh%NeY
641 call prof_rapstart( 'hevi_cal_vi_matbnd', 3)
642
643 call vi_construct_matbnd_uv( pmatbnd_uv(:,:,:), & ! (out)
644 kl_uv, ku_uv, nz_1d_uv, & ! (in)
645 prog_vars(:,:,:,ke_xy), & ! (in)
646 dens_hyd_z(:,:,ke_xy), pres_hyd_z(:,:,ke_xy), & ! (in)
647 g13_z(:,:,ke_xy), g23_z(:,:,ke_xy), gsqrtv_z(:,:,ke_xy), & ! (in)
648 alph(:,:,ke_xy), & ! (in)
649 rtot_z(:,:,ke_xy), cptot_ov_cvtot(:,:,ke_xy), & ! (in)
650 dz, lift, intrpmat_vpordm1, & ! (in)
651 impl_fac, dt, & ! (in)
652 lmesh, elem, nz(:,:,ke_xy), vmapm, vmapp, ke_xy, 1 ) ! (in)
653
654 call prof_rapend( 'hevi_cal_vi_matbnd', 3)
655
656 call prof_rapstart( 'hevi_cal_vi_lin', 3)
657 !$omp parallel private(ij, v, ke_z, info_uv, ColMask)
658 !$omp do
659 do ij=1, elem%Nnode_h1D**2
660 call linalgebra_solvelineq_bndmat( pmatbnd_uv(:,:,ij), b1d_uv(:,:,:,ij,ke_xy), ipiv_uv(:,ij), nz_1d_uv, kl_uv, ku_uv, 2, vi_use_lapack_flag )
661
662 colmask(:) = elem%Colmask(:,ij)
663 do ke_z=1, lmesh%NeZ
664 prog_vars(colmask(:),ke_z,momx_vid,ke_xy) = prog_vars(colmask(:),ke_z,momx_vid,ke_xy) + b1d_uv(:,ke_z,1,ij,ke_xy)
665 prog_vars(colmask(:),ke_z,momy_vid,ke_xy) = prog_vars(colmask(:),ke_z,momy_vid,ke_xy) + b1d_uv(:,ke_z,2,ij,ke_xy)
666 end do
667 end do ! for ij
668 !$omp end do
669 !$omp end parallel
670 call prof_rapend( 'hevi_cal_vi_lin', 3)
671
672 end do ! for ke_xy
673
674 call prof_rapstart( 'hevi_cal_vi_ax', 3)
675 call vi_eval_ax( &
676 dens_dt(:,:), momz_dt(:,:), rhot_dt(:,:), & ! (out, dummy)
677 alph(:,:,:), & ! (in)
678 prog_vars, prog_vars0, & ! (in)
679 ddens_, momx_, momy_, momz_, drhot_, & ! (in)
680 dens_hyd_z, pres_hyd_z, & ! (in)
681 rtot_z, cptot_ov_cvtot, & ! (in)
682 dz, lift, intrpmat_vpordm1, & ! (in)
683 gnnm_z, g13_z, g23_z, gsqrtv_z, & ! (in)
684 impl_fac, dt, & ! (in)
685 lmesh, elem, nz, vmapm, vmapp, & ! (in)
686 b1d(:,:,:,:,:) ) ! (out)
687 call prof_rapend( 'hevi_cal_vi_ax', 3)
688
689 do ke_xy=1, lmesh%NeX * lmesh%NeY
690 call prof_rapstart( 'hevi_cal_vi_matbnd', 3)
691 call vi_construct_matbnd( pmatbnd(:,:,:), & ! (out)
692 kl, ku, nz_1d, & ! (in)
693 prog_vars(:,:,:,ke_xy), & ! (in)
694 dens_hyd_z(:,:,ke_xy), pres_hyd_z(:,:,ke_xy), & ! (in)
695 g13_z(:,:,ke_xy), g23_z(:,:,ke_xy), gsqrtv_z(:,:,ke_xy), & ! (in)
696 alph(:,:,ke_xy), & ! (in)
697 rtot_z(:,:,ke_xy), cptot_ov_cvtot(:,:,ke_xy), & ! (in)
698 dz, lift, intrpmat_vpordm1, & ! (in)
699 impl_fac, dt, & ! (in)
700 lmesh, elem, nz(:,:,ke_xy), vmapm, vmapp, ke_xy, 1 ) ! (in)
701
702 call prof_rapend( 'hevi_cal_vi_matbnd', 3)
703
704 call prof_rapstart( 'hevi_cal_vi_lin', 3)
705 !$omp parallel private(ij, v, ke_z, info, ColMask)
706 !$omp do
707 do ij=1, elem%Nnode_h1D**2
708 call linalgebra_solvelineq_bndmat( pmatbnd(:,:,ij), b1d(:,:,:,ij,ke_xy), ipiv(:,ij), nz_1d, kl, ku, 1, vi_use_lapack_flag )
709 colmask(:) = elem%Colmask(:,ij)
710 do ke_z=1, lmesh%NeZ
711 prog_vars(colmask(:),ke_z,dens_vid,ke_xy) = prog_vars(colmask(:),ke_z,dens_vid,ke_xy) + b1d(1,:,ke_z,ij,ke_xy)
712 prog_vars(colmask(:),ke_z,momz_vid,ke_xy) = prog_vars(colmask(:),ke_z,momz_vid,ke_xy) + b1d(2,:,ke_z,ij,ke_xy)
713 prog_vars(colmask(:),ke_z,rhot_vid,ke_xy) = prog_vars(colmask(:),ke_z,rhot_vid,ke_xy) + b1d(3,:,ke_z,ij,ke_xy)
714 end do
715 end do ! for ij
716 !$omp end do
717 !$omp end parallel
718 call prof_rapend( 'hevi_cal_vi_lin', 3)
719
720 end do ! for ke_xy
721 end do ! itr nlin
722
723 call prof_rapend( 'hevi_cal_vi_itr', 3)
724 end if
725
726
727 call prof_rapstart( 'hevi_cal_vi_retrun_var', 3)
728 if ( abs(impl_fac) > 0.0_rp) then
729 !$omp parallel do collapse(2) private(ke_xy, ke_z, ke)
730 do ke_xy=1, lmesh%NeX * lmesh%NeY
731 do ke_z=1, lmesh%NeZ
732 ke = ke_xy + (ke_z-1)*lmesh%NeX*lmesh%NeY
733 dens_dt(:,ke) = ( prog_vars(:,ke_z,dens_vid,ke_xy) - ddens_(:,ke) ) / impl_fac
734 momx_dt(:,ke) = ( prog_vars(:,ke_z,momx_vid,ke_xy) - momx_(:,ke) ) / impl_fac
735 momy_dt(:,ke) = ( prog_vars(:,ke_z,momy_vid,ke_xy) - momy_(:,ke) ) / impl_fac
736 momz_dt(:,ke) = ( prog_vars(:,ke_z,momz_vid,ke_xy) - momz_(:,ke) ) / impl_fac
737 rhot_dt(:,ke) = ( prog_vars(:,ke_z,rhot_vid,ke_xy) - drhot_(:,ke) ) / impl_fac
738 end do
739 end do
740 else
741 call vi_eval_ax_uv( &
742 momx_dt(:,:), momy_dt(:,:), & ! (out)
743 alph(:,:,:), & ! (out, dummy)
744 prog_vars, prog_vars0, & ! (in)
745 ddens_, momx_, momy_, momz_, drhot_, & ! (in)
746 dens_hyd_z, pres_hyd_z, & ! (in)
747 rtot_z, cptot_ov_cvtot, & ! (in)
748 dz, lift, intrpmat_vpordm1, & ! (in)
749 gnnm_z, g13_z, g23_z, gsqrtv_z, & ! (in)
750 impl_fac, dt, & ! (in)
751 lmesh, elem, nz, vmapm, vmapp ) ! (in)
752
753 call vi_eval_ax( &
754 dens_dt(:,:), momz_dt(:,:), rhot_dt(:,:), & ! (out)
755 alph(:,:,:), & ! (in, dummy)
756 prog_vars, prog_vars0, & ! (in)
757 ddens_, momx_, momy_, momz_, drhot_, & ! (in)
758 dens_hyd_z, pres_hyd_z, & ! (in)
759 rtot_z, cptot_ov_cvtot, & ! (in)
760 dz, lift, intrpmat_vpordm1, & ! (in)
761 gnnm_z, g13_z, g23_z, gsqrtv_z, & ! (in)
762 impl_fac, dt, & ! (in)
763 lmesh, elem, nz, vmapm, vmapp ) ! (in)
764 end if
765 call prof_rapend( 'hevi_cal_vi_retrun_var', 3)
766
767 return
subroutine, public atm_dyn_dgm_nonhydro3d_rhot_hevi_common_eval_ax(dens_t, momz_t, rhot_t, alph, prog_vars, prog_vars0, ddens00, momx00, momy00, momz00, drhot00, dens_hyd, pres_hyd, rtot, cptot_ov_cvtot, dz, lift, intrpmat_vpordm1, gnnm, g13, g23, gsqrtv, impl_fac, dt, lmesh, elem, nz, vmapm, vmapp, b1d_ij)
subroutine, public atm_dyn_dgm_nonhydro3d_rhot_hevi_common_eval_ax_uv(momx_t, momy_t, alph, prog_vars, prog_vars0, ddens00, momx00, momy00, momz00, drhot00, dens_hyd, pres_hyd, rtot, cptot_ov_cvtot, dz, lift, intrpmat_vpordm1, gnnm, g13, g23, gsqrtv, impl_fac, dt, lmesh, elem, nz, vmapm, vmapp, b1d_ij_uv)
subroutine, public atm_dyn_dgm_nonhydro3d_rhot_hevi_common_construct_matbnd(pmatbnd, kl, ku, nz_1d, prog_vars0, dens_hyd, pres_hyd, g13, g23, gsqrtv, alph, rtot, cptot_ov_cvtot, dz, lift, intrpmat_vpordm1, impl_fac, dt, lmesh, elem, nz, vmapm, vmapp, ke_x, ke_y)
subroutine, public atm_dyn_dgm_nonhydro3d_rhot_hevi_common_construct_matbnd_uv(pmatbnd_uv, kl_uv, ku_uv, nz_1d_uv, prog_vars0, dens_hyd, pres_hyd, g13, g23, gsqrtv, alph, rtot, cptot_ov_cvtot, dz, lift, intrpmat_vpordm1, impl_fac, dt, lmesh, elem, nz, vmapm, vmapp, ke_x, ke_y)
Module common / Linear algebra.
subroutine, public linalgebra_solvelineq_bndmat(a, b, ipiv, n, kl, ku, nrhs, use_lapack)
Calculate the solution of linear equations with a band matrix.

References scale_atm_dyn_dgm_nonhydro3d_rhot_hevi_common::atm_dyn_dgm_nonhydro3d_rhot_hevi_common_construct_matbnd(), scale_atm_dyn_dgm_nonhydro3d_rhot_hevi_common::atm_dyn_dgm_nonhydro3d_rhot_hevi_common_construct_matbnd_uv(), scale_atm_dyn_dgm_nonhydro3d_rhot_hevi_common::atm_dyn_dgm_nonhydro3d_rhot_hevi_common_eval_ax(), scale_atm_dyn_dgm_nonhydro3d_rhot_hevi_common::atm_dyn_dgm_nonhydro3d_rhot_hevi_common_eval_ax_uv(), scale_atm_dyn_dgm_nonhydro3d_common::intrpmat_vpordm1, scale_linalgebra::linalgebra_solvelineq_bndmat(), and scale_atm_dyn_dgm_nonhydro3d_rhot_hevi_common::vi_use_lapack_flag.

◆ atm_dyn_dgm_nonhydro3d_rhot_hevi_cal_vi()

subroutine, public scale_atm_dyn_dgm_nonhydro3d_rhot_hevi::atm_dyn_dgm_nonhydro3d_rhot_hevi_cal_vi ( real(rp), dimension(elem%np,lmesh%nea), intent(out) dens_dt,
real(rp), dimension(elem%np,lmesh%nea), intent(out) momx_dt,
real(rp), dimension(elem%np,lmesh%nea), intent(out) momy_dt,
real(rp), dimension(elem%np,lmesh%nea), intent(out) momz_dt,
real(rp), dimension(elem%np,lmesh%nea), intent(out) rhot_dt,
real(rp), dimension(elem%np,lmesh%nea), intent(in) ddens_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) momx_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) momy_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) momz_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) drhot_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) dens_hyd,
real(rp), dimension(elem%np,lmesh%nea), intent(in) pres_hyd,
real(rp), dimension(elem%np,lmesh%nea), intent(in) ddens0_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) momx0_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) momy0_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) momz0_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) drhot0_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) rtot,
real(rp), dimension(elem%np,lmesh%nea), intent(in) cvtot,
real(rp), dimension(elem%np,lmesh%nea), intent(in) cptot,
class(elementoperationbase3d), intent(in) element3d_operation,
class(sparsemat), intent(in) dz,
class(sparsemat), intent(in) lift,
real(rp), intent(in) impl_fac,
real(rp), intent(in) dt,
class(localmesh3d), intent(in) lmesh,
class(elementbase3d), intent(in) elem,
class(localmesh2d), intent(in) lmesh2d,
class(elementbase2d), intent(in) elem2d )

Definition at line 772 of file scale_atm_dyn_dgm_nonhydro3d_rhot_hevi.F90.

780
788
789 implicit none
790
791 class(LocalMesh3D), intent(in) :: lmesh
792 class(ElementBase3D), intent(in) :: elem
793 class(LocalMesh2D), intent(in) :: lmesh2D
794 class(ElementBase2D), intent(in) :: elem2D
795 real(RP), intent(out) :: DENS_dt(elem%Np,lmesh%NeA)
796 real(RP), intent(out) :: MOMX_dt(elem%Np,lmesh%NeA)
797 real(RP), intent(out) :: MOMY_dt(elem%Np,lmesh%NeA)
798 real(RP), intent(out) :: MOMZ_dt(elem%Np,lmesh%NeA)
799 real(RP), intent(out) :: RHOT_dt(elem%Np,lmesh%NeA)
800 real(RP), intent(in) :: DDENS_(elem%Np,lmesh%NeA)
801 real(RP), intent(in) :: MOMX_(elem%Np,lmesh%NeA)
802 real(RP), intent(in) :: MOMY_(elem%Np,lmesh%NeA)
803 real(RP), intent(in) :: MOMZ_(elem%Np,lmesh%NeA)
804 real(RP), intent(in) :: DRHOT_(elem%Np,lmesh%NeA)
805 real(RP), intent(in) :: DENS_hyd(elem%Np,lmesh%NeA)
806 real(RP), intent(in) :: PRES_hyd(elem%Np,lmesh%NeA)
807 real(RP), intent(in) :: DDENS0_(elem%Np,lmesh%NeA)
808 real(RP), intent(in) :: MOMX0_(elem%Np,lmesh%NeA)
809 real(RP), intent(in) :: MOMY0_(elem%Np,lmesh%NeA)
810 real(RP), intent(in) :: MOMZ0_(elem%Np,lmesh%NeA)
811 real(RP), intent(in) :: DRHOT0_(elem%Np,lmesh%NeA)
812 real(RP), intent(in) :: Rtot(elem%Np,lmesh%NeA)
813 real(RP), intent(in) :: CVtot(elem%Np,lmesh%NeA)
814 real(RP), intent(in) :: CPtot(elem%Np,lmesh%NeA)
815 class(ElementOperationBase3D), intent(in) :: element3D_operation
816 class(SparseMat), intent(in) :: Dz, Lift
817 real(RP), intent(in) :: impl_fac
818 real(RP), intent(in) :: dt
819
820 real(RP) :: PROG_VARS (elem%Np,lmesh%NeX*lmesh%NeY,lmesh%NeZ,PRGVAR_NUM)
821 real(RP) :: PROG_VARS0(elem%Np,lmesh%NeX*lmesh%NeY,lmesh%NeZ,PRGVAR_NUM)
822 real(RP) :: alph(elem%NfpTot,lmesh%Ne)
823 real(RP) :: GsqrtV(elem%Np,lmesh%Ne)
824
825 integer :: vmapM(elem%NfpTot,lmesh%Ne)
826 integer :: vmapP(elem%NfpTot,lmesh%Ne)
827 integer :: ke_xy, ke_z, ke, ke2d
828 integer :: p
829 integer :: itr_nlin
830 integer :: ij, i, j, im, jm
831 logical :: is_converged
832
833 real(RP), allocatable :: b_uv(:,:,:,:)
834 real(RP), allocatable :: b (:,:,:)
835
836 real(RP) :: DENS(elem%Np,lmesh%Ne)
837 real(RP) :: W(elem%Np,lmesh%Ne)
838 real(RP) :: WT(elem%Np,lmesh%Ne)
839 real(RP) :: POT(elem%Np,lmesh%Ne)
840 real(RP) :: DPDRHOT(elem%Np,lmesh%Ne)
841 !------------------------------------------------------------------------
842
843 call prof_rapstart( 'hevi_cal_vi_prep', 3)
844
845 call lmesh%GetVmapZ3D( vmapm, vmapp ) ! (out)
847 elem%Nnode_h1D ) ! (in)
848
849 !-
850 !$omp parallel private( ke_xy, ke_z, ke, ke2D, p )
851 !$omp do collapse(2)
852 do ke_z =1, lmesh%NeZ
853 do ke_xy=1, lmesh%NeX * lmesh%NeY
854 ke = ke_xy + (ke_z-1)*lmesh%NeX*lmesh%NeY
855 ke2d = lmesh%EMap3Dto2D(ke)
856
857 prog_vars(:,ke_xy,ke_z,dens_vid) = ddens0_(:,ke)
858 prog_vars(:,ke_xy,ke_z,momx_vid) = momx0_(:,ke)
859 prog_vars(:,ke_xy,ke_z,momy_vid) = momy0_(:,ke)
860 prog_vars(:,ke_xy,ke_z,momz_vid) = momz0_(:,ke)
861 prog_vars(:,ke_xy,ke_z,rhot_vid) = drhot0_(:,ke)
862
863 do p=1, elem%Np
864 gsqrtv(p,ke) = lmesh%Gsqrt(p,ke) / lmesh%GsqrtH(elem%IndexH2Dto3D(p),ke2d)
865 end do
866 end do
867 end do
868 !$omp end do
869 !$omp workshare
870 prog_vars0(:,:,:,:) = prog_vars(:,:,:,:)
871 !$omp end workshare
872 !$omp end parallel
873 call prof_rapend( 'hevi_cal_vi_prep', 3)
874
875 !--
876
877 if ( abs(impl_fac) > 0.0_rp ) then
878 call prof_rapstart( 'hevi_cal_vi_itr', 3)
879
880 allocate( b_uv(im*elem%Nnode_v,2,jm,lmesh%Ne) )
881 allocate( b(im*3*elem%Nnode_v,jm,lmesh%Ne) )
882
883 ! G = (q^n+1 - q^n*) + impl_fac * A(q^n+1) = 0
884 ! dG/dq^n+1 del[q] = - G(q^n*)
885 do itr_nlin = 1, 1
886 call prof_rapstart( 'hevi_cal_vi_ax_uv', 3)
887 call vi_eval_ax_uv( &
888 momx_dt(:,:), momy_dt(:,:), alph(:,:), & ! (out)
889 prog_vars, prog_vars0, & ! (in)
890 ddens_, momx_, momy_, momz_, drhot_, & ! (in)
891 dens_hyd, pres_hyd, & ! (in)
892 rtot, cptot, cvtot, gsqrtv, & ! (in)
893 impl_fac, dt, & ! (in)
894 lmesh, elem, vmapm, vmapp, & ! (in)
895 element3d_operation, & ! (in)
896 im, jm, & ! (in)
897 b_uv ) ! (out)
898 call prof_rapend( 'hevi_cal_vi_ax_uv', 3)
899
900 call vi_solve_uv( prog_vars, & ! (inout)
901 b_uv, alph, & ! (in)
902 gsqrtv, & ! (in)
903 impl_fac, dt, & ! (in)
904 vmapm, vmapp, lmesh, elem, im, jm ) ! (in)
905
906 !----
907
908 call prof_rapstart( 'hevi_cal_vi_ax', 3)
909 call vi_eval_ax( &
910 dens_dt(:,:), momz_dt(:,:), rhot_dt(:,:), & ! (out)
911 alph(:,:), & ! (in)
912 prog_vars, prog_vars0, ddens_, momx_, momy_, momz_, drhot_, & ! (in)
913 dens_hyd, pres_hyd, rtot, cptot, cvtot, gsqrtv, & ! (in)
914 impl_fac, dt, lmesh, elem, vmapm, vmapp, & ! (in)
915 element3d_operation, & ! (in)
916 im, jm, & ! (in)
917 b, dens, w, wt, pot, dpdrhot ) ! (out)
918 call prof_rapend( 'hevi_cal_vi_ax', 3)
919
920 call vi_solve( prog_vars, & ! (inout)
921 b, dens, w, wt, pot, dpdrhot, alph, & ! (in)
922 gsqrtv, intrpmat_vpordm1, & ! (in)
923 impl_fac, dt, & ! (in)
924 vmapm, vmapp, lmesh, elem, im, jm ) ! (in)
925 end do ! itr nlin
926
927 call prof_rapend( 'hevi_cal_vi_itr', 3)
928 end if
929
930 call prof_rapstart( 'hevi_cal_vi_retrun_var', 3)
931 if ( abs(impl_fac) > 0.0_rp) then
932 !$omp parallel do collapse(2) private(ke_xy, ke_z, ke)
933 do ke_z=1, lmesh%NeZ
934 do ke_xy=1, lmesh%NeX * lmesh%NeY
935 ke = ke_xy + (ke_z-1)*lmesh%NeX*lmesh%NeY
936 dens_dt(:,ke) = ( prog_vars(:,ke_xy,ke_z,dens_vid) - ddens_(:,ke) ) / impl_fac
937 momx_dt(:,ke) = ( prog_vars(:,ke_xy,ke_z,momx_vid) - momx_(:,ke) ) / impl_fac
938 momy_dt(:,ke) = ( prog_vars(:,ke_xy,ke_z,momy_vid) - momy_(:,ke) ) / impl_fac
939 momz_dt(:,ke) = ( prog_vars(:,ke_xy,ke_z,momz_vid) - momz_(:,ke) ) / impl_fac
940 rhot_dt(:,ke) = ( prog_vars(:,ke_xy,ke_z,rhot_vid) - drhot_(:,ke) ) / impl_fac
941 end do
942 end do
943 else
944 call vi_eval_ax_uv( &
945 momx_dt(:,:), momy_dt(:,:), alph(:,:), & ! (out)
946 prog_vars, prog_vars0, & ! (in)
947 ddens_, momx_, momy_, momz_, drhot_, & ! (in)
948 dens_hyd, pres_hyd, & ! (in)
949 rtot, cptot, cvtot, gsqrtv, & ! (in)
950 impl_fac, dt, & ! (in)
951 lmesh, elem, vmapm, vmapp, & ! (in)
952 element3d_operation, im, jm ) ! (in)
953
954 call vi_eval_ax( &
955 dens_dt(:,:), momz_dt(:,:), rhot_dt(:,:), & ! (out, dummy)
956 alph(:,:), & ! (in)
957 prog_vars, prog_vars0, ddens_, momx_, momy_, momz_, drhot_, & ! (in)
958 dens_hyd, pres_hyd, rtot, cptot, cvtot, gsqrtv, & ! (in)
959 impl_fac, dt, lmesh, elem, vmapm, vmapp, & ! (in)
960 element3d_operation, im, jm ) ! (in)
961 end if
962 call prof_rapend( 'hevi_cal_vi_retrun_var', 3)
963
964 return
module FElib / Fluid dyn solver / Atmosphere / HEVI / Common
subroutine, public atm_dyn_dgm_hevi_common_linalgebra_get_param(im, jm, nnode_h1d)
module FElib / Fluid dyn solver / Atmosphere / Nonhydrostatic model / HEVI / Common
subroutine, public atm_dyn_dgm_nonhydro3d_rhot_hevi_common_solve(prog_vars, b, dens, w, wt, pot, dpdrhot, alph, gsqrtv, intrpmat_vpordm1, impl_fac, dt, vmapm, vmapp, lmesh, elem, im, jm)
subroutine, public atm_dyn_dgm_nonhydro3d_rhot_hevi_common_eval_ax(dens_t, momz_t, rhot_t, alph, prog_vars, prog_vars0, ddens00, momx00, momy00, momz00, drhot00, dens_hyd, pres_hyd, rtot, cptot, cvtot, gsqrtv, impl_fac, dt, lmesh, elem, vmapm, vmapp, element3d_operation, im, jm, b, dens, w, wt, pot, dpdrhot)
subroutine, public atm_dyn_dgm_nonhydro3d_rhot_hevi_common_solve_uv(prog_vars, b_uv, alph, gsqrtv, impl_fac, dt, vmapm, vmapp, lmesh, elem, im, jm)
subroutine, public atm_dyn_dgm_nonhydro3d_rhot_hevi_common_eval_ax_uv(momx_t, momy_t, alph, prog_vars, prog_vars0, ddens00, momx00, momy00, momz00, drhot00, dens_hyd, pres_hyd, rtot, cptot, cvtot, gsqrtv, impl_fac, dt, lmesh, elem, vmapm, vmapp, element3d_operation, im, jm, b1d_ij_uv)

References scale_atm_dyn_dgm_hevi_common_linalgebra::atm_dyn_dgm_hevi_common_linalgebra_get_param(), scale_atm_dyn_dgm_nonhydro3d_rhot_hevi_common_2::atm_dyn_dgm_nonhydro3d_rhot_hevi_common_eval_ax(), scale_atm_dyn_dgm_nonhydro3d_rhot_hevi_common_2::atm_dyn_dgm_nonhydro3d_rhot_hevi_common_eval_ax_uv(), scale_atm_dyn_dgm_nonhydro3d_rhot_hevi_common_2::atm_dyn_dgm_nonhydro3d_rhot_hevi_common_solve(), scale_atm_dyn_dgm_nonhydro3d_rhot_hevi_common_2::atm_dyn_dgm_nonhydro3d_rhot_hevi_common_solve_uv(), and scale_atm_dyn_dgm_nonhydro3d_common::intrpmat_vpordm1.