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

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

Functions/Subroutines

subroutine, public atm_dyn_dgm_nonhydro2d_init (mesh)
 
subroutine, public atm_dyn_dgm_nonhydro2d_prepair_expfilter (elem, etac, alpha, ord)
 
subroutine, public atm_dyn_dgm_nonhydro2d_final ()
 
subroutine, public atm_dyn_dgm_nonhydro2d_filter_prgvar (ddens_, momx_, momz_, drhot_, lmesh, elem)
 
subroutine, public atm_dyn_dgm_nonhydro2d_cal_tend (dens_dt, momx_dt, momz_dt, rhot_dt, ddens_, momx_, momz_, drhot_, dens_hyd, pres_hyd, gxu_, gzu_, gxw_, gzw_, gxpt_, gzpt_, visccoef_h, visccoef_v, diffcoef_h, diffcoef_v, dx, dz, sx, sz, lift, lmesh, elem)
 
subroutine, public atm_dyn_dgm_nonhydro2d_cal_grad_diffvars (gxu_, gzu_, gxw_, gzw_, gxpt_, gzpt_, ddens_, momx_, momz_, drhot_, dens_hyd, pres_hyd, dx, dz, lift, lmesh, elem)
 

Detailed Description

module FElib / Fluid dyn solver / Atmosphere / 2D nonhydrostatic model / HEVE

Description
HEVE DGM scheme for 2D Atmospheric Dynamical process. The governing equations is a fully compressibile nonhydrostic equations, which consist of mass, momentum, and thermodynamics equations.
Author
Yuta Kawai, Team SCALE

Function/Subroutine Documentation

◆ atm_dyn_dgm_nonhydro2d_init()

subroutine, public scale_atm_dyn_dgm_nonhydro2d::atm_dyn_dgm_nonhydro2d_init ( class(meshbase2d), intent(in), target mesh)

Definition at line 84 of file scale_atm_dyn_dgm_nonhydro2d.F90.

85
86 implicit none
87 class(MeshBase2D), intent(in), target :: mesh
88
89
90 integer :: p1, p_
91 real(RP), allocatable :: invV_VPOrdM1(:,:)
92
93 type(ElementBase2D), pointer :: elem
94 !--------------------------------------------
95
96 elem => mesh%refElem2D
97
98 allocate( invv_vpordm1(elem%Np,elem%Np) )
99 invv_vpordm1(:,:) = elem%invV
100 do p1=1, elem%PolyOrder+1
101 p_ = p1 + elem%PolyOrder*(elem%PolyOrder + 1)
102 invv_vpordm1(p_,:) = 0.0_rp
103 end do
104
105 allocate( intrpmat_vpordm1(elem%Np,elem%Np) )
106 intrpmat_vpordm1(:,:) = matmul(elem%V, invv_vpordm1)
107
108 return

◆ atm_dyn_dgm_nonhydro2d_prepair_expfilter()

subroutine, public scale_atm_dyn_dgm_nonhydro2d::atm_dyn_dgm_nonhydro2d_prepair_expfilter ( class(elementbase2d), intent(in) elem,
real(rp), intent(in) etac,
real(rp), intent(in) alpha,
integer, intent(in) ord )

Definition at line 111 of file scale_atm_dyn_dgm_nonhydro2d.F90.

114
115 implicit none
116 class(elementbase2D), intent(in) :: elem
117 real(RP), intent(in) :: etac
118 real(RP), intent(in) :: alpha
119 integer, intent(in) :: ord
120
121 real(RP) :: filter1D(elem%Nfp)
122 real(RP) :: eta
123 integer :: p1, p2
124 integer :: l
125 !----------------------------------------------------
126
127 filter1d(:) = 1.0_rp
128 do p1=1, elem%Nfp
129 eta = dble(p1-1)/dble(elem%PolyOrder)
130 if ( eta > etac .and. p1 /= 1) then
131 filter1d(p1) = exp( - alpha*( ((eta - etac)/(1.0_rp - etac))**ord ))
132 end if
133 end do
134
135 allocate( filtermat(elem%Np,elem%Np) )
136 filtermat(:,:) = 0.0_rp
137 do p2=1, elem%Nfp
138 do p1=1, elem%Nfp
139 l = p1 + (p2-1)*elem%Nfp
140 filtermat(l,l) = filter1d(p1) * filter1d(p2)
141 end do
142 end do
143 filtermat(:,:) = matmul(filtermat, elem%invV)
144 filtermat(:,:) = matmul(elem%V, filtermat)
145
146 return

