Package 'ccrtm'

Title: Coupled Chain Radiative Transfer Models
Description: A set of radiative transfer models to quantitatively describe the absorption, reflectance and transmission of solar energy in vegetation, and model remotely sensed spectral signatures of vegetation at distinct spatial scales (leaf,canopy and stand). The main principle behind ccrtm is that many radiative transfer models can form a coupled chain, basically models that feed into each other in a linked chain (from leaf, to canopy, to stand, to atmosphere). It allows the simulation of spectral datasets in the solar spectrum (400-2500nm) using leaf models as PROSPECT5, 5b, and D which can be coupled with canopy models as 'FLIM', 'SAIL','SAIL2' and 'INFORM' for top of canopy reflectance, and with atmospheric models such as 'SKYL'and 'SMAC' for calculation of top of the atmosphere reflectance. All models can run in forward mode, and a selection are invertable (backward simulations) if provided with spectral data. Jacquemoud et al 2008 provides a comprehensive overview of these and other models <doi:10.1016/j.rse.2008.01.026>.
Authors: Marco D. Visser [aut, cre]
Maintainer: Marco D. Visser <[email protected]>
License: GPL (>= 2)
Version: 0.4.1
Built: 2025-02-15 12:06:43 UTC
Source: https://github.com/marcodvisser/ccrtm

Help Index


Backward implementation (inversion) of coupled Radiative Transfer Models.

Description

Backward implementation (inversion) of coupled Radiative Transfer Models.

Usage

bRTM(fm = rho ~ prospect5, data = NULL)

Arguments

fm

A formula specifying which rtm to run (see details).

data

(measured) reflectance spectra. Expected are reflectance values between 0 and 1 for wavelength between 400 and 2400 at 1 nm steps. The range 400 to 2400 is based on the largest common range in most leaf spectral datasets - and hence is a range that can be generated by most spectrometers.

Details

Formula: In general the form of the formula specifies the both the model and the data supplied (transmittance or reflectance), however, currently only reflectance data is expected (transmission not included yet).

Models: At current the following radiative transfer models are implemented for backward / inversion mode

Example of a formula Model
rho ~ prospect5 prospect5
rho ~ prospectd prospectd

Inversion is rapid, and based on emulation of prospect models by a multivariate neural net (MANN) and a partial least squares regression (PLSR) model. The two methods are selected as the performance of NN or PLSR differ for each inverted parameter - with one method outperforming the other depending on the parameter. The predictions are then combined using a linear Bayesian mixing model that weights the NN and PLSR prediction for each parameter - and includes an estimate of model inversion uncertainty.

Model inversion uncertainty estimates the 95% credible intervals under which the parameter will fall compared to a perfect inversion. Model inversion uncertainty arises due to parameter identifiability issues, and does not reflect the uncertainty in the data. Uncertainty in the data should be estimated with replicate measurements and standard statistical methods (not implemented).

Questions and requests can be made on the ccrtm github page.

Value

a list of inverted parameters and their 95% CI

Examples

## get reflectance for a single leaf on simulated spectra

## make a parameter list
parameters<-list(prospectd=c(N=3,Cab=40,Car=15,Cw=0.01,Cm=0.025,Canth=26,Cbrown=4))

## simulate spectra at the inversion requirements 
ref <- fRTM(rho~prospectd,pars=parameters,wl=400:2400)

## reorder with replicates measurements over rows, and make into matrix
refdata<-t(as.matrix(ref))

fit<-bRTM(rho~prospectd,data=refdata)
summary(fit)

## compare fit with simulation on log-scale so all parameter are visible
plot(parameters$prospectd,fit$mu,xlab="expected",ylab="inverted",pch=16,log="xy")

## add uncertainty
segments(parameters$prospectd,fit$lower.ci,
parameters$prospectd,fit$upper.ci,lwd=2)

## 1 to 1 line
abline(0,1)

## Inversion for multiple leaf spectra

## using lower-level vectorized prospect function
set.seed(1234)
## we simulate all spectra at once
nsim<-300 ## number of leaves

## random leaf parameters
parmat<-cbind(N=runif(nsim,1,6),
Cab=runif(nsim,5,80),
Car=runif(nsim,1,40),
Cw=runif(nsim,0.001,.02),
Cm=runif(nsim,0.002,0.03)+0.01,
Canth=runif(nsim,0,6),
Cbrown=runif(nsim,0,4)
)

## simulate with the lower level prospect for rapid simulation
## of many leaves
ref<-ccrtm:::.prospectdv(parmat)[[1]][,1:2001] ## subset to 400:2400 wl

## invert the simulations
fit<-bRTM(rho~prospectd,data=ref)
summary(fit)

## check inversion performace for LMA
plot(parmat[,"Cm"],fit$mu[,"Cm"],xlab="expected",ylab="inverted",pch=16)

## add uncertainty
segments(parmat[,"Cm"],fit$lower.ci[,"Cm"],
         parmat[,"Cm"],fit$upper.ci[,"Cm"],lwd=2)
abline(0,1)

## replace the simulated ref with measured reflectance over wavelengths 400:2400
## to invert for spectrometer data

Leaf inclination distribution function Ellipsoidal distribution function

Description

Leaf inclination distribution function Ellipsoidal distribution function

Usage

cambell(ala, tx1, tx2)

Arguments

ala

average leaf angle parameter

tx1

angle in degree

tx2

angle in degree

Value

angle fraction value


ccrtm: Coupled Chain Radiative Transfer Models.

Description

