CanopyLayers
Structures
Land.CanopyLayers.Canopy4RT
— Typemutable struct Canopy4RT
A canopy struct for the radiation transfer module
Fields
nLayer::Int64
: Number of canopy layers
LAI::AbstractFloat
: Leaf Area Index
Ω::AbstractFloat
: Clumping factor
clump_a::AbstractFloat
: Structure factor a
clump_b::AbstractFloat
: Structure factor b
leaf_width::AbstractFloat
: Leaf width
hc::AbstractFloat
: Vegetation height
LIDFa::AbstractFloat
: Leaf Inclination
LIDFb::AbstractFloat
: Variation in leaf inclination
hot::AbstractFloat
: HotSpot parameter (still need to check!)
height::AbstractFloat
: Canopy height [m]
z0m::AbstractFloat
: Canopy roughness [m]
z0h::AbstractFloat
: Tree roughtnes [m]
d::AbstractFloat
: Canopy displacement height [m]
Cd::AbstractFloat
: m/sqrt(s) turbulent transfer coefficient
litab::Vector{FT} where FT<:AbstractFloat
: List of mean inclination angles [°]
litab_bnd::Matrix{FT} where FT<:AbstractFloat
: List of inclination angle boundaries [°]
lazitab::Vector{FT} where FT<:AbstractFloat
: List of mean azimuth angles [°]
cos_ttlo::Vector{FT} where FT<:AbstractFloat
: Cosine of lazitab
cos_philo::Vector{FT} where FT<:AbstractFloat
: Cosine of lazitab - raa (relative azimuth angle), update with time
cos_ttli::Vector{FT} where FT<:AbstractFloat
: Cosine of litab
sin_ttli::Vector{FT} where FT<:AbstractFloat
: Sine of litab
vol_scatt::Vector{FT} where FT<:AbstractFloat
: Cache for volome scatter function
lidf::Vector{FT} where FT<:AbstractFloat
: Inclination angles weight distribution
xl::Vector{FT} where FT<:AbstractFloat
: List of level location (level = layer + 1)
dx::AbstractFloat
: 1/nLayer
nAzi::Int64
: Number of azimuth angles
nIncl::Int64
: Number of inclination angles
Land.CanopyLayers.CanopyOpticals
— Typemutable struct CanopyOpticals{FT}
A struct for canopy optical properties
Fields
nAzi::Int64
: Number of azimuth angles
nIncl::Int64
: Number of inclination agles
nLayer::Int64
: Number of canopy layers
nWL::Int64
: Number of wave lengths
sdb::Any
: Solar -> Diffuse backscatter weight
sdf::Any
: Solar -> Diffuse forward scatter weight
dob::Any
: Diffuse -> Directional backscatter weight
dof::Any
: Diffuse -> Directional forward scatter weight
ddb::Any
: Diffuse -> Diffuse backscatter weight
ddf::Any
: Diffuse -> Diffuse forward scatter weight
ks::Any
: Solar beam extinction coefficient weight
ko::Any
: Outgoing beam extinction coefficient weight
bf::Any
: ?
sob::Any
: Weight of specular2directional backscatter coefficient
sof::Any
: Weight of specular2directional forward coefficient
Ps::Vector
: Probability of directly viewing a leaf in solar direction
Po::Vector
: Probability of directly viewing a leaf in viewing direction
Pso::Vector
: Bi-directional probability of directly viewing a leaf (solar->canopy->viewing)
fs::Matrix
: conversion factor fs to compute irradiance on inclined leaf
absfs::Matrix
: abs(fs)
absfsfo::Matrix
: abs(fs*fo)
fsfo::Matrix
: fs*fo
fo::Matrix
: conversion factor fo for angle towards observer (not sun like fs)
absfo::Matrix
: abs(fo)
cosΘ_l::Matrix
: Cosine of leaf azimuths
cos2Θ_l::Matrix
: cos of leaf azimuth sqared
sigb::Matrix
: diffuse backscatter scattering coefficient for diffuse incidence
sigf::Matrix
: diffuse forward scattering coefficient for diffuse incidence
sb::Matrix
: diffuse backscatter scattering coefficient for specular incidence
sf::Matrix
: diffuse forward scattering coefficient for specular incidence
vb::Matrix
: directional backscatter scattering coefficient for diffuse incidence
vf::Matrix
: directional forward scattering coefficient for diffuse incidence
w::Matrix
: bidirectional scattering coefficent (directional-directional)
a::Matrix
: attenuation
Xsd::Matrix
: Effective layer transmittance (direct->diffuse)
Xdd::Matrix
: Effective layer transmittance (diffuse->diffuse)
R_sd::Matrix
: Effective layer reflectance (direct->diffuse)
R_dd::Matrix
: Effective layer reflectance (diffuse->diffuse)
Es_::Matrix
: Solar direct radiation per layer)
Land.CanopyLayers.CanopyRads
— Typemutable struct CanopyRads{FT}
A struct for canopy radiation information
Fields
nAzi::Int64
: Number of azimuth angles
nIncl::Int64
: Number of inclination agles
nLayer::Int64
: Number of canopy layers
nLevel::Int64
: Number of canopy levels
nWL::Int64
: Number of wave lengths
nWLF::Int64
: Number of wave lengths for SIF
intEout::Any
: Integrated TOC outgoing flux [W m⁻²]
incomingPAR::Any
: Incident spectrally integrated total PAR [mol m⁻² s⁻¹]
incomingPAR_direct::Any
: Incident spectrally integrated direct PAR [mol m⁻² s⁻¹]
incomingPAR_diffuse::Any
: Incident spectrally integrated diffuse PAR [mol m⁻² s⁻¹]
RnSoil_diffuse::Any
: Net radiation of shaded soil [W m⁻²]
RnSoil_direct::Any
: Net Short-wave radiation of sunlit soil [W m⁻²]
RnSoil::Any
: Net Short-wave radiation of soil (shaded + sunlit) [W m⁻²]
RnSoilLW::Any
: Net long-wave radiation of soil (shaded + sunlit) [W m⁻²]
absPAR_shade::Vector
: Net PAR of shaded leaves [mol m⁻² s⁻¹]
absPAR_shadeCab::Vector
: Net PAR by Cab+Car of shaded leaves [moles m⁻² s⁻¹]
intNetSW_sunlit::Vector
: Spectrally integrated net absorbed direct radiation in each layer per leaf area [W m⁻²]
intNetSW_shade::Vector
: Spectrally integrated net absorbed diffuse radiation in each layer per leaf area [W m⁻²]
intNetLW_sunlit::Vector
: Spectrally integrated net absorbed direct radiation in each layer per leaf area [W m⁻²]
intNetLW_shade::Vector
: Spectrally integrated net absorbed diffuse radiation in each layer per leaf area [W m⁻²]
T_sun::Vector
: Leaf temperature (sunlit) [K]
T_shade::Vector
: Leaf temperature (shaded) [K]
ϕ_shade::Vector
: Fluorescence yield for shaded leaves
H_shade::Vector
: Sensible Heat flux H of shaded leaves [W m⁻²]
LE_shade::Vector
: Latent Heat flux LE of shaded leaves [W m⁻²]
NPQ_shade::Vector
: NPQ of shaded leaves
GPP_shade::Vector
: GPP of shaded leaves [μmol m⁻² s⁻¹]
gs_shade::Vector
: gs of shaded leaves [mol m⁻² s⁻¹]
ψl_shade::Vector
: Leaf water potential of shaded leaves [MPa]
Cc_shade::Vector
: Cc of shaded leaves [µmol/mol]
Pi_shade::Vector
: internal CO₂ concentration of shaded leaves [µmol/mol]
Lo::Vector
: Short-wave TOC outgoing radiance in observation direction [mW m⁻² nm⁻¹ sr⁻¹]
Eout::Vector
: Short-wave TOC outgoing radiation [mW m⁻² nm⁻¹]
alb_obs::Vector
: Short-wave Albedo in viewing direction
alb_direct::Vector
: Short-wave Albedo for direct incoming radiation
alb_diffuse::Vector
: Short-wave Albedo for diffuse incoming radiation
E_up::Matrix
: Upwelling diffuse short-wave radiation within canopy [mW m⁻² nm⁻¹]
E_down::Matrix
: Downwelling diffuse short-wave radiation within canopy [mW m⁻² nm⁻¹]
netSW_sunlit::Matrix
: Net absorbed direct radiation in each layer [mW m⁻² nm⁻¹]
netSW_shade::Matrix
: net absorbed diffuse radiation in each layer [mW m⁻² nm⁻¹]
absPAR_sun::Array{FT, 3} where FT
: net PAR of sunlit leaves [mol m⁻² s⁻¹]
absPAR_sunCab::Array{FT, 3} where FT
: net PAR by Cab+Car of sunlit leaves [mol m⁻² s⁻¹]
T_sun3D::Array{FT, 3} where FT
: Leaf temperature (sunlit) [K]
ϕ_sun::Array{FT, 3} where FT
: Fluorescence yield for sunlit leaves
H_sun::Array{FT, 3} where FT
: Sensible Heat flux H of sunlit leaves [W m⁻²]
LE_sun::Array{FT, 3} where FT
: Latent Heat flux LE of sunlit leaves [W m⁻²]
NPQ_sun::Array{FT, 3} where FT
: NPQ of sunlit leaves
GPP_sun::Array{FT, 3} where FT
: GPP of sunlit leaves [μmol m⁻² s⁻¹]
gs_sun::Array{FT, 3} where FT
: gs of sunlit leaves [mol m⁻² s⁻¹]
ψl_sun::Array{FT, 3} where FT
: Leaf water potential of sunlit leaves [MPa]
Cc_sun::Array{FT, 3} where FT
: Cc of sunlit leaves [µmol/mol]
Pi_sun::Array{FT, 3} where FT
: Internal CO₂ concentration of sunlit leaves [µmol/mol]
SIF_hemi::Vector
: Hemispheric total outgoing SIF flux [mW m⁻² nm⁻¹]
)
SIF_obs::Vector
: Observer-direction outgoing SIF radiance (mW m⁻² nm⁻¹ sr⁻¹))
SIF_obs_sunlit::Vector
: Observer-direction outgoing SIF radiance, sunlit leaves (mW m⁻² nm⁻¹ sr⁻¹)
SIF_obs_shaded::Vector
: Observer-direction outgoing SIF radiance, shaded leaves (mW m⁻² nm⁻¹ sr⁻¹)
SIF_obs_scattered::Vector
: Observer-direction outgoing SIF radiance, scattered (mW m⁻² nm⁻¹ sr⁻¹)
SIF_obs_soil::Vector
: Observer-direction outgoing SIF radiance, soil-reflected (mW m⁻² nm⁻¹ sr⁻¹)
SIF_sum::Vector
: Total SIF sum of layer sources [mW m⁻² nm⁻¹]
)
Land.CanopyLayers.IncomingRadiation
— Typemutable struct IncomingRadiation{FT}
Incoming radiation information
Fields
E_direct::Vector
: Direct incoming radiation [mW m⁻² nm⁻¹]
E_diffuse::Vector
: Diffuse incoming radiation [mW m⁻² nm⁻¹]
Land.CanopyLayers.LeafBios
— Typemutable struct LeafBios{FT}
A struct of leaf biological parameters
Fields
nWL::Int64
: Number of wave length
nWLE::Int64
: Number of wave length for excitation
nWLF::Int64
: Number of wave length for SIF
N::Any
: Leaf structure parameter
Cab::Any
: Chlorophyll a+b content [µg cm⁻²]
Car::Any
: Carotenoid content [µg cm⁻²]
Ant::Any
: Anthocynanin content [µg cm⁻²]
Cs::Any
: Senescent material fraction
Cw::Any
: Equivalent water thickness [cm]
Cm::Any
: Dry matter content (dry leaf mass per unit area) [g cm⁻²]
Cx::Any
: Fractionation between Zeaxanthin and Violaxanthin in Car (1=all Zeaxanthin) (-)
fqe::Any
: Leaf fluorescence efficiency (Fo standard)
ρ_LW::Any
: Broadband thermal reflectance (-)
τ_LW::Any
: Broadband thermal transmission (-)
ρ_SW::Vector
: Shortwave leaf reflectance
τ_SW::Vector
: Shortwave leaf transmission
α_SW::Vector
: Shortwave absorption
kChlrel::Vector
: Relative absorbtion by Chlorophyll+Car
kChlrel_old::Vector
: Relative absorbtion by Chlorophyll
Mb::Matrix
: Fluorescence excitation matrix backwards
Mf::Matrix
: Fluorescence excitation matrix forwards
ndub::Int64
: Doubling adding layers
Land.CanopyLayers.LeafOpticals
— Typemutable struct LeafOpticals{FT}
Struct for leaf optical properties
Fields
nr::Vector
Km::Vector
Kab::Vector
Kant::Vector
Kcar::Vector
Kw::Vector
KBrown::Vector
phi::Vector
KcaV::Vector
KcaZ::Vector
lambda::Vector
: Wave length [nm]
, same as WL
in WaveLengths
`
Land.CanopyLayers.RTDimensions
— Typemutable struct RTDimensions
Struct that stores matrix dimension information
Fields
nAzi::Int64
: Number of azimuth angles
nIncl::Int64
: Number of inclination agles
nLayer::Int64
: Number of canopy layers
nLevel::Int64
: Number of canopy layer boundaries nLayer+1
nPAR::Int64
: Number of PAR wave lengths
nWL::Int64
: Number of wave lengths
nWLE::Int64
: Number of wave length for excitation
nWLF::Int64
: Number of wave lengths for SIF
Land.CanopyLayers.SoilOpticals
— Typemutable struct SoilOpticals{FT}
A struct of soil optical parameters
Fields
T::Any
: Soil surface temperature
color::Int64
: Soil color class
ρ_NIR::Any
: Shortwave albedo for NIR
ρ_PAR::Any
: Shortwave albedo for PAR
ρ_SW::Vector
: Shortwave albedo that matches WL
from WaveLengths
ρ_SW_SIF::Vector
: Shortwave albedo that matches WLF
from WaveLengths
ε_SW::Vector
: Shortwave absorption that equals 1 - ρ_SW
SW_mat_raw_4::Matrix
: Shortwave albedo matrix from 4 bands, wavelengths are 400:10:2500 nm
SW_mat_raw_2::Matrix
: Shortwave albedo matrix from 2 bands, wavelengths are 400:10:2500 nm
SW_mat_4::Matrix
: Shortwave albedo matrix from 4 bands with WL
from WaveLengths
SW_mat_2::Matrix
: Shortwave albedo matrix from 2 bands with WL
from WaveLengths
SW_vec_4::Vector
: Shortwave albedo weight from 4 bands
SW_vec_2::Vector
: Shortwave albedo weight from 2 bands
ρ_SW_raw::Vector
: Shortwave albedo, wavelengths are 400:10:2500 nm
ρ_LW::Vector
: Longtwave albedo
dry_NIR::Any
: Mean value for day band 1 in NIR region
dry_PAR::Any
: Mean value for day band 1 in PAR region
wet_NIR::Any
: Mean value for day band 1 in NIR region
wet_PAR::Any
: Mean value for day band 1 in PAR region
Land.CanopyLayers.SolarAngles
— Typestruct SolarAngles{FT}
Struct for observation and solar angles
Fields
hza::Any
: Hill zenith angle [°]
, hill slope angle
haa::Any
: Hill azimuth angle [°]
, 0 for North, 180 for south
saa::Any
: Solar azimuth angle [°]
, a function of time
sza::Any
: Solar zenith angle [°]
, a function of lat and time
vaa::Any
: Viewing azimuth angle [°]
vza::Any
: Viewing zenith angle [°]
raa::Any
: Relative azimuth angle [°]
, difference between saa and vaa
Land.CanopyLayers.WaveLengths
— Typemutable struct WaveLengths{FT}
Struct for pre-set wave length parameters
Fields
minwlPAR::Any
: Minimal WL for PAR [nm]
maxwlPAR::Any
: Maximal WL for PAR [nm]
minwlNIR::Any
: Minimal WL for NIR [nm]
maxwlNIR::Any
: Maximal WL for NIR [nm]
minwle::Any
: Minimal WL for SIF excitation [nm]
maxwle::Any
: Maximal WL for SIF excitation [nm]
minwlf::Any
: Minimal WL for SIF emission/fluorescence [nm]
maxwlf::Any
: Maximal WL for SIF emission/fluorescence [nm]
sWL::Vector
: Standard wave length [nm]
dWL::Vector
: Differential wavelength
optis::Land.CanopyLayers.LeafOpticals
: Leaf optical parameter set
WL::Vector
: Wave length [nm]
iWLE::Vector{Int64}
: Index of WLE in WL
iWLF::Vector{Int64}
: Index of WLF in WL
iPAR::Vector{Int64}
: index of wlPAR in WL
iNIR::Vector{Int64}
: index of wlNIR in WL
WLE::Vector
: excitation wave length [nm]
WLF::Vector
: Fluorescence wave length [nm]
WL_iPAR::Vector
: Wave length for PAR
dWL_iPAR::Vector
: Differential wave length for PAR
dWL_iWLE::Vector
: Differential wave length for iWLE
nPAR::Int64
: Length of WL_iPAR
nWL::Int64
: Length of WL
nWLE::Int64
: length of WLE
nWLF::Int64
: length of WLF
Caches
Land.CanopyLayers.CFCache
— Typemutable struct CFCache{FT}
Cache to speed canopy_fluxes!
by pre-allocating arrays
Fields
abs_wave::Vector
: absorbed energy from wave lengths
absfs_lidf::Vector
: absfs' * lidf [nAzi]
E_all::Vector
: wave length energy [same as dWL]
E_iPAR::Vector
: wave length energy [same as iPAR]
lPs::Vector
: lPs [nLayer]
kChlrel::Vector
: kChlrel [same as iPAR]
PAR_diff::Vector
: diffusive PAR [same as iPAR]
PAR_diffCab::Vector
: diffusive PAR for photosynthesis [same as iPAR]
PAR_dir::Vector
: direct PAR [same as iPAR]
PAR_dirCab::Vector
: diffusive PAR for photosynthesis [same as iPAR]
Land.CanopyLayers.create_cf_cache
— Functioncreate_cf_cache(FT, rt_dim::RTDimensions)
Create a CFCache
type struct, given
FT
Floating number typert_dim
RTDimensions
type struct
Land.CanopyLayers.CGCache
— Typemutable struct CGCache{FT}
Cache to speed canopy_geometry!
by pre-allocating arrays
Fields
_Co::Vector
: cos_ttli .* cos(vza) dim: nIncl
_Cs::Vector
: cos_ttli .* cos(sza) dim: nIncl
_So::Vector
: sin_ttli .* sin(vza) dim: nIncl
_Ss::Vector
: sin_ttli .* sin(sza) dim: nIncl
_1s::Matrix
: maxtrix filled with 1 dim: (1, nAzi)
_2d::Matrix
: 2D array to speed up _cds and _cdo dim: (nIncl, nAzi)
_cdo::Matrix
: Co * _1s .+ _So * cosphilo' dim: (nIncl, nAzi)
_cds::Matrix
: Cs * _1s .+ _Ss * costtlo' dim: (nIncl, nAzi)
Land.CanopyLayers.create_cg_cache
— Functioncreate_cg_cache(FT, rt_dim::RTDimensions)
Create a CGCache
type struct, given
FT
Floating number typert_dim
RTDimensions
type struct
Land.CanopyLayers.SFCache
— Typemutable struct SFCache{FT}
Cache to speed SIF_fluxes!
by pre-allocating arrays
Fields
M⁻_sun::Vector
M⁺_sun::Vector
wfEs::Vector
sfEs::Vector
sbEs::Vector
M⁺⁻::Vector
M⁺⁺::Vector
M⁻⁺::Vector
M⁻⁻::Vector
sun_dwl_iWlE::Vector
tmp_dwl_iWlE::Vector
ϕ_cosΘ_lidf::Vector
vfEplu_shade::Vector
vbEmin_shade::Vector
vfEplu_sun::Vector
vbEmin_sun::Vector
sigfEmin_shade::Vector
sigbEmin_shade::Vector
sigfEmin_sun::Vector
sigbEmin_sun::Vector
sigfEplu_shade::Vector
sigbEplu_shade::Vector
sigfEplu_sun::Vector
sigbEplu_sun::Vector
zeroB::Vector
tmp_1d_nWlF::Vector
tmp_1d_nLayer::Vector
dnorm::Vector
τ_dd::Matrix
: transmission of diffusive light?
ρ_dd::Matrix
: extinction of diffuse light?
Xdd::Matrix
Rdd::Matrix
Y::Matrix
U::Matrix
S⁻::Matrix
S⁺::Matrix
piLs::Matrix
piLd::Matrix
Fsmin::Matrix
Fsplu::Matrix
Fdmin::Matrix
Fdplu::Matrix
Femo::Matrix
M⁺::Matrix
M⁻::Matrix
ϕ_cosΘ::Matrix
F⁻::Matrix
F⁺::Matrix
net_diffuse::Matrix
tmp_2d_nWlF_nLayer::Matrix
tmp_2d_nWlF_nLayer_2::Matrix
Land.CanopyLayers.create_sf_cache
— Functioncreate_sf_cache(FT, rt_dim::RTDimensions)
Create a SFCache
type struct, given
FT
Floating number typert_dim
RTDimensions
type struct
Land.CanopyLayers.SWCache
— Typemutable struct SWCache{FT}
Cache to speed short_wave!
by pre-allocating arrays
Fields
dnorm::Vector
: dnorm?
piLo::Vector
: pi * Lo
piLoc::Vector
: pi * Lo from canopy
piLos::Vector
: pi * Lo from soil
piLoc2::Matrix
: pi * Lo from canopy 2D matrix
ρ_dd::Matrix
: extinction of diffuse light?
ρ_sd::Matrix
: extinction of direct light?
τ_dd::Matrix
: transmission of diffusive light?
τ_sd::Matrix
: transmission of direct light?
Land.CanopyLayers.create_sw_cache
— Functioncreate_sw_cache(FT, rt_dim::RTDimensions)
Create a CGCache
type struct, given
FT
Floating number typert_dim
RTDimensions
type struct
Land.CanopyLayers.RTCache
— Typemutable struct RTCache{FT}
Collection of caches to speed up RT module
Fields
cf_con::Land.CanopyLayers.CFCache
: CFCache
type cache
cg_con::Land.CanopyLayers.CGCache
: CGCache
type cache
sf_con::Land.CanopyLayers.SFCache
: SFCache
type cache
sw_con::Land.CanopyLayers.SWCache
: SWCache
type cache
Land.CanopyLayers.create_rt_cache
— Functioncreate_rt_cache(FT, rt_dim::RTDimensions)
Create an RTCache
, given
FT
Floating number typert_dim
RTDimensions
type struct
Initialization of Structures
Land.CanopyLayers.create_canopy_opticals
— Functioncreate_canopy_opticals( FT, rt_dim::RTDimensions)
Create a CanopyOpticals
struct, given
FT
Floating number typert_dim
RTDimensions
type struct
Land.CanopyLayers.create_canopy_rads
— Functioncreate_canopy_rads(FT, rt_dim::RTDimensions)
Create a CanopyRads
struct, given
FT
Floating number typert_dim
RTDimensions
type struct
Land.CanopyLayers.create_canopy_rt
— Functioncreate_canopy_rt(FT; nLayer::Int = 20, LAI::Number = FT(3))
Create Canopy4RT
, given
FT
Floating number typenLayer
Number of canopy layersLAI
Leaf area index
Land.CanopyLayers.create_incoming_radiation
— Functioncreate_incoming_radiation(
wls::WaveLengths{FT},
wlfn::String = FILE_SUN
) where {FT<:AbstractFloat}
Create an AbstractIncomingRadiation
struct, given
wls
WaveLengths
type structwlfn
File that saves incoming wave information
Land.CanopyLayers.create_leaf_bios
— Functioncreate_leaf_bios(FT, rt_dim::RTDimensions)
Create a LeafBios
type struct, given
FT
Floating number typert_dim
RTDimensions
type struct
Land.CanopyLayers.create_leaf_opticals
— Functioncreate_leaf_opticals(
sWL::Array{FT,1},
file::String = OPTI_2021
) where {FT<:AbstractFloat}
Create an AbstractLeafOptiPara
struct, given
sWL
Standard wave lengthopfn
File that saves optical parameters
Land.CanopyLayers.create_rt_dims
— Functioncreate_rt_dims(
can::Canopy4RT{FT},
wls::WaveLengths{FT}
) where {FT<:AbstractFloat}
Create RTDimensions
, given
can
Canopy4RT
type structwls
WaveLengths
type struct
Land.CanopyLayers.create_wave_length
— Functioncreate_wave_length(FT, sWLs = [collect(400.0:10.0: 650.1); collect(655.0: 5.0: 770.1); collect(780.0:25.0:2400.1)];
max_NIR::Number = 2500, max_PAR::Number = 700, min_NIR::Number = 700, min_PAR::Number = 400, opti_file::String = OPTI_2021)
Create WaveLengths
type struct, given
FT
Floating number typesWLs
Shortwave wavelength binsmax_NIR
Maximal NIR wavelengthmax_PAR
Maximal PAR wavelengthmin_NIR
Minimal NIR wavelengthmin_PAR
Minimal PAR wavelengthopti_file
Input reference optical file path
Land.CanopyLayers.initialize_rt_module
— Functioninitialize_rt_module(FT; nLayer::Int = 20, LAI::Number = FT(3))
Initialize the RT module and return the sturctures, given
FT
Floating number typenLayer
Number of canopy layersLAI
Leaf area index
This function initializes and returns
angles
SolarAngles
can
Canopy4RT
can_opt
CanopyOpticals
can_rad
CanopyRads
in_rad
IncomingRadiation
leaves
Array{LeafBios
,1}rt_con
RTCache
rt_dim
RTDimensions
soil
SoilOpticals
wls
WaveLengths
Big Leaf Model
Land.CanopyLayers.big_leaf_partition
— Functionbig_leaf_partition(
lai::FT,
zenith::FT,
r_all::FT,
r_dir::FT = FT(0.8)
) where {FT<:AbstractFloat}
Partition the big-leaf canopy into sunlit and shaded layers, given
lai
Leaf area indexzenith
Zenith angle in degreer_all
Total radiation in[W m⁻²]
r_dir
Direct radiation partition inr_all
The function returns
ratio
ratio of sunlit leaves out of all leavesq_slm
Mean sunlit layer PARq_shm
Mean shaded layer PARe_sl
Mean sunlit layer absorbed total energye_sh
Mean shaded layer absorbed total energy
SCOPE Model
Land.CanopyLayers.canopy_fluxes!
— Functioncanopy_fluxes!(
can::Canopy4RT{FT},
can_opt::CanopyOpticals{FT},
can_rad::CanopyRads{FT},
in_rad::IncomingRadiation{FT},
soil::SoilOpticals{FT},
leaves::Array{LeafBios{FT},1},
wls::WaveLengths{FT},
rt_con::RTCache{FT}
) where {FT<:AbstractFloat}
Computes a variety of integrated fluxes from the spectrally resolved computations in the short-wave Canopy RT (e.g. absorbed soil radiation, absorbed direct and diffuse PAR by layer (and angles for direct), net direct and diffuse energy balance per layer), given
can
Canopy4RT
type structcan_opt
CanopyOpticals
type structcan_rad
CanopyRads
type structin_rad
IncomingRadiation
type structsoil
SoilOpticals
type structleaves
Array ofLeafBios
type structwls
WaveLengths
type structrt_con
RTCache
type cache
Land.CanopyLayers.canopy_geometry!
— Functioncanopy_geometry!(
can::Canopy4RT{FT},
angles::SolarAngles{FT},
can_opt::CanopyOpticals{FT},
rt_con::RTCache{FT}
) where {FT<:AbstractFloat}
Computes canopy optical properties (extinction coefficients for direct and diffuse light) based on the SAIL model. Most important input parameters are leaf inclination and azimuth distribution functions and sun-sensor geometry. Canopy clumping Ω is implemented as in Pinty et al (2015), given
can
Canopy4RT
type structangles
SolarAngles
type structcan_opt
CanopyOpticals
type structrt_con
RTCache
type cache
Land.CanopyLayers.canopy_matrices!
— Functioncanopy_matrices!(
leaves::Array{LeafBios{FT},1},
can_opt::CanopyOpticals{FT}
) where {FT<:AbstractFloat}
Compute scattering coefficient matrices for direct and diffuse light given geometry dependent overall extinction coefficients and pigment dependent leaf reflectance and transmission (computed via fluspect). This function has to be called before short_wave!
can be used.
leaves
Array ofLeafBios
type structcan_opt
CanopyOpticals
type struct
Land.CanopyLayers.clumping_factor!
— Functionclumping_factor!(
can::Canopy4RT{FT},
angles::SolarAngles{FT}
) where {FT<:AbstractFloat}
Calculate the clumping factor, given
can
Canopy4RT
type structangles
SolarAngles
type struct
Land.CanopyLayers.diffusive_S
— Functiondiffusive_S(τ_dd::Array{FT},
ρ_dd::Array{FT},
S⁻::Array{FT},
S⁺::Array{FT},
boundary_top::Array{FT},
boundary_bottom::Array{FT},
rsoil::Array{FT}
) where {FT<:AbstractFloat}
Computes 2-stream diffusive radiation transport (used for thermal and SIF) given:
τ_dd
A 2D Array with layer reflectancesρ_dd
A 2D Array with layer transmissionsS⁻
A 2D Array with layer source terms in the downwelling directionS⁺
A 2D Array with layer source terms in the upwelling directionboundary_top
A 1D array with downwelling radiation at the top (top of canopy)boundary_bottom
A 1D array with upwnwelling radiation at the bottom (soil)rsoil
A 1D array with soil reflectance
Land.CanopyLayers.diffusive_S!
— Functiondiffusive_S!(
sf_con::SFCache{FT},
soil::SoilOpticals{FT},
rt_dim::RTDimensions
) where {FT<:AbstractFloat}
Computes 2-stream diffusive radiation transport (used for thermal and SIF), given
sf_con
SFCache
type cachesoil
SoilOpticals
type structrt_dim
RTDimensions
type struct
Land.CanopyLayers.fluspect!
— Functionfluspect!(leaf::LeafBios{FT},
wls::WaveLengths{FT};
APAR_car::Bool = true
) where {FT<:AbstractFloat}
Computes leaf optical properties (reflectance and transittance) based on pigment concentrations. Also computes Fluorescence excitation matrices. Mostly based on PROSPECT-D for leaf reflectance/transmission and FluSpec for fluorescence.
leaf
LeafBios
type structwls
WaveLengths
type structAPAR_car
If true, include Car absorption in APAR for photosynthesis
Land.CanopyLayers.short_wave!
— Functionshort_wave!(can::Canopy4RT{FT},
can_opt::CanopyOpticals{FT},
can_rad::CanopyRads{FT},
in_rad::IncomingRadiation{FT},
soil::SoilOpticals{FT},
rt_con::RTCache{FT}
) where {FT<:AbstractFloat}
Simulate the short wave radiation through the canopy, given
can
Canopy4RT
type structcan_opt
CanopyOpticals
type structcan_rad
CanopyRads
type structin_rad
IncomingRadiation
type structsoil
SoilOpticals
type structrt_con
RTCache
type cache
Land.CanopyLayers.SIF_fluxes!
— FunctionSIF_fluxes!(leaves::Array{LeafBios{FT},1},
can_opt::CanopyOpticals{FT},
can_rad::CanopyRads{FT},
can::Canopy4RT{FT},
soil::SoilOpticals{FT},
wls::WaveLengths{FT},
rt_con::RTCache{FT},
rt_dim::RTDimensions;
photon::Bool = true
) where {FT<:AbstractFloat}
Computes 2-stream diffusive radiation transport for SIF radiation (calls [diffusive_S!
] internally). Layer reflectance and transmission is computed from SW optical properties, layer sources from absorbed light and SIF efficiencies. Boundary conditions are zero SIF incoming from atmosphere or soil.
leaves
Array ofLeafBios
type structcan_opt
CanopyOpticals
type structcan_rad
CanopyRads
type structcan
Canopy4RT
type structsoil
SoilOpticals
type structwls
WaveLengths
type structrt_con
RTCache
type cachert_dim
RTDimensions
type structphoton
If true, use photon unit in the matrix conversion
SIF_fluxes!(leaf::LeafBios{FT},
in_rad::IncomingRadiation{FT},
wls::WaveLengths{FT},
rt_con::RTCache{FT},
fqe::FT = FT(0.01);
photon::Bool = true
) where {FT<:AbstractFloat}
Leaf level SIF, given
leaf
LeafBios
type structin_rad
IncomingRadiation
type structwls
WaveLengths
type structrt_con
RTCache
type cachefqe
Fluorescence quantum yield (default at 1%)photon
If true, use photon unit in the matrix conversion
Note that in_rad
assumes direct light with zenith angle of 0, and a zenith angle correction needs to be made before passing it to this function. The up- and down-ward SIF are stored in sf_con
as M⁻_sun
and M⁺_sun
.
Land.CanopyLayers.thermal_fluxes!
— Functionthermal_fluxes!(
leaves::Array{LeafBios{FT},1},
can_opt::CanopyOpticals{FT},
can_rad::CanopyRads{FT},
can::Canopy4RT{FT},
soil::SoilOpticals{FT},
incLW::Array{FT},
wls::WaveLengths{FT}
) where {FT<:AbstractFloat}
Computes 2-stream diffusive radiation transport for thermal radiation (calls [diffusive_S
] internally). Layer reflectance and transmission is computed from LW optical properties, layer sources from temperature and Planck law, boundary conditions from the atmosphere and soil emissivity and temperature. Currently only uses Stefan Boltzmann law to compute spectrally integrated LW but can be easily adjusted to be spectrally resolved.
leaves
Array ofLeafBios
type structcan_opt
CanopyOpticals
type structcan_rad
CanopyRads
type structcan
Canopy4RT
type structsoil
SoilOpticals
type structincLW
A 1D array with incoming long-wave radiationwls
WaveLengths
type struct
Indicies
Land.CanopyLayers.BLUE
— FunctionBLUE(can_rad::CanopyRads{FT},
wls::WaveLengths{FT}
) where {FT<:AbstractFloat}
Return the BLUE @ 469 nm, given
can_rad
CanopyRads
type structwls
WaveLengths
type struct
Land.CanopyLayers.EVI
— FunctionEVI(can_rad::CanopyRads{FT},
wls::WaveLengths{FT}
) where {FT<:AbstractFloat}
Return the EVI, given
can_rad
CanopyRads
type structwls
WaveLengths
type struct
Land.CanopyLayers.EVI2
— FunctionEVI2(can_rad::CanopyRads{FT},
wls::WaveLengths{FT}
) where {FT<:AbstractFloat}
Return the EVI2, given
can_rad
CanopyRads
type structwls
WaveLengths
type struct
Land.CanopyLayers.LSWI
— FunctionLSWI(can_rad::CanopyRads{FT},
wls::WaveLengths{FT}
) where {FT<:AbstractFloat}
Return the LSWI, given
can_rad
CanopyRads
type structwls
WaveLengths
type struct
Land.CanopyLayers.NDVI
— FunctionNDVI(can_rad::CanopyRads{FT},
wls::WaveLengths{FT}
) where {FT<:AbstractFloat}
Return the NDVI, given
can_rad
CanopyRads
type structwls
WaveLengths
type struct
Land.CanopyLayers.NIR
— FunctionNIR(can_rad::CanopyRads{FT},
wls::WaveLengths{FT}
) where {FT<:AbstractFloat}
Return the NIR @ 858.5 nm, given
can_rad
CanopyRads
type structwls
WaveLengths
type struct
Land.CanopyLayers.NIRv
— FunctionNIRv(can_rad::CanopyRads{FT},
wls::WaveLengths{FT}
) where {FT<:AbstractFloat}
Return the NIRv, given
can_rad
CanopyRads
type structwls
WaveLengths
type struct
Land.CanopyLayers.RED
— FunctionRED(can_rad::CanopyRads{FT},
wls::WaveLengths{FT}
) where {FT<:AbstractFloat}
Return the RED @ 645 nm, given
can_rad
CanopyRads
type structwls
WaveLengths
type struct
Land.CanopyLayers.REF_WL
— FunctionREF_WL(wls::WaveLengths{FT},
can_rad::CanopyRads{FT}
wls::WaveLengths{FT},
twl::FT
) where {FT<:AbstractFloat}
Return the Reflectance, given
can_rad
CanopyRads
type structwls
WaveLengths
type structtwl
Target wave length in nm
Land.CanopyLayers.SIF_740
— FunctionSIF_740(can_rad::CanopyRads{FT},
wls::WaveLengths{FT}
) where {FT<:AbstractFloat}
Return the SIF @ 740 nm, given
can_rad
CanopyRads
type structwls
WaveLengths
type struct
Land.CanopyLayers.SIF_757
— FunctionSIF_757(can_rad::CanopyRads{FT},
wls::WaveLengths{FT};
oco::Int = 2
) where {FT<:AbstractFloat}
Return the SIF @ 758.68 (OCO2) or 758.77 (OCO3) nm, given
can_rad
CanopyRads
type structwls
WaveLengths
type structoco
Integer to indentify OCO2 or OCO3
Land.CanopyLayers.SIF_771
— FunctionSIF_771(can_rad::CanopyRads{FT},
wls::WaveLengths{FT};
oco::Int = 2
) where {FT<:AbstractFloat}
Return the SIF @ 769.94 (OCO2) or 770.005 (OCO3) nm, given
can_rad
CanopyRads
type structwls
WaveLengths
type structoco
Integer to indentify OCO2 or OCO3
Land.CanopyLayers.SIF_WL
— FunctionSIF_WL(wls::WaveLengths{FT},
can_rad::CanopyRads{FT}
wls::WaveLengths{FT},
twl::FT
) where {FT<:AbstractFloat}
Return the SIF, given
can_rad
CanopyRads
type structwls
WaveLengths
type structtwl
Target SIF wave length in nm
Land.CanopyLayers.SWIR
— FunctionSWIR(can_rad::CanopyRads{FT},
wls::WaveLengths{FT}
) where {FT<:AbstractFloat}
Return the SWIR @ 2130 nm, given
can_rad
CanopyRads
type structwls
WaveLengths
type struct
Utils
Land.CanopyLayers.calctav
— Functioncalctav(α::FT, nr::FT) where {FT<:AbstractFloat}
Computes transmission of isotropic radiation across an interface between two dielectrics (Stern F., 1964; Allen W.A., 1973)). From calctav.m in PROSPECT-D
α
angle of incidencenr
Index of refraction
Land.CanopyLayers.dcum
— Functiondcum(a::FT, b::FT, t::FT) where {FT<:AbstractFloat}
TODO Add function description
Land.CanopyLayers.dladgen
— Functiondladgen(a::FT, b::FT, litab_bnd::Array{FT,2}) where {FT<:AbstractFloat}
TODO Calculate the freqency of WHAT?
Land.CanopyLayers.e2phot
— Functione2phot(λ::Array{FT}, E::Array{FT}) where {FT<:AbstractFloat}
Calculates the number of moles of photons, given
λ
An array of wave length in[nm]
, converted to[m]
by _FACE
Joules of energy
Land.CanopyLayers.e2phot!
— Functione2phot!(
λ::Array{FT,1},
E::Array{FT,1},
cache::Array{FT,1}
) where {FT<:AbstractFloat}
Calculates the number of moles of photons, given
λ
An array of wave length in[nm]
, converted to[m]
by _FACE
Joules of energycache
Cache to avoid memory allocations
Land.CanopyLayers.psofunction
— Functionpsofunction(K::FT,
k::FT,
Ω::FT,
LAI::FT,
q::FT,
dso::FT,
xl::FT
) where {FT<:AbstractFloat}
TODO explain the variables
Return the probability of observing a sunlit leaf at depth xl
(pso
, see eq 31 in vdT 2009), given
xl
Leaf depth in the canopy
Land.CanopyLayers.volscatt!
— Functionvolscatt!(cache::Array{FT,1},
sza::FT,
vza::FT,
raa::FT,
ttl::FT
) where {FT<:AbstractFloat}
Calculate interception parameters (chi_s
and chi_s
) and leaf reflectance multiplier (frho
) and transmittance multiplier (ftau
), given
cache
Array cache for resultssza
Solar zenith anglevza
Viewing zenith angleraa
Relative azimuth anglettl
Leaf inclination angle