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

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

Functions/Subroutines

subroutine, public atm_dyn_dgm_globalnonhydro3d_rhot_hevi_init (mesh)
subroutine, public atm_dyn_dgm_globalnonhydro3d_rhot_hevi_final ()
subroutine, public atm_dyn_dgm_globalnonhydro3d_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_globalnonhydro3d_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_globalnonhydro3d_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_globalnonhydro3d_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 / Global nonhydrostatic model / HEVI

Description
HEVI DGM scheme for Global 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, Xuanzhengbo Ren, and Team SCALE

Function/Subroutine Documentation

◆ atm_dyn_dgm_globalnonhydro3d_rhot_hevi_init()

subroutine, public scale_atm_dyn_dgm_globalnonhydro3d_rhot_hevi::atm_dyn_dgm_globalnonhydro3d_rhot_hevi_init ( class(meshbase3d), intent(in) mesh)

Definition at line 77 of file scale_atm_dyn_dgm_globalnonhydro3d_rhot_hevi.F90.

80 implicit none
81 class(MeshBase3D), intent(in) :: mesh
82 !--------------------------------------------
83
84 call atm_dyn_dgm_nonhydro3d_common_init( mesh )
86
87 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_globalnonhydro3d_rhot_hevi_final()

subroutine, public scale_atm_dyn_dgm_globalnonhydro3d_rhot_hevi::atm_dyn_dgm_globalnonhydro3d_rhot_hevi_final

◆ atm_dyn_dgm_globalnonhydro3d_rhot_hevi_cal_tend_asis()

subroutine, public scale_atm_dyn_dgm_globalnonhydro3d_rhot_hevi::atm_dyn_dgm_globalnonhydro3d_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 106 of file scale_atm_dyn_dgm_globalnonhydro3d_rhot_hevi.F90.