A collection of radiative transfer models that can form a coupled chain to model radiative transfer across multiple spatial scales from leaf to canopy to stand.

Details

Currently implemented models that can be coupled:

  • 1 = PROSPECT 5, 5B and D

  • 2 = FOURSAIL, and FOURSAIL2

  • 3 = FLIM

Currently being tested / or to be implemented models

  • 1 = LIBERTY, PROCOSINE

  • 2 = INFORM*

*available as lower-level library (see ccrtm github page).

Generating predictions

To generate model prediction the typical approach is to use the fRTM function and apply a formula that specifies the both the expected output (left hand) and the different models you would like to couple to generate the output (right hand).

At current the following radiative transfer models, and corresponding formula, are given in the next table

Example of a formula Model
rho ~ prospect5 prospect5
rho ~ prospectd prospectd
rho ~ prospect5 + foursail PROSAIL
rho ~ prospect5 + foursail PROSAIL
rho ~ prospectd + foursail2 PROSAIL
rho ~ prospectd + prospect5 + foursail2 PROSAIL2
rho ~ prospectd + foursail2b PROSAIL2b
rho ~ prospectd + foursail + flim + sky INFORM*
rho ~ prospectd + foursail + flim INFORM*

*INFORM is currently restricted to a lower-level function only. See the cctrm github readme page on how to use it.

In the above examples, prospectd can be replaced with prospect5

  • or vice versa - if so desired.

Also, note that the formula "rho ~ prospectd + foursail2" is the same as "rho ~ prospectd + prospectd + foursail2" and both expect a names list of 3 parameter vectors.

See the help files for details on each right hand component. For instance, ?foursail provides more elaboration on the 4SAIL model and gives an example for lower-level implementations of each component model. See also ?flim, ?foursail2, ?foursail2b, ?prospect5, and ?prospectd.

Tranmission can also be returned if specified in the left-hand component of the formula:

Formula Model
rho + tau ~ prospect5 prospect5
rho + tau ~ prospectd + foursail 4SAIL

The examples above indicate that the users wishes to predict transmission next to reflectance. More specifically, The first returns leaf reflectance and transmission while the second returns 4 components of canopy reflectance and canopy transmission in the solar and viewing direction.

PROSAIL and diffuse and direct light

In contrast to some FORTRAN and MATLAB implementation, the sky light model is not implemented by default in ccrtm. This is because it is not a standard component of 4SAIL. In addition, this would affect limit the application of other more realistic atmospheric models. You can apply it by using ?skyl on model output obtained from fRTM (see example in ?skyl). The sky light model is implemented for the INFORM model as per March 2022.

Parameters

Parameters are given as input to fRTM as a named list. See ?getDefaults for examples on how to structure parameters. Individual models can consulted on each parameter (e.g. see ?foursail2b or ?prospect5).

Leaf inclination models.

Canopy models such as foursail use leaf inclination models. In cctrm four inclination models are implemented. see ?lidf for more details.

Author(s)

Marco D. Visser


Leaf inclination distribution function cummulative lagden function from Wout Verhoef's dissertation Extended here for any angle

Description

Leaf inclination distribution function cummulative lagden function from Wout Verhoef's dissertation Extended here for any angle

Usage

cdcum(a, b, theta)

Arguments

a

parameter

b

parameter

theta

angle in degrees

Value

angle fraction value


Function to check and return parameters

Description

Function to check and return parameters

Usage

checkPars(pars, fm, ordN)

refractive index and specific absorption coefficients for PROSPECT 5

Description

see http://teledetection.ipgp.jussieu.fr/prosail/ for more details on the data.

Usage

data(prospect5)

details

data_prospect5 (february, 25th 2008) The dataset contains the following labels (columns):

  • 1 = wavelength (nm)

  • 2 = refractive index of leaf material ( or the ratio of the velocity of light in a vacuum to its velocity in "leaf medium").

  • 3 = specific absorption coefficient of chlorophyll (a+b) (cm2.microg-1)

  • 4 = specific absorption coefficient of carotenoids (cm2.microg-1)

  • 5 = specific absorption coefficient of brown pigments (arbitrary units)

  • 6 = specific absorption coefficient of water (cm-1)

  • 7 = specific absorption coefficient of dry matter (g.cm-1)

  • 8 = direct light

  • 9 = diffuse light

  • 10 = dry soil

  • 11 = wet soil

references

Feret et al. (2008), PROSPECT-4 and 5: Advances in the Leaf Optical Properties Model Separating Photosynthetic Pigments, Remote Sensing of Environment


refractive index and specific absorption coefficients for PROSPECT D

Description

see http://teledetection.ipgp.jussieu.fr/prosail/ for more details on the data.

Usage

data(prospectd)

details

data_prospect5 (february, 25th 2008) The dataset contains the following labels (columns):

  • 1 = wavelength (nm)

  • 2 = refractive index of leaf material ( or the ratio of the velocity of light in a vacuum to its velocity in "leaf medium").

  • 3 = specific absorption coefficient of chlorophyll (a+b) (cm2.microg-1)

  • 4 = specific absorption coefficient of carotenoids (cm2.microg-1)

  • 5 = specific absorption coefficient of brown pigments (arbitrary units)

  • 6 = specific absorption coefficient of water (cm-1)

  • 7 = specific absorption coefficient of dry matter (g.cm-1)

  • 8 = direct light

  • 9 = diffuse light

  • 10 = dry soil

  • 11 = wet soil

references

Feret et al. (2008), PROSPECT-4 and 5: Advances in the Leaf Optical Properties Model Separating Photosynthetic Pigments, Remote Sensing of Environment


