55 real(RP),
private :: ATMOS_PHY_TB_DNS_NU = 1.512e-5_rp
57 real(RP),
private :: ATMOS_PHY_TB_DNS_MU = 1.512e-5_rp
70 character(len=*),
intent(in) :: type_tb
71 integer,
intent(out) :: i_tke_out
75 if(
io_l )
write(
io_fid_log,*)
'++++++ Module[Turbulence Tracer] / Categ[ATMOS PHYSICS] / Origin[SCALElib]' 76 if(
io_l )
write(
io_fid_log,*)
'*** Tracers for Deardorff (1980) 1.5th TKE Model' 78 if ( type_tb /=
'DNS' )
then 79 write(*,*)
'xxx ATMOS_PHY_TB_TYPE is not DNS. Check!' 94 real(RP),
intent(in) :: cdz(
ka)
95 real(RP),
intent(in) :: cdx(
ia)
96 real(RP),
intent(in) :: cdy(
ja)
97 real(RP),
intent(in) :: cz (
ka,
ia,
ja)
99 namelist / param_atmos_phy_tb_dns / &
100 atmos_phy_tb_dns_nu, &
107 if(
io_l )
write(
io_fid_log,*)
'++++++ Module[Turbulence] / Categ[ATMOS PHYSICS] / Origin[SCALElib]' 112 read(
io_fid_conf,nml=param_atmos_phy_tb_dns,iostat=ierr)
114 if(
io_l )
write(
io_fid_log,*)
'*** Not found namelist. Default used.' 115 elseif( ierr > 0 )
then 116 write(*,*)
'xxx Not appropriate names in namelist PARAM_ATMOS_PHY_TB_DNS. Check!' 126 qflx_sgs_MOMZ, qflx_sgs_MOMX, qflx_sgs_MOMY, &
127 qflx_sgs_rhot, qflx_sgs_rhoq, &
128 RHOQ_t, nu, Ri, Pr, &
129 MOMZ, MOMX, MOMY, RHOT, DENS, QTRC, N2, &
130 SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_Q, &
131 GSQRT, J13G, J23G, J33G, MAPF, dt )
150 real(RP),
intent(out) :: qflx_sgs_momz(
ka,
ia,
ja,3)
151 real(RP),
intent(out) :: qflx_sgs_momx(
ka,
ia,
ja,3)
152 real(RP),
intent(out) :: qflx_sgs_momy(
ka,
ia,
ja,3)
153 real(RP),
intent(out) :: qflx_sgs_rhot(
ka,
ia,
ja,3)
154 real(RP),
intent(out) :: qflx_sgs_rhoq(
ka,
ia,
ja,3,
qa)
156 real(RP),
intent(inout) :: rhoq_t(
ka,
ia,
ja,
qa)
158 real(RP),
intent(out) :: nu(
ka,
ia,
ja)
159 real(RP),
intent(out) :: ri(
ka,
ia,
ja)
160 real(RP),
intent(out) :: pr(
ka,
ia,
ja)
162 real(RP),
intent(in) :: momz(
ka,
ia,
ja)
163 real(RP),
intent(in) :: momx(
ka,
ia,
ja)
164 real(RP),
intent(in) :: momy(
ka,
ia,
ja)
165 real(RP),
intent(in) :: rhot(
ka,
ia,
ja)
166 real(RP),
intent(in) :: dens(
ka,
ia,
ja)
167 real(RP),
intent(in) :: qtrc(
ka,
ia,
ja,
qa)
168 real(RP),
intent(in) :: n2(
ka,
ia,
ja)
170 real(RP),
intent(in) :: sflx_mw(
ia,
ja)
171 real(RP),
intent(in) :: sflx_mu(
ia,
ja)
172 real(RP),
intent(in) :: sflx_mv(
ia,
ja)
173 real(RP),
intent(in) :: sflx_sh(
ia,
ja)
174 real(RP),
intent(in) :: sflx_q (
ia,
ja,
qa)
176 real(RP),
intent(in) :: gsqrt (
ka,
ia,
ja,7)
177 real(RP),
intent(in) :: j13g (
ka,
ia,
ja,7)
178 real(RP),
intent(in) :: j23g (
ka,
ia,
ja,7)
179 real(RP),
intent(in) :: j33g
180 real(RP),
intent(in) :: mapf (
ia,
ja,2,4)
181 real(DP),
intent(in) :: dt
183 real(RP) :: pott(
ka,
ia,
ja)
188 integer :: k, i, j, iq
192 qflx_sgs_momz(:,:,:,:) = undef
193 qflx_sgs_momx(:,:,:,:) = undef
194 qflx_sgs_momy(:,:,:,:) = undef
195 qflx_sgs_rhot(:,:,:,:) = undef
196 qflx_sgs_rhoq(:,:,:,:,:) = undef
209 pott(k,i,j) = rhot(k,i,j) / dens(k,i,j)
225 qflx_sgs_momz(k,i,j,
zdir) = -atmos_phy_tb_dns_nu * ( momz(k,i,j)-momz(k-1,i,j) ) * rcdz(k)
232 qflx_sgs_momz(
ks,i,j,
zdir) = 0.0_rp
233 qflx_sgs_momz(
ke,i,j,
zdir) = 0.0_rp
241 qflx_sgs_momz(k,i,j,
xdir) = -atmos_phy_tb_dns_nu * ( momz(k,i+1,j)-momz(k,i,j) ) * rfdx(i) * mapf(i,j,1,
i_xy)
250 qflx_sgs_momz(k,i,j,
ydir) = -atmos_phy_tb_dns_nu * ( momz(k,i,j+1)-momz(k,i,j) ) * rfdy(j) * mapf(i,j,2,
i_xy)
260 qflx_sgs_momx(k,i,j,
zdir) = -atmos_phy_tb_dns_nu * ( momx(k+1,i,j)-momx(k,i,j) ) * rfdz(k)
267 qflx_sgs_momx(
ks-1,i,j,
zdir) = 0.0_rp
268 qflx_sgs_momx(
ke ,i,j,
zdir) = 0.0_rp
276 qflx_sgs_momx(k,i,j,
xdir) = -atmos_phy_tb_dns_nu * ( momx(k,i,j)-momx(k,i-1,j) ) * rcdx(i) * mapf(i,j,1,
i_uy)
285 qflx_sgs_momx(k,i,j,
ydir) = -atmos_phy_tb_dns_nu * ( momx(k,i,j+1)-momx(k,i,j) ) * rfdy(j) * mapf(i,j,2,
i_uy)
296 qflx_sgs_momy(k,i,j,
zdir) = -atmos_phy_tb_dns_nu * ( momy(k+1,i,j)-momy(k,i,j) ) * rfdz(k)
303 qflx_sgs_momy(
ks-1,i,j,
zdir) = 0.0_rp
304 qflx_sgs_momy(
ke ,i,j,
zdir) = 0.0_rp
312 qflx_sgs_momy(k,i,j,
xdir) = -atmos_phy_tb_dns_nu * ( momy(k,i+1,j)-momy(k,i,j) ) * rfdx(i) * mapf(i,j,1,
i_xv)
321 qflx_sgs_momy(k,i,j,
ydir) = -atmos_phy_tb_dns_nu * ( momy(k,i,j)-momy(k,i,j-1) ) * rcdy(j) * mapf(i,j,2,
i_xv)
332 qflx_sgs_rhot(k,i,j,
zdir) = -0.5_rp * ( dens(k+1,i,j)+dens(k,i,j) ) &
333 * atmos_phy_tb_dns_mu * ( pott(k+1,i,j)-pott(k,i,j) ) * rfdz(k)
340 qflx_sgs_rhot(
ks-1,i,j,
zdir) = 0.0_rp
341 qflx_sgs_rhot(
ke ,i,j,
zdir) = 0.0_rp
349 qflx_sgs_rhot(k,i,j,
xdir) = -0.5_rp * ( dens(k,i+1,j)+dens(k,i,j) ) &
350 * atmos_phy_tb_dns_mu * ( pott(k,i+1,j)-pott(k,i,j) ) * rfdx(i) * mapf(i,j,1,
i_xy)
359 qflx_sgs_rhot(k,i,j,
ydir) = -0.5_rp * ( dens(k,i,j+1)+dens(k,i,j) ) &
360 * atmos_phy_tb_dns_mu * ( pott(k,i,j+1)-pott(k,i,j) ) * rfdy(j) * mapf(i,j,2,
i_xy)
382 qflx_sgs_rhoq(k,i,j,
zdir,iq) = -0.5_rp * ( dens(k+1,i,j)+dens(k,i,j) ) &
383 * atmos_phy_tb_dns_mu * ( qtrc(k+1,i,j,iq)-qtrc(k,i,j,iq) ) * rfdz(k)
389 qflx_sgs_rhoq(
ks-1,i,j,
zdir,iq) = 0.0_rp
390 qflx_sgs_rhoq(
ke ,i,j,
zdir,iq) = 0.0_rp
398 qflx_sgs_rhoq(k,i,j,
xdir,iq) = -0.5_rp * ( dens(k,i+1,j)+dens(k,i,j) ) &
399 * atmos_phy_tb_dns_mu * ( qtrc(k,i+1,j,iq)-qtrc(k,i,j,iq) ) * rfdx(i) * mapf(i,j,1,
i_xy)
408 qflx_sgs_rhoq(k,i,j,
ydir,iq) = -0.5_rp * ( dens(k,i,j+1)+dens(k,i,j) ) &
409 * atmos_phy_tb_dns_mu * ( qtrc(k,i,j+1,iq)-qtrc(k,i,j,iq) ) * rfdy(j) * mapf(i,j,2,
i_xy)
integer, public is
start point of inner domain: x, local
integer, public je
end point of inner domain: y, local
real(rp), dimension(:), allocatable, public grid_rcdy
reciprocal of center-dy
subroutine, public prc_mpistop
Abort MPI.
subroutine, public atmos_phy_tb_dns(qflx_sgs_MOMZ, qflx_sgs_MOMX, qflx_sgs_MOMY, qflx_sgs_rhot, qflx_sgs_rhoq, RHOQ_t, nu, Ri, Pr, MOMZ, MOMX, MOMY, RHOT, DENS, QTRC, N2, SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_Q, GSQRT, J13G, J23G, J33G, MAPF, dt)
integer, public iblock
block size for cache blocking: x
logical, public io_l
output log or not? (this process)
real(rp), dimension(:), allocatable, public grid_rcdx
reciprocal of center-dx
integer, parameter, public zdir
integer, parameter, public ydir
integer, public ke
end point of inner domain: z, local
integer, parameter, public xdir
logical, dimension(qa_max), public tracer_advc
real(rp), dimension(:), allocatable, public grid_rfdy
reciprocal of face-dy
subroutine, public check(current_line, v)
Undefined value checker.
real(rp), dimension(:), allocatable, public grid_rcdz
reciprocal of center-dz
real(rp), public const_undef
logical, public io_nml
output log or not? (for namelist, this process)
subroutine, public atmos_phy_tb_dns_setup(CDZ, CDX, CDY, CZ)
integer, public ia
of whole cells: x, local, with HALO
integer, public ka
of whole cells: z, local, with HALO
integer, public jblock
block size for cache blocking: y
real(rp), public const_grav
standard acceleration of gravity [m/s2]
subroutine, public atmos_phy_tb_dns_config(TYPE_TB, I_TKE_out)
Config.
integer, public js
start point of inner domain: y, local
integer, parameter, public const_undef2
undefined value (INT2)
module ATMOSPHERE / Physics Turbulence
integer, public ks
start point of inner domain: z, local
integer, public ie
end point of inner domain: x, local
real(rp), dimension(:), allocatable, public grid_rfdx
reciprocal of face-dx
real(rp), dimension(:), allocatable, public grid_rfdz
reciprocal of face-dz
integer, public io_fid_conf
Config file ID.
integer, public io_fid_log
Log file ID.
integer, public io_fid_nml
Log file ID (only for output namelist)
integer, public ja
of whole cells: y, local, with HALO