GFS Operational Physics Documentation  gsm/branches/DTC/phys-doc-all phys-doc-all R82971

This module includes NCEP's modifications of the rrtmg-sw radiation code from AER. More...

Detailed Description

This module includes NCEP's modifications of the rrtmg-sw radiation code from AER.

The SW radiation model in the current NOAA Environmental Modeling System (NEMS) was adapted from the RRTM radiation model developed by AER Inc. (Clough et al., 2005 [14]; Mlawer et al., 1997 [45]). It contains 14 spectral bands spanning a spectral wavenumber range of \(50000-820 cm^{-1}\) (corresponding to a wavelength range \(0.2-12.2\mu m\)), each spectral band focuses on a specific set of atmospheric absorbing species as shown in Table 1. To achieve great computation efficiency while at the same time to maintain a high degree of accuracy, the RRTM radiation model employs a corrected-k distribution method (i.e. mapping the highly spectral changing absorption coefficient, k, into a monotonic and smooth varying cumulative probability function, g). In the RRTM-SW, there are 16 unevenly distributed g points for each of the 14 bands for a total of 224 g points. The GCM version of the code (RRTMG-SW) uses a reduced number (various between 2 to 16) of g points for each of the bands that totals to 112 instead of the full set of 224. To get high quality for the scheme, many advanced techniques are used in RRTM such as carefully selecting the band structure to handle various major (key-species) and minor absorbers; deriving a binary parameter for a paired key molecular species in the same domain; and using two pressure regions (dividing level is at about 96mb) for optimal treatment of various species, etc.

Table 1. RRTMG-SW spectral bands and the corresponding absorbing species

Band #Wavenumber Range Lower Atm (Key)Lower Atm (Minor)Mid/Up Atm (Key)Mid/Up Atm (Minor)
16 2600-3250 H2O,CH4 CH4
17 3250-4000 H2O,CO2 H2O,CO2
18 4000-4650 H2O,CH4 CH4
19 4650-5150 H2O,CO2 CO2
20 5150-6150 H2O CH4 H2O CH4
21 6150-7700 H2O,CO2 H2O,CO2
22 7700-8050 H2O,O2 O2
23 8050-12850 H2O
24 12850-16000 H2O,O2 O3 O2 O3
25 16000-22650 H2O O3 O3
26 22650-29000
27 29000-38000 O3 O3
28 38000-50000 O3,O2 O3,O2
29 820-2600 H2O CO2 CO2 H2O


scattering due to clouds greatly complicate the SW radiative transfer computations. To balance the trade-off between computation and speed, RRTMG-SW uses a two-stream approximation method with a delta-function adjustment. Several variations of the delta-two method are included in the radiation transfer code; each holds its own strength and shortcomings (King and Harshvadhan, 1986 [37] ; \(R\ddot{a}is\ddot{a}nen\),2002 [51] ; Barker et al., 2015 [4]). The default (the same in operation runs) selection (iswmode=2) activates the Practical Improved Flux Method (PIFM) by Zdunkowski et al.(1980) [62] . In dealing with a column of cloudy atmosphere, two approaches are included in the RRTMG-SW. One is the commonly used treatment that sees each of the cloud contaminated layers as independent, partially and uniformly filled slabs. Cloud inhomogeneity within and the nature coherence among adjacent cloud layers are largely ignored to reduce the overwhelm complexities associated with scattering process. The effective layer reflectance and transmittance are weighted mean according to cloud fraction. The approach may overestimate cloud effect, especially for multi-layered cloud system associated with deep convection. In NEMS radiation code, to mitigate this shortcoming without increase computation cost, the cloud contaminated column is divided into two parts based on the column's total cloud coverage (a maximum-random overlapping is used in the operational models) to form a cloud free part and an overcast part. Layered clouds are then normalized by the total cloud amount before going through radiative transfer calculations. Fluxes from the cloud-free part and cloudy part are combined together to obtain the final result.
On the other hand, the Monte-Carlo Independent Column Approximation (McICA) (Pincus et al.,2003 [49] ; \(R\ddot{a}is\ddot{a}nen\) and Barker, 2004 [50]), provides a simple and effective way to solve cloud overlapping issue without increasing computational burden. The method is based on the concept of an ICA scheme that divides each grid column into a large number of sub-columns, and statistically redistributes layered clouds (under an assumed overlapping condition, such as the maximum-random method) into the sub-columns (i.e. at any layer it will be either clear or overcast). Thus the grid domain averaged flux under ICA scheme can be expressed as:

\[ \overline{F}=\frac{1}{N}\sum_{n=1}^N F_{n} =\frac{1}{N}\sum_{n=1}^N\sum_{k=1}^K F_{n,k} \]

Where \(N\) is the number of total sub-columns, and \(K\) is the number of spectral terms in integration. \(F_{n}\) is flux obtained in the \(n^{th}\) sub-column, that is the summation of total of \(K\) spectral corresponding fluxes, \(F_{n,k}\) . The double integrations (summations) make ICA impractical for GCM applications. The McICA method is to divide a model grid into \(K\) sub-columns and randomly to pair a sub-column's cloud profile with one of the radiative spectral intervals (e.g. the g-point in RRTM). The double summations will then be reduced to only one:

\[ \overline{F}=\frac{1}{N}\sum_{n=1}^N\sum_{k=1}^K F_{n,k} \approx\overline{F}=\sum_{k=1}^K F_{S_{k},k} \]

The RRTM-SW package includes three files:

Author
Eli J. Mlawer, emlaw.nosp@m.er@a.nosp@m.er.co.nosp@m.m
Jennifer S. Delamere, jdela.nosp@m.mer@.nosp@m.aer.c.nosp@m.om
Michael J. Iacono, miaco.nosp@m.no@a.nosp@m.er.co.nosp@m.m
Shepard A. Clough
Version
NCEP SW v5.1 Nov 2012 -RRTMG-SW v3.8

The authors wish to acknowledge the contributions of the following people: Steven J. Taubman, Karen Cady-Pereira, Patrick D. Brown, Ronald E. Farren, Luke Chen, Robert Bergstrom.

Collaboration diagram for module_radsw_main:

Modules

 module_radsw_kgbnn
 

Modules

module  module_radsw_parameters
 This module is for specifying the band structures and program parameters used by the RRTMG-SW scheme.
 

Functions/Subroutines

subroutine, public module_radsw_main::rswinit (me)
 This subroutine initializes non-varying module variables, conversion factors, and look-up tables. More...
 
subroutine module_radsw_main::mcica_subcol (cldf, nlay, ipseed, lcloudy )
 This subroutine computes the sub-colum cloud profile flag array. More...
 
subroutine module_radsw_main::setcoef (pavel, tavel, h2ovmr, nlay, nlp1, laytrop, jp, jt, jt1, fac00, fac01, fac10, fac11, selffac, selffrac, indself, forfac, forfrac, indfor )
 This subroutine computes various coefficients needed in radiative transfer calculation. More...
 
subroutine module_radsw_main::spcvrtm (ssolar, cosz, sntz, albbm, albdf, sfluxzen, cldfmc, cf1, cf0, taug, taur, tauae, ssaae, asyae, taucw, ssacw, asycw, nlay, nlp1, fxupc, fxdnc, fxup0, fxdn0, ftoauc, ftoau0, ftoadc, fsfcuc, fsfcu0, fsfcdc, fsfcd0, sfbmc, sfdfc, sfbm0, sfdf0, suvbfc, suvbf0 )
 This subroutine computes the shortwave radiative fluxes using two-stream method of h. barder and mcica,the monte-carlo independent column approximation, for the representation of sub-grid cloud variability (i.e. cloud overlap). More...
 

constant values

real(kind=kind_phys), parameter module_radsw_main::eps = 1.0e-6
 
real(kind=kind_phys), parameter module_radsw_main::oneminus = 1.0 - eps
 
real(kind=kind_phys), parameter module_radsw_main::bpade = 1.0/0.278
 pade approx constant
 
real(kind=kind_phys), parameter module_radsw_main::stpfac = 296.0/1013.0
 
real(kind=kind_phys), parameter module_radsw_main::ftiny = 1.0e-12
 
real(kind=kind_phys), parameter module_radsw_main::s0 = 1368.22
 internal solar constant
 
real(kind=kind_phys), parameter module_radsw_main::f_zero = 0.0
 
real(kind=kind_phys), parameter module_radsw_main::f_one = 1.0
 

atomic weights for conversion from mass to volume mixing ratios

real(kind=kind_phys), parameter module_radsw_main::amdw = con_amd/con_amw
 
real(kind=kind_phys), parameter module_radsw_main::amdo3 = con_amd/con_amo3
 

band indices

integer, dimension(nblow:nbhgh) module_radsw_main::nspa
 
integer, dimension(nblow:nbhgh) module_radsw_main::nspb
 
integer, dimension(nblow:nbhgh) module_radsw_main::idxsfc
 band index for sfc flux
 
integer, dimension(nblow:nbhgh) module_radsw_main::idxebc
 band index for cld prop
 
integer, parameter module_radsw_main::nuvb = 27
 uv-b band index
 

logical flags for optional output fields

logical module_radsw_main::lhswb = .false.
 
logical module_radsw_main::lhsw0 = .false.
 
logical module_radsw_main::lflxprf = .false.
 
logical module_radsw_main::lfdncmp = .false.
 
real(kind=kind_phys), dimension(0:ntbmx) module_radsw_main::exp_tbl
 those data will be set up only once by "rswinit"
 
real(kind=kind_phys) module_radsw_main::heatfac
 the factor for heating rates (in k/day, or k/sec set by subroutine 'rswinit')
 
integer, parameter module_radsw_main::ipsdsw0 = 1
 initial permutation seed used for sub-column cloud scheme
 
subroutine, public module_radsw_main::swrad (plyr, plvl, tlyr, tlvl, qlyr, olyr, gasvmr, clouds, icseed, aerosols, sfcalb, cosz, solcon, NDAY, idxday, npts, nlay, nlp1, lprnt, hswc, topflx, sfcflx, HSW0, HSWB, FLXPRF, FDNCMP )
 This subroutine is the main SW radiation routine. More...
 
subroutine module_radsw_main::cldprop (cfrac, cliqp, reliq, cicep, reice, cdat1, cdat2, cdat3, cdat4, cf1, nlay, ipseed, taucw, ssacw, asycw, cldfrc, cldfmc )
 This subroutine computes the cloud optical properties for each cloudy layer and g-point interval. More...
 
subroutine module_radsw_main::spcvrtc (ssolar, cosz, sntz, albbm, albdf, sfluxzen, cldfrc, cf1, cf0, taug, taur, tauae, ssaae, asyae, taucw, ssacw, asycw, nlay, nlp1, fxupc, fxdnc, fxup0, fxdn0, ftoauc, ftoau0, ftoadc, fsfcuc, fsfcu0, fsfcdc, fsfcd0, sfbmc, sfdfc, sfbm0, sfdf0, suvbfc, suvbf0 )
 This subroutine computes the shortwave radiative fluxes using two-stream method. More...
 
subroutine module_radsw_main::vrtqdr (zrefb, zrefd, ztrab, ztrad, zldbt, ztdbt, NLAY, NLP1, zfu, zfd )
 This subroutine is called by spcvrtc() and spcvrtm(), and computes the upward and downward radiation fluxes. More...
 
subroutine module_radsw_main::taumol (colamt, colmol, fac00, fac01, fac10, fac11, jp, jt, jt1, laytrop, forfac, forfrac, indfor, selffac, selffrac, indself, nlay, sfluxzen, taug, taur )
 This subroutine calculates optical depths for gaseous absorption and rayleigh scattering
subroutine called taumol## (## = 16-29) More...
 
subroutine taumol16
 The subroutine computes the optical depth in band 16: 2600-3250 cm-1 (low - h2o,ch4; high - ch4)
 

Function/Subroutine Documentation

subroutine module_radsw_main::cldprop ( real (kind=kind_phys), dimension(nlay), intent(in)  cfrac,
real (kind=kind_phys), dimension(nlay), intent(in)  cliqp,
real (kind=kind_phys), dimension(nlay), intent(in)  reliq,
real (kind=kind_phys), dimension(nlay), intent(in)  cicep,
real (kind=kind_phys), dimension(nlay), intent(in)  reice,
real (kind=kind_phys), dimension(nlay), intent(in)  cdat1,
real (kind=kind_phys), dimension(nlay), intent(in)  cdat2,
real (kind=kind_phys), dimension(nlay), intent(in)  cdat3,
real (kind=kind_phys), dimension(nlay), intent(in)  cdat4,
real (kind=kind_phys), intent(in)  cf1,
integer, intent(in)  nlay,
integer, intent(in)  ipseed,
real (kind=kind_phys), dimension(nlay,nbdsw), intent(out)  taucw,
real (kind=kind_phys), dimension(nlay,nbdsw), intent(out)  ssacw,
real (kind=kind_phys), dimension(nlay,nbdsw), intent(out)  asycw,
real (kind=kind_phys), dimension(nlay), intent(out)  cldfrc,
real (kind=kind_phys), dimension(nlay,ngptsw), intent(out)  cldfmc 
)
private

This subroutine computes the cloud optical properties for each cloudy layer and g-point interval.

Parameters
cfraclayer cloud fraction
for physparam::iswcliq > 0 (prognostic cloud scheme) - - -
cliqplayer in-cloud liq water path ( \(g/m^2\))
reliqmean eff radius for liq cloud (micron)
ciceplayer in-cloud ice water path ( \(g/m^2\))
reicemean eff radius for ice cloud (micron)
cdat1layer rain drop water path ( \(g/m^2\))
cdat2effective radius for rain drop (micron)
cdat3layer snow flake water path( \(g/m^2\))
cdat4mean eff radius for snow flake(micron)
for physparam::iswcliq = 0 (diagnostic cloud scheme) - - -
cliqpnot used
cicepnot used
reliqnot used
reicenot used
cdat1layer cloud optical depth
cdat2layer cloud single scattering albedo
cdat3layer cloud asymmetry factor
cdat4optional use
cf1effective total cloud cover at surface
nlayvertical layer number
ipseedpermutation seed for generating random numbers (isubcsw>0)
taucwcloud optical depth, w/o delta scaled
ssacwweighted cloud single scattering albedo (ssa = ssacw / taucw)
asycwweighted cloud asymmetry factor (asy = asycw / ssacw)
cldfrccloud fraction of grid mean value
cldfmccloud fraction for each sub-column

General Algorithm

  1. Compute cloud radiative properties for a cloudy column.
    • Compute optical properties for rain and snow.
      For rain: tauran/ssaran/asyran
      For snow: tausnw/ssasnw/asysnw
    • Calculation of absorption coefficients due to water clouds
      For water clouds: tauliq/ssaliq/asyliq
    • Calculation of absorption coefficients due to ice clouds
      For ice clouds: tauice/ssaice/asyice
    • For Prognostic cloud scheme: sum up the cloud optical property:
      \( taucw=tauliq+tauice+tauran+tausnw \)
      \( ssacw=ssaliq+ssaice+ssaran+ssasnw \)
      \( asycw=asyliq+asyice+asyran+asysnw \)
  2. if physparam::isubcsw > 0, call mcica_subcol() to distribute cloud properties to each g-point.

Definition at line 1564 of file radsw_main.f.

References module_radsw_cldprtb::a0r, module_radsw_cldprtb::a0s, module_radsw_cldprtb::abari, module_radsw_cldprtb::asyice2, module_radsw_cldprtb::asyice3, module_radsw_cldprtb::asyliq1, module_radsw_cldprtb::b0r, module_radsw_cldprtb::b0s, module_radsw_cldprtb::bbari, module_radsw_cldprtb::c0r, module_radsw_cldprtb::c0s, module_radsw_cldprtb::cbari, module_radsw_cldprtb::dbari, module_radsw_cldprtb::ebari, module_radsw_cldprtb::extice2, module_radsw_cldprtb::extice3, module_radsw_cldprtb::extliq1, module_radsw_cldprtb::fbari, idxebc, physparam::isubcsw, physparam::iswcice, physparam::iswcliq, mcica_subcol(), module_radsw_parameters::nbhgh, module_radsw_parameters::nblow, module_radsw_parameters::ngptsw, module_radsw_cldprtb::ssaice2, module_radsw_cldprtb::ssaice3, and module_radsw_cldprtb::ssaliq1.