d = stand density (d) cd = crown diameter (cd) h = mean crown height (h) lai = leaf area index (lai) alpha = light extinction coefficient (alpha) tts = Solar zenith angle (tts) tto = Observer zenith angle (tto) psi = Sun-sensor azimuth angle (psi)

Description

d = stand density (d) cd = crown diameter (cd) h = mean crown height (h) lai = leaf area index (lai) alpha = light extinction coefficient (alpha) tts = Solar zenith angle (tts) tto = Observer zenith angle (tto) psi = Sun-sensor azimuth angle (psi)

Usage

## S3 method for class 'flim'
defaults(x, simple = TRUE)

eigen decomposition for PROSPECT5

Description

data reduction used on simulated data from PROSPECT5 (for NN and PLSR)


eigen decomposition for PROSPECTD

Description

data reduction used on simulated data from PROSPECTD (for NN and PLSR)


Forest Light Interaction Model (FLIM)

Description

The FLIM model was first described by Rosema et al (1992). In FLIM forests are assumed a discontinous mix of tree crowns and gaps. Reflectance is modelled in terms of the probabilty to observe either a gap (background) or a tree crown. Both gaps and crowns may be shaded.

Usage

flim(Rc, Rg, To = NULL, Ts = NULL, params)

Arguments

Rc

Canopy reflectance at infinite depth

Rg

soil/background reflectance

To

transmission in viewing direction

Ts

transmission in sun direction

params

a named vector of parameters:

  • 1 = stand density (d) (1)

  • 2 = crown diameter (cd) (1)

  • 3 = mean crown height (h)

  • 4 = leaf area index (lai) (2)

  • 5 = light extinction coefficient (alpha) (2)

  • 6 = Solar zenith angle (tts)

  • 7 = Observer zenith angle (tto)

  • 8 = Sun-sensor azimuth angle (psi)

area

area of stand (m2)

Details

(1) Parameters are confounded (d & cd), confounded parameters pairs cannot be inversely estimated, one of the pairs should be set to 1 - or given explicitly. (2) required if T0, Ts are not given.

Value

a list with reflectance, and the fractions of shaded and sunexplosed crowns, shaded and sun exposed open space.

References

Rosema, A., Verhoef, W., Noorbergen, H., Borgesius, J.J. (1992). A new forest light interaction model in support of forest monitoring. Remote Sens. Environ. 42, 23-41.

Examples

parvec<- c(alpha = 0.5,lai=5,cd=15,h=30,d=10,tto=10,tts=20,psi=30)
flim(0.1,0.1,params=parvec)

Optimized R implementation of foursail (4SAIL)

Description

The foursail (or 4SAIL) radiative transfer model is commonly used to simulate bidirectional reflectance distribution functions within vegetation canopies. Foursail (4SAIL) refers to "Scattering by Arbitrary Inclined Leaves" in a 4-stream model. The four-streams represents the scattering and absorption of upward, downward and two directional radiative fluxes with four linear differential equations in a 1-D canopy. The model was initially developed by Verhoef (1984), who extended work by Suits (1971) 4-steam model.

Usage

foursail(rho, tau, bgr, param)

Arguments

rho

input leaf reflectance from 400-2500nm (can be measured or modeled)

tau

input leaf transmittance from 400-2500nm (can be measured or modeled)

bgr

background reflectance. Usual input is soil reflectance spectra from 400-2500nm (can be measured or modeled)

param

A named vector of SAIL parameter values (note: program ignores case):

  • 1 = Leaf angle distribution function parameter a (LIDFa)

  • 2 = Leaf angle distribution function parameter b (LIDFb)

  • 3 = Leaf angle distribution function type (see ?lidf)

  • 4 = Leaf area index (LAI)

  • 5 = Hot spot effect parameter (hspot)

  • 6 = Solar zenith angle (tts)

  • 7 = Observer zenith angle (tto)

  • 8 = Sun-sensor azimuth angle (psi)

Value

spectra matrixwith 4 reflectance factors and canopy transmission for wavelengths 400 to 2500nm:

  • 1 = bi-hemispherical reflectance (rddt). White-sky albedo: the reflectance of the canopy under diffuse illumination. The BRDF integrated over all viewing and illumination directions.

  • 2 = hemispherical directional reflectance (rsdt). Black-sky albedo: reflectance of a surface under direct light without a diffuse component. It is the integral of the BRDF over all viewing directions.

  • 3 = directional hemispherical reflectance (rdot). Diffuse reflectance in the vieweing direction.

  • 4 = bi-directional reflectance (rsot). The ratio of reflected radiance in the viewing direction to the incoming radiant flux in the solar direction.

  • 5 = Canopy transmission of diffuse light through the canopy (taud).

  • 6 = transmission of direct light through the canopy in the solar direction (taus).

  • 7 = transmission of direct light through the canopy in the sensor/viewing direction (tauo).

References

Suits, G.H., 1971. The calculation of the directional reflectance of a vegetative canopy. Remote Sens. Environ. 2, 117-125.

Verhoef, W. (1984). Light scattering by leaf layers with application to canopy reflectance modeling: The SAIL model. Remote Sens. Environ. 16, 125-141.

Examples

## lower-level implementation example
## see ?fRTM for the typical mode of simulation
## e.g. fRTM(rho~prospectd+foursail)

## 1) get parameters
params<-getDefaults("foursail")

##  getDefaults(rho~prospectd+foursail) will also work
pars<-params$foursail

## ensure the vector is named
names(pars) <- names(params$foursail)

