FE-Project
All Classes Namespaces Files Functions Variables Pages
Functions/Subroutines
scale_atm_phy_tb_dgm_dns Module Reference

module FElib / Atmosphere / Physics turbulence More...

Functions/Subroutines

subroutine, public atm_phy_tb_dgm_dns_init (mesh)
 
subroutine, public atm_phy_tb_dgm_dns_final ()
 
subroutine, public atm_phy_tb_dgm_dns_cal_grad (t11, t12, t13, t21, t22, t23, t31, t32, t33, df1, df2, df3, tke, nu, kh, ddens_, momx_, momy_, momz_, drhot_, dens_hyd, pres_hyd, pres, pt, dx, dy, dz, sx, sy, sz, lift, lmesh, elem, lmesh2d, elem2d, is_bound)
 Calculate parameterized stress tensor and eddy heat flux with turbulent model.
 

Detailed Description

module FElib / Atmosphere / Physics turbulence

Description
Turbulnce process for DNS
Author
Team SCALE
Yuta Kawai, Team SCALE
NAMELIST
  • PARAM_ATMOS_PHY_TB_DGM_DNS
    nametypedefault valuecomment
    DNS_NU real(RP) 1.512E-5_RP kinematic viscosity coefficient [m2/s] for air at 20degC
    DNS_MU real(RP) 1.8E-5_RP molecular diffusive coefficient [m2/s] for air at 20degC

History Output
No history output

Function/Subroutine Documentation

◆ atm_phy_tb_dgm_dns_init()

subroutine, public scale_atm_phy_tb_dgm_dns::atm_phy_tb_dgm_dns_init ( class(meshbase3d), intent(in) mesh)

Definition at line 71 of file scale_atm_phy_tb_dgm_dns.F90.

72 implicit none
73 class(MeshBase3D), intent(in) :: mesh
74
75 logical :: consistent_tke = .true.
76
77 namelist / param_atmos_phy_tb_dgm_dns / &
78 dns_nu, dns_mu
79
80 integer :: ierr
81 !--------------------------------------------------------------------
82
83 log_newline
84 log_info("ATMOS_PHY_TB_dgm_dns_setup",*) 'Setup'
85 log_info("ATMOS_PHY_TB_dgm_dns_setup",*) 'Eddy Viscocity Model for DNS'
86
87 !--- read namelist
88 rewind(io_fid_conf)
89 read(io_fid_conf,nml=param_atmos_phy_tb_dgm_dns,iostat=ierr)
90 if( ierr < 0 ) then !--- missing
91 log_info("ATMOS_PHY_TB_dgm_dns_setup",*) 'Not found namelist. Default used.'
92 elseif( ierr > 0 ) then !--- fatal error
93 log_error("ATMOS_PHY_TB_dgm_dns_setup",*) 'Not appropriate names in namelist PARAM_ATMOS_PHY_TB_DGM_DNS. Check!'
94 call prc_abort
95 endif
96 log_nml(param_atmos_phy_tb_dgm_dns)
97
98 return

◆ atm_phy_tb_dgm_dns_final()

subroutine, public scale_atm_phy_tb_dgm_dns::atm_phy_tb_dgm_dns_final

Definition at line 102 of file scale_atm_phy_tb_dgm_dns.F90.

103 implicit none
104 !--------------------------------------------------------------------
105
106 return

◆ atm_phy_tb_dgm_dns_cal_grad()