Referenced by swrad().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine module_radsw_main::mcica_subcol ( real (kind=kind_phys), dimension(nlay), intent(in)  cldf,
integer, intent(in)  nlay,
integer, intent(in)  ipseed,
logical, dimension(nlay,ngptsw), intent(out)  lcloudy 
)
private

This subroutine computes the sub-colum cloud profile flag array.

Parameters
cldflayer cloud fraction
nlaynumber of model vertical layers
ipseedpermute seed for random num generator
lcloudysub-colum cloud profile flag array

Definition at line 1924 of file radsw_main.f.

References physparam::iovrsw, and module_radsw_parameters::ngptsw.

Referenced by cldprop().

Here is the caller graph for this function:

subroutine, public module_radsw_main::rswinit ( integer, intent(in)  me)

This subroutine initializes non-varying module variables, conversion factors, and look-up tables.

Parameters
meprint control for parallel process

Definition at line 1375 of file radsw_main.f.

References bpade, physcons::con_cp, physcons::con_g, exp_tbl, heatfac, physparam::icldflg, physparam::iovrsw, physparam::isubcsw, physparam::iswcliq, physparam::iswmode, physparam::iswrate, physparam::iswrgas, and module_radsw_parameters::ntbmx.

Referenced by module_radiation_driver::radinit().

Here is the caller graph for this function:

subroutine module_radsw_main::setcoef ( real (kind=kind_phys), dimension(:), intent(in)  pavel,
real (kind=kind_phys), dimension(:), intent(in)  tavel,
real (kind=kind_phys), dimension(:), intent(in)  h2ovmr,
integer, intent(in)  nlay,
integer, intent(in)  nlp1,
integer, intent(out)  laytrop,
integer, dimension(nlay), intent(out)  jp,
integer, dimension(nlay), intent(out)  jt,
integer, dimension(nlay), intent(out)  jt1,
real (kind=kind_phys), dimension(nlay), intent(out)  fac00,
real (kind=kind_phys), dimension(nlay), intent(out)  fac01,
real (kind=kind_phys), dimension(nlay), intent(out)  fac10,
real (kind=kind_phys), dimension(nlay), intent(out)  fac11,
real (kind=kind_phys), dimension(nlay), intent(out)  selffac,
real (kind=kind_phys), dimension(nlay), intent(out)  selffrac,
integer, dimension(nlay), intent(out)  indself,
real (kind=kind_phys), dimension(nlay), intent(out)  forfac,
real (kind=kind_phys), dimension(nlay), intent(out)  forfrac,
integer, dimension(nlay), intent(out)  indfor 
)
private

This subroutine computes various coefficients needed in radiative transfer calculation.

Parameters
pavellayer pressure (mb)
tavellayer temperature (k)
h2ovmrlayer w.v. volumn mixing ratio (kg/kg)
nlaytotal number of vertical layers
nlp1total number of vertical levels
laytroptropopause layer index (unitless)
jpindices of lower reference pressure
jt,jt1indices of lower reference temperatures at levels of jp and jp+1
facijfactors mltiply the reference ks,i,j=0/1 for lower/higher of the 2 appropriate temperature and altitudes.
selffacscale factor for w. v. self-continuum equals (w.v. density)/(atmospheric density at 296k and 1013 mb)
seffracfactor for temperature interpolation of reference w.v. self-continuum data
indselfindex of lower ref temp for selffac
forfacscale factor for w. v. foreign-continuum
forfracfactor for temperature interpolation of reference w.v. foreign-continuum data
indforindex of lower ref temp for forfac

Definition at line 2107 of file radsw_main.f.

References module_radsw_ref::preflog, and module_radsw_ref::tref.

Referenced by swrad().

Here is the caller graph for this function:

subroutine module_radsw_main::spcvrtc ( real (kind=kind_phys), intent(in)  ssolar,
real (kind=kind_phys), intent(in)  cosz,
real (kind=kind_phys), intent(in)  sntz,
real (kind=kind_phys), dimension(2), intent(in)  albbm,
real (kind=kind_phys), dimension(2), intent(in)  albdf,
real (kind=kind_phys), dimension(ngptsw), intent(in)  sfluxzen,
real (kind=kind_phys), dimension(nlay), intent(in)  cldfrc,
real (kind=kind_phys), intent(in)  cf1,
real (kind=kind_phys), intent(in)  cf0,
real (kind=kind_phys), dimension(nlay,ngptsw), intent(in)  taug,
real (kind=kind_phys), dimension(nlay,ngptsw), intent(in)  taur,
real (kind=kind_phys), dimension(nlay,nbdsw), intent(in)  tauae,
real (kind=kind_phys), dimension(nlay,nbdsw), intent(in)  ssaae,
real (kind=kind_phys), dimension(nlay,nbdsw), intent(in)  asyae,
real (kind=kind_phys), dimension(nlay,nbdsw), intent(in)  taucw,
real (kind=kind_phys), dimension(nlay,nbdsw), intent(in)  ssacw,
real (kind=kind_phys), dimension(nlay,nbdsw), intent(in)  asycw,
integer, intent(in)  nlay,
integer, intent(in)  nlp1,
real (kind=kind_phys), dimension(nlp1,nbdsw), intent(out)  fxupc,
real (kind=kind_phys), dimension(nlp1,nbdsw), intent(out)  fxdnc,
real (kind=kind_phys), dimension(nlp1,nbdsw), intent(out)  fxup0,
real (kind=kind_phys), dimension(nlp1,nbdsw), intent(out)  fxdn0,
real (kind=kind_phys), intent(out)  ftoauc,
real (kind=kind_phys), intent(out)  ftoau0,
real (kind=kind_phys), intent(out)  ftoadc,
real (kind=kind_phys), intent(out)  fsfcuc,
real (kind=kind_phys), intent(out)  fsfcu0,
real (kind=kind_phys), intent(out)  fsfcdc,
real (kind=kind_phys), intent(out)  fsfcd0,
real (kind=kind_phys), dimension(2), intent(out)  sfbmc,
real (kind=kind_phys), dimension(2), intent(out)  sfdfc,
real (kind=kind_phys), dimension(2), intent(out)  sfbm0,
real (kind=kind_phys), dimension(2), intent(out)  sfdf0,
real (kind=kind_phys), intent(out)  suvbfc,
real (kind=kind_phys), intent(out)  suvbf0 
)
private