## 2) get leaf reflectance and transmission 
rt<-fRTM(rho+tau~prospectd)

## 3) get soil reflectance to model background reflectance
data(soil)

## a linear mixture soil model
bgRef<- pars["psoil"]*soil[,"drySoil"] + (1-pars["psoil"])*soil[,"wetSoil"]


## 4) run 4SAIL
result<-foursail(rt[,"rho"],rt[,"tau"],bgRef,pars)
head(result)

R implementation of the foursail2 model with 2 canopy layers.

Description

The foursail2 model is a two layer implementation of the foursail model described in Verhoef and Bach (2007). Layers are assumed identical in particle inclination and hotspot, but may differ in the relative density and types of particles (see foursail2b for a layer specific inclination angle). In comparison to foursail, the background (soil), can now be non-Lambertain, having it own 4-stream BDRF (not implemented here but may be input by the user). There are two types of particles, generalized to primary and secondary (originally termed "green" and "brown" particles). The realtive abundance of the secondary particle in the top canopy is regulated by the dissociation paramerter.The model 4SAIL2 combines with prospect, libery or procosine for the reflectance and transmittance of the particles, and with the the foursail or Hapke elements for the background reflectance. If run alone, these require direct inputs which could be measured leaf reflectance.

Usage

foursail2(
  rhoA,
  tauA,
  rhoB = NULL,
  tauB = NULL,
  bgr,
  rsobgr = NULL,
  rdobgr = NULL,
  rsdbgr = NULL,
  rddbgr = NULL,
  param
)

Arguments

rhoA

primary particle reflectance from 400-2500nm (can be measured or modeled)

tauA

primary particle transmittance from 400-2500nm (can be measured or modeled)

rhoB

secondary particle reflectance from 400-2500nm (can be measured or modeled)

tauB

secondary particle reflectance from 400-2500nm (can be measured or modeled)

bgr

background reflectance. Usual input is soil reflectance spectra from 400-2500nm (can be measured or modeled)

rsobgr

: background bidirectional reflectance (rso)

rdobgr

: background directional hemispherical reflectance (rdo)

rsdbgr

: background hemispherical directional reflectance (rsd)

rddbgr

: background bi-hemispherical diffuse reflectance (rdd)

param

A named vector of 4SAIL2 parameter values (note: program ignores case):

  • 1 = Leaf angle distribution function parameter a (LIDFa)

  • 2 = Leaf angle distribution function parameter b (LIDFb)

  • 3 = Leaf angle distribution function type (TypeLidf, see ?lidf)

  • 4 = Total Leaf Area Index (LAI), including primary and secondary particles (brown and green leafs).

  • 5 = fraction secondary particles ("brown leaf fraction", fb)

  • 6 = Canopy dissociation factor for secondary particles ("diss")

  • 7 = Hot spot effect parameter (hspot). Often defined as the ratio of mean leaf width and canopy height.

  • 7 = vertical crown coverage fraction (Cv), models clumping in combination with parameter zeta.

  • 7 = tree shape factor (zeta), defined as the ratio of crown diameter and height.

  • 6 = Solar zenith angle (tts)

  • 7 = Observer zenith angle (tto)

  • 8 = Sun-sensor azimuth angle (psi)

Value

spectra matrixwith 4 reflectance factors and canopy transmission for wavelengths 400 to 2500nm:

  • 1 = bi-hemispherical reflectance (rddt). White-sky albedo: the reflectance of the canopy under diffuse illumination. The BRDF integrated over all viewing and illumination directions. Diffuse reflectance for diffuse incidence.

  • 2 = hemispherical directional reflectance (rsdt). Black-sky albedpo: reflectance of a surface under direct light without a diffuse component. It is the integral of the BRDF over all viewing directions. Diffuse reflectance for direct solar incidence.

  • 3 = directional hemispherical reflectance (rdot). Diffuse reflectance in the vieweing direction.

  • 4 = bi-directional reflectance (rsot). The ratio of reflected radiance in the viewing direction to the incoming radiant flux in the solar direction.

References

Verhoef, W., Bach, H. (2007). Coupled soil-leaf-canopy and atmosphere radiative transfer modeling to simulate hyperspectral multi-angular surface reflectance and TOA radiance data. Remote Sens. Environ. 109, 166-182.

Examples

## see ?foursail for lower-level implementations
fRTM(rho~prospect5+foursail2)

R implementation of the foursail2 model with 2 canopy layers.

Description

The foursail2b model is a two layer implementation of the foursail model described in Zhang et al (2005). Layers are assumed identical in hotspot, but may differ in the relative density, inclination and types of particles. In comparison to foursail, the background (soil), can now be non-Lambertain, having it own 4-stream BDRF (not implemented here but may be input by the user). There are two types of particles, generalized to primary and secondary (originally termed "green" and "brown" particles). The realtive abundance of the secondary particle in the top canopy is regulated by the dissociation paramerter. The model 4SAIL2 combines with prospect, libery or procosine for the reflectance and transmittance of the particles, and with the the foursail or Hapke elements for the background reflectance. If run alone, these require direct inputs which could be measured leaf reflectance.

Usage

foursail2b(
  rhoA,
  tauA,
  rhoB = NULL,
  tauB = NULL,
  bgr,
  rsobgr = NULL,
  rdobgr = NULL,
  rsdbgr = NULL,
  rddbgr = NULL,
  param
)

Arguments

rhoA

primary particle reflectance from 400-2500nm (can be measured or modeled)

tauA

primary particle transmittance from 400-2500nm (can be measured or modeled)

rhoB