subroutine, public scale_atm_phy_tb_dgm_dns::atm_phy_tb_dgm_dns_cal_grad ( real(rp), dimension(elem%np,lmesh%nea), intent(out) t11,
real(rp), dimension(elem%np,lmesh%nea), intent(out) t12,
real(rp), dimension(elem%np,lmesh%nea), intent(out) t13,
real(rp), dimension(elem%np,lmesh%nea), intent(out) t21,
real(rp), dimension(elem%np,lmesh%nea), intent(out) t22,
real(rp), dimension(elem%np,lmesh%nea), intent(out) t23,
real(rp), dimension(elem%np,lmesh%nea), intent(out) t31,
real(rp), dimension(elem%np,lmesh%nea), intent(out) t32,
real(rp), dimension(elem%np,lmesh%nea), intent(out) t33,
real(rp), dimension(elem%np,lmesh%nea), intent(out) df1,
real(rp), dimension(elem%np,lmesh%nea), intent(out) df2,
real(rp), dimension(elem%np,lmesh%nea), intent(out) df3,
real(rp), dimension(elem%np,lmesh%nea), intent(out) tke,
real(rp), dimension(elem%np,lmesh%nea), intent(out) nu,
real(rp), dimension(elem%np,lmesh%nea), intent(out) kh,
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) pres,
real(rp), dimension(elem%np,lmesh%nea), intent(in) pt,
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,
logical, dimension(elem%nfptot,lmesh%ne), intent(in) is_bound )

Calculate parameterized stress tensor and eddy heat flux with turbulent model.

Parameters
[out]t11(1,1) component of stress tensor
[out]t12(1,2) component of stress tensor
[out]t13(1,3) component of stress tensor
[out]t21(2,1) component of stress tensor
[out]t22(2,2) component of stress tensor
[out]t23(2,3) component of stress tensor
[out]t31(3,1) component of stress tensor
[out]t32(3,2) component of stress tensor
[out]t33(3,3) component of stress tensor
[out]df1Diffusive heat flux in x1 direction / density
[out]df2Diffusive heat flux in x2 direction / density
[out]df3Diffusive heat flux in x3 direction / density
[out]tkeParameterized turbulent kinetic energy
[out]nuEddy viscosity
[out]khEddy diffusivity
[in]ddens_Density perturbation
[in]momx_Momentum in x1 direction
[in]momy_Momentum in x2 direction
[in]momz_Momentum in x3 direction
[in]drhot_Density x potential temperature perturbation
[in]dens_hydReference pressure in hydrostatic balance
[in]pres_hydReference density in hydrostatic balance
[in]presPressure
[in]ptPotential temperature
[in]dzDifferential matrix managed by sparse matrix type
[in]szStiffness matrix managed by sparse matrix type
[in]liftLifting matrix managed by sparse matrix type
[in]is_boundFlag whether nodes are located at domain boundaries

Definition at line 112 of file scale_atm_phy_tb_dgm_dns.F90.

