GFS Operational Physics Documentation  Revision: 81451

This module computes cloud related quantities for radiation computations. More...

Detailed Description

This module computes cloud related quantities for radiation computations.

Knowledge of cloud properties and their vertical structure is important for meteorological studies due to their impact on both the Earth's radiation budget and adiabatic heating within the atmosphere. Cloud properties in the US National Oceanic and Atmospheric Administration National Centers for Environmental Prediction Global Forecast System (GFS) include (i) cloud liquid/ice water path; (ii) the fraction of clouds; (iii) effective radius of water/ice droplet:

Version
NCEP-Radiation_clouds v5.1 Nov 2012

This module has three externally accessible subroutines:

and two internally accessable only subroutines:

General Algorithm

  1. Cloud Liquid/Ice Water Path (LWP,IWP)
    We define the fraction of liquid and ice cloud as:
    Fraction of ice cloud (F): \(F=(273.16K-T)/20\)
    LWP = total cloud condensate path X (1-F)
    IWP = total clod condensate path X F
  2. GFS Cloud Fraction
    The cloud fraction in a given grid box of the GFS model is computed using the parameterization scheme of Xu and Randall(1996) [54] :

    \[ \sigma =RH^{k_{1}}\left[1-exp\left(-\frac{k_{2}q_{l}}{\left[\left(1-RH\right)q_{s}\right]^{k_{3}}}\right)\right] \]

    Where \(RH\) is relative humidity, \(q_{l}\) is the cloud condensate, \(q_{s}\) is saturation specific humidity, \(k_{1}(=0.25)\), \(k_{2}(=100)\), \(k_{3}(=0.49)\) are the empirical parameters. The cloud condensate is partitioned into cloud water and ice in radiation based on temperature. Cloud drop effective radius ranges 5-10 microns over land depending on temperature. Ice crystal radius is function of ice water content (Heymsfield and McFarquhar (1996) [25]). Maximum-randomly cloud overlapping is used in both long-wave radiation and short-wave radiation. Convective clouds are not considered in radiation.
  3. The parameterization of effective radius of water/ice droplet ( \(r_{e}\))
    Two methods has been used to parameterize cloud properties in the GFS model. The first method makes use of a diagnostic cloud scheme, in which cloud properties are determined based on model-predicted temperature, pressure, and boundary layer circulation from Harshvardhan et al. (1989) [23] . The diagnostic scheme is now replaced with a prognostic scheme that uses cloud condensate information instead (NCEP Office Note 441).
    For the parameterization of effective radius, \(r_{ew}\), of water droplet, we fix \(r_{ew}\) to a value of \(10\mu m\) over the oceans. Over the land, \(\) is defined as:

    \[ r_{ew} = 5+5\times F \]

    Thus, the effective radius of cloud water droplets will reach to a minimum values of \(5\mu m\) when F=0, and to a maximum value of \(10\mu m\) when the ice fraction is increasing.
    For ice clouds, following Heymsfield and McFarquhar (1996) [25], we have made the effective ice droplet radius to be an empirical function of ice water concentration (IWC) and environmental temperature as:

    \[ r_{ei}=\begin{cases}(1250/9.917)IWC^{0.109} & T <-50^0C \\(1250/9.337)IWC^{0.080} & -50^0C \leq T<-40^0C\\(1250/9.208)IWC^{0.055} & -40^0C\leq T <-30^0C\\(1250/9.387)IWC^{0.031} & -30^0C \leq T\end{cases} \]

    where IWC and IWP satisfy:

    \[ IWP_{\triangle Z}=\int_{\triangle Z} IWCdZ \]

Collaboration diagram for module_radiation_clouds:

Variables

integer, parameter, public module_radiation_clouds::nf_clds = 9
 number of fields in cloud array
 
integer, parameter, public module_radiation_clouds::nk_clds = 3
 number of cloud vertical domains
 
real(kind=kind_phys), dimension(nk_clds+1, 2), save module_radiation_clouds::ptopc
 pressure limits of cloud domain interfaces (low,mid,high) in mb (0.1kPa)
 
real(kind=kind_phys), parameter module_radiation_clouds::climit = 0.001
 
real(kind=kind_phys), parameter module_radiation_clouds::climit2 =0.05
 
real(kind=kind_phys), parameter module_radiation_clouds::ovcst = 1.0 - 1.0e-8
 
real(kind=kind_phys), parameter module_radiation_clouds::reliq_def = 10.0
 default liq radius to 10 micron
 
