Go to the documentation of this file.
76 log_info(
"ATMOS_PHY_CP_driver_setup",*)
'Setup'
88 case (
'KF-JMAPPLIB' )
92 log_error(
"ATMOS_PHY_CP_driver_setup",*)
'ATMOS_PHY_CP_TYPE (', trim(
atmos_phy_cp_type),
') is invalid. Check!'
97 log_info(
"ATMOS_PHY_CP_driver_setup",*)
'this component is never called.'
115 log_info(
"ATMOS_PHY_CP_driver_finalize",*)
'Finalize'
203 sflx_snow => atmos_phy_cp_sflx_snow, &
204 sflx_engi => atmos_phy_cp_sflx_engi, &
205 cloudtop => atmos_phy_cp_cloudtop, &
206 cloudbase => atmos_phy_cp_cloudbase, &
207 cldfrac_dp => atmos_phy_cp_cldfrac_dp, &
208 cldfrac_sh => atmos_phy_cp_cldfrac_sh, &
209 w0mean => atmos_phy_cp_w0mean, &
210 kf_nca => atmos_phy_cp_kf_nca
218 logical,
intent(in) :: update_flag
222 real(
rp) :: sflx_prec(
ia,
ja)
224 integer ::
k, i, j, iq
239 if ( update_flag )
then
248 u(:,:,:),
v(:,:,:), &
256 rhoqv_t_cp(:,:,:), rhohyd_t_cp(:,:,:,:), &
258 sflx_rain(:,:), sflx_snow(:,:), &
260 cloudtop(:,:), cloudbase(:,:), &
261 cldfrac_dp(:,:,:), cldfrac_sh(:,:,:) )
268 sflx_prec(i,j) = sflx_rain(i,j) + sflx_snow(i,j)
273 case (
'KF-JMAPPLIB' )
279 sflx_prec(i,j) = sflx_rain(i,j) + sflx_snow(i,j)
287 u(:,:,:),
v(:,:,:),
w(:,:,:), &
291 qc(:,:,:),
qi(:,:,:), &
292 us(:,:), pblh(:,:), sflx_buoy(:,:), &
293 cz(:,:,:), fz(:,:,:), &
296 rhoqv_t_cp(:,:,:), rhohyd_t_cp(:,:,:,:), &
297 w0mean(:,:,:), kf_nca(:,:), &
298 sflx_rain(:,:), sflx_snow(:,:), &
300 cloudtop(:,:), cloudbase(:,:) )
307 sflx_engi(i,j) = sflx_rain(i,j) * cv_water *
temp(
ks,i,j) &
308 + sflx_snow(i,j) * ( cv_ice *
temp(
ks,i,j) - lhf )
320 mflx_cloudbase(i,j) = 0.0_rp
327 call file_history_in( w0mean(:,:,:),
'w0mean',
'running mean vertical wind velocity',
'kg/m2/s', fill_halo=.true. )
328 call file_history_in( mflx_cloudbase(:,:),
'CBMFX',
'cloud base mass flux',
'kg/m2/s', fill_halo=.true. )
329 call file_history_in( sflx_rain(:,:),
'RAIN_CP',
'surface rain rate by CP',
'kg/m2/s', fill_halo=.true. )
330 call file_history_in( sflx_snow(:,:),
'SNOW_CP',
'surface snow rate by CP',
'kg/m2/s', fill_halo=.true. )
331 call file_history_in( sflx_prec(:,:),
'PREC_CP',
'surface precipitation rate by CP',
'kg/m2/s', fill_halo=.true. )
332 call file_history_in( cloudtop(:,:),
'CUMHGT',
'CP cloud top height',
'm', fill_halo=.true. )
333 call file_history_in( cloudbase(:,:),
'CUBASE',
'CP cloud base height',
'm', fill_halo=.true. )
334 call file_history_in( cldfrac_dp(:,:,:),
'CUMFRC_DP',
'CP cloud fraction (deep)',
'1', fill_halo=.true. )
335 call file_history_in( cldfrac_sh(:,:,:),
'CUMFRC_SH',
'CP cloud fraction (shallow)',
'1', fill_halo=.true. )
336 call file_history_in( kf_nca(:,:),
'kf_nca',
'advection or cumulus convection timescale for KF',
's', fill_halo=.true. )
338 call file_history_in( dens_t_cp(:,:,:),
'DENS_t_CP',
'tendency DENS in CP',
'kg/m3/s' , fill_halo=.true. )
339 call file_history_in( rhot_t_cp(:,:,:),
'RHOT_t_CP',
'tendency RHOT in CP',
'K*kg/m3/s', fill_halo=.true. )
341 call file_history_in( rhoqv_t_cp(:,:,:),
'QV_t_CP',
'tendency rho*QV in CP',
'kg/m3/s', fill_halo=.true. )
343 call file_history_in( rhohyd_t_cp(:,:,:,iq), trim(hyd_name(iq))//
'_t_CP', &
344 'tendency rho*'//trim(hyd_name(iq))//
' in CP',
'kg/m3/s', fill_halo=.true. )
356 dens_t(
k,i,j) = dens_t(
k,i,j) + dens_t_cp(
k,i,j)
357 rhot_t(
k,i,j) = rhot_t(
k,i,j) + rhot_t_cp(
k,i,j)
366 rhoqv_t_cp(:,:,:), rhohyd_t_cp(:,:,:,:), &
375 rhoq_t(
k,i,j,iq) = rhoq_t(
k,i,j,iq) + rhoq_t_cp(
k,i,j,iq)
384 dens_t_cp(:,:,:),
'DENS_t_CP', &
385 atmos_grid_cartesc_real_vol(:,:,:), &
386 atmos_grid_cartesc_real_totvol )
388 rhot_t_cp(:,:,:),
'RHOT_t_CP', &
389 atmos_grid_cartesc_real_vol(:,:,:), &
390 atmos_grid_cartesc_real_totvol )
394 rhoq_t_cp(:,:,:,iq), trim(
tracer_name(iq))//
'_t_CP', &
395 atmos_grid_cartesc_real_vol(:,:,:), &
396 atmos_grid_cartesc_real_totvol )
real(rp), dimension(:,:,:), pointer, public momz_av
real(rp), dimension(:,:), allocatable, public atmos_phy_cp_mflx_cloudbase
module Atmosphere / Physics Cloud Microphysics
real(rp), dimension(:,:,:), allocatable, public atmos_phy_cp_rhoqv_t
integer, public ke
end point of inner domain: z, local
real(rp), dimension(:,:,:,:), allocatable, public atmos_phy_cp_rhohyd_t
subroutine, public prc_abort
Abort Process.
subroutine, public atmos_phy_cp_kf_setup(KA, KS, KE, IA, IS, IE, JA, JS, JE, CZ, AREA, WARMRAIN_in)
Setup initial setup for Kain-Fritsch Cumulus Parameterization.
real(rp), dimension(:,:,:), pointer, public momx_av
module atmosphere / physics / cumulus / kf-jmapplib
real(rp), public atmos_grid_cartesc_real_totvolzxv
total volume (zxv, local) [m3]
real(rp), dimension(:,:,:,:), allocatable, public rhoq_tp
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_cz
geopotential height [m] (zxy)
real(rp), dimension(:,:,:), allocatable, target, public pott
character(len=h_short), dimension(n_hyd), parameter, public hyd_name
real(rp), dimension(:,:,:), allocatable, public atmos_phy_cp_rhot_t
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_volwxy
control volume (wxy) [m3]
real(rp), dimension(:,:,:,:), pointer, public qtrc_av
module atmosphere / physics / cumulus / Common
real(rp), dimension(:,:), allocatable, public atmos_phy_bl_zi
real(rp), dimension(:,:,:), allocatable, target, public qdry
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_ustar
real(rp), dimension(:,:), allocatable, public atmos_phy_bl_sflx_buoy
subroutine, public atmos_phy_cp_kf_tendency(KA, KS, KE, IA, IS, IE, JA, JS, JE, DENS, U, V, RHOT, TEMP, PRES, QDRY, QV_in, w0avg, FZ, KF_DTSECD, DENS_t, RHOT_t, RHOQV_t, RHOQ_t, nca, SFLX_rain, SFLX_snow, SFLX_engi, cloudtop, cloudbase, cldfrac_dp, cldfrac_sh)
ATMOS_PHY_CP_kf calculate Kain-Fritsch Cumulus Parameterization.
real(rp), public atmos_grid_cartesc_real_totvolwxy
total volume (wxy, local) [m3]
real(rp), dimension(:,:,:), pointer, public rhot_av
module atmosphere / hydrometeor
module ATMOSPHERIC Surface Variables
module atmosphere / physics / cumulus / Kain-Fritsch
real(rp), dimension(:,:,:), allocatable, target, public rhot
module Atmosphere GRID CartesC Real(real space)
real(rp), dimension(:,:,:,:), allocatable, target, public qtrc
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_volzuy
control volume (zuy) [m3]
module atmosphere / physics / cloud microphysics
real(rp), dimension(:,:), allocatable, public atmos_phy_cp_sflx_rain
real(dp), public time_dtsec_atmos_phy_cp
time interval of physics(cumulus ) [sec]
integer, parameter, public rp
integer, public ie
end point of inner domain: x, local
subroutine, public atmos_phy_cp_kf_jmapplib_tendency(KA, KS, KE, IA, IS, IE, JA, JS, JE, DENS, U, V, W, TEMP, POTT, PRES, EXNER, QDRY, QV, QC, QI, us, PBLH, SFLX_BUOY, CZ, FZ, DENS_t, RHOT_t, RHOQV_t, RHOQ_t, w0avg, nca, SFLX_rain, SFLX_snow, SFLX_prec, cloudtop, cloudbase)
ATMOS_PHY_CP_KF_JMAPPLIB_tendency calculate tendency by the virtical eddy viscosity.
module atmosphere / physics / PBL
real(rp), dimension(:,:,:), allocatable, target, public dens
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_volzxv
control volume (zxv) [m3]
module atmosphere / grid / cartesC index
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_vol
control volume (zxy) [m3]
real(rp), dimension(:,:,:), allocatable, target, public momz
real(rp), dimension(:,:,:), pointer, public momy_av
real(rp), dimension(:,:,:), allocatable, target, public v
character(len=h_short), dimension(qa_max), public tracer_name
real(rp), dimension(:,:,:), allocatable, target, public w
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_area
horizontal area ( xy, normal z) [m2]
real(rp), dimension(:,:,:), allocatable, public momz_tp
real(rp), dimension(:,:,:), allocatable, target, public momx
real(rp), dimension(:,:,:), allocatable, target, public exner
integer, public is
start point of inner domain: x, local
logical, public atmos_sw_phy_cp
subroutine, public atmos_phy_cp_kf_finalize
finalize
real(rp), dimension(:,:,:), allocatable, target, public temp
real(rp), dimension(:,:,:), allocatable, public dens_tp
real(rp), public atmos_grid_cartesc_real_totvolzuy
total volume (zuy, local) [m3]
real(rp), dimension(:,:,:), allocatable, target, public momy
real(rp), dimension(:,:,:), allocatable, pointer, target, public qv
module ATMOSPHERE / Physics Cumulus
character(len=h_short), public atmos_phy_cp_type
real(rp), dimension(:,:,:), allocatable, target, public pres
real(rp), dimension(:,:,:), pointer, public dens_av
real(rp), dimension(:,:,:), allocatable, target, public u
integer, public ks
start point of inner domain: z, local
logical, public atmos_hydrometeor_ice_phase
subroutine, public atmos_phy_cp_common_wmean(KA, KS, KE, IA, IS, IE, JA, JS, JE, W, TIME_DTSEC, CP_DTSEC, W0_mean)
ATMOS_PHY_CP_wmean running mean vertical wind velocity comment for W0 imported from WRF.
subroutine, public atmos_phy_cp_driver_setup
Setup.
logical, public statistics_checktotal
calc&report variable totals to logfile?
real(dp), public time_dtsec
time interval of model [sec]
module ATMOSPHERIC Variables
real(rp), dimension(:,:,:), pointer, public qi
real(rp), dimension(:,:,:), allocatable, public rhot_tp
integer, public js
start point of inner domain: y, local
subroutine, public atmos_phy_cp_driver_finalize
finalize
real(rp), dimension(:,:,:), allocatable, public atmos_phy_cp_dens_t
real(rp), dimension(:,:,:), allocatable, public momy_tp
subroutine, public atmos_phy_cp_driver_calc_tendency(update_flag)
Driver.
real(rp), public lhf
latent heat of fusion for use [J/kg]
real(rp), dimension(:,:,:), allocatable, public momx_tp
subroutine, public atmos_phy_cp_kf_jmapplib_setup(KA, KS, KE, IA, JA, dx, dt)
ATMOS_PHY_CP_KF_JMAPPLIB_setup Setup.
real(rp), dimension(:,:,:), pointer, public qc
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_fz
geopotential height [m] (wxy)
real(rp), public atmos_grid_cartesc_real_totvol
total volume (zxy, local) [m3]
module atmosphere / grid / cartesC
subroutine, public atmos_phy_cp_common_setup
Setup.
subroutine, public atmos_phy_mp_driver_qhyd2qtrc(KA, KS, KE, IA, IS, IE, JA, JS, JE, QV, QHYD, QTRC, QNUM)
module Atmosphere / Physics Cumulus
integer, public je
end point of inner domain: y, local
integer, parameter, public n_hyd
real(rp), public cv_water
CV for water [J/kg/K].
real(rp), public cv_ice
CV for ice [J/kg/K].