FE-Project
Loading...
Searching...
No Matches
Functions/Subroutines
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, 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(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.

98
101
102 implicit none
103
104 class(LocalMesh3D), intent(in) :: lmesh
105 class(ElementBase3D), intent(in) :: elem
106 class(LocalMesh2D), intent(in) :: lmesh2D
107 class(ElementBase2D), intent(in) :: elem2D
108 class(ElementOperationBase3D), intent(in) :: element3D_operation
109 type(SparseMat), intent(in) :: Dx, Dy, Dz, Sx, Sy, Sz, Lift
110 real(RP), intent(out) :: DENS_dt(elem%Np,lmesh%NeA)
111 real(RP), intent(out) :: MOMX_dt(elem%Np,lmesh%NeA)
112 real(RP), intent(out) :: MOMY_dt(elem%Np,lmesh%NeA)
113 real(RP), intent(out) :: MOMZ_dt(elem%Np,lmesh%NeA)
114 real(RP), intent(out) :: ETOT_dt(elem%Np,lmesh%NeA)
115 real(RP), intent(in) :: DDENS_(elem%Np,lmesh%NeA)
116 real(RP), intent(in) :: MOMX_(elem%Np,lmesh%NeA)
117 real(RP), intent(in) :: MOMY_(elem%Np,lmesh%NeA)
118 real(RP), intent(in) :: MOMZ_(elem%Np,lmesh%NeA)
119 real(RP), intent(in) :: ETOT_(elem%Np,lmesh%NeA)
120 real(RP), intent(in) :: DPRES_(elem%Np,lmesh%NeA)
121 real(RP), intent(in) :: DENS_hyd(elem%Np,lmesh%NeA)
122 real(RP), intent(in) :: PRES_hyd(elem%Np,lmesh%NeA)
123 real(RP), intent(in) :: PRES_hyd_ref(elem%Np,lmesh%NeA)
124 real(RP), intent(in) :: CORIOLIS(elem2D%Np,lmesh2D%NeA)
125 real(RP), intent(in) :: Rtot(elem%Np,lmesh%NeA)
126 real(RP), intent(in) :: CVtot(elem%Np,lmesh%NeA)
127 real(RP), intent(in) :: CPtot(elem%Np,lmesh%NeA)
128 real(RP), intent(in) :: DPhydDx(elem%Np,lmesh%NeA)
129 real(RP), intent(in) :: DPhydDy(elem%Np,lmesh%NeA)
130
131 real(RP) :: Fx(elem%Np), Fy(elem%Np), Fz(elem%Np), LiftDelFlx(elem%Np)
132 real(RP) :: DPRES_hyd(elem%Np), GradPhyd_x(elem%Np), GradPhyd_y(elem%Np)
133 real(RP) :: del_flux(elem%NfpTot,lmesh%Ne,PRGVAR_NUM)
134 real(RP) :: del_flux_hyd(elem%NfpTot,lmesh%Ne,2)
135 real(RP) :: rdens_(elem%Np), u_(elem%Np), v_(elem%Np), w_(elem%Np), wt_(elem%Np)
136 real(RP) :: Cori(elem%Np)
137 real(RP) :: drho(elem%Np)
138 real(RP) :: GsqrtV(elem%Np), RGsqrtV(elem%Np)
139
140 integer :: ke, ke2d
141
142 real(RP) :: gamm, rgamm
143 real(RP) :: rP0
144 real(RP) :: RovP0, P0ovR
145
146 real(RP) :: enthalpy_(elem%Np)
147 !------------------------------------------------------------------------
148
149 call prof_rapstart('cal_dyn_tend_bndflux', 3)
150 call get_ebnd_flux( &
151 del_flux, del_flux_hyd, & ! (out)
152 ddens_, momx_, momy_, momz_, etot_, dpres_, dens_hyd, pres_hyd, & ! (in)
153 rtot, cvtot, cptot, & ! (in)
154 lmesh%Gsqrt, lmesh%GI3(:,:,1), lmesh%GI3(:,:,2), lmesh%zlev(:,:), & ! (in)
155 lmesh%normal_fn(:,:,1), lmesh%normal_fn(:,:,2), lmesh%normal_fn(:,:,3), & ! (in)
156 lmesh%vmapM, lmesh%vmapP, & ! (in)
157 lmesh, elem, lmesh2d, elem2d ) ! (in)
158 call prof_rapend('cal_dyn_tend_bndflux', 3)
159
160 !-----
161 call prof_rapstart('cal_dyn_tend_interior', 3)
162 gamm = cpdry / cvdry
163 rgamm = cvdry / cpdry
164 rp0 = 1.0_rp / pres00
165 rovp0 = rdry * rp0
166 p0ovr = pres00 / rdry
167
168 !$omp parallel do private( ke, ke2d, Cori, &
169 !$omp rdens_, u_, v_, w_, wt_, &
170 !$omp enthalpy_, drho, &
171 !$omp 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 rdens_(:) = 1.0_rp / (ddens_(:,ke) + dens_hyd(:,ke))
184 u_(:) = momx_(:,ke) * rdens_(:)
185 v_(:) = momy_(:,ke) * rdens_(:)
186 w_(:) = momz_(:,ke) * rdens_(:)
187 wt_(:) = w_(:) * rgsqrtv(:) + lmesh%GI3(:,ke,1) * u_(:) + lmesh%GI3(:,ke,2) * v_(:)
188
189 ! DPRES_(:) = ( CPtot(:,ke) / CVtot(:,ke) - 1.0_RP ) &
190 ! * ( ETOT_(:,ke) - Grav * ( DDENS_(:,ke) + DENS_hyd(:,ke) ) * lmesh%zlev(:,ke) &
191 ! - 0.5_RP * ( MOMX_(:,ke) * u_(:) + MOMY_(:,ke) * v_(:) + MOMZ_(:,ke) * w_(:) ) ) &
192 ! - PRES_hyd(:,ke)
193
194 enthalpy_(:) = etot_(:,ke) + pres_hyd(:,ke) + dpres_(:,ke)
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 !-- ETOT
273 call sparsemat_matmul(dx, lmesh%Gsqrt(:,ke) * u_(:) * enthalpy_(:), fx)
274 call sparsemat_matmul(dy, lmesh%Gsqrt(:,ke) * v_(:) * enthalpy_(:), fy)
275 call sparsemat_matmul(dz, lmesh%Gsqrt(:,ke) * wt_(:) * enthalpy_(:), fz)
276 call sparsemat_matmul(lift, lmesh%Fscale(:,ke) * del_flux(:,ke,etot_vid), liftdelflx)
277
278 etot_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 / 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.