113
116 use scale_const, only: &
117 ohm => const_ohm
118 implicit none
119
120 class(LocalMesh3D), intent(in) :: lmesh
121 class(ElementBase3D), intent(in) :: elem
122 class(LocalMesh2D), intent(in) :: lmesh2D
123 class(ElementBase2D), intent(in) :: elem2D
124 class(ElementOperationBase3D), intent(in) :: element3D_operation
125 type(SparseMat), intent(in) :: Dx, Dy, Dz, Sx, Sy, Sz, Lift
126 real(RP), intent(out) :: DENS_dt(elem%Np,lmesh%NeA)
127 real(RP), intent(out) :: MOMX_dt(elem%Np,lmesh%NeA)
128 real(RP), intent(out) :: MOMY_dt(elem%Np,lmesh%NeA)
129 real(RP), intent(out) :: MOMZ_dt(elem%Np,lmesh%NeA)
130 real(RP), intent(out) :: RHOT_dt(elem%Np,lmesh%NeA)
131 real(RP), intent(in) :: DDENS_(elem%Np,lmesh%NeA)
132 real(RP), intent(in) :: MOMX_(elem%Np,lmesh%NeA)
133 real(RP), intent(in) :: MOMY_(elem%Np,lmesh%NeA)
134 real(RP), intent(in) :: MOMZ_(elem%Np,lmesh%NeA)
135 real(RP), intent(in) :: DRHOT_(elem%Np,lmesh%NeA)
136 real(RP), intent(in) :: DPRES_(elem%Np,lmesh%NeA)
137 real(RP), intent(in) :: DENS_hyd(elem%Np,lmesh%NeA)
138 real(RP), intent(in) :: PRES_hyd(elem%Np,lmesh%NeA)
139 real(RP), intent(in) :: PRES_hyd_ref(elem%Np,lmesh%NeA)
140 real(RP), intent(in) :: THERM_hyd(elem%Np,lmesh%NeA)
141 real(RP), intent(in) :: CORIOLIS(elem2D%Np,lmesh2D%NeA)
142 real(RP), intent(in) :: Rtot (elem%Np,lmesh%NeA)
143 real(RP), intent(in) :: CVtot(elem%Np,lmesh%NeA)
144 real(RP), intent(in) :: CPtot(elem%Np,lmesh%NeA)
145 real(RP), intent(in) :: DPhydDx(elem%Np,lmesh%NeA)
146 real(RP), intent(in) :: DPhydDy(elem%Np,lmesh%NeA)
147
148 real(RP) :: Fx(elem%Np), Fy(elem%Np), Fz(elem%Np), LiftDelFlx(elem%Np)
149 real(RP) :: DPRES_hyd(elem%Np), GradPhyd_x(elem%Np), GradPhyd_y(elem%Np)
150 real(RP) :: del_flux(elem%NfpTot,lmesh%Ne,PRGVAR_NUM)
151 real(RP) :: del_flux_hyd(elem%NfpTot,lmesh%Ne,2)
152 real(RP) :: RHOT_(elem%Np)
153 real(RP) :: rdens_(elem%Np), u_(elem%Np), v_(elem%Np), w_(elem%Np), wt_(elem%Np)
154
155 real(RP) :: G11(elem%Np), G12(elem%Np), G22(elem%Np), Rgam2(elem%Np)
156 real(RP) :: GsqrtV(elem%Np), RGsqrtV(elem%Np)
157 real(RP) :: X2D(elem2D%Np,lmesh2D%Ne), Y2D(elem2D%Np,lmesh2D%Ne)
158 real(RP) :: X(elem%Np), Y(elem%Np), twoOVdel2(elem%Np)
159 real(RP) :: CORI(elem%Np,2)
160 logical :: is_panel1to4
161 real(RP) :: s
162
163 integer :: ke, ke2d
164
165 real(RP) :: gamm, rgamm
166 real(RP) :: rP0
167 real(RP) :: RovP0, P0ovR
168 !------------------------------------------------------------------------
169
170 call prof_rapstart('cal_dyn_tend_bndflux', 3)
171 call get_ebnd_flux( &
172 del_flux, del_flux_hyd, & ! (out)
173 ddens_, momx_, momy_, momz_, drhot_, dpres_, dens_hyd, pres_hyd, & ! (in)
174 rtot, cvtot, cptot, & ! (in)
175 lmesh%Gsqrt, lmesh%GIJ(:,:,1,1), lmesh%GIJ(:,:,1,2), lmesh%GIJ(:,:,2,2), & ! (in)
176 lmesh%GsqrtH, lmesh%GI3(:,:,1), lmesh%GI3(:,:,2), & ! (in)
177 lmesh%normal_fn(:,:,1), lmesh%normal_fn(:,:,2), lmesh%normal_fn(:,:,3), & ! (in)
178 lmesh%vmapM, lmesh%vmapP, elem%IndexH2Dto3D_bnd, & ! (in)
179 lmesh, elem, lmesh2d, elem2d ) ! (in)
180 call prof_rapend('cal_dyn_tend_bndflux', 3)
181
182 !-----
183 call prof_rapstart('cal_dyn_tend_interior', 3)
184 gamm = cpdry / cvdry
185 rgamm = cvdry / cpdry
186 rp0 = 1.0_rp / pres00
187 rovp0 = rdry * rp0
188 p0ovr = pres00 / rdry
189
190 s = 1.0_rp
191 is_panel1to4 = .true.
192 if ( lmesh%panelID == 5 ) then
193 is_panel1to4 = .false.
194 else if ( lmesh%panelID == 6 ) then
195 is_panel1to4 = .false.
196 s = - 1.0_rp
197 end if
198
199 !$omp parallel private( &
200 !$omp RHOT_, rdens_, u_, v_, w_, wt_, &
201 !$omp Fx, Fy, Fz, LiftDelFlx, &
202 !$omp DPRES_hyd, GradPhyd_x, GradPhyd_y, &
203 !$omp G11, G12, G22, Rgam2, GsqrtV, RGsqrtV, &
204 !$omp X, Y, twoOVdel2, &
205 !$omp CORI, ke, ke2D )
206
207 !$omp do
208 do ke2d = lmesh2d%NeS, lmesh2d%NeE
209 x2d(:,ke2d) = tan(lmesh2d%pos_en(:,ke2d,1))
210 y2d(:,ke2d) = tan(lmesh2d%pos_en(:,ke2d,2))
211 end do
212
213 !$omp do
214 do ke = lmesh%NeS, lmesh%NeE
215 !--
216 ke2d = lmesh%EMap3Dto2D(ke)
217 rgam2(:) = 1.0_rp / lmesh%gam(:,ke)**2
218 g11(:) = lmesh%GIJ(elem%IndexH2Dto3D,ke2d,1,1) * rgam2(:)
219 g12(:) = lmesh%GIJ(elem%IndexH2Dto3D,ke2d,1,2) * rgam2(:)
220 g22(:) = lmesh%GIJ(elem%IndexH2Dto3D,ke2d,2,2) * rgam2(:)
221 gsqrtv(:) = lmesh%Gsqrt(:,ke) * rgam2(:) / lmesh%GsqrtH(elem%IndexH2Dto3D,ke2d)
222 rgsqrtv(:) = 1.0_rp / gsqrtv(:)
223
224 !--
225 rhot_(:) = p0ovr * ( pres_hyd(:,ke) * rp0 )**rgamm + drhot_(:,ke)
226 ! DPRES_(:) = PRES00 * ( Rtot(:,ke) * rP0 * RHOT_(:) )**( CPtot(:,ke) / CVtot(:,ke) ) &
227 ! - PRES_hyd(:,ke)
228
229 rdens_(:) = 1.0_rp / ( ddens_(:,ke) + dens_hyd(:,ke) )
230 u_(:) = momx_(:,ke) * rdens_(:)
231 v_(:) = momy_(:,ke) * rdens_(:)
232 w_(:) = momz_(:,ke) * rdens_(:)
233 wt_(:) = w_(:) * rgsqrtv(:) + lmesh%GI3(:,ke,1) * u_(:) + lmesh%GI3(:,ke,2) * v_(:)
234
235 x(:) = x2d(elem%IndexH2Dto3D,ke2d)
236 y(:) = y2d(elem%IndexH2Dto3D,ke2d)
237 twoovdel2(:) = 2.0_rp / ( 1.0_rp + x(:)**2 + y(:)**2 )
238
239 cori(:,1) = s * ohm * twoovdel2(:) * ( - x(:) * y(:) * momx_(:,ke) + (1.0_rp + y(:)**2) * momy_(:,ke) )
240 cori(:,2) = s * ohm * twoovdel2(:) * ( - (1.0_rp + x(:)**2) * momx_(:,ke) + x(:) * y(:) * momy_(:,ke) )
241 if ( is_panel1to4 ) then
242 cori(:,1) = s * y(:) * cori(:,1)
243 cori(:,2) = s * y(:) * cori(:,2)
244 end if
245
246 !-- Gradient hydrostatic pressure
247
248 dpres_hyd(:) = pres_hyd(:,ke) - pres_hyd_ref(:,ke)
249
250 call sparsemat_matmul(dx, gsqrtv(:) * dpres_hyd(:), fx)
251 call sparsemat_matmul(dz, gsqrtv(:) * lmesh%GI3(:,ke,1) * dpres_hyd(:), fz)
252 call sparsemat_matmul(lift, lmesh%Fscale(:,ke) * del_flux_hyd(:,ke,1), liftdelflx)
253 gradphyd_x(:) = lmesh%Escale(:,ke,1,1) * fx(:) &
254 + lmesh%Escale(:,ke,3,3) * fz(:) &
255 + liftdelflx(:)
256
257 call sparsemat_matmul(dy, gsqrtv(:) * dpres_hyd(:), fy)
258 call sparsemat_matmul(dz, gsqrtv(:) * lmesh%GI3(:,ke,2) * dpres_hyd(:), fz)
259 call sparsemat_matmul(lift, lmesh%Fscale(:,ke) * del_flux_hyd(:,ke,2), liftdelflx)
260 gradphyd_y(:) = lmesh%Escale(:,ke,2,2) * fy(:) &
261 + lmesh%Escale(:,ke,3,3) * fz(:) &
262 + liftdelflx(:)
263
264 !-- DENS
265 call sparsemat_matmul(dx, lmesh%Gsqrt(:,ke) * momx_(:,ke), fx)
266 call sparsemat_matmul(dy, lmesh%Gsqrt(:,ke) * momy_(:,ke), fy)
267 call sparsemat_matmul(lift, lmesh%Fscale(:,ke) * del_flux(:,ke,dens_vid), liftdelflx)
268
269 dens_dt(:,ke) = - ( &
270 lmesh%Escale(:,ke,1,1) * fx(:) &
271 + lmesh%Escale(:,ke,2,2) * fy(:) &
272 + liftdelflx(:) ) / lmesh%Gsqrt(:,ke)
273
274 !-- MOMX
275 call sparsemat_matmul(dx, lmesh%Gsqrt(:,ke) * ( u_(:) * momx_(:,ke) + g11(:) * dpres_(:,ke) ), fx)
276 call sparsemat_matmul(dy, lmesh%Gsqrt(:,ke) * ( v_(:) * momx_(:,ke) + g12(:) * dpres_(:,ke) ), fy)
277 call sparsemat_matmul(dz, lmesh%Gsqrt(:,ke) * ( wt_(:) * momx_(:,ke) &
278 + ( lmesh%GI3(:,ke,1) * g11(:) + lmesh%GI3(:,ke,2) * g12(:) ) * dpres_(:,ke) ), fz)
279 call sparsemat_matmul(lift, lmesh%Fscale(:,ke) * del_flux(:,ke,momx_vid), liftdelflx)
280
281 momx_dt(:,ke) = &
282 - ( lmesh%Escale(:,ke,1,1) * fx(:) &
283 + lmesh%Escale(:,ke,2,2) * fy(:) &
284 + lmesh%Escale(:,ke,3,3) * fz(:) &
285 + liftdelflx(:) ) / lmesh%Gsqrt(:,ke) &
286 + twoovdel2(:) * y(:) * &
287 ( - x(:) * y(:) * u_(:) + (1.0_rp + y(:)**2) * v_(:) ) * momx_(:,ke) &
288 - ( g11(:) * gradphyd_x(:) + g12(:) * gradphyd_y(:) ) * rgsqrtv(:) &
289 + cori(:,1)
290
291 !-- MOMY
292 call sparsemat_matmul(dx, lmesh%Gsqrt(:,ke) * ( u_(:) * momy_(:,ke) + g12(:) * dpres_(:,ke) ), fx)
293 call sparsemat_matmul(dy, lmesh%Gsqrt(:,ke) * ( v_(:) * momy_(:,ke) + g22(:) * dpres_(:,ke) ), fy)
294 call sparsemat_matmul(dz, lmesh%Gsqrt(:,ke) * ( wt_(:) * momy_(:,ke) &
295 + ( lmesh%GI3(:,ke,1) * g12(:) + lmesh%GI3(:,ke,2) * g22(:) ) * dpres_(:,ke) ), fz)
296 call sparsemat_matmul(lift, lmesh%Fscale(:,ke) * del_flux(:,ke,momy_vid), liftdelflx)
297
298 momy_dt(:,ke) = &
299 - ( lmesh%Escale(:,ke,1,1) * fx(:) &
300 + lmesh%Escale(:,ke,2,2) * fy(:) &
301 + lmesh%Escale(:,ke,3,3) * fz(:) &
302 + liftdelflx(:) ) / lmesh%Gsqrt(:,ke) &
303 + twoovdel2(:) * x(:) * &
304 ( (1.0_rp + x(:)**2) * u_(:) - x(:) * y(:) * v_(:) ) * momy_(:,ke) &
305 - ( g12(:) * gradphyd_x(:) + g22(:) * gradphyd_y(:) ) * rgsqrtv(:) &
306 + cori(:,2)
307
308 !-- MOMZ
309 call sparsemat_matmul(dx, lmesh%Gsqrt(:,ke) * u_(:) * momz_(:,ke), fx)
310 call sparsemat_matmul(dy, lmesh%Gsqrt(:,ke) * v_(:) * momz_(:,ke), fy)
311 call sparsemat_matmul(dz, lmesh%Gsqrt(:,ke) * wt_(:) * momz_(:,ke), fz)
312 call sparsemat_matmul(lift, lmesh%Fscale(:,ke) * del_flux(:,ke,momz_vid), liftdelflx)
313
314 momz_dt(:,ke) = &
315 - ( lmesh%Escale(:,ke,1,1) * fx(:) &
316 + lmesh%Escale(:,ke,2,2) * fy(:) &
317 + lmesh%Escale(:,ke,3,3) * fz(:) &
318 + liftdelflx(:) ) / lmesh%Gsqrt(:,ke)
319
320 !-- RHOT
321 call sparsemat_matmul(dx, lmesh%Gsqrt(:,ke) * u_(:) * rhot_(:), fx)
322 call sparsemat_matmul(dy, lmesh%Gsqrt(:,ke) * v_(:) * rhot_(:), fy)
323 call sparsemat_matmul(lift, lmesh%Fscale(:,ke) * del_flux(:,ke,rhot_vid), liftdelflx)
324
325 rhot_dt(:,ke) = &
326 - ( lmesh%Escale(:,ke,1,1) * fx(:) &
327 + lmesh%Escale(:,ke,2,2) * fy(:) &
328 + liftdelflx(:) ) / lmesh%Gsqrt(:,ke)
329 end do
330 !$omp end parallel
331 call prof_rapend('cal_dyn_tend_interior', 3)
332
333 return
module FElib / Fluid dyn solver / Atmosphere / Nonhydrostatic model / HEVI / Numflux
subroutine, public atm_dyn_dgm_nonhydro3d_rhot_hevi_numflux_get_generalhvc_asis(del_flux, del_flux_hyd, ddens_, momx_, momy_, momz_, drhot_, dpres_, dens_hyd, pres_hyd, rtot, cvtot, cptot, gsqrt, g11, g12, g22, gsqrth, g13, g23, nx, ny, nz, vmapm, vmapp, im2dto3d, lmesh, elem, lmesh2d, elem2d)

