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

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

Functions/Subroutines

subroutine, public atm_dyn_dgm_nonhydro3d_etot_hevi_init (mesh)
 
subroutine, public atm_dyn_dgm_nonhydro3d_etot_hevi_final ()
 
subroutine, public atm_dyn_dgm_nonhydro3d_etot_hevi_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)
 
subroutine, public atm_dyn_dgm_nonhydro3d_etot_hevi_cal_vi (dens_dt, momx_dt, momy_dt, momz_dt, etot_dt, ddens_, momx_, momy_, momz_, etot_, dens_hyd, pres_hyd, ddens0_, momx0_, momy0_, momz0_, etot0_, rtot, cvtot, cptot, element3d_operation, dz, lift, impl_fac, dt, lmesh, elem, lmesh2d, elem2d)
 

Detailed Description

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

Description
HEVI DGM scheme for Atmospheric dynamical process The governing equations is a fully compressibile nonhydrostic 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_hevi_init()

subroutine, public scale_atm_dyn_dgm_nonhydro3d_etot_hevi::atm_dyn_dgm_nonhydro3d_etot_hevi_init ( class(meshbase3d), intent(in) mesh)

Definition at line 73 of file scale_atm_dyn_dgm_nonhydro3d_etot_hevi.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_nonhydro3d_etot_hevi_final()

subroutine, public scale_atm_dyn_dgm_nonhydro3d_etot_hevi::atm_dyn_dgm_nonhydro3d_etot_hevi_final

Definition at line 84 of file scale_atm_dyn_dgm_nonhydro3d_etot_hevi.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_etot_hevi_cal_tend()

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

101
104 use, intrinsic :: ieee_arithmetic, only: isnan => ieee_is_nan
105 implicit none
106
107 class(LocalMesh3D), intent(in) :: lmesh
108 class(ElementBase3D), intent(in) :: elem
109 class(LocalMesh2D), intent(in) :: lmesh2D
110 class(ElementBase2D), intent(in) :: elem2D
111 class(ElementOperationBase3D), intent(in) :: element3D_operation
112 type(SparseMat), intent(in) :: Dx, Dy, Dz, Sx, Sy, Sz, Lift
113 real(RP), intent(out) :: DENS_dt(elem%Np,lmesh%NeA)
114 real(RP), intent(out) :: MOMX_dt(elem%Np,lmesh%NeA)
115 real(RP), intent(out) :: MOMY_dt(elem%Np,lmesh%NeA)
116 real(RP), intent(out) :: MOMZ_dt(elem%Np,lmesh%NeA)
117 real(RP), intent(out) :: ETOT_dt(elem%Np,lmesh%NeA)
118 real(RP), intent(in) :: DDENS_(elem%Np,lmesh%NeA)
119 real(RP), intent(in) :: MOMX_(elem%Np,lmesh%NeA)
120 real(RP), intent(in) :: MOMY_(elem%Np,lmesh%NeA)
121 real(RP), intent(in) :: MOMZ_(elem%Np,lmesh%NeA)
122 real(RP), intent(in) :: ETOT_(elem%Np,lmesh%NeA)
123 real(RP), intent(in) :: DPRES_(elem%Np,lmesh%NeA)
124 real(RP), intent(in) :: DENS_hyd(elem%Np,lmesh%NeA)
125 real(RP), intent(in) :: PRES_hyd(elem%Np,lmesh%NeA)
126 real(RP), intent(in) :: PRES_hyd_ref(elem%Np,lmesh%NeA)
127 real(RP), intent(in) :: CORIOLIS(elem2D%Np,lmesh2D%NeA)
128 real(RP), intent(in) :: Rtot(elem%Np,lmesh%NeA)
129 real(RP), intent(in) :: CVtot(elem%Np,lmesh%NeA)
130 real(RP), intent(in) :: CPtot(elem%Np,lmesh%NeA)
131 real(RP), intent(in) :: DPhydDx(elem%Np,lmesh%NeA)
132 real(RP), intent(in) :: DPhydDy(elem%Np,lmesh%NeA)
133
134 real(RP) :: Fx(elem%Np), Fy(elem%Np), Fz(elem%Np), LiftDelFlx(elem%Np)
135 real(RP) :: DPRES_hyd(elem%Np), GradPhyd_x(elem%Np), GradPhyd_y(elem%Np)
136 real(RP) :: del_flux(elem%NfpTot,lmesh%Ne,PRGVAR_NUM)
137 real(RP) :: del_flux_hyd(elem%NfpTot,lmesh%Ne,2)
138 real(RP) :: rdens_(elem%Np), u_(elem%Np), v_(elem%Np), w_(elem%Np), wt_(elem%Np)
139 real(RP) :: GsqrtV(elem%Np), RGsqrtV(elem%Np)
140 real(RP) :: Cori(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 integer :: p
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_, etot_, dpres_, dens_hyd, pres_hyd, & ! (in)
157 rtot, cvtot, cptot, & ! (in)
158 lmesh%Gsqrt, lmesh%GI3(:,:,1), lmesh%GI3(:,:,2), lmesh%zlev(:,:), & ! (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( ke2d, Cori, &
173 !$omp rdens_, u_, v_, w_, wt_, &
174 !$omp enthalpy_, &
175 !$omp DPRES_hyd, GradPhyd_x, GradPhyd_y, &
176 !$omp GsqrtV, RGsqrtV, &
177 !$omp Fx, Fy, Fz, LiftDelFlx )
178 do ke = lmesh%NeS, lmesh%NeE
179 !--
180 ke2d = lmesh%EMap3Dto2D(ke)
181 cori(:) = coriolis(elem%IndexH2Dto3D(:),ke2d)
182
183 gsqrtv(:) = lmesh%Gsqrt(:,ke) / lmesh%GsqrtH(elem%IndexH2Dto3D,ke2d)
184 rgsqrtv(:) = 1.0_rp / gsqrtv(:)
185
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 ! DPRES_(:) = ( CPtot(:,ke) / CVtot(:,ke) - 1.0_RP ) &
194 ! * ( ETOT_(:,ke) - Grav * ( DDENS_(:,ke) + DENS_hyd(:,ke) ) * lmesh%zlev(:,ke) &
195 ! - 0.5_RP * ( MOMX_(:,ke) * u_(:) + MOMY_(:,ke) * v_(:) + MOMZ_(:,ke) * w_(:) ) ) &
196 ! - PRES_hyd(:,ke)
197
198 enthalpy_(:) = etot_(:,ke) + pres_hyd(:,ke) + dpres_(:,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(lift, lmesh%Fscale(:,ke) * del_flux(:,ke,dens_vid), liftdelflx)
222
223 dens_dt(:,ke) = - ( &
224 lmesh%Escale(:,ke,1,1) * fx(:) &
225 + lmesh%Escale(:,ke,2,2) * fy(:) &
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), fz)
262 call sparsemat_matmul(lift, lmesh%Fscale(:,ke) * del_flux(:,ke,momz_vid), liftdelflx)
263
264 momz_dt(:,ke) = &
265 - ( lmesh%Escale(:,ke,1,1) * fx(:) &
266 + lmesh%Escale(:,ke,2,2) * fy(:) &
267 + lmesh%Escale(:,ke,3,3) * fz(:) &
268 + liftdelflx(:) ) / lmesh%Gsqrt(:,ke)
269
270 !-- ETOT
271 call sparsemat_matmul(dx, lmesh%Gsqrt(:,ke) * u_(:) * enthalpy_(:), fx)
272 call sparsemat_matmul(dy, lmesh%Gsqrt(:,ke) * v_(:) * enthalpy_(:), fy)
273 call sparsemat_matmul(lift, lmesh%Fscale(:,ke) * del_flux(:,ke,etot_vid), liftdelflx)
274
275 etot_dt(:,ke) = &
276 - ( lmesh%Escale(:,ke,1,1) * fx(:) &
277 + lmesh%Escale(:,ke,2,2) * fy(:) &
278 + liftdelflx(:) ) / lmesh%Gsqrt(:,ke)
279
280 end do
281 call prof_rapend( 'cal_dyn_tend_interior', 3)
282
283 return
module FElib / Fluid dyn solver / Atmosphere / Nonhydrostatic model / HEVI / Numflux
subroutine, public atm_dyn_dgm_nonhydro3d_etot_hevi_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_hevi_numflux::atm_dyn_dgm_nonhydro3d_etot_hevi_numflux_get_generalvc().

◆ atm_dyn_dgm_nonhydro3d_etot_hevi_cal_vi()

subroutine, public scale_atm_dyn_dgm_nonhydro3d_etot_hevi::atm_dyn_dgm_nonhydro3d_etot_hevi_cal_vi ( 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) dens_hyd,
real(rp), dimension(elem%np,lmesh%nea), intent(in) pres_hyd,
real(rp), dimension(elem%np,lmesh%nea), intent(in) ddens0_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) momx0_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) momy0_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) momz0_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) etot0_,
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,
class(elementoperationbase3d), intent(in) element3d_operation,
class(sparsemat), intent(in) dz,
class(sparsemat), intent(in) lift,
real(rp), intent(in) impl_fac,
real(rp), intent(in) dt,
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_etot_hevi.F90.

297
303
304 implicit none
305
306 class(LocalMesh3D), intent(in) :: lmesh
307 class(ElementBase3D), intent(in) :: elem
308 class(LocalMesh2D), intent(in) :: lmesh2D
309 class(ElementBase2D), intent(in) :: elem2D
310 real(RP), intent(out) :: DENS_dt(elem%Np,lmesh%NeA)
311 real(RP), intent(out) :: MOMX_dt(elem%Np,lmesh%NeA)
312 real(RP), intent(out) :: MOMY_dt(elem%Np,lmesh%NeA)
313 real(RP), intent(out) :: MOMZ_dt(elem%Np,lmesh%NeA)
314 real(RP), intent(out) :: ETOT_dt(elem%Np,lmesh%NeA)
315 real(RP), intent(in) :: DDENS_(elem%Np,lmesh%NeA)
316 real(RP), intent(in) :: MOMX_(elem%Np,lmesh%NeA)
317 real(RP), intent(in) :: MOMY_(elem%Np,lmesh%NeA)
318 real(RP), intent(in) :: MOMZ_(elem%Np,lmesh%NeA)
319 real(RP), intent(in) :: ETOT_(elem%Np,lmesh%NeA)
320 real(RP), intent(in) :: DENS_hyd(elem%Np,lmesh%NeA)
321 real(RP), intent(in) :: PRES_hyd(elem%Np,lmesh%NeA)
322 real(RP), intent(in) :: DDENS0_(elem%Np,lmesh%NeA)
323 real(RP), intent(in) :: MOMX0_(elem%Np,lmesh%NeA)
324 real(RP), intent(in) :: MOMY0_(elem%Np,lmesh%NeA)
325 real(RP), intent(in) :: MOMZ0_(elem%Np,lmesh%NeA)
326 real(RP), intent(in) :: ETOT0_(elem%Np,lmesh%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 class(ElementOperationBase3D), intent(in) :: element3D_operation
331 class(SparseMat), intent(in) :: Dz, Lift
332 real(RP), intent(in) :: impl_fac
333 real(RP), intent(in) :: dt
334
335 real(RP) :: PROG_VARS (elem%Np,lmesh%NeZ,PRGVAR_NUM,lmesh%NeX*lmesh%NeY)
336 real(RP) :: DPRES (elem%Np,lmesh%NeZ,lmesh%NeX*lmesh%NeY)
337 real(RP) :: PROG_VARS0(elem%Np,lmesh%NeZ,PRGVAR_NUM,lmesh%NeX*lmesh%NeY)
338 real(RP) :: DPRES0 (elem%Np,lmesh%NeZ,lmesh%NeX*lmesh%NeY)
339 real(RP) :: b1D(3,elem%Nnode_v,lmesh%NeZ,elem%Nnode_h1D**2,lmesh%NeX*lmesh%NeY)
340 real(RP) :: GeoPot (elem%Np,lmesh%NeZ,lmesh%NeX*lmesh%NeY)
341 real(RP) :: KinHovDENS(elem%Np,lmesh%NeZ,lmesh%NeX*lmesh%NeY)
342 integer :: ipiv(elem%Nnode_v*3*lmesh%NeZ,elem%Nnode_h1D**2)
343 real(RP) :: b1D_uv(elem%Nnode_v,lmesh%NeZ,2,elem%Nnode_h1D**2,lmesh%NeX*lmesh%NeY)
344 integer :: ipiv_uv(elem%Nnode_v*1*lmesh%NeZ,elem%Nnode_h1D**2)
345 real(RP) :: alph(elem%NfpTot,lmesh%NeZ,lmesh%NeX*lmesh%NeY)
346 real(RP) :: Rtot_z(elem%Np,lmesh%NeZ,lmesh%NeX*lmesh%NeY)
347 real(RP) :: CPtot_ov_CVtot(elem%Np,lmesh%NeZ,lmesh%NeX*lmesh%NeY)
348 real(RP) :: DENS_hyd_z(elem%Np,lmesh%NeZ,lmesh%NeX*lmesh%NeY)
349 real(RP) :: PRES_hyd_z(elem%Np,lmesh%NeZ,lmesh%NeX*lmesh%NeY)
350 real(RP) :: GnnM_z(elem%Np,lmesh%NeZ,lmesh%NeX*lmesh%NeY)
351 real(RP) :: G13_z(elem%Np,lmesh%NeZ,lmesh%NeX*lmesh%NeY)
352 real(RP) :: G23_z(elem%Np,lmesh%NeZ,lmesh%NeX*lmesh%NeY)
353 real(RP) :: GsqrtV_z(elem%Np,lmesh%NeZ,lmesh%NeX*lmesh%NeY)
354 real(RP) :: nz(elem%NfpTot,lmesh%NeZ,lmesh%NeX*lmesh%NeY)
355 integer :: vmapM(elem%NfpTot,lmesh%NeZ)
356 integer :: vmapP(elem%NfpTot,lmesh%NeZ)
357 integer :: ColMask(elem%Nnode_v)
358 integer :: ke_xy, ke_z, ke, ke2D, v
359 integer :: itr_nlin
360 integer :: kl, ku, nz_1D
361 integer :: kl_uv, ku_uv, nz_1D_uv
362 integer :: ij, info
363 logical :: is_converged
364
365 real(RP), allocatable :: PmatBnd(:,:,:)
366 real(RP), allocatable :: PmatBnd_uv(:,:,:)
367
368 real(RP) :: DENS_(elem%Np)
369 !------------------------------------------------------------------------
370
371 call prof_rapstart( 'hevi_cal_vi_prep', 3)
372
373 nz_1d = elem%Nnode_v * 3 * lmesh%NeZ
374 kl = ( elem%Nnode_v + 1 ) * 3 - 1
375 ku = kl
376 nz_1d_uv = elem%Nnode_v * 1 * lmesh%NeZ
377 kl_uv = elem%Nnode_v
378 ku_uv = kl_uv
379 allocate( pmatbnd(2*kl+ku+1,nz_1d,elem%Nnode_h1D**2) )
380 allocate( pmatbnd_uv(2*kl_uv+ku_uv+1,nz_1d_uv,elem%Nnode_h1D**2) )
381
382 call lmesh%GetVmapZ1D( vmapm, vmapp ) ! (out)
383
384 !-
385
386 !$omp parallel private( ke_xy, ke_z, ke, ke2D, DENS_ )
387 !$omp do collapse(2)
388 do ke_xy=1, lmesh%NeX*lmesh%NeY
389 do ke_z=1, lmesh%NeZ
390 ke = ke_xy + (ke_z-1)*lmesh%NeX*lmesh%NeY
391 ke2d = lmesh%EMap3Dto2D(ke)
392
393 prog_vars(:,ke_z,dens_vid,ke_xy) = ddens0_(:,ke)
394 prog_vars(:,ke_z,momx_vid,ke_xy) = momx0_(:,ke)
395 prog_vars(:,ke_z,momy_vid,ke_xy) = momy0_(:,ke)
396 prog_vars(:,ke_z,momz_vid,ke_xy) = momz0_(:,ke)
397 prog_vars(:,ke_z,etot_vid,ke_xy) = etot0_(:,ke)
398
399 dens_hyd_z(:,ke_z,ke_xy) = dens_hyd(:,ke)
400 pres_hyd_z(:,ke_z,ke_xy) = pres_hyd(:,ke)
401 geopot(:,ke_z,ke_xy) = grav * lmesh%zlev(:,ke)
402
403 dens_(:) = dens_hyd(:,ke) + ddens0_(:,ke)
404 kinhovdens(:,ke_z,ke_xy) = 0.5_rp * ( &
405 momx0_(:,ke) * ( lmesh%G_ij(elem%IndexH2Dto3D,ke2d,1,1) * momx0_(:,ke) + lmesh%G_ij(elem%IndexH2Dto3D,ke2d,2,1) * momy0_(:,ke) ) &
406 + momy0_(:,ke) * ( lmesh%G_ij(elem%IndexH2Dto3D,ke2d,2,1) * momx0_(:,ke) + lmesh%G_ij(elem%IndexH2Dto3D,ke2d,2,2) * momy0_(:,ke) ) &
407 ) / dens_(:)**2
408
409
410 rtot_z(:,ke_z,ke_xy) = rtot(:,ke)
411 cptot_ov_cvtot(:,ke_z,ke_xy) = cptot(:,ke) / cvtot(:,ke)
412
413 dpres(:,ke_z,ke_xy) = &
414 ( cptot_ov_cvtot(:,ke_z,ke_xy) - 1.0_rp ) &
415 * ( etot0_(:,ke) - ( dens_(:) * ( kinhovdens(:,ke_z,ke_xy) + geopot(:,ke_z,ke_xy) ) + 0.5_rp * momz0_(:,ke)**2 / dens_(:) ) ) &
416 - pres_hyd(:,ke)
417
418 nz(:,ke_z,ke_xy) = lmesh%normal_fn(:,ke,3)
419 g13_z(:,ke_z,ke_xy) = lmesh%GI3(:,ke,1)
420 g23_z(:,ke_z,ke_xy) = lmesh%GI3(:,ke,2)
421 gsqrtv_z(:,ke_z,ke_xy) = lmesh%Gsqrt(:,ke) / lmesh%GsqrtH(elem%IndexH2Dto3D,ke2d)
422
423 gnnm_z(:,ke_z,ke_xy) = ( 1.0_rp / gsqrtv_z(:,ke_z,ke_xy)**2 &
424 + g13_z(:,ke_z,ke_xy)**2 + g23_z(:,ke_z,ke_xy) )
425 end do
426 end do
427 !$omp end do
428 !$omp workshare
429 prog_vars0(:,:,:,:) = prog_vars(:,:,:,:)
430 dpres0(:,:,:) = dpres(:,:,:)
431 !$omp end workshare
432 !$omp end parallel
433
434 call prof_rapend( 'hevi_cal_vi_prep', 3)
435
436 !--
437
438 if ( abs(impl_fac) > 0.0_rp ) then
439 call prof_rapstart( 'hevi_cal_vi_itr', 3)
440
441 ! G = (q^n+1 - q^n*) + impl_fac * A(q^n+1) = 0
442 ! dG/dq^n+1 del[q] = - G(q^n*)
443 do itr_nlin = 1, 1
444 call prof_rapstart( 'hevi_cal_vi_ax', 3)
445
446 call vi_eval_ax_uv( &
447 momx_dt(:,:), momy_dt(:,:), alph(:,:,:), & ! (out, dummy)
448 prog_vars, dpres, prog_vars0, dpres0, & ! (in)
449 ddens_, momx_, momy_, momz_, etot_, & ! (in)
450 dens_hyd_z, pres_hyd_z, & ! (in)
451 rtot_z, cptot_ov_cvtot, & ! (in)
452 dz, lift, intrpmat_vpordm1, & ! (in)
453 gnnm_z, g13_z, g23_z, gsqrtv_z, & ! (in)
454 impl_fac, dt, & ! (in)
455 lmesh, elem, nz, vmapm, vmapp, & ! (in)
456 b1d_uv(:,:,:,:,:) ) ! (out)
457
458 call prof_rapend( 'hevi_cal_vi_ax', 3)
459
460 do ke_xy=1, lmesh%NeX * lmesh%NeY
461 call prof_rapstart( 'hevi_cal_vi_matbnd', 3)
462 call vi_construct_matbnd_uv( pmatbnd_uv(:,:,:), & ! (out)
463 kl_uv, ku_uv, nz_1d_uv, & ! (in)
464 prog_vars(:,:,:,ke_xy), kinhovdens(:,:,ke_xy), & ! (in)
465 dens_hyd_z(:,:,ke_xy), pres_hyd_z(:,:,ke_xy), & ! (in)
466 g13_z(:,:,ke_xy), g23_z(:,:,ke_xy), gsqrtv_z(:,:,ke_xy), & ! (in)
467 alph(:,:,ke_xy), & ! (in)
468 rtot_z(:,:,ke_xy), cptot_ov_cvtot(:,:,ke_xy), & ! (in)
469 geopot(:,:,ke_xy), & ! (in)
470 dz, lift, intrpmat_vpordm1, & ! (in)
471 impl_fac, dt, & ! (in)
472 lmesh, elem, nz(:,:,ke_xy), vmapm, vmapp, ke_xy, 1 ) ! (in)
473 call prof_rapend( 'hevi_cal_vi_matbnd', 3)
474
475 call prof_rapstart( 'hevi_cal_vi_lin', 3)
476 !$omp parallel private(ij, v, ke_z, info, ColMask)
477 !$omp do
478 do ij=1, elem%Nnode_h1D**2
479 call dgbsv( nz_1d_uv, kl_uv, ku_uv, 2, pmatbnd_uv(:,:,ij), 2*kl_uv+ku_uv+1, ipiv_uv(:,ij), b1d_uv(:,:,:,ij,ke_xy), nz_1d_uv, info)
480
481 colmask(:) = elem%Colmask(:,ij)
482 do ke_z=1, lmesh%NeZ
483 prog_vars(colmask(:),ke_z,momx_vid,ke_xy) = prog_vars(colmask(:),ke_z,momx_vid,ke_xy) + b1d_uv(:,ke_z,1,ij,ke_xy)
484 prog_vars(colmask(:),ke_z,momy_vid,ke_xy) = prog_vars(colmask(:),ke_z,momy_vid,ke_xy) + b1d_uv(:,ke_z,2,ij,ke_xy)
485 end do
486 end do ! for ij
487 !$omp end do
488 !$omp end parallel
489 call prof_rapend( 'hevi_cal_vi_lin', 3)
490
491 end do ! for ke_xy
492
493 call prof_rapstart( 'hevi_cal_vi_ax', 3)
494 call vi_eval_ax( &
495 dens_dt(:,:), momz_dt(:,:), etot_dt(:,:), & ! (out, dummy)
496 alph(:,:,:), & ! (in)
497 prog_vars, dpres, prog_vars0, dpres0, & ! (in)
498 ddens_, momx_, momy_, momz_, etot_, & ! (in)
499 dens_hyd_z, pres_hyd_z, & ! (in)
500 rtot_z, cptot_ov_cvtot, & ! (in)
501 dz, lift, intrpmat_vpordm1, & ! (in)
502 gnnm_z, g13_z, g23_z, gsqrtv_z, & ! (in)
503 impl_fac, dt, & ! (in)
504 lmesh, elem, nz, vmapm, vmapp, & ! (in)
505 b1d(:,:,:,:,:) ) ! (out)
506 call prof_rapend( 'hevi_cal_vi_ax', 3)
507
508 do ke_xy=1, lmesh%NeX * lmesh%NeY
509 call prof_rapstart( 'hevi_cal_vi_matbnd', 3)
510 call vi_construct_matbnd( pmatbnd(:,:,:), & ! (out)
511 kl, ku, nz_1d, & ! (in)
512 prog_vars(:,:,:,ke_xy), kinhovdens(:,:,ke_xy), & ! (in)
513 dens_hyd_z(:,:,ke_xy), pres_hyd_z(:,:,ke_xy), & ! (in)
514 g13_z(:,:,ke_xy), g23_z(:,:,ke_xy), gsqrtv_z(:,:,ke_xy), & ! (in)
515 alph(:,:,ke_xy), & ! (in)
516 rtot_z(:,:,ke_xy), cptot_ov_cvtot(:,:,ke_xy), & ! (in)
517 geopot(:,:,ke_xy), & ! (in)
518 dz, lift, intrpmat_vpordm1, & ! (in)
519 impl_fac, dt, & ! (in)
520 lmesh, elem, nz(:,:,ke_xy), vmapm, vmapp, ke_xy, 1 ) ! (in)
521 call prof_rapend( 'hevi_cal_vi_matbnd', 3)
522
523 call prof_rapstart( 'hevi_cal_vi_lin', 3)
524 !$omp parallel private(ij, v, ke_z, info, ColMask, DENS_)
525 !$omp do
526 do ij=1, elem%Nnode_h1D**2
527 call dgbsv( nz_1d, kl, ku, 1, pmatbnd(:,:,ij), 2*kl+ku+1, ipiv(:,ij), b1d(:,:,:,ij,ke_xy), nz_1d, info)
528
529 colmask(:) = elem%Colmask(:,ij)
530 do ke_z=1, lmesh%NeZ
531 prog_vars(colmask(:),ke_z,dens_vid,ke_xy) = prog_vars(colmask(:),ke_z,dens_vid,ke_xy) + b1d(1,:,ke_z,ij,ke_xy)
532 prog_vars(colmask(:),ke_z,momz_vid,ke_xy) = prog_vars(colmask(:),ke_z,momz_vid,ke_xy) + b1d(2,:,ke_z,ij,ke_xy)
533 prog_vars(colmask(:),ke_z,etot_vid,ke_xy) = prog_vars(colmask(:),ke_z,etot_vid,ke_xy) + b1d(3,:,ke_z,ij,ke_xy)
534 end do
535 end do ! for ij
536 !$omp end do
537 !$omp do
538 do ke_z=1, lmesh%NeZ
539 dens_(:) = dens_hyd_z(:,ke_z,ke_xy) + prog_vars(:,ke_z,dens_vid,ke_xy)
540 dpres(:,ke_z,ke_xy) = &
541 ( cptot_ov_cvtot(:,ke_z,ke_xy) - 1.0_rp ) &
542 * ( prog_vars(:,ke_z,etot_vid,ke_xy) &
543 - ( dens_(:) * ( kinhovdens(:,ke_z,ke_xy) + geopot(:,ke_z,ke_xy) ) + 0.5_rp * prog_vars(:,ke_z,momz_vid,ke_xy)**2 / dens_(:) ) ) &
544 - pres_hyd_z(:,ke_z,ke_xy)
545 end do
546 !$omp end parallel
547 call prof_rapend( 'hevi_cal_vi_lin', 3)
548
549 end do ! for ke_xy
550 end do ! itr nlin
551
552 call prof_rapend( 'hevi_cal_vi_itr', 3)
553 end if
554
555 call prof_rapstart( 'hevi_cal_vi_retrun_var', 3)
556 if ( abs(impl_fac) > 0.0_rp) then
557 !$omp parallel do collapse(2) private(ke_xy, ke_z, ke)
558 do ke_xy=1, lmesh%NeX * lmesh%NeY
559 do ke_z=1, lmesh%NeZ
560 ke = ke_xy + (ke_z-1)*lmesh%NeX*lmesh%NeY
561 dens_dt(:,ke) = ( prog_vars(:,ke_z,dens_vid,ke_xy) - ddens_(:,ke) ) / impl_fac
562 momx_dt(:,ke) = ( prog_vars(:,ke_z,momx_vid,ke_xy) - momx_(:,ke) ) / impl_fac
563 momy_dt(:,ke) = ( prog_vars(:,ke_z,momy_vid,ke_xy) - momy_(:,ke) ) / impl_fac
564 momz_dt(:,ke) = ( prog_vars(:,ke_z,momz_vid,ke_xy) - momz_(:,ke) ) / impl_fac
565 etot_dt(:,ke) = ( prog_vars(:,ke_z,etot_vid,ke_xy) - etot_(:,ke) ) / impl_fac
566 end do
567 end do
568 else
569 call vi_eval_ax_uv( &
570 momx_dt(:,:), momy_dt(:,:), alph(:,:,:), & ! (out)
571 prog_vars, dpres, prog_vars0, dpres0, & ! (in)
572 ddens_, momx_, momy_, momz_, etot_, & ! (in)
573 dens_hyd_z, pres_hyd_z, & ! (in)
574 rtot_z, cptot_ov_cvtot, & ! (in)
575 dz, lift, intrpmat_vpordm1, & ! (in)
576 gnnm_z, g13_z, g23_z, gsqrtv_z, & ! (in)
577 impl_fac, dt, & ! (in)
578 lmesh, elem, nz, vmapm, vmapp ) ! (in)
579
580 call vi_eval_ax( &
581 dens_dt(:,:), momz_dt(:,:), etot_dt(:,:), & ! (out)
582 alph(:,:,:), & ! (in)
583 prog_vars, dpres, prog_vars0, dpres0, & ! (in)
584 ddens_, momx_, momy_, momz_, etot_, & ! (in)
585 dens_hyd_z, pres_hyd_z, & ! (in)
586 rtot_z, cptot_ov_cvtot, & ! (in)
587 dz, lift, intrpmat_vpordm1, & ! (in)
588 gnnm_z, g13_z, g23_z, gsqrtv_z, & ! (in)
589 impl_fac, dt, & ! (in)
590 lmesh, elem, nz, vmapm, vmapp ) ! (in)
591 end if
592 call prof_rapend( 'hevi_cal_vi_retrun_var', 3)
593
594 return
module FElib / Fluid dyn solver / Atmosphere / Nonhydrostatic model / HEVI / Common
subroutine, public atm_dyn_dgm_nonhydro3d_etot_hevi_common_construct_matbnd(pmatbnd, kl, ku, nz_1d, prog_vars0, kinhovdens00, dens_hyd, pres_hyd, g13, g23, gsqrtv, alph, rtot, cptot_ov_cvtot, geopot, dz, lift, intrpmat_vpordm1, impl_fac, dt, lmesh, elem, nz, vmapm, vmapp, ke_x, ke_y)
subroutine, public atm_dyn_dgm_nonhydro3d_etot_hevi_common_eval_ax(dens_t, momz_t, etot_t, alph, prog_vars, dpres, prog_vars0, dpres0, ddens00, momx00, momy00, momz00, entot00, dens_hyd, pres_hyd, rtot, cptot_ov_cvtot, dz, lift, intrpmat_vpordm1, gnnm, g13, g23, gsqrtv, impl_fac, dt, lmesh, elem, nz, vmapm, vmapp, b1d_ij)
subroutine, public atm_dyn_dgm_nonhydro3d_etot_hevi_common_eval_ax_uv(momx_t, momy_t, alph, prog_vars, dpres, prog_vars0, dpres0, ddens00, momx00, momy00, momz00, entot00, dens_hyd, pres_hyd, rtot, cptot_ov_cvtot, dz, lift, intrpmat_vpordm1, gnnm, g13, g23, gsqrtv, impl_fac, dt, lmesh, elem, nz, vmapm, vmapp, b1d_ij_uv)
subroutine, public atm_dyn_dgm_nonhydro3d_etot_hevi_common_construct_matbnd_uv(pmatbnd_uv, kl_uv, ku_uv, nz_1d_uv, prog_vars0, kinhovdens00, dens_hyd, pres_hyd, g13, g23, gsqrtv, alph, rtot, cptot_ov_cvtot, geopot, dz, lift, intrpmat_vpordm1, impl_fac, dt, lmesh, elem, nz, vmapm, vmapp, ke_x, ke_y)

References scale_atm_dyn_dgm_nonhydro3d_etot_hevi_common::atm_dyn_dgm_nonhydro3d_etot_hevi_common_construct_matbnd(), scale_atm_dyn_dgm_nonhydro3d_etot_hevi_common::atm_dyn_dgm_nonhydro3d_etot_hevi_common_construct_matbnd_uv(), scale_atm_dyn_dgm_nonhydro3d_etot_hevi_common::atm_dyn_dgm_nonhydro3d_etot_hevi_common_eval_ax(), scale_atm_dyn_dgm_nonhydro3d_etot_hevi_common::atm_dyn_dgm_nonhydro3d_etot_hevi_common_eval_ax_uv(), and scale_atm_dyn_dgm_nonhydro3d_common::intrpmat_vpordm1.