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.jl

Code docs:

Types

CanopyOptics.LeafOpticalPropertiesType
struct PigmentOpticalProperties{FT}

A struct which stores important absorption cross sections of pigments, water, etc

Fields

  • λ

    Wavelength [length]

  • nᵣ

    Refractive index of leaf material

  • Kcab

    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.LeafProspectProPropertiesType
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 units

  • Cw

    Equivalent water thickness [cm], typical about 0.002-0.015

  • Cm

    Dry matter content (dry leaf mass per unit area) [g cm⁻²], typical about 0.003-0.016

  • Cprot

    protein content [g/cm]

  • Ccbc

    Carbone-based constituents content in [g/cm⁻²] (cellulose, lignin, sugars...)

CanopyOptics.dirVectorType
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.compute_Z_matricesMethod
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 model BiLambertianCanopyScattering, uses R,T,nQuad from that model.
  • μ::Array{FT,1}: Quadrature points ∈ [0,1]
  • LD a AbstractLeafDistribution struct 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.createLeafOpticalStructMethod
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 sections

Examples

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.createLeafOpticalStructMethod
createLeafOpticalStruct()
As in createLeafOpticalStruct(λ_bnds) but reads in the in Prospect-PRO database at original resolution (400-2500nm in 1nm steps)
CanopyOptics.dielectricMethod

dielectric(mod::SoilMW, T::FT,f::FT)

Computes the complex dieletric of soil (Ulaby & Long book)

Arguments

  • mod an SoilMW type struct (includes sandfrac, clayfrac, water_frac, ρ)
  • T Temperature in [⁰K]
  • f Frequency in [GHz]

Examples

julia> w = SoilMW(sand_frac=0.2,clay_frac=0.2, water_frac = 0.3, ρ=1.7 )     # Create struct for salty seawater
julia> dielectric(w,283.0,10.0)
53.45306674215052 + 38.068642430090044im
CanopyOptics.dielectricMethod

dielectric(mod::LiquidPureWater, T::FT,f::FT)

Computes the complex dieletric of liquid pure walter (Ulaby & Long book)

Arguments

  • mod an LiquidSaltWater type struct
  • T Temperature in [⁰K]
  • f Frequency in [GHz]

Examples

julia> w = CanopyOptics.LiquidPureWater()     # Create struct for salty seawater
julia> CanopyOptics.dielectric(w,283.0,10.0)
53.45306674215052 + 38.068642430090044im
CanopyOptics.dielectricMethod

dielectric(mod::LiquidSaltWater, T::FT,f::FT)

Computes the complex dieletric of liquid salty walter (Ulaby & Long book)

Arguments

  • mod an LiquidSaltWater type struct, provide Salinity in PSU
  • T Temperature in [⁰K]
  • f Frequency in [GHz]
  • S Salinity in [PSU] comes from the mod LiquidSaltWater struct

Examples

julia> w = CanopyOptics.LiquidSaltWater(S=10.0)     # Create struct for salty seawater
julia> CanopyOptics.dielectric(w,283.0,10.0)
51.79254073931596 + 38.32304044382495im
CanopyOptics.dielectricMethod

dielectric(mod::PureIce, T::FT,f::FT)

Computes the complex dieletric of liquid pure walter (Ulaby & Long book)

Arguments

  • mod an PureIce type struct
  • T Temperature in [⁰K]
  • f Frequency in [GHz]

Examples

julia> w = CanopyOptics.PureIce()     # Create struct for salty seawater
julia> CanopyOptics.dielectric(w,283.0,10.0)
53.45306674215052 + 38.068642430090044im
CanopyOptics.prospectMethod
prospect(leaf::LeafProspectProProperties{FT},
                optis) where {FT<:AbstractFloat}

Computes leaf optical properties (reflectance and transittance) based on PROSPECT-PRO

Arguments

Examples

julia-repl julia> opti = createLeafOpticalStruct((400.0:5:2400)*u"nm"); julia> leaf = LeafProspectProProperties{Float64}(Ccab=30.0); julia> T,R = prospect(leaf,opti);`

CanopyOptics.wood_backwardMethod

Calculation of average backscatter cross-section of a finite length cylinder, exact series solution (KARAM, FUNG, AND ANTAR, 1988)

CanopyOptics.wood_forwardMethod

Calculation of average forward scattering amplitudes of a finite length cylinder, exact series solution