References scale_atm_dyn_dgm_nonhydro3d_rhot_hevi_numflux::atm_dyn_dgm_nonhydro3d_rhot_hevi_numflux_get_generalhvc_asis().

◆ atm_dyn_dgm_globalnonhydro3d_rhot_hevi_cal_tend()

subroutine, public scale_atm_dyn_dgm_globalnonhydro3d_rhot_hevi::atm_dyn_dgm_globalnonhydro3d_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 337 of file scale_atm_dyn_dgm_globalnonhydro3d_rhot_hevi.F90.

344
347 use scale_const, only: &
348 ohm => const_ohm
349 implicit none
350
351 class(LocalMesh3D), intent(in) :: lmesh
352 class(ElementBase3D), intent(in) :: elem
353 class(LocalMesh2D), intent(in) :: lmesh2D
354 class(ElementBase2D), intent(in) :: elem2D
355 class(ElementOperationBase3D), intent(in) :: element3D_operation
356 type(SparseMat), intent(in) :: Dx, Dy, Dz, Sx, Sy, Sz, Lift
357 real(RP), intent(out) :: DENS_dt(elem%Np,lmesh%NeA)
358 real(RP), intent(out) :: MOMX_dt(elem%Np,lmesh%NeA)
359 real(RP), intent(out) :: MOMY_dt(elem%Np,lmesh%NeA)
360 real(RP), intent(out) :: MOMZ_dt(elem%Np,lmesh%NeA)
361 real(RP), intent(out) :: RHOT_dt(elem%Np,lmesh%NeA)
362 real(RP), intent(in) :: DDENS_(elem%Np,lmesh%NeA)
363 real(RP), intent(in) :: MOMX_(elem%Np,lmesh%NeA)
364 real(RP), intent(in) :: MOMY_(elem%Np,lmesh%NeA)
365 real(RP), intent(in) :: MOMZ_(elem%Np,lmesh%NeA)
366 real(RP), intent(in) :: DRHOT_(elem%Np,lmesh%NeA)
367 real(RP), intent(in) :: DPRES_(elem%Np,lmesh%NeA)
368 real(RP), intent(in) :: DENS_hyd(elem%Np,lmesh%NeA)
369 real(RP), intent(in) :: PRES_hyd(elem%Np,lmesh%NeA)
370 real(RP), intent(in) :: PRES_hyd_ref(elem%Np,lmesh%NeA)
371 real(RP), intent(in) :: THERM_hyd(elem%Np,lmesh%NeA)
372 real(RP), intent(in) :: CORIOLIS(elem2D%Np,lmesh2D%NeA)
373 real(RP), intent(in) :: Rtot (elem%Np,lmesh%NeA)
374 real(RP), intent(in) :: CVtot(elem%Np,lmesh%NeA)
375 real(RP), intent(in) :: CPtot(elem%Np,lmesh%NeA)
376 real(RP), intent(in) :: DPhydDx(elem%Np,lmesh%NeA)
377 real(RP), intent(in) :: DPhydDy(elem%Np,lmesh%NeA)
378
379 real(RP) :: Flux(elem%Np,3,5), DFlux(elem%Np,4,5)
380 real(RP) :: del_flux(elem%NfpTot,PRGVAR_NUM,lmesh%Ne)
381 real(RP) :: u_, v_, w_, pt_
382 real(RP) :: RHOT_(elem%Np)
383 real(RP) :: RDENS_(elem%Np), GsqrtV(elem%Np), RGsqrtV(elem%Np), RGsqrt(elem%Np)
384 real(RP) :: Gsqrt_, GsqrtDPRES_, E11, E22, E33
385
386 real(RP) :: G11(elem%Np), G12(elem%Np), G22(elem%Np), Rgam2(elem%Np)
387 real(RP) :: X2D(elem2D%Np,lmesh2D%Ne), Y2D(elem2D%Np,lmesh2D%Ne)
388 real(RP) :: X(elem%Np), Y(elem%Np), twoOVdel2(elem%Np)
389 real(RP) :: CORI(elem%Np,2)
390 logical :: is_panel1to4
391 real(RP) :: s
392
393 integer :: ke, ke2d
394 integer :: p
395
396 real(RP) :: gamm, rgamm
397 real(RP) :: rP0
398 real(RP) :: RovP0, P0ovR
399 !------------------------------------------------------------------------
400
401 call prof_rapstart('cal_dyn_tend_bndflux', 3)
402 call get_ebnd_flux( &
403 del_flux, & ! (out)
404 ddens_, momx_, momy_, momz_, drhot_, & ! (in)
405 dpres_, dens_hyd, pres_hyd, therm_hyd, rtot, cvtot, cptot, & ! (in)
406 lmesh%Gsqrt, lmesh%GIJ(:,:,1,1), lmesh%GIJ(:,:,1,2), lmesh%GIJ(:,:,2,2), & ! (in)
407 lmesh%GsqrtH, lmesh%gam, lmesh%GI3(:,:,1), lmesh%GI3(:,:,2), & ! (in)
408 lmesh%normal_fn(:,:,1), lmesh%normal_fn(:,:,2), lmesh%normal_fn(:,:,3), & ! (in)
409 lmesh%vmapM, lmesh%vmapP, elem%IndexH2Dto3D_bnd, & ! (in)
410 lmesh, elem, lmesh2d, elem2d ) ! (in)
411 call prof_rapend('cal_dyn_tend_bndflux', 3)
412
413 !-----
414 call prof_rapstart('cal_dyn_tend_interior', 3)
415 gamm = cpdry / cvdry
416 rgamm = cvdry / cpdry
417 rp0 = 1.0_rp / pres00
418 rovp0 = rdry * rp0
419 p0ovr = pres00 / rdry
420
421 s = 1.0_rp
422 is_panel1to4 = .true.
423 if ( lmesh%panelID == 5 ) then
424 is_panel1to4 = .false.
425 else if ( lmesh%panelID == 6 ) then
426 is_panel1to4 = .false.
427 s = - 1.0_rp
428 end if
429
430 !$omp parallel private( ke, ke2d, p, &
431 !$omp u_, v_, w_, pt_, &
432 !$omp CORI, GsqrtDPRES_, &
433 !$omp RDENS_, RGsqrt, GsqrtV, RGsqrtV, Gsqrt_, Rgam2, &
434 !$omp G11, G12, G22, X, Y, twoOVdel2, &
435 !$omp E11, E22, E33, DFlux, Flux )
436
437 !$omp do
438 do ke2d = lmesh2d%NeS, lmesh2d%NeE
439 x2d(:,ke2d) = tan(lmesh2d%pos_en(:,ke2d,1))
440 y2d(:,ke2d) = tan(lmesh2d%pos_en(:,ke2d,2))
441 end do
442
443 !$omp do
444 do ke = lmesh%NeS, lmesh%NeE
445 !--
446 ke2d = lmesh%EMap3Dto2D(ke)
447 do p=1, elem%Np
448 rgam2(p) = 1.0_rp / lmesh%gam(p,ke)**2
449 g11(p) = lmesh%GIJ(elem%IndexH2Dto3D(p),ke2d,1,1) * rgam2(p)
450 g12(p) = lmesh%GIJ(elem%IndexH2Dto3D(p),ke2d,1,2) * rgam2(p)
451 g22(p) = lmesh%GIJ(elem%IndexH2Dto3D(p),ke2d,2,2) * rgam2(p)
452 end do
453 do p=1, elem%Np
454 gsqrtv(p) = lmesh%Gsqrt(p,ke) * rgam2(p) / lmesh%GsqrtH(elem%IndexH2Dto3D(p),ke2d)
455 rgsqrtv(p) = 1.0_rp / gsqrtv(p)
456 rgsqrt(p) = 1.0_rp / lmesh%Gsqrt(p,ke)
457 rdens_(p) = 1.0_rp / ( ddens_(p,ke) + dens_hyd(p,ke) )
458 end do
459
460 !--
461 do p=1, elem%Np
462 gsqrt_ = lmesh%Gsqrt(p,ke)
463 flux(p,1,dens_vid) = gsqrt_ * momx_(p,ke)
464 flux(p,2,dens_vid) = gsqrt_ * momy_(p,ke)
465 flux(p,3,dens_vid) = gsqrt_ * ( &
466 momz_(p,ke) * rgsqrtv(p) &
467 + lmesh%GI3(p,ke,1) * momx_(p,ke) &
468 + lmesh%GI3(p,ke,2) * momy_(p,ke) )
469 end do
470 do p=1, elem%Np
471 pt_ = ( therm_hyd(p,ke) + drhot_(p,ke) ) * rdens_(p)
472 flux(p,1,rhot_vid) = flux(p,1,dens_vid) * pt_
473 flux(p,2,rhot_vid) = flux(p,2,dens_vid) * pt_
474! Flux(p,3,RHOT_VID) = 0.0_RP ! for HEVI
475
476 w_ = momz_(p,ke) * rdens_(p)
477 flux(p,1,momz_vid) = flux(p,1,dens_vid) * w_
478 flux(p,2,momz_vid) = flux(p,2,dens_vid) * w_
479 flux(p,3,momz_vid) = flux(p,3,dens_vid) * w_
480 end do
481 do p=1, elem%Np
482 gsqrtdpres_ = lmesh%Gsqrt(p,ke) * dpres_(p,ke)
483
484 u_ = momx_(p,ke) * rdens_(p)
485 v_ = momy_(p,ke) * rdens_(p)
486
487 flux(p,1,momx_vid) = flux(p,1,dens_vid) * u_ + g11(p) * gsqrtdpres_
488 flux(p,2,momx_vid) = flux(p,2,dens_vid) * u_ + g12(p) * gsqrtdpres_
489 flux(p,3,momx_vid) = flux(p,3,dens_vid) * u_ + gsqrtdpres_ * ( g11(p) * lmesh%GI3(p,ke,1) + g12(p) * lmesh%GI3(p,ke,2) )
490
491 flux(p,1,momy_vid) = flux(p,1,dens_vid) * v_ + g12(p) * gsqrtdpres_
492 flux(p,2,momy_vid) = flux(p,2,dens_vid) * v_ + g22(p) * gsqrtdpres_
493 flux(p,3,momy_vid) = flux(p,3,dens_vid) * v_ + gsqrtdpres_ * ( g12(p) * lmesh%GI3(p,ke,1) + g22(p) * lmesh%GI3(p,ke,2) )
494 end do
495
496 call element3d_operation%Div_var5( &
497 flux, del_flux(:,:,ke), & ! (in)
498 dflux ) ! (out)
499
500 do p=1, elem%Np
501 e11 = lmesh%Escale(p,ke,1,1)
502 e22 = lmesh%Escale(p,ke,2,2)
503 e33 = lmesh%Escale(p,ke,3,3)
504
505 dens_dt(p,ke) = - ( &
506 e11 * dflux(p,1,dens_vid) &
507 + e22 * dflux(p,2,dens_vid) &
508 + dflux(p,4,dens_vid) ) * rgsqrt(p)
509
510 rhot_dt(p,ke) = - ( &
511 e11 * dflux(p,1,rhot_vid) &
512 + e22 * dflux(p,2,rhot_vid) &
513 + dflux(p,4,rhot_vid) ) * rgsqrt(p)
514 end do
515
516 do p=1, elem%Np
517 e11 = lmesh%Escale(p,ke,1,1)
518 e22 = lmesh%Escale(p,ke,2,2)
519 e33 = lmesh%Escale(p,ke,3,3)
520
521 momz_dt(p,ke) = - ( &
522 e11 * dflux(p,1,momz_vid) &
523 + e22 * dflux(p,2,momz_vid) &
524 + e33 * dflux(p,3,momz_vid) &
525 + dflux(p,4,momz_vid) ) * rgsqrt(p)
526 end do
527
528 !--
529 do p=1, elem%Np
530 x(p) = x2d(elem%IndexH2Dto3D(p),ke2d)
531 y(p) = y2d(elem%IndexH2Dto3D(p),ke2d)
532 twoovdel2(p) = 2.0_rp / ( 1.0_rp + x(p)**2 + y(p)**2 )
533 end do
534
535 do p=1, elem%Np
536 cori(p,1) = s * ohm * twoovdel2(p) * ( - x(p) * y(p) * momx_(p,ke) + ( 1.0_rp + y(p)**2 ) * momy_(p,ke) )
537 cori(p,2) = s * ohm * twoovdel2(p) * ( - ( 1.0_rp + x(p)**2 ) * momx_(p,ke) + x(p) * y(p) * momy_(p,ke) )
538 end do
539 if ( is_panel1to4 ) then
540 do p=1, elem%Np
541 cori(p,1) = s * y(p) * cori(p,1)
542 cori(p,2) = s * y(p) * cori(p,2)
543 end do
544 end if
545
546 do p=1, elem%Np
547 u_ = momx_(p,ke) * rdens_(p)
548 v_ = momy_(p,ke) * rdens_(p)
549
550 momx_dt(p,ke) = - ( g11(p) * dphyddx(p,ke) + g12(p) * dphyddy(p,ke) ) &
551 - twoovdel2(p) * y(p) * &
552 ( x(p) * y(p) * u_ - (1.0_rp + y(p)**2) * v_ ) * momx_(p,ke) &
553 + cori(p,1)
554
555 momy_dt(p,ke) = - ( g12(p) * dphyddx(p,ke) + g22(p) * dphyddy(p,ke) ) &
556 - twoovdel2(p) * x(p) * &
557 ( - (1.0_rp + x(p)**2) * u_ + x(p) * y(p) * v_ ) * momy_(p,ke) &
558 + cori(p,2)
559 end do
560
561 do p=1, elem%Np
562 e11 = lmesh%Escale(p,ke,1,1)
563 e22 = lmesh%Escale(p,ke,2,2)
564 e33 = lmesh%Escale(p,ke,3,3)
565
566 momx_dt(p,ke) = momx_dt(p,ke) - ( &
567 e11 * dflux(p,1,momx_vid) &
568 + e22 * dflux(p,2,momx_vid) &
569 + e33 * dflux(p,3,momx_vid) &
570 + dflux(p,4,momx_vid) ) * rgsqrt(p)
571
572 momy_dt(p,ke) = momy_dt(p,ke) - ( &
573 e11 * dflux(p,1,momy_vid) &
574 + e22 * dflux(p,2,momy_vid) &
575 + e33 * dflux(p,3,momy_vid) &
576 + dflux(p,4,momy_vid) ) * rgsqrt(p)
577 end do
578 end do
579 !$omp end parallel
580 call prof_rapend('cal_dyn_tend_interior', 3)
581
582 return
subroutine, public atm_dyn_dgm_nonhydro3d_rhot_hevi_numflux_get_generalhvc(del_flux, ddens_, momx_, momy_, momz_, drhot_, dpres, dens_hyd, pres_hyd, therm_hyd, rtot, cvtot, cptot, gsqrt, g11, g12, g22, gsqrth, gam, g13, g23, nx, ny, nz, vmapm, vmapp, im2dto3d, lmesh, elem, lmesh2d, elem2d)