This subroutine computes the shortwave radiative fluxes using two-stream method.

Parameters
ssolarincoming solar flux at top
coszcosine solar zenith angle
sntzsecant solar zenith angle
albbmsurface albedo for direct beam radiation
albdfsurface albedo for diffused radiation
sfluxzenspectral distribution of incoming solar flux
cldfrclayer cloud fraction
cf1>0: cloudy sky, otherwise: clear sky
cf0=1-cf1
taugspectral optical depth for gases
tauroptical depth for rayleigh scattering
tauaeaerosols optical depth
ssaaeaerosols single scattering albedo
asyaeaerosols asymmetry factor
taucwweighted cloud optical depth
ssacwweighted cloud single scat albedo
asycwweighted cloud asymmetry factor
nlay,nlp1number of layers/levels
fxupctot sky upward flux
fxdnctot sky downward flux
fxup0clr sky upward flux
fxdn0clr sky downward flux
ftoauctot sky toa upwd flux
ftoau0clr sky toa upwd flux
ftoadctoa downward (incoming) solar flux
fsfcuctot sky sfc upwd flux
fsfcu0clr sky sfc upwd flux
fsfcdctot sky sfc dnwd flux
fsfcd0clr sky sfc dnwd flux
sfbmctot sky sfc dnwd beam flux (nir/uv+vis)
sfdfctot sky sfc dnwd diff flux (nir/uv+vis)
sfbm0clr sky sfc dnwd beam flux (nir/uv+vis)
sfdf0clr sky sfc dnwd diff flux (nir/uv+vis)
suvbfctot sky sfc dnwd uv-b flux
suvbf0clr sky sfc dnwd uv-b flux

General Algorithm

  1. Initialize output fluxes.
  2. Compute clear-sky optical parameters, layer reflectance and transmittance.
    • Set up toa direct beam and surface values (beam and diff)
    • Delta scaling for clear-sky condition
    • General two-stream expressions for physparam::iswmode
    • Compute homogeneous reflectance and transmittance for both conservative and non-conservative scattering
    • Pre-delta-scaling clear and cloudy direct beam transmittance
    • Call swflux() to compute the upward and downward radiation fluxes
  3. Compute total sky optical parameters, layer reflectance and transmittance.
    • Set up toa direct beam and surface values (beam and diff)
    • Delta scaling for total-sky condition
    • General two-stream expressions for physparam::iswmode
    • Compute homogeneous reflectance and transmittance for conservative scattering and non-conservative scattering
    • Pre-delta-scaling clear and cloudy direct beam transmittance
    • Call swflux() to compute the upward and downward radiation fluxes
  4. Process and save outputs.

Definition at line 2300 of file radsw_main.f.

References bpade, exp_tbl, idxsfc, physparam::iswmode, module_radsw_parameters::nblow, module_radsw_parameters::ngb, module_radsw_parameters::ngptsw, module_radsw_parameters::ntbmx, nuvb, and vrtqdr().

Referenced by swrad().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine module_radsw_main::spcvrtm ( real (kind=kind_phys), intent(in)  ssolar,
real (kind=kind_phys), intent(in)  cosz,
real (kind=kind_phys), intent(in)  sntz,
real (kind=kind_phys), dimension(2), intent(in)  albbm,
real (kind=kind_phys), dimension(2), intent(in)  albdf,
real (kind=kind_phys), dimension(ngptsw), intent(in)  sfluxzen,
real (kind=kind_phys), dimension(nlay,ngptsw), intent(in)  cldfmc,
real (kind=kind_phys), intent(in)  cf1,
real (kind=kind_phys), intent(in)  cf0,
real (kind=kind_phys), dimension(nlay,ngptsw), intent(in)  taug,
real (kind=kind_phys), dimension(nlay,ngptsw), intent(in)  taur,
real (kind=kind_phys), dimension(nlay,nbdsw), intent(in)  tauae,
real (kind=kind_phys), dimension(nlay,nbdsw), intent(in)  ssaae,
real (kind=kind_phys), dimension(nlay,nbdsw), intent(in)  asyae,
real (kind=kind_phys), dimension(nlay,nbdsw), intent(in)  taucw,
real (kind=kind_phys), dimension(nlay,nbdsw), intent(in)  ssacw,
real (kind=kind_phys), dimension(nlay,nbdsw), intent(in)  asycw,
integer, intent(in)  nlay,
integer, intent(in)  nlp1,
real (kind=kind_phys), dimension(nlp1,nbdsw), intent(out)  fxupc,
real (kind=kind_phys), dimension(nlp1,nbdsw), intent(out)  fxdnc,
real (kind=kind_phys), dimension(nlp1,nbdsw), intent(out)  fxup0,
real (kind=kind_phys), dimension(nlp1,nbdsw), intent(out)  fxdn0,
real (kind=kind_phys), intent(out)  ftoauc,
real (kind=kind_phys), intent(out)  ftoau0,
real (kind=kind_phys), intent(out)  ftoadc,
real (kind=kind_phys), intent(out)  fsfcuc,
real (kind=kind_phys), intent(out)  fsfcu0,
real (kind=kind_phys), intent(out)  fsfcdc,
real (kind=kind_phys), intent(out)  fsfcd0,
real (kind=kind_phys), dimension(2), intent(out)  sfbmc,
real (kind=kind_phys), dimension(2), intent(out)  sfdfc,
real (kind=kind_phys), dimension(2), intent(out)  sfbm0,
real (kind=kind_phys), dimension(2), intent(out)  sfdf0,
real (kind=kind_phys), intent(out)  suvbfc,
real (kind=kind_phys), intent(out)  suvbf0 
)
private

