CanopyOptics.jl
A package to compute canopy scattering properties
Package Features
- Use leaf angle distributions to compute bi-Lambertian canopy scattering matrices
- Compute closed-form cosine Fourier moments for the bi-Lambertian canopy kernel
- Compose diffuse and specular leaf-surface scattering models additively
- Request scalar or Stokes-expanded canopy Z matrices with
npol - Combine leaf and wood canopy components with separate angle distributions
- Differentiate canopy Z matrices with ForwardDiff leaf optical parameters
- Compute polarized specular leaf-surface reflection terms
- Compute leaf reflectance and transmittance based on Prospect-PRO
Z-matrix convention
Canopy scattering matrices use Z[i_out, j_in]: rows are outgoing streams and columns are incoming streams. Z⁺⁺ is same-sign transmission and Z⁻⁺ is sign-change reflection. The leaf single-scattering albedo ϖ = R + T is not folded into Z; layer solvers multiply by ϖ separately. For conservative bi-Lambertian leaves, the m = 0 column integral of Z⁺⁺ + Z⁻⁺ is therefore approximately 2 under the quadrature weights.
Specular components include their Fresnel/roughness strength in the returned Z contribution. If a downstream RT solver multiplies canopy Z by a separate single-scattering albedo, it must provide a matching effective albedo when using specular or diffuse+specular composite leaves.
Installation
The latest release of CanopyOptics can be installed from the Julia REPL prompt with
julia> ]add https://github.com/RemoteSensingTools/CanopyOptics.jlCode docs:
Types
CanopyOptics.AbstractCanopyScatteringType — Type
Abstract Type for canopy scattering
CanopyOptics.AbstractClumping — Type
AbstractClumping{FT}Angular clumping model for canopy interception.
Clumping modifies the effective projected area used in gap probability and extinction calculations:
\[G_{eff}(μ) = Ω(μ) G(μ).\]
It does not change the leaf or wood scattering matrix Z; it changes how often canopy elements are encountered along a path.
CanopyOptics.AbstractHotSpot — Type
AbstractHotSpot{FT}Model for bidirectional hotspot corrections to canopy gap probability.
Hotspot terms are propagation corrections, not leaf scattering matrices: they modify the joint probability that the solar and viewing paths see the same canopy gaps. They should therefore be applied where direct-beam source terms are assembled, not inside compute_Z_matrices.
CanopyOptics.AbstractLUTWoodReflectance — Type
AbstractLUTWoodReflectance{FT}Base type for wood reflectance models backed by tabulated spectral data.
CanopyOptics.AbstractLeafDistribution — Type
Abstract Type for leaf distributions
CanopyOptics.AbstractWoodReflectance — Type
AbstractWoodReflectance{FT}Spectral reflectance model for bark or woody canopy elements.
Wood reflectance is kept separate from angular scattering. A reflectance model answers "what is the bark albedo at this wavelength or wavenumber?", while a canopy scattering model answers "how does that woody surface redirect light?".
CanopyOptics.BiLambertianCanopyScattering — Type
Model for bi-lambertian canopy leaf scattering
CanopyOptics.CanopyComponent — Type
CanopyComponent(scatterer, LAD, area_index; clumping = nothing)
CanopyComponent(; scatterer, LAD, area_index, clumping = nothing)One canopy scattering population.
scatterer is a canopy scattering model such as BiLambertianCanopyScattering or LambertianWoodCanopyScattering, LAD is that population's element angle distribution, and area_index is its one-sided area index: LAI for leaves, WAI/BAI for woody material. Optional per-component clumping is used by bulk_G, but mixed Z-matrix normalization keeps using the unclumped Ross G convention.
CanopyOptics.CanopyQuadrature — Type
CanopyQuadrature(; n_leaf = 64, n_azimuth = 48)Integration controls used when projecting canopy scattering models onto Fourier Z matrices.
n_leaf controls the leaf-inclination quadrature used by the bi-Lambertian kernel. n_azimuth controls the outgoing-azimuth quadrature used by the specular kernel and legacy brute-force reference paths.
CanopyOptics.CompositeCanopyScattering — Type
CompositeCanopyScattering(components...)Additive canopy scattering model. Components are evaluated independently and their Z matrices are summed, so a leaf can carry both diffuse bi-Lambertian and specular surface terms.
CanopyOptics.ConstantClumping — Type
ConstantClumping(; Ω = 1)
ConstantClumping(; Ω₀ = 1)Direction-independent clumping index.
Values Ω < 1 reduce effective extinction relative to a random leaf distribution. Ω = 1 preserves the unclumped Beer-Lambert limit. Values above one represent more regular spacing and are allowed for completeness.
CanopyOptics.ConstantWoodReflectance — Type
ConstantWoodReflectance(R)
ConstantWoodReflectance(; R = 0.25)Wavelength-independent bark/wood reflectance.
CanopyOptics.EmpiricalDirectionalClumping — Type
EmpiricalDirectionalClumping(; Ω₀ = 0.7, c = 2, e = 2)Simple empirical angular clumping model:
\[Ω(θ) = \frac{Ω_0}{Ω_0 + (1 - Ω_0)\exp(-c θ^e)}, \qquad θ = \arccos μ.\]
The default shape gives Ω(0) = Ω₀ at nadir and approaches one toward the horizon, where clumping is visually less apparent along the slant path.
This is a compact directional parameterization, not the full Chen-Leblanc row/crown gap-size model. ChenLeblancClumping remains as a compatibility alias for older code that used the previous name.
CanopyOptics.KuuskHotSpot — Type
KuuskHotSpot(; h = 0.1)
KuuskHotSpot{FT}(; h = FT(0.1))Kuusk-style hotspot correction with dimensionless size parameter h, roughly leaf size divided by canopy height. Larger h broadens the hotspot. h <= 0 disables the correction and returns the independent-gap probability.
CanopyOptics.LUTWoodReflectance — Type
LUTWoodReflectance(grid, R; grid_unit = :nm, extrapolation = :clamp)Piecewise-linear bark/wood reflectance lookup table.
grid and R must have the same length, with at least two points. The constructor sorts the table by spectral coordinate. grid_unit is :nm for wavelength in nanometers or :cm_inv for wavenumber. extrapolation = :clamp uses the nearest endpoint outside the tabulated range; :error throws.
CanopyOptics.LambertianWoodCanopyScattering — Type
LambertianWoodCanopyScattering(reflectance)
LambertianWoodCanopyScattering(R)
LambertianWoodCanopyScattering(; R = 0.25, reflectance = nothing)Opaque Lambertian woody element: reflection only, no transmission.
This is the angular-scattering counterpart to AbstractWoodReflectance. For scalar Z-matrix assembly it delegates to BiLambertianCanopyScattering(R, 0), preserving the existing canopy normalization convention.
CanopyOptics.LeafDistribution — Type
struct LeafDistribution{FT<:Real}A struct that defines the leaf angular distribution in radians (from 0->π/2; here scaled to [0,1])
Fields
LD: Julia Univariate Distribution from Distributions.jsscaling: Scaling factor to normalize distribution (here mostly 2/π as Beta distribution is from [0,1])
CanopyOptics.LeafOpticalProperties — Type
struct PigmentOpticalProperties{FT}A struct which stores important absorption cross sections of pigments, water, etc
Fields
λ: Wavelength[length]nᵣ: Refractive index of leaf materialKcab: specific absorption coefficient of chlorophyll (a+b)[cm² μg⁻¹]Kcar: specific absorption coefficient of carotenoids[cm² μg⁻¹]Kant: specific absorption coefficient of Anthocyanins[cm² nmol⁻¹]Kb: specific absorption coefficient of brown pigments (arbitrary units)Kw: specific absorption coefficient of water[cm⁻¹]Km: specific absorption coefficient of dry matter[cm² g⁻¹]Kp: specific absorption coefficient of proteins[cm² g⁻¹]Kcbc: specific absorption coefficient of carbon based constituents[cm² g⁻¹]
CanopyOptics.LeafProspectProProperties — Type
struct LeafProperties{FT}A struct which stores important variables of leaf chemistry and structure
Fields
N: Leaf structure parameter [0-3]Ccab: Chlorophyll a+b content[µg cm⁻²]Ccar: Carotenoid content[µg cm⁻²]Canth: Anthocynanin content[nmol cm⁻²]Cbrown: Brown pigments content in arbitrary unitsCw: Equivalent water thickness[cm], typical about 0.002-0.015Cm: Dry matter content (dry leaf mass per unit area)[g cm⁻²], typical about 0.003-0.016Cprot: protein content[g/cm]Ccbc: Carbone-based constituents content in[g/cm⁻²](cellulose, lignin, sugars...)
CanopyOptics.LiquidPureWater — Type
Pure liquid water
CanopyOptics.LiquidSaltWater — Type
Salty liquid water
CanopyOptics.MixedCanopy — Type
MixedCanopy(components...)
MixedCanopy((components...,))Area-weighted canopy made from multiple scattering populations.
Each component may have its own scatterer, angle distribution, area index, and clumping model. Z matrices are mixed column-by-column using the incoming-stream projected area weights
\[w_c(μ_j) = \frac{AI_c G_c(μ_j)} {\sum_d AI_d G_d(μ_j)}.\]
The optional component clumping models are intentionally not included in this Z normalization; use bulk_G with clumped = true for propagation.
CanopyOptics.NoClumping — Type
NoClumping()
NoClumping{FT}()Default clumping model with Ω(μ) = 1.
CanopyOptics.NoHotSpot — Type
NoHotSpot()
NoHotSpot{FT}()Default model with independent sun and view gaps. The joint gap probability is exp(-(k_s + k_o)L).
CanopyOptics.PolynomialWoodReflectance — Type
PolynomialWoodReflectance(coeffs; grid_unit = :nm, x_offset = 0, x_scale = 1)Polynomial bark/wood reflectance model.
Coefficients use ascending order:
\[R(x) = c_0 + c_1 \hat{x} + c_2 \hat{x}^2 + \cdots, \qquad \hat{x} = (x - x_{offset}) / x_{scale}.\]
Use x_offset and x_scale to fit on a normalized wavelength or wavenumber coordinate instead of raw nanometers.
CanopyOptics.PureIce — Type
Pure Ice
CanopyOptics.SpecularCanopyScattering — Type
Model for specular canopy leaf scattering
CanopyOptics.dirVector — Type
dirVector{FT}Struct for spherical coordinate directions in θ (elevation angle) and ϕ (azimuth angle)
Fields
θϕ
CanopyOptics.dirVector_μ — Type
dirVector_μ{FT}Struct for spherical coordinate directions in θ (elevation angle) and ϕ (azimuth angle)"
Fields
μϕ
Functions
CanopyOptics.G — Method
G(μ::Array{FT}, LD::AbstractLeafDistribution; nLeg=20)Returns the integrated projection of leaf area in the direction of μ, assumes azimuthally uniform distribution and a LD distribution for leaf polar angle θ. This function is often referred to as the function O(B) (Goudriaan 1977) or G(Ζ) (Ross 1975,1981), see Bonan modeling book, eqs. 14.21-14.26.
Arguments
μan array of cos(θ) (directions [0,1])LDanAbstractLeafDistributiontype struct, includes a leaf distribution functionnLegan optional parameter for the number of legendre polynomials to integrate over the leaf distribution (default=20)
Examples
julia> μ,w = CanopyOptics.gauleg(10,0.0,1.0); # Create 10 quadrature points in μ
julia> LD = CanopyOptics.spherical_leaves() # Create a default spherical leaf distribution
julia> G = CanopyOptics.G(μ, LD) # Compute G(μ)
10-element Vector{Float64}:
0.5002522783000879
0.5002715115149204
0.5003537989277846
0.5004432798701134
0.5005134448870893
0.5003026448466977
0.4999186257540982
0.4994511190721635
0.49907252201082375
0.49936166823681594CanopyOptics.G2 — Method
G2(μ::AbstractArray{FT}, LD::AbstractLeafDistribution; nLeg=40)Reference implementation of G that integrates in leaf-normal cosine μ_L instead of leaf inclination angle θ_L.
This path is retained for cross-checking the production θ_L quadrature in G. New code should normally call G.
CanopyOptics.abs_components — Method
Given a complex number, return the complex number with absolute value of both parts
CanopyOptics.afsal — Method
Calculate forward scattering amplitudes
CanopyOptics.asal — Method
Calculate scattering coefficients for a thin disk
CanopyOptics.beta_leaves — Method
beta_leaves(α::FT, β::FT) where FTCreates a LeafDistribution from an arbitrary Beta(α, β) distribution over [0, π/2]. The parameters α and β correspond to ν and μ in the RAMI benchmark convention. Use planophile_leaves, spherical_leaves, etc. for standard named distributions.
CanopyOptics.bfG — Method
bfG(μ::Array{FT}, LD::AbstractLeafDistribution; nLeg=20)Brute-force two-angle reference for the Ross projection factor G(μ).
This directly integrates abs(Ω ⋅ Ω_L) over leaf inclination and leaf azimuth. It is useful for tests or debugging quadrature changes, but it is intentionally slower than G.
CanopyOptics.bulk_G — Method
bulk_G(canopy, μ; clumped = true)Area-index-weighted canopy projected area:
\[G_{bulk}(μ) = \sum_c AI_c G_c(μ),\]
or \sum_c AI_c Ω_c(μ)G_c(μ) when clumped = true.
CanopyOptics.canopy_extinction — Method
canopy_extinction(G_value, μ)Convert projected area G(μ) into the directional extinction coefficient k = G(μ) / μ used by Beer-Lambert canopy gap probabilities.
CanopyOptics.clumping_index — Method
clumping_index(model, μ)Evaluate Ω(μ) for one direction cosine or an array of direction cosines.
CanopyOptics.component_G — Method
component_G(component, μ; clumped = false)Projected area for one CanopyComponent. Set clumped = true to apply the component's clumping model.
CanopyOptics.compute_Z_matrices — Method
compute_Z_matrices(canopy::MixedCanopy, μ, m; quadrature, npol = 1)
compute_Z_matrices(canopy::MixedCanopy, μ, ignored_LD, m; quadrature, npol = 1)Compute canopy Z matrices for a mixture of components with their own angle distributions and area indices. The ignored_LD compatibility form exists so callers that already carry (scatterer, μ, LAD, m) can pass a MixedCanopy as the scatterer without changing their call shape.
CanopyOptics.compute_Z_matrices — Method
compute_Z_matrices(mod::AbstractCanopyScatteringType,
μ, LD, m_range::AbstractUnitRange; quadrature, npol = 1)Compute a stack of canopy Z matrices for Fourier moments in m_range.
The scalar default returns arrays with shape (length(μ), length(μ), length(m_range)). Passing npol = 3 or 4 returns Stokes-expanded arrays with the first two dimensions multiplied by npol. In both cases the package convention is Z[i_out, j_in, k], where k = m - first(m_range) + 1.
CanopyOptics.compute_Z_matrices — Method
compute_Z_matrices(mod::BiLambertianCanopyScattering,
μ::AbstractVector, LD::AbstractLeafDistribution,
m::Integer; quadrature = CanopyQuadrature(), npol = 1) -> (Z⁺⁺, Z⁻⁺)Compute one cosine Fourier moment of the bi-Lambertian canopy scattering matrices. Rows are outgoing streams and columns are incoming streams: Z[i_out, j_in].
npol = 1 returns scalar Stokes-I matrices. npol = 3 and npol = 4 expand the scalar diffuse result into vSmartMOM-style Stokes blocks with only I→I populated.
This method uses the closed-form Fourier moments from compute_Z_matrices_aniso_analytic. For all moments in one call, pass a range such as 0:m_max.
CanopyOptics.compute_Z_matrices — Method
compute_Z_matrices(mod::CompositeCanopyScattering, μ, LD, m; quadrature, npol = 1)Compute one Fourier moment for an additive canopy-scattering model.
Each component is evaluated independently and the resulting (Z⁺⁺, Z⁻⁺) matrices are summed. Diffuse components use the vSmartMOM convention where scalar single-scattering albedo is divided out of Z; specular components carry their Fresnel/roughness strength in their kernel. Consumers that multiply by a separate canopy single-scattering albedo must account for that when using specular components in a full RT layer.
CanopyOptics.compute_Z_matrices — Method
compute_Z_matrices(mod::LambertianWoodCanopyScattering,
μ, LD, m; spectral_coordinate = nothing, grid_unit = :nm,
npol = 1)Compute one Fourier moment for opaque Lambertian wood.
The wood reflectance is evaluated, then the angular Z construction delegates to BiLambertianCanopyScattering(R, 0). For non-constant wood reflectance models, pass spectral_coordinate in the unit specified by grid_unit.
CanopyOptics.compute_Z_matrices — Method
compute_Z_matrices(mod::SpecularCanopyScattering,
μ::AbstractVector, LD::AbstractLeafDistribution,
m::Integer; quadrature = CanopyQuadrature())Convenience method for non-Array vector inputs. It materializes μ with collect and delegates to the specular azimuth-integration method.
CanopyOptics.compute_Z_matrices_aniso_analytic — Method
compute_Z_matrices_aniso_analytic(mod::BiLambertianCanopyScattering,
μ::AbstractVector{FT},
LD::AbstractLeafDistribution,
m_max::Int;
quadrature = CanopyQuadrature(),
npol = 1) -> (Z⁺⁺, Z⁻⁺)Compute all scalar Fourier moments m = 0:m_max of the bi-Lambertian canopy phase matrices without azimuthal quadrature.
The scalar default returns arrays of shape (length(μ), length(μ), m_max + 1). Passing npol = 3 or npol = 4 expands those scalar moments into Stokes-block matrices of shape (npol*length(μ), npol*length(μ), m_max + 1). Because bi-Lambertian scattering is unpolarized here, only the I→I entry of each stream block is populated.
All cases use the vSmartMOM convention
\[Z[i_{out}, j_{in}, m+1],\]
where Z⁺⁺ is same-sign transmission and Z⁻⁺ is sign-change reflection. The single-scattering albedo ϖ = R + T is divided out, so for conservative leaves the m = 0 column integral of Z⁺⁺ + Z⁻⁺ targets ≈ 2 under the same quadrature and G(μ) convention used by vSmartMOM.
Origin trace
Shultis and Myneni (1988), Eq. 45 writes the azimuthally averaged Lambertian canopy kernels as products of one-direction projected-area functions. For each leaf inclination this implementation generalizes that factorization from m = 0 to arbitrary cosine Fourier order by using the closed-form moments returned by _one_sided_projection_moments.
For signed incoming/outgoing directions, let (Pᵢ, Nᵢ) and (Pₒ, Nₒ) be the positive-face and negative-face projection moments. The circular-convolution theorem gives the per-leaf Fourier coefficients
\[Ψ^+_m = P_{i,m}P_{o,m} + N_{i,m}N_{o,m}, \qquad Ψ^-_m = P_{i,m}N_{o,m} + N_{i,m}P_{o,m}.\]
The code multiplies these products by 2 before applying the existing f_0 = 2, f_{m>0} = 4 factors. That leading 2 is not leaf physics; it converts the 1/(2π) Fourier coefficients above to the half-range cosine moment used by the historical compute_Z_matrices_aniso implementation: (2/π)∫_0^π Γ(Δϕ)cos(mΔϕ)dΔϕ. This is the normalization expected by vSmartMOM's elemental kernels, where ϖ is applied outside Z.
CanopyOptics.compute_Z_matrices_aniso_analytic — Method
compute_Z_matrices_aniso_analytic(mod::CompositeCanopyScattering,
μ, LD, m_max; quadrature)Compute all moments 0:m_max for a composite canopy-scattering model by delegating to the stacked compute_Z_matrices API.
CanopyOptics.compute_reflection — Method
compute_reflection(mod::SpecularCanopyScattering, Ωⁱⁿ::dirVector_μ{FT}, Ωᵒᵘᵗ::dirVector_μ{FT},
LD::AbstractLeafDistribution) where FTReturns the specular bidirectional scattering coefficient for direction pair (Ωⁱⁿ, Ωᵒᵘᵗ), following Knyazikhin & Marshak, Discrete Ordinates Method for Photon Transport in Leaf Canopies, Eq. 2.39:
$f_s(Ω' \to Ω) = \frac{1}{8}\,g_L(θ^*)\,K(κ, α^*)\,F_r(n_r, α^*)$
where:
θ*= polar angle of the specular leaf normal (fromgetSpecularΩ)α*= incidence half-angle =arccos(Ωⁱⁿ ⋅ Ωᵒᵘᵗ) / 2K(κ, α*)= Nilson–Kuusk roughness factor (fromK)Fᵣ(nᵣ, α*)= unpolarized Fresnel reflectance (fromFᵣ)
For Stokes-vector propagation, use compute_reflection_mueller, which constructs the Fresnel Mueller matrix and rotates it into the scattering-plane basis used by the canopy Fourier Z matrices.
Note: only reflection is currently modelled; specular transmission is not yet implemented.
CanopyOptics.compute_reflection — Method
compute_reflection(mod::SpecularCanopyScattering,
Ωⁱⁿ::dirVector{FT}, Ωᵒᵘᵗ::dirVector{FT}, LD) where FTAngle-based variant of the specular leaf bidirectional scattering coefficient.
This uses the same Knyazikhin-Marshak specular-leaf model as the dirVector_μ method below, but accepts polar-angle direction vectors. The dirVector_μ method is the one used by Z-matrix assembly; this method is retained for interactive angle-grid checks.
CanopyOptics.compute_reflection_mueller — Method
compute_reflection_mueller(mod, Ωin, Ωout, LD, npol)Specular leaf bidirectional scattering Mueller matrix for one direction pair.
The scalar I→I element reduces to compute_reflection. For npol = 3 or 4, the Fresnel Mueller matrix is rotated between the scattering plane and the local leaf incidence/reflection planes before the Nilson-Kuusk roughness and leaf-normal PDF factors are applied.
CanopyOptics.createLeafOpticalStruct — Method
createLeafOpticalStruct(λ_bnds)
Loads in the PROSPECT-PRO database of pigments (and other) absorption cross section in leaves, returns a [`LeafOpticalProperties`](@ref) type struct with spectral units attached.Arguments
- `λ_bnds` an array (with or without a spectral grid unit) that defines the upper and lower limits over which to average the absorption cross sectionsExamples
julia> using Unitful # Need to include Unitful package
julia> opti = createLeafOpticalStruct((400.0:5:2400)*u"nm"); # in nm
julia> opti = createLeafOpticalStruct((0.4:0.1:2.4)*u"μm"); # in μm
julia> opti = CanopyOptics.createLeafOpticalStruct((10000.0:100:25000.0)u"1/cm"); # in wavenumber (cm⁻¹)CanopyOptics.createLeafOpticalStruct — Method
createLeafOpticalStruct()
As in createLeafOpticalStruct(λ_bnds) but reads in the in Prospect-PRO database at original resolution (400-2500nm in 1nm steps)CanopyOptics.dielectric — Method
dielectric(mod::SoilMW, T::FT,f::FT)
Computes the complex dielectric constant of moist soil (Ulaby & Long book, Eqs. 4.66–4.70).
Arguments
modaSoilMWtype struct (fields:sand_frac,clay_frac,mᵥ,ρ)TTemperature in[K]fFrequency in[GHz]
Examples
julia> w = SoilMW(sand_frac=0.2, clay_frac=0.2, mᵥ=0.3, ρ=1.7)
julia> dielectric(w, 283.0, 10.0)CanopyOptics.dielectric — Method
dielectric(mod::LiquidPureWater, T::FT,f::FT)
Computes the complex dieletric of liquid pure walter (Ulaby & Long book)
Arguments
modanLiquidSaltWatertype structTTemperature in[⁰K]fFrequency in[GHz]
Examples
julia> w = CanopyOptics.LiquidPureWater() # Create struct for salty seawater
julia> CanopyOptics.dielectric(w,283.0,10.0)
53.45306674215052 + 38.068642430090044imCanopyOptics.dielectric — Method
dielectric(mod::LiquidSaltWater, T::FT,f::FT)
Computes the complex dieletric of liquid salty walter (Ulaby & Long book)
Arguments
modanLiquidSaltWatertype struct, provide Salinity in PSUTTemperature in[⁰K]fFrequency in[GHz]SSalinity in[PSU]comes from themodLiquidSaltWaterstruct
Examples
julia> w = CanopyOptics.LiquidSaltWater(S=10.0) # Create struct for salty seawater
julia> CanopyOptics.dielectric(w,283.0,10.0)
51.79254073931596 + 38.32304044382495imCanopyOptics.dielectric — Method
dielectric(mod::PureIce, T::FT,f::FT)
Computes the complex dielectric constant of pure ice (Ulaby & Long book, Eqs. 4.22–4.25). Valid temperature range: 233–273.15 K.
Arguments
modaPureIcetype structTTemperature in[K], must be ∈ [233, 273.15]fFrequency in[GHz]
Examples
julia> w = CanopyOptics.PureIce()
julia> CanopyOptics.dielectric(w, 253.0, 10.0) # T = -20°C, f = 10 GHz
# Real part ≈ 3.17 (ice has low imaginary dielectric compared to liquid water)CanopyOptics.effective_G — Method
effective_G(clumping, G_value, μ)Apply clumping to a Ross projection factor: G_eff = Ω(μ) G.
This helper is intended for propagation/extinction code. Z-matrix construction should generally keep using the unclumped leaf-angle G normalization.
CanopyOptics.erectophile_leaves — Function
erectophile_leaves(FT=Float64)Creates a LeafDistribution following a erectophile (mostly vertical) distribution using a Beta distribution
Standard values for θₗ (63.24) and s (18.5) from Bonan https://doi.org/10.1017/9781107339217.003, Table 2.1
CanopyOptics.flat_leaves — Function
flat_leaves(FT=Float64)Creates a LeafDistribution representing nearly horizontal (flat) leaves using a strongly planophile Beta(1, 100) distribution, which peaks sharply near θ_L ≈ 0.
CanopyOptics.hotspot_correction — Method
hotspot_correction(model, k_s, k_o, μ_s, μ_o, dϕ, L)Return the multiplicative correction C_hs applied to the independent joint gap probability. k_s and k_o are directional extinction coefficients, μ_s and μ_o are positive cosines for sun and view paths, dϕ is relative azimuth in radians, and L is cumulative LAI or area index.
CanopyOptics.hotspot_separation — Method
hotspot_separation(μ_s, μ_o, dϕ, h)Angular separation parameter used by the Kuusk hotspot correction:
\[α = \frac{1}{h}\sqrt{\tan^2 θ_s + \tan^2 θ_o - 2\tan θ_s\tan θ_o\cos Δϕ}.\]
CanopyOptics.joint_gap_probability — Method
joint_gap_probability(model, k_s, k_o, μ_s, μ_o, dϕ, L)Return the full correlated bidirectional gap probability
\[P_{so}(L) = \exp[-(k_s + k_o)L] C_{hs}(L).\]
Use NoHotSpot for the independent-path Beer-Lambert result.
Hotspot corrections belong in the direct-beam source assembly of the canopy RT consumer, not in compute_Z_matrices. A typical per-layer skeleton is:
G_s = G([μ_s], LAD)[1]
G_o = G([μ_o], LAD)[1]
k_s = canopy_extinction(G_s, μ_s)
k_o = canopy_extinction(G_o, μ_o)
P_so = joint_gap_probability(hotspot, k_s, k_o, μ_s, μ_o, dϕ, L_cumulative)
P_o = exp(-k_o * L_cumulative)
direct_source_scale = P_so / P_oThe consumer then applies direct_source_scale to the direct solar source term before multiplying by the canopy scattering kernel.
CanopyOptics.plagiophile_leaves — Function
plagiophile_leaves(FT=Float64)Creates a LeafDistribution following a plagiophile (mostly 45degrees) distribution using a Beta distribution
Standard values for θₗ (45.0) and s (16.27) from Bonan https://doi.org/10.1017/9781107339217.003, Table 2.1
CanopyOptics.planophile_leaves — Function
planophile_leaves(FT=Float64)Creates a LeafDistribution following a planophile (mostly horizontal) distribution using a Beta distribution
Standard values for θₗ (26.76) and s (18.51) from Bonan https://doi.org/10.1017/9781107339217.003, Table 2.1
CanopyOptics.planophile_leaves2 — Function
planophile_leaves2(FT=Float64)Creates a LeafDistribution following a planophile distribution using Beta(1, 2.768), an analytically convenient approximation to Bunnik's planophile distribution. See Shultis & Myneni (1988) Eq. (59): gL(θL) = (2/π)(1 + cos 2θ_L).
CanopyOptics.prospect — Method
prospect(leaf::LeafProspectProProperties{FT},
optis) where {FT<:AbstractFloat}Computes leaf optical properties (reflectance and transittance) based on PROSPECT-PRO
Arguments
leafLeafProspectProPropertiestype struct which provides leaf compositionoptisLeafOpticalPropertiestype struct, which provides absorption cross sections and spectral grid
Examples
julia> opti = createLeafOpticalStruct((400.0:5:2400)*u"nm");
julia> leaf = LeafProspectProProperties{Float64}(Ccab=30.0);
julia> T,R = prospect(leaf,opti);CanopyOptics.spherical_leaves — Function
spherical_leaves(FT=Float64)Creates a LeafDistribution following a spherical (distributed as facets on a sphere) distribution using a Beta distribution
Standard values for θₗ (57.3) and s (21.55) from Bonan https://doi.org/10.1017/9781107339217.003, Table 2.1
CanopyOptics.uniform_leaves — Function
uniform_leaves(FT=Float64)Creates a LeafDistribution following a uniform distribution (all angles equally likely)
CanopyOptics.wood_backward — Method
Calculation of average backscatter cross-section of a finite length cylinder, exact series solution (KARAM, FUNG, AND ANTAR, 1988)
CanopyOptics.wood_forward — Method
Calculation of average forward scattering amplitudes of a finite length cylinder, exact series solution
CanopyOptics.wood_reflectance — Method
wood_reflectance(model, x; grid_unit = :nm)
wood_reflectance(model, grid; grid_unit = :nm)Evaluate a wood reflectance model at one spectral coordinate or a vector of coordinates. grid_unit describes the input coordinate unit.
CanopyOptics.βparameters — Method
βparameters(θₗ, s)Returns Beta distribution parameters (α, β) from the mean leaf inclination angle θₗ and its standard deviation s (both in radians), following Bonan Eqs. 2.9–2.10:
$α = \frac{1 - (s^2 + θₗ^2)/(θₗ π/2)}{(s^2 + θₗ^2)/θₗ^2 - 1}$
$β = \left(\frac{π/2}{θₗ} - 1\right) α$
The resulting Beta(α, β) distribution is defined on [0, 1], which maps to leaf zenith angles [0, π/2] via θ_L = (π/2) u for u ~ Beta(α, β).