secondary particle reflectance from 400-2500nm (can be measured or modeled)

tauB

secondary particle reflectance from 400-2500nm (can be measured or modeled)

bgr

background reflectance. Usual input is soil reflectance spectra from 400-2500nm (can be measured or modeled)

rsobgr

: background bidirectional reflectance (rso)

rdobgr

: background directional hemispherical reflectance (rdo)

rsdbgr

: background hemispherical directional reflectance (rsd)

rddbgr

: background bi-hemispherical diffuse reflectance (rdd)

param

A named vector of 4SAIL2 parameter values (note: program ignores case):

  • 1 = Mean leaf angle for first (top) layer (LIDFa)

  • 2 = Mean leaf angle for second (bottom) layer (LIDFb)

  • 3 = Leaf angle distribution function type (ignored, only value 2 allow)

  • 4 = Total Leaf Area Index (LAI), including primary and secondary particles (brown and green leafs).

  • 5 = fraction secondary particles ("brown leaf fraction", fb)

  • 6 = Canopy dissociation factor for secondary particles ("diss")

  • 7 = Hot spot effect parameter (hspot). Often defined as the ratio of mean leaf width and canopy height.

  • 7 = vertical crown coverage fraction (Cv), models clumping in combination with parameter zeta.

  • 7 = tree shape factor (zeta), defined as the ratio of crown diameter and height.

  • 6 = Solar zenith angle (tts)

  • 7 = Observer zenith angle (tto)

  • 8 = Sun-sensor azimuth angle (psi)

Details

Leaf inclination angles: leaf angles in 4SAIL2b are set for each layer and only the Cambell leaf angle distribution model is allowed. This means that each layer has a single parameter that defines leaf angles.

Value

spectra matrixwith 4 reflectance factors and canopy transmission for wavelengths 400 to 2500nm:

  • 1 = bi-hemispherical reflectance (rddt). White-sky albedo: the reflectance of the canopy under diffuse illumination. The BRDF integrated over all viewing and illumination directions. Diffuse reflectance for diffuse incidence.

  • 2 = hemispherical directional reflectance (rsdt). Black-sky albedpo: reflectance of a surface under direct light without a diffuse component. It is the integral of the BRDF over all viewing directions. Diffuse reflectance for direct solar incidence.

  • 3 = directional hemispherical reflectance (rdot). Diffuse reflectance in the vieweing direction.

  • 4 = bi-directional reflectance (rsot). The ratio of reflected radiance in the viewing direction to the incoming radiant flux in the solar direction.

References

Zhang, Q., Xiao, X., Braswell, B., Linder, E., Baret, F., Moore, B. (2005). Estimating light absorption by chlorophyll, leaf and canopy in a deciduous broadleaf forest using MODIS data and a radiative transfer model. Remote Sens. Environ. 99, 357-371.

Examples

## see ?foursail for lower-level implementations
fRTM(rho~prospectd+foursail2b)

Forward implementation of coupled Radiative Transfer Models.

Description

Forward implementation of coupled Radiative Transfer Models.

Usage

fRTM(fm = rho + tau ~ prospect5 + foursail, pars = NULL, wl = 400:2500)

Arguments

fm

A formula specifying which rtm to run (see details).

pars

a named list of named parameter vectors for all models. The parameter list for a model call as rho ~ prospect + foursail, therefore, contains two vectors: the first with parameters for prospect and the second with parameters for foursail. See ?getDefaults for an example of a parameter list. If left empty default parameters are generated.

wl

wavelengths (in nm) add only if certain wavelengths are required as output. Input is expected to integers between 400 and 2500, or will be forced to be an integer. Integers outside the 400:2500 range will not be returned.

Details

Formula: In general the form of the formula specifies the both the expected output (left hand) and the different models you would like to couple to generate the output (right hand).

Models: At current the following radiative transfer models are implemented

Example of a formula Model
rho ~ prospect5 prospect5
rho ~ prospectd prospectd
rho ~ prospect5 + foursail PROSAIL
rho ~ prospect5 + foursail PROSAIL
rho ~ prospectd + foursail2 PROSAIL
rho ~ prospectd + prospect5 + foursail2 PROSAIL2*
rho ~ prospectd + foursail2b PROSAIL2b*
rho ~ prospectd + foursail + flim + sky INFORM**
rho ~ prospectd + foursail + flim INFORM**
  • Note that the formula "rho ~ prospectd + foursail2" is the same as "rho ~ prospectd + prospectd + foursail2" and both expect a names list of 3 parameter vectors (leaf type 1, leaf type 2, and the canopy parameters).

** INFORM is currently restricted to a lower-level function only. See the cctrm github readme page on how to use it.

In the above examples with additive components, prospectd can be replaced with prospect5 - or vice versa - if so desired. See the help files for details on each right hand component. For instance, ?foursail provides more elaboration on the 4SAIL model and gives an example for lower-level implementations of each component model.

Tranmission can also be returned if specified in the left-hand component of the formula:

Formula Model
rho + tau ~ prospect5 prospect5
rho + tau ~ prospectd + foursail 4SAIL

The examples above indicate that the users wishes to predict transmission next to reflectance. More specifically, The first returns leaf reflectance and transmission while the second returns 4 components of canopy reflectance and canopy transmission in the solar and viewing direction.

More details are given in ?cctrm.

Questions and requests can be made on the ccrtm github page.

Value

spectra matrix with reflectance (and transmission, depending on the formula inputs). See seperate model helpfiles for details.

Examples

## setup graphics for plots
oldpar<-par()
par(mfrow=c(3,2))