This subroutine computes the shortwave radiative fluxes using two-stream method of h. barder and mcica,the monte-carlo independent column approximation, for the representation of sub-grid cloud variability (i.e. cloud overlap).

Parameters
ssolarincoming solar flux at top
coszcosine solar zenith angle
sntzsecant solar zenith angle
albbmsurface albedo for direct beam radiation
albdfsurface albedo for diffused radiation
sfluxzenspectral distribution of incoming solar flux
cldfmclayer cloud fraction for g-point
cf1>0: cloudy sky, otherwise: clear sky
cf0=1-cf1
taugspectral optical depth for gases
tauroptical depth for rayleigh scattering
tauaeaerosols optical depth
ssaaeaerosols single scattering albedo
asyaeaerosols asymmetry factor
taucwweighted cloud optical depth
ssacwweighted cloud single scat albedo
asycwweighted cloud asymmetry factor
nlay,nlp1number of layers/levels
fxupctot sky upward flux
fxdnctot sky downward flux
fxup0clr sky upward flux
fxdn0clr sky downward flux
ftoauctot sky toa upwd flux
ftoau0clr sky toa upwd flux
ftoadctoa downward (incoming) solar flux
fsfcuctot sky sfc upwd flux
fsfcu0clr sky sfc upwd flux
fsfcdctot sky sfc dnwd flux
fsfcd0clr sky sfc dnwd flux
sfbmctot sky sfc dnwd beam flux (nir/uv+vis)
sfdfctot sky sfc dnwd diff flux (nir/uv+vis)
sfbm0clr sky sfc dnwd beam flux (nir/uv+vis)
sfdf0clr sky sfc dnwd diff flux (nir/uv+vis)
suvbfctot sky sfc dnwd uv-b flux
suvbf0clr sky sfc dnwd uv-b flux
  1. Initialize output fluxes.
  2. Compute clear-sky optical parameters, layer reflectance and transmittance.
    • Set up toa direct beam and surface values (beam and diff)
    • Delta scaling for clear-sky condition
    • General two-stream expressions for physparam::iswmode
    • Compute homogeneous reflectance and transmittance for both conservative and non-conservative scattering
    • Pre-delta-scaling clear and cloudy direct beam transmittance
    • Call swflux() to compute the upward and downward radiation fluxes
  3. Compute total sky optical parameters, layer reflectance and transmittance.
    • Set up toa direct beam and surface values (beam and diff)
    • Delta scaling for total-sky condition
    • General two-stream expressions for physparam::iswmode
    • Compute homogeneous reflectance and transmittance for conservative scattering and non-conservative scattering
    • Pre-delta-scaling clear and cloudy direct beam transmittance
    • Call swflux() to compute the upward and downward radiation fluxes
  4. Process and save outputs.

Definition at line 3060 of file radsw_main.f.

References bpade, exp_tbl, idxsfc, physparam::iswmode, module_radsw_parameters::nblow, module_radsw_parameters::ngb, module_radsw_parameters::ngptsw, module_radsw_parameters::ntbmx, nuvb, and vrtqdr().

Referenced by swrad().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine, public module_radsw_main::swrad ( real (kind=kind_phys), dimension(npts,nlay), intent(in)  plyr,
real (kind=kind_phys), dimension(npts,nlp1), intent(in)  plvl,
real (kind=kind_phys), dimension(npts,nlay), intent(in)  tlyr,
real (kind=kind_phys), dimension(npts,nlp1), intent(in)  tlvl,
real (kind=kind_phys), dimension(npts,nlay), intent(in)  qlyr,
real (kind=kind_phys), dimension(npts,nlay), intent(in)  olyr,
real (kind=kind_phys), dimension(npts,nlay,9), intent(in)  gasvmr,
real (kind=kind_phys), dimension(npts,nlay,9), intent(in)  clouds,
integer, dimension(:), intent(in)  icseed,
real (kind=kind_phys), dimension(npts,nlay,nbdsw,3), intent(in)  aerosols,
real (kind=kind_phys), dimension(npts,4), intent(in)  sfcalb,
real (kind=kind_phys), dimension(npts), intent(in)  cosz,
real (kind=kind_phys), intent(in)  solcon,
integer, intent(in)  NDAY,
integer, dimension(:), intent(in)  idxday,
integer, intent(in)  npts,
integer, intent(in)  nlay,
integer, intent(in)  nlp1,
logical, intent(in)  lprnt,
real (kind=kind_phys), dimension(npts,nlay), intent(out)  hswc,
type (topfsw_type), dimension(npts), intent(out)  topflx,
type (sfcfsw_type), dimension(npts), intent(out)  sfcflx,
real (kind=kind_phys), dimension(npts,nlay), intent(out), optional  HSW0,
real (kind=kind_phys), dimension(npts,nlay,nbdsw), intent(out), optional  HSWB,
type (profsw_type), dimension(npts,nlp1), intent(out), optional  FLXPRF,
type (cmpfsw_type), dimension(npts), intent(out), optional  FDNCMP 
)

