19 #include "inc_openmp.h" 52 integer,
private :: atmos_phy_sf_flg_mom_flux = 0
56 real(RP),
private,
parameter :: atmos_phy_sf_u_maxm = 100.0_rp
57 real(RP),
private :: atmos_phy_sf_u_minm = 0.0_rp
58 real(RP),
private,
parameter :: atmos_phy_sf_cm_max = 2.5e-3_rp
59 real(RP),
private :: atmos_phy_sf_cm_min = 1.0e-5_rp
61 real(RP),
private :: atmos_phy_sf_const_ustar = 0.25_rp
62 real(RP),
private :: atmos_phy_sf_const_cm = 0.0011_rp
63 real(RP),
private :: atmos_phy_sf_const_sh = 15.0_rp
64 real(RP),
private :: atmos_phy_sf_const_lh = 115.0_rp
66 logical,
private :: atmos_phy_sf_flg_sh_diurnal = .false.
67 real(RP),
private :: atmos_phy_sf_const_freq = 24.0_rp
78 character(len=*),
intent(in) :: ATMOS_PHY_SF_TYPE
80 namelist / param_atmos_phy_sf_const / &
81 atmos_phy_sf_flg_mom_flux, &
82 atmos_phy_sf_u_minm, &
83 atmos_phy_sf_cm_min, &
84 atmos_phy_sf_const_ustar, &
85 atmos_phy_sf_const_cm, &
86 atmos_phy_sf_const_sh, &
87 atmos_phy_sf_const_lh, &
88 atmos_phy_sf_flg_sh_diurnal, &
89 atmos_phy_sf_const_freq
95 if(
io_l )
write(
io_fid_log,*)
'++++++ Module[SURFACE FLUX] / Categ[ATMOS PHYSICS] / Origin[SCALElib]' 98 if ( atmos_phy_sf_type /=
'CONST' )
then 99 write(*,*)
'xxx ATMOS_PHY_SF_TYPE is not CONST. Check!' 105 read(
io_fid_conf,nml=param_atmos_phy_sf_const,iostat=ierr)
107 if(
io_l )
write(
io_fid_log,*)
'*** Not found namelist. Default used.' 108 elseif( ierr > 0 )
then 109 write(*,*)
'xxx Not appropriate names in namelist PARAM_ATMOS_PHY_SF_CONST. Check!' 120 ATM_TEMP, ATM_PRES, ATM_W, ATM_U, ATM_V, &
124 SFC_DENS, SFC_PRES, &
125 SFLX_LW_dn, SFLX_SW_dn, &
126 SFC_TEMP, SFC_albedo, &
127 SFC_Z0M, SFC_Z0H, SFC_Z0E, &
128 SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_LH, &
136 atmos_thermodyn_templhv
141 real(RP),
intent(in) :: ATM_TEMP (
ia,
ja)
142 real(RP),
intent(in) :: ATM_PRES (
ia,
ja)
143 real(RP),
intent(in) :: ATM_W (
ia,
ja)
144 real(RP),
intent(in) :: ATM_U (
ia,
ja)
145 real(RP),
intent(in) :: ATM_V (
ia,
ja)
146 real(RP),
intent(in) :: ATM_DENS (
ia,
ja)
147 real(RP),
intent(in) :: ATM_QTRC (
ia,
ja,
qa)
148 real(RP),
intent(in) :: ATM_Z1 (
ia,
ja)
149 real(DP),
intent(in) :: dt
150 real(RP),
intent(in) :: SFC_DENS (
ia,
ja)
151 real(RP),
intent(in) :: SFC_PRES (
ia,
ja)
152 real(RP),
intent(in) :: SFLX_LW_dn(
ia,
ja)
153 real(RP),
intent(in) :: SFLX_SW_dn(
ia,
ja)
154 real(RP),
intent(in) :: SFC_TEMP (
ia,
ja)
155 real(RP),
intent(in) :: SFC_albedo(
ia,
ja,2)
156 real(RP),
intent(inout) :: SFC_Z0M (
ia,
ja)
157 real(RP),
intent(inout) :: SFC_Z0H (
ia,
ja)
158 real(RP),
intent(inout) :: SFC_Z0E (
ia,
ja)
159 real(RP),
intent(out) :: SFLX_MW (
ia,
ja)
160 real(RP),
intent(out) :: SFLX_MU (
ia,
ja)
161 real(RP),
intent(out) :: SFLX_MV (
ia,
ja)
162 real(RP),
intent(out) :: SFLX_SH (
ia,
ja)
163 real(RP),
intent(out) :: SFLX_LH (
ia,
ja)
164 real(RP),
intent(out) :: SFLX_QTRC (
ia,
ja,
qa)
165 real(RP),
intent(out) :: U10 (
ia,
ja)
166 real(RP),
intent(out) :: V10 (
ia,
ja)
167 real(RP),
intent(out) :: T2 (
ia,
ja)
168 real(RP),
intent(out) :: Q2 (
ia,
ja)
170 real(RP) :: ATM_Uabs(
ia,
ja)
172 real(RP) :: Cm(
ia,
ja)
175 real(RP) :: modulation
177 real(RP) :: LHV(
ia,
ja)
181 if(
io_l )
write(
io_fid_log,*)
'*** Physics step: Surface flux(const)' 185 atm_uabs(i,j) = sqrt( atm_w(i,j)*atm_w(i,j) &
186 + atm_u(i,j)*atm_u(i,j) &
187 + atm_v(i,j)*atm_v(i,j) )
191 if ( atmos_phy_sf_flg_mom_flux == 0 )
then 194 cm(i,j) = atmos_phy_sf_const_cm
197 elseif( atmos_phy_sf_flg_mom_flux == 1 )
then 200 cm(i,j) = ( atmos_phy_sf_const_ustar / atm_uabs(i,j) )**2
201 cm(i,j) = min( max( cm(i,j), atmos_phy_sf_cm_min ), atmos_phy_sf_cm_max )
210 uabs_lim = min( max( atm_uabs(i,j), atmos_phy_sf_u_minm ), atmos_phy_sf_u_maxm )
212 sflx_mw(i,j) = -cm(i,j) * uabs_lim * sfc_dens(i,j) * atm_w(i,j)
213 sflx_mu(i,j) = -cm(i,j) * uabs_lim * sfc_dens(i,j) * atm_u(i,j)
214 sflx_mv(i,j) = -cm(i,j) * uabs_lim * sfc_dens(i,j) * atm_v(i,j)
220 if ( atmos_phy_sf_flg_sh_diurnal )
then 221 modulation = sin( 2.0_rp * pi *
time_nowsec / 3600.0_rp / atmos_phy_sf_const_freq )
228 sflx_sh(i,j) = atmos_phy_sf_const_sh * modulation
229 sflx_lh(i,j) = atmos_phy_sf_const_lh
234 call atmos_thermodyn_templhv( lhv, atm_temp )
236 sflx_qtrc(:,:,:) = 0.0_rp
239 sflx_qtrc(i,j,
i_qv) = sflx_lh(i,j) / lhv(i,j)
247 r10 = 10.0_rp / atm_z1(i,j)
249 u10(i,j) = r10 * atm_u(i,j)
250 v10(i,j) = r10 * atm_v(i,j)
256 t2(i,j) = atm_temp(i,j)
257 q2(i,j) = atm_qtrc(i,j,
i_qv)
integer, public is
start point of inner domain: x, local
integer, public je
end point of inner domain: y, local
subroutine, public prc_mpistop
Abort MPI.
real(dp), public time_nowsec
subday part of current time [sec]
logical, public io_l
output log or not? (this process)
subroutine, public atmos_phy_sf_const(ATM_TEMP, ATM_PRES, ATM_W, ATM_U, ATM_V, ATM_DENS, ATM_QTRC, ATM_Z1, dt, SFC_DENS, SFC_PRES, SFLX_LW_dn, SFLX_SW_dn, SFC_TEMP, SFC_albedo, SFC_Z0M, SFC_Z0H, SFC_Z0E, SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_LH, SFLX_QTRC, U10, V10, T2, Q2)
Constant flux.
integer, public ia
of x whole cells (local, with HALO)
integer, public js
start point of inner domain: y, local
integer, public ie
end point of inner domain: x, local
module ATMOSPHERE / Thermodynamics
logical, public io_lnml
output log or not? (for namelist, this process)
module ATMOSPHERE / Physics Surface fluxes
real(rp), public const_pi
pi
integer, public io_fid_conf
Config file ID.
integer, public io_fid_log
Log file ID.
subroutine, public atmos_phy_sf_const_setup(ATMOS_PHY_SF_TYPE)
Setup.
integer, public ja
of y whole cells (local, with HALO)