GMTB Common Community Physics Package (CCPP) Scientific Documentation  Version 1.0
GFS radsw Main

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

Detailed Description

The SW radiation model in the current NOAA Environmental Modeling System (NEMS) was adapted from the RRTM radiation model developed by AER Inc. ([16]; [63]). 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 ([49]; [73]; [6]). The default (the same in operation runs) selection (iswmode=2) activates the Practical Improved Flux Method (PIFM) by [86] . 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) ([71]; [72]), 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.

Argument Table

local_name standard_name long_name units rank type kind intent optional
plyr air_pressure_at_layer_for_radiation_in_hPa air pressure layer hPa 2 real kind_phys in F
plvl air_pressure_at_interface_for_radiation_in_hPa air pressure level hPa 2 real kind_phys in F
tlyr air_temperature_at_layer_for_radiation air temperature layer K 2 real kind_phys in F
tlvl air_temperature_at_interface_for_radiation air temperature level K 2 real kind_phys in F
qlyr water_vapor_specific_humidity_at_layer_for_radiation specific humidity layer kg kg-1 2 real kind_phys in F
olyr ozone_concentration_at_layer_for_radiation ozone concentration layer kg kg-1 2 real kind_phys in F
gasvmr_co2 volume_mixing_ratio_co2 volume mixing ratio co2 kg kg-1 2 real kind_phys in F
gasvmr_n2o volume_mixing_ratio_n2o volume mixing ratio no2 kg kg-1 2 real kind_phys in F
gasvmr_ch4 volume_mixing_ratio_ch4 volume mixing ratio ch4 kg kg-1 2 real kind_phys in F
gasvmr_o2 volume_mixing_ratio_o2 volume mixing ratio o2 kg kg-1 2 real kind_phys in F
gasvmr_co volume_mixing_ratio_co volume mixing ratio co kg kg-1 2 real kind_phys in F
gasvmr_cfc11 volume_mixing_ratio_cfc11 volume mixing ratio cfc11 kg kg-1 2 real kind_phys in F
gasvmr_cfc12 volume_mixing_ratio_cfc12 volume mixing ratio cfc12 kg kg-1 2 real kind_phys in F
gasvmr_cfc22 volume_mixing_ratio_cfc22 volume mixing ratio cfc22 kg kg-1 2 real kind_phys in F
gasvmr_ccl4 volume_mixing_ratio_ccl4 volume mixing ratio ccl4 kg kg-1 2 real kind_phys in F
icseed seed_random_numbers_sw seed for random number generation for shortwave radiation none 1 integer in F
aeraod aerosol_optical_depth_for_shortwave_bands_01-16 aerosol optical depth for shortwave bands 01-16 none 3 real kind_phys in F
aerssa aerosol_single_scattering_albedo_for_shortwave_bands_01-16 aerosol single scattering albedo for shortwave bands 01-16 frac 3 real kind_phys in F
aerasy aerosol_asymmetry_parameter_for_shortwave_bands_01-16 aerosol asymmetry paramter for shortwave bands 01-16 none 3 real kind_phys in F
sfcalb_nir_dir surface_albedo_due_to_near_IR_direct surface albedo due to near IR direct beam frac 1 real kind_phys in F
sfcalb_nir_dif surface_albedo_due_to_near_IR_diffused surface albedo due to near IR diffused beam frac 1 real kind_phys in F
sfcalb_uvis_dir surface_albedo_due_to_UV_and_VIS_direct surface albedo due to UV+VIS direct beam frac 1 real kind_phys in F
sfcalb_uvis_dif surface_albedo_due_to_UV_and_VIS_diffused surface albedo due to UV+VIS diffused beam frac 1 real kind_phys in F
cosz cosine_of_zenith_angle cosine of the solar zenit angle none 1 real kind_phys in F
solcon solar_constant solar constant W m-2 0 real kind_phys in F
nday daytime_points_dimension daytime points dimension count 0 integer in F
idxday daytime_points daytime points index 1 integer in F
npts horizontal_loop_extent horizontal dimension count 0 integer in F
nlay adjusted_vertical_layer_dimension_for_radiation number of vertical layers for radiation count 0 integer in F
nlp1 adjusted_vertical_level_dimension_for_radiation number of vertical levels for radiation count 0 integer in F
lprnt flag_print flag to print flag 0 logical in F
cld_cf total_cloud_fraction total cloud fraction frac 2 real kind_phys in F
lsswr flag_to_calc_sw flag to calculate SW irradiances flag 0 logical in F
hswc tendency_of_air_temperature_due_to_shortwave_heating_on_radiation_time_step shortwave total sky heating rate K s-1 2 real kind_phys inout F
topflx sw_fluxes_top_atmosphere shortwave total sky fluxes at the top of the atm W m-2 1 topfsw_type inout F
sfcflx sw_fluxes_sfc shortwave total sky fluxes at the Earth surface W m-2 1 sfcfsw_type inout F
hsw0 tendency_of_air_temperature_due_to_shortwave_heating_assuming_clear_sky_on_radiation_time_step shortwave clear sky heating rate K s-1 2 real kind_phys inout T
hswb sw_heating_rate_spectral shortwave total sky heating rate (spectral) K s-1 3 real kind_phys inout T
flxprf sw_fluxes sw fluxes total sky / csk and up / down at levels W m-2 2 profsw_type inout T
fdncmp components_of_surface_downward_shortwave_fluxes derived type for special components of surface downward shortwave fluxes W m-2 1 cmpfsw_type inout T
cld_lwp cloud_liquid_water_path cloud liquid water path g m-2 2 real kind_phys in T
cld_ref_liq mean_effective_radius_for_liquid_cloud mean effective radius for liquid cloud micron 2 real kind_phys in T
cld_iwp cloud_ice_water_path cloud ice water path g m-2 2 real kind_phys in T
cld_ref_ice mean_effective_radius_for_ice_cloud mean effective radius for ice cloud micron 2 real kind_phys in T
cld_rwp cloud_rain_water_path cloud rain water path g m-2 2 real kind_phys in T
cld_ref_rain mean_effective_radius_for_rain_drop mean effective radius for rain drop micron 2 real kind_phys in T
cld_swp cloud_snow_water_path cloud snow water path g m-2 2 real kind_phys in T
cld_ref_snow mean_effective_radius_for_snow_flake mean effective radius for snow flake micron 2 real kind_phys in T
cld_od cloud_optical_depth cloud optical depth none 2 real kind_phys in T
cld_ssa cloud_single_scattering_albedo cloud single scattering albedo frac 2 real kind_phys in T
cld_asy cloud_asymmetry_parameter cloud asymmetry parameter none 2 real kind_phys in T
errmsg error_message error message for error handling in CCPP none 0 character len=* out F
errflg error_flag error flag for error handling in CCPP flag 0 integer out F

