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

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

Functions/Subroutines

subroutine, public atm_dyn_dgm_nonhydro3d_etot_heve_init (mesh)
subroutine, public atm_dyn_dgm_nonhydro3d_etot_heve_final ()
subroutine, public atm_dyn_dgm_nonhydro3d_etot_heve_cal_tend (dens_dt, momx_dt, momy_dt, momz_dt, etot_dt, ddens_, momx_, momy_, momz_, etot_, 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 (total energy conservation) equations.
Author
Yuta Kawai, Team SCALE

Function/Subroutine Documentation

◆ atm_dyn_dgm_nonhydro3d_etot_heve_init()

subroutine, public scale_atm_dyn_dgm_nonhydro3d_etot_heve::atm_dyn_dgm_nonhydro3d_etot_heve_init ( class(meshbase3d), intent(in) mesh)

Definition at line 71 of file scale_atm_dyn_dgm_nonhydro3d_etot_heve.F90.

72 implicit none
73 class(MeshBase3D), intent(in) :: mesh
74 !--------------------------------------------
75
76 call atm_dyn_dgm_nonhydro3d_common_init( mesh )
77 return

References scale_atm_dyn_dgm_nonhydro3d_common::atm_dyn_dgm_nonhydro3d_common_init().

◆ atm_dyn_dgm_nonhydro3d_etot_heve_final()

subroutine, public scale_atm_dyn_dgm_nonhydro3d_etot_heve::atm_dyn_dgm_nonhydro3d_etot_heve_final

Definition at line 81 of file scale_atm_dyn_dgm_nonhydro3d_etot_heve.F90.

82 implicit none
83 !--------------------------------------------
84
85 call atm_dyn_dgm_nonhydro3d_common_final()
86 return

References scale_atm_dyn_dgm_nonhydro3d_common::atm_dyn_dgm_nonhydro3d_common_final().

◆ atm_dyn_dgm_nonhydro3d_etot_heve_cal_tend()

subroutine, public scale_atm_dyn_dgm_nonhydro3d_etot_heve::atm_dyn_dgm_nonhydro3d_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) etot_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(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 92 of file scale_atm_dyn_dgm_nonhydro3d_etot_heve.F90.

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

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