18 #include "inc_openmp.h" 51 real(RP),
private,
parameter :: atmos_phy_sf_u_maxm = 100.0_rp
52 real(RP),
private,
parameter :: atmos_phy_sf_u_maxh = 100.0_rp
53 real(RP),
private,
parameter :: atmos_phy_sf_u_maxe = 100.0_rp
54 real(RP),
private :: atmos_phy_sf_u_minm = 0.0_rp
55 real(RP),
private :: atmos_phy_sf_u_minh = 0.0_rp
56 real(RP),
private :: atmos_phy_sf_u_mine = 0.0_rp
69 character(len=*),
intent(in) :: ATMOS_PHY_SF_TYPE
71 namelist / param_atmos_phy_sf_bulk / &
72 atmos_phy_sf_u_minm, &
73 atmos_phy_sf_u_minh, &
80 if(
io_l )
write(
io_fid_log,*)
'++++++ Module[SURFACE FLUX] / Categ[ATMOS PHYSICS] / Origin[SCALElib]' 83 if ( atmos_phy_sf_type /=
'BULK' )
then 84 write(*,*)
'xxx ATMOS_PHY_SF_TYPE is not BULK. Check!' 90 read(
io_fid_conf,nml=param_atmos_phy_sf_bulk,iostat=ierr)
92 if(
io_l )
write(
io_fid_log,*)
'*** Not found namelist. Default used.' 93 elseif( ierr > 0 )
then 94 write(*,*)
'xxx Not appropriate names in namelist PARAM_ATMOS_PHY_SF_BULK. Check!' 99 call sf_bulkcoef_setup
107 ATM_TEMP, ATM_PRES, ATM_W, ATM_U, ATM_V, &
111 SFC_DENS, SFC_PRES, &
112 SFLX_LW_dn, SFLX_SW_dn, &
113 SFC_TEMP, SFC_albedo, SFC_beta, &
114 SFC_Z0M, SFC_Z0H, SFC_Z0E, &
115 SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_LH, &
126 saturation_pres2qsat_all => atmos_saturation_pres2qsat_all
128 atmos_thermodyn_templhv
135 real(RP),
intent(in) :: ATM_TEMP (
ia,
ja)
136 real(RP),
intent(in) :: ATM_PRES (
ia,
ja)
137 real(RP),
intent(in) :: ATM_W (
ia,
ja)
138 real(RP),
intent(in) :: ATM_U (
ia,
ja)
139 real(RP),
intent(in) :: ATM_V (
ia,
ja)
140 real(RP),
intent(in) :: ATM_DENS (
ia,
ja)
141 real(RP),
intent(in) :: ATM_QTRC (
ia,
ja,
qa)
142 real(RP),
intent(in) :: ATM_Z1 (
ia,
ja)
143 real(DP),
intent(in) :: dt
144 real(RP),
intent(in) :: SFC_DENS (
ia,
ja)
145 real(RP),
intent(in) :: SFC_PRES (
ia,
ja)
146 real(RP),
intent(in) :: SFLX_LW_dn(
ia,
ja)
147 real(RP),
intent(in) :: SFLX_SW_dn(
ia,
ja)
148 real(RP),
intent(in) :: SFC_TEMP (
ia,
ja)
149 real(RP),
intent(in) :: SFC_albedo(
ia,
ja,2)
150 real(RP),
intent(in) :: SFC_beta (
ia,
ja)
151 real(RP),
intent(inout) :: SFC_Z0M (
ia,
ja)
152 real(RP),
intent(inout) :: SFC_Z0H (
ia,
ja)
153 real(RP),
intent(inout) :: SFC_Z0E (
ia,
ja)
154 real(RP),
intent(out) :: SFLX_MW (
ia,
ja)
155 real(RP),
intent(out) :: SFLX_MU (
ia,
ja)
156 real(RP),
intent(out) :: SFLX_MV (
ia,
ja)
157 real(RP),
intent(out) :: SFLX_SH (
ia,
ja)
158 real(RP),
intent(out) :: SFLX_LH (
ia,
ja)
159 real(RP),
intent(out) :: SFLX_QTRC (
ia,
ja,
qa)
160 real(RP),
intent(out) :: U10 (
ia,
ja)
161 real(RP),
intent(out) :: V10 (
ia,
ja)
162 real(RP),
intent(out) :: T2 (
ia,
ja)
163 real(RP),
intent(out) :: Q2 (
ia,
ja)
165 real(RP) :: SFC_Z0M_t(
ia,
ja)
166 real(RP) :: SFC_Z0H_t(
ia,
ja)
167 real(RP) :: SFC_Z0E_t(
ia,
ja)
169 real(RP) :: SFC_QSAT(
ia,
ja)
171 real(RP) :: PBL(
ia,
ja)
172 real(RP) :: LHV (
ia,
ja)
181 if(
io_l )
write(
io_fid_log,*)
'*** Physics step: Surface flux(bulk)' 194 sfc_z0m(:,:) = sfc_z0m(:,:) + sfc_z0m_t(:,:) * dt
195 sfc_z0h(:,:) = sfc_z0h(:,:) + sfc_z0h_t(:,:) * dt
196 sfc_z0e(:,:) = sfc_z0e(:,:) + sfc_z0e_t(:,:) * dt
198 call saturation_pres2qsat_all( sfc_qsat(:,:), &
202 call atmos_thermodyn_templhv( lhv, atm_temp )
204 sflx_qtrc(:,:,:) = 0.0_rp
218 atm_qtrc(i,j,
i_qv), &
229 sflx_mw(i,j) = -atm_dens(i,j) * ustar**2 / uabs * atm_w(i,j)
230 sflx_mu(i,j) = -atm_dens(i,j) * ustar**2 / uabs * atm_u(i,j)
231 sflx_mv(i,j) = -atm_dens(i,j) * ustar**2 / uabs * atm_v(i,j)
234 sflx_sh(i,j) = -cpdry * atm_dens(i,j) * ustar * tstar
235 sflx_lh(i,j) = -lhv(i,j) * atm_dens(i,j) * ustar * qstar * sfc_beta(i,j)
238 sflx_qtrc(i,j,
i_qv) = sflx_lh(i,j) / lhv(i,j)
246 u10(i,j) = atm_u(i,j) *
log( 10.0_rp / sfc_z0m(i,j) ) /
log( atm_z1(i,j) / sfc_z0m(i,j) )
247 v10(i,j) = atm_v(i,j) *
log( 10.0_rp / sfc_z0m(i,j) ) /
log( atm_z1(i,j) / sfc_z0m(i,j) )
248 t2(i,j) = sfc_temp(i,j) + ( atm_temp(i,j) - sfc_temp(i,j) ) &
249 * (
log( 2.0_rp / sfc_z0m(i,j) ) *
log( 2.0_rp / sfc_z0h(i,j) ) ) &
250 / (
log( atm_z1(i,j) / sfc_z0m(i,j) ) *
log( atm_z1(i,j) / sfc_z0h(i,j) ) )
251 q2(i,j) = sfc_qsat(i,j) + ( atm_qtrc(i,j,
i_qv) - sfc_qsat(i,j) ) &
252 * (
log( 2.0_rp / sfc_z0m(i,j) ) *
log( 2.0_rp / sfc_z0e(i,j) ) ) &
253 / (
log( atm_z1(i,j) / sfc_z0m(i,j) ) *
log( atm_z1(i,j) / sfc_z0e(i,j) ) )
integer, public is
start point of inner domain: x, local
integer, public je
end point of inner domain: y, local
subroutine, public atmos_phy_sf_bulk_setup(ATMOS_PHY_SF_TYPE)
Setup.
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
module ATMOSPHERE / Saturation adjustment
subroutine, public prc_mpistop
Abort MPI.
logical, public io_l
output log or not? (this process)
subroutine, public atmos_phy_sf_bulkcoef_setup
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
subroutine, public atmos_phy_sf_bulk(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_beta, SFC_Z0M, SFC_Z0H, SFC_Z0E, SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_LH, SFLX_QTRC, U10, V10, T2, Q2)
Calculate surface flux.
integer, public ia
of x whole cells (local, with HALO)
procedure(bc), pointer, public bulkflux
integer, public js
start point of inner domain: y, local
module ATMOSPHERE / Physics Surface fluxes
procedure(rl), pointer, public roughness
procedure(bc), pointer, public atmos_phy_sf_bulkcoef
subroutine, public log(type, message)
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 Surface roughness length
integer, public io_fid_conf
Config file ID.
module ATMOSPHERE / Physics Surface bulk coefficient
integer, public io_fid_log
Log file ID.
integer, public ja
of y whole cells (local, with HALO)