RRTMG Shortwave Radiation Scheme General Algorithm

Modules

module  module_radsw_parameters
 This module is for specifying the band structures and program parameters used by the RRTMG-SW scheme.
 
module  module_radsw_ref
 This module contains the reference pressures (in logarithm form) at 59 vertical levels (TOA is omitted), and the mid-latitude summer (MLS) standard temperature profile for the 59 pressure layers that are used to establish pre calculated transmission tables.
 
module  module_radsw_cldprtb
 This module contains cloud radiative property coefficients.
 
module  module_radsw_sflux
 This module contains various indexes and coefficients for SW spectral bands, as well as the spectral distribution of solar flux. The values of spectral solar flux are derived based on a prescribed solar constant ( \(1368.22 W/m^2\)). Scaling will be applied for the actual inputted solar constant value.
 
module  module_radsw_kgb16
 This module sets up absorption coefficients for band 16: 2600-3250 cm-1 (low - h2o, ch4; high - ch4)
 
module  module_radsw_kgb17
 This module sets up absorption coeffients for band 17: 3250-4000 cm-1 (low - h2o,co2; high - h2o,co2)
 
module  module_radsw_kgb18
 This module sets up absorption coeffients for band 18: 4000-4650 cm-1 (low - h2o,ch4; high - ch4)
 
module  module_radsw_kgb19
 This module sets up absorption coeffients for band 19: 4650-5150 cm-1 (low - h2o,co2; high - co2)
 
module  module_radsw_kgb20
 This module sets up absorption coeffients for band 20: 5150-6150 cm-1 (low - h2o; high - h2o)
 
module  module_radsw_kgb21
 This module sets up absorption coeffients for band 21: 6150-7700 cm-1 (low - h2o,co2; high - h2o,co2)
 
module  module_radsw_kgb22
 This module sets up absorption coeffients for band 22: 7700-8050 cm-1 (low - h2o, o2; high - o2)
 
module  module_radsw_kgb23
 This module sets up absorption coeffients for band 23: 8050-12850 cm-1 (low - h2o; high - nothing)
 
module  module_radsw_kgb24
 This module sets up absorption coeffients for band 24: 12850-16000 cm-1 (low - h2o, o2; high - o2)
 
module  module_radsw_kgb25
 This module sets up absorption coeffients for band 25: 16000-22650 cm-1 (low - h2o; high - nothing)
 
