54 real(RP),
private :: atmos_phy_tb_dns_nu = 1.512e-5_rp
56 real(RP),
private :: atmos_phy_tb_dns_mu = 1.512e-5_rp
69 character(len=*),
intent(in) :: TYPE_TB
71 real(RP),
intent(in) :: CDZ(
ka)
72 real(RP),
intent(in) :: CDX(
ia)
73 real(RP),
intent(in) :: CDY(
ja)
74 real(RP),
intent(in) :: CZ (
ka,
ia,
ja)
76 namelist / param_atmos_phy_tb_dns / &
77 atmos_phy_tb_dns_nu, &
85 if ( type_tb /=
'DNS' )
then 86 write(*,*)
'xxx ATMOS_PHY_TB_TYPE is not DNS. Check!' 92 read(
io_fid_conf,nml=param_atmos_phy_tb_dns,iostat=ierr)
94 if(
io_l )
write(
io_fid_log,*)
'*** Not found namelist. Default used.' 95 elseif( ierr > 0 )
then 96 write(*,*)
'xxx Not appropriate names in namelist PARAM_ATMOS_PHY_TB_DNS. Check!' 106 qflx_sgs_MOMZ, qflx_sgs_MOMX, qflx_sgs_MOMY, &
107 qflx_sgs_rhot, qflx_sgs_rhoq, &
108 tke, tke_t, nu, Ri, Pr, N2, &
109 MOMZ, MOMX, MOMY, RHOT, DENS, QTRC, &
110 SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_QV, &
111 GSQRT, J13G, J23G, J33G, MAPF, dt )
130 real(RP),
intent(out) :: qflx_sgs_MOMZ(
ka,
ia,
ja,3)
131 real(RP),
intent(out) :: qflx_sgs_MOMX(
ka,
ia,
ja,3)
132 real(RP),
intent(out) :: qflx_sgs_MOMY(
ka,
ia,
ja,3)
133 real(RP),
intent(out) :: qflx_sgs_rhot(
ka,
ia,
ja,3)
134 real(RP),
intent(out) :: qflx_sgs_rhoq(
ka,
ia,
ja,3,
qa)
136 real(RP),
intent(inout) :: tke(
ka,
ia,
ja)
137 real(RP),
intent(out) :: tke_t(
ka,
ia,
ja)
138 real(RP),
intent(out) :: nu (
ka,
ia,
ja)
139 real(RP),
intent(out) :: Ri (
ka,
ia,
ja)
140 real(RP),
intent(out) :: Pr (
ka,
ia,
ja)
141 real(RP),
intent(out) :: N2 (
ka,
ia,
ja)
143 real(RP),
intent(in) :: MOMZ(
ka,
ia,
ja)
144 real(RP),
intent(in) :: MOMX(
ka,
ia,
ja)
145 real(RP),
intent(in) :: MOMY(
ka,
ia,
ja)
146 real(RP),
intent(in) :: RHOT(
ka,
ia,
ja)
147 real(RP),
intent(in) :: DENS(
ka,
ia,
ja)
148 real(RP),
intent(in) :: QTRC(
ka,
ia,
ja,
qa)
150 real(RP),
intent(in) :: SFLX_MW(
ia,
ja)
151 real(RP),
intent(in) :: SFLX_MU(
ia,
ja)
152 real(RP),
intent(in) :: SFLX_MV(
ia,
ja)
153 real(RP),
intent(in) :: SFLX_SH(
ia,
ja)
154 real(RP),
intent(in) :: SFLX_QV(
ia,
ja)
156 real(RP),
intent(in) :: GSQRT (
ka,
ia,
ja,7)
157 real(RP),
intent(in) :: J13G (
ka,
ia,
ja,7)
158 real(RP),
intent(in) :: J23G (
ka,
ia,
ja,7)
159 real(RP),
intent(in) :: J33G
160 real(RP),
intent(in) :: MAPF (
ia,
ja,2,4)
161 real(DP),
intent(in) :: dt
163 real(RP) :: POTT(
ka,
ia,
ja)
168 integer :: k, i, j, iq
172 qflx_sgs_momz(:,:,:,:) = undef
173 qflx_sgs_momx(:,:,:,:) = undef
174 qflx_sgs_momy(:,:,:,:) = undef
175 qflx_sgs_rhot(:,:,:,:) = undef
176 qflx_sgs_rhoq(:,:,:,:,:) = undef
181 tke_t(:,:,:) = 0.0_rp
192 pott(k,i,j) = rhot(k,i,j) / dens(k,i,j)
208 qflx_sgs_momz(k,i,j,
zdir) = -atmos_phy_tb_dns_nu * ( momz(k,i,j)-momz(k-1,i,j) ) * rcdz(k)
215 qflx_sgs_momz(
ks,i,j,
zdir) = 0.0_rp
216 qflx_sgs_momz(
ke,i,j,
zdir) = 0.0_rp
224 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)
233 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)
243 qflx_sgs_momx(k,i,j,
zdir) = -atmos_phy_tb_dns_nu * ( momx(k+1,i,j)-momx(k,i,j) ) * rfdz(k)
250 qflx_sgs_momx(
ks-1,i,j,
zdir) = 0.0_rp
251 qflx_sgs_momx(
ke ,i,j,
zdir) = 0.0_rp
259 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)
268 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)
279 qflx_sgs_momy(k,i,j,
zdir) = -atmos_phy_tb_dns_nu * ( momy(k+1,i,j)-momy(k,i,j) ) * rfdz(k)
286 qflx_sgs_momy(
ks-1,i,j,
zdir) = 0.0_rp
287 qflx_sgs_momy(
ke ,i,j,
zdir) = 0.0_rp
295 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)
304 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)
315 qflx_sgs_rhot(k,i,j,
zdir) = -0.5_rp * ( dens(k+1,i,j)+dens(k,i,j) ) &
316 * atmos_phy_tb_dns_mu * ( pott(k+1,i,j)-pott(k,i,j) ) * rfdz(k)
323 qflx_sgs_rhot(
ks-1,i,j,
zdir) = 0.0_rp
324 qflx_sgs_rhot(
ke ,i,j,
zdir) = 0.0_rp
332 qflx_sgs_rhot(k,i,j,
xdir) = -0.5_rp * ( dens(k,i+1,j)+dens(k,i,j) ) &
333 * atmos_phy_tb_dns_mu * ( pott(k,i+1,j)-pott(k,i,j) ) * rfdx(i) * mapf(i,j,1,
i_xy)
342 qflx_sgs_rhot(k,i,j,
ydir) = -0.5_rp * ( dens(k,i,j+1)+dens(k,i,j) ) &
343 * atmos_phy_tb_dns_mu * ( pott(k,i,j+1)-pott(k,i,j) ) * rfdy(j) * mapf(i,j,2,
i_xy)
363 qflx_sgs_rhoq(k,i,j,
zdir,iq) = -0.5_rp * ( dens(k+1,i,j)+dens(k,i,j) ) &
364 * atmos_phy_tb_dns_mu * ( qtrc(k+1,i,j,iq)-qtrc(k,i,j,iq) ) * rfdz(k)
370 qflx_sgs_rhoq(
ks-1,i,j,
zdir,iq) = 0.0_rp
371 qflx_sgs_rhoq(
ke ,i,j,
zdir,iq) = 0.0_rp
379 qflx_sgs_rhoq(k,i,j,
xdir,iq) = -0.5_rp * ( dens(k,i+1,j)+dens(k,i,j) ) &
380 * atmos_phy_tb_dns_mu * ( qtrc(k,i+1,j,iq)-qtrc(k,i,j,iq) ) * rfdx(i) * mapf(i,j,1,
i_xy)
389 qflx_sgs_rhoq(k,i,j,
ydir,iq) = -0.5_rp * ( dens(k,i,j+1)+dens(k,i,j) ) &
390 * 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.
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
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
subroutine, public atmos_phy_tb_dns(qflx_sgs_MOMZ, qflx_sgs_MOMX, qflx_sgs_MOMY, qflx_sgs_rhot, qflx_sgs_rhoq, tke, tke_t, nu, Ri, Pr, N2, MOMZ, MOMX, MOMY, RHOT, DENS, QTRC, SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_QV, GSQRT, J13G, J23G, J33G, MAPF, dt)
real(rp), public const_undef
integer, public ia
of x whole cells (local, with HALO)
integer, public ka
of z whole cells (local, with HALO)
subroutine, public atmos_phy_tb_dns_setup(TYPE_TB, CDZ, CDX, CDY, CZ)
integer, public jblock
block size for cache blocking: y
real(rp), public const_grav
standard acceleration of gravity [m/s2]
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 ja
of y whole cells (local, with HALO)