CoreRT API
CoreRT owns the adding-doubling solver, the RT model containers, and the linearized/Jacobian-facing model layout.
Forward and Linearized Modes
vSmartMOM.FwdMode Type
FwdMode()Marker selecting the forward radiative-transfer model construction path.
sourcevSmartMOM.LinMode Type
LinMode()Marker selecting the linearized model construction path used for Jacobians.
sourcevSmartMOM.CoreRT.model_from_parameters_lin Function
model_from_parameters_lin(params)Convenience alias for model_from_parameters(LinMode(), params). Returns (model, lin_model).
vSmartMOM.CoreRT.rt_run_lin Function
rt_run_lin(model, lin_model, NAer, NGas, NSurf; i_band=1)Convenience alias for the linearized rt_run overload. Equivalent to rt_run(model, lin_model, NAer, NGas, NSurf; i_band).
Model Types
vSmartMOM.CoreRT.RTModel Type
RTModel{ARCH, FT} <: AbstractRTModel{ARCH, FT}Central model state for vSmartMOM radiative transfer.
Only two type parameters:
ARCH: compute architecture (CPUorGPU)FT: floating-point precision (Float32orFloat64)
All physics sub-components are organized hierarchically:
solver: RT solver configuration (polarization, quadrature, truncation)geometry: observation geometry (SZA, VZA, VAZ, observer altitude)quad_points: precomputed quadrature points and weightsatmosphere: atmospheric state (profile + spectral bands)optics: all precomputed optical properties (Rayleigh, aerosol, absorption)surfaces: per-band surface BRDF models
Fields
architecture: Compute architecture (CPU/GPU)solver: RT solver configurationnumerics: Numerical-knob parameters (doubling threshold, etc.) — seeRTNumericalParametersgeometry: Observation geometryquad_points: Quadrature points and weightsatmosphere: Atmospheric state (profile + spectral bands)optics: Precomputed optical propertiessurfaces: Surface models, one per spectral bandsources: Source-term scene (v0.6 source-term refactor). Defaults to a singleSolarBeam— i.e. the historical unit-Stokes-I direct solar illumination — and can be overridden atrt_runcall time via thesources=kwarg. Mirrors the abstract-typedsurfacesslot.
vSmartMOM.CoreRT.SolverConfig Type
SolverConfig{FT, PT, QT}Immutable RT solver configuration: polarization mode, quadrature scheme, Fourier truncation, and Legendre truncation. Fixed after model construction; never differentiated.
Fields
polarization_type: Type of polarization (Stokes_I / IQU / IQUV)quadrature_type: Quadrature type (RadauQuad / GaussLegQuad)m_max_bands: Per-band Fourier loop bound, order semantics: loop runsm = 0:m_max_bands[iBand]. Equalsn_fourier_moments_bands .- 1.n_fourier_moments_bands: Per-band Fourier moment count:n_fourier_moments_bands[iBand] = m_max_bands[iBand] + 1. Provided so consumers don't have to remember the count↔order convention.l_max: Per-band max truncated Legendre indexl_trunc: Legendre truncation order (user-specified)Δ_angle: Exclusion angle for forward peak [deg]depol: Depolarization factoruse_component_traits: Phase C flag: whentrue(current default), per-bandm_max_bandsare derived fromcomponent_m_max(c, ctx)traits across active components — Cox-Munk / RPV / RossLi / canopy run to their fulluser_l_capinstead of the historical half-truncated aggregator. Flip tofalseto fall back to the legacymin(ceil((l_max+1)/2), params.max_m)aggregator (bit-equal to Phase B).
vSmartMOM.CoreRT.Atmosphere Type
Atmosphere{FT, VMR}The atmospheric column: profile data plus spectral grid definitions.
Fields
profile: Atmospheric profile (T, p, q, vmr, vcd, etc.)spec_bands: Spectral bands — Vector of wavenumber grids, one per band
vSmartMOM.CoreRT.RayleighScattering Type
RayleighScattering{FT, GC}Precomputed Rayleigh and Cabannes (elastic) scattering properties. Derived from the depolarization ratio; fixed for a given spectral band.
Fields
greek_rayleigh: Greek coefficients for total Rayleigh scattering (single set or per-band)greek_cabannes: Greek coefficients for Cabannes (pure elastic) Rayleigh per bandϖ_Cabannes: Pure elastic (Cabannes) fraction of Rayleigh scattering per band
vSmartMOM.CoreRT.AerosolState Type
AerosolState{FT, AO}Per-band aerosol scattering optics and optical depth profiles. Primary differentiable state for aerosol retrievals.
Fields
aerosol_optics: Truncated aerosol optics: aerosol_optics[iBand][iAer]τ_aer: Aerosol optical depth profiles:τ_aer[iBand][iAer, iLayer](or[iAer, nSpec, iLayer]for lin)
vSmartMOM.CoreRT.Optics Type
Optics{FT, RS, AS}Container for all precomputed optical properties that feed into the RT solver. Groups Rayleigh scattering, aerosol scattering, absorption optical depths, and Rayleigh optical depths.
AD boundary: aerosols and τ_abs hold differentiable state; rayleigh and τ_rayl are typically fixed.
Fields
rayleigh: Rayleigh/Cabannes scattering propertiesaerosols: Aerosol scattering optics and optical depthsτ_abs: Absorption optical depth: τ_abs[iBand][nSpec × nLayers]τ_rayl: Rayleigh optical depth: τ_rayl[iBand][nSpec × nLayers]
vSmartMOM.CoreRT.OpticsLin Type
OpticsLin{FT}Linearized (Jacobian) counterpart of Optics. Each field stores derivatives of the corresponding forward-model field with respect to the physical state vector.
Fields
τ̇_abs: ∂τ_abs/∂x per band: Vector of arrays [NGas × nSpec × nLayers]τ̇_aer: ∂τ_aer/∂x per band: Vector of arrays [NAer × 7 × nSpec × nLayers]lin_aerosol_optics: Linearized aerosol optics per band per aerosol
vSmartMOM.CoreRT.vSmartMOM_Parameters Type
vSmartMOM_Parameters{FT<:Real}Top-level container for all user-specified model parameters before any derived quantities are computed. Groups radiative-transfer settings (spectral bands, BRDF, quadrature, polarization), observation geometry (SZA, VZA, VAZ), the atmospheric profile (T, p, q), and optional absorption/scattering sub-parameter structs.
Fields
spec_bands: Spectral bands — Vector of wavenumber grids, one per bandbrdf: Surface (Bidirectional Reflectance Distribution Function)quadrature_type: Quadrature type for RT streams (RadauQuad/GaussLegQuad)polarization_type: Type of polarization (I/IQ/IQU/IQUV)max_m: Hard cutoff for maximum number of Fourier moments to loop overΔ_angle: Exclusion angle for forward peak in degrees (legacy — seetruncation)l_trunc: Truncation length for legendre terms (scalar for now, can donBandlater)truncation: Phase-function truncation method. Defaults toδBGE(l_trunc, Δ_angle)for backward compatibility; set explicitly toNoTruncation()for canopy-only or smooth-phase-function runs (Sanghavi & Stephens 2015 §2 — thef_tr → 0limit). TheΔ_anglefield above is then ignored.depol: Depolarization factornumerics: Numerical-knob parameters (doubling threshold, etc.) — defaults toRTNumericalParameters{FT}()if not specified in YAML.float_type: Float type to use in the RT (Float64/Float32)architecture: Architecture to use for calculations (CPU/GPU)sza: Solar zenith angle [deg]vza: Viewing zenith angles [deg]vaz: Viewing azimuthal angles [deg]obs_alt: Altitude of observer [Pa]T: Temperature Profile [K]p: Pressure Profile [hPa]q: Specific humidity profileprofile_reduction_n: Length of profile reductionabsorption_params: Optional struct that holds all absorption-related parametersscattering_params: Optional struct that holds all aerosol scattering-related parametersnstreams: Phase D — primary user-facing resolution knob: weighted streams per hemisphere. Public contract:stream_l_cap = 2·nstreams - 1. Whennothing, the parser populatedmax_m/l_truncfrom a legacy YAML; new-schema configs set this and may omitmax_m/l_trunc.m_max_override: Phase D — explicit per-band Fourier loop bound (order, optional). When set, clamps the trait-aggregator output. Lower hard cap on top ofstream_l_cap.stream_l_cap: Phase D — derived projection cap:2·nstreams - 1for new schema, orl_truncfor legacy. Internally consumed by the trait aggregator and the truncation-resolver. Always populated; defaults tol_truncuntilnstreamsis explicitly set.legacy_l_cap_override: Phase D — legacyl_truncvalue retained verbatim from YAML when the user explicitly set it.nothingfor new-schema configs that derivestream_l_capfromnstreamsinstead.
vSmartMOM.CoreRT.AtmosphericProfile Type
AtmosphericProfile{FT, VMR}Vertical atmospheric state on a pressure grid, ordered from TOA to BOA.
Stores temperature, pressure (full levels and half levels), humidity, dry and wet vertical column densities, and volume mixing ratios for trace gases. Constructed internally by model_from_parameters from the raw arrays in vSmartMOM_Parameters.
p_full has N levels; p_half has N-1 layer-boundary pressures. T, q, vmr_h2o, vcd_dry, vcd_h2o, and Δz are layer quantities of length N.
Fields
T: Temperature Profilep_full: Pressure Profile (Full)q: Specific humidity profilep_half: Pressure Levelsvmr_h2o: H2O Volume Mixing Ratio Profilevcd_dry: Vertical Column Density (Dry)vcd_h2o: Vertical Column Density (H2O)vmr: Volume Mixing Ratio of Constituent GasesΔz: Layer height (meters)
vSmartMOM.CoreRT.QuadPoints Type
QuadPoints{FT}Gauss or Gauss-Radau quadrature points and weights used for the angular discretisation of the radiative transfer equation. Also stores the cosine of the solar zenith angle (μ₀) and its index within the quadrature grid.
Fields
μ₀: μ₀, cos(SZA)iμ₀: Index in quadrature points with suniμ₀Nstart: Index in quadrature points with sun (in qp_μN)qp_μ: Quadrature pointswt_μ: Weights of quadrature pointsqp_μN: Quadrature points (repeated for polarizations)wt_μN: Weights of quadrature points (repeated for polarizations)Nquad: Total number of quadrature points (weighted streams + zero-weight SZA/VZA output nodes)Nstreams: Number of weighted streams per hemisphere (count of nonzero weights). Public contract: stream_l_cap = 2·Nstreams - 1.
vSmartMOM.CoreRT.CompositeLayer Type
CompositeLayer{FT} <: AbstractLayerAccumulated (composite) layer matrices produced by the interaction step of the adding/doubling RT method. Stores the reflectance (R), transmission (T), and source (J) matrices for the combined atmospheric slab from the top of atmosphere down to the current layer.
Sign convention: − = outgoing (upward, decreasing τ), + = incoming (downward, increasing τ).
Fields
R⁻⁺::AbstractArray{FT,3}: reflectance from incoming (+) to outgoing (−)R⁺⁻::AbstractArray{FT,3}: reflectance from outgoing (−) to incoming (+)T⁺⁺::AbstractArray{FT,3}: transmission in the incoming (+) directionT⁻⁻::AbstractArray{FT,3}: transmission in the outgoing (−) directionJ₀⁺::AbstractArray{FT,3}: source in the incoming (+) directionJ₀⁻::AbstractArray{FT,3}: source in the outgoing (−) direction
vSmartMOM.CoreRT.AddedLayer Type
AddedLayer{FT} <: AbstractLayerSingle homogeneous layer matrices produced by the elemental/doubling steps of the RT solver. After doubling, these represent a full atmospheric layer that is subsequently combined with the composite layer via the interaction step.
Lower-case field names (r, t, j) distinguish added-layer quantities from composite-layer quantities (R, T, J).
Fields
r⁻⁺::AbstractArray{FT,3}: reflectance from incoming (+) to outgoing (−)t⁺⁺::AbstractArray{FT,3}: transmission in the incoming (+) directionr⁺⁻::AbstractArray{FT,3}: reflectance from outgoing (−) to incoming (+)t⁻⁻::AbstractArray{FT,3}: transmission in the outgoing (−) directionj₀⁺::AbstractArray{FT,3}: source in the incoming (+) directionj₀⁻::AbstractArray{FT,3}: source in the outgoing (−) directiontemp1,temp2: pre-allocated workspace arrays to avoid allocationstemp1_ptr,temp2_ptr: CUDA pointer arrays (ignored on CPU)
vSmartMOM.CoreRT.RTNumericalParameters Type
RTNumericalParameters{FT}Centralised home for tunable numerical knobs that were previously hardcoded in the RT kernels. Lives once on vSmartMOM_Parameters / RTModel, threaded through to the kernels where needed. Designed to grow as more constants are lifted (Glint Legendre expansion cap, etc.).
Fields
dτ_max_threshold: Doubling resolution: cap on the elemental-layer optical depthdτ_max = threshold · μ_min, whereμ_minis the smallest TRUE (positive-weight) quadrature stream cosine. Smaller values produce finer initial layers (more doublings → more accurate single-scatter at the cost of more roundoff through doubling). Default0.001is conservative legacy; raising to~0.01reduces doubling iterations.dτ_min_floor: Absolute floor ondτ_max(and hence ondτ_initial), expressed as a multiple ofeps(FT). Prevents the doubling discretisation from collapsing below FT precision regardless of geometry. Default1024·eps(FT)≈ 1.2e-4 for Float32 / 2.3e-13 for Float64 — the F32 value capsndoublaround 13–14 for τ_total≤1 (well within Float32 representability); the F64 value is far below any sensiblethreshold·μ_minso it never activates in practice.blas_threads: BLAS thread cap applied at everyrt_runinvocation.nothingmeans leave the BLAS thread setting alone (use whateverLinearAlgebra.BLAS.get_num_threads()returned at session start). Why cap: the batched-GEMM call sites incpu_batched.jland the elemental/doubling/interaction kernels operate on small matrices (NSTREAMS·n_stokes ≈ 12-48 per side) tiled across a wide spectral batch. Multi-threaded BLAS coordination dominates the work at that matrix size and serializes against the spectral-batch parallelism. The VLIDORT baseline harness empirically picked 8 as the sweet spot on a 128-thread machine;1is the most defensive choice when a single process owns all spectral points and threading happens at the batched outer level. Set via thenumerics.blas_threadsYAML key, or programmatically withRTNumericalParameters{FT}(blas_threads = 8).verbose: Verbose output flag. Whenfalse(default)rt_run/rt_run_linsuppress the per-callprint_timer()dump at the end of the run. Settrue(YAML keynumerics.verbose: trueor programmaticallyRTNumericalParameters{FT}(verbose = true)) to bring back the timing tree — useful for profiling, noisy for production loops. Per-bandComputing profile for X...andAOD at band X...messages were demoted to@debugand are silent unless you setENV["JULIA_DEBUG"] = "vSmartMOM".
Source-term types (
SolarBeam,BlackbodySource,SurfaceSIF,SourceSet,NoSource, AD-mode traits,prepare_source/surface_source_contribute!dispatchers) are documented on the Source terms page.
Jacobian Parameter Layout
vSmartMOM.CoreRT.ParameterLayout Type
ParameterLayoutDescribes the ordering of physical parameters in the Jacobian derivative dimension.
Instead of hardcoding 7*NAer + NGas + NSurf throughout the codebase, all index arithmetic goes through this struct. Each aerosol carries aerosol_params sub-parameters (currently 7: τ_ref, nᵣ, nᵢ, rₘ, σ_g, p₀, σ_p), followed by one slot per gas VMR and one per surface parameter.
Example
layout = ParameterLayout(n_aerosols=1, n_gases=2, n_surface=1)
aerosol_range(layout, 1) # 1:7
gas_range(layout) # 8:9
surface_range(layout) # 10:10
n_total(layout) # 10vSmartMOM.CoreRT.surface_index Function
Index of a specific surface parameter (1-based within the surface block).
sourcevSmartMOM.CoreRT.n_layer_params Function
Number of layer-level parameters (aerosol + gas, excluding surface and canopy).
sourcevSmartMOM.CoreRT.canopy_range Function
Index range for canopy parameters (LAI, leaf_R, leaf_T, ...).
source