real(kind=kind_phys), parameter module_radiation_clouds::reice_def = 50.0
 default ice radius to 50 micron
 
real(kind=kind_phys), parameter module_radiation_clouds::rrain_def = 1000.0
 default rain radius to 1000 micron
 
real(kind=kind_phys), parameter module_radiation_clouds::rsnow_def = 250.0
 default snow radius to 250 micron
 
integer, parameter module_radiation_clouds::nbin =100
 rh in one percent interval
 
integer, parameter module_radiation_clouds::nlon =2
 =1,2 for eastern and western hemispheres
 
integer, parameter module_radiation_clouds::nlat =4
 =1,4 for 60n-30n,30n-equ,equ-30s,30s-60s
 
integer, parameter module_radiation_clouds::mcld =4
 =1,4 for bl,low,mid,hi cld type
 
integer, parameter module_radiation_clouds::nseal =2
 =1,2 for land,sea
 
real(kind=kind_phys), parameter module_radiation_clouds::cldssa_def = 0.99
 default cld single scat albedo
 
real(kind=kind_phys), parameter module_radiation_clouds::cldasy_def = 0.84
 default cld asymmetry factor
 
real(kind=kind_phys), dimension(3) module_radiation_clouds::xlabdy
 lat bndry between tuning regions
 
real(kind=kind_phys), dimension(3) module_radiation_clouds::xlobdy
 lon bndry between tuning regions
 
real(kind=kind_phys), parameter module_radiation_clouds::xlim =5.0
 +/- xlim for transition
 
real(kind=kind_phys), parameter module_radiation_clouds::vvcld1 = 0.0003e0
 low cloud vertical velocity adjustment boundaries in mb/sec
 
real(kind=kind_phys), parameter module_radiation_clouds::vvcld2 =-0.0005e0
 low cloud vertical velocity adjustment boundaries in mb/sec
 
real(kind=kind_phys), dimension(nbin, nlon, nlat, mcld, nseal) module_radiation_clouds::rhcl
 tuned relative humidity relation table for diagnostic cloud scheme
 
integer module_radiation_clouds::llyr = 2
 upper limit of boundary layer clouds
 
integer module_radiation_clouds::iovr = 1
 maximum-random cloud overlapping method
 
subroutine, public module_radiation_clouds::cld_init
 This subroutine is an initialization program for cloud-radiation calculations and sets up boundary layer cloud top. More...
 
subroutine, public module_radiation_clouds::progcld1
 This subroutine computes cloud related quantities using zhao/moorthi's prognostic cloud microphysics scheme. More...
 
subroutine, public module_radiation_clouds::progcld2
 This subroutine computes cloud related quantities using ferrier's prognostic cloud microphysics scheme. More...
 
subroutine, public module_radiation_clouds::progcld3
 This subroutine computes cloud related quantities using zhao/moorthi's prognostic cloud microphysics scheme + pdfcld. More...
 
subroutine, public module_radiation_clouds::diagcld1
 This subroutine computes cloud fractions for radiation calculations. More...
 
subroutine module_radiation_clouds::gethml
 This subroutine computes high, mid, low, total, and boundary cloud fractions and cloud top/bottom layer indices for model diagnostic output. The three cloud domain boundaries are defined by ptopc. The cloud overlapping method is defined by control flag 'iovr', which is also used by LW and SW radiation programs. More...
 
subroutine module_radiation_clouds::rhtable
 cld-rh relations obtained from mitchell-hahn procedure.
 

Function/Subroutine Documentation

subroutine, public module_radiation_clouds::cld_init ( )

This subroutine is an initialization program for cloud-radiation calculations and sets up boundary layer cloud top.

Parameters
simodel vertical sigma layer interface
NLAYvertical layer number
meprint control flag

General Algorithm

  1. Call rhtable() to set up tuned relative humidity table.
  2. Compute the top of BL cld (llyr), which is the topmost non cld(low) layer for stratiform (at or above lowest 0.1 of the atmosphere).

Definition at line 332 of file radiation_clouds.f.

References physparam::icldflg, physparam::icmphys, iovr, physparam::iovrlw, physparam::iovrsw, physparam::ivflip, llyr, and rhtable().

Referenced by module_radiation_driver::radinit().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine, public module_radiation_clouds::diagcld1 ( )

This subroutine computes cloud fractions for radiation calculations.

