module ATMOSPHERE / Physics Surface fluxes
More...
module ATMOSPHERE / Physics Surface fluxes
- Description
- Flux from/to bottom boundary of atmosphere (surface)
- Author
- Team SCALE
- NAMELIST
-
- History Output
name | description | unit | variable |
SFLX_{TRACER_NAME} | surface {TRACER_NAME} flux;
{TRACER_NAME} depends on the physics schemes, e.g., QV, QC, QR. | kg/m2/s | SFLX_QTRC |
GHFLX | ground heat flux (downward) | W/m2 | SFLX_GH |
LHFLX | latent heat flux | W/m2 | SFLX_LH |
MSLP | mean sea-level pressure | Pa | MSLP |
MUFLX | u-momentum flux | kg/m/s2 | SFLX_MU |
MVFLX | v-momentum flux | kg/m/s2 | SFLX_MV |
MWFLX | w-momentum flux | kg/m/s2 | SFLX_MW |
Q2 | 2m specific humidity | kg/kg | Q2 |
Qstar | moisuter scale | kg/kg | Qstar |
RLmo | inverse of Obukhov length | 1/m | RLmo |
SFC_ALB_IR_dif | surface albedo (IR; diffuse) | 1 | SFC_albedo |
SFC_ALB_IR_dir | surface albedo (IR; direct | 1 | SFC_albedo |
SFC_ALB_NIR_dif | surface albedo (NIR; diffuse | 1 | SFC_albedo |
SFC_ALB_NIR_dir | surface albedo (NIR; direct | 1 | SFC_albedo |
SFC_ALB_VIS_dif | surface albedo (VIS; diffuse | 1 | SFC_albedo |
SFC_ALB_VIS_dir | surface albedo (VIS; direct | 1 | SFC_albedo |
SFC_DENS | surface atmospheric density | kg/m3 | SFC_DENS |
SFC_PRES | surface atmospheric pressure | Pa | SFC_PRES |
SFC_TEMP | surface skin temperature | K | SFC_TEMP |
SFC_Z0E | roughness length (vapor) | m | SFC_Z0E |
SFC_Z0H | roughness length (heat) | m | SFC_Z0H |
SFC_Z0M | roughness length (momentum) | m | SFC_Z0M |
SFLX_ENGI | ground internal energy flux (merged) | W/m2 | SFLX_ENGI |
SHFLX | sensible heat flux | W/m2 | SFLX_SH |
T2 | 2m air temperature | K | T2 |
Tstar | temperature scale | K | Tstar |
U10 | 10m x-wind | m/s | U10 |
U10m | 10m eastward wind | m/s | U10m |
Uabs10 | 10m absolute wind | m/s | Uabs10 |
Ustar | friction velocity | m/s | Ustar |
V10 | 10m y-wind | m/s | V10 |
V10m | 10m northward wind | m/s | V10m |
Wstar | convective velocity scale | m/s | Wstar |
◆ atmos_phy_sf_driver_setup()
subroutine, public mod_atmos_phy_sf_driver::atmos_phy_sf_driver_setup |
Setup.
Definition at line 50 of file mod_atmos_phy_sf_driver.F90.
83 log_info(
"ATMOS_PHY_SF_driver_setup",*)
'Setup'
88 log_info(
"ATMOS_PHY_SF_driver_setup",*)
'Coupler is enabled.'
97 log_error(
"ATMOS_PHY_SF_driver_setup",*)
'invalid Surface flux type(', trim(
atmos_phy_sf_type),
'). CHECK!'
104 log_info(
"ATMOS_PHY_SF_driver_setup",*)
'this component is never called.'
105 log_info(
"ATMOS_PHY_SF_driver_setup",*)
'surface fluxes are set to zero.'
106 sflx_mw(:,:) = 0.0_rp
107 sflx_mu(:,:) = 0.0_rp
108 sflx_mv(:,:) = 0.0_rp
109 sflx_sh(:,:) = 0.0_rp
110 sflx_lh(:,:) = 0.0_rp
111 sflx_shex(:,:) = 0.0_rp
112 sflx_qvex(:,:) = 0.0_rp
118 log_info(
"ATMOS_PHY_SF_driver_setup",*)
'SFC_TEMP, SFC_albedo is set in ATMOS_PHY_SF_vars.'
122 sflx_qtrc(:,:,:) = 0.0_rp
123 sflx_engi(:,:) = 0.0_rp
References scale_atmos_phy_sf_bulk::atmos_phy_sf_bulk_setup(), scale_atmos_phy_sf_const::atmos_phy_sf_const_setup(), 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_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_engi, 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_mw, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_qtrc, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_qvex, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_sh, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_shex, mod_atmos_phy_sf_vars::atmos_phy_sf_tstar, mod_atmos_admin::atmos_phy_sf_type, mod_atmos_phy_sf_vars::atmos_phy_sf_ustar, mod_atmos_phy_sf_vars::atmos_phy_sf_wstar, mod_atmos_admin::atmos_sw_phy_sf, mod_cpl_admin::cpl_sw, and scale_prc::prc_abort().
Referenced by mod_atmos_driver::atmos_driver_setup().
◆ atmos_phy_sf_driver_calc_tendency()
subroutine, public mod_atmos_phy_sf_driver::atmos_phy_sf_driver_calc_tendency |
( |
logical, intent(in) |
update_flag | ) |
|
calculation tendency
Definition at line 131 of file mod_atmos_phy_sf_driver.F90.
159 bulkflux_diagnose_scales
220 logical,
intent(in) :: update_flag
222 real(RP) :: ATM_W (IA,JA)
223 real(RP) :: ATM_U (IA,JA)
224 real(RP) :: ATM_V (IA,JA)
225 real(RP) :: ATM_DENS(IA,JA)
226 real(RP) :: ATM_TEMP(IA,JA)
227 real(RP) :: ATM_PRES(IA,JA)
228 real(RP) :: ATM_QV (IA,JA)
229 real(RP) :: SFLX_SH2(IA,JA)
230 real(RP) :: SFLX_QV (IA,JA)
231 real(RP) :: CP_t, CV_t
239 if ( update_flag )
then
242 call bottom_estimate( ka, ks, ke, ia, is, ie, ja, js, je, &
246 sfc_dens(:,:), sfc_pres(:,:) )
251 atm_dens(i,j) =
dens(ks,i,j)
252 atm_temp(i,j) =
temp(ks,i,j)
253 atm_pres(i,j) =
pres(ks,i,j)
262 sflx_sh2(i,j) = sflx_sh(i,j) - sflx_shex(i,j)
270 sflx_qv(i,j) = 0.0_rp
277 sflx_qv(i,j) = sflx_qtrc(i,j,
i_qv) - sflx_qvex(i,j)
282 call bulkflux_diagnose_scales( ia, is, ie, ja, js, je, &
283 sflx_mw(:,:), sflx_mu(:,:), sflx_mv(:,:), &
284 sflx_sh2(:,:), sflx_qv(:,:), &
285 sfc_dens(:,:), sfc_temp(:,:), pbl_zi(:,:), &
286 ustar(:,:), tstar(:,:), qstar(:,:), &
287 wstar(:,:), rlmo(:,:) )
293 atm_u(i,j) =
u(ks,i,j)
294 atm_v(i,j) =
v(ks,i,j)
295 atm_w(i,j) = atm_u(i,j) * tansl_x(i,j) + atm_v(i,j) * tansl_y(i,j)
296 atm_qv(i,j) =
qv(ks,i,j)
304 atm_w(:,:), atm_u(:,:), atm_v(:,:), &
305 atm_temp(:,:), atm_pres(:,:), atm_qv(:,:), &
306 sfc_dens(:,:), sfc_temp(:,:), sfc_pres(:,:), &
307 sfc_z0m(:,:), sfc_z0h(:,:), sfc_z0e(:,:), &
308 pbl_zi(:,:), z1(:,:), &
309 sflx_mw(:,:), sflx_mu(:,:), sflx_mv(:,:), &
310 sflx_sh(:,:), sflx_lh(:,:), sflx_qv(:,:), &
311 ustar(:,:), tstar(:,:), qstar(:,:), &
314 u10(:,:), v10(:,:), t2(:,:), q2(:,:) )
319 atm_w(:,:), atm_u(:,:), atm_v(:,:), sfc_temp(:,:), &
320 z1(:,:), sfc_dens(:,:), &
321 sflx_mw(:,:), sflx_mu(:,:), sflx_mv(:,:), &
322 sflx_sh(:,:), sflx_lh(:,:), sflx_qv(:,:), &
328 t2(:,:) = atm_temp(:,:)
329 q2(:,:) = atm_qv(:,:)
334 sflx_qtrc(:,:,
i_qv) = sflx_qv(:,:)
335 sflx_engi(:,:) = sflx_qv(:,:) * ( tracer_cv(
i_qv) * sfc_temp(:,:) +
lhv )
347 rdz = 1.0_rp / ( fz(ks,i,j) - fz(ks-1,i,j) )
348 momz_t_sf(i,j) = sflx_mw(i,j) / ( cz(ks+1,i,j) - cz(ks,i,j) )
349 rhou_t_sf(i,j) = sflx_mu(i,j) * rdz
350 rhov_t_sf(i,j) = sflx_mv(i,j) * rdz
358 rdz = 1.0_rp / ( fz(ks,i,j) - fz(ks-1,i,j) )
359 dens_t_sf(i,j) = 0.0_rp
362 engi_t = sflx_engi(i,j) * rdz
364 work = sflx_qtrc(i,j,iq) * rdz
366 rhoq_t_sf(i,j,iq) = work
367 dens_t_sf(i,j) = dens_t_sf(i,j) + work * tracer_mass(iq)
368 cp_t = cp_t + work * tracer_cp(iq)
369 cv_t = cv_t + work * tracer_cv(iq)
370 engi_t = engi_t - tracer_engi0(iq) * work
372 cp_t = ( cp_t -
cptot(ks,i,j) * dens_t_sf(i,j) ) / atm_dens(i,j)
373 cv_t = ( cv_t -
cvtot(ks,i,j) * dens_t_sf(i,j) ) / atm_dens(i,j)
375 rhoh_sf(i,j) = sflx_sh(i,j) * rdz + engi_t &
376 - ( cp_t + log( atm_pres(i,j) / pre00 ) * (
cvtot(ks,i,j) /
cptot(ks,i,j) * cp_t - cv_t ) ) * atm_dens(i,j) * atm_temp(i,j)
385 momz_t(ks,i,j) = momz_t(ks,i,j) + momz_t_sf(i,j)
386 rhou_t(ks,i,j) = rhou_t(ks,i,j) + rhou_t_sf(i,j)
387 rhov_t(ks,i,j) = rhov_t(ks,i,j) + rhov_t_sf(i,j)
388 rhoh(ks,i,j) = rhoh(ks,i,j) + rhoh_sf(i,j)
389 dens_t(ks,i,j) = dens_t(ks,i,j) + dens_t_sf(i,j)
398 rhoq_t(ks,i,j,iq) = rhoq_t(ks,i,j,iq) + rhoq_t_sf(i,j,iq)
409 call statistics_total( ia, is, ie, ja, js, je, &
410 sflx_qtrc(:,:,iq),
'SFLX_'//trim(tracer_name(iq)), &
416 call statistics_total( ia, is, ie, ja, js, je, &
417 dens_t_sf(:,:),
'DENS_t_SF', &
420 call statistics_total( ia, is, ie, ja, js, je, &
421 momz_t_sf(:,:),
'MOMZ_t_SF', &
424 call statistics_total( ia, is, ie, ja, js, je, &
425 rhou_t_sf(:,:),
'RHOU_t_SF', &
428 call statistics_total( ia, is, ie, ja, js, je, &
429 rhov_t_sf(:,:),
'RHOV_t_SF', &
432 call statistics_total( ia, is, ie, ja, js, je, &
433 rhoh_sf(:,:),
'RHOH_SF', &
438 call statistics_total( ia, is, ie, ja, js, je, &
439 rhoq_t_sf(:,:,iq), trim(tracer_name(iq))//
'_t_SF', &
References scale_atmos_bottom::atmos_bottom_estimate(), scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_area, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_cz, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_fz, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_totarea, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_z1, scale_atmos_hydrometeor::atmos_hydrometeor_dry, mod_atmos_phy_bl_vars::atmos_phy_bl_zi, mod_atmos_phy_rd_vars::atmos_phy_rd_sflx_lw_dn, mod_atmos_phy_rd_vars::atmos_phy_rd_sflx_sw_dn, scale_atmos_phy_sf_bulk::atmos_phy_sf_bulk_flux(), scale_atmos_phy_sf_const::atmos_phy_sf_const_flux(), mod_atmos_phy_sf_vars::atmos_phy_sf_dens_t, mod_atmos_phy_sf_vars::atmos_phy_sf_momz_t, 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_rhoh, mod_atmos_phy_sf_vars::atmos_phy_sf_rhoq_t, mod_atmos_phy_sf_vars::atmos_phy_sf_rhou_t, mod_atmos_phy_sf_vars::atmos_phy_sf_rhov_t, 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_engi, 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_mw, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_qtrc, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_qvex, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_sh, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_shex, mod_atmos_phy_sf_vars::atmos_phy_sf_t2, mod_atmos_phy_sf_vars::atmos_phy_sf_tstar, mod_atmos_admin::atmos_phy_sf_type, 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, scale_const::const_pre00, scale_const::const_undef, mod_cpl_admin::cpl_sw, mod_atmos_vars::cptot, mod_atmos_vars::cvtot, mod_atmos_vars::dens, mod_atmos_vars::dens_av, mod_atmos_vars::dens_tp, history_output(), scale_atmos_hydrometeor::i_qv, 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_atmos_hydrometeor::lhv, mod_atmos_vars::momz_tp, mod_atmos_vars::pott, mod_atmos_vars::pres, scale_tracer::qa, mod_atmos_vars::qv, mod_atmos_vars::rhoh_p, mod_atmos_vars::rhoq_tp, mod_atmos_vars::rhot, mod_atmos_vars::rhot_av, mod_atmos_vars::rhou_tp, mod_atmos_vars::rhov_tp, scale_statistics::statistics_checktotal, mod_atmos_vars::temp, scale_time::time_dtsec_atmos_phy_sf, scale_topography::topography_tansl_x, scale_topography::topography_tansl_y, scale_tracer::tracer_cp, scale_tracer::tracer_cv, scale_tracer::tracer_engi0, 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().
◆ history_output()
subroutine mod_atmos_phy_sf_driver::history_output |
Definition at line 451 of file mod_atmos_phy_sf_driver.F90.
456 barometric_law_mslp => atmos_hydrostatic_barometric_law_mslp
495 real(RP) :: MSLP (IA,JA)
496 real(RP) :: Uabs10(IA,JA)
497 real(RP) :: U10m (IA,JA)
498 real(RP) :: V10m (IA,JA)
507 uabs10(i,j) = sqrt( u10(i,j)**2 + v10(i,j)**2 )
508 u10m(i,j) = u10(i,j) * rotc(i,j,1) - v10(i,j) * rotc(i,j,2)
509 v10m(i,j) = u10(i,j) * rotc(i,j,2) + v10(i,j) * rotc(i,j,1)
514 call barometric_law_mslp( ka, ks, ke, ia, is, ie, ja, js, je, &
519 call file_history_in( sfc_dens(:,:),
'SFC_DENS',
'surface atmospheric density',
'kg/m3' )
520 call file_history_in( sfc_pres(:,:),
'SFC_PRES',
'surface atmospheric pressure',
'Pa' )
521 call file_history_in( sfc_temp(:,:),
'SFC_TEMP',
'surface skin temperature',
'K' )
522 call file_history_in( sfc_albedo(:,:,i_r_direct ,i_r_ir ),
'SFC_ALB_IR_dir' ,
'surface albedo (IR; direct',
'1' , fill_halo=.true. )
523 call file_history_in( sfc_albedo(:,:,i_r_diffuse,i_r_ir ),
'SFC_ALB_IR_dif' ,
'surface albedo (IR; diffuse)',
'1' , fill_halo=.true. )
524 call file_history_in( sfc_albedo(:,:,i_r_direct ,i_r_nir),
'SFC_ALB_NIR_dir',
'surface albedo (NIR; direct',
'1' , fill_halo=.true. )
525 call file_history_in( sfc_albedo(:,:,i_r_diffuse,i_r_nir),
'SFC_ALB_NIR_dif',
'surface albedo (NIR; diffuse',
'1' , fill_halo=.true. )
526 call file_history_in( sfc_albedo(:,:,i_r_direct ,i_r_vis),
'SFC_ALB_VIS_dir',
'surface albedo (VIS; direct',
'1' , fill_halo=.true. )
527 call file_history_in( sfc_albedo(:,:,i_r_diffuse,i_r_vis),
'SFC_ALB_VIS_dif',
'surface albedo (VIS; diffuse',
'1' , fill_halo=.true. )
528 call file_history_in( sfc_z0m(:,:),
'SFC_Z0M',
'roughness length (momentum)',
'm' , fill_halo=.true. )
529 call file_history_in( sfc_z0h(:,:),
'SFC_Z0H',
'roughness length (heat)',
'm' , fill_halo=.true. )
530 call file_history_in( sfc_z0e(:,:),
'SFC_Z0E',
'roughness length (vapor)',
'm' , fill_halo=.true. )
531 call file_history_in( sflx_mw(:,:),
'MWFLX',
'w-momentum flux',
'kg/m/s2' )
532 call file_history_in( sflx_mu(:,:),
'MUFLX',
'u-momentum flux',
'kg/m/s2' )
533 call file_history_in( sflx_mv(:,:),
'MVFLX',
'v-momentum flux',
'kg/m/s2' )
534 call file_history_in( sflx_sh(:,:),
'SHFLX',
'sensible heat flux',
'W/m2' , fill_halo=.true. )
535 call file_history_in( sflx_lh(:,:),
'LHFLX',
'latent heat flux',
'W/m2' , fill_halo=.true. )
536 call file_history_in( sflx_gh(:,:),
'GHFLX',
'ground heat flux (downward)',
'W/m2' , fill_halo=.true. )
538 call file_history_in( sflx_qtrc(:,:,iq),
'SFLX_'//trim(tracer_name(iq)), &
539 'surface '//trim(tracer_name(iq))//
' flux', &
540 'kg/m2/s' , fill_halo=.true. )
542 call file_history_in( sflx_engi(:,:),
'SFLX_ENGI',
'ground internal energy flux (merged)',
'W/m2' , fill_halo=.true. )
544 call file_history_in( ustar(:,:),
'Ustar',
'friction velocity',
'm/s' , fill_halo=.true. )
545 call file_history_in( tstar(:,:),
'Tstar',
'temperature scale',
'K' , fill_halo=.true. )
546 call file_history_in( qstar(:,:),
'Qstar',
'moisuter scale',
'kg/kg', fill_halo=.true. )
547 call file_history_in( wstar(:,:),
'Wstar',
'convective velocity scale',
'm/s', fill_halo=.true. )
548 call file_history_in( rlmo(:,:),
'RLmo',
'inverse of Obukhov length',
'1/m' , fill_halo=.true. )
549 call file_history_in( uabs10(:,:),
'Uabs10',
'10m absolute wind',
'm/s' , fill_halo=.true. )
550 call file_history_in( u10(:,:),
'U10',
'10m x-wind',
'm/s' , fill_halo=.true. )
551 call file_history_in( v10(:,:),
'V10',
'10m y-wind',
'm/s' , fill_halo=.true. )
552 call file_history_in( u10m(:,:),
'U10m',
'10m eastward wind',
'm/s' , fill_halo=.true. )
553 call file_history_in( v10m(:,:),
'V10m',
'10m northward wind',
'm/s' , fill_halo=.true. )
554 call file_history_in( t2(:,:),
'T2 ',
'2m air temperature',
'K' , fill_halo=.true. )
555 call file_history_in( q2(:,:),
'Q2 ',
'2m specific humidity',
'kg/kg', fill_halo=.true. )
556 call file_history_in( mslp(:,:),
'MSLP',
'mean sea-level pressure',
'Pa' , fill_halo=.true., standard_name=
'air_pressure_at_mean_sea_level' )
References scale_atmos_grid_cartesc_metric::atmos_grid_cartesc_metric_rotc, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_cz, scale_atmos_hydrometeor::atmos_hydrometeor_dry, 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_albedo, 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_engi, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_gh, 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_mw, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_qtrc, 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, scale_const::const_undef, scale_atmos_hydrometeor::i_qv, scale_cpl_sfc_index::i_r_diffuse, scale_cpl_sfc_index::i_r_direct, scale_cpl_sfc_index::i_r_ir, scale_cpl_sfc_index::i_r_nir, scale_cpl_sfc_index::i_r_vis, scale_atmos_grid_cartesc_index::ie, scale_atmos_grid_cartesc_index::is, 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, mod_atmos_vars::pres, scale_tracer::qa, mod_atmos_vars::qv, mod_atmos_vars::temp, and scale_tracer::tracer_name.
Referenced by atmos_phy_sf_driver_calc_tendency().
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0m
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_dens
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
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_z0e
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)
real(rp), dimension(:,:,:), allocatable, public rhov_tp
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_qstar
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_rhoh
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_ustar
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_wstar
real(rp), public atmos_grid_cartesc_real_totarea
total area (xy, local) [m2]
real(rp), dimension(:,:,:), allocatable, target, public cptot
module Atmosphere Grid CartesianC metirc
real(rp), dimension(:,:,:), pointer, public rhot_av
module atmosphere / hydrometeor
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_mv
module Atmosphere / Physics Radiation
module ATMOSPHERIC Surface Variables
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_mw
logical, public atmos_hydrometeor_dry
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
module atmosphere / physics / surface / bulk
subroutine, public atmos_phy_sf_bulk_setup
Setup.
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_qvex
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_u10
real(rp), dimension(:,:,:), allocatable, public rhou_tp
character(len=h_short), public atmos_phy_sf_type
real(rp), dimension(:,:), allocatable, public atmos_phy_rd_sflx_lw_dn
module atmosphere / physics / PBL
real(rp), dimension(:,:), allocatable, public atmos_phy_rd_sflx_sw_dn
real(rp), dimension(:,:,:), allocatable, target, public dens
real(rp), dimension(:,:,:), allocatable, public atmos_phy_sf_rhoq_t
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_shex
real(rp), dimension(:,:,:,:), allocatable, public atmos_phy_sf_sfc_albedo
real(rp), public lhv
latent heat of vaporization for use [J/kg]
real(dp), public time_dtsec_atmos_phy_sf
time interval of physics(surface flux) [sec]
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
real(rp), dimension(:,:), allocatable, public atmos_phy_bl_zi
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
subroutine, public atmos_phy_sf_const_setup
Setup.
real(rp), dimension(:,:,:), allocatable, target, public temp
real(rp), dimension(:,:,:), allocatable, public dens_tp
real(rp), dimension(:,:,:), allocatable, pointer, target, public qv
real(rp), dimension(:,:,:), allocatable, target, public atmos_phy_sf_sflx_qtrc
real(rp), dimension(:,:,:), allocatable, target, public pres
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_engi
real(rp), dimension(:,:,:), pointer, public dens_av
real(rp), dimension(:,:,:), allocatable, target, public u
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_t2
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_lh
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_dens_t
module atmosphere / hydrostatic barance
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_pres
real(rp), dimension(:,:,:), allocatable, public rhoh_p
logical, public statistics_checktotal
calc&report variable totals to logfile?
module ATMOSPHERIC Variables
module atmosphere / bottom boundary extrapolation
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_sh
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_q2
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_fz
geopotential height [m] (wxy)
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_gh
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_rhou_t
real(rp), public const_undef
subroutine, public atmos_phy_sf_const_flux(IA, IS, IE, JA, JS, JE, ATM_W, ATM_U, ATM_V, ATM_TEMP, ATM_Z1, SFC_DENS, SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_LH, SFLX_QV, U10, V10)
Constant flux.
real(rp), public const_pre00
pressure reference [Pa]
real(rp), dimension(:,:,:), allocatable, target, public cvtot
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)
module atmosphere / physics / surface / const
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_tstar
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_momz_t
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_rhov_t
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_metric_rotc
rotation coefficient
logical, public atmos_sw_phy_sf