CCPP Scidoc for SRW v2.1.0  SRW v2.1.0
Common Community Physics Package Developed at DTC

◆ satmedmfvdifq_run()

subroutine satmedmfvdifq::satmedmfvdifq_run ( integer, intent(in)  im,
integer, intent(in)  km,
logical, intent(in)  progsigma,
integer, intent(in)  ntrac,
integer, intent(in)  ntcw,
integer, intent(in)  ntrw,
integer, intent(in)  ntiw,
integer, intent(in)  ntke,
real(kind=kind_phys), intent(in)  grav,
real(kind=kind_phys), intent(in)  rd,
real(kind=kind_phys), intent(in)  cp,
real(kind=kind_phys), intent(in)  rv,
real(kind=kind_phys), intent(in)  hvap,
real(kind=kind_phys), intent(in)  hfus,
real(kind=kind_phys), intent(in)  fv,
real(kind=kind_phys), intent(in)  eps,
real(kind=kind_phys), intent(in)  epsm1,
real(kind=kind_phys), dimension(:,:), intent(inout)  dv,
real(kind=kind_phys), dimension(:,:), intent(inout)  du,
real(kind=kind_phys), dimension(:,:), intent(inout)  tdt,
real(kind=kind_phys), dimension(:,:,:), intent(inout)  rtg,
real(kind=kind_phys), dimension(:,:), intent(inout)  tmf,
real(kind=kind_phys), dimension(:,:), intent(in)  u1,
real(kind=kind_phys), dimension(:,:), intent(in)  v1,
real(kind=kind_phys), dimension(:,:), intent(in)  t1,
real(kind=kind_phys), dimension(:,:,:), intent(in)  q1,
real(kind=kind_phys), dimension(:,:), intent(in)  swh,
real(kind=kind_phys), dimension(:,:), intent(in)  hlw,
real(kind=kind_phys), dimension(:), intent(in)  xmu,
real(kind=kind_phys), dimension(:), intent(in)  garea,
real(kind=kind_phys), dimension(:), intent(in)  zvfun,
real(kind=kind_phys), dimension(:), intent(in)  psk,
real(kind=kind_phys), dimension(:), intent(in)  rbsoil,
real(kind=kind_phys), dimension(:), intent(in)  zorl,
real(kind=kind_phys), dimension(:), intent(in)  u10m,
real(kind=kind_phys), dimension(:), intent(in)  v10m,
real(kind=kind_phys), dimension(:), intent(in)  fm,
real(kind=kind_phys), dimension(:), intent(in)  fh,
real(kind=kind_phys), dimension(:), intent(in)  tsea,
real(kind=kind_phys), dimension(:), intent(in)  heat,
real(kind=kind_phys), dimension(:), intent(in)  evap,
real(kind=kind_phys), dimension(:), intent(in)  stress,
real(kind=kind_phys), dimension(:), intent(in)  spd1,
integer, dimension(:), intent(out)  kpbl,
real(kind=kind_phys), dimension(:,:), intent(in)  prsi,
real(kind=kind_phys), dimension(:,:), intent(in)  del,
real(kind=kind_phys), dimension(:,:), intent(in)  prsl,
real(kind=kind_phys), dimension(:,:), intent(in)  prslk,
real(kind=kind_phys), dimension(:,:), intent(in)  phii,
real(kind=kind_phys), dimension(:,:), intent(in)  phil,
real(kind=kind_phys), intent(in)  delt,
logical, intent(in)  dspheat,
real(kind=kind_phys), dimension(:), intent(out)  dusfc,
real(kind=kind_phys), dimension(:), intent(out)  dvsfc,
real(kind=kind_phys), dimension(:), intent(out)  dtsfc,
real(kind=kind_phys), dimension(:), intent(out)  dqsfc,
real(kind=kind_phys), dimension(:), intent(out)  hpbl,
real(kind=kind_phys), dimension(:,:), intent(out)  dkt,
real(kind=kind_phys), dimension(:,:), intent(out)  dku,
integer, dimension(:), intent(in)  kinver,
real(kind=kind_phys), intent(in)  xkzm_m,
real(kind=kind_phys), intent(in)  xkzm_h,
real(kind=kind_phys), intent(in)  xkzm_s,
real(kind=kind_phys), intent(in)  dspfac,
real(kind=kind_phys), intent(in)  bl_upfr,
real(kind=kind_phys), intent(in)  bl_dnfr,
real(kind=kind_phys), intent(in)  rlmx,
real(kind=kind_phys), intent(in)  elmx,
integer, intent(in)  sfc_rlm,
integer, intent(in)  ntqv,
real(kind=kind_phys), dimension(:,:,:), intent(inout)  dtend,
integer, dimension(:,:), intent(in)  dtidx,
integer, intent(in)  index_of_temperature,
integer, intent(in)  index_of_x_wind,
integer, intent(in)  index_of_y_wind,
integer, intent(in)  index_of_process_pbl,
logical, intent(in)  gen_tend,
logical, intent(in)  ldiag3d,
character(len=*), intent(out)  errmsg,
integer, intent(out)  errflg 
)