Parameters
plyr(IX,NLAY), model layer mean pressure in mb (100Pa)
plvl(IX,NLP1), model level pressure in mb (100Pa)
tlyr(IX,NLAY), model layer mean temperature in K
rhly(IX,NLAY), layer relative humidity
vvel(IX,NLAY), layer mean vertical velocity in mb/sec
cv(IX), fraction of convective cloud
cvt,cvb(IX), conv cloud top/bottom pressure in mb
xlat(IX), grid latitude in radians, default to pi/2 -> -pi/2 range, otherwise see in-line comment
xlon(IX), grid longitude in radians, ok for both 0->2pi or -pi -> +pi ranges
slmsk(IX), sea/land mask array (sea:0,land:1,sea-ice:2)
IXhorizontal dimention
NLAY,NLP1vertical layer/level dimensions
clouds(IX,NLAY,NF_CLDS), cloud profiles
(:,:,1) - layer total cloud fraction
(:,:,2) - layer cloud optical depth
(:,:,3) - layer cloud single scattering albedo
(:,:,4) - layer cloud asymmetry factor
clds(IX,5), fraction of clouds for low, mid, hi, tot, bl
mtop(IX,3), vertical indices for low, mid, hi cloud tops
mbot(IX,3), vertical indices for low, mid, hi cloud bases

General Algorithm

  1. Get rh-cld relation for this lat.
  2. Linear transition between latitudinal regions.
  3. Get rh-cld relationship for each grid point, interpolating longitudinally between regions if necessary.
  4. Find top pressure for each cloud domain.
  5. Compute stratiform cloud optical depth.
  6. Calculate potential temperature and its lapse rate.
  7. Diagnostic method to find cloud amount cldtot, cldcnv.
    • Find the lowest low cloud top sigma level, computed for each lat cause domain definition changes with latitude.
    • Find the stratosphere cut off layer for high cloud (about 250mb). It is assumed to be above the layerwith dthdp less than -0.25 in the high cloud domain.
    • Put convective cloud into 'cldcnv', no merge at this point.
    • Calculate stratiform cloud and put into array 'cldtot' using the cld-rh relationship from table look-up, where tables obtained using k.mitchell frequency distribution tuning.
    • Compute vertical velocity adjustment on low clouds.
    • Calculate diagnostic cloud optical depth.
  8. Call gethml(), to compute low, mid, high, total, and boundary layer cloud fractions and cloud top/bottom layer indices for low, mid, and high clouds.

Definition at line 1855 of file radiation_clouds.f.

References cldasy_def, cldssa_def, physcons::con_pi, physcons::con_ttp, gethml(), physparam::ivflip, llyr, mcld, nbin, nlon, nseal, ptopc, rhcl, vvcld1, vvcld2, xlabdy, xlim, and xlobdy.

Referenced by module_radiation_driver::grrad().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine module_radiation_clouds::gethml ( )
private

This subroutine computes high, mid, low, total, and boundary cloud fractions and cloud top/bottom layer indices for model diagnostic output. The three cloud domain boundaries are defined by ptopc. The cloud overlapping method is defined by control flag 'iovr', which is also used by LW and SW radiation programs.

Parameters
plyr(IX,NLAY), model layer mean pressure in mb (100Pa)
ptop1(IX,4), pressure limits of cloud domain interfaces (sfc,low,mid,high) in mb (100Pa)
cldtot(IX,NLAY), total or stratiform cloud profile in fraction
cldcnv(IX,NLAY), convective cloud (for diagnostic scheme only)
IXhorizontal dimension
NLAYvertical layer dimensions
clds(IX,5), fraction of clouds for low, mid, hi, tot, bl
mtop(IX,3),vertical indices for low, mid, hi cloud tops
mbot(IX,3),vertical indices for low, mid, hi cloud bases

Detailed Algorithm

  1. Calculate total and BL cloud fractions (maximum-random cloud overlapping is operational).
  2. Calculte high, mid, low cloud fractions and vertical indices of cloud tops/bases.

Definition at line 2472 of file radiation_clouds.f.

References iovr, physparam::ivflip, llyr, and nk_clds.

Referenced by diagcld1(), progcld1(), progcld2(), and progcld3().

Here is the caller graph for this function:

subroutine, public module_radiation_clouds::progcld1 ( )

This subroutine computes cloud related quantities using zhao/moorthi's prognostic cloud microphysics scheme.