module  module_radsw_kgb26
 This module sets up absorption coeffients for band 26: 22650-29000 cm-1 (low - nothing; high - nothing)
 
module  module_radsw_kgb27
 This module sets up absorption coeffients for band 27: 29000-38000 cm-1 (low - o3; high - o3)
 
module  module_radsw_kgb28
 This module sets up absorption coeffients for band 28: 38000-50000 cm-1 (low - o3,o2; high - o3,o2)
 
module  module_radsw_kgb29
 This module sets up absorption coeffients for band 29: 820-2600 cm-1 (low - h2o; high - co2)
 

Functions/Subroutines

subroutine, public rrtmg_sw::rrtmg_sw_run (plyr, plvl, tlyr, tlvl, qlyr, olyr, gasvmr_co2, gasvmr_n2o, gasvmr_ch4, gasvmr_o2, gasvmr_co, gasvmr_cfc11, gasvmr_cfc12, gasvmr_cfc22, gasvmr_ccl4, icseed, aeraod, aerssa, aerasy, sfcalb_nir_dir, sfcalb_nir_dif, sfcalb_uvis_dir, sfcalb_uvis_dif, cosz, solcon, NDAY, idxday, npts, nlay, nlp1, lprnt, cld_cf, lsswr, hswc, topflx, sfcflx, HSW0, HSWB, FLXPRF, FDNCMP, cld_lwp, cld_ref_liq, cld_iwp, cld_ref_ice, cld_rwp, cld_ref_rain, cld_swp, cld_ref_snow, cld_od, cld_ssa, cld_asy, errmsg, errflg
 
subroutine, public rrtmg_sw::rswinit (me)
 This subroutine initializes non-varying module variables, conversion factors, and look-up tables. More...
 
subroutine rrtmg_sw::mcica_subcol (cldf, nlay, ipseed, lcloudy )
 This subroutine computes the sub-colum cloud profile flag array. More...
 
subroutine rrtmg_sw::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 rrtmg_sw::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...
 
subroutine taumol17
 The subroutine computes the optical depth in band 17: 3250-4000 cm-1 (low - h2o,co2; high - h2o,co2)
 
subroutine taumol18
 The subroutine computes the optical depth in band 18: 4000-4650 cm-1 (low - h2o,ch4; high - ch4)
 
subroutine taumol19
 The subroutine computes the optical depth in band 19: 4650-5150 cm-1 (low - h2o,co2; high - co2)
 
subroutine taumol20
 The subroutine computes the optical depth in band 20: 5150-6150 cm-1 (low - h2o; high - h2o)
 
subroutine taumol21
 The subroutine computes the optical depth in band 21: 6150-7700 cm-1 (low - h2o,co2; high - h2o,co2)
 
subroutine taumol22
 The subroutine computes the optical depth in band 22: 7700-8050 cm-1 (low - h2o,o2; high - o2)
 
subroutine taumol23
 The subroutine computes the optical depth in band 23: 8050-12850 cm-1 (low - h2o; high - nothing)
 
subroutine taumol24
 The subroutine computes the optical depth in band 24: 12850-16000 cm-1 (low - h2o,o2; high - o2)
 
subroutine taumol25
 The subroutine computes the optical depth in band 25: 16000-22650 cm-1 (low - h2o; high - nothing)
 
subroutine taumol26
 The subroutine computes the optical depth in band 26: 22650-29000 cm-1 (low - nothing; high - nothing)
 
subroutine taumol27
 The subroutine computes the optical depth in band 27: 29000-38000 cm-1 (low - o3; high - o3)
 
subroutine taumol28
 The subroutine computes the optical depth in band 28: 38000-50000 cm-1 (low - o3,o2; high - o3,o2)
 
subroutine taumol29
 The subroutine computes the optical depth in band 29: 820-2600 cm-1 (low - h2o; high - co2)
 
subroutine rrtmg_sw::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 rrtmg_sw::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 rrtmg_sw::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 rrtmg_sw::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 rrtmg_sw::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
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.

References rrtmg_sw::mcica_subcol().

Referenced by rrtmg_sw::rrtmg_sw_run().

Here is the call graph for this function:

subroutine rrtmg_sw::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
Parameters
cldflayer cloud fraction
nlaynumber of model vertical layers
ipseedpermute seed for random num generator
lcloudysub-colum cloud profile flag array

Referenced by rrtmg_sw::cldprop().

subroutine, public rrtmg_sw::rrtmg_sw_run ( 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), intent(in)  gasvmr_co2,
real(kind=kind_phys), dimension(npts,nlay), intent(in)  gasvmr_n2o,
real(kind=kind_phys), dimension(npts,nlay), intent(in)  gasvmr_ch4,
real(kind=kind_phys), dimension(npts,nlay), intent(in)  gasvmr_o2,
real(kind=kind_phys), dimension(npts,nlay), intent(in)  gasvmr_co,
real(kind=kind_phys), dimension(npts,nlay), intent(in)  gasvmr_cfc11,
real(kind=kind_phys), dimension(npts,nlay), intent(in)  gasvmr_cfc12,
real(kind=kind_phys), dimension(npts,nlay), intent(in)  gasvmr_cfc22,
real(kind=kind_phys), dimension(npts,nlay), intent(in)  gasvmr_ccl4,
integer, dimension(:), intent(in)  icseed,
real(kind=kind_phys), dimension(npts,nlay,nbdsw), intent(in)  aeraod,
real(kind=kind_phys), dimension(npts,nlay,nbdsw), intent(in)  aerssa,
real(kind=kind_phys), dimension(npts,nlay,nbdsw), intent(in)  aerasy,
real (kind=kind_phys), dimension(npts), intent(in)  sfcalb_nir_dir,
real (kind=kind_phys), dimension(npts), intent(in)  sfcalb_nir_dif,
real (kind=kind_phys), dimension(npts), intent(in)  sfcalb_uvis_dir,
real (kind=kind_phys), dimension(npts), intent(in)  sfcalb_uvis_dif,
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(in)  cld_cf,
logical, intent(in)  lsswr,
real (kind=kind_phys), dimension(npts,nlay), intent(inout)  hswc,
type (topfsw_type), dimension(npts), intent(inout)  topflx,
type (sfcfsw_type), dimension(npts), intent(inout)  sfcflx,
real (kind=kind_phys), dimension(npts,nlay), intent(inout), optional  HSW0,
real (kind=kind_phys), dimension(npts,nlay,nbdsw), intent(inout), optional  HSWB,
type (profsw_type), dimension(npts,nlp1), intent(inout), optional  FLXPRF,
type (cmpfsw_type), dimension(npts), intent(inout), optional  FDNCMP,
real (kind=kind_phys), dimension(npts,nlay), intent(in), optional  cld_lwp,
real (kind=kind_phys), dimension(npts,nlay), intent(in), optional  cld_ref_liq,
real (kind=kind_phys), dimension(npts,nlay), intent(in), optional  cld_iwp,
real (kind=kind_phys), dimension(npts,nlay), intent(in), optional  cld_ref_ice,
real (kind=kind_phys), dimension(npts,nlay), intent(in), optional  cld_rwp,
real (kind=kind_phys), dimension(npts,nlay), intent(in), optional  cld_ref_rain,
real (kind=kind_phys), dimension(npts,nlay), intent(in), optional  cld_swp,
real (kind=kind_phys), dimension(npts,nlay), intent(in), optional  cld_ref_snow,
real (kind=kind_phys), dimension(npts,nlay), intent(in), optional  cld_od,
real (kind=kind_phys), dimension(npts,nlay), intent(in), optional  cld_ssa,
real (kind=kind_phys), dimension(npts,nlay), intent(in), optional  cld_asy,
character(len=*), intent(out)  errmsg,
integer, intent(out)  errflg 
)
  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.