Argument Table

satmedmfvdifq_run argument table
local_namestandard_namelong_nameunitstypedimensionskindintent
imhorizontal_loop_extenthorizontal loop extentcountinteger()in
kmvertical_layer_dimensionvertical layer dimensioncountinteger()in
progsigmado_prognostic_updraft_area_fractionflag for prognostic sigma in cumuls schemeflaglogical()in
ntracnumber_of_vertical_diffusion_tracersnumber of tracers to diffuse verticallycountinteger()in
ntcwindex_for_liquid_cloud_condensate_vertical_diffusion_tracertracer index for cloud condensate (or liquid water)indexinteger()in
ntrwindex_for_rain_water_vertical_diffusion_tracertracer index for rain water in the vertically diffused tracer arrayindexinteger()in
ntiwindex_for_ice_cloud_condensate_vertical_diffusion_tracertracer index for ice water in the vertically diffused tracer arrayindexinteger()in
ntkeindex_for_turbulent_kinetic_energy_vertical_diffusion_tracerindex for turbulent kinetic energy in the vertically diffused tracer arrayindexinteger()in
gravgravitational_accelerationgravitational accelerationm s-2real()kind_physin
rdgas_constant_of_dry_airideal gas constant for dry airJ kg-1 K-1real()kind_physin
cpspecific_heat_of_dry_air_at_constant_pressurespecific heat of dry air at constant pressureJ kg-1 K-1real()kind_physin
rvgas_constant_water_vaporideal gas constant for water vaporJ kg-1 K-1real()kind_physin
hvaplatent_heat_of_vaporization_of_water_at_0clatent heat of evaporation/sublimationJ kg-1real()kind_physin
hfuslatent_heat_of_fusion_of_water_at_0clatent heat of fusionJ kg-1real()kind_physin
fvratio_of_vapor_to_dry_air_gas_constants_minus_one(rv/rd) - 1 (rv = ideal gas constant for water vapor)nonereal()kind_physin
epsratio_of_dry_air_to_water_vapor_gas_constantsrd/rvnonereal()kind_physin
epsm1ratio_of_dry_air_to_water_vapor_gas_constants_minus_one(rd/rv) - 1nonereal()kind_physin
dvprocess_split_cumulative_tendency_of_y_windupdated tendency of the y windm s-2real(ccpp_constant_one:horizontal_loop_extent, ccpp_constant_one:vertical_layer_dimension)kind_physinout
duprocess_split_cumulative_tendency_of_x_windupdated tendency of the x windm s-2real(ccpp_constant_one:horizontal_loop_extent, ccpp_constant_one:vertical_layer_dimension)kind_physinout
tdtprocess_split_cumulative_tendency_of_air_temperatureupdated tendency of the temperatureK s-1real(ccpp_constant_one:horizontal_loop_extent, ccpp_constant_one:vertical_layer_dimension)kind_physinout
rtgtendency_of_vertically_diffused_tracer_concentrationupdated tendency of the tracers due to vertical diffusion in PBL schemekg kg-1 s-1real(ccpp_constant_one:horizontal_loop_extent, ccpp_constant_one:vertical_layer_dimension, ccpp_constant_one:number_of_vertical_diffusion_tracers)kind_physinout
tmfinstantaneous_tendency_of_specific_humidity_due_to_pblinstantaneous_tendency_of_specific_humidity_due_to_PBLkg kg-1 s-1real(ccpp_constant_one:horizontal_loop_extent, ccpp_constant_one:vertical_layer_dimension)kind_physinout
u1x_windx component of layer windm s-1real(ccpp_constant_one:horizontal_loop_extent, ccpp_constant_one:vertical_layer_dimension)kind_physin
v1y_windy component of layer windm s-1real(ccpp_constant_one:horizontal_loop_extent, ccpp_constant_one:vertical_layer_dimension)kind_physin
t1air_temperaturelayer mean air temperatureKreal(ccpp_constant_one:horizontal_loop_extent, ccpp_constant_one:vertical_layer_dimension)kind_physin
q1vertically_diffused_tracer_concentrationtracer concentration diffused by PBL schemekg kg-1real(ccpp_constant_one:horizontal_loop_extent, ccpp_constant_one:vertical_layer_dimension, ccpp_constant_one:number_of_vertical_diffusion_tracers)kind_physin
swhtendency_of_air_temperature_due_to_shortwave_heating_on_radiation_timesteptotal sky shortwave heating rateK s-1real(ccpp_constant_one:horizontal_loop_extent, ccpp_constant_one:vertical_layer_dimension)kind_physin
hlwtendency_of_air_temperature_due_to_longwave_heating_on_radiation_timesteptotal sky longwave heating rateK s-1real(ccpp_constant_one:horizontal_loop_extent, ccpp_constant_one:vertical_layer_dimension)kind_physin
xmuzenith_angle_temporal_adjustment_factor_for_shortwave_fluxeszenith angle temporal adjustment factor for shortwavenonereal(ccpp_constant_one:horizontal_loop_extent)kind_physin
gareacell_areaarea of the grid cellm2real(ccpp_constant_one:horizontal_loop_extent)kind_physin
zvfunfunction_of_surface_roughness_length_and_green_vegetation_fractionfunction of surface roughness length and green vegetation fractionnonereal(ccpp_constant_one:horizontal_loop_extent)kind_physin
psksurface_dimensionless_exner_functiondimensionless Exner function at the surface interfacenonereal(ccpp_constant_one:horizontal_loop_extent)kind_physin
rbsoilbulk_richardson_number_at_lowest_model_levelbulk Richardson number at the surfacenonereal(ccpp_constant_one:horizontal_loop_extent)kind_physin
zorlsurface_roughness_lengthsurface roughness length in cmcmreal(ccpp_constant_one:horizontal_loop_extent)kind_physin
u10mx_wind_at_10mx component of wind at 10 mm s-1real(ccpp_constant_one:horizontal_loop_extent)kind_physin
v10my_wind_at_10my component of wind at 10 mm s-1real(ccpp_constant_one:horizontal_loop_extent)kind_physin
fmmonin_obukhov_similarity_function_for_momentumMonin-Obukhov similarity function for momentumnonereal(ccpp_constant_one:horizontal_loop_extent)kind_physin
fhmonin_obukhov_similarity_function_for_heatMonin-Obukhov similarity function for heatnonereal(ccpp_constant_one:horizontal_loop_extent)kind_physin
tseasurface_skin_temperaturesurface skin temperatureKreal(ccpp_constant_one:horizontal_loop_extent)kind_physin
heatkinematic_surface_upward_sensible_heat_flux_reduced_by_surface_roughness_and_vegetationkinematic surface upward sensible heat flux reduced by surface roughness and vegetationK m s-1real(ccpp_constant_one:horizontal_loop_extent)kind_physin
evapsurface_upward_specific_humidity_fluxkinematic surface upward latent heat fluxkg kg-1 m s-1real(ccpp_constant_one:horizontal_loop_extent)kind_physin
stresssurface_wind_stresssurface wind stressm2 s-2real(ccpp_constant_one:horizontal_loop_extent)kind_physin
spd1wind_speed_at_lowest_model_layerwind speed at lowest model levelm s-1real(ccpp_constant_one:horizontal_loop_extent)kind_physin
kpblvertical_index_at_top_of_atmosphere_boundary_layerPBL top model level indexindexinteger(ccpp_constant_one:horizontal_loop_extent)out
prsiair_pressure_at_interfaceair pressure at model layer interfacesPareal(ccpp_constant_one:horizontal_loop_extent, ccpp_constant_one:vertical_interface_dimension)kind_physin
delair_pressure_difference_between_midlayerspres(k) - pres(k+1)Pareal(ccpp_constant_one:horizontal_loop_extent, ccpp_constant_one:vertical_layer_dimension)kind_physin
prslair_pressuremean layer pressurePareal(ccpp_constant_one:horizontal_loop_extent, ccpp_constant_one:vertical_layer_dimension)kind_physin
prslkdimensionless_exner_functionExner function at layersnonereal(ccpp_constant_one:horizontal_loop_extent, ccpp_constant_one:vertical_layer_dimension)kind_physin
phiigeopotential_at_interfacegeopotential at model layer interfacesm2 s-2real(ccpp_constant_one:horizontal_loop_extent, ccpp_constant_one:vertical_interface_dimension)kind_physin
philgeopotentialgeopotential at model layer centersm2 s-2real(ccpp_constant_one:horizontal_loop_extent, ccpp_constant_one:vertical_layer_dimension)kind_physin
delttimestep_for_physicstime step for physicssreal()kind_physin
dspheatflag_tke_dissipation_heatingflag for using TKE dissipation heatingflaglogical()in
dusfcinstantaneous_surface_x_momentum_fluxx momentum fluxPareal(ccpp_constant_one:horizontal_loop_extent)kind_physout
dvsfcinstantaneous_surface_y_momentum_fluxy momentum fluxPareal(ccpp_constant_one:horizontal_loop_extent)kind_physout
dtsfcinstantaneous_surface_upward_sensible_heat_fluxsurface upward sensible heat fluxW m-2real(ccpp_constant_one:horizontal_loop_extent)kind_physout
dqsfcinstantaneous_surface_upward_latent_heat_fluxsurface upward latent heat fluxW m-2real(ccpp_constant_one:horizontal_loop_extent)kind_physout
hpblatmosphere_boundary_layer_thicknessPBL thicknessmreal(ccpp_constant_one:horizontal_loop_extent)kind_physout
dktatmosphere_heat_diffusivityatmospheric heat diffusivitym2 s-1real(ccpp_constant_one:horizontal_loop_extent, ccpp_constant_one:vertical_layer_dimension)kind_physout
dkuatmosphere_momentum_diffusivityatmospheric momentum diffusivitym2 s-1real(ccpp_constant_one:horizontal_loop_extent, ccpp_constant_one:vertical_layer_dimension)kind_physout
kinverindex_of_highest_temperature_inversionindex of highest temperature inversionindexinteger(ccpp_constant_one:horizontal_loop_extent)in
xkzm_matmosphere_momentum_diffusivity_due_to_backgroundbackground value of momentum diffusivitym2 s-1real()kind_physin
xkzm_hatmosphere_heat_diffusivity_due_to_backgroundbackground value of heat diffusivitym2 s-1real()kind_physin
xkzm_ssigma_pressure_threshold_at_upper_extent_of_background_diffusionsigma level threshold for background diffusivitynonereal()kind_physin
dspfacmultiplicative_tuning_parameter_for_tke_dissipative_heatingtke dissipative heating factornonereal()kind_physin
bl_upfrupdraft_area_fraction_in_scale_aware_tke_moist_edmf_pbl_schemeupdraft fraction in boundary layer mass flux schemenonereal()kind_physin
bl_dnfrdowndraft_area_fraction_in_scale_aware_tke_moist_edmf_pbl_schemedowndraft fraction in boundary layer mass flux schemenonereal()kind_physin
rlmxmaximum_allowed_mixing_length_in_boundary_layer_mass_flux_schememaximum allowed mixing length in boundary layer mass flux schememreal()kind_physin
elmxmaximum_allowed_dissipation_mixing_length_in_boundary_layer_mass_flux_schememaximum allowed dissipation mixing length in boundary layer mass flux schememreal()kind_physin
sfc_rlmchoice_of_near_surface_mixing_length_in_boundary_layer_mass_flux_schemechoice of near surface mixing length in boundary layer mass flux schemenoneinteger()in
ntqvindex_of_specific_humidity_in_tracer_concentration_arraytracer index for water vapor (specific humidity)indexinteger()in
dtendcumulative_change_of_state_variablesdiagnostic tendencies for state variablesmixedreal(ccpp_constant_one:horizontal_loop_extent, ccpp_constant_one:vertical_layer_dimension, ccpp_constant_one:cumulative_change_of_state_variables_outer_index_max)kind_physin
dtidxcumulative_change_of_state_variables_outer_indexindex of state-variable and process in last dimension of diagnostic tendencies array AKA cumulative_change_indexindexinteger(ccpp_constant_one:number_of_tracers_plus_one_hundred, ccpp_constant_one:number_of_cumulative_change_processes)in
index_of_temperatureindex_of_temperature_in_cumulative_change_indexindex of temperature in first dimension of array cumulative change indexindexinteger()in
index_of_x_windindex_of_x_wind_in_cumulative_change_indexindex of x-wind in first dimension of array cumulative change indexindexinteger()in
index_of_y_windindex_of_y_wind_in_cumulative_change_indexindex of x-wind in first dimension of array cumulative change indexindexinteger()in
index_of_process_pblindex_of_subgrid_scale_vertical_mixing_process_in_cumulative_change_indexindex of subgrid scale vertical mixing process in second dimension of array cumulative change indexindexinteger()in
gen_tendflag_for_generic_tendency_due_to_planetary_boundary_layertrue if GFS_PBL_generic should calculate tendenciesflaglogical()in
ldiag3dflag_for_diagnostics_3dflag for 3d diagnostic fieldsflaglogical()in
errmsgccpp_error_messageerror message for error handling in CCPPnonecharacter()len=*out
errflgccpp_error_codeerror code for error handling in CCPP1integer()out