120
121 implicit none
122
123 class(LocalMesh3D), intent(in) :: lmesh
124 class(ElementBase3D), intent(in) :: elem
125 class(LocalMesh2D), intent(in) :: lmesh2D
126 class(ElementBase2D), intent(in) :: elem2D
127 real(RP), intent(out) :: T11(elem%Np,lmesh%NeA)
128 real(RP), intent(out) :: T12(elem%Np,lmesh%NeA)
129 real(RP), intent(out) :: T13(elem%Np,lmesh%NeA)
130 real(RP), intent(out) :: T21(elem%Np,lmesh%NeA)
131 real(RP), intent(out) :: T22(elem%Np,lmesh%NeA)
132 real(RP), intent(out) :: T23(elem%Np,lmesh%NeA)
133 real(RP), intent(out) :: T31(elem%Np,lmesh%NeA)
134 real(RP), intent(out) :: T32(elem%Np,lmesh%NeA)
135 real(RP), intent(out) :: T33(elem%Np,lmesh%NeA)
136 real(RP), intent(out) :: DF1(elem%Np,lmesh%NeA)
137 real(RP), intent(out) :: DF2(elem%Np,lmesh%NeA)
138 real(RP), intent(out) :: DF3(elem%Np,lmesh%NeA)
139 real(RP), intent(out) :: TKE(elem%Np,lmesh%NeA)
140 real(RP), intent(out) :: Nu(elem%Np,lmesh%NeA)
141 real(RP), intent(out) :: Kh(elem%Np,lmesh%NeA)
142 real(RP), intent(in) :: DDENS_(elem%Np,lmesh%NeA)
143 real(RP), intent(in) :: MOMX_ (elem%Np,lmesh%NeA)
144 real(RP), intent(in) :: MOMY_ (elem%Np,lmesh%NeA)
145 real(RP), intent(in) :: MOMZ_ (elem%Np,lmesh%NeA)
146 real(RP), intent(in) :: DRHOT_(elem%Np,lmesh%NeA)
147 real(RP), intent(in) :: DENS_hyd(elem%Np,lmesh%NeA)
148 real(RP), intent(in) :: PRES_hyd(elem%Np,lmesh%NeA)
149 real(RP), intent(in) :: PRES(elem%Np,lmesh%NeA)
150 real(RP), intent(in) :: PT(elem%Np,lmesh%NeA)
151 type(SparseMat), intent(in) :: Dx, Dy, Dz
152 type(SparseMat), intent(in) :: Sx, Sy, Sz
153 type(SparseMat), intent(in) :: Lift
154 logical, intent(in) :: is_bound(elem%NfpTot,lmesh%Ne)
155
156 real(RP) :: Fx(elem%Np), Fy(elem%Np), Fz(elem%Np), LiftDelFlx(elem%Np)
157 real(RP) :: DENS(elem%Np), RDENS(elem%Np), RHOT(elem%Np), Q(elem%Np)
158 real(RP) :: DdensDxi(elem%Np,3)
159 real(RP) :: DVelDxi(elem%Np,3,3)
160 real(RP) :: del_flux_rho (elem%NfpTot,lmesh%Ne,3)
161 real(RP) :: del_flux_mom (elem%NfpTot,lmesh%Ne,3,3)
162 real(RP) :: del_flux_rhot(elem%NfpTot,lmesh%Ne,3)
163
164 real(RP) :: S11(elem%Np), S12(elem%Np), S22(elem%Np), S23(elem%Np), S31(elem%Np), S33(elem%Np)
165 real(RP) :: SkkOvThree
166 real(RP) :: coef
167
168 integer :: ke
169 integer :: p
170 !--------------------------------------------------------------------
171
172 call cal_del_flux_grad( del_flux_rho, del_flux_mom, del_flux_rhot, & ! (out)
173 ddens_, momx_, momy_, momz_, drhot_, dens_hyd, pres_hyd, pt, & ! (in)
174 lmesh%normal_fn(:,:,1), lmesh%normal_fn(:,:,2), lmesh%normal_fn(:,:,3), & ! (in)
175 lmesh%vmapM, lmesh%vmapP, & ! (in)
176 lmesh, elem, is_bound ) ! (in)
177
178 !$omp parallel do private( &
179 !$omp Fx, Fy, Fz, LiftDelFlx, &
180 !$omp DENS, RHOT, RDENS, Q, DdensDxi, DVelDxi, &
181 !$omp S11, S12, S22, S23, S31, S33, &
182 !$omp coef, SkkOvThree, &
183 !$omp p )
184 do ke=lmesh%NeS, lmesh%NeE
185 !---
186 dens(:) = dens_hyd(:,ke) + ddens_(:,ke)
187 rdens(:) = 1.0_rp / dens(:)
188 rhot(:) = dens(:) * pt(:,ke)
189
190 ! gradient of density
191 call sparsemat_matmul( dx, dens, fx )
192 call sparsemat_matmul( lift, lmesh%Fscale(:,ke) * del_flux_rho(:,ke,1), liftdelflx )
193 ddensdxi(:,1) = lmesh%Escale(:,ke,1,1) * fx(:) + liftdelflx(:)
194
195 call sparsemat_matmul( dy, dens, fy )
196 call sparsemat_matmul( lift, lmesh%Fscale(:,ke) * del_flux_rho(:,ke,2), liftdelflx )
197 ddensdxi(:,2) = lmesh%Escale(:,ke,2,2) * fy(:) + liftdelflx(:)
198
199 call sparsemat_matmul( dz, dens, fz )
200 call sparsemat_matmul( lift, lmesh%Fscale(:,ke) * del_flux_rho(:,ke,3), liftdelflx )
201 ddensdxi(:,3) = lmesh%Escale(:,ke,3,3) * fz(:) + liftdelflx(:)
202
203 ! gradient of u
204 q(:) = momx_(:,ke) * rdens(:)
205
206 call sparsemat_matmul( dx, momx_(:,ke), fx )
207 call sparsemat_matmul( lift, lmesh%Fscale(:,ke) * del_flux_mom(:,ke,1,1), liftdelflx )
208 dveldxi(:,1,1) = ( lmesh%Escale(:,ke,1,1) * fx(:) + liftdelflx(:) - q(:) * ddensdxi(:,1) ) * rdens(:)
209
210 call sparsemat_matmul( dy, momx_(:,ke), fy )
211 call sparsemat_matmul( lift, lmesh%Fscale(:,ke) * del_flux_mom(:,ke,2,1), liftdelflx )
212 dveldxi(:,2,1) = ( lmesh%Escale(:,ke,2,2) * fy(:) + liftdelflx(:) - q(:) * ddensdxi(:,2) ) * rdens(:)
213
214 call sparsemat_matmul( dz, momx_(:,ke), fz )
215 call sparsemat_matmul( lift, lmesh%Fscale(:,ke) * del_flux_mom(:,ke,3,1), liftdelflx )
216 dveldxi(:,3,1) = ( lmesh%Escale(:,ke,3,3) * fz(:) + liftdelflx(:) - q(:) * ddensdxi(:,3) ) * rdens(:)
217
218 ! gradient of v
219 q(:) = momy_(:,ke) * rdens(:)
220
221 call sparsemat_matmul( dx, momy_(:,ke), fx )
222 call sparsemat_matmul( lift, lmesh%Fscale(:,ke) * del_flux_mom(:,ke,1,2), liftdelflx )
223 dveldxi(:,1,2) = ( lmesh%Escale(:,ke,1,1) * fx(:) + liftdelflx(:) - q(:) * ddensdxi(:,1) ) * rdens(:)
224
225 call sparsemat_matmul( dy, momy_(:,ke), fy )
226 call sparsemat_matmul( lift, lmesh%Fscale(:,ke) * del_flux_mom(:,ke,2,2), liftdelflx )
227 dveldxi(:,2,2) = ( lmesh%Escale(:,ke,2,2) * fy(:) + liftdelflx(:) - q(:) * ddensdxi(:,2) ) * rdens(:)
228
229 call sparsemat_matmul( dz, momy_(:,ke), fz )
230 call sparsemat_matmul( lift, lmesh%Fscale(:,ke) * del_flux_mom(:,ke,3,2), liftdelflx )
231 dveldxi(:,3,2) = ( lmesh%Escale(:,ke,3,3) * fz(:) + liftdelflx(:) - q(:) * ddensdxi(:,3) ) * rdens(:)
232
233 ! gradient of w
234 q(:) = momz_(:,ke) * rdens(:)
235
236 call sparsemat_matmul( dx, momz_(:,ke), fx )
237 call sparsemat_matmul( lift, lmesh%Fscale(:,ke) * del_flux_mom(:,ke,1,3), liftdelflx )
238 dveldxi(:,1,3) = ( lmesh%Escale(:,ke,1,1) * fx(:) + liftdelflx(:) - q(:) * ddensdxi(:,1) ) * rdens(:)
239
240 call sparsemat_matmul( dy, momz_(:,ke), fy )
241 call sparsemat_matmul( lift, lmesh%Fscale(:,ke) * del_flux_mom(:,ke,2,3), liftdelflx )
242 dveldxi(:,2,3) = ( lmesh%Escale(:,ke,2,2) * fy(:) + liftdelflx(:) - q(:) * ddensdxi(:,2) ) * rdens(:)
243
244 call sparsemat_matmul( dz, momz_(:,ke), fz )
245 call sparsemat_matmul( lift, lmesh%Fscale(:,ke) * del_flux_mom(:,ke,3,3), liftdelflx )
246 dveldxi(:,3,3) = ( lmesh%Escale(:,ke,3,3) * fz(:) + liftdelflx(:) - q(:) * ddensdxi(:,3) ) * rdens(:)
247
248 ! gradient of pt
249 q(:) = rhot(:) * rdens(:)
250
251 call sparsemat_matmul( dx, rhot, fx )
252 call sparsemat_matmul( lift, lmesh%Fscale(:,ke) * del_flux_rhot(:,ke,1), liftdelflx )
253 df1(:,ke) = ( lmesh%Escale(:,ke,1,1) * fx(:) + liftdelflx(:) - q(:) * ddensdxi(:,1) ) * rdens(:)
254
255 call sparsemat_matmul( dy, rhot, fy )
256 call sparsemat_matmul( lift, lmesh%Fscale(:,ke) * del_flux_rhot(:,ke,2), liftdelflx )
257 df2(:,ke) = ( lmesh%Escale(:,ke,2,2) * fy(:) + liftdelflx(:) - q(:) * ddensdxi(:,2) ) * rdens(:)
258
259 call sparsemat_matmul( dz, rhot, fz )
260 call sparsemat_matmul( lift, lmesh%Fscale(:,ke) * del_flux_rhot(:,ke,3), liftdelflx )
261 df3(:,ke) = ( lmesh%Escale(:,ke,3,3) * fz(:) + liftdelflx(:) - q(:) * ddensdxi(:,3) ) * rdens(:)
262
263
264 ! Calculate the component of strain velocity tensor
265 s11(:) = dveldxi(:,1,1)
266 s12(:) = 0.5_rp * ( dveldxi(:,1,2) + dveldxi(:,2,1) )
267 s22(:) = dveldxi(:,2,2)
268 s23(:) = 0.5_rp * ( dveldxi(:,2,3) + dveldxi(:,3,2) )
269 s31(:) = 0.5_rp * ( dveldxi(:,1,3) + dveldxi(:,3,1) )
270 s33(:) = dveldxi(:,3,3)
271
272 ! Caclulate eddy viscosity & eddy diffusivity
273
274 nu(:,ke) = dns_nu
275 kh(:,ke) = dns_mu
276
277 !---
278
279 do p=1, elem%Np
280 skkovthree = ( s11(p) + s22(p) + s33(p) ) * oneoverthree
281 coef = 2.0_rp * nu(p,ke)
282
283 t11(p,ke) = dens(p) * ( coef * ( s11(p) - skkovthree ) )
284 t12(p,ke) = dens(p) * coef * s12(p)
285 t13(p,ke) = dens(p) * coef * s31(p)
286
287 t21(p,ke) = dens(p) * coef * s12(p)
288 t22(p,ke) = dens(p) * ( coef * ( s22(p) - skkovthree ) )
289 t23(p,ke) = dens(p) * coef * s23(p)
290
291 t31(p,ke) = dens(p) * coef * s31(p)
292 t32(p,ke) = dens(p) * coef * s23(p)
293 t33(p,ke) = dens(p) * ( coef * ( s33(p) - skkovthree ) )
294 end do
295
296 df1(:,ke) = kh(:,ke) * df1(:,ke)
297 df2(:,ke) = kh(:,ke) * df2(:,ke)
298 df3(:,ke) = kh(:,ke) * df3(:,ke)
299 end do
300
301 return