References rrtmg_sw::cldprop(), rrtmg_sw::setcoef(), rrtmg_sw::spcvrtc(), rrtmg_sw::spcvrtm(), and rrtmg_sw::taumol().

Here is the call graph for this function:

subroutine, public rrtmg_sw::rswinit ( integer, intent(in)  me)
Parameters
meprint control for parallel process
subroutine rrtmg_sw::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
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

Referenced by rrtmg_sw::rrtmg_sw_run().

subroutine rrtmg_sw::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
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.

References rrtmg_sw::vrtqdr().

Referenced by rrtmg_sw::rrtmg_sw_run().

Here is the call graph for this function:

subroutine rrtmg_sw::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
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.

References rrtmg_sw::vrtqdr().

Referenced by rrtmg_sw::rrtmg_sw_run().

Here is the call graph for this function:

subroutine rrtmg_sw::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
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

References taumol16(), taumol17(), taumol18(), taumol19(), taumol20(), taumol21(), taumol22(), taumol23(), taumol24(), taumol25(), taumol26(), taumol27(), taumol28(), and taumol29().

Referenced by rrtmg_sw::rrtmg_sw_run().

Here is the call graph for this function:

subroutine rrtmg_sw::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
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.

Referenced by rrtmg_sw::spcvrtc(), and rrtmg_sw::spcvrtm().