GFS satmedmfvdifq General Algorithm

satmedmfvdifq_run() computes subgrid vertical turbulence mixing using the scale-aware TKE-based moist eddy-diffusion mass-flux (EDMF) parameterization of Han and Bretherton (2019) [78] .

  1. The local turbulent mixing is represented by an eddy-diffusivity scheme which is a function of a prognostic TKE.
  2. For the convective boundary layer, nonlocal transport by large eddies (mfpbltq.f), is represented using a mass flux approach (Siebesma et al.(2007) [175] ).
  3. A mass-flux approach is also used to represent the stratocumulus-top-induced turbulence (mfscuq.f).

GFS satmedmfvdifq Detailed Algorithm

Compute preliminary variables from input arguments

  • Compute physical height of the layer centers and interfaces from the geopotential height (zi and zl)
  • Initialize variables needed for prognostic cumulus closure
  • Compute horizontal grid size (gdx)
  • Initialize tke value at vertical layer centers and interfaces from tracer (tke and tkeh)
  • Compute reciprocal of \( \Delta z \) (rdzt)
  • Compute reciprocal of pressure (tx1, tx2)
  • Compute minimum turbulent mixing length (rlmnz)
  • Compute background vertical diffusivities for scalars and momentum (xkzo and xkzmo)
  • set background diffusivities with xkzm_h & xkzm_m for gdx >= xkgdx and as a function of horizontal grid size for gdx < xkgdx
    xkzm_hx = xkzm_h * (gdx / xkgdx)
    xkzm_mx = xkzm_m * (gdx / xkgdx)
  • Some output variables and logical flags are initialized
  • Compute \(\theta\)(theta), and \(q_l\)(qlx), \(\theta_e\)(thetae), \(\theta_v\)(thvx), \(\theta_{l,v}\) (thlvx) including ice water
  • Compute an empirical cloud fraction based on Xu and Randall (1996) [193]
  • Compute buoyancy modified by clouds
  • Initialize diffusion coefficients to 0 and calculate the total radiative heating rate (dku, dkt, radx)
  • Compute stable/unstable PBL flag (pblflg) based on the total surface energy flux (false if the total surface energy flux is into the surface)

