FE-Project
Loading...
Searching...
No Matches
Functions/Subroutines
scale_atm_dyn_dgm_globalnonhydro3d_etot_heve Module Reference

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

Functions/Subroutines

subroutine, public atm_dyn_dgm_globalnonhydro3d_etot_heve_init (mesh)
 
subroutine, public atm_dyn_dgm_globalnonhydro3d_etot_heve_final ()
 
subroutine, public atm_dyn_dgm_globalnonhydro3d_etot_heve_cal_tend (dens_dt, momx_dt, momy_dt, momz_dt, entot_dt, ddens_, momx_, momy_, momz_, etot_, 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 / Global nonhydrostatic model / HEVE

Description
HEVE DGM scheme for Global Atmospheric Dynamical process. The governing equations is a fully compressible nonhydrostatic equations, which consist of mass, momentum, and thermodynamics (total energy conservation) equations.
Author
Yuta Kawai, Team SCALE

Function/Subroutine Documentation

◆ atm_dyn_dgm_globalnonhydro3d_etot_heve_init()

subroutine, public scale_atm_dyn_dgm_globalnonhydro3d_etot_heve::atm_dyn_dgm_globalnonhydro3d_etot_heve_init ( class(meshbase3d), intent(in) mesh)

Definition at line 73 of file scale_atm_dyn_dgm_globalnonhydro3d_etot_heve.F90.

74
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_globalnonhydro3d_etot_heve_final()

subroutine, public scale_atm_dyn_dgm_globalnonhydro3d_etot_heve::atm_dyn_dgm_globalnonhydro3d_etot_heve_final

Definition at line 84 of file scale_atm_dyn_dgm_globalnonhydro3d_etot_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_globalnonhydro3d_etot_heve_cal_tend()

subroutine, public scale_atm_dyn_dgm_globalnonhydro3d_etot_heve::atm_dyn_dgm_globalnonhydro3d_etot_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) entot_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) etot_,
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 95 of file scale_atm_dyn_dgm_globalnonhydro3d_etot_heve.F90.

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

References scale_atm_dyn_dgm_nonhydro3d_etot_heve_numflux::atm_dyn_dgm_nonhydro3d_etot_heve_numflux_get_generalhvc(), and scale_atm_dyn_dgm_nonhydro3d_common::intrpmat_vpordm1.