CanopyOptics.AbstractCanopyScatteringType
— TypeAbstract Type for canopy scattering
CanopyOptics.AbstractLeafDistribution
— TypeAbstract Type for leaf distributions
CanopyOptics.AbstractLeafProperties
— TypeAbstract Type for leaf composition
CanopyOptics.AbstractMaterial
— TypeAbstract Material Type (can be all)
CanopyOptics.AbstractProspectProperties
— TypeAbstract type for Prospect model versions
CanopyOptics.AbstractSoil
— TypeAbstract soil type
CanopyOptics.AbstractWater
— TypeAbstract water type
CanopyOptics.BiLambertianCanopyScattering
— TypeModel for bi-lambertian canopy leaf scattering
CanopyOptics.Leaf
— TypeAll leaf-related parameters
CanopyOptics.LeafDistribution
— Typestruct LeafDistribution{FT<:AbstractFloat}
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
— Typestruct 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
— Typestruct 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
— TypePure liquid water
CanopyOptics.LiquidSaltWater
— TypeSalty liquid water
CanopyOptics.PureIce
— TypePure Ice
CanopyOptics.SoilMW
— TypeSoil MW properties
CanopyOptics.SpecularCanopyScattering
— TypeModel for specular canopy leaf scattering
CanopyOptics.Wood
— TypeAll wood-related parameters (branch or trunk)
CanopyOptics.dirVector
— TypedirVector{FT}
Struct for spherical coordinate directions in θ (elevation angle) and ϕ (azimuth angle)
Fields
θ
ϕ
CanopyOptics.dirVector_μ
— TypedirVector_μ{FT}
Struct for spherical coordinate directions in θ (elevation angle) and ϕ (azimuth angle)"
Fields
μ
ϕ
CanopyOptics.A
— MethodIntegrated projection of leaf area for a single leaf inclination of θₗ, assumes azimuthally uniform distribution
CanopyOptics.Asm
— MethodIntegrated projection of leaf area for a single leaf inclination of θₗ, assumes azimuthally uniform distribution
CanopyOptics.Fᵣ
— MethodFresnel reflection, needs pol_type soon too!
CanopyOptics.G
— MethodG(μ::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])LD
anAbstractLeafDistribution
type struct, includes a leaf distribution functionnLeg
an 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.49936166823681594
CanopyOptics.H
— MethodEq 46 in Shultis and Myneni
CanopyOptics.K
— MethodThe reduction factor proposed by Nilson and Kuusk, κ ≈ 0.1-0.3, returns exp(-κ * tan(abs(α))
CanopyOptics.abs_components
— MethodGiven a complex number, return the complex number with absolute value of both parts
CanopyOptics.abs_imag_only
— MethodGiven a complex number, return the complex number with same real part and absolute value of imaginary part
CanopyOptics.afsal
— MethodCalculate forward scattering amplitudes
CanopyOptics.asal
— MethodCalculate scattering coefficients for a thin disk
CanopyOptics.bfG
— MethodBrute Force G calculation (for testing
CanopyOptics.calctav
— Methodcalctav(α::FT, nr::FT) where {FT<:AbstractFloat}
Computes transmission of isotropic radiation across a dielectric surface (Stern F., 1964; Allen W.A.,Appl. Opt., 3(1):111-113 1973)). From calctav.m in PROSPECT-D
Arguments
α
angle of incidence [degrees]nr
Index of refraction
CanopyOptics.compute_Z_matrices
— Methodcompute_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]LD
aAbstractLeafDistribution
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.compute_lambertian_Γ
— Methodcompute_lambertian_Γ(μ::Array{FT,1},μꜛ::Array{FT,1}, r,t, LD::AbstractLeafDistribution; nLeg = 20)
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.
Arguments
μ::Array{FT,1}
: Quadrature points incoming direction (cos(θ))μꜛ::Array{FT,1}
: Quadrature points outgoing direction (cos(θ))r
: Leaf lambertian reflectancet
: Leaf lambertian transmittanceLD
aAbstractLeafDistribution
struct that describes the leaf angular distribution function.nLeg = 20
: number of quadrature points used for integration over all leaf angles (default is 20).
CanopyOptics.compute_Ψ
— MethodEq 45 in Shultis and Myneni, fixed grid in μ for both μ and μ'
CanopyOptics.compute_Ψ
— MethodEq 45 in Shultis and Myneni, fixed grid in μ for both μ and μ'
CanopyOptics.createLeafOpticalStruct
— MethodcreateLeafOpticalStruct(λ_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.createLeafOpticalStruct
— MethodcreateLeafOpticalStruct()
As in createLeafOpticalStruct(λ_bnds) but reads in the in Prospect-PRO database at original resolution (400-2500nm in 1nm steps)
CanopyOptics.cylinder
— MethodCompute the coefficients of the series expansion of cylinder scattering amplitude.
Equation 24 in Electromagnetic Wave Scattering from Some Vegetation Samples (KARAM, FUNG, AND ANTAR, 1988) http://www2.geog.ucl.ac.uk/~plewis/kuusk/karam.pdf
CanopyOptics.dielectric
— Methoddielectric(mod::SoilMW, T::FT,f::FT)
Computes the complex dieletric of soil (Ulaby & Long book)
Arguments
mod
anSoilMW
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.dielectric
— Methoddielectric(mod::LiquidPureWater, T::FT,f::FT)
Computes the complex dieletric of liquid pure walter (Ulaby & Long book)
Arguments
mod
anLiquidSaltWater
type structT
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.dielectric
— Methoddielectric(mod::LiquidSaltWater, T::FT,f::FT)
Computes the complex dieletric of liquid salty walter (Ulaby & Long book)
Arguments
mod
anLiquidSaltWater
type struct, provide Salinity in PSUT
Temperature in[⁰K]
f
Frequency in[GHz]
S
Salinity in[PSU]
comes from themod
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.dielectric
— Methoddielectric(mod::PureIce, T::FT,f::FT)
Computes the complex dieletric of liquid pure walter (Ulaby & Long book)
Arguments
mod
anPureIce
type structT
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.erectophile_leaves
— Functionerectophile_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.funcm
— Method3.1.2 in SAATCHI, MCDONALD
CanopyOptics.grdoh
— MethodOh et al, 1992 semi-emprical model for rough surface scattering
Equation 3.1.5 in Coherent Effects in Microwave Backscattering Models for Forest Canopies (SAATCHI, MCDONALD) https://www.researchgate.net/publication/3202927Semi-empiricalmodelofthe_ ensemble-averageddifferentialMuellermatrixformicrowavebackscatteringfrombaresoilsurfaces
CanopyOptics.plagiophile_leaves
— Functionplagiophile_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
— Functionplanophile_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.prospect
— Methodprospect(leaf::LeafProspectProProperties{FT},
optis) where {FT<:AbstractFloat}
Computes leaf optical properties (reflectance and transittance) based on PROSPECT-PRO
Arguments
leaf
LeafProspectProProperties
type struct which provides leaf compositionoptis
LeafOpticalProperties
type 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.scattering
— MethodCalculate scattering cross-section of a finite length cylinder with specified orientation
CanopyOptics.scattering_amplitude
— MethodCalculate bistatic scattering amplitude in the cylinder coordinate
Equations 26, 27 in Electromagnetic Wave Scattering from Some Vegetation Samples (KARAM, FUNG, AND ANTAR, 1988) http://www2.geog.ucl.ac.uk/~plewis/kuusk/karam.pdf
CanopyOptics.scattering_model_lambertian
— MethodEq 37 in Shultis and Myneni
CanopyOptics.spherical_leaves
— Functionspherical_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
— Functionuniform_leaves(FT=Float64)
Creates a LeafDistribution
following a uniform distribution (all angles equally likely)
CanopyOptics.unit2nm
— MethodConvert an input array with Spectral units (or none) to a grid in nm
CanopyOptics.wood_backward
— MethodCalculation of average backscatter cross-section of a finite length cylinder, exact series solution (KARAM, FUNG, AND ANTAR, 1988)
CanopyOptics.wood_forward
— MethodCalculation of average forward scattering amplitudes of a finite length cylinder, exact series solution
CanopyOptics.zeta
— MethodCalculate zeta function using bessel functions
Equation 26 in Electromagnetic Wave Scattering from Some Vegetation Samples (KARAM, FUNG, AND ANTAR, 1988) http://www2.geog.ucl.ac.uk/~plewis/kuusk/karam.pdf
CanopyOptics.βparameters
— MethodGet Beta Distribution parameters using standard tabulated values θₗ,s, Eq. 2.9 and 2.10 in Bonan