This subroutine is the main SW radiation routine.

Parameters
plyrmodel layer mean pressure in mb
plvlmodel level pressure in mb
tlyrmodel layer mean temperature in K
tlvlmodel level temperature in K (not in use)
qlyrlayer specific humidity in gm/gm
olyrlayer ozone concentration in gm/gm
gasvmratmospheric constent gases
(:,:,1) - co2 volume mixing ratio
(:,:,2) - n2o volume mixing ratio
(:,:,3) - ch4 volume mixing ratio
(:,:,4) - o2 volume mixing ratio
(:,:,5) - co volume mixing ratio (not used)
(:,:,6) - cfc11 volume mixing ratio (not used)
(:,:,7) - cfc12 volume mixing ratio (not used)
(:,:,8) - cfc22 volume mixing ratio (not used)
(:,:,9) - ccl4 volume mixing ratio (not used)
cloudscloud profile
(:,:,1) - layer total cloud fraction
(:,:,2) - layer in-cloud liq water path ( \(g/m^2\))
(:,:,3) - mean eff radius for liq cloud (micron)
(:,:,4) - layer in-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)
icseedauxiliary special cloud related array.
aerosolsaerosol optical properties
(:,:,:,1) - optical depth
(:,:,:,2) - single scattering albedo
(:,:,:,3) - asymmetry parameter
sfcalbsurface albedo in fraction
(:,1) - near ir direct beam albedo
(:,2) - near ir diffused albedo
(:,3) - uv+vis direct beam albedo
(:,4) - uv+vis diffused albedo
coszcosine of solar zenith angle
solconsolar constant ( \(W/m^2\))
NDAYnum of daytime points
idxdayindex array for daytime points
nptsnumber of horizontal points
nlay,nlp1vertical layer/lavel numbers
lprntlogical check print flag
hswctotal sky heating rates (k/sec or k/day)
topflxradiation fluxes at toa ( \(W/m^2\)), components:
upfxc - total sky upward flux at toa
dnflx - total sky downward flux at toa
upfx0 - clear sky upward flux at toa
sfcflxradiation fluxes at sfc ( \(W/m^2\)), components:
upfxc - total sky upward flux at sfc
dnfxc - total sky downward flux at sfc
upfx0 - clear sky upward flux at sfc
dnfx0 - clear sky downward flux at sfc
hswbspectral band total sky heating rates
hsw0clear sky heating rates (k/sec or k/day)
flxprflevel radiation fluxes ( \( W/m^2 \)), components:
dnfxc - total sky downward flux at interface
upfxc - total sky upward flux at interface
dnfx0 - clear sky downward flux at interface
upfx0 - clear sky upward flux at interface
fdncmpsurface downward fluxes ( \(W/m^2\)), components:
uvbfc - total sky downward uv-b flux at sfc
uvbf0 - clear sky downward uv-b flux at sfc
nirbm - downward surface nir direct beam flux
nirdf - downward surface nir diffused flux
visbm - downward surface uv+vis direct beam flux
visdf - downward surface uv+vis diffused flux

General Algorithm

  1. Compute solar constant adjustment factor (s0fac) according to solcon.
  2. Initial output arrays (and optional) as zero.
  3. Change random number seed value for each radiation invocation (isubcsw =1 or 2).
  4. Prepare surface albedo: bm,df - dir,dif; 1,2 - nir,uvv.
  5. Prepare atmospheric profile for use in rrtm.
  6. Set absorber and gas column amount, convert from volume mixing ratio to molec/cm2 based on coldry (scaled to 1.0e-20)
    • colamt(nlay,maxgas):column amounts of absorbing gases 1 to maxgas are for h2o,co2,o3,n2o,ch4,o2,co, respectively ( \( mol/cm^2 \))
  7. Read aerosol optical properties from 'aerosols'.
  8. Read cloud optical properties from 'clouds'.
  9. Compute fractions of clear sky view:
    • random overlapping
    • max/ran overlapping
    • maximum overlapping
  10. For cloudy sky column, call cldprop() to compute the cloud optical properties for each cloudy layer.
  11. Call setcoef() to compute various coefficients needed in radiative transfer calculations.
  12. Call taumol() to calculate optical depths for gaseous absorption and rayleigh scattering
  13. Call the 2-stream radiation transfer model:
    • if physparam::isubcsw .le.0, using standard cloud scheme, call spcvrtc().
    • if physparam::isubcsw .gt.0, using mcica cloud scheme, call spcvrtm().
  14. Save outputs.

Definition at line 592 of file radsw_main.f.