## get reflectance for a leaf
ref <- fRTM(rho~prospect5)
plot(ref,main="Prospect 5")

## get reflectance and transmission for a leaf
reftrans <- fRTM(rho+tau~prospect5)
plot(reftrans,main="Prospect 5")

## get reflectance for a single layered canopy
ref <- fRTM(rho~prospect5+foursail)
plot(ref,main="Prospect 5 + 4SAIL")

## get reflectance for a 2 layered canopy with two leaf types
ref <- fRTM(rho~prospectd+prospect5+foursail2)
plot(ref,main="Prospect D + Prospect 5  + 4SAIL2")

## edit the parameters: sparse vegetation LAI
parlist <- getDefaults(rho~prospectd+prospect5+foursail2)
parlist$foursail2["LAI"] <- 0.05

## update reflectance
ref <- fRTM(rho~prospect5+prospectd+foursail2,parlist)
plot(ref,main="LAI=0.05")

## change leaf area index to dense vegetation
parlist$foursail2["LAI"]<-8.5

## update reflectance
ref <- fRTM(rho~prospect5+prospectd+foursail2,parlist)
plot(ref,main="LAI=8.5")

par(oldpar)

S3- methods for Generate defaults settings and parameters for all supported models. See ?ccrtm for details.

Description

S3- methods for Generate defaults settings and parameters for all supported models. See ?ccrtm for details.

Usage

getDefaults(model = NULL, ...)

Arguments

model

a ccrtm formula (rho ~ prospectd) or character vector of modelnames (e.g. "prospect5")

...

not used.

Value

a data.frame with default model parameters


invert a requested RTM (internal function)

Description

List of aliases: prospect5, prospectd

Usage

invertRTM(pars)

Arguments

pars

the required parameters (vector or list), and newdata

modReq

model request object built in bRTM

Value

prediction from the requested model


Kullback-Lieber divergence function D(spec1 || spec2) = sum(spec1 * log(spec1 / spec2))

Description

Kullback-Lieber divergence function D(spec1 || spec2) = sum(spec1 * log(spec1 / spec2))

Usage

KLd(spec1, spec2)

Arguments

spec1

spectral signal 1

spec2

spectral signal 2 at identical wavelengths

Value

the KL divergence between the vector inputs


Leaf inclination distribution function models s3 method for calling leaf models.

Description

Leaf inclination distribution function models s3 method for calling leaf models.

Usage

lidf(pars)

Arguments

pars

a parameter vector c(angles, LIDFa, LIDFb) with a class lidf.modelnumber. Models include:

  • 1 = Dlagden distribution (1, lidf.1)

  • 2 = Ellipsoid (Campebll) distribution (2, lidf.2)

  • 3 = Beta distribution (3, lidf.3)

  • 4 = One parameter beta distribution (4, lidf.4)

Models 1 and 2 are the standard models from the SAIL model. Two parameter models use parameters LIDFa and LIDFb, while single parameter models use only LIDFa (ignoring any supplied LIDFb).

More information on the Dlagden and Ellipsoid parameter is given in Verhoef, W. (1998),theory of radiative transfer models applied in optical remote sensing of vegetation canopies (PhD thesis).

The beta distribution is the typical beta distribution as often implemented (as in dbeta(x,LIDFa, LIDFb)). Where x is a value between 0 and 90, that gives the angular density over 0 and 90 degrees (rescaled to 0 and 1).

The one parameter beta distribution is given by LIDFa*x^(LIDFa-1). Where x is a value between 0 and 90, that given the angular density over 0 and 90 degrees (rescaled to 0 and 1).

Value

a vector of proportions for each leaf angle calculated from each leaf inclination model.


Bayesian fitted weight matrix for PROSPECT5

Description

Weight coefficients for neural network and plsr predictions .


Bayesian fitted weight matrix for PROSPECTD

Description

Weight coefficients for neural network and plsr predictions .


fitted weight matrix for PROSPECT5

Description

Weight matrices for a fit neural network on simulated data from PROSPECT5.


fitted weight matrix for PROSPECTD

Description

Weight matrices for a fit neural network on simulated data from PROSPECTD.


Plot RTM return spectra vs. wavelength

Description

Plot RTM return spectra vs. wavelength

Usage

## S3 method for class 'rtm.spectra'
plot(x, ...)

Arguments

x

predictions from an RTM

...

additional plot arguments

Value

plots to the device a ccrtm standard spectra plot based on the function call returned from fRTM.


fitted PLSR for PROSPECT5

Description

A partial least squares model fit on simulated data from PROSPECT5.


fitted PLSR for PROSPECTD

Description

A partial least squares model fit on simulated data from PROSPECTD.


RTM inversion

Description

RTM inversion

Usage

## S3 method for class 'rtm.inversion'
print(x, ...)

Arguments

x

predictions from an RTM

...

additional plot arguments

Value

prints the inverted parameters


RTM generic print function

Description

RTM generic print function

Usage

## S3 method for class 'rtm.spectra'
print(x, ...)

Arguments

x

predictions from an RTM

...

additional plot arguments

Value

prints the standard information from a simulated ccrtm spectra plot


PROSPECT model version 5 and 5B

Description

The PROSPECT5(b) leaf reflectance model. The model was implemented based on Jacquemoud and Ustin (2019), and is further described in detail in Feret et al (2008). PROSPECT models use the plate models developed in Allen (1969) and Stokes (1862). Set Cbrown to 0 for prospect version 5.

Usage

prospect5(param)

Arguments

param