References scale_atm_dyn_dgm_nonhydro3d_rhot_hevi_numflux::atm_dyn_dgm_nonhydro3d_rhot_hevi_numflux_get_generalhvc().

◆ atm_dyn_dgm_globalnonhydro3d_rhot_hevi_cal_vi_asis()

subroutine, public scale_atm_dyn_dgm_globalnonhydro3d_rhot_hevi::atm_dyn_dgm_globalnonhydro3d_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 586 of file scale_atm_dyn_dgm_globalnonhydro3d_rhot_hevi.F90.

594
602
603 implicit none
604
605 class(LocalMesh3D), intent(in) :: lmesh
606 class(ElementBase3D), intent(in) :: elem
607 class(LocalMesh2D), intent(in) :: lmesh2D
608 class(ElementBase2D), intent(in) :: elem2D
609 real(RP), intent(out) :: DENS_dt(elem%Np,lmesh%NeA)
610 real(RP), intent(out) :: MOMX_dt(elem%Np,lmesh%NeA)
611 real(RP), intent(out) :: MOMY_dt(elem%Np,lmesh%NeA)
612 real(RP), intent(out) :: MOMZ_dt(elem%Np,lmesh%NeA)
613 real(RP), intent(out) :: RHOT_dt(elem%Np,lmesh%NeA)
614 real(RP), intent(in) :: DDENS_(elem%Np,lmesh%NeA)
615 real(RP), intent(in) :: MOMX_(elem%Np,lmesh%NeA)
616 real(RP), intent(in) :: MOMY_(elem%Np,lmesh%NeA)
617 real(RP), intent(in) :: MOMZ_(elem%Np,lmesh%NeA)
618 real(RP), intent(in) :: DRHOT_(elem%Np,lmesh%NeA)
619 real(RP), intent(in) :: DENS_hyd(elem%Np,lmesh%NeA)
620 real(RP), intent(in) :: PRES_hyd(elem%Np,lmesh%NeA)
621 real(RP), intent(in) :: DDENS0_(elem%Np,lmesh%NeA)
622 real(RP), intent(in) :: MOMX0_(elem%Np,lmesh%NeA)
623 real(RP), intent(in) :: MOMY0_(elem%Np,lmesh%NeA)
624 real(RP), intent(in) :: MOMZ0_(elem%Np,lmesh%NeA)
625 real(RP), intent(in) :: DRHOT0_(elem%Np,lmesh%NeA)
626 real(RP), intent(in) :: Rtot(elem%Np,lmesh%NeA)
627 real(RP), intent(in) :: CVtot(elem%Np,lmesh%NeA)
628 real(RP), intent(in) :: CPtot(elem%Np,lmesh%NeA)
629 class(ElementOperationBase3D), intent(in) :: element3D_operation
630 class(SparseMat), intent(in) :: Dz, Lift
631 real(RP), intent(in) :: impl_fac
632 real(RP), intent(in) :: dt
633
634 real(RP) :: PROG_VARS (elem%Np,lmesh%NeZ,PRGVAR_NUM,lmesh%NeX*lmesh%NeY)
635 real(RP) :: PROG_VARS0(elem%Np,lmesh%NeZ,PRGVAR_NUM,lmesh%NeX*lmesh%NeY)
636 real(RP) :: b1D(3,elem%Nnode_v,lmesh%NeZ,elem%Nnode_h1D**2,lmesh%NeX*lmesh%NeY)
637 integer :: ipiv(elem%Nnode_v*3*lmesh%NeZ,elem%Nnode_h1D**2)
638 real(RP) :: b1D_uv(elem%Nnode_v,lmesh%NeZ,2,elem%Nnode_h1D**2,lmesh%NeX*lmesh%NeY)
639 integer :: ipiv_uv(elem%Nnode_v*1*lmesh%NeZ,elem%Nnode_h1D**2)
640 real(RP) :: alph(elem%NfpTot,lmesh%NeZ,lmesh%NeX*lmesh%NeY)
641 real(RP) :: Rtot_z(elem%Np,lmesh%NeZ,lmesh%NeX*lmesh%NeY)
642 real(RP) :: CPtot_ov_CVtot(elem%Np,lmesh%NeZ,lmesh%NeX*lmesh%NeY)
643 real(RP) :: DENS_hyd_z(elem%Np,lmesh%NeZ,lmesh%NeX*lmesh%NeY)
644 real(RP) :: PRES_hyd_z(elem%Np,lmesh%NeZ,lmesh%NeX*lmesh%NeY)
645 real(RP) :: GnnM_z(elem%Np,lmesh%NeZ,lmesh%NeX*lmesh%NeY)
646 real(RP) :: G13_z(elem%Np,lmesh%NeZ,lmesh%NeX*lmesh%NeY)
647 real(RP) :: G23_z(elem%Np,lmesh%NeZ,lmesh%NeX*lmesh%NeY)
648 real(RP) :: GsqrtV_z(elem%Np,lmesh%NeZ,lmesh%NeX*lmesh%NeY)
649 real(RP) :: nz(elem%NfpTot,lmesh%NeZ,lmesh%NeX*lmesh%NeY)
650 integer :: vmapM(elem%NfpTot,lmesh%NeZ)
651 integer :: vmapP(elem%NfpTot,lmesh%NeZ)
652 integer :: ColMask(elem%Nnode_v)
653 integer :: ke_xy, ke_z, ke, ke2D, v
654 integer :: itr_nlin
655 integer :: kl, ku, nz_1D
656 integer :: kl_uv, ku_uv, nz_1D_uv
657 integer :: ij
658 logical :: is_converged
659
660 real(RP), allocatable :: PmatBnd(:,:,:)
661 real(RP), allocatable :: PmatBnd_uv(:,:,:)
662 integer :: info, info_uv
663 !------------------------------------------------------------------------
664
665 call prof_rapstart( 'hevi_cal_vi_prep', 3)
666
667 nz_1d = elem%Nnode_v * 3 * lmesh%NeZ
668 kl = ( elem%Nnode_v + 1 ) * 3 - 1
669 ku = kl
670 nz_1d_uv = elem%Nnode_v * 1 * lmesh%NeZ
671 kl_uv = elem%Nnode_v
672 ku_uv = kl_uv
673 allocate( pmatbnd(2*kl+ku+1,nz_1d,elem%Nnode_h1D**2) )
674 allocate( pmatbnd_uv(2*kl_uv+ku_uv+1,nz_1d_uv,elem%Nnode_h1D**2) )
675
676 call lmesh%GetVmapZ1D( vmapm, vmapp ) ! (out)
677
678 !-
679
680 !$omp parallel private( ke_xy, ke_z, ke, ke2D )
681 !$omp do collapse(2)
682 do ke_xy=1, lmesh%NeX*lmesh%NeY
683 do ke_z=1, lmesh%NeZ
684 ke = ke_xy + (ke_z-1)*lmesh%NeX*lmesh%NeY
685 ke2d = lmesh%EMap3Dto2D(ke)
686
687 prog_vars(:,ke_z,dens_vid,ke_xy) = ddens0_(:,ke)
688 prog_vars(:,ke_z,momx_vid,ke_xy) = momx0_(:,ke)
689 prog_vars(:,ke_z,momy_vid,ke_xy) = momy0_(:,ke)
690 prog_vars(:,ke_z,momz_vid,ke_xy) = momz0_(:,ke)
691 prog_vars(:,ke_z,rhot_vid,ke_xy) = drhot0_(:,ke)
692
693 dens_hyd_z(:,ke_z,ke_xy) = dens_hyd(:,ke)
694 pres_hyd_z(:,ke_z,ke_xy) = pres_hyd(:,ke)
695
696 rtot_z(:,ke_z,ke_xy) = rtot(:,ke)
697 cptot_ov_cvtot(:,ke_z,ke_xy) = cptot(:,ke) / cvtot(:,ke)
698
699 nz(:,ke_z,ke_xy) = lmesh%normal_fn(:,ke,3)
700 g13_z(:,ke_z,ke_xy) = lmesh%GI3(:,ke,1)
701 g23_z(:,ke_z,ke_xy) = lmesh%GI3(:,ke,2)
702 gsqrtv_z(:,ke_z,ke_xy) = lmesh%Gsqrt(:,ke) / lmesh%GsqrtH(elem%IndexH2Dto3D,ke2d)
703
704 gnnm_z(:,ke_z,ke_xy) = ( &
705 1.0_rp / gsqrtv_z(:,ke_z,ke_xy)**2 &
706 + g13_z(:,ke_z,ke_xy) * ( lmesh%GIJ(elem%IndexH2Dto3D,ke2d,1,1) * g13_z(:,ke_z,ke_xy) &
707 + lmesh%GIJ(elem%IndexH2Dto3D,ke2d,1,2) * g23_z(:,ke_z,ke_xy) ) &
708 + g23_z(:,ke_z,ke_xy) * ( lmesh%GIJ(elem%IndexH2Dto3D,ke2d,1,2) * g13_z(:,ke_z,ke_xy) &
709 + lmesh%GIJ(elem%IndexH2Dto3D,ke2d,2,2) * g23_z(:,ke_z,ke_xy) ) )
710 end do
711 end do
712 !$omp end do
713 !$omp workshare
714 prog_vars0(:,:,:,:) = prog_vars(:,:,:,:)
715 !$omp end workshare
716 !$omp end parallel
717
718 call prof_rapend( 'hevi_cal_vi_prep', 3)
719
720 !--
721
722 if ( abs(impl_fac) > 0.0_rp ) then
723 call prof_rapstart( 'hevi_cal_vi_itr', 3)
724
725 ! G = (q^n+1 - q^n*) + impl_fac * A(q^n+1) = 0
726 ! dG/dq^n+1 del[q] = - G(q^n*)
727 do itr_nlin = 1, 1
728 call prof_rapstart( 'hevi_cal_vi_ax_uv', 3)
729
730 call vi_eval_ax_uv( &
731 momx_dt(:,:), momy_dt(:,:), alph(:,:,:), & ! (out)
732 prog_vars, prog_vars0, & ! (in)
733 ddens_, momx_, momy_, momz_, drhot_, & ! (in)
734 dens_hyd_z, pres_hyd_z, & ! (in)
735 rtot_z, cptot_ov_cvtot, & ! (in)
736 dz, lift, intrpmat_vpordm1, & ! (in)
737 gnnm_z, g13_z, g23_z, gsqrtv_z, & ! (in)
738 impl_fac, dt, & ! (in)
739 lmesh, elem, nz, vmapm, vmapp, & ! (in)
740 b1d_uv(:,:,:,:,:) ) ! (out)
741
742 call prof_rapend( 'hevi_cal_vi_ax_uv', 3)
743
744 do ke_xy=1, lmesh%NeX * lmesh%NeY
745 call prof_rapstart( 'hevi_cal_vi_matbnd_uv', 3)
746
747 call vi_construct_matbnd_uv( pmatbnd_uv(:,:,:), & ! (out)
748 kl_uv, ku_uv, nz_1d_uv, & ! (in)
749 prog_vars(:,:,:,ke_xy), & ! (in)
750 dens_hyd_z(:,:,ke_xy), pres_hyd_z(:,:,ke_xy), & ! (in)
751 g13_z(:,:,ke_xy), g23_z(:,:,ke_xy), gsqrtv_z(:,:,ke_xy), & ! (in)
752 alph(:,:,ke_xy), & ! (in)
753 rtot_z(:,:,ke_xy), cptot_ov_cvtot(:,:,ke_xy), & ! (in)
754 dz, lift, intrpmat_vpordm1, & ! (in)
755 impl_fac, dt, & ! (in)
756 lmesh, elem, nz(:,:,ke_xy), vmapm, vmapp, ke_xy, 1 ) ! (in)
757
758 call prof_rapend( 'hevi_cal_vi_matbnd_uv', 3)
759
760 call prof_rapstart( 'hevi_cal_vi_lin_uv', 3)
761 !$omp parallel private(ij, v, ke_z, info_uv, ColMask)
762 !$omp do
763 do ij=1, elem%Nnode_h1D**2
764 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 )
765
766 colmask(:) = elem%Colmask(:,ij)
767 do ke_z=1, lmesh%NeZ
768 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)
769 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)
770 end do
771 end do ! for ij
772 !$omp end do
773 !$omp end parallel
774 call prof_rapend( 'hevi_cal_vi_lin_uv', 3)
775 end do ! for ke_xy
776
777 call prof_rapstart( 'hevi_cal_vi_ax', 3)
778 call vi_eval_ax( &
779 dens_dt(:,:), momz_dt(:,:), rhot_dt(:,:), & ! (out, dummy)
780 alph(:,:,:), & ! (in)
781 prog_vars, prog_vars0, & ! (in)
782 ddens_, momx_, momy_, momz_, drhot_, & ! (in)
783 dens_hyd_z, pres_hyd_z, & ! (in)
784 rtot_z, cptot_ov_cvtot, & ! (in)
785 dz, lift, intrpmat_vpordm1, & ! (in)
786 gnnm_z, g13_z, g23_z, gsqrtv_z, & ! (in)
787 impl_fac, dt, & ! (in)
788 lmesh, elem, nz, vmapm, vmapp, & ! (in)
789 b1d(:,:,:,:,:) ) ! (out)
790 call prof_rapend( 'hevi_cal_vi_ax', 3)
791
792 do ke_xy=1, lmesh%NeX * lmesh%NeY
793 call prof_rapstart( 'hevi_cal_vi_matbnd', 3)
794 call vi_construct_matbnd( pmatbnd(:,:,:), & ! (out)
795 kl, ku, nz_1d, & ! (in)
796 prog_vars(:,:,:,ke_xy), & ! (in)
797 dens_hyd_z(:,:,ke_xy), pres_hyd_z(:,:,ke_xy), & ! (in)
798 g13_z(:,:,ke_xy), g23_z(:,:,ke_xy), gsqrtv_z(:,:,ke_xy), & ! (in)
799 alph(:,:,ke_xy), & ! (in)
800 rtot_z(:,:,ke_xy), cptot_ov_cvtot(:,:,ke_xy), & ! (in)
801 dz, lift, intrpmat_vpordm1, & ! (in)
802 impl_fac, dt, & ! (in)
803 lmesh, elem, nz(:,:,ke_xy), vmapm, vmapp, ke_xy, 1 ) ! (in)
804
805 call prof_rapend( 'hevi_cal_vi_matbnd', 3)
806
807 call prof_rapstart( 'hevi_cal_vi_lin', 3)
808 !$omp parallel private(ij, v, ke_z, info, ColMask)
809 !$omp do
810 do ij=1, elem%Nnode_h1D**2
811 call linalgebra_solvelineq_bndmat( pmatbnd(:,:,ij), b1d(:,:,:,ij,ke_xy), ipiv(:,ij), nz_1d, kl, ku, 1, vi_use_lapack_flag )
812
813 colmask(:) = elem%Colmask(:,ij)
814 do ke_z=1, lmesh%NeZ
815 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)
816 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)
817 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)
818 end do
819 end do ! for ij
820 !$omp end do
821 !$omp end parallel
822 call prof_rapend( 'hevi_cal_vi_lin', 3)
823
824 end do ! for ke_xy
825 end do ! itr nlin
826
827 call prof_rapend( 'hevi_cal_vi_itr', 3)
828 end if
829
830 call prof_rapstart( 'hevi_cal_vi_retrun_var', 3)
831 if ( abs(impl_fac) > 0.0_rp) then
832 !$omp parallel do collapse(2) private(ke_xy, ke_z, ke)
833 do ke_xy=1, lmesh%NeX * lmesh%NeY
834 do ke_z=1, lmesh%NeZ
835 ke = ke_xy + (ke_z-1)*lmesh%NeX*lmesh%NeY
836 dens_dt(:,ke) = ( prog_vars(:,ke_z,dens_vid,ke_xy) - ddens_(:,ke) ) / impl_fac
837 momx_dt(:,ke) = ( prog_vars(:,ke_z,momx_vid,ke_xy) - momx_(:,ke) ) / impl_fac
838 momy_dt(:,ke) = ( prog_vars(:,ke_z,momy_vid,ke_xy) - momy_(:,ke) ) / impl_fac
839 momz_dt(:,ke) = ( prog_vars(:,ke_z,momz_vid,ke_xy) - momz_(:,ke) ) / impl_fac
840 rhot_dt(:,ke) = ( prog_vars(:,ke_z,rhot_vid,ke_xy) - drhot_(:,ke) ) / impl_fac
841 end do
842 end do
843 else
844 call vi_eval_ax_uv( &
845 momx_dt(:,:), momy_dt(:,:), & ! (out)
846 alph(:,:,:), & ! (out, dummy)
847 prog_vars, prog_vars0, & ! (in)
848 ddens_, momx_, momy_, momz_, drhot_, & ! (in)
849 dens_hyd_z, pres_hyd_z, & ! (in)
850 rtot_z, cptot_ov_cvtot, & ! (in)
851 dz, lift, intrpmat_vpordm1, & ! (in)
852 gnnm_z, g13_z, g23_z, gsqrtv_z, & ! (in)
853 impl_fac, dt, & ! (in)
854 lmesh, elem, nz, vmapm, vmapp ) ! (in)
855
856 call vi_eval_ax( &
857 dens_dt(:,:), momz_dt(:,:), rhot_dt(:,:), & ! (out)
858 alph(:,:,:), & ! (in, dummy)
859 prog_vars, prog_vars0, & ! (in)
860 ddens_, momx_, momy_, momz_, drhot_, & ! (in)
861 dens_hyd_z, pres_hyd_z, & ! (in)
862 rtot_z, cptot_ov_cvtot, & ! (in)
863 dz, lift, intrpmat_vpordm1, & ! (in)
864 gnnm_z, g13_z, g23_z, gsqrtv_z, & ! (in)
865 impl_fac, dt, & ! (in)
866 lmesh, elem, nz, vmapm, vmapp ) ! (in)
867 end if
868 call prof_rapend( 'hevi_cal_vi_retrun_var', 3)
869 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_globalnonhydro3d_rhot_hevi_cal_vi()

