module atmosphere / physics / PBL
More...
module atmosphere / physics / PBL
- Description
- Planetary boundary layer turbulence
- Author
- Team SCALE
- NAMELIST
-
- History Output
name | description | unit | variable |
Kh_BL | eddy diffusion | m2/s | Kh |
Nu_BL | eddy viscosity | m2/s | Nu |
QL_BL | liquid water content in partial condensation | kg/kg | QL |
RHOT_t_BL | RHOT tendency (BL) | K.kg/m3/s | RHOT_t_BL |
RHOU_t_BL | MOMX tendency (BL) | kg/m2/s2 | RHOU_t_BL |
RHOV_t_BL | MOMY tendency (BL) | kg/m2/s2 | RHOV_t_BL |
Zi_BL | depth of the boundary layer | m | Zi |
cldfrac_BL | cloud fraction in partial condensation | 1 | cldfrac |
{TRACER_NAME}_t_BL | RHO*{TRACER_NAME} tendency (BL);
{TRACER_NAME} depends on the physics schemes, e.g., QV, QC, QR. | kg/m3/s | RHOQ_t_BL |
◆ atmos_phy_bl_driver_tracer_setup()
subroutine, public mod_atmos_phy_bl_driver::atmos_phy_bl_driver_tracer_setup |
Config.
Definition at line 52 of file mod_atmos_phy_bl_driver.F90.
76 log_info(
"ATMOS_PHY_BL_driver_tracer_setup",*)
'Setup'
89 case (
'MYNN-JMAPPLIB' )
98 log_error(
"ATMOS_PHY_BL_driver_tracer_setup",*)
'ATMOS_PHY_BL_TYPE is invalid: ', trim(
atmos_phy_bl_type)
References scale_atmos_phy_bl_mynn::atmos_phy_bl_mynn_desc, scale_atmos_phy_bl_mynn_jmapplib::atmos_phy_bl_mynn_jmapplib_desc, scale_atmos_phy_bl_mynn_jmapplib::atmos_phy_bl_mynn_jmapplib_name, scale_atmos_phy_bl_mynn_jmapplib::atmos_phy_bl_mynn_jmapplib_ntracer, scale_atmos_phy_bl_mynn_jmapplib::atmos_phy_bl_mynn_jmapplib_units, scale_atmos_phy_bl_mynn::atmos_phy_bl_mynn_name, scale_atmos_phy_bl_mynn::atmos_phy_bl_mynn_ntracer, scale_atmos_phy_bl_mynn::atmos_phy_bl_mynn_tracer_setup(), scale_atmos_phy_bl_mynn::atmos_phy_bl_mynn_units, mod_atmos_admin::atmos_phy_bl_type, mod_atmos_admin::atmos_sw_phy_bl, scale_prc::prc_abort(), mod_atmos_phy_bl_vars::qe, mod_atmos_phy_bl_vars::qs, and scale_tracer::tracer_regist().
Referenced by mod_atmos_driver::atmos_driver_tracer_setup().
◆ atmos_phy_bl_driver_setup()
subroutine, public mod_atmos_phy_bl_driver::atmos_phy_bl_driver_setup |
Setup.
Definition at line 109 of file mod_atmos_phy_bl_driver.F90.
130 log_info(
"ATMOS_PHY_BL_driver_setup",*)
'Setup'
138 case (
'MYNN-JMAPPLIB' )
145 log_info(
"ATMOS_PHY_BL_driver_setup",*)
'this component is never called.'
References scale_atmos_grid_cartesc::atmos_grid_cartesc_cz, scale_atmos_phy_bl_mynn_jmapplib::atmos_phy_bl_mynn_jmapplib_setup(), scale_atmos_phy_bl_mynn::atmos_phy_bl_mynn_setup(), mod_atmos_phy_bl_vars::atmos_phy_bl_sflx_buoy, mod_atmos_admin::atmos_phy_bl_type, mod_atmos_phy_bl_vars::atmos_phy_bl_zi, mod_atmos_admin::atmos_sw_phy_bl, scale_bulkflux::bulkflux_type, scale_atmos_grid_cartesc::dx, scale_atmos_grid_cartesc_index::ka, scale_atmos_grid_cartesc_index::ke, scale_atmos_grid_cartesc_index::ks, and scale_time::time_dtsec_atmos_phy_bl.
Referenced by mod_atmos_driver::atmos_driver_setup().
◆ atmos_phy_bl_driver_finalize()
subroutine, public mod_atmos_phy_bl_driver::atmos_phy_bl_driver_finalize |
◆ atmos_phy_bl_driver_mkinit()
subroutine, public mod_atmos_phy_bl_driver::atmos_phy_bl_driver_mkinit |
( |
real(rp), intent(in) |
TKE_CONST | ) |
|
make initial state
Definition at line 191 of file mod_atmos_phy_bl_driver.F90.
263 atmos_vars_get_diagnostic
265 real(RP),
intent(in) :: TKE_CONST
267 real(RP) :: N2 (KA,IA,JA)
268 real(RP) :: POTL(KA,IA,JA)
269 real(RP) :: POTV(KA,IA,JA)
271 real(RP) :: ATM_TEMP(IA,JA)
272 real(RP) :: ATM_PRES(IA,JA)
273 real(RP) :: ATM_U (IA,JA)
274 real(RP) :: ATM_V (IA,JA)
275 real(RP) :: ATM_W (IA,JA)
276 real(RP) :: ATM_QV (IA,JA)
278 real(RP) :: QW(KA,IA,JA)
280 integer :: k, i, j, iq
287 if (
qtrc(ks,is,je,
qs) == undef )
then
291 call comm_vars8(
dens, 1)
292 call comm_vars8(
momz, 2)
293 call comm_vars8(
momx, 3)
294 call comm_vars8(
momy, 4)
295 call comm_vars8(
rhot, 5)
296 call comm_vars8(
qv, 6)
297 call comm_vars8(
qc, 7)
298 call comm_vars8(
qi, 8)
299 call comm_wait(
dens, 1)
300 call comm_wait(
momz, 2)
301 call comm_wait(
momx, 3)
302 call comm_wait(
momy, 4)
303 call comm_wait(
rhot, 5)
304 call comm_wait(
qv, 6)
305 call comm_wait(
qc, 7)
306 call comm_wait(
qi, 8)
309 call atmos_vars_get_diagnostic(
"N2", n2 )
310 call atmos_vars_get_diagnostic(
"POTL", potl )
311 call atmos_vars_get_diagnostic(
"POTV", potv )
313 call bottom_estimate( ka, ks, ke, ia, isb, ieb, ja, jsb, jeb, &
317 sfc_dens(:,:), sfc_pres(:,:) )
323 atm_temp(i,j) =
temp(ks,i,j)
324 atm_pres(i,j) =
pres(ks,i,j)
325 atm_u(i,j) =
u(ks,i,j)
326 atm_v(i,j) =
v(ks,i,j)
327 atm_w(i,j) = atm_u(i,j) * tansl_x(i,j) + atm_v(i,j) * tansl_y(i,j)
328 atm_qv(i,j) =
qv(ks,i,j)
334 atm_w(:,:), atm_u(:,:), atm_v(:,:), &
335 atm_temp(:,:), atm_pres(:,:), atm_qv(:,:), &
336 sfc_dens(:,:), sfc_temp(:,:), sfc_pres(:,:), &
337 sfc_z0m(:,:), sfc_z0h(:,:), sfc_z0e(:,:), &
339 sflx_mw(:,:), sflx_mu(:,:), sflx_mv(:,:), &
340 sflx_sh(:,:), sflx_lh(:,:), sflx_qv(:,:), &
341 ustar(:,:), tstar(:,:), qstar(:,:), &
344 u10(:,:), v10(:,:), t2(:,:), q2(:,:) )
350 qw(k,i,j) =
qv(k,i,j) +
qc(k,i,j) +
qi(k,i,j)
357 ka, ks, ke, ia, isb, ieb, ja, jsb, jeb, &
359 dens(:,:,:),
u(:,:,:),
v(:,:,:),
w(:,:,:),
pott(:,:,:), &
361 qdry(:,:,:),
qv(:,:,:), qw(:,:,:), &
362 potl(:,:,:), potv(:,:,:), &
364 sflx_mu(:,:), sflx_mv(:,:), sflx_sh(:,:), sflx_qv(:,:), &
365 ustar(:,:), tstar(:,:), qstar(:,:), rlmo(:,:), &
367 cz(:,:,:), fz(:,:,:), f2h(:,:,:,:), &
374 case (
'MYNN-JMAPPLIB' )
380 if (
qtrc(k,i,j,
qs) == undef )
then
381 qtrc(k,i,j,
qs) = tke_const
393 if (
qtrc(k,i,j,iq) == undef )
then
394 qtrc(k,i,j,iq) = 0.0_rp
References scale_atmos_bottom::atmos_bottom_estimate(), scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_cz, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_f2h, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_fz, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_z1, scale_atmos_phy_bl_mynn::atmos_phy_bl_mynn_mkinit(), mod_atmos_admin::atmos_phy_bl_type, mod_atmos_phy_bl_vars::atmos_phy_bl_zi, scale_atmos_phy_sf_bulk::atmos_phy_sf_bulk_flux(), mod_atmos_phy_sf_vars::atmos_phy_sf_q2, mod_atmos_phy_sf_vars::atmos_phy_sf_qstar, mod_atmos_phy_sf_vars::atmos_phy_sf_rlmo, mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_dens, mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_pres, mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_temp, mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_z0e, mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_z0h, mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_z0m, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_lh, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_mu, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_mv, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_qv, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_sh, mod_atmos_phy_sf_vars::atmos_phy_sf_t2, mod_atmos_phy_sf_vars::atmos_phy_sf_tstar, mod_atmos_phy_sf_vars::atmos_phy_sf_u10, mod_atmos_phy_sf_vars::atmos_phy_sf_ustar, mod_atmos_phy_sf_vars::atmos_phy_sf_v10, mod_atmos_phy_sf_vars::atmos_phy_sf_wstar, mod_atmos_admin::atmos_sw_phy_bl, mod_atmos_vars::atmos_vars_calc_diagnostics(), scale_bulkflux::bulkflux_type, scale_const::const_undef, mod_atmos_vars::dens, mod_atmos_vars::exner, scale_atmos_grid_cartesc_index::ia, 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::jsb, scale_tracer::k, scale_atmos_grid_cartesc_index::ka, scale_atmos_grid_cartesc_index::ke, scale_atmos_grid_cartesc_index::ks, scale_landuse::landuse_frac_land, mod_atmos_vars::momx, mod_atmos_vars::momy, mod_atmos_vars::momz, mod_atmos_vars::pott, mod_atmos_vars::pres, mod_atmos_vars::qc, mod_atmos_vars::qdry, mod_atmos_phy_bl_vars::qe, mod_atmos_vars::qi, mod_atmos_phy_bl_vars::qs, mod_atmos_vars::qtrc, mod_atmos_vars::qv, mod_atmos_vars::rhot, mod_atmos_vars::temp, scale_time::time_dtsec_atmos_phy_bl, scale_topography::topography_tansl_x, scale_topography::topography_tansl_y, mod_atmos_vars::u, mod_atmos_vars::v, and mod_atmos_vars::w.
Referenced by mod_mkinit::tke_setup().
◆ atmos_phy_bl_driver_calc_tendency()
subroutine, public mod_atmos_phy_bl_driver::atmos_phy_bl_driver_calc_tendency |
( |
logical, intent(in) |
update_flag | ) |
|
calculate tendency
Definition at line 409 of file mod_atmos_phy_bl_driver.F90.
453 atmos_vars_get_diagnostic
479 logical,
intent(in) :: update_flag
481 real(RP) :: Nu(KA,IA,JA)
482 real(RP) :: Kh(KA,IA,JA)
484 real(RP) :: QW(KA,IA,JA)
486 real(RP) :: N2 (KA,IA,JA)
487 real(RP) :: POTL(KA,IA,JA)
488 real(RP) :: POTV(KA,IA,JA)
490 real(RP),
pointer :: RHOQV_t(:,:,:)
492 integer :: k, i, j, iq
495 if ( update_flag )
then
500 rhoq_t_bl(:,:,:,:) = 0.0_rp
505 call atmos_vars_get_diagnostic(
"N2", n2 )
506 call atmos_vars_get_diagnostic(
"POTL", potl )
507 call atmos_vars_get_diagnostic(
"POTV", potv )
512 qw(k,i,j) =
qv(k,i,j) +
qc(k,i,j) +
qi(k,i,j)
518 rhoqv_t => rhoq_t_bl(:,:,:,
i_qv)
520 allocate( rhoqv_t(ka,ia,ja) )
524 ka, ks, ke, ia, is, ie, ja, js, je, &
525 dens(:,:,:),
u(:,:,:),
v(:,:,:),
w(:,:,:), &
528 qdry(:,:,:),
qv(:,:,:), qw(:,:,:), &
529 potl(:,:,:), potv(:,:,:), &
531 sflx_mu(:,:), sflx_mv(:,:), sflx_sh(:,:), sflx_qv(:,:), &
532 ustar(:,:), tstar(:,:), qstar(:,:), rlmo(:,:), &
534 cz(:,:,:), fz(:,:,:), f2h(:,:,:,:), dt_bl, &
536 rhou_t_bl(:,:,:), rhov_t_bl(:,:,:), rhot_t_bl(:,:,:), &
537 rhoqv_t(:,:,:), rhoq_t_bl(:,:,:,
qs:
qe), &
538 nu(:,:,:), kh(:,:,:), &
539 ql(:,:,:), cldfrac(:,:,:), zi(:,:), sflx_buoy(:,:) )
540 if (
i_qv <= 0 )
then
542 deallocate( rhoqv_t )
544 case (
'MYNN-JMAPPLIB' )
546 rhoqv_t => rhoq_t_bl(:,:,:,
i_qv)
548 allocate( rhoqv_t(ka,ia,ja) )
553 ka, ks, ke, ia, is, ie, ja, js, je, &
554 dens(:,:,:),
u(:,:,:),
v(:,:,:), &
557 qdry(:,:,:),
qv(:,:,:),
qc(:,:,:),
qi(:,:,:), &
558 sfc_dens(:,:), sfc_pres(:,:), &
559 sflx_mu(:,:), sflx_mv(:,:), sflx_sh(:,:), sflx_qv(:,:), &
560 ustar(:,:), rlmo(:,:), &
561 cz(:,:,:), fz(:,:,:), f2h(:,:,:,:), dt_bl, &
562 rhou_t_bl(:,:,:), rhov_t_bl(:,:,:), rhot_t_bl(:,:,:), &
563 rhoqv_t(:,:,:), rhoq_t_bl(:,:,:,
qs:
qe), &
564 nu(:,:,:), kh(:,:,:), &
565 zi = zi(:,:), sflx_buoy = sflx_buoy(:,:) )
567 if (
i_qv <= 0 )
then
569 deallocate( rhoqv_t )
575 if ( ( .not. tracer_advc(iq) ) .or. iq==
i_qv .or. (iq>=
qs .and. iq<=
qe) ) cycle
577 ka, ks, ke, ia, is, ie, ja, js, je, &
582 cz(:,:,:), fz(:,:,:), f2h(:,:,:,:), &
583 dt_bl, tracer_name(iq), &
584 rhoq_t_bl(:,:,:,iq) )
588 call file_history_in( nu(:,:,:),
'Nu_BL',
'eddy viscosity',
'm2/s', fill_halo=.true., dim_type=
"ZHXY" )
589 call file_history_in( kh(:,:,:),
'Kh_BL',
'eddy diffusion',
'm2/s', fill_halo=.true., dim_type=
"ZHXY" )
591 call file_history_in( ql(:,:,:),
'QL_BL',
'liquid water content in partial condensation',
'kg/kg', fill_halo=.true. )
592 call file_history_in( cldfrac(:,:,:),
'cldfrac_BL',
'cloud fraction in partial condensation',
'1', fill_halo=.true. )
594 call file_history_in( zi(:,:),
'Zi_BL',
'depth of the boundary layer',
'm', fill_halo=.true. )
596 call file_history_in( rhou_t_bl(:,:,:),
'RHOU_t_BL',
'MOMX tendency (BL)',
'kg/m2/s2', fill_halo=.true. )
597 call file_history_in( rhov_t_bl(:,:,:),
'RHOV_t_BL',
'MOMY tendency (BL)',
'kg/m2/s2', fill_halo=.true. )
598 call file_history_in( rhot_t_bl(:,:,:),
'RHOT_t_BL',
'RHOT tendency (BL)',
'K.kg/m3/s', fill_halo=.true. )
601 if ( .not. tracer_advc(iq) ) cycle
602 call file_history_in( rhoq_t_bl(:,:,:,iq), trim(tracer_name(iq))//
'_t_BL', &
603 'RHO*'//trim(tracer_name(iq))//
' tendency (BL)',
'kg/m3/s', fill_halo=.true. )
607 call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
608 rhou_t_bl(:,:,:),
'RHOU_t_BL', &
611 call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
612 rhov_t_bl(:,:,:),
'RHOV_t_BL', &
615 call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
616 rhot_t_bl(:,:,:),
'RHOT_t_BL', &
619 call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
620 nu(:,:,:),
'Nu_BL', &
623 call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
624 kh(:,:,:),
'Kh_BL', &
629 if ( .not. tracer_advc(iq) ) cycle
630 call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
631 rhoq_t_bl(:,:,:,iq), trim(tracer_name(iq))//
'_t_BL', &
647 rhou_t(k,i,j) = rhou_t(k,i,j) + rhou_t_bl(k,i,j)
648 rhov_t(k,i,j) = rhov_t(k,i,j) + rhov_t_bl(k,i,j)
649 rhot_t(k,i,j) = rhot_t(k,i,j) + rhot_t_bl(k,i,j)
659 if ( .not. tracer_advc(iq) ) cycle
666 if ( tracer_advc(iq) ) &
668 rhoq_t(k,i,j,iq) = rhoq_t(k,i,j,iq) + rhoq_t_bl(k,i,j,iq)
References scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_cz, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_f2h, 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_vol, mod_atmos_phy_bl_vars::atmos_phy_bl_cldfrac, mod_atmos_phy_bl_vars::atmos_phy_bl_mix_tracers, scale_atmos_phy_bl_mynn_jmapplib::atmos_phy_bl_mynn_jmapplib_tendency(), scale_atmos_phy_bl_mynn::atmos_phy_bl_mynn_tendency(), mod_atmos_phy_bl_vars::atmos_phy_bl_ql, mod_atmos_phy_bl_vars::atmos_phy_bl_rhoq_t, mod_atmos_phy_bl_vars::atmos_phy_bl_rhot_t, mod_atmos_phy_bl_vars::atmos_phy_bl_rhou_t, mod_atmos_phy_bl_vars::atmos_phy_bl_rhov_t, mod_atmos_phy_bl_vars::atmos_phy_bl_sflx_buoy, scale_atmos_phy_bl_common::atmos_phy_bl_tendency_tracer(), mod_atmos_admin::atmos_phy_bl_type, mod_atmos_phy_bl_vars::atmos_phy_bl_zi, mod_atmos_phy_sf_vars::atmos_phy_sf_qstar, mod_atmos_phy_sf_vars::atmos_phy_sf_rlmo, mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_dens, mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_pres, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_mu, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_mv, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_qtrc, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_qv, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_sh, mod_atmos_phy_sf_vars::atmos_phy_sf_tstar, mod_atmos_phy_sf_vars::atmos_phy_sf_ustar, scale_bulkflux::bulkflux_type, mod_atmos_vars::dens, mod_atmos_vars::dens_av, mod_atmos_vars::exner, scale_atmos_hydrometeor::i_qv, 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, scale_landuse::landuse_frac_land, mod_atmos_vars::pott, mod_atmos_vars::pres, scale_tracer::qa, mod_atmos_vars::qc, mod_atmos_vars::qdry, mod_atmos_phy_bl_vars::qe, mod_atmos_vars::qi, mod_atmos_phy_bl_vars::qs, mod_atmos_vars::qtrc, mod_atmos_vars::qtrc_av, mod_atmos_vars::qv, mod_atmos_vars::rhoq_tp, mod_atmos_vars::rhot_tp, mod_atmos_vars::rhou_tp, mod_atmos_vars::rhov_tp, scale_statistics::statistics_checktotal, scale_time::time_dtsec_atmos_phy_bl, scale_tracer::tracer_advc, scale_tracer::tracer_mass, 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(), and mod_atmos_driver::atmos_driver_calc_tendency_from_sflux().
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0m
character(len=h_short), public atmos_phy_bl_type
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_dens
character(len=h_long), dimension(:), allocatable, public atmos_phy_bl_mynn_desc
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_mu
subroutine, public prc_abort
Abort Process.
subroutine, public atmos_phy_sf_bulk_flux(IA, IS, IE, JA, JS, JE, ATM_W, ATM_U, ATM_V, ATM_TEMP, ATM_PRES, ATM_QV, SFC_DENS, SFC_TEMP, SFC_PRES, SFC_Z0M, SFC_Z0H, SFC_Z0E, PBL, ATM_Z1, SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_LH, SFLX_QV, Ustar, Tstar, Qstar, Wstar, RLmo, U10, V10, T2, Q2)
Calculate surface flux.
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, public atmos_phy_sf_v10
real(rp), dimension(:,:,:), allocatable, target, public pott
character(len=h_long), dimension(1), public atmos_phy_bl_mynn_jmapplib_desc
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0e
real(rp), dimension(:,:,:,:), pointer, public qtrc_av
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_z1
Height of the lowermost grid from surface (cell center) [m].
real(rp), dimension(:,:), allocatable, public topography_tansl_y
tan(slope_y)
character(len=h_short), dimension(:), allocatable, public atmos_phy_bl_mynn_units
real(rp), dimension(:,:,:), allocatable, public rhov_tp
real(rp), dimension(:,:), allocatable, public atmos_phy_bl_zi
module atmosphere / physics / pbl / mynn
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_qstar
real(rp), dimension(:,:,:), allocatable, target, public qdry
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_ustar
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_wstar
subroutine, public atmos_phy_bl_mynn_mkinit(KA, KS, KE, IA, IS, IE, JA, JS, JE, PROG, DENS, U, V, W, POTT, PRES, EXNER, N2, QDRY, QV, Qw, POTL, POTV, SFC_DENS, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_QV, us, ts, qs, RLmo, frac_land, CZ, FZ, F2H, BULKFLUX_type)
ATMOS_PHY_BL_MYNN_mkinit initialize TKE.
real(rp), dimension(:,:), allocatable, public atmos_phy_bl_sflx_buoy
real(dp), public time_dtsec_atmos_phy_bl
time interval of physics(pbl ) [sec]
real(rp), dimension(:,:,:), allocatable, public atmos_phy_bl_rhou_t
character(len=h_short), dimension(1), public atmos_phy_bl_mynn_jmapplib_units
module atmosphere / hydrometeor
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_mv
real(rp), dimension(:,:,:,:), allocatable, target, public atmos_phy_bl_rhoq_t
module ATMOSPHERIC Surface Variables
real(rp), dimension(:,:,:), allocatable, target, public rhot
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0h
module Atmosphere GRID CartesC Real(real space)
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_temp
real(rp), dimension(:,:,:,:), allocatable, target, public qtrc
module atmosphere / physics / surface / bulk
subroutine, public atmos_phy_bl_mynn_tendency(KA, KS, KE, IA, IS, IE, JA, JS, JE, DENS, U, V, W, POTT, PROG, PRES, EXNER, N2, QDRY, QV, Qw, POTL, POTV, SFC_DENS, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_QV, us, ts, qs, RLmo, frac_land, CZ, FZ, F2H, dt_DP, BULKFLUX_type, RHOU_t, RHOV_t, RHOT_t, RHOQV_t, RPROG_t, Nu, Kh, Qlp, cldfrac, Zi, SFLX_BUOY)
ATMOS_PHY_BL_MYNN_tendency calculate tendency by the vertical eddy viscosity.
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_u10
real(rp), dimension(:,:,:), allocatable, public atmos_phy_bl_rhov_t
real(rp), dimension(:,:,:), allocatable, public rhou_tp
module atmosphere / physics / PBL
real(rp), dimension(:,:,:), allocatable, target, public dens
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_vol
control volume (zxy) [m3]
subroutine, public atmos_phy_bl_mynn_jmapplib_setup(KA, KS, KE, CZ, dt, PBL_MAX, SHCU_MAX)
ATMOS_PHY_BL_MYNN_JMAPPLIB_setup Setup.
real(rp), dimension(:,:), allocatable, public landuse_frac_land
land fraction
real(rp), dimension(:,:,:), allocatable, target, public momz
subroutine, public atmos_phy_bl_mynn_finalize
real(rp), dimension(:,:), allocatable, public topography_tansl_x
tan(slope_x)
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_rlmo
real(rp), dimension(:,:,:), allocatable, target, public v
subroutine, public atmos_phy_bl_tendency_tracer(KA, KS, KE, IA, IS, IE, JA, JS, JE, DENS, QTRC, SFLX_Q, Kh, MASS, CZ, FZ, F2H, DDT, TRACER_NAME, RHOQ_t)
ATMOS_PHY_BL_tendency_tracer calculate tendency of tracers by the eddy viscosity.
real(rp), dimension(:,:,:), allocatable, target, public w
real(rp), dimension(:,:,:), allocatable, target, public momx
real(rp), dimension(:,:,:), allocatable, target, public exner
real(rp), dimension(:,:,:), allocatable, target, public temp
real(rp), dimension(:,:), pointer, public atmos_phy_sf_sflx_qv
real(rp), dimension(:,:,:), allocatable, target, public momy
character(len=h_short), public bulkflux_type
real(rp), dimension(:,:,:), allocatable, pointer, target, public qv
real(rp), dimension(:,:,:), allocatable, target, public atmos_phy_sf_sflx_qtrc
subroutine, public tracer_regist(QS, NQ, NAME, DESC, UNIT, CV, CP, R, ENGI0, ADVC, MASS)
Regist tracer.
integer, parameter, public atmos_phy_bl_mynn_jmapplib_ntracer
real(rp), dimension(:,:,:), allocatable, target, public pres
real(rp), dimension(:,:,:), pointer, public dens_av
real(rp), dimension(:,:,:), allocatable, target, public u
real(rp), dimension(:,:,:,:), allocatable, public atmos_grid_cartesc_real_f2h
coefficient for interpolation from full to half levels
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_t2
character(len=h_short), dimension(1), public atmos_phy_bl_mynn_jmapplib_name
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_lh
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_pres
logical, public statistics_checktotal
calc&report variable totals to logfile?
module ATMOSPHERIC Variables
module atmosphere / bottom boundary extrapolation
real(rp), dimension(:,:,:), pointer, public qi
real(rp), dimension(:,:,:), allocatable, public rhot_tp
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_sh
module atmosphere / physics / pbl / mynn-jmapplib
logical, public atmos_sw_phy_bl
integer, public atmos_phy_bl_mynn_ntracer
subroutine, public atmos_vars_calc_diagnostics
Calc diagnostic variables.
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_q2
real(rp), dimension(:,:,:), pointer, public qc
subroutine, public atmos_phy_bl_mynn_jmapplib_tendency(KA, KS, KE, IA, IS, IE, JA, JS, JE, DENS, U, V, POTT, PROG, PRES, EXNER, QDRY, QV, QC, QI, SFC_DENS, SFC_PRES, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_QV, us, RLmo, CZ, FZ, F2H, dt, RHOU_t, RHOV_t, RHOT_t, RHOQV_t, RPROG_t, Nu, Kh, Zi, SFLX_BUOY)
ATMOS_PHY_BL_MYNN_JMAPPLIB_tendency calculate tendency by the virtical eddy viscosity.
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]
real(rp), dimension(:,:,:), allocatable, public atmos_phy_bl_cldfrac
module atmosphere / grid / cartesC
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cz
center coordinate [m]: z, local
real(rp), dimension(:,:,:), allocatable, public atmos_phy_bl_rhot_t
real(rp), dimension(:,:,:), allocatable, public atmos_phy_bl_ql
subroutine, public atmos_phy_bl_mynn_jmapplib_finalize
real(rp), public const_undef
character(len=h_short), dimension(:), allocatable, public atmos_phy_bl_mynn_name
subroutine, public atmos_bottom_estimate(KA, KS, KE, IA, IS, IE, JA, JS, JE, DENS, PRES, QV, SFC_TEMP, FZ, SFC_DENS, SFC_PRES)
Calc bottom boundary of atmosphere (just above surface)
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_tstar
logical, public atmos_phy_bl_mix_tracers
module atmosphere / physics / pbl / common
subroutine, public atmos_phy_bl_mynn_setup(BULKFLUX_type, dx, TKE_MIN, PBL_MAX)
ATMOS_PHY_BL_MYNN_setup Setup.
subroutine, public atmos_phy_bl_mynn_tracer_setup
ATMOS_PHY_BL_MYNN_tracer_setup Tracer Setup.