◆ atm_dyn_dgm_nonhydro2d_final()

subroutine, public scale_atm_dyn_dgm_nonhydro2d::atm_dyn_dgm_nonhydro2d_final

Definition at line 149 of file scale_atm_dyn_dgm_nonhydro2d.F90.

150 implicit none
151 !--------------------------------------------
152
153 if( allocated(intrpmat_vpordm1) ) deallocate( intrpmat_vpordm1 )
154 if( allocated(filtermat) ) deallocate( filtermat )
155
156 return

◆ atm_dyn_dgm_nonhydro2d_filter_prgvar()

subroutine, public scale_atm_dyn_dgm_nonhydro2d::atm_dyn_dgm_nonhydro2d_filter_prgvar ( real(rp), dimension(elem%np,lmesh%nea), intent(inout) ddens_,
real(rp), dimension(elem%np,lmesh%nea), intent(inout) momx_,
real(rp), dimension(elem%np,lmesh%nea), intent(inout) momz_,
real(rp), dimension(elem%np,lmesh%nea), intent(inout) drhot_,
class(localmesh2d), intent(in) lmesh,
class(elementbase2d), intent(in) elem )

Definition at line 161 of file scale_atm_dyn_dgm_nonhydro2d.F90.

163 implicit none
164
165 class(LocalMesh2D), intent(in) :: lmesh
166 class(elementbase2D), intent(in) :: elem
167 real(RP), intent(inout) :: DDENS_(elem%Np,lmesh%NeA)
168 real(RP), intent(inout) :: MOMX_(elem%Np,lmesh%NeA)
169 real(RP), intent(inout) :: MOMZ_(elem%Np,lmesh%NeA)
170 real(RP), intent(inout) :: DRHOT_(elem%Np,lmesh%NeA)
171
172 integer :: k
173 !------------------------------------
174
175 do k=1, lmesh%Ne
176 ddens_(:,k) = matmul(filtermat,ddens_(:,k))
177 momx_(:,k) = matmul(filtermat,momx_(:,k))
178 momz_(:,k) = matmul(filtermat,momz_(:,k))
179 drhot_(:,k) = matmul(filtermat,drhot_(:,k))
180 end do
181
182 return

◆ atm_dyn_dgm_nonhydro2d_cal_tend()

subroutine, public scale_atm_dyn_dgm_nonhydro2d::atm_dyn_dgm_nonhydro2d_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) 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) 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) gxu_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) gzu_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) gxw_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) gzw_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) gxpt_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) gzpt_,
real(rp), intent(in) visccoef_h,
real(rp), intent(in) visccoef_v,
real(rp), intent(in) diffcoef_h,
real(rp), intent(in) diffcoef_v,
type(sparsemat), intent(in) dx,
type(sparsemat), intent(in) dz,
type(sparsemat), intent(in) sx,
type(sparsemat), intent(in) sz,
type(sparsemat), intent(in) lift,
class(localmesh2d), intent(in) lmesh,
class(elementbase2d), intent(in) elem )

Definition at line 185 of file scale_atm_dyn_dgm_nonhydro2d.F90.

