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

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

Functions/Subroutines

subroutine, public atm_dyn_dgm_nonhydro3d_rhot_heve_init (mesh)
subroutine, public atm_dyn_dgm_nonhydro3d_rhot_heve_final ()
subroutine, public atm_dyn_dgm_nonhydro3d_rhot_heve_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_heve_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_heve_cal_tend_cco (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)

Detailed Description

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

Description
HEVE 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, Xuanzhengbo Ren, and Team SCALE

Function/Subroutine Documentation

◆ atm_dyn_dgm_nonhydro3d_rhot_heve_init()

subroutine, public scale_atm_dyn_dgm_nonhydro3d_rhot_heve::atm_dyn_dgm_nonhydro3d_rhot_heve_init ( class(meshbase3d), intent(in) mesh)

Definition at line 74 of file scale_atm_dyn_dgm_nonhydro3d_rhot_heve.F90.

75 implicit none
76 class(MeshBase3D), intent(in) :: mesh
77 !--------------------------------------------
78
79 call atm_dyn_dgm_nonhydro3d_common_init( mesh )
80 return

References scale_atm_dyn_dgm_nonhydro3d_common::atm_dyn_dgm_nonhydro3d_common_init().

◆ atm_dyn_dgm_nonhydro3d_rhot_heve_final()

subroutine, public scale_atm_dyn_dgm_nonhydro3d_rhot_heve::atm_dyn_dgm_nonhydro3d_rhot_heve_final

Definition at line 84 of file scale_atm_dyn_dgm_nonhydro3d_rhot_heve.F90.

85 implicit none
86 !--------------------------------------------
87
88 call atm_dyn_dgm_nonhydro3d_common_final()
89 return

References scale_atm_dyn_dgm_nonhydro3d_common::atm_dyn_dgm_nonhydro3d_common_final().

◆ atm_dyn_dgm_nonhydro3d_rhot_heve_cal_tend_asis()

subroutine, public scale_atm_dyn_dgm_nonhydro3d_rhot_heve::atm_dyn_dgm_nonhydro3d_rhot_heve_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 95 of file scale_atm_dyn_dgm_nonhydro3d_rhot_heve.F90.

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

◆ atm_dyn_dgm_nonhydro3d_rhot_heve_cal_tend()

subroutine, public scale_atm_dyn_dgm_nonhydro3d_rhot_heve::atm_dyn_dgm_nonhydro3d_rhot_heve_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 292 of file scale_atm_dyn_dgm_nonhydro3d_rhot_heve.F90.

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

◆ atm_dyn_dgm_nonhydro3d_rhot_heve_cal_tend_cco()

subroutine, public scale_atm_dyn_dgm_nonhydro3d_rhot_heve::atm_dyn_dgm_nonhydro3d_rhot_heve_cal_tend_cco ( 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 493 of file scale_atm_dyn_dgm_nonhydro3d_rhot_heve.F90.

500
503
504 implicit none
505
506 class(LocalMesh3D), intent(in) :: lmesh
507 class(ElementBase3D), intent(in) :: elem
508 class(LocalMesh2D), intent(in) :: lmesh2D
509 class(ElementBase2D), intent(in) :: elem2D
510 class(ElementOperationBase3D), intent(in) :: element3D_operation
511 type(SparseMat), intent(in) :: Dx, Dy, Dz, Sx, Sy, Sz, Lift
512 real(RP), intent(out) :: DENS_dt(elem%Np,lmesh%NeA)
513 real(RP), intent(out) :: MOMX_dt(elem%Np,lmesh%NeA)
514 real(RP), intent(out) :: MOMY_dt(elem%Np,lmesh%NeA)
515 real(RP), intent(out) :: MOMZ_dt(elem%Np,lmesh%NeA)
516 real(RP), intent(out) :: RHOT_dt(elem%Np,lmesh%NeA)
517 real(RP), intent(in) :: DDENS_(elem%Np,lmesh%NeA)
518 real(RP), intent(in) :: MOMX_(elem%Np,lmesh%NeA)
519 real(RP), intent(in) :: MOMY_(elem%Np,lmesh%NeA)
520 real(RP), intent(in) :: MOMZ_(elem%Np,lmesh%NeA)
521 real(RP), intent(in) :: DRHOT_(elem%Np,lmesh%NeA)
522 real(RP), intent(in) :: DPRES_(elem%Np,lmesh%NeA)
523 real(RP), intent(in) :: DENS_hyd(elem%Np,lmesh%NeA)
524 real(RP), intent(in) :: PRES_hyd(elem%Np,lmesh%NeA)
525 real(RP), intent(in) :: PRES_hyd_ref(elem%Np,lmesh%NeA)
526 real(RP), intent(in) :: THERM_hyd(elem%Np,lmesh%NeA)
527 real(RP), intent(in) :: CORIOLIS(elem2D%Np,lmesh2D%NeA)
528 real(RP), intent(in) :: Rtot(elem%Np,lmesh%NeA)
529 real(RP), intent(in) :: CVtot(elem%Np,lmesh%NeA)
530 real(RP), intent(in) :: CPtot(elem%Np,lmesh%NeA)
531 real(RP), intent(in) :: DPhydDx(elem%Np,lmesh%NeA)
532 real(RP), intent(in) :: DPhydDy(elem%Np,lmesh%NeA)
533
534 real(RP) :: Flux(elem%Np,3,5), DFlux(elem%Np,3,5)
535 real(RP) :: u_, v_, w_, pt_
536 real(RP) :: cor
537 real(RP) :: drho(elem%Np)
538 real(RP) :: RDENS_(elem%Np), GsqrtV(elem%Np), RGsqrtV(elem%Np), RGsqrt(elem%Np)
539 real(RP) :: Gsqrt_, GsqrtDPRES_, E11, E22, E33
540
541 integer :: ke, ke2d
542 integer :: p
543
544 real(RP) :: gamm, rgamm
545 real(RP) :: rP0
546 real(RP) :: RovP0, P0ovR
547 !------------------------------------------------------------------------
548
549 call prof_rapstart('cal_dyn_tend_interior', 3)
550 gamm = cpdry / cvdry
551 rgamm = cvdry / cpdry
552 rp0 = 1.0_rp / pres00
553 rovp0 = rdry * rp0
554 p0ovr = pres00 / rdry
555
556 !$omp parallel do private( ke, ke2d, p, cor, &
557 !$omp u_, v_, w_, pt_, &
558 !$omp drho, &
559 !$omp RDENS_, RGsqrt, GsqrtV, RGsqrtV, Gsqrt_, GsqrtDPRES_, &
560 !$omp E11, E22, E33, DFlux, Flux )
561 do ke = lmesh%NeS, lmesh%NeE
562 !--
563 ke2d = lmesh%EMap3Dto2D(ke)
564
565 do p=1, elem%Np
566 gsqrtv(p) = lmesh%Gsqrt(p,ke) / lmesh%GsqrtH(elem%IndexH2Dto3D(p),ke2d)
567 rgsqrtv(p) = 1.0_rp / gsqrtv(p)
568 rgsqrt(p) = 1.0_rp / lmesh%Gsqrt(p,ke)
569 rdens_(p) = 1.0_rp / (ddens_(p,ke) + dens_hyd(p,ke))
570 end do
571
572 !--
573
574 do p=1, elem%Np
575 gsqrt_ = lmesh%Gsqrt(p,ke)
576 flux(p,1,dens_vid) = gsqrt_ * momx_(p,ke)
577 flux(p,2,dens_vid) = gsqrt_ * momy_(p,ke)
578 flux(p,3,dens_vid) = gsqrt_ * ( &
579 momz_(p,ke) * rgsqrtv(p) &
580 + lmesh%GI3(p,ke,1) * momx_(p,ke) &
581 + lmesh%GI3(p,ke,2) * momy_(p,ke) )
582 end do
583 do p=1, elem%Np
584 pt_ = ( therm_hyd(p,ke) + drhot_(p,ke) ) * rdens_(p)
585
586 flux(p,1,rhot_vid) = flux(p,1,dens_vid) * pt_
587 flux(p,2,rhot_vid) = flux(p,2,dens_vid) * pt_
588 flux(p,3,rhot_vid) = flux(p,3,dens_vid) * pt_
589
590 w_ = momz_(p,ke) * rdens_(p)
591 flux(p,1,momz_vid) = flux(p,1,dens_vid) * w_
592 flux(p,2,momz_vid) = flux(p,2,dens_vid) * w_
593 flux(p,3,momz_vid) = flux(p,3,dens_vid) * w_ + lmesh%Gsqrt(p,ke) * dpres_(p,ke) * rgsqrtv(p)
594 end do
595 do p=1, elem%Np
596 gsqrtdpres_ = lmesh%Gsqrt(p,ke) * dpres_(p,ke)
597
598 u_ = momx_(p,ke) * rdens_(p)
599 flux(p,1,momx_vid) = flux(p,1,dens_vid) * u_ + gsqrtdpres_
600 flux(p,2,momx_vid) = flux(p,2,dens_vid) * u_
601 flux(p,3,momx_vid) = flux(p,3,dens_vid) * u_ + gsqrtdpres_ * lmesh%GI3(p,ke,1)
602
603 v_ = momy_(p,ke) * rdens_(p)
604 flux(p,1,momy_vid) = flux(p,1,dens_vid) * v_
605 flux(p,2,momy_vid) = flux(p,2,dens_vid) * v_ + gsqrtdpres_
606 flux(p,3,momy_vid) = flux(p,3,dens_vid) * v_ + gsqrtdpres_ * lmesh%GI3(p,ke,2)
607 end do
608
609 call element3d_operation%Div_var5_2( &
610 flux, & ! (in)
611 dflux ) ! (out)
612
613 do p=1, elem%Np
614 e11 = lmesh%Escale(p,ke,1,1)
615 e22 = lmesh%Escale(p,ke,2,2)
616 e33 = lmesh%Escale(p,ke,3,3)
617
618 dens_dt(p,ke) = - ( &
619 e11 * dflux(p,1,dens_vid) &
620 + e22 * dflux(p,2,dens_vid) &
621 + e33 * dflux(p,3,dens_vid) &
622 ) * rgsqrt(p)
623
624 rhot_dt(p,ke) = - ( &
625 e11 * dflux(p,1,rhot_vid) &
626 + e22 * dflux(p,2,rhot_vid) &
627 + e33 * dflux(p,3,rhot_vid) &
628 ) * rgsqrt(p)
629 end do
630
631 call element3d_operation%VFilterPM1( ddens_(:,ke), & ! (in)
632 drho ) ! (out)
633
634 do p=1, elem%Np
635 e11 = lmesh%Escale(p,ke,1,1)
636 e22 = lmesh%Escale(p,ke,2,2)
637 e33 = lmesh%Escale(p,ke,3,3)
638
639 momz_dt(p,ke) = - ( &
640 e11 * dflux(p,1,momz_vid) &
641 + e22 * dflux(p,2,momz_vid) &
642 + e33 * dflux(p,3,momz_vid) &
643 ) * rgsqrt(p) &
644 - grav * drho(p)
645 end do
646
647 !--
648 do p=1, elem%Np
649 cor = coriolis(elem%IndexH2Dto3D(p),ke2d)
650 momx_dt(p,ke) = - dphyddx(p,ke) &
651 + cor * momy_(p,ke)
652 momy_dt(p,ke) = - dphyddy(p,ke) &
653 - cor * momx_(p,ke)
654 end do
655
656 do p=1, elem%Np
657 e11 = lmesh%Escale(p,ke,1,1)
658 e22 = lmesh%Escale(p,ke,2,2)
659 e33 = lmesh%Escale(p,ke,3,3)
660
661 momx_dt(p,ke) = momx_dt(p,ke) - ( &
662 e11 * dflux(p,1,momx_vid) &
663 + e22 * dflux(p,2,momx_vid) &
664 + e33 * dflux(p,3,momx_vid) &
665 ) * rgsqrt(p)
666
667 momy_dt(p,ke) = momy_dt(p,ke) - ( &
668 e11 * dflux(p,1,momy_vid) &
669 + e22 * dflux(p,2,momy_vid) &
670 + e33 * dflux(p,3,momy_vid) &
671 ) * rgsqrt(p)
672 end do
673 end do
674
675 call prof_rapend('cal_dyn_tend_interior', 3)
676
677 return

References scale_atm_dyn_dgm_nonhydro3d_rhot_heve_numflux::atm_dyn_dgm_nonhydro3d_rhot_heve_numflux_get_generalvc().