FE-Project
Loading...
Searching...
No Matches
Functions/Subroutines
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, 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, 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 73 of file scale_atm_dyn_dgm_nonhydro3d_rhot_heve.F90.

74 implicit none
75 class(MeshBase3D), intent(in) :: mesh
76 !--------------------------------------------
77
78 call atm_dyn_dgm_nonhydro3d_common_init( mesh )
79 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 83 of file scale_atm_dyn_dgm_nonhydro3d_rhot_heve.F90.

84 implicit none
85 !--------------------------------------------
86
87 call atm_dyn_dgm_nonhydro3d_common_final()
88 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(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 94 of file scale_atm_dyn_dgm_nonhydro3d_rhot_heve.F90.

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

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