Calculate the PBL height

The calculation of the boundary layer height follows Troen and Mahrt (1986) [185] section 3. The approach is to find the level in the column where a modified bulk Richardson number exceeds a critical value.

  • Compute critical bulk Richardson number ( \(Rb_{cr}\)) (crb)
  • For the unstable PBL, crb is a constant (0.25)
  • For the stable boundary layer (SBL), \(Rb_{cr}\) varies with the surface Rossby number, \(R_{0}\), as given by Vickers and Mahrt (2004) [188]

    \[ Rb_{cr}=0.16(10^{-7}R_{0})^{-0.18} \]

    \[ R_{0}=\frac{U_{10}}{f_{0}z_{0}} \]

    where \(U_{10}\) is the wind speed at 10m above the ground surface, \(f_0\) the Coriolis parameter, and \(z_{0}\) the surface roughness length. To avoid too much variation, we restrict \(Rb_{cr}\) to vary within the range of 0.15~0.35

Compute \(\frac{\Delta t}{\Delta z}\) , \(u_*\)

  • Compute buoyancy \(\frac{\partial \theta_v}{\partial z}\) (bf) and the wind shear squared (shr2)
  • Given the thermal's properties and the critical Richardson number, a loop is executed to find the first level above the surface (kpblx) where the modified Richardson number is greater than the critical Richardson number, using equation 10a from Troen and Mahrt (1996) [185] (also equation 8 from Hong and Pan (1996) [92]):
  • Once the level is found, some linear interpolation is performed to find the exact height of the boundary layer top (where \(R_{i} > Rb_{cr}\)) and the PBL height (hpbl and kpbl) and the PBL top index are saved.

Compute Monin-Obukhov similarity parameters

  • Calculate the Monin-Obukhov nondimensional stability paramter, commonly referred to as \(\zeta\) using the following equation from Businger et al.(1971) [32] (eqn 28):

    \[ \zeta = Ri_{sfc}\frac{F_m^2}{F_h} = \frac{z}{L} \]

    where \(F_m\) and \(F_h\) are surface Monin-Obukhov stability functions calculated in sfc_diff.f and \(L\) is the Obukhov length.

Calculate the nondimensional gradients of momentum and temperature ( \(\phi_m\) (phim) and \(\phi_h\)(phih)) are calculated using eqns 5 and 6 from Hong and Pan (1996) [92] depending on the surface layer stability:

  • For the unstable and neutral conditions:

    \[ \phi_m=(1-16\frac{0.1h}{L})^{-1/4} \phi_h=(1-16\frac{0.1h}{L})^{-1/2} \]

  • For the stable regime

    \[ \phi_m=\phi_t=(1+5\frac{0.1h}{L}) \]

The \(z/L\) (zol) is used as the stability criterion for the PBL.Currently, strong unstable (convective) PBL for \(z/L < -0.02\) and weakly and moderately unstable PBL for \(0>z/L>-0.02\)

  • Compute the velocity scale \(w_s\) (wscale) (eqn 22 of Han et al. 2019). It is represented by the value scaled at the top of the surface layer:

    \[ w_s=(u_*^3+7\alpha\kappa w_*^3)^{1/3} \]

    where \(u_*\) (ustar) is the surface friction velocity, \(\alpha\) is the ratio of the surface layer height to the PBL height (specified as sfcfrac =0.1), \(\kappa =0.4\) is the von Karman constant, and \(w_*\) is the convective velocity scale defined as eqn23 of Han et al.(2019):

    \[ w_{*}=[(g/T)\overline{(w'\theta_v^{'})}_0h]^{1/3} \]

The counter-gradient terms for temperature and humidity are calculated.

  • Equation 4 of Hong and Pan (1996) [92] and are used to calculate the "scaled virtual temperature excess near the surface" (equation 9 in Hong and Pan (1996) [92]) for use in the mass-flux algorithm.

Determine whether stratocumulus layers exist and compute quantities

  • Starting at the PBL top and going downward, if the level is less than 2.5 km and \(q_l\geq q_{lcr}\) then set kcld = k (find the cloud top index in the PBL. If no cloud water above the threshold is hound, scuflg is set to F.

Starting at the PBL top and going downward, if the level is less than the cloud top, find the level of the minimum radiative heating rate wihin the cloud. If the level of the minimum is the lowest model level or the minimum radiative heating rate is positive, then set scuflg to F.

Compute components for mass flux mixing by large thermals

  • If the PBL is convective, the updraft properties are initialized to be the same as the state variables.

Call mfpbltq(), which is an EDMF parameterization (Siebesma et al.(2007) [175]) to take into account nonlocal transport by large eddies. For details of the mfpbltq subroutine, step into its documentation ::mfpbltq

  • Call mfscuq(), which is a new mass-flux parameterization for stratocumulus-top-induced turbulence mixing. For details of the mfscuq subroutine, step into its documentation ::mfscuq

Compute Prandtl number \f$P_r\f$ (prn) and exchange coefficient varying with height

Compute an asymtotic mixing length

  • Following Bougeault and Lacarrere(1989), the characteristic length scale ( \(l_2\)) (eqn 10 in Han et al.(2019) [78]) is given by:

    \[ l_2=min(l_{up},l_{down}) \]

    and dissipation length scale \(l_d\) is given by:

    \[ l_d=(l_{up}l_{down})^{1/2} \]

    where \(l_{up}\) and \(l_{down}\) are the distances that a parcel having an initial TKE can travel upward and downward before being stopped by buoyancy effects.
  • Compute the surface layer length scale ( \(l_1\)) following Nakanishi (2001) [140] (eqn 9 of Han et al.(2019) [78])
  • If sfc_rlm=1, use zk for elm within surface layer

Compute eddy diffusivities

Compute TKE.

  • Compute a minimum TKE deduced from background diffusivity for momentum.

Compute buoyancy and shear productions of TKE

  • First predict tke due to tke production & dissipation(diss)
  • Compute updraft & downdraft properties for TKE
  • Compute tridiagonal matrix elements for turbulent kinetic energy
  • Call tridit() to solve tridiagonal problem for TKE
  • Recover the tendency of tke

Compute tridiagonal matrix elements for heat and moisture

  • Call tridin() to solve tridiagonal problem for heat and moisture
  • Recover the tendencies of heat and moisture

Add TKE dissipative heating to temperature tendency

Compute tridiagonal matrix elements for momentum

  • Call tridi2() to solve tridiagonal problem for momentum
  • Recover the tendencies of momentum

Save PBL height for diagnostic purpose

References mfpbltq_mod::mfpbltq(), mfscuq_mod::mfscuq(), tridi_mod::tridi2(), tridi_mod::tridin(), and tridi_mod::tridit().

Here is the call graph for this function: