CanopyOptics.jl
A package to compute canopy scattering properties
Package Features
- Use leaf angle distributions to compute bi-lambertian area scattering matrices
- Compute specular reflection
- Compute leaf reflectance and transmittance based on Prospect-PRO
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.BiLambertianCanopyScattering — Type
Model for bi-lambertian canopy leaf scattering
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.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.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.compute_Z_matrices — Method
compute_Z_matrices(mod::BiLambertianCanopyScattering, μ::Array{FT,1}, LD::AbstractLeafDistribution, m::Int)Computes the single scattering Z matrices (𝐙⁺⁺ for same incoming and outgoing sign of μ, 𝐙⁻⁺ for a change in direction). Internally computes the azimuthally-averaged area scattering transfer function following Shultis and Myneni (https://doi.org/10.1016/0022-4073(88)90079-9), Eq 43::
$Γ(μ' -> μ) = \int_0^1 dμ_L g_L(μ_L)[t_L Ψ⁺(μ, μ', μ_L) + r_L Ψ⁻(μ, μ', μ_L)]$
assuming an azimuthally uniform leaf angle distribution. Normalized Γ as 𝐙 = 4Γ/(ϖ⋅G(μ)). Returns 𝐙⁺⁺, 𝐙⁻⁺
Arguments
mod: A bilambertian canopy scattering modelBiLambertianCanopyScattering, uses R,T,nQuad from that model.μ::Array{FT,1}: Quadrature points ∈ [0,1]LDaAbstractLeafDistributionstruct that describes the leaf angular distribution function.m: Fourier moment (for azimuthally uniform leave distributions such as here, only m=0 returns non-zero matrices)
CanopyOptics.compute_Z_matrices — Method
compute_Z_matrices(mod::SpecularCanopyScattering, μ::Array{FT,1},
LD::AbstractLeafDistribution, m::Int) where FTComputes the Fourier-m component of the single-scattering phase matrices (𝐙⁺⁺, 𝐙⁻⁺) for a specular leaf surface by integrating compute_reflection over the azimuthal quadrature grid.
𝐙⁺⁺[i,j]: same-hemisphere scattering (μ>0 → μ>0, forward scatter)𝐙⁻⁺[i,j]: opposite-hemisphere scattering (μ>0 → μ<0, backscatter)
Azimuth integration uses nQuad Gauss-Legendre points over [0, 2π] from mod; Fourier weights are cos(m ϕ).
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 full Stokes-vector propagation (vSmartMOM.jl), use fresnel_components to obtain r_s, r_p and construct the 4×4 Mueller reflection matrix directly.
Note: only reflection is currently modelled; specular transmission is not yet implemented.
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.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.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