References cldprop(), physcons::con_amd, physcons::con_amw, physcons::con_avgd, physcons::con_g, heatfac, physparam::iovrsw, ipsdsw0, physparam::isubcsw, physparam::iswcliq, physparam::iswrgas, physparam::ivflip, s0, setcoef(), spcvrtc(), spcvrtm(), and taumol().

Referenced by module_radiation_driver::grrad().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine module_radsw_main::taumol ( real (kind=kind_phys), dimension(nlay,maxgas), intent(in)  colamt,
real (kind=kind_phys), dimension(nlay), intent(in)  colmol,
real (kind=kind_phys), dimension(nlay), intent(in)  fac00,
real (kind=kind_phys), dimension(nlay), intent(in)  fac01,
real (kind=kind_phys), dimension(nlay), intent(in)  fac10,
real (kind=kind_phys), dimension(nlay), intent(in)  fac11,
integer, dimension(nlay), intent(in)  jp,
integer, dimension(nlay), intent(in)  jt,
integer, dimension(nlay), intent(in)  jt1,
integer, intent(in)  laytrop,
real (kind=kind_phys), dimension(nlay), intent(in)  forfac,
real (kind=kind_phys), dimension(nlay), intent(in)  forfrac,
integer, dimension(nlay), intent(in)  indfor,
real (kind=kind_phys), dimension(nlay), intent(in)  selffac,
real (kind=kind_phys), dimension(nlay), intent(in)  selffrac,
integer, dimension(nlay), intent(in)  indself,
integer, intent(in)  nlay,
real (kind=kind_phys), dimension(ngptsw), intent(out)  sfluxzen,
real (kind=kind_phys), dimension(nlay,ngptsw), intent(out)  taug,
real (kind=kind_phys), dimension(nlay,ngptsw), intent(out)  taur 
)
private

This subroutine calculates optical depths for gaseous absorption and rayleigh scattering
subroutine called taumol## (## = 16-29)

Parameters
colamtcolumn amounts of absorbing gases the index are for h2o, co2, o3, n2o, ch4, and o2, respectively \((mol/cm^2)\)
colmoltotal column amount (dry air+water vapor)
facijfor each layer, these are factors that are needed to compute the interpolation factors that multiply the appropriate reference k-values. a value of 0/1 for i,j indicates that the corresponding factor multiplies reference k-value for the lower/higher of the two appropriate temperatures, and altitudes, respectively.
jpthe index of the lower (in altitude) of the two appropriate ref pressure levels needed for interpolation.
jt,jt1the indices of the lower of the two approp ref temperatures needed for interpolation (for pressure levels jp and jp+1, respectively)
laytroptropopause layer index
forfacscale factor needed to foreign-continuum.
forfracfactor needed for temperature interpolation
indforindex of the lower of the two appropriate reference temperatures needed for foreign-continuum interpolation
selffacscale factor needed to h2o self-continuum.
selffracfactor needed for temperature interpolation of reference h2o self-continuum data
indselfindex of the lower of the two appropriate reference temperatures needed for the self-continuum interpolation
nlaynumber of vertical layers
sfluxzenspectral distribution of incoming solar flux
taugspectral optical depth for gases
tauropt depth for rayleigh scattering

General Algorithm

Definition at line 3898 of file radsw_main.f.

References module_radsw_sflux::ibx, module_radsw_sflux::ix1, module_radsw_sflux::ix2, module_radsw_sflux::layreffr, module_radsw_parameters::nbhgh, module_radsw_parameters::nblow, module_radsw_parameters::ngs, module_radsw_sflux::sfluxref01, module_radsw_sflux::sfluxref02, module_radsw_sflux::sfluxref03, module_radsw_sflux::specwt, module_radsw_sflux::strrat, taumol16(), taumol17(), taumol18(), taumol19(), taumol20(), taumol21(), taumol22(), taumol23(), taumol24(), taumol25(), taumol26(), taumol27(), taumol28(), and taumol29().

Referenced by swrad().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine module_radsw_main::vrtqdr ( real (kind=kind_phys), dimension(nlp1), intent(in)  zrefb,
real (kind=kind_phys), dimension(nlp1), intent(in)  zrefd,
real (kind=kind_phys), dimension(nlp1), intent(in)  ztrab,
real (kind=kind_phys), dimension(nlp1), intent(in)  ztrad,
real (kind=kind_phys), dimension(nlp1), intent(in)  zldbt,
real (kind=kind_phys), dimension(nlp1), intent(in)  ztdbt,
integer, intent(in)  NLAY,
integer, intent(in)  NLP1,
real (kind=kind_phys), dimension(nlp1), intent(out)  zfu,
real (kind=kind_phys), dimension(nlp1), intent(out)  zfd 
)
private

This subroutine is called by spcvrtc() and spcvrtm(), and computes the upward and downward radiation fluxes.

Parameters
zrefblayer direct beam reflectivity
zrefdlayer diffuse reflectivity
ztrablayer direct beam transmissivity
ztradlayer diffuse transmissivity
zldbtlayer mean beam transmittance
ztdbttotal beam transmittance at levels
NLAY,NLP1number of layers/levels
zfuupward flux at layer interface
zfddownward flux at layer interface

General Algorithm

  1. Link lowest layer with surface.
  2. Pass from bottom to top.
  3. Upper boundary conditions
  4. Pass from top to bottom
  5. Up and down-welling fluxes at levels.

Definition at line 3762 of file radsw_main.f.

Referenced by spcvrtc(), and spcvrtm().

Here is the caller graph for this function: