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

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

Functions/Subroutines

subroutine, public atm_dyn_dgm_nonhydro3d_rhot_hevi_init (mesh)
 
subroutine, public atm_dyn_dgm_nonhydro3d_rhot_hevi_final ()
 
subroutine, public atm_dyn_dgm_nonhydro3d_rhot_hevi_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)
 
subroutine, public atm_dyn_dgm_nonhydro3d_rhot_hevi_cal_vi (dens_dt, momx_dt, momy_dt, momz_dt, rhot_dt, ddens_, momx_, momy_, momz_, drhot_, dens_hyd, pres_hyd, ddens0_, momx0_, momy0_, momz0_, drhot0_, 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 compressible nonhydrostatic equations, which consist of mass, momentum, and thermodynamics (density * potential temperature conservation) equations.
Author
Yuta Kawai, Team SCALE

Function/Subroutine Documentation

◆ atm_dyn_dgm_nonhydro3d_rhot_hevi_init()

subroutine, public scale_atm_dyn_dgm_nonhydro3d_rhot_hevi::atm_dyn_dgm_nonhydro3d_rhot_hevi_init ( class(meshbase3d), intent(in) mesh)

Definition at line 73 of file scale_atm_dyn_dgm_nonhydro3d_rhot_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_rhot_hevi_final()

subroutine, public scale_atm_dyn_dgm_nonhydro3d_rhot_hevi::atm_dyn_dgm_nonhydro3d_rhot_hevi_final

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

subroutine, public scale_atm_dyn_dgm_nonhydro3d_rhot_hevi::atm_dyn_dgm_nonhydro3d_rhot_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) 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 95 of file scale_atm_dyn_dgm_nonhydro3d_rhot_hevi.F90.

101
104
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) :: RHOT_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) :: DRHOT_(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) :: RHOT_(elem%Np)
139 real(RP) :: rdens_(elem%Np), u_(elem%Np), v_(elem%Np), w_(elem%Np), wt_(elem%Np)
140 real(RP) :: GsqrtV(elem%Np), RGsqrtV(elem%Np)
141 real(RP) :: Cori(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( ke2d, Cori, &
170 !$omp RHOT_, rdens_, u_, v_, w_, wt_, &
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 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 !-- Gradient hydrostatic pressure
194
195 dpres_hyd(:) = pres_hyd(:,ke) - pres_hyd_ref(:,ke)
196
197 call sparsemat_matmul(dx, gsqrtv(:) * dpres_hyd(:), fx)
198 call sparsemat_matmul(dz, gsqrtv(:) * lmesh%GI3(:,ke,1) * dpres_hyd(:), fz)
199 call sparsemat_matmul(lift, lmesh%Fscale(:,ke) * del_flux_hyd(:,ke,1), liftdelflx)
200 gradphyd_x(:) = lmesh%Escale(:,ke,1,1) * fx(:) &
201 + lmesh%Escale(:,ke,3,3) * fz(:) &
202 + liftdelflx(:)
203
204 call sparsemat_matmul(dy, gsqrtv(:) * dpres_hyd(:), fy)
205 call sparsemat_matmul(dz, gsqrtv(:) * lmesh%GI3(:,ke,2) * dpres_hyd(:), fz)
206 call sparsemat_matmul(lift, lmesh%Fscale(:,ke) * del_flux_hyd(:,ke,2), liftdelflx)
207 gradphyd_y(:) = lmesh%Escale(:,ke,2,2) * fy(:) &
208 + lmesh%Escale(:,ke,3,3) * fz(:) &
209 + liftdelflx(:)
210 !-- DENS
211 call sparsemat_matmul(dx, lmesh%Gsqrt(:,ke) * momx_(:,ke), fx)
212 call sparsemat_matmul(dy, lmesh%Gsqrt(:,ke) * momy_(:,ke), fy)
213 call sparsemat_matmul(lift, lmesh%Fscale(:,ke) * del_flux(:,ke,dens_vid), liftdelflx)
214
215 dens_dt(:,ke) = - ( &
216 lmesh%Escale(:,ke,1,1) * fx(:) &
217 + lmesh%Escale(:,ke,2,2) * fy(:) &
218 + liftdelflx(:) ) / lmesh%Gsqrt(:,ke)
219
220 !-- MOMX
221 call sparsemat_matmul(dx, lmesh%Gsqrt(:,ke) * ( u_(:) * momx_(:,ke) + dpres_(:,ke) ), fx)
222 call sparsemat_matmul(dy, lmesh%Gsqrt(:,ke) * v_(:) * momx_(:,ke) , fy)
223 call sparsemat_matmul(dz, lmesh%Gsqrt(:,ke) * ( wt_(:) * momx_(:,ke) &
224 + lmesh%GI3(:,ke,1) * dpres_(:,ke) ), fz)
225 call sparsemat_matmul(lift, lmesh%Fscale(:,ke) * del_flux(:,ke,momx_vid), liftdelflx)
226
227 momx_dt(:,ke) = &
228 - ( lmesh%Escale(:,ke,1,1) * fx(:) &
229 + lmesh%Escale(:,ke,2,2) * fy(:) &
230 + lmesh%Escale(:,ke,3,3) * fz(:) &
231 + liftdelflx(:) ) / lmesh%Gsqrt(:,ke) &
232 - gradphyd_x(:) * rgsqrtv(:) &
233 + cori(:) * momy_(:,ke)
234
235 !-- MOMY
236 call sparsemat_matmul(dx, lmesh%Gsqrt(:,ke) * u_(:) * momy_(:,ke) , fx)
237 call sparsemat_matmul(dy, lmesh%Gsqrt(:,ke) * ( v_(:) * momy_(:,ke) + dpres_(:,ke) ), fy)
238 call sparsemat_matmul(dz, lmesh%Gsqrt(:,ke) * ( wt_(:) * momy_(:,ke) &
239 + lmesh%GI3(:,ke,2) * dpres_(:,ke) ), fz)
240 call sparsemat_matmul(lift, lmesh%Fscale(:,ke) * del_flux(:,ke,momy_vid), liftdelflx)
241
242 momy_dt(:,ke) = &
243 - ( lmesh%Escale(:,ke,1,1) * fx(:) &
244 + lmesh%Escale(:,ke,2,2) * fy(:) &
245 + lmesh%Escale(:,ke,3,3) * fz(:) &
246 + liftdelflx(:) ) / lmesh%Gsqrt(:,ke) &
247 - gradphyd_y(:) * rgsqrtv(:) &
248 - cori(:) * momx_(:,ke)
249
250 !-- MOMZ
251 call sparsemat_matmul(dx, lmesh%Gsqrt(:,ke) * u_(:) * momz_(:,ke), fx)
252 call sparsemat_matmul(dy, lmesh%Gsqrt(:,ke) * v_(:) * momz_(:,ke), fy)
253 call sparsemat_matmul(dz, lmesh%Gsqrt(:,ke) * wt_(:) * momz_(:,ke), fz)
254 call sparsemat_matmul(lift, lmesh%Fscale(:,ke) * del_flux(:,ke,momz_vid), liftdelflx)
255
256 momz_dt(:,ke) = &
257 - ( lmesh%Escale(:,ke,1,1) * fx(:) &
258 + lmesh%Escale(:,ke,2,2) * fy(:) &
259 + lmesh%Escale(:,ke,3,3) * fz(:) &
260 + liftdelflx(:) ) / lmesh%Gsqrt(:,ke)
261
262 !-- RHOT
263 call sparsemat_matmul(dx, lmesh%Gsqrt(:,ke) * u_(:) * rhot_(:), fx)
264 call sparsemat_matmul(dy, lmesh%Gsqrt(:,ke) * v_(:) * rhot_(:), fy)
265 call sparsemat_matmul(lift, lmesh%Fscale(:,ke) * del_flux(:,ke,rhot_vid), liftdelflx)
266
267 rhot_dt(:,ke) = &
268 - ( lmesh%Escale(:,ke,1,1) * fx(:) &
269 + lmesh%Escale(:,ke,2,2) * fy(:) &
270 + liftdelflx(:) ) / lmesh%Gsqrt(:,ke)
271 end do
272 call prof_rapend( 'cal_dyn_tend_interior', 3)
273
274 return
module FElib / Fluid dyn solver / Atmosphere / Nonhydrostatic model / HEVI / Numflux
subroutine, public atm_dyn_dgm_nonhydro3d_rhot_hevi_numflux_get_generalvc(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_hevi_numflux::atm_dyn_dgm_nonhydro3d_rhot_hevi_numflux_get_generalvc().

◆ atm_dyn_dgm_nonhydro3d_rhot_hevi_cal_vi()

subroutine, public scale_atm_dyn_dgm_nonhydro3d_rhot_hevi::atm_dyn_dgm_nonhydro3d_rhot_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) 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) 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) drhot0_,
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 280 of file scale_atm_dyn_dgm_nonhydro3d_rhot_hevi.F90.

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

References scale_atm_dyn_dgm_nonhydro3d_rhot_hevi_common::atm_dyn_dgm_nonhydro3d_rhot_hevi_common_construct_matbnd(), scale_atm_dyn_dgm_nonhydro3d_rhot_hevi_common::atm_dyn_dgm_nonhydro3d_rhot_hevi_common_construct_matbnd_uv(), scale_atm_dyn_dgm_nonhydro3d_rhot_hevi_common::atm_dyn_dgm_nonhydro3d_rhot_hevi_common_eval_ax(), scale_atm_dyn_dgm_nonhydro3d_rhot_hevi_common::atm_dyn_dgm_nonhydro3d_rhot_hevi_common_eval_ax_uv(), and scale_atm_dyn_dgm_nonhydro3d_common::intrpmat_vpordm1.