module ATMOSPHERE / Physics Cumulus
More...
module ATMOSPHERE / Physics Cumulus
- Description
- Cumulus parameterization driver
- Author
- Team SCALE
- NAMELIST
-
- History Output
name | description | unit | variable |
CBMFX | cloud base mass flux | kg/m2/s | MFLX_cloudbase |
CUBASE | CP cloud base height | m | cloudbase |
CUMFRC_DP | CP cloud fraction (deep) | 1 | cldfrac_dp |
CUMFRC_SH | CP cloud fraction (shallow) | 1 | cldfrac_sh |
CUMHGT | CP cloud top height | m | cloudtop |
DENS_t_CP | tendency DENS in CP | kg/m3/s | DENS_t_CP |
PREC_CP | surface precipitation rate by CP | kg/m2/s | SFLX_prec |
QV_t_CP | tendency rho*QV in CP | kg/m3/s | RHOQV_t_CP |
RAIN_CP | surface rain rate by CP | kg/m2/s | SFLX_rain |
RHOT_t_CP | tendency RHOT in CP | K*kg/m3/s | RHOT_t_CP |
SNOW_CP | surface snow rate by CP | kg/m2/s | SFLX_snow |
kf_nca | advection or cumulus convection timescale for KF | s | kf_nca |
{HYD_NAME}_t_CP | tendency rho*{HYD_NAME} in CP;
{HYD_NAME} is QC, QR, QI, QS, QG, QH. | kg/m3/s | RHOHYD_t_CP |
w0mean | running mean vertical wind velocity | kg/m2/s | w0mean |
◆ atmos_phy_cp_driver_setup()
subroutine, public mod_atmos_phy_cp_driver::atmos_phy_cp_driver_setup |
Setup.
Definition at line 49 of file mod_atmos_phy_cp_driver.F90.
72 log_info(
"ATMOS_PHY_CP_driver_setup",*)
'Setup'
85 log_error(
"ATMOS_PHY_CP_driver_setup",*)
'ATMOS_PHY_CP_TYPE (', trim(
atmos_phy_cp_type),
') is invalid. Check!'
90 log_info(
"ATMOS_PHY_CP_driver_setup",*)
'this component is never called.'
References scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_area, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_cz, scale_atmos_hydrometeor::atmos_hydrometeor_ice_phase, scale_atmos_phy_cp_common::atmos_phy_cp_common_setup(), scale_atmos_phy_cp_kf::atmos_phy_cp_kf_setup(), mod_atmos_admin::atmos_phy_cp_type, mod_atmos_admin::atmos_sw_phy_cp, scale_atmos_grid_cartesc_index::ia, scale_atmos_grid_cartesc_index::ie, scale_atmos_grid_cartesc_index::is, scale_atmos_grid_cartesc_index::ja, scale_atmos_grid_cartesc_index::je, scale_atmos_grid_cartesc_index::js, scale_atmos_grid_cartesc_index::ka, scale_atmos_grid_cartesc_index::ke, scale_atmos_grid_cartesc_index::ks, scale_prc::prc_abort(), scale_time::time_dtsec, and scale_time::time_dtsec_atmos_phy_cp.
Referenced by mod_atmos_driver::atmos_driver_setup().
◆ atmos_phy_cp_driver_calc_tendency()
subroutine, public mod_atmos_phy_cp_driver::atmos_phy_cp_driver_calc_tendency |
( |
logical, intent(in) |
update_flag | ) |
|
Driver.
Definition at line 99 of file mod_atmos_phy_cp_driver.F90.
163 sflx_snow => atmos_phy_cp_sflx_snow, &
164 sflx_engi => atmos_phy_cp_sflx_engi, &
165 cloudtop => atmos_phy_cp_cloudtop, &
166 cloudbase => atmos_phy_cp_cloudbase, &
167 cldfrac_dp => atmos_phy_cp_cldfrac_dp, &
168 cldfrac_sh => atmos_phy_cp_cldfrac_sh, &
169 w0mean => atmos_phy_cp_w0mean, &
170 kf_nca => atmos_phy_cp_kf_nca
173 logical,
intent(in) :: update_flag
175 real(RP) :: RHOQ_t_CP(KA,IA,JA,QS_MP:QE_MP)
177 real(RP) :: SFLX_prec(IA,JA)
179 integer :: k, i, j, iq
187 call file_history_in( w0mean(:,:,:),
'w0mean',
'running mean vertical wind velocity',
'kg/m2/s', fill_halo=.true. )
189 if ( update_flag )
then
194 u(:,:,:),
v(:,:,:), &
202 rhoqv_t_cp(:,:,:), rhohyd_t_cp(:,:,:,:), &
204 sflx_rain(:,:), sflx_snow(:,:), &
206 cloudtop(:,:), cloudbase(:,:), &
207 cldfrac_dp(:,:,:), cldfrac_sh(:,:,:) )
215 mflx_cloudbase(i,j) = 0.0_rp
224 sflx_prec(i,j) = sflx_rain(i,j) + sflx_snow(i,j)
228 call file_history_in( mflx_cloudbase(:,:),
'CBMFX',
'cloud base mass flux',
'kg/m2/s', fill_halo=.true. )
229 call file_history_in( sflx_rain(:,:),
'RAIN_CP',
'surface rain rate by CP',
'kg/m2/s', fill_halo=.true. )
230 call file_history_in( sflx_snow(:,:),
'SNOW_CP',
'surface snow rate by CP',
'kg/m2/s', fill_halo=.true. )
231 call file_history_in( sflx_prec(:,:),
'PREC_CP',
'surface precipitation rate by CP',
'kg/m2/s', fill_halo=.true. )
232 call file_history_in( cloudtop(:,:),
'CUMHGT',
'CP cloud top height',
'm', fill_halo=.true. )
233 call file_history_in( cloudbase(:,:),
'CUBASE',
'CP cloud base height',
'm', fill_halo=.true. )
234 call file_history_in( cldfrac_dp(:,:,:),
'CUMFRC_DP',
'CP cloud fraction (deep)',
'1', fill_halo=.true. )
235 call file_history_in( cldfrac_sh(:,:,:),
'CUMFRC_SH',
'CP cloud fraction (shallow)',
'1', fill_halo=.true. )
236 call file_history_in( kf_nca(:,:),
'kf_nca',
'advection or cumulus convection timescale for KF',
's', fill_halo=.true. )
238 call file_history_in( dens_t_cp(:,:,:),
'DENS_t_CP',
'tendency DENS in CP',
'kg/m3/s' , fill_halo=.true. )
239 call file_history_in( rhot_t_cp(:,:,:),
'RHOT_t_CP',
'tendency RHOT in CP',
'K*kg/m3/s', fill_halo=.true. )
241 call file_history_in( rhoqv_t_cp(:,:,:),
'QV_t_CP',
'tendency rho*QV in CP',
'kg/m3/s', fill_halo=.true. )
243 call file_history_in( rhohyd_t_cp(:,:,:,iq), trim(hyd_name(iq))//
'_t_CP', &
244 'tendency rho*'//trim(hyd_name(iq))//
' in CP',
'kg/m3/s', fill_halo=.true. )
253 dens_t(k,i,j) = dens_t(k,i,j) + dens_t_cp(k,i,j)
254 rhot_t(k,i,j) = rhot_t(k,i,j) + rhot_t_cp(k,i,j)
260 rhoqv_t_cp(:,:,:), rhohyd_t_cp(:,:,:,:), &
268 rhoq_t(k,i,j,iq) = rhoq_t(k,i,j,iq) + rhoq_t_cp(k,i,j,iq)
275 call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
276 dens_t_cp(:,:,:),
'DENS_t_CP', &
277 atmos_grid_cartesc_real_vol(:,:,:), &
278 atmos_grid_cartesc_real_totvol )
279 call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
280 rhot_t_cp(:,:,:),
'RHOT_t_CP', &
281 atmos_grid_cartesc_real_vol(:,:,:), &
282 atmos_grid_cartesc_real_totvol )
285 call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
286 rhoq_t_cp(:,:,:,iq), trim(
tracer_name(iq))//
'_t_CP', &
287 atmos_grid_cartesc_real_vol(:,:,:), &
288 atmos_grid_cartesc_real_totvol )
References scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_fz, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_totvol, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_totvolwxy, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_totvolzuy, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_totvolzxv, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_vol, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_volwxy, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_volzuy, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_volzxv, scale_atmos_phy_cp_common::atmos_phy_cp_common_wmean(), mod_atmos_phy_cp_vars::atmos_phy_cp_dens_t, scale_atmos_phy_cp_kf::atmos_phy_cp_kf_tendency(), mod_atmos_phy_cp_vars::atmos_phy_cp_mflx_cloudbase, mod_atmos_phy_cp_vars::atmos_phy_cp_rhohyd_t, mod_atmos_phy_cp_vars::atmos_phy_cp_rhoqv_t, mod_atmos_phy_cp_vars::atmos_phy_cp_rhot_t, mod_atmos_phy_cp_vars::atmos_phy_cp_sflx_rain, mod_atmos_admin::atmos_phy_cp_type, mod_atmos_phy_mp_driver::atmos_phy_mp_driver_qhyd2qtrc(), scale_const::const_pre00, scale_atmos_hydrometeor::cv_water, mod_atmos_vars::dens, mod_atmos_vars::dens_av, mod_atmos_vars::dens_tp, scale_atmos_hydrometeor::hyd_name, scale_atmos_grid_cartesc_index::ia, scale_atmos_grid_cartesc_index::ie, scale_atmos_grid_cartesc_index::ieb, scale_atmos_grid_cartesc_index::is, scale_atmos_grid_cartesc_index::isb, scale_atmos_grid_cartesc_index::ja, scale_atmos_grid_cartesc_index::je, scale_atmos_grid_cartesc_index::jeb, scale_atmos_grid_cartesc_index::js, scale_atmos_grid_cartesc_index::jsb, scale_tracer::k, scale_atmos_grid_cartesc_index::ka, scale_atmos_grid_cartesc_index::ke, scale_atmos_grid_cartesc_index::ks, mod_atmos_vars::momx, mod_atmos_vars::momx_av, mod_atmos_vars::momx_tp, mod_atmos_vars::momy, mod_atmos_vars::momy_av, mod_atmos_vars::momy_tp, mod_atmos_vars::momz, mod_atmos_vars::momz_av, mod_atmos_vars::momz_tp, scale_atmos_hydrometeor::n_hyd, mod_atmos_vars::pres, mod_atmos_vars::qdry, mod_atmos_phy_mp_vars::qe_mp, mod_atmos_phy_mp_vars::qs_mp, mod_atmos_vars::qtrc, mod_atmos_vars::qtrc_av, mod_atmos_vars::qv, mod_atmos_vars::rhoq_tp, mod_atmos_vars::rhot, mod_atmos_vars::rhot_av, mod_atmos_vars::rhot_tp, scale_statistics::statistics_checktotal, mod_atmos_vars::temp, scale_time::time_dtsec, scale_time::time_dtsec_atmos_phy_cp, scale_tracer::tracer_name, mod_atmos_vars::u, mod_atmos_vars::v, and mod_atmos_vars::w.
Referenced by mod_atmos_driver::atmos_driver_calc_tendency().
real(rp), dimension(:,:,:), pointer, public momz_av
real(rp), dimension(:,:), allocatable, public atmos_phy_cp_mflx_cloudbase
module Atmosphere / Physics Cloud Microphysics
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
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)
character(len=h_short), dimension(n_hyd), parameter, public hyd_name
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, target, public qdry
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
real(rp), dimension(:,:,:), allocatable, public atmos_phy_cp_rhoqv_t
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(dp), public time_dtsec_atmos_phy_cp
time interval of physics(cumulus ) [sec]
real(rp), dimension(:,:,:), allocatable, target, public dens
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_volzxv
control volume (zxv) [m3]
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, public atmos_phy_cp_rhot_t
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
logical, public atmos_sw_phy_cp
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
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
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.
logical, public statistics_checktotal
calc&report variable totals to logfile?
real(rp), dimension(:,:,:,:), allocatable, public atmos_phy_cp_rhohyd_t
real(dp), public time_dtsec
time interval of model [sec]
module ATMOSPHERIC Variables
real(rp), dimension(:,:,:), allocatable, public rhot_tp
real(rp), dimension(:,:), allocatable, public atmos_phy_cp_sflx_rain
real(rp), dimension(:,:,:), allocatable, public momy_tp
real(rp), dimension(:,:,:), allocatable, public momx_tp
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]
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
real(rp), public const_pre00
pressure reference [Pa]
integer, parameter, public n_hyd
real(rp), dimension(:,:,:), allocatable, public atmos_phy_cp_dens_t
real(rp), public cv_water
CV for water [J/kg/K].