Optics
RRTMGP.Optics.AbstractOpticalProps
— TypeAbstractOpticalProps
Optical properties for one scalar and two stream calculations.
RRTMGP.Optics.OneScalar
— TypeOneScalar{FTA2D,AD} <: AbstractOpticalProps
Single scalar approximation for optical depth, used in calculations accounting for extinction and emission
Fields
layerdata
: storage for optical thicknessτ
: view into optical depth
RRTMGP.Optics.TwoStream
— TypeTwoStream{FTA2D} <: AbstractOpticalProps
Two stream approximation for optical properties, used in calculations accounting for extinction and emission
Fields
layerdata
: storage for optical depth, single scattering albedo and asymmerty parameterτ
: view into optical depthssa
: view into single scattering albedog
: view into asymmetry parameter
RRTMGP.Optics.compute_col_gas!
— Functioncompute_col_gas!(
context,
p_lev,
col_dry,
param_set,
vmr_h2o,
lat,
max_threads = Int(256),
)
This function computes the column amounts of dry or moist air.
RRTMGP.Optics.compute_relative_humidity!
— Functioncompute_relative_humidity!(
device::ClimaComms.AbstractCPUDevice,
rh::AbstractArray{FT, 2},
p_lay::AbstractArray{FT, 2},
t_lay::AbstractArray{FT, 2},
param_set::RP.ARP,
vmr_h2o::Union{AbstractArray{FT, 2}, Nothing} = nothing,
) where {FT}
This function computes the relative humidity.
RRTMGP.Optics.compute_optical_props!
— Functioncompute_optical_props!(
op::AbstractOpticalProps,
as::AtmosphericState{FT},
sf::AbstractSourceLW{FT},
gcol::Int,
igpt::Int,
lkp::LookUpLW{FT},
lkp_cld::Union{LookUpCld,PadeCld,Nothing} = nothing,
lkp_aero::Union{LookUpAerosolMerra, Nothing} = nothing,
) where {FT<:AbstractFloat}
Computes optical properties for the longwave problem.
compute_optical_props!(
op::AbstractOpticalProps,
as::AtmosphericState,
gcol::Int,
igpt::Int,
lkp::LookUpSW,
lkp_cld::Union{LookUpCld,PadeCld,Nothing} = nothing,
)
Computes optical properties for the shortwave problem.
compute_optical_props!(
op::AbstractOpticalProps,
as::GrayAtmosphericState,
sf::AbstractSourceLW,
gcol::Int,
)
Computes optical properties for the longwave gray radiation problem.
compute_optical_props!(
op::AbstractOpticalProps,
as::GrayAtmosphericState,
gcol::Int,
)
Computes optical properties for the shortwave gray radiation problem.
RRTMGP.Optics.compute_col_gas_kernel!
— Functioncompute_col_gas_kernel!(
col_gas,
p_lev,
mol_m_dry,
mol_m_h2o,
avogadro
helmert1,
vmr_h2o
lat,
glay, gcol,
)
This function computes the column amounts of dry or moist air.
RRTMGP.Optics.compute_relative_humidity_kernel!
— Functioncomputerelativehumiditykernel!( rh::AbstractArray{FT, 2}, play::AbstractArray{FT, 2}, tlay::AbstractArray{FT, 2}, vmrh2o::AbstractArray{FT, 2}, mwd::FT, tref::FT, qlay_min::FT, glay::Int, gcol::Int, ) where {FT}
This function computes the relative humidity.
RRTMGP.Optics.compute_interp_frac_temp
— Functioncompute_interp_frac_temp(Δ_t_ref, n_t_ref, t_ref, t_lay)
compute interpolation fraction for temperature.
RRTMGP.Optics.compute_interp_frac_press
— Functioncompute_interp_frac_press(
lkp::AbstractLookUp,
p_lay,
tropo,
glay,
gcol,
)
Compute interpolation fraction for pressure.
RRTMGP.Optics.compute_interp_frac_η
— Functioncompute_interp_frac_η(
lkp::AbstractLookUp,
vmr,
tropo,
jtemp,
ibnd,
glay,
gcol,
)
Compute interpolation fraction for binary species parameter.
RRTMGP.Optics.compute_gas_optics
— Functioncompute_gas_optics(
lkp::Union{LookUpLW, LookUpSW},
vmr,
col_dry,
igpt,
ibnd,
p_lay,
t_lay,
glay, gcol,
) where {FT<:AbstractFloat}
Compute optical thickness, single scattering albedo, and asymmetry parameter.
RRTMGP.Optics.compute_τ_minor
— Functioncompute_τ_minor(
lkp::AbstractLookUp,
vmr,
vmr_h2o::FT,
col_dry,
p_lay::FT,
t_lay::FT,
jtemp::Int,
ftemp::FT,
jη1::Int,
jη2::Int,
fη1::FT,
fη2::FT,
igpt,
ibnd,
glay,
gcol,
) where {FT<:AbstractFloat}
Compute optical thickness contributions from minor gases.
RRTMGP.Optics.compute_τ_rayleigh
— Functioncompute_τ_rayleigh(
lkp::LookUpSW,
col_dry::FT,
vmr_h2o::FT,
jtemp::Int,
ftemp::FT,
jη1::Int,
jη2::Int,
fη1::FT,
fη2::FT,
igpt::Int,
) where {FT<:AbstractFloat}
Compute Rayleigh scattering optical depths for shortwave problem
RRTMGP.Optics.build_cloud_mask!
— Functionbuild_cloud_mask!(cld_mask, cld_frac, ::MaxRandomOverlap)
Builds McICA-sampled cloud mask from cloud fraction data for maximum-random overlap
Reference: https://github.com/AER-RC/RRTMG_SW/
RRTMGP.Optics.add_cloud_optics_2stream!
— Functionadd_cloud_optics_2stream!(
τ,
ssa,
g,
cld_mask,
cld_r_eff_liq,
cld_r_eff_ice,
cld_path_liq,
cld_path_ice,
ice_rgh,
lkp_cld,
ibnd;
delta_scaling = false,
)
This function computes the TwoStream clouds optics properties and adds them to the TwoStream gas optics properties.
RRTMGP.Optics.compute_lookup_cld_liq_props
— Functioncompute_lookup_cld_liq_props(
nsize_liq,
radliq_lwr,
radliq_upr,
lut_extliq,
lut_ssaliq,
lut_asyliq,
re_liq,
cld_path_liq,
)
This function computes the TwoStream
cloud liquid properties using the LookUpTable
method.
RRTMGP.Optics.compute_lookup_cld_ice_props
— Functioncompute_lookup_cld_ice_props(
nsize_ice,
radice_lwr,
radice_upr,
lut_extice,
lut_ssaice,
lut_asyice,
re_ice,
cld_path_ice,
)
This function computes the TwoStream
cloud ice properties using the LookUpTable
method.
RRTMGP.Optics.compute_pade_cld_props
— Functioncompute_pade_cld_props(
lkp_cld::PadeCld,
re_liq::FT,
re_ice::FT,
ice_rgh::Int,
cld_path_liq::FT,
cld_path_ice::FT,
ibnd::UInt8,
) where {FT}
This function computes the TwoSteam
cloud optics properties using the pade method.
RRTMGP.Optics.pade_eval
— Functionpade_eval(
ibnd,
re,
irad,
m,
n,
pade_coeffs,
irgh::Union{Int,Nothing} = nothing,
)
Evaluate Pade approximant of order [m/n]
RRTMGP.Optics.add_aerosol_optics_1scalar!
— Functionadd_aerosol_optics_1scalar!(
τ,
aero_mask,
aero_size,
aero_mass,
rel_hum,
lkp_aero,
ibnd,
)
This function computes the OneScalar
aerosol optics properties and adds them to the exising OneScalar
optics properties.
RRTMGP.Optics.add_aerosol_optics_2stream!
— Functionadd_aerosol_optics_2stream!(
τ,
ssa,
g,
aero_mask,
aero_size,
aero_mass,
rel_hum,
lkp_aero,
ibnd;
delta_scaling = false,
)
This function computes the TwoStream
aerosol optics properties and adds them to the exising TwoStream
optics properties.
RRTMGP.Optics.compute_lookup_aerosol
— Functioncompute_lookup_aerosol(lkp_aero, ibnd::Int, aero_mass, aero_size, rh::FT, glay) where {FT}
Compute cumulative aerosol optical properties for various aerosol particles given the aerosol lookup table (lkp_aero
), band number (ibnd
), aerosol mass (aeromass
), aerosol size (aerosize
), and relative humidity (rh
).
RRTMGP.Optics.compute_lookup_dust_props
— Functioncompute_lookup_dust_props(lkp_aero, ibnd, aeromass::FT, aerosize::FT) where {FT}
Compute aerosol optical properties for dust particles given the aerosol lookup table (lkp_aero
), band number (ibnd
), aerosol mass (aeromass
) and aerosol size (aerosize
).
RRTMGP.Optics.compute_lookup_sea_salt_props
— Functioncompute_lookup_sea_salt_props(lkp_aero, ibnd, aeromass::FT, aerosize::FT, rh::FT) where {FT}
Compute aerosol optical properties for sea salt particles given the aerosol lookup table (lkp_aero
), band number (ibnd
), aerosol mass (aeromass
), aerosol size (aerosize
), and relative humidity (rh
).
RRTMGP.Optics.compute_lookup_sulfate_props
— Functioncompute_lookup_sulfate_props(lkp_aero, ibnd, aeromass::FT, rh::FT) where {FT}
Compute aerosol optical properties for sulfate particles given the aerosol lookup table (lkp_aero
), band number (ibnd
), aerosol mass (aeromass
), and relative humidity (rh
).
RRTMGP.Optics.compute_lookup_black_carbon_props
— Functioncompute_lookup_black_carbon_props(lkp_aero, ibnd, aeromass::FT) where {FT}
Compute aerosol optical properties for hydrophobic black carbon particles given the aerosol lookup table (lkp_aero
), band number (ibnd
), and aerosol mass (aeromass
).
compute_lookup_black_carbon_props(lkp_aero, ibnd, aeromass::FT, rh::FT) where {FT}
Compute aerosol optical properties for hydrophilic black carbon particles given the aerosol lookup table (lkp_aero
), band number (ibnd
), aerosol mass (aeromass
), and relative humidity (rh
).
RRTMGP.Optics.compute_lookup_organic_carbon_props
— Functioncompute_lookup_organic_carbon_props(lkp_aero, ibnd, aeromass)
Compute aerosol optical properties for hydrophobic organic carbon particles given the aerosol lookup table (lkp_aero
), band number (ibnd
), and aerosol mass (aeromass
).
compute_lookup_organic_carbon_props(lkp_aero, ibnd, aeromass::FT, rh::FT) where {FT}
Compute aerosol optical properties for hydrophilic organic carbon particles given the aerosol lookup table (lkp_aero
), band number (ibnd
), aerosol mass (aeromass
), and relative humidity (rh
).
RRTMGP.Optics.locate_merra_size_bin
— Functionlocate_merra_size_bin(size_bin_limits, aerosize)
Locate merra bin number for a given aerosol size (aerosize
).
RRTMGP.Optics.compute_gray_optical_thickness_lw
— Functioncompute_gray_optical_thickness_lw(
params::GrayOpticalThicknessSchneider2004{FT},
p0,
Δp,
p,
lat,
) where {FT<:AbstractFloat}
This functions calculates the optical thickness based on pressure and lapse rate for a gray atmosphere. See Schneider 2004, J. Atmos. Sci. (2004) 61 (12): 1317–1340. DOI: https://doi.org/10.1175/1520-0469(2004)061<1317:TTATTS>2.0.CO;2
compute_gray_optical_thickness_lw(
params::GrayOpticalThicknessOGorman2008{FT},
p0,
Δp,
p,
lat,
) where {FT}
This functions calculates the optical thickness based on pressure and lapse rate for a gray atmosphere. See O'Gorman 2008, Journal of Climate Vol 21, Page(s): 3815–3832. DOI: https://doi.org/10.1175/2007JCLI2065.1
RRTMGP.Optics.compute_gray_optical_thickness_sw
— Functioncompute_gray_optical_thickness_sw(
params::GrayOpticalThicknessOGorman2008{FT},
p0,
Δp,
p,
rest...,
) where {FT}
This functions calculates the optical thickness based on pressure for a gray atmosphere. See O'Gorman 2008, Journal of Climate Vol 21, Page(s): 3815–3832. DOI: https://doi.org/10.1175/2007JCLI2065.1
RRTMGP.Optics.interp1d_equispaced
— Functioninterp1d_equispaced(xi::FT, x, y) where {FT<:AbstractFloat}
Perform 1D linear interpolation. This function assumes x
to be uniformly spaced.
RRTMGP.Optics.interp2d
— Functioninterp2d(
fη1::FT,
fη2::FT,
ftemp::FT,
coeff::FTA2D,
jη1::Int,
jη2::Int,
jtemp::Int,
) where {FT,FTA2D}
Perform 2D linear interpolation.
fminor[1, 1] = (FT(1) - fη1) * (1-ftemp)
fminor[2, 1] = fη1 * (1-ftemp)
fminor[1, 2] = (FT(1) - fη2) * ftemp
fminor[2, 2] = fη2 * ftemp
RRTMGP.Optics.interp3d
— Functioninterp3d(
jη1::Int,
jη2::Int,
fη1::FT,
fη2::FT,
jtemp::Int,
ftemp::FT,
jpresst::Int,
fpress::FT,
coeff::FTA3D,
s1::FT = FT(1),
s2::FT = FT(1),
) where {FT<:AbstractFloat,FTA4D<:AbstractArray{FT,4}}
Perform 3D linear interpolation.
fmajor[1, 1, itemp] = (FT(1) - fpress) * fminor[1, itemp]
fmajor[2, 1, itemp] = (FT(1) - fpress) * fminor[2, itemp]
fmajor[1, 2, itemp] = fpress * fminor[1, itemp]
fmajor[2, 2, itemp] = fpress * fminor[2, itemp]
where,
fminor[1, itemp=1] = (FT(1) - fη1) * (1-ftemp)
fminor[2, itemp=1] = fη1 * (1-ftemp)
fminor[1, itemp=2] = (FT(1) - fη2) * ftemp
fminor[2, itemp=2] = fη2 * ftemp
RRTMGP.Optics.increment_2stream
— Functionincrement_2stream!(τ1::FT, ssa1::FT, g1::FT, τ2::FT, ssa2::FT, g2::FT) where {FT}
Increment TwoStream optical properties τ1
, ssa1
and g1
with τ2
, ssa2
and g2
. Here τ
is the optical thickness, ssa
is the single-scattering albedo, and g
is the symmetry parameter.
RRTMGP.Optics.delta_scale
— Functiondelta_scale(τ, ssa, g)
delta-scale two stream optical properties.
RRTMGP.Optics.loc_lower
— Functionloc_lower(xi, Δx, n, x)
Return the location of the left (lower) point of the interval in which xi
is located in vector x
. This function assumes Δx
is uniform.
function loc_lower(xi, x)
Return the location of the left (lower) point of the interval in which xi
is located in vector x
.
RRTMGP.Optics.interp1d_loc_factor
— Functionfunction interp1d_loc_factor(xi::FT, x::AbstractArray{FT, 1}) where {FT <: AbstractFloat}
Computes the weights for linear interpolation. This works with non-uniformly spaced x
.