Parameters
plyr(IX,NLAY), model layer mean pressure in mb (100Pa)
plvl(IX,NLP1), model level pressure in mb (100Pa)
tlyr(IX,NLAY), model layer mean temperature in K
tvly(IX,NLAY), model layer virtual temperature in K
qlyr(IX,NLAY), layer specific humidity in gm/gm
qstl(IX,NLAY), layer saturate humidity in gm/gm
rhly(IX,NLAY), layer relative humidity \( (=qlyr/qstl) \)
clw(IX,NLAY), layer cloud condensate amount
xlat(IX), grid latitude in radians, default to pi/2 -> -pi/2 range, otherwise see in-line comment
xlon(IX), grid longitude in radians (not used)
slmsk(IX), sea/land mask array (sea:0,land:1,sea-ice:2)
IXhorizontal dimention
NLAY,NLP1vertical layer/level dimensions
clouds(IX,NLAY,NF_CLDS), cloud profiles
(:,:,1) - layer total cloud fraction
(:,:,2) - layer cloud liq water path \((g/m^2)\)
(:,:,3) - mean eff radius for liq cloud (micron)
(:,:,4) - layer cloud ice water path \((g/m^2)\)
(:,:,5) - mean eff radius for ice cloud (micron)
(:,:,6) - layer rain drop water path (not assigned)
(:,:,7) - mean eff radius for rain drop (micron)
(:,:,8) - layer snow flake water path (not assigned)
(:,:,9) - mean eff radius for snow flake (micron)
*** fu's scheme need to be normalized by snow density \( (g/m^3/1.0e6)\)
clds(IX,5), fraction of clouds for low, mid, hi, tot, bl
mtop(IX,3), vertical indices for low, mid, hi cloud tops
mbot(IX,3), vertical indices for low, mid, hi cloud bases

General Algorithm

  1. Find top pressure for each cloud domain for given latitude.
  2. Compute cloud liquid/ice condensate path in \( g/m^2 \) .
  3. Compute effective liquid cloud droplet radius over land.
  4. Calculate layer cloud fraction.
  5. Compute effective ice cloud droplet radius following Heymsfield and McFarquhar (1996) [25].
  6. Call gethml() to compute low,mid,high,total, and boundary layer cloud fractions and clouds top/bottom layer indices for low, mid, and high clouds.

Definition at line 482 of file radiation_clouds.f.

References physcons::con_pi, physcons::con_ttp, gethml(), physparam::ivflip, physparam::lcnorm, physparam::lcrick, nf_clds, ptopc, reice_def, reliq_def, rrain_def, and rsnow_def.

Referenced by module_radiation_driver::grrad().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine, public module_radiation_clouds::progcld2 ( )

This subroutine computes cloud related quantities using ferrier's prognostic cloud microphysics scheme.

Parameters
plyr(IX,NLAY), model layer mean pressure in mb (100Pa)
plvl(IX,NLP1), model level pressure in mb (100Pa)
tlyr(IX,NLAY), model layer mean temperature in K
tvly(IX,NLAY), model layer virtual temperature in K
qlyr(IX,NLAY), layer specific humidity in gm/gm
qstl(IX,NLAY), layer saturate humidity in gm/gm
rhly(IX,NLAY), layer relative humidity (=qlyr/qstl)
clw(IX,NLAY), layer cloud condensate amount
f_ice(IX,NLAY), fraction of layer cloud ice (ferrier micro-phys)
f_rain(IX,NLAY), fraction of layer rain water (ferrier micro-phys)
r_rime(IX,NLAY), mass ratio of total ice to unrimed ice (>=1)
flgmin(IX), minimum large ice fraction
xlat(IX), grid latitude in radians, default to pi/2 -> -pi/2 range, otherwise see in-line comment
xlon(IX), grid longitude in radians (not used)
slmsk(IX), sea/land mask array (sea:0,land:1,sea-ice:2)
IXhorizontal dimention
NLAY,NLP1vertical layer/level dimensions
clouds(IX,NLAY,NF_CLDS), cloud profiles
(:,:,1) - layer total cloud fraction
(:,:,2) - layer cloud liq water path \((g/m^2)\)
(:,:,3) - mean eff radius for liq cloud (micron)
(:,:,4) - layer cloud ice water path \((g/m^2)\)
(:,:,5) - mean eff radius for ice cloud (micron)
(:,:,6) - layer rain drop water path \((g/m^2)\)
(:,:,7) - mean eff radius for rain drop (micron)
(:,:,8) - layer snow flake water path \((g/m^2)\)
(:,:,9) - mean eff radius for snow flake (micron)
*** fu's scheme need to be normalized by snow density (g/m**3/1.0e6)
clds(IX,5), fraction of clouds for low, mid, hi, tot, bl
mtop(IX,3), vertical indices for low, mid, hi cloud tops
mbot(IX,3), vertical indices for low, mid, hi cloud bases

