CCPP Scientific Documentation
v5.0.0
module_radiation_clouds Module Reference

This module computes cloud related quantities for radiation computations.

Functions/Subroutines

subroutine, public cld_init (si, NLAY, imp_physics, me)
 This subroutine is an initialization program for cloud-radiation calculations and sets up boundary layer cloud top. More...
 
subroutine, public progcld1 (plyr, plvl, tlyr, tvly, qlyr, qstl, rhly, clw, xlat, xlon, slmsk, dz, delp, IX, NLAY, NLP1, uni_cld, lmfshal, lmfdeep2, cldcov, effrl, effri, effrr, effrs, effr_in, clouds, clds, mtop, mbot, de_lgth )
 This subroutine computes cloud related quantities using zhao/moorthi's prognostic cloud microphysics scheme. More...
 
subroutine, public progcld2 (plyr, plvl, tlyr, tvly, qlyr, qstl, rhly, clw, xlat, xlon, slmsk, dz, delp, f_ice, f_rain, r_rime, flgmin, IX, NLAY, NLP1, lmfshal, lmfdeep2, clouds, clds, mtop, mbot, de_lgth )
 This subroutine computes cloud related quantities using Ferrier's prognostic cloud microphysics scheme. More...
 
subroutine, public progcld3 (plyr, plvl, tlyr, tvly, qlyr, qstl, rhly, clw, cnvw, cnvc, xlat, xlon, slmsk, dz, delp, ix, nlay, nlp1, deltaq, sup, kdt, me, clouds, clds, mtop, mbot, de_lgth )
 This subroutine computes cloud related quantities using zhao/moorthi's prognostic cloud microphysics scheme + pdfcld. More...
 
subroutine, public progcld4 (plyr, plvl, tlyr, tvly, qlyr, qstl, rhly, clw, cnvw, cnvc, xlat, xlon, slmsk, cldtot, dz, delp, IX, NLAY, NLP1, clouds, clds, mtop, mbot, de_lgth )
 This subroutine computes cloud related quantities using GFDL Lin MP prognostic cloud microphysics scheme. More...
 
subroutine, public progcld4o (plyr, plvl, tlyr, tvly, qlyr, qstl, rhly, clw, xlat, xlon, slmsk, dz, delp, ntrac, ntcw, ntiw, ntrw, ntsw, ntgl, ntclamt, IX, NLAY, NLP1, clouds, clds, mtop, mbot, de_lgth )
 This subroutine computes cloud related quantities using GFDL Lin MP prognostic cloud microphysics scheme. Moist species from MP are fed into the corresponding arrays for calculation of cloud fractions. More...
 
subroutine, public progcld5 (plyr, plvl, tlyr, qlyr, qstl, rhly, clw, xlat, xlon, slmsk, dz, delp, ntrac, ntcw, ntiw, ntrw, ntsw, ntgl, IX, NLAY, NLP1, uni_cld, lmfshal, lmfdeep2, cldcov, re_cloud, re_ice, re_snow, clouds, clds, mtop, mbot, de_lgth )
 This subroutine computes cloud related quantities using Thompson/WSM6 cloud microphysics scheme. More...
 
subroutine, public progclduni (plyr, plvl, tlyr, tvly, ccnd, ncnd, xlat, xlon, slmsk, dz, delp, IX, NLAY, NLP1, cldtot, effrl, effri, effrr, effrs, effr_in, clouds, clds, mtop, mbot, de_lgth )
 This subroutine computes cloud related quantities using for unified cloud microphysics scheme. More...
 
subroutine, public gethml (plyr, ptop1, cldtot, cldcnv, dz, de_lgth, IX, NLAY, clds, mtop, mbot )
 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, public get_alpha_dcorr (nCol, nLev, lat, con_pi, deltaZ, de_lgth, cloud_overlap_param)
 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, public get_alpha_exp (nlon, nlay, dzlay, iovrlp, latdeg, juldat, yearlen, cldf, alpha)
 This program derives the exponential transition, alpha, from maximum to random overlap needed to define the fractional cloud vertical correlation for the exponential (EXP, iovrlp=4) or the exponential-random (ER, iovrlp=5) cloud overlap options for RRTMGP. For exponential, the transition from maximum to random with distance through model layers occurs without regard to the configuration of clear and cloudy layers. For the ER method, each block of adjacent cloudy layers is treated with a separate transition from maximum to random, and blocks of cloudy layers separated by one or more clear layers are correlated randomly. /param nlon : number of model longitude points /param nlay : vertical layer dimension /param dzlay(nlon,nlay) : distance between the center of model layers /param iovrlp : cloud overlap method : 0 = random : 1 = maximum-random : 2 = maximum : 3 = decorrelation (NOAA/Hou) : 4 = exponential (AER) : 5 = exponential-random (AER) /param latdeg(nlon) : latitude (in degrees 90 -> -90) /param juldat : day of the year (fractional julian day) /param yearlen : current length of the year (365/366 days) /param cldf(nlon,nlay) : cloud fraction /param idcor : decorrelation length method : 0 = constant value (AER; decorr_con) : 1 = latitude and day of year varying value (AER; Oreopoulos, et al., 2012) /param decorr_con : decorrelation length constant. More...
 

Variables

character(40), parameter vtagcld ='NCEP-Radiation_clouds v5.1 Nov 2012 '
 
real(kind=kind_phys), parameter gfac =1.0e5/con_g
 
real(kind=kind_phys), parameter gord =con_g/con_rd
 
integer, parameter, public nf_clds = 9
 number of fields in cloud array More...
 
integer, parameter, public nk_clds = 3
 number of cloud vertical domains More...
 
real(kind=kind_phys), dimension(nk_clds+1, 2), save ptopc
 pressure limits of cloud domain interfaces (low, mid, high) in mb (0.1kPa) More...
 
real(kind=kind_phys), parameter climit = 0.001
 
real(kind=kind_phys), parameter climit2 =0.05
 
real(kind=kind_phys), parameter ovcst = 1.0 - 1.0e-8
 
real(kind=kind_phys), parameter reliq_def = 10.0
 default liq radius to 10 micron More...
 
real(kind=kind_phys), parameter reice_def = 50.0
 default ice radius to 50 micron More...
 
real(kind=kind_phys), parameter rrain_def = 1000.0
 default rain radius to 1000 micron More...
 
real(kind=kind_phys), parameter rsnow_def = 250.0
 default snow radius to 250 micron More...
 
real(kind=kind_phys), parameter cldssa_def = 0.99
 default cld single scat albedo More...
 
real(kind=kind_phys), parameter cldasy_def = 0.84
 default cld asymmetry factor More...
 
integer llyr = 2
 upper limit of boundary layer clouds More...
 
integer iovr = 1
 maximum-random cloud overlapping method More...