13 #include "inc_openmp.h" 81 if(
io_l )
write(
io_fid_log,*)
'++++++ Module[DRIVER] / Categ[ATMOS PHY_SF] / Origin[SCALE-RM]' 94 if(
io_l )
write(
io_fid_log,*)
'*** this component is never called.' 95 if(
io_l )
write(
io_fid_log,*)
'*** surface fluxes are set to zero.' 100 sflx_lh(:,:) = 0.0_rp
101 sflx_qtrc(:,:,:) = 0.0_rp
103 if(
io_l )
write(
io_fid_log,*)
'*** SFC_TEMP, SFC_albedo is set in ATMOS_PHY_SF_vars.' 164 barometric_law_mslp => atmos_hydrostatic_barometric_law_mslp
166 thermodyn_qd => atmos_thermodyn_qd, &
167 thermodyn_cp => atmos_thermodyn_cp
217 logical,
intent(in) :: update_flag
219 real(RP) :: Uabs10(
ia,
ja)
220 real(RP) :: MSLP (
ia,
ja)
233 if ( update_flag )
then 236 call bottom_estimate(
dens(:,:,:), &
278 uabs10(i,j) = sqrt( u10(i,j)**2 + v10(i,j)**2 )
282 call barometric_law_mslp( mslp(:,:), sfc_pres(:,:), t2(:,:),
topo_zsfc(:,:) )
284 call hist_in( sfc_dens(:,:),
'SFC_DENS',
'surface atmospheric density',
'kg/m3' )
285 call hist_in( sfc_pres(:,:),
'SFC_PRES',
'surface atmospheric pressure',
'Pa' )
286 call hist_in( sfc_temp(:,:),
'SFC_TEMP',
'surface skin temperature (merged)',
'K' )
287 call hist_in( sfc_albedo(:,:,i_lw),
'SFC_ALB_LW',
'surface albedo (longwave,merged)',
'0-1' , nohalo=.true. )
288 call hist_in( sfc_albedo(:,:,i_sw),
'SFC_ALB_SW',
'surface albedo (shortwave,merged)',
'0-1' , nohalo=.true. )
289 call hist_in( sfc_z0m(:,:),
'SFC_Z0M',
'roughness length (momentum)',
'm' , nohalo=.true. )
290 call hist_in( sfc_z0h(:,:),
'SFC_Z0H',
'roughness length (heat)',
'm' , nohalo=.true. )
291 call hist_in( sfc_z0e(:,:),
'SFC_Z0E',
'roughness length (vapor)',
'm' , nohalo=.true. )
292 call hist_in( sflx_mw(:,:),
'MWFLX',
'w-momentum flux (merged)',
'kg/m2/s' )
293 call hist_in( sflx_mu(:,:),
'MUFLX',
'u-momentum flux (merged)',
'kg/m2/s' )
294 call hist_in( sflx_mv(:,:),
'MVFLX',
'v-momentum flux (merged)',
'kg/m2/s' )
295 call hist_in( sflx_sh(:,:),
'SHFLX',
'sensible heat flux (merged)',
'W/m2' , nohalo=.true. )
296 call hist_in( sflx_lh(:,:),
'LHFLX',
'latent heat flux (merged)',
'W/m2' , nohalo=.true. )
297 call hist_in( sflx_gh(:,:),
'GHFLX',
'ground heat flux (merged)',
'W/m2' , nohalo=.true. )
298 call hist_in( uabs10(:,:),
'Uabs10',
'10m absolute wind',
'm/s' , nohalo=.true. )
299 call hist_in( u10(:,:),
'U10',
'10m x-wind',
'm/s' , nohalo=.true. )
300 call hist_in( v10(:,:),
'V10',
'10m y-wind',
'm/s' , nohalo=.true. )
301 call hist_in( t2(:,:),
'T2 ',
'2m air temperature',
'K' , nohalo=.true. )
302 call hist_in( q2(:,:),
'Q2 ',
'2m specific humidity',
'kg/kg' , nohalo=.true. )
303 call hist_in( mslp(:,:),
'MSLP',
'mean sea-level pressure',
'Pa' )
305 call comm_vars8( sflx_mu(:,:), 1 )
306 call comm_vars8( sflx_mv(:,:), 2 )
307 call comm_wait ( sflx_mu(:,:), 1 )
308 call comm_wait ( sflx_mv(:,:), 2 )
313 momz_t_sf(i,j) = sflx_mw(i,j) * rfdz(
ks) / gsqrt(
ks,i,j,
i_xyw)
320 momx_t_sf(i,j) = 0.5_rp * ( sflx_mu(i,j) + sflx_mu(i+1,j) ) * rcdz(
ks) / gsqrt(
ks,i,j,
i_uyz)
327 momy_t_sf(i,j) = 0.5_rp * ( sflx_mv(i,j) + sflx_mv(i,j+1) ) * rcdz(
ks) / gsqrt(
ks,i,j,
i_xvz)
334 call thermodyn_qd( qdry, q )
335 call thermodyn_cp( cptot, q, qdry )
336 rtot = rdry * qdry + rvap * q(
i_qv)
337 rhot_t_sf(i,j) = ( sflx_sh(i,j) * rcdz(
ks) / ( cptot * gsqrt(
ks,i,j,
i_xyz) ) ) &
342 dens_t_sf(:,:) = 0.0_rp
346 work = sflx_qtrc(i,j,iq) * rcdz(
ks) / gsqrt(
ks,i,j,
i_xyz)
347 dens_t_sf(i,j) = dens_t_sf(i,j) + work
348 rhoq_t_sf(i,j,iq) = work
355 rhot_t_sf(i,j) = rhot_t_sf(i,j) + dens_t_sf(i,j) *
rhot(
ks,i,j) /
dens(
ks,i,j)
363 dens_t(
ks,i,j) = dens_t(
ks,i,j) + dens_t_sf(i,j)
364 momz_t(
ks,i,j) = momz_t(
ks,i,j) + momz_t_sf(i,j)
365 momx_t(
ks,i,j) = momx_t(
ks,i,j) + momx_t_sf(i,j)
366 momy_t(
ks,i,j) = momy_t(
ks,i,j) + momy_t_sf(i,j)
367 rhot_t(
ks,i,j) = rhot_t(
ks,i,j) + rhot_t_sf(i,j)
374 rhoq_t(
ks,i,j,iq) = rhoq_t(
ks,i,j,iq) + rhoq_t_sf(i,j,iq)
380 call stat_total( total, dens_t_sf(:,:),
'DENS_t_SF' )
381 call stat_total( total, momz_t_sf(:,:),
'MOMZ_t_SF' )
382 call stat_total( total, momx_t_sf(:,:),
'MOMX_t_SF' )
383 call stat_total( total, momy_t_sf(:,:),
'MOMY_t_SF' )
384 call stat_total( total, rhot_t_sf(:,:),
'RHOT_t_SF' )
387 call stat_total( total, rhoq_t_sf(:,:,iq), trim(
aq_name(iq))//
'_t_SF' )
real(rp), dimension(:,:,:), allocatable, public dens_tp
integer, public is
start point of inner domain: x, local
subroutine, public atmos_phy_sf_setup(SF_TYPE)
logical, public statistics_checktotal
calc&report variable totals to logfile?
integer, public je
end point of inner domain: y, local
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0e
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_momz_t
integer, public const_i_lw
long-wave radiation index
real(rp), dimension(:,:,:), allocatable, public atmos_phy_sf_sflx_qtrc
subroutine, public prc_mpistop
Abort MPI.
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_mw
real(rp), dimension(:,:,:), allocatable, public atmos_phy_sf_rhoq_t
real(rp), dimension(:,:,:), allocatable, target, public rhot
real(rp), dimension(:,:,:), allocatable, public momy_tp
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_v10
logical, public io_l
output log or not? (this process)
module ATMOSPHERE / Bottom boundary treatment
real(rp), dimension(:,:,:), allocatable, public pres
subroutine, public atmos_phy_sf_driver_resume
Resume.
module ATMOSPHERIC Variables
module ATMOSPHERE / Physics Surface fluxes
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_rhot_t
real(rp), dimension(:,:,:), allocatable, public rhot_tp
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_temp
real(rp), dimension(:,:,:), allocatable, public real_cz
geopotential height [m] (cell center)
real(rp), dimension(:,:,:), allocatable, target, public dens
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_mv
real(rp), dimension(:), allocatable, public grid_rcdz
reciprocal of center-dz
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
real(rp), dimension(:,:), allocatable, public atmos_phy_rd_sflx_lw_dn
module Atmosphere / Physics Radiation
character(len=h_short), public atmos_phy_sf_type
module ATMOSPHERIC Surface Variables
integer, public ia
of x whole cells (local, with HALO)
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_u10
real(rp), dimension(:,:,:), allocatable, public temp
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_momx_t
real(rp), dimension(:,:), allocatable, public real_z1
Height of the lowermost grid from surface (cell center) [m].
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_dens_t
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0h
real(rp), dimension(:,:), allocatable, public atmos_phy_rd_sflx_sw_dn
procedure(sf), pointer, public atmos_phy_sf
real(rp), public const_pre00
pressure reference [Pa]
real(rp), dimension(:,:,:,:), allocatable, public rhoq_tp
real(rp), dimension(:,:,:), allocatable, public atmos_phy_sf_sfc_albedo
character(len=h_short), dimension(:), allocatable, public aq_name
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_lh
subroutine, public atmos_phy_sf_driver_setup
Setup.
module ATMOSPHERE / Hydrostatic barance
integer, public js
start point of inner domain: y, local
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_pres
real(dp), public time_dtsec_atmos_phy_sf
time interval of physics(surface flux) [sec]
logical, public atmos_sw_phy_sf
module ATMOSPHERE / Physics Surface fluxes
real(rp), parameter, public const_rvap
specific gas constant (water vapor) [J/kg/K]
subroutine, public atmos_bottom_estimate(DENS, PRES, CZ, Zsfc, Z1, SFC_DENS, SFC_PRES)
Calc bottom boundary of atmosphere (just above surface)
integer, public ks
start point of inner domain: z, local
real(rp), dimension(:,:,:), allocatable, public momx_tp
subroutine, public prof_rapstart(rapname_base, level)
Start raptime.
real(rp), dimension(:,:,:), allocatable, public momz_tp
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_sh
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_momy_t
real(rp), dimension(:,:,:), allocatable, public v
integer, public ie
end point of inner domain: x, local
real(rp), dimension(:,:,:,:), allocatable, public gtrans_gsqrt
transformation metrics from Z to Xi, {G}^1/2
real(rp), dimension(:,:,:), allocatable, public w
module ATMOSPHERE / Thermodynamics
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_q2
integer, public const_i_sw
short-wave radiation index
real(rp), dimension(:,:), allocatable, public topo_zsfc
absolute ground height [m]
real(rp), dimension(:,:,:), allocatable, public u
real(rp), dimension(:), allocatable, public grid_rfdz
reciprocal of face-dz
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0m
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_t2
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_dens
integer, public io_fid_log
Log file ID.
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_gh
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_mu
subroutine, public prof_rapend(rapname_base, level)
Save raptime.
real(rp), dimension(:,:,:,:), allocatable, target, public qtrc
subroutine, public atmos_phy_sf_driver(update_flag)
Driver.
integer, public ja
of y whole cells (local, with HALO)