Low-Level Kernel Modules
The documentation below provides details for the low-level modules that wrap the underlying Fortran code. These modules are generated by pybind11 from C headers of RTE-RRTMGP.
The low-level modules provide a direct interface to the Fortran functions in the RTE+RRTMGP package. These modules are used by the high-level modules to provide a more user-friendly interface to the RTE+RRTMGP package. See Python Interface for more information about the high-level modules.
See the RTE-RRTMGP repository on GitHub and the RTE-RRTMGP documentation for more information about these functions.
pyrte_rrtmgp.kernels.rte module
Kernel functions for RTE.
- pyrte_rrtmgp.kernels.rte.delta_scale_2str(nlay: int, ngpt: int, tau: ndarray[tuple[Any, ...], dtype[float64]], ssa: ndarray[tuple[Any, ...], dtype[float64]], g: ndarray[tuple[Any, ...], dtype[float64]]) ndarray[tuple[Any, ...], dtype[float64]]
Apply the delta-scaling transformation to two-stream radiative properties.
- Parameters:
nlay – Number of layers
ngpt – Number of g-points
tau – Optical depth array to be modified (ncol, nlay, ngpt)
ssa – Single scattering albedo array to be modified (ncol, nlay, ngpt)
g – Asymmetry parameter array to be modified (ncol, nlay, ngpt)
- Returns:
Modified tau, ssa, and g arrays
- Return type:
tuple
- pyrte_rrtmgp.kernels.rte.delta_scale_2str_f(nlay: int, ngpt: int, tau: ndarray[tuple[Any, ...], dtype[float64]], ssa: ndarray[tuple[Any, ...], dtype[float64]], g: ndarray[tuple[Any, ...], dtype[float64]], f: ndarray[tuple[Any, ...], dtype[float64]]) ndarray[tuple[Any, ...], dtype[float64]]
Apply the delta-scaling to two-stream with forward scattering fraction.
- Parameters:
nlay – Number of layers
ngpt – Number of g-points
tau – Optical depth array to be modified (ncol, nlay, ngpt)
ssa – Single scattering albedo array to be modified (ncol, nlay, ngpt)
g – Asymmetry parameter array to be modified (ncol, nlay, ngpt)
f – Forward scattering fraction array (ncol, nlay, ngpt)
- Returns:
Modified tau, ssa, and g arrays
- Return type:
tuple
- pyrte_rrtmgp.kernels.rte.inc_1scalar_by_1scalar_bybnd(nlay: int, ngpt: int, tau_inout: ndarray[tuple[Any, ...], dtype[float64]], tau_in: ndarray[tuple[Any, ...], dtype[float64]], nbnd: int, band_lims_gpoint: ndarray[tuple[Any, ...], dtype[int32]]) ndarray[tuple[Any, ...], dtype[float64]]
Increment one set of scalar optical properties with another set by band.
- Parameters:
nlay – Number of layers
ngpt – Number of g-points
tau_inout – Optical depth array to be modified (ncol, nlay, ngpt)
tau_in – Input optical depth array (ncol, nlay, nbnd)
nbnd – Number of bands
band_lims_gpoint – Band limits for g-points (2, nbnd)
- pyrte_rrtmgp.kernels.rte.inc_1scalar_by_2stream_bybnd(nlay: int, ngpt: int, tau_inout: ndarray[tuple[Any, ...], dtype[float64]], tau_in: ndarray[tuple[Any, ...], dtype[float64]], ssa_in: ndarray[tuple[Any, ...], dtype[float64]], nbnd: int, band_lims_gpoint: ndarray[tuple[Any, ...], dtype[int32]]) None
Increment scalar optical properties with 2-stream properties by band.
- Parameters:
nlay – Number of layers
ngpt – Number of g-points
tau_inout – Optical depth array to be modified (ncol, nlay, ngpt)
tau_in – Input optical depth array (ncol, nlay, nbnd)
ssa_in – Input single scattering albedo array (ncol, nlay, nbnd)
nbnd – Number of bands
band_lims_gpoint – Band limits for g-points (2, nbnd)
- pyrte_rrtmgp.kernels.rte.inc_2stream_by_1scalar_bybnd(nlay: int, ngpt: int, tau_inout: ndarray[tuple[Any, ...], dtype[float64]], ssa_inout: ndarray[tuple[Any, ...], dtype[float64]], tau_in: ndarray[tuple[Any, ...], dtype[float64]], nbnd: int, band_lims_gpoint: ndarray[tuple[Any, ...], dtype[int32]]) ndarray[tuple[Any, ...], dtype[float64]]
Increment 2-stream optical properties with scalar properties by band.
- Parameters:
nlay – Number of layers
ngpt – Number of g-points
tau_inout – Optical depth array to be modified (ncol, nlay, ngpt)
ssa_inout – Single scattering albedo array to be modified (ncol, nlay, ngpt)
tau_in – Input optical depth array (ncol, nlay, nbnd)
nbnd – Number of bands
band_lims_gpoint – Band limits for g-points (2, nbnd)
- pyrte_rrtmgp.kernels.rte.inc_2stream_by_2stream_bybnd(nlay: int, ngpt: int, tau_inout: ndarray[tuple[Any, ...], dtype[float64]], ssa_inout: ndarray[tuple[Any, ...], dtype[float64]], g_inout: ndarray[tuple[Any, ...], dtype[float64]], tau_in: ndarray[tuple[Any, ...], dtype[float64]], ssa_in: ndarray[tuple[Any, ...], dtype[float64]], g_in: ndarray[tuple[Any, ...], dtype[float64]], nbnd: int, band_lims_gpoint: ndarray[tuple[Any, ...], dtype[int32]]) ndarray[tuple[Any, ...], dtype[float64]]
Increment one set of 2-stream optical properties with another by band.
- Parameters:
nlay – Number of layers
ngpt – Number of g-points
tau_inout – Optical depth array to be modified (ncol, nlay, ngpt)
ssa_inout – Single scattering albedo array to be modified (ncol, nlay, ngpt)
g_inout – Asymmetry parameter array to be modified (ncol, nlay, ngpt)
tau_in – Input optical depth array (ncol, nlay, nbnd)
ssa_in – Input single scattering albedo array (ncol, nlay, nbnd)
g_in – Input asymmetry parameter array (ncol, nlay, nbnd)
nbnd – Number of bands
band_lims_gpoint – Band limits for g-points (2, nbnd)
- pyrte_rrtmgp.kernels.rte.increment_1scalar_by_1scalar(nlay: int, ngpt: int, tau_inout: ndarray[tuple[Any, ...], dtype[float64]], tau_in: ndarray[tuple[Any, ...], dtype[float64]]) None
Increment one set of optical properties with another set (scalar by scalar).
- Parameters:
nlay – Number of layers
ngpt – Number of g-points
tau_inout – Optical depth array to be modified (ncol, nlay, ngpt)
tau_in – Input optical depth array (ncol, nlay, ngpt)
- pyrte_rrtmgp.kernels.rte.increment_1scalar_by_2stream(nlay: int, ngpt: int, tau_inout: ndarray[tuple[Any, ...], dtype[float64]], tau_in: ndarray[tuple[Any, ...], dtype[float64]], ssa_in: ndarray[tuple[Any, ...], dtype[float64]]) None
Increment scalar optical properties with 2-stream properties.
- Parameters:
nlay – Number of layers
ngpt – Number of g-points
tau_inout – Optical depth array to be modified (ncol, nlay, ngpt)
tau_in – Input optical depth array (ncol, nlay, ngpt)
ssa_in – Input single scattering albedo array (ncol, nlay, ngpt)
- pyrte_rrtmgp.kernels.rte.increment_2stream_by_1scalar(nlay: int, ngpt: int, tau_inout: ndarray[tuple[Any, ...], dtype[float64]], ssa_inout: ndarray[tuple[Any, ...], dtype[float64]], tau_in: ndarray[tuple[Any, ...], dtype[float64]]) None
Increment 2-stream optical properties with scalar properties.
- Parameters:
nlay – Number of layers
ngpt – Number of g-points
tau_inout – Optical depth array to be modified (ncol, nlay, ngpt)
ssa_inout – Single scattering albedo array to be modified (ncol, nlay, ngpt)
tau_in – Input optical depth array (ncol, nlay, ngpt)
- pyrte_rrtmgp.kernels.rte.increment_2stream_by_2stream(nlay: int, ngpt: int, tau_inout: ndarray[tuple[Any, ...], dtype[float64]], ssa_inout: ndarray[tuple[Any, ...], dtype[float64]], g_inout: ndarray[tuple[Any, ...], dtype[float64]], tau_in: ndarray[tuple[Any, ...], dtype[float64]], ssa_in: ndarray[tuple[Any, ...], dtype[float64]], g_in: ndarray[tuple[Any, ...], dtype[float64]]) None
Increment one set of 2-stream optical properties with another.
- Parameters:
nlay – Number of layers
ngpt – Number of g-points
tau_inout – Optical depth array to be modified (ncol, nlay, ngpt)
ssa_inout – Single scattering albedo array to be modified (ncol, nlay, ngpt)
g_inout – Asymmetry parameter array to be modified (ncol, nlay, ngpt)
tau_in – Input optical depth array (ncol, nlay, ngpt)
ssa_in – Input single scattering albedo array (ncol, nlay, ngpt)
g_in – Input asymmetry parameter array (ncol, nlay, ngpt)
- pyrte_rrtmgp.kernels.rte.lw_solver(nlay: int, ngpt: int, ds: ndarray[tuple[Any, ...], dtype[float64]], weights: ndarray[tuple[Any, ...], dtype[float64]], tau: ndarray[tuple[Any, ...], dtype[float64]], ssa: ndarray[tuple[Any, ...], dtype[float64]], g: ndarray[tuple[Any, ...], dtype[float64]], lay_source: ndarray[tuple[Any, ...], dtype[float64]], lev_source: ndarray[tuple[Any, ...], dtype[float64]], sfc_emis: ndarray[tuple[Any, ...], dtype[float64]], sfc_src: ndarray[tuple[Any, ...], dtype[float64]], sfc_src_jac: ndarray[tuple[Any, ...], dtype[float64]], inc_flux: ndarray[tuple[Any, ...], dtype[float64]], top_at_1: bool = True, nmus: int = 1, do_broadband: bool = True, do_Jacobians: bool = False, do_rescaling: bool = False) Tuple[ndarray[tuple[Any, ...], dtype[float64]], ndarray[tuple[Any, ...], dtype[float64]], ndarray[tuple[Any, ...], dtype[float64]], ndarray[tuple[Any, ...], dtype[float64]]]
Perform longwave radiation transfer calculations without scattering.
- This function solves the longwave radiative transfer equation for absorption
and emission.
- Parameters:
nlay – Number of layers
ngpt – Number of g-points
ds – Integration weights with shape (ncol, ngpt, n_quad_angs)
weights – Gaussian quadrature weights with shape (n_quad_angs,)
tau – Optical depths with shape (ncol, nlay, ngpt)
ssa – Single scattering albedos with shape (ncol, nlay, ngpt)
g – Asymmetry parameters with shape (ncol, nlay, ngpt)
lay_source – Layer source terms with shape (ncol, nlay, ngpt)
lev_source – Level source terms with shape (ncol, nlay+1, ngpt)
sfc_emis – Surface emissivities with shape (ncol, ngpt) or (ncol,)
sfc_src – Surface source terms with shape (ncol, ngpt)
sfc_src_jac – Surface source Jacobians with shape (ncol, nlay+1)
inc_flux – Incident fluxes with shape (ncol, ngpt)
top_at_1 – Whether the top of the atmosphere is at index 1
nmus – Number of quadrature points (1-4)
do_broadband – Whether to compute broadband fluxes
do_Jacobians – Whether to compute Jacobians
do_rescaling – Whether to perform flux rescaling
- Returns:
flux_up_jac: Upward flux Jacobians with shape (ncol, nlay+1) broadband_up: Upward broadband fluxes with shape (ncol, nlay+1) broadband_dn: Downward broadband fluxes with shape (ncol, nlay+1) flux_up: Upward fluxes with shape (ncol, nlay+1, ngpt) flux_dn: Downward fluxes with shape (ncol, nlay+1, ngpt)
- Return type:
Tuple containing
- pyrte_rrtmgp.kernels.rte.lw_solver_2stream(ncol: int, nlay: int, ngpt: int, tau: ndarray[tuple[Any, ...], dtype[float64]], ssa: ndarray[tuple[Any, ...], dtype[float64]], g: ndarray[tuple[Any, ...], dtype[float64]], lay_source: ndarray[tuple[Any, ...], dtype[float64]], lev_source: ndarray[tuple[Any, ...], dtype[float64]], sfc_emis: ndarray[tuple[Any, ...], dtype[float64]], sfc_src: ndarray[tuple[Any, ...], dtype[float64]], inc_flux: ndarray[tuple[Any, ...], dtype[float64]], top_at_1: bool = True) Tuple[ndarray[tuple[Any, ...], dtype[float64]], ndarray[tuple[Any, ...], dtype[float64]]]
Solve the longwave radiative transfer equation using the 2-stream approximation.
This function implements the two-stream approximation for longwave radiative transfer, accounting for both absorption and scattering processes.
- Parameters:
ncol – Number of columns
nlay – Number of layers
ngpt – Number of g-points
tau – Optical depths with shape (ncol, nlay, ngpt)
ssa – Single-scattering albedos with shape (ncol, nlay, ngpt)
g – Asymmetry parameters with shape (ncol, nlay, ngpt)
lay_source – Layer source terms with shape (ncol, nlay, ngpt)
lev_source – Level source terms with shape (ncol, nlay+1, ngpt)
sfc_emis – Surface emissivities with shape (ncol, ngpt) or (ncol,)
sfc_src – Surface source terms with shape (ncol, ngpt)
inc_flux – Incident fluxes with shape (ncol, ngpt)
top_at_1 – Whether the top of the atmosphere is at index 1
- Returns:
flux_up: Upward fluxes with shape (ncol, nlay+1, ngpt) flux_dn: Downward fluxes with shape (ncol, nlay+1, ngpt)
- Return type:
Tuple containing
- pyrte_rrtmgp.kernels.rte.sw_solver_2stream(nlay: int, ngpt: int, tau: ndarray[tuple[Any, ...], dtype[float64]], ssa: ndarray[tuple[Any, ...], dtype[float64]], g: ndarray[tuple[Any, ...], dtype[float64]], mu0: ndarray[tuple[Any, ...], dtype[float64]], sfc_alb_dir: ndarray[tuple[Any, ...], dtype[float64]], sfc_alb_dif: ndarray[tuple[Any, ...], dtype[float64]], inc_flux_dir: ndarray[tuple[Any, ...], dtype[float64]], inc_flux_dif: ndarray[tuple[Any, ...], dtype[float64]], top_at_1: bool = True, has_dif_bc: bool = False, do_broadband: bool = True) Tuple[ndarray[tuple[Any, ...], dtype[float64]], ndarray[tuple[Any, ...], dtype[float64]], ndarray[tuple[Any, ...], dtype[float64]], ndarray[tuple[Any, ...], dtype[float64]], ndarray[tuple[Any, ...], dtype[float64]], ndarray[tuple[Any, ...], dtype[float64]]]
Perform shortwave radiation transfer calculations using the 2-stream approximation.
This function implements the two-stream approximation for shortwave radiative transfer, computing direct, diffuse upward and downward fluxes, as well as optional broadband fluxes.
- Parameters:
nlay – Number of layers
ngpt – Number of g-points
tau – Optical depths with shape (ncol, nlay, ngpt)
ssa – Single scattering albedos with shape (ncol, nlay, ngpt)
g – Asymmetry parameters with shape (ncol, nlay, ngpt)
mu0 – Cosine of solar zenith angles with shape (ncol, ngpt)
sfc_alb_dir – Direct surface albedos with shape (ncol, ngpt) or (ncol,)
sfc_alb_dif – Diffuse surface albedos with shape (ncol, ngpt) or (ncol,)
inc_flux_dir – Direct incident fluxes with shape (ncol, ngpt)
inc_flux_dif – Diffuse incident fluxes with shape (ncol, ngpt)
top_at_1 – Whether the top of the atmosphere is at index 1
has_dif_bc – Whether the boundary condition includes diffuse fluxes
do_broadband – Whether to compute broadband fluxes
- Returns:
flux_up: Upward fluxes with shape (ncol, nlay+1, ngpt) flux_dn: Downward fluxes with shape (ncol, nlay+1, ngpt) flux_dir: Direct fluxes with shape (ncol, nlay+1, ngpt) broadband_up: Broadband upward fluxes with shape (ncol, nlay+1) broadband_dn: Broadband downward fluxes with shape (ncol, nlay+1) broadband_dir: Broadband direct fluxes with shape (ncol, nlay+1)
- Return type:
Tuple containing
- pyrte_rrtmgp.kernels.rte.sw_solver_noscat(ncol: int, nlay: int, ngpt: int, tau: ndarray[tuple[Any, ...], dtype[float64]], mu0: ndarray[tuple[Any, ...], dtype[float64]], inc_flux_dir: ndarray[tuple[Any, ...], dtype[float64]], top_at_1: bool = True) ndarray[tuple[Any, ...], dtype[float64]]
Perform shortwave radiation transfer calculations without scattering.
This function solves the shortwave radiative transfer equation in the absence of scattering, computing direct beam fluxes only.
- Parameters:
tau – Optical depths with shape (ncol, nlay, ngpt)
mu0 – Cosine of solar zenith angles with shape (ncol, nlay)
inc_flux_dir – Direct beam incident fluxes with shape (ncol, ngpt)
top_at_1 – Whether the top of the atmosphere is at index 1
- Returns:
Direct-beam fluxes with shape (ncol, nlay+1, ngpt)
pyrte_rrtmgp.kernels.rrtmgp module
Kernel functions for RRTMGP.
- pyrte_rrtmgp.kernels.rrtmgp.compute_cld_from_table(nlay: int, ngpt: int, mask: ndarray[tuple[Any, ...], dtype[bool]], lwp: ndarray[tuple[Any, ...], dtype[float64]], re: ndarray[tuple[Any, ...], dtype[float64]], nsteps: int, step_size: float, offset: float, tau_table: ndarray[tuple[Any, ...], dtype[float64]], ssa_table: ndarray[tuple[Any, ...], dtype[float64]], asy_table: ndarray[tuple[Any, ...], dtype[float64]]) Tuple[ndarray[tuple[Any, ...], dtype[float64]], ndarray[tuple[Any, ...], dtype[float64]], ndarray[tuple[Any, ...], dtype[float64]]]
Compute cloud optical properties using lookup tables.
This function computes cloud optical properties (optical depth, single scattering albedo, and asymmetry parameter) using pre-computed lookup tables. The tables are indexed by liquid water path (lwp) and effective radius (re).
- Parameters:
nlay – Number of atmospheric layers
ngpt – Number of g-points
mask – Cloud mask with shape (ncol, nlay)
lwp – Liquid water path with shape (ncol, nlay)
re – Effective radius with shape (ncol, nlay)
nsteps – Number of steps in the lookup tables
step_size – Step size for the lookup tables
offset – Offset for the lookup tables
tau_table – Optical depth lookup table with shape (nsteps, ngpt)
ssa_table – Single scattering albedo lookup table with shape (nsteps, ngpt)
asy_table – Asymmetry parameter lookup table with shape (nsteps, ngpt)
- Returns:
tau: Cloud optical depth with shape (ncol, nlay, ngpt)
taussa: Product of tau and single scattering albedo with shape (ncol, nlay, ngpt)
taussag: Product of taussa and asymmetry parameter with shape (ncol, nlay, ngpt)
- Return type:
Tuple containing
- pyrte_rrtmgp.kernels.rrtmgp.compute_planck_source(nlay: int, nbnd: int, ngpt: int, nflav: int, neta: int, npres: int, ntemp: int, nPlanckTemp: int, tlay: ndarray[tuple[Any, ...], dtype[float64]], tlev: ndarray[tuple[Any, ...], dtype[float64]], tsfc: ndarray[tuple[Any, ...], dtype[float64]], top_at_1: bool, fmajor: ndarray[tuple[Any, ...], dtype[float64]], jeta: ndarray[tuple[Any, ...], dtype[int32]], tropo: ndarray[tuple[Any, ...], dtype[bool]], jtemp: ndarray[tuple[Any, ...], dtype[int32]], jpress: ndarray[tuple[Any, ...], dtype[int32]], band_lims_gpt: ndarray[tuple[Any, ...], dtype[int32]], pfracin: ndarray[tuple[Any, ...], dtype[float64]], temp_ref_min: float, temp_ref_max: float, totplnk: ndarray[tuple[Any, ...], dtype[float64]], gpoint_flavor: ndarray[tuple[Any, ...], dtype[int32]]) Tuple[ndarray[tuple[Any, ...], dtype[float64]], ndarray[tuple[Any, ...], dtype[float64]], ndarray[tuple[Any, ...], dtype[float64]], ndarray[tuple[Any, ...], dtype[float64]]]
Compute the Planck source function for radiative transfer calculations.
This function calculates the Planck blackbody emission source terms needed for longwave radiative transfer calculations.
- Parameters:
nlay – Number of atmospheric layers
nbnd – Number of spectral bands
ngpt – Number of g-points
nflav – Number of gas flavors
neta – Number of eta points
npres – Number of pressure points
ntemp – Number of temperature points
nPlanckTemp – Number of temperatures for Planck function
tlay – Layer temperatures with shape (ncol, nlay)
tlev – Level temperatures with shape (ncol, nlay+1)
tsfc – Surface temperatures with shape (ncol,)
top_at_1 – Whether the top of the atmosphere is at index 1
fmajor – Major gas interpolation weights with shape (2, 2, 2, ncol, nlay, nflav)
jeta – Eta interpolation indices with shape (2, ncol, nlay, nflav)
tropo – Troposphere mask with shape (ncol, nlay)
jtemp – Temperature interpolation indices with shape (ncol, nlay)
jpress – Pressure interpolation indices with shape (ncol, nlay)
gpoint_bands – g-point for each band (ngpt)
band_lims_gpt – Band limits in g-point space with shape (2, nbnd)
pfracin – Planck fractions with shape (ntemp, neta, npres+1, ngpt)
temp_ref_min – Minimum reference temperature
temp_ref_max – Maximum reference temperature
totplnk – Total Planck function by band with shape (nPlanckTemp, nbnd)
gpoint_flavor – G-point flavors with shape (2, ngpt)
- Returns:
sfc_src: Surface emission with shape (ncol, ngpt)
lay_src: Layer emission with shape (ncol, nlay, ngpt)
lev_src: Level emission with shape (ncol, nlay+1, ngpt)
sfc_src_jac: Surface emission Jacobian with shape (ncol, ngpt)
- Return type:
Tuple containing
- pyrte_rrtmgp.kernels.rrtmgp.compute_tau_absorption(nlay: int, nbnd: int, ngpt: int, ngas: int, nflav: int, neta: int, npres: int, ntemp: int, nminorlower: int, nminorklower: int, nminorupper: int, nminorkupper: int, idx_h2o: int, gpoint_flavor: ndarray[tuple[Any, ...], dtype[int32]], band_lims_gpt: ndarray[tuple[Any, ...], dtype[int32]], kmajor: ndarray[tuple[Any, ...], dtype[float64]], kminor_lower: ndarray[tuple[Any, ...], dtype[float64]], kminor_upper: ndarray[tuple[Any, ...], dtype[float64]], minor_limits_gpt_lower: ndarray[tuple[Any, ...], dtype[int32]], minor_limits_gpt_upper: ndarray[tuple[Any, ...], dtype[int32]], minor_scales_with_density_lower: ndarray[tuple[Any, ...], dtype[bool]], minor_scales_with_density_upper: ndarray[tuple[Any, ...], dtype[bool]], scale_by_complement_lower: ndarray[tuple[Any, ...], dtype[bool]], scale_by_complement_upper: ndarray[tuple[Any, ...], dtype[bool]], idx_minor_lower: ndarray[tuple[Any, ...], dtype[int32]], idx_minor_upper: ndarray[tuple[Any, ...], dtype[int32]], idx_minor_scaling_lower: ndarray[tuple[Any, ...], dtype[int32]], idx_minor_scaling_upper: ndarray[tuple[Any, ...], dtype[int32]], kminor_start_lower: ndarray[tuple[Any, ...], dtype[int32]], kminor_start_upper: ndarray[tuple[Any, ...], dtype[int32]], tropo: ndarray[tuple[Any, ...], dtype[bool]], col_mix: ndarray[tuple[Any, ...], dtype[float64]], fmajor: ndarray[tuple[Any, ...], dtype[float64]], fminor: ndarray[tuple[Any, ...], dtype[float64]], play: ndarray[tuple[Any, ...], dtype[float64]], tlay: ndarray[tuple[Any, ...], dtype[float64]], col_gas: ndarray[tuple[Any, ...], dtype[float64]], jeta: ndarray[tuple[Any, ...], dtype[int32]], jtemp: ndarray[tuple[Any, ...], dtype[int32]], jpress: ndarray[tuple[Any, ...], dtype[int32]]) ndarray[tuple[Any, ...], dtype[float64]]
Compute the absorption optical depth for atmospheric profiles.
This function calculates the total absorption optical depth by combining contributions from major and minor gas species in both the upper and lower atmosphere.
- Parameters:
nlay – Number of atmospheric layers
nbnd – Number of spectral bands
ngpt – Number of g-points
ngas – Number of gases
nflav – Number of gas flavors
neta – Number of eta points
npres – Number of pressure points
ntemp – Number of temperature points
nminorlower – Number of minor species in lower atmosphere
nminorklower – Number of minor absorption coefficients in lower atmosphere
nminorupper – Number of minor species in upper atmosphere
nminorkupper – Number of minor absorption coefficients in upper atmosphere
idx_h2o – Index of water vapor
gpoint_flavor – G-point flavors with shape (2, ngpt)
band_lims_gpt – Band limits in g-point space with shape (2, nbnd)
kmajor – Major gas absorption coefficients
kminor_lower – Minor gas absorption coefficients for lower atmosphere
kminor_upper – Minor gas absorption coefficients for upper atmosphere
minor_limits_gpt_lower – G-point limits for minor gases in lower atmosphere
minor_limits_gpt_upper – G-point limits for minor gases in upper atmosphere
minor_scales_with_density_lower – Density scaling flags for lower atmosphere
minor_scales_with_density_upper – Density scaling flags for upper atmosphere
scale_by_complement_lower – Complement scaling flags for lower atmosphere
scale_by_complement_upper – Complement scaling flags for upper atmosphere
idx_minor_lower – Minor gas indices for lower atmosphere
idx_minor_upper – Minor gas indices for upper atmosphere
idx_minor_scaling_lower – Minor gas scaling indices for lower atmosphere
idx_minor_scaling_upper – Minor gas scaling indices for upper atmosphere
kminor_start_lower – Starting indices for minor gases in lower atmosphere
kminor_start_upper – Starting indices for minor gases in upper atmosphere
tropo – Troposphere mask with shape (ncol, nlay)
col_mix – Gas mixing ratios with shape (2, ncol, nlay, nflav)
fmajor – Major gas interpolation weights
fminor – Minor gas interpolation weights
play – Layer pressures with shape (ncol, nlay)
tlay – Layer temperatures with shape (ncol, nlay)
col_gas – Gas concentrations with shape (ncol, nlay, ngas)
jeta – Eta interpolation indices
jtemp – Temperature interpolation indices
jpress – Pressure interpolation indices
- Returns:
Absorption optical depth with shape (ncol, nlay, ngpt)
- pyrte_rrtmgp.kernels.rrtmgp.compute_tau_rayleigh(nlay: int, nbnd: int, ngpt: int, ngas: int, nflav: int, neta: int, ntemp: int, gpoint_flavor: ndarray[tuple[Any, ...], dtype[int32]], band_lims_gpt: ndarray[tuple[Any, ...], dtype[int32]], krayl: ndarray[tuple[Any, ...], dtype[float64]], idx_h2o: int, col_dry: ndarray[tuple[Any, ...], dtype[float64]], col_gas: ndarray[tuple[Any, ...], dtype[float64]], fminor: ndarray[tuple[Any, ...], dtype[float64]], jeta: ndarray[tuple[Any, ...], dtype[int32]], tropo: ndarray[tuple[Any, ...], dtype[bool]], jtemp: ndarray[tuple[Any, ...], dtype[int32]]) ndarray[tuple[Any, ...], dtype[float64]]
Compute Rayleigh scattering optical depth.
This function calculates the optical depth due to Rayleigh scattering by air molecules.
- Parameters:
nlay – Number of atmospheric layers
nbnd – Number of spectral bands
ngpt – Number of g-points
ngas – Number of gases
nflav – Number of gas flavors
neta – Number of eta points
ntemp – Number of temperature points
gpoint_flavor – G-point flavors with shape (2, ngpt)
band_lims_gpt – Band limits in g-point space with shape (2, nbnd)
krayl – Rayleigh scattering coefficients with shape (ntemp, neta, ngpt, 2)
idx_h2o – Index of water vapor
col_dry – Dry air column amounts with shape (ncol, nlay)
col_gas – Gas concentrations with shape (ncol, nlay, ngas + 1)
fminor – Minor gas interpolation weights
jeta – Eta interpolation indices
tropo – Troposphere mask with shape (ncol, nlay)
jtemp – Temperature interpolation indices
- Returns:
Rayleigh scattering optical depth with shape (ncol, nlay, ngpt)
- pyrte_rrtmgp.kernels.rrtmgp.interpolation(nlay: int, ngas: int, nflav: int, neta: int, npres: int, ntemp: int, flavor: ndarray[tuple[Any, ...], dtype[int32]], press_ref: ndarray[tuple[Any, ...], dtype[float64]], temp_ref: ndarray[tuple[Any, ...], dtype[float64]], press_ref_trop: float, vmr_ref: ndarray[tuple[Any, ...], dtype[float64]], play: ndarray[tuple[Any, ...], dtype[float64]], tlay: ndarray[tuple[Any, ...], dtype[float64]], col_gas: ndarray[tuple[Any, ...], dtype[float64]]) Tuple[ndarray[tuple[Any, ...], dtype[int32]], ndarray[tuple[Any, ...], dtype[float64]], ndarray[tuple[Any, ...], dtype[float64]], ndarray[tuple[Any, ...], dtype[float64]], ndarray[tuple[Any, ...], dtype[bool]], ndarray[tuple[Any, ...], dtype[int32]], ndarray[tuple[Any, ...], dtype[int32]]]
Interpolate the RRTMGP coefficients to the current atmospheric state.
This function performs interpolation of gas optics coefficients based on the current atmospheric temperature and pressure profiles.
- Parameters:
nlay – Number of atmospheric layers
ngas – Number of gases
nflav – Number of gas flavors
neta – Number of mixing fraction points
npres – Number of reference pressure grid points
ntemp – Number of reference temperature grid points
flavor – Index into vmr_ref of major gases for each flavor with shape (nflav,)
press_ref – Reference pressure grid with shape (npres,)
temp_ref – Reference temperature grid with shape (ntemp,)
press_ref_trop – Reference pressure at the tropopause
vmr_ref – Reference volume mixing ratios with shape (ngas,)
play – Layer pressures with shape (ncol, nlay)
tlay – Layer temperatures with shape (ncol, nlay)
col_gas – Gas concentrations with shape (ncol, nlay, ngas)
- Returns:
jtemp: Temperature interpolation indices with shape (ncol, nlay)
fmajor: Major gas interpolation fractions with shape (2, 2, 2, ncol, nlay, nflav)
fminor: Minor gas interpolation fractions with shape (2, 2, ncol, nlay, nflav)
col_mix: Mixing fractions with shape (2, ncol, nlay, nflav)
tropo: Boolean mask for troposphere with shape (ncol, nlay)
jeta: Binary species interpolation indices with shape (2, ncol, nlay, nflav)
jpress: Pressure interpolation indices with shape (ncol, nlay)
- Return type:
Tuple containing