191
192 implicit none
193
194 class(LocalMesh2D), intent(in) :: lmesh
195 class(elementbase2D), intent(in) :: elem
196 type(SparseMat), intent(in) :: Dx, Dz, Sx, Sz, Lift
197 real(RP), intent(out) :: DENS_dt(elem%Np,lmesh%NeA)
198 real(RP), intent(out) :: MOMX_dt(elem%Np,lmesh%NeA)
199 real(RP), intent(out) :: MOMZ_dt(elem%Np,lmesh%NeA)
200 real(RP), intent(out) :: RHOT_dt(elem%Np,lmesh%NeA)
201 real(RP), intent(in) :: DDENS_(elem%Np,lmesh%NeA)
202 real(RP), intent(in) :: MOMX_(elem%Np,lmesh%NeA)
203 real(RP), intent(in) :: MOMZ_(elem%Np,lmesh%NeA)
204 real(RP), intent(in) :: DRHOT_(elem%Np,lmesh%NeA)
205 real(RP), intent(in) :: DENS_hyd(elem%Np,lmesh%NeA)
206 real(RP), intent(in) :: PRES_hyd(elem%Np,lmesh%NeA)
207 real(RP), intent(in) :: GxU_(elem%Np,lmesh%NeA)
208 real(RP), intent(in) :: GzU_(elem%Np,lmesh%NeA)
209 real(RP), intent(in) :: GxW_(elem%Np,lmesh%NeA)
210 real(RP), intent(in) :: GzW_(elem%Np,lmesh%NeA)
211 real(RP), intent(in) :: GxPT_(elem%Np,lmesh%NeA)
212 real(RP), intent(in) :: GzPT_(elem%Np,lmesh%NeA)
213 real(RP), intent(in) :: viscCoef_h
214 real(RP), intent(in) :: viscCoef_v
215 real(RP), intent(in) :: diffCoef_h
216 real(RP), intent(in) :: diffCoef_v
217
218 real(RP) :: Fx(elem%Np), Fz(elem%Np), LiftDelFlx(elem%Np)
219 real(RP) :: del_flux(elem%NfpTot,lmesh%Ne,PROG_VARS_NUM)
220 real(RP) :: dens_(elem%Np), RHOT_(elem%Np), dpres_(elem%Np)
221 real(RP) :: pres_(elem%Np), u_(elem%Np), w_(elem%Np)
222
223 integer :: k
224 !------------------------------------------------------------------------
225
226 call prof_rapstart( 'cal_dyn_tend_bndflux', 2)
227 call cal_del_flux_dyn( del_flux, & ! (out)
228 ddens_, momx_, momz_, drhot_, dens_hyd, pres_hyd, & ! (in)
229 gxu_, gzu_, gxw_, gzw_, gxpt_, gzpt_, & ! (in)
230 visccoef_h, visccoef_v, & ! (in)
231 diffcoef_h, diffcoef_v, & ! (in)
232 lmesh%normal_fn(:,:,1), lmesh%normal_fn(:,:,2), & ! (in)
233 lmesh%vmapM, lmesh%vmapP, & ! (in)
234 lmesh, elem ) ! (in)
235 call prof_rapend( 'cal_dyn_tend_bndflux', 2)
236
237 !-----
238 call prof_rapstart( 'cal_dyn_tend_interior', 2)
239 do k = lmesh%NeS, lmesh%NeE
240 !--
241
242 rhot_(:) = pres00/rdry * (pres_hyd(:,k)/pres00)**(cvdry/cpdry) + drhot_(:,k)
243 pres_(:) = pres00 * (rdry*rhot_(:)/pres00)**(cpdry/cvdry)
244 dpres_(:) = pres_(:) - pres_hyd(:,k)
245 dens_(:) = ddens_(:,k) + dens_hyd(:,k)
246
247 u_(:) = momx_(:,k)/dens_(:)
248 w_(:) = momz_(:,k)/dens_(:)
249
250 !-- DENS
251 call sparsemat_matmul(dx, momx_(:,k), fx)
252 call sparsemat_matmul(dz, momz_(:,k), fz)
253 call sparsemat_matmul(lift, lmesh%Fscale(:,k)*del_flux(:,k,vars_ddens_id), liftdelflx)
254
255 dens_dt(:,k) = - ( &
256 lmesh%Escale(:,k,1,1) * fx(:) &
257 + lmesh%Escale(:,k,2,2) * fz(:) &
258 + liftdelflx(:) )
259
260 !-- MOMX
261 call sparsemat_matmul(dx, u_(:)*momx_(:,k) + dpres_(:) - visccoef_h*dens_(:)*gxu_(:,k), fx)
262 call sparsemat_matmul(dz, w_(:)*momx_(:,k) - visccoef_v*dens_(:)*gzu_(:,k) ,fz)
263 call sparsemat_matmul(lift, lmesh%Fscale(:,k)*del_flux(:,k,vars_momx_id), liftdelflx)
264
265 momx_dt(:,k) = - ( &
266 lmesh%Escale(:,k,1,1) * fx(:) &
267 + lmesh%Escale(:,k,2,2) * fz(:) &
268 + liftdelflx(:) )
269
270 !-- MOMZ
271 call sparsemat_matmul(dx, u_(:)*momz_(:,k) - visccoef_h*dens_(:)*gxw_(:,k), fx)
272 call sparsemat_matmul(dz, w_(:)*momz_(:,k) + dpres_(:) - visccoef_v*dens_(:)*gzw_(:,k), fz)
273 call sparsemat_matmul(lift, lmesh%Fscale(:,k)*del_flux(:,k,vars_momz_id), liftdelflx)
274
275 momz_dt(:,k) = - ( &
276 lmesh%Escale(:,k,1,1) * fx(:) &
277 + lmesh%Escale(:,k,2,2) * fz(:) &
278 + liftdelflx(:) ) &
279 - matmul(intrpmat_vpordm1, ddens_(:,k)) * grav
280 !- DDENS_(:,k)*Grav
281
282
283 !-- RHOT
284 call sparsemat_matmul(dx, u_(:)*rhot_(:) - diffcoef_h*dens_(:)*gxpt_(:,k), fx)
285 call sparsemat_matmul(dz, w_(:)*rhot_(:) - diffcoef_v*dens_(:)*gzpt_(:,k), fz)
286 call sparsemat_matmul(lift, lmesh%Fscale(:,k)*del_flux(:,k,vars_drhot_id), liftdelflx)
287
288 rhot_dt(:,k) = - ( &
289 lmesh%Escale(:,k,1,1) * fx(:) &
290 + lmesh%Escale(:,k,2,2) * fz(:) &
291 + liftdelflx )
292
293 end do
294 call prof_rapend( 'cal_dyn_tend_interior', 2)
295
296 return