General Algorithm

  1. Find top pressure (ptopc) for each cloud domain for given latitude.
  2. Seperate cloud condensate into liquid, ice, and rain types, and save the liquid+ice condensate in array clw2 for later calculation of cloud fraction.
  3. Call module_microphysics::rsipath2(), in Ferrier's scheme, to compute layer's cloud liquid, ice, rain, and snow water condensate path and the partical effective radius for liquid droplet, rain drop, and snow flake.
  4. Calculate layer cloud fraction.
  5. Calculate effective ice cloud droplet radius.
  6. Call gethml(), to compute low, mid, high, total, and boundary layer cloud fractions and clouds top/bottom layer indices for low, mid, and high clouds.

Definition at line 934 of file radiation_clouds.f.

References physcons::con_g, physcons::con_pi, physcons::con_rd, physcons::con_t0c, physcons::con_ttp, gethml(), physparam::ivflip, physparam::lcnorm, physparam::lcrick, physparam::lnoprec, ptopc, reice_def, reliq_def, rrain_def, and rsnow_def.

Referenced by module_radiation_driver::grrad().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine, public module_radiation_clouds::progcld3 ( )

This subroutine computes cloud related quantities using zhao/moorthi's prognostic cloud microphysics scheme + pdfcld.

Parameters
plyr(ix,nlay), model layer mean pressure in mb (100pa)
plvl(ix,nlp1), model level pressure in mb (100pa)
tlyr(ix,nlay), model layer mean temperature in K
tvly(ix,nlay), model layer virtual temperature in K
qlyr(ix,nlay), layer specific humidity in gm/gm
qstl(ix,nlay), layer saturate humidity in gm/gm
rhly(ix,nlay), layer relative humidity (=qlyr/qstl)
clw(ix,nlay), layer cloud condensate amount
xlat(ix), grid latitude in radians, default to pi/2 -> -pi/2 range, otherwise see in-line comment
xlon(ix), grid longitude in radians (not used)
slmsk(ix), sea/land mask array (sea:0,land:1,sea-ice:2)
ixhorizontal dimention
nlay,nlp1vertical layer/level dimensions
deltaq(ix,nlay), half total water distribution width
supsupersaturation
kdt
meprint control flag
clouds(ix,nlay,nf_clds), cloud profiles
(:,:,1) - layer total cloud fraction
(:,:,2) - layer cloud liq water path (g/m**2)
(:,:,3) - mean eff radius for liq cloud (micron)
(:,:,4) - layer cloud ice water path (g/m**2)
(:,:,5) - mean eff radius for ice cloud (micron)
(:,:,6) - layer rain drop water path not assigned
(:,:,7) - mean eff radius for rain drop (micron)
(:,:,8) - layer snow flake water path not assigned
(:,:,9) - mean eff radius for snow flake(micron)
clds(ix,5), fraction of clouds for low, mid, hi, tot, bl
mtop(ix,3), vertical indices for low, mid, hi cloud tops
mbot(ix,3), vertical indices for low, mid, hi cloud bases

General Algorithm

  1. Find top pressure (ptopc) for each cloud domain for given latitude.
  2. Calculate liquid/ice condensate path in \( g/m^2 \)
  3. Calculate effective liquid cloud droplet radius over land.
  4. Calculate layer cloud fraction.
  5. Calculate effective ice cloud droplet radius following Heymsfield and McFarquhar (1996) [25].
  6. Call gethml() to compute low,mid,high,total, and boundary layer cloud fractions and clouds top/bottom layer indices for low, mid, and high clouds.

Definition at line 1433 of file radiation_clouds.f.

References physcons::con_pi, physcons::con_thgni, physcons::con_ttp, gethml(), physparam::ivflip, physparam::lcnorm, physparam::lcrick, nf_clds, ptopc, reice_def, reliq_def, rrain_def, and rsnow_def.

Referenced by module_radiation_driver::grrad().

Here is the call graph for this function:

Here is the caller graph for this function: