The Mass-Flux shallow convection scheme parameterizes the effect of shallow convection on the environment much like the Simplified Arakawa-Schubert Deep Convection scheme with a few key modifications. Perhaps most importantly, no quasi-equilibrium assumption is necessary since the shallow cloud base mass flux is parameterized from the surface buoyancy flux. Further, there are no convective downdrafts, the entrainment rate is greater than for deep convection, and the shallow convection is limited to not extend over the level where \(p=0.7p_{sfc}\). More...
The Mass-Flux shallow convection scheme parameterizes the effect of shallow convection on the environment much like the Simplified Arakawa-Schubert Deep Convection scheme with a few key modifications. Perhaps most importantly, no quasi-equilibrium assumption is necessary since the shallow cloud base mass flux is parameterized from the surface buoyancy flux. Further, there are no convective downdrafts, the entrainment rate is greater than for deep convection, and the shallow convection is limited to not extend over the level where \(p=0.7p_{sfc}\).
This scheme was designed to replace the previous eddy-diffusivity approach to shallow convection with a mass-flux based approach as it is used for deep convection. Differences between the shallow and deep SAS schemes are presented in Han and Pan (2011) [23] . Like the deep scheme, it uses the working concepts put forth in Arakawa and Schubert (1974) [2] but includes modifications and simplifications from Grell (1993) [21] such as only one cloud type (the deepest possible, up to \(p=0.7p_{sfc}\)), rather than a spectrum based on cloud top heights or assumed entrainment rates, although it assumes no convective downdrafts. It contains many modifications associated with deep scheme as discussed in Han and Pan (2011) [23] , including the calculation of cloud top, a greater CFL-criterion-based maximum cloud base mass flux, and the inclusion of convective overshooting.
This space is reserved for a description of how this scheme uses information from other scheme types and/or how information calculated in this scheme is used in other scheme types.
Files | |
file | shalcnv.f |
Contains the entire SAS shallow convection scheme. | |
subroutine | shalcnv (im, ix, km, jcap, delt, delp, prslp, psp, phil, ql, q1, t1, u1, v1, rn, kbot, ktop, kcnv, islimsk, dot, ncloud, hpbl, heat, evap, ud_mf, dt_mf, cnvw, cnvc) |
This subroutine contains the entirety of the SAS shallow convection scheme. More... | |
subroutine shalcnv | ( | integer | im, |
integer | ix, | ||
integer | km, | ||
integer | jcap, | ||
real(kind=kind_phys) | delt, | ||
real(kind=kind_phys), dimension(ix,km) | delp, | ||
real(kind=kind_phys), dimension(ix,km) | prslp, | ||
real(kind=kind_phys), dimension(im) | psp, | ||
real(kind=kind_phys), dimension(ix,km) | phil, | ||
real(kind=kind_phys), dimension(ix,km,2) | ql, | ||
real(kind=kind_phys), dimension(ix,km) | q1, | ||
real(kind=kind_phys), dimension(ix,km) | t1, | ||
real(kind=kind_phys), dimension(ix,km) | u1, | ||
real(kind=kind_phys), dimension(ix,km) | v1, | ||
real(kind=kind_phys), dimension(im) | rn, | ||
integer, dimension(im) | kbot, | ||
integer, dimension(im) | ktop, | ||
integer, dimension(im) | kcnv, | ||
integer, dimension(im), intent(in) | islimsk, | ||
real(kind=kind_phys), dimension(ix,km) | dot, | ||
integer | ncloud, | ||
real(kind=kind_phys), dimension(im) | hpbl, | ||
real(kind=kind_phys), dimension(im) | heat, | ||
real(kind=kind_phys), dimension(im) | evap, | ||
real(kind=kind_phys), dimension(im,km) | ud_mf, | ||
real(kind=kind_phys), dimension(im,km) | dt_mf, | ||
real(kind=kind_phys), dimension(ix,km) | cnvw, | ||
real(kind=kind_phys), dimension(ix,km) | cnvc | ||
) |
This subroutine contains the entirety of the SAS shallow convection scheme.
This routine follows the Simplified Arakawa-Schubert Deep Convection scheme quite closely, although it can be interpreted as only having the "static" and "feedback" control portions, since the "dynamic" control is not necessary to find the cloud base mass flux. The algorithm is simplified from SAS deep convection by excluding convective downdrafts and being confined to operate below \(p=0.7p_{sfc}\). Also, entrainment is both simpler and stronger in magnitude compared to the deep scheme.
[in] | im | number of used points |
[in] | ix | horizontal dimension |
[in] | km | vertical layer dimension |
[in] | jcap | number of spectral wave trancation |
[in] | delt | physics time step in seconds |
[in] | delp | pressure difference between level k and k+1 (Pa) |
[in] | prslp | mean layer presure (Pa) |
[in] | psp | surface pressure (Pa) |
[in] | phil | layer geopotential ( \(m^s/s^2\)) |
[in,out] | ql | cloud water or ice (kg/kg) |
[in,out] | q1 | updated tracers (kg/kg) |
[in,out] | t1 | updated temperature (K) |
[in,out] | u1 | updated zonal wind ( \(m s^{-1}\)) |
[in,out] | v1 | updated meridional wind ( \(m s^{-1}\)) |
[out] | rn | convective rain (m) |
[out] | kbot | index for cloud base |
[out] | ktop | index for cloud top |
[out] | kcnv | flag to denote deep convection (0=no, 1=yes) |
[in] | islimsk | sea/land/ice mask (=0/1/2) |
[in] | dot | layer mean vertical velocity (Pa/s) |
[in] | ncloud | number of cloud species |
[in] | hpbl | PBL height (m) |
[in] | heat | surface sensible heat flux (K m/s) |
[in] | evap | surface latent heat flux (kg/kg m/s) |
[out] | ud_mf | updraft mass flux multiplied by time step ( \(kg/m^2\)) |
[out] | dt_mf | ud_mf at cloud top ( \(kg/m^2\)) |
[out] | cnvw | convective cloud water (kg/kg) |
[out] | cnvc | convective cloud cover (unitless) |
\[ \overline{w'\theta_v'}=\overline{w'\theta'}+\left(\frac{R_v}{R_d}-1\right)T_0\overline{w'q'} \]
where \(\overline{w'\theta'}\) is the surface sensible heat flux, \(\overline{w'q'}\) is the surface latent heat flux, \(R_v\) is the gas constant for water vapor, \(R_d\) is the gas constant for dry air, and \(T_0\) is a reference temperature.\[ \frac{1}{\eta}\frac{\partial \eta}{\partial z} = \lambda_e - \lambda_d \]
where \(\eta\) is the normalized mass flux, \(\lambda_e\) is the entrainment rate and \(\lambda_d\) is the detrainment rate. The normalized mass flux increases upward below the cloud base and decreases upward above.\[ A_u=\int_{z_0}^{z_t}\frac{g}{c_pT(z)}\frac{\eta}{1 + \gamma}[h(z)-h^*(z)]dz \]
(discretized according to Grell (1993) [21] equation B.10 using B.2 and B.3 of Arakawa and Schubert (1974) [2] and assuming \(\eta=1\)) where \(A_u\) is the updraft cloud work function, \(z_0\) and \(z_t\) are cloud base and cloud top, respectively, \(\gamma = \frac{L}{c_p}\left(\frac{\partial \overline{q_s}}{\partial T}\right)_p\) and other quantities are previously defined.\[ M_c = 0.03\rho w_* \]
where \(M_c\) is the cloud base mass flux, \(\rho\) is the air density, and \(w_*=\left(\frac{g}{T_0}\overline{w'\theta_v'}h\right)^{1/3}\) with \(h\) the PBL height and other quantities have been defined previously.Definition at line 58 of file shalcnv.f.
References physcons::con_cliq, physcons::con_cp, physcons::con_cvap, physcons::con_g, physcons::con_hvap, physcons::con_rd, physcons::con_rv, and physcons::con_t0c.
Referenced by gbphys().