◆ atm_dyn_dgm_nonhydro2d_cal_grad_diffvars()

subroutine, public scale_atm_dyn_dgm_nonhydro2d::atm_dyn_dgm_nonhydro2d_cal_grad_diffvars ( real(rp), dimension(elem%np,lmesh%nea), intent(out) gxu_,
real(rp), dimension(elem%np,lmesh%nea), intent(out) gzu_,
real(rp), dimension(elem%np,lmesh%nea), intent(out) gxw_,
real(rp), dimension(elem%np,lmesh%nea), intent(out) gzw_,
real(rp), dimension(elem%np,lmesh%nea), intent(out) gxpt_,
real(rp), dimension(elem%np,lmesh%nea), intent(out) gzpt_,
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) 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,
type(sparsemat), intent(in) dx,
type(sparsemat), intent(in) dz,
type(sparsemat), intent(in) lift,
class(localmesh2d), intent(in) lmesh,
class(elementbase2d), intent(in) elem )

Definition at line 416 of file scale_atm_dyn_dgm_nonhydro2d.F90.

419
420 implicit none
421
422 class(LocalMesh2D), intent(in) :: lmesh
423 class(elementbase2D), intent(in) :: elem
424 type(SparseMat), intent(in) :: Dx, Dz, Lift
425 real(RP), intent(out) :: GxU_(elem%Np,lmesh%NeA)
426 real(RP), intent(out) :: GzU_(elem%Np,lmesh%NeA)
427 real(RP), intent(out) :: GxW_(elem%Np,lmesh%NeA)
428 real(RP), intent(out) :: GzW_(elem%Np,lmesh%NeA)
429 real(RP), intent(out) :: GxPT_(elem%Np,lmesh%NeA)
430 real(RP), intent(out) :: GzPT_(elem%Np,lmesh%NeA)
431 real(RP), intent(in) :: DDENS_(elem%Np,lmesh%NeA)
432 real(RP), intent(in) :: MOMX_(elem%Np,lmesh%NeA)
433 real(RP), intent(in) :: MOMZ_(elem%Np,lmesh%NeA)
434 real(RP), intent(in) :: DRHOT_(elem%Np,lmesh%NeA)
435 real(RP), intent(in) :: DENS_hyd(elem%Np,lmesh%NeA)
436 real(RP), intent(in) :: PRES_hyd(elem%Np,lmesh%NeA)
437
438 integer :: k
439 real(RP) :: DENS_(elem%Np), U_(elem%Np), W_(elem%Np), DTHETA_(elem%Np), RHOT_(elem%Np), RHOT_hyd(elem%Np)
440 real(RP) :: Fx(elem%Np), Fz(elem%Np), LiftDelFlx(elem%Np)
441 real(RP) :: del_flux(elem%NfpTot,lmesh%Ne,AUX_DIFFVARS_NUM)
442 !------------------------------------------------------------------------------
443
444 call cal_del_graddiffvar( del_flux, & ! (out)
445 ddens_, momx_, momz_, drhot_, dens_hyd, pres_hyd, & ! (in)
446 lmesh%normal_fn(:,:,1), lmesh%normal_fn(:,:,2), & ! (in)
447 lmesh%vmapM, lmesh%vmapP, & ! (in)
448 lmesh, elem ) ! (in)
449
450 do k=lmesh%NeS, lmesh%NeE
451 dens_(:) = ddens_(:,k) + dens_hyd(:,k)
452 rhot_hyd(:) = pres00/rdry * (pres_hyd(:,k)/pres00)**(cvdry/cpdry)
453 rhot_(:) = rhot_hyd(:) + drhot_(:,k)
454
455 u_(:) = momx_(:,k)/dens_(:)
456 w_(:) = momz_(:,k)/dens_(:)
457 dtheta_(:) = rhot_(:)/dens_(:) - rhot_hyd(:)/dens_hyd(:,k)
458
459 call sparsemat_matmul(dx, u_, fx)
460 call sparsemat_matmul(lift, lmesh%Fscale(:,k)*del_flux(:,k,vars_gxu_id), liftdelflx)
461 gxu_(:,k) = lmesh%Escale(:,k,1,1)*fx(:) + liftdelflx(:)
462
463 call sparsemat_matmul(dz, u_, fz)
464 call sparsemat_matmul(lift, lmesh%Fscale(:,k)*del_flux(:,k,vars_gzu_id), liftdelflx)
465 gzu_(:,k) = lmesh%Escale(:,k,2,2)*fz(:) + liftdelflx(:)
466
467 call sparsemat_matmul(dx, w_, fx)
468 call sparsemat_matmul(lift, lmesh%Fscale(:,k)*del_flux(:,k,vars_gxw_id), liftdelflx)
469 gxw_(:,k) = lmesh%Escale(:,k,1,1)*fx(:) + liftdelflx(:)
470
471 call sparsemat_matmul(dz, w_, fz)
472 call sparsemat_matmul(lift, lmesh%Fscale(:,k)*del_flux(:,k,vars_gzw_id), liftdelflx)
473 gzw_(:,k) = lmesh%Escale(:,k,2,2)*fz(:) + liftdelflx(:)
474
475 call sparsemat_matmul(dx, dtheta_, fx)
476 call sparsemat_matmul(lift, lmesh%Fscale(:,k)*del_flux(:,k,vars_gxpt_id), liftdelflx)
477 gxpt_(:,k) = lmesh%Escale(:,k,1,1)*fx(:) + liftdelflx(:)
478
479 call sparsemat_matmul(dz, dtheta_, fz)
480 call sparsemat_matmul(lift, lmesh%Fscale(:,k)*del_flux(:,k,vars_gzpt_id), liftdelflx)
481 gzpt_(:,k) = lmesh%Escale(:,k,2,2)*fz(:) + liftdelflx(:)
482 end do
483
484 return