subroutine, public scale_atm_dyn_dgm_globalnonhydro3d_rhot_hevi::atm_dyn_dgm_globalnonhydro3d_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 873 of file scale_atm_dyn_dgm_globalnonhydro3d_rhot_hevi.F90.

881
889
890 implicit none
891
892 class(LocalMesh3D), intent(in) :: lmesh
893 class(ElementBase3D), intent(in) :: elem
894 class(LocalMesh2D), intent(in) :: lmesh2D
895 class(ElementBase2D), intent(in) :: elem2D
896 real(RP), intent(out) :: DENS_dt(elem%Np,lmesh%NeA)
897 real(RP), intent(out) :: MOMX_dt(elem%Np,lmesh%NeA)
898 real(RP), intent(out) :: MOMY_dt(elem%Np,lmesh%NeA)
899 real(RP), intent(out) :: MOMZ_dt(elem%Np,lmesh%NeA)
900 real(RP), intent(out) :: RHOT_dt(elem%Np,lmesh%NeA)
901 real(RP), intent(in) :: DDENS_(elem%Np,lmesh%NeA)
902 real(RP), intent(in) :: MOMX_(elem%Np,lmesh%NeA)
903 real(RP), intent(in) :: MOMY_(elem%Np,lmesh%NeA)
904 real(RP), intent(in) :: MOMZ_(elem%Np,lmesh%NeA)
905 real(RP), intent(in) :: DRHOT_(elem%Np,lmesh%NeA)
906 real(RP), intent(in) :: DENS_hyd(elem%Np,lmesh%NeA)
907 real(RP), intent(in) :: PRES_hyd(elem%Np,lmesh%NeA)
908 real(RP), intent(in) :: DDENS0_(elem%Np,lmesh%NeA)
909 real(RP), intent(in) :: MOMX0_(elem%Np,lmesh%NeA)
910 real(RP), intent(in) :: MOMY0_(elem%Np,lmesh%NeA)
911 real(RP), intent(in) :: MOMZ0_(elem%Np,lmesh%NeA)
912 real(RP), intent(in) :: DRHOT0_(elem%Np,lmesh%NeA)
913 real(RP), intent(in) :: Rtot(elem%Np,lmesh%NeA)
914 real(RP), intent(in) :: CVtot(elem%Np,lmesh%NeA)
915 real(RP), intent(in) :: CPtot(elem%Np,lmesh%NeA)
916 class(ElementOperationBase3D), intent(in) :: element3D_operation
917 class(SparseMat), intent(in) :: Dz, Lift
918 real(RP), intent(in) :: impl_fac
919 real(RP), intent(in) :: dt
920
921 real(RP) :: PROG_VARS (elem%Np,lmesh%NeX*lmesh%NeY,lmesh%NeZ,PRGVAR_NUM)
922 real(RP) :: PROG_VARS0(elem%Np,lmesh%NeX*lmesh%NeY,lmesh%NeZ,PRGVAR_NUM)
923 real(RP) :: alph(elem%NfpTot,lmesh%Ne)
924 real(RP) :: GsqrtV(elem%Np,lmesh%Ne)
925
926 integer :: vmapM(elem%NfpTot,lmesh%Ne)
927 integer :: vmapP(elem%NfpTot,lmesh%Ne)
928 integer :: ke_xy, ke_z, ke, ke2d
929 integer :: p
930 integer :: itr_nlin
931 integer :: ij, i, j, im, jm
932 logical :: is_converged
933
934 real(RP), allocatable :: b_uv(:,:,:,:)
935 real(RP), allocatable :: b (:,:,:)
936
937 real(RP) :: DENS(elem%Np,lmesh%Ne)
938 real(RP) :: W(elem%Np,lmesh%Ne)
939 real(RP) :: WT(elem%Np,lmesh%Ne)
940 real(RP) :: POT(elem%Np,lmesh%Ne)
941 real(RP) :: DPDRHOT(elem%Np,lmesh%Ne)
942 !------------------------------------------------------------------------
943
944 call prof_rapstart( 'hevi_cal_vi_prep', 3)
945
946 call lmesh%GetVmapZ3D( vmapm, vmapp ) ! (out)
948 elem%Nnode_h1D ) ! (in)
949
950 !-
951 !$omp parallel private( ke_xy, ke_z, ke, ke2D, p )
952 !$omp do collapse(2)
953 do ke_z =1, lmesh%NeZ
954 do ke_xy=1, lmesh%NeX * lmesh%NeY
955 ke = ke_xy + (ke_z-1)*lmesh%NeX*lmesh%NeY
956 ke2d = lmesh%EMap3Dto2D(ke)
957
958 prog_vars(:,ke_xy,ke_z,dens_vid) = ddens0_(:,ke)
959 prog_vars(:,ke_xy,ke_z,momx_vid) = momx0_(:,ke)
960 prog_vars(:,ke_xy,ke_z,momy_vid) = momy0_(:,ke)
961 prog_vars(:,ke_xy,ke_z,momz_vid) = momz0_(:,ke)
962 prog_vars(:,ke_xy,ke_z,rhot_vid) = drhot0_(:,ke)
963
964 do p=1, elem%Np
965 gsqrtv(p,ke) = lmesh%Gsqrt(p,ke) / ( lmesh%gam(p,ke)**2 * lmesh%GsqrtH(elem%IndexH2Dto3D(p),ke2d) )
966 end do
967 end do
968 end do
969 !$omp end do
970 !$omp workshare
971 prog_vars0(:,:,:,:) = prog_vars(:,:,:,:)
972 !$omp end workshare
973 !$omp end parallel
974 call prof_rapend( 'hevi_cal_vi_prep', 3)
975
976 !--
977
978 if ( abs(impl_fac) > 0.0_rp ) then
979 call prof_rapstart( 'hevi_cal_vi_itr', 3)
980
981 allocate( b_uv(im*elem%Nnode_v,2,jm,lmesh%Ne) )
982 allocate( b(im*3*elem%Nnode_v,jm,lmesh%Ne) )
983
984 ! G = (q^n+1 - q^n*) + impl_fac * A(q^n+1) = 0
985 ! dG/dq^n+1 del[q] = - G(q^n*)
986 do itr_nlin = 1, 1
987 call prof_rapstart( 'hevi_cal_vi_ax_uv', 3)
988 call vi_eval_ax_uv( &
989 momx_dt(:,:), momy_dt(:,:), alph(:,:), & ! (out)
990 prog_vars, prog_vars0, & ! (in)
991 ddens_, momx_, momy_, momz_, drhot_, & ! (in)
992 dens_hyd, pres_hyd, & ! (in)
993 rtot, cptot, cvtot, gsqrtv, & ! (in)
994 impl_fac, dt, & ! (in)
995 lmesh, elem, vmapm, vmapp, & ! (in)
996 element3d_operation, & ! (in)
997 im, jm, & ! (in)
998 b_uv ) ! (out)
999 call prof_rapend( 'hevi_cal_vi_ax_uv', 3)
1000
1001 call vi_solve_uv( prog_vars, & ! (inout)
1002 b_uv, alph, & ! (in)
1003 gsqrtv, & ! (in)
1004 impl_fac, dt, & ! (in)
1005 vmapm, vmapp, lmesh, elem, im, jm ) ! (in)
1006
1007 !----
1008
1009 call prof_rapstart( 'hevi_cal_vi_ax', 3)
1010 call vi_eval_ax( &
1011 dens_dt(:,:), momz_dt(:,:), rhot_dt(:,:), & ! (out)
1012 alph(:,:), & ! (in)
1013 prog_vars, prog_vars0, ddens_, momx_, momy_, momz_, drhot_, & ! (in)
1014 dens_hyd, pres_hyd, rtot, cptot, cvtot, gsqrtv, & ! (in)
1015 impl_fac, dt, lmesh, elem, vmapm, vmapp, & ! (in)
1016 element3d_operation, & ! (in)
1017 im, jm, & ! (in)
1018 b, dens, w, wt, pot, dpdrhot ) ! (out)
1019 call prof_rapend( 'hevi_cal_vi_ax', 3)
1020
1021 call vi_solve( prog_vars, & ! (inout)
1022 b, dens, w, wt, pot, dpdrhot, alph, & ! (in)
1023 gsqrtv, intrpmat_vpordm1, & ! (in)
1024 impl_fac, dt, & ! (in)
1025 vmapm, vmapp, lmesh, elem, im, jm ) ! (in)
1026 end do ! itr nlin
1027
1028 call prof_rapend( 'hevi_cal_vi_itr', 3)
1029 end if
1030
1031 call prof_rapstart( 'hevi_cal_vi_retrun_var', 3)
1032 if ( abs(impl_fac) > 0.0_rp) then
1033 !$omp parallel do collapse(2) private(ke_xy, ke_z, ke)
1034 do ke_z=1, lmesh%NeZ
1035 do ke_xy=1, lmesh%NeX * lmesh%NeY
1036 ke = ke_xy + (ke_z-1)*lmesh%NeX*lmesh%NeY
1037 dens_dt(:,ke) = ( prog_vars(:,ke_xy,ke_z,dens_vid) - ddens_(:,ke) ) / impl_fac
1038 momx_dt(:,ke) = ( prog_vars(:,ke_xy,ke_z,momx_vid) - momx_(:,ke) ) / impl_fac
1039 momy_dt(:,ke) = ( prog_vars(:,ke_xy,ke_z,momy_vid) - momy_(:,ke) ) / impl_fac
1040 momz_dt(:,ke) = ( prog_vars(:,ke_xy,ke_z,momz_vid) - momz_(:,ke) ) / impl_fac
1041 rhot_dt(:,ke) = ( prog_vars(:,ke_xy,ke_z,rhot_vid) - drhot_(:,ke) ) / impl_fac
1042 end do
1043 end do
1044 else
1045 call vi_eval_ax_uv( &
1046 momx_dt(:,:), momy_dt(:,:), alph(:,:), & ! (out)
1047 prog_vars, prog_vars0, & ! (in)
1048 ddens_, momx_, momy_, momz_, drhot_, & ! (in)
1049 dens_hyd, pres_hyd, & ! (in)
1050 rtot, cptot, cvtot, gsqrtv, & ! (in)
1051 impl_fac, dt, & ! (in)
1052 lmesh, elem, vmapm, vmapp, & ! (in)
1053 element3d_operation, im, jm ) ! (in)
1054
1055 call vi_eval_ax( &
1056 dens_dt(:,:), momz_dt(:,:), rhot_dt(:,:), & ! (out, dummy)
1057 alph(:,:), & ! (in)
1058 prog_vars, prog_vars0, ddens_, momx_, momy_, momz_, drhot_, & ! (in)
1059 dens_hyd, pres_hyd, rtot, cptot, cvtot, gsqrtv, & ! (in)
1060 impl_fac, dt, lmesh, elem, vmapm, vmapp, & ! (in)
1061 element3d_operation, im, jm ) ! (in)
1062 end if
1063 call prof_rapend( 'hevi_cal_vi_retrun_var', 3)
1064
1065 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.