A named vector of PROSPECT parameters (note: program ignores case):

  • 1 = leaf structure parameter (N)

  • 2 = chlorophyll a+b content in ug/cm2 (Cab)

  • 3 = carotenoids content in ug/cm2 (Car)

  • 4 = brown pigments content in arbitrary units (Cbrown)

  • 5 = equivalent water thickness in g/cm2 (Cw)

  • 6 = leaf dry matter content in g/cm2 - lma - (Cm)

Value

spectra matrix with leaf reflectance and transmission for wavelengths 400 to 2500nm:

  • 1 = leaf reflectance (rho)

  • 2 = leaf transmission (tau)

References

Jacquemoud, S., and Ustin, S. (2019). Leaf optical properties. Cambridge University Press.

Feret, J.B., Francois, C., Asner, G.P., Gitelson, A.A., Martin, R.E., Bidel, L.P.R., Ustin, S.L., le Maire, G., Jacquemoud, S. (2008), PROSPECT-4 and 5: Advances in the leaf optical properties model separating photosynthetic pigments. Remote Sens. Environ. 112, 3030-3043.

Allen W.A., Gausman H.W., Richardson A.J., Thomas J.R. (1969), Interaction of isotropic ligth with a compact plant leaf, Journal of the Optical Society of American, 59:1376-1379.

Stokes G.G. (1862), On the intensity of the light reflected from or transmitted through a pile of plates, Proceedings of the Royal Society of London, 11:545-556.

Examples

## see ?fRTM for the typical mode of simulation
## e.g. fRTM(rho~prospect5)

## 1) get parameters
defaultpars<-getDefaults(rho~prospect5)
## getDefaults("prospect5") will also work

## 2) get leaf reflectance and transmission
rt<-fRTM(rho+tau~prospect5,defaultpars)

## lower-level implementation example
## Alternatively implement directly
mypars<-c("N"=1,"Cab"=35,"Car"=20,"Cbrown"=3,"Cw"=0.01,"Cm"=0.01)
prospect5(mypars)

PROSPECT model version D

Description

The PROSPECTD leaf reflectance model. The model was implemented based on Jacquemoud and Ustin (2019), and is further described in detail in Feret et al (2017). PROSPECT models use the plate models developed in Allen (1969) and Stokes (1862).

Usage

prospectd(param)

Arguments

param

A named vector of PROSPECT parameters (note: program ignores case):

  • 1 = leaf structure parameter (N)

  • 2 = chlorophyll a+b content in ug/cm2 (Cab)

  • 3 = carotenoids content in ug/cm2 (Car)

  • 4 = Leaf anthocyanin content (ug/cm2) (Canth)

  • 5 = brown pigments content in arbitrary units (Cbrown)

  • 6 = equivalent water thickness in g/cm2 (Cw)

  • 7 = leaf dry matter content in g/cm2 - lma - (Cm)

Value

spectra matrix with leaf reflectance and transmission for wavelengths 400 to 2500nm:

  • 1 = leaf reflectance (rho)

  • 2 = leaf transmission (tau)

References

Jacquemoud, S., and Ustin, S. (2019). Leaf optical properties. Cambridge University Press.

Feret, J.B., Gitelson, A.A., Noble, S.D., Jacquemoud, S. (2017). PROSPECT-D: Towards modeling leaf optical properties through a complete lifecycle. Remote Sens. Environ. 193, 204-215.

Allen W.A., Gausman H.W., Richardson A.J., Thomas J.R. (1969), Interaction of isotropic ligth with a compact plant leaf, Journal of the Optical Society of American, 59:1376-1379.

Stokes G.G. (1862), On the intensity of the light reflected from or transmitted through a pile of plates, Proceedings of the Royal Society of London, 11:545-556.

Examples

## see ?fRTM for the typical mode of simulation
## e.g. fRTM(rho~prospectd) 

## 1) get parameters
defaultpars<-getDefaults(rho~prospectd) 
## getDefaults("prospectd") will also work

## 2) get leaf reflectance and transmission 
rt<-fRTM(rho+tau~prospectd,defaultpars)

## lower-level implementation example
## Alternatively implement directly (case ignored for parameters)
mypars<-c("N"=1,"Cab"=35,"Car"=20,"Canth"=15,"Cbrown"=3,"Cw"=0.01,"Cm"=0.01)
prospectd(mypars)

R implementation of foursail (pure R)

Description

The pure R version of foursail is included in the package as an easy way to review the code, and to check more optimized versions against. Model originally developed by Wout Verhoef.

Usage

r_foursail(rho, tau, bgr, param)

Arguments

rho

input leaf reflectance from 400-2500nm (can be measured or modeled)

tau

input leaf transmittance from 400-2500nm (can be measured or modeled)

bgr

background reflectance. Usual input is soil reflectance spectra from 400-2500nm (can be measured or modeled)

param

A named vector of SAIL parameter values (note: program ignores case):

  • 1 = Leaf angle distribution function parameter a (LIDFa)

  • 2 = Leaf angle distribution function parameter b (LIDFb)

  • 3 = Leaf angle distribution function type (see ?lidfFun)

  • 4 = Leaf area index (LAI)

  • 5 = Hot spot effect parameter (hspot) - The foliage hot spot size parameter is equal to the ratio of the correlation length of leaf projections in the horizontal plane and the canopy height (Verhoef & Bach 2007).

  • 6 = Solar zenith angle (tts)

  • 7 = Observer zenith angle (tto)

  • 8 = Sun-sensor azimuth angle (psi)

Value

spectra matrixwith 4 reflectance factors and canopy transmission for wavelengths 400 to 2500nm:

  • 1 = bi-hemispherical reflectance (rddt). White-sky albedo: the reflectance of the canopy under diffuse illumination. The BRDF integrated over all viewing and illumination directions.

  • 2 = hemispherical directional reflectance (rsdt). Black-sky albedo: reflectance of a surface under direct light without a diffuse component. It is the integral of the BRDF over all viewing directions.

  • 3 = directional hemispherical reflectance (rdot). Diffuse reflectance in the vieweing direction.

  • 4 = bi-directional reflectance (rsot). The ratio of reflected radiance in the viewing direction to the incoming radiant flux in the solar direction.

  • 5 = Canopy transmission of diffuse light through the canopy (taud).

  • 6 = transmission of direct light through the canopy (taus).

Author(s)

Marco D. Visser (R implementation)


run a requested RTM (internal function)

Description

List of aliases: prospect5, prospectd, prosail5, prosaild, prosail2_55,prosail2_dd, prosail2_5d, prosail2_d5, rtm.inform5, rtm.informd

Usage

runRTM(pars)

Arguments

pars

the required parameters (vector or list)

modReq

model request object built in fRTM

Value

prediction from the requested model


The SAIL BDRF function

Description

The SAIL BDRF function

Usage

sail_BDRF(
  w,
  lai,
  sumint,
  tsstoo,
  rsoil,
  rdd,
  tdd,
  tsd,
  rsd,
  tdo,
  rdo,
  tss,
  too,
  rsod
)

Arguments

w

goemeric reflectance parameter

lai

leaf area index

sumint

exp int

tsstoo

Bi-directional gap fraction

rsoil

background reflectance

rdd

Bi-hemispherical reflectance over all in & outgoing directions (white-sky albedo).

tdd

Bi-hemispherical transmittance (diffuse transmittance in all directions)

tsd

Directional hemispherical transmittance for solar flux

rsd

Directional hemispherical reflectance for solar flux (diffuse solar angle)

tdo

Directional hemispherical transmittance (diffuse in viewing direction)

rdo

Directional hemispherical reflectance in viewing direction

tss

Direct transmittance of solar flux

too

Direct transmittance in viewing direction

rsod

Multi scattering contribution

Value

spectra matrixwith 4 reflectance factors and canopy transmission for wavelengths 400 to 2500nm:

  • 1 = bi-hemispherical reflectance (rddt). White-sky albedo: the reflectance of the canopy under diffuse illumination. The BRDF integrated over all viewing and illumination directions.

  • 2 = hemispherical directional reflectance (rsdt). Black-sky albedo: reflectance of a surface under direct light without a diffuse component. It is the integral of the BRDF over all viewing directions.

  • 3 = directional hemispherical reflectance (rdot). Diffuse reflectance in the vieweing direction.

  • 4 = bi-directional reflectance (rsot). The ratio of reflected radiance in the viewing direction to the incoming radiant flux in the solar direction.

  • 5 = Canopy transmission of diffuse light through the canopy (taud).

  • 6 = transmission of direct light through the canopy in the solar direction (taus).

  • 7 = transmission of direct light through the canopy in the sensor/viewing direction (tauo).


Sky light model

Description

Simple skyl atmospheric model.

Usage

skyl(rddt, rsdt, rdot, rsot, Es, Ed, tts, skyl = NULL)

Arguments

rddt

Bi-hemispherical reflectance

rsdt

Directional-hemispherical reflectance for solar incident flux

rdot

Hemispherical-directional reflectance in viewing direction

rsot

Bi-directional reflectance factor

Es

Solar flux

Ed

Diffuse flux

tts

solar angle

skyl

diffuse fraction, if NULL skyl is estimated using the tts (solar angle).

Details

The version implemented here can also include a dependence of the sun zenith angle after Danner et al. (2019) who build on recommendations from Francois et al. (2002).

Value

a list with hemispherical and directional reflectance.

References

Francois, C., Ottle, C., Olioso, A., Prevot, L., Bruguier, N., Ducros, Y.(2002). Conversion of 400-1100 nm vegetation albedo measurements into total shortwave broadband albedo using a canopy radiative transfer model. Agronomie 22, 611-618.

Danner M, Berger K, Wocher M, Mauser W, Hank T. Fitted PROSAIL Parameterization of Leaf Inclinations, Water Content and Brown Pigment Content for Winter Wheat and Maize Canopies. Remote Sensing. 2019; 11(10):1150.

Examples

data(solar)
rt<-fRTM(rho~prospect5+foursail)
skyl(rt[,"rddt"],rt[,"rsdt"],rt[,"rdot"],rt[,"rsot"],
Es=solar[,1],Ed=solar[,2],tts=45,skyl=NULL)

soil reflectance

Description

soil reflectance

Usage

data(soil)

details

  • 1 = wet soil

  • 2 = dry soil

references

Feret et al. (2008), PROSPECT-4 and 5: Advances in the Leaf Optical Properties Model Separating Photosynthetic Pigments, Remote Sensing of Environment


direct and diffuse light

Description

direct and diffuse light

Usage

data(solar)

details

  • 1 = direct light

  • 2 = diffuse light

references

Feret et al. (2008), PROSPECT-4 and 5: Advances in the Leaf Optical Properties Model Separating Photosynthetic Pigments, Remote Sensing of Environment


RTM inversion summary

Description

RTM inversion summary

Usage

## S3 method for class 'rtm.inversion'
summary(x, ...)

Arguments

x

predictions from an RTM

...

additional plot arguments

Value

summarizes the inverted parameters