Title: | PEcAn Functions Used for Radiative Transfer Modeling |
---|---|
Description: | Functions for performing forward runs and inversions of radiative transfer models (RTMs). Inversions can be performed using maximum likelihood, or more complex hierarchical Bayesian methods. Underlying numerical analyses are optimized for speed using Fortran code. |
Authors: | Mike Dietze [aut], Shawn Serbin [aut], Alexey Shiklomanov [aut, cre], University of Illinois, NCSA [cph] |
Maintainer: | Alexey Shiklomanov <[email protected]> |
License: | BSD_3_clause + file LICENSE |
Version: | 1.7.3.9000 |
Built: | 2024-11-20 22:25:48 UTC |
Source: | https://github.com/PecanProject/pecan |
Select spectra
## S3 method for class 'spectra' spectra[i, j, drop = FALSE]
## S3 method for class 'spectra' spectra[i, j, drop = FALSE]
spectra |
Object of class |
i , j
|
indices specifying elements to extract or replace.
See |
drop |
Coerce result to the lowest dimension possible?
See |
Select spectra by wavelength
## S3 method for class 'spectra' spectra[[wavelength, j]]
## S3 method for class 'spectra' spectra[[wavelength, j]]
spectra |
Object of class |
wavelength |
Wavelength vector to select |
j |
index specifying elements to extract or replace.
See |
Assign values to spectra by wavelength
## S3 replacement method for class 'spectra' spectra[[wavelength, j]] <- value
## S3 replacement method for class 'spectra' spectra[[wavelength, j]] <- value
spectra |
Object of class |
wavelength |
Wavelength vector to select |
j |
index specifying elements to extract or replace.
See |
value |
Vector or matrix of values to assign |
Burn-in and thinning of MCMC samples
burnin.thin( samples, target = 5000, burnin.ratio = 2, auto = TRUE, burnin = NULL, thin = NULL )
burnin.thin( samples, target = 5000, burnin.ratio = 2, auto = TRUE, burnin = NULL, thin = NULL )
samples |
Matrix of MCMC samples |
target |
Target number of samples (default = 5000). Only applicable if auto=TRUE. |
burnin.ratio |
Fraction of samples to burn-in; i.e. 2 means to remove first 1/2 of samples, 3 means 1/3, etc. (default = 2). Only applicable if auto=TRUE. |
auto |
Whether or not to perform automatic burnin and thin based on target number of samples. |
burnin |
Number of samples to discard as burnin (auto must be FALSE) |
thin |
Thinning interval (auto must be FALSE) |
Combine spectra by wavelength
## S3 method for class 'spectra' cbind(...)
## S3 method for class 'spectra' cbind(...)
... |
Spectra to combine |
Uses Gelman multivariate Gelman-Rubin diagnostic to check if multiple MCMC chains have converged
check.convergence(jags_out, threshold = 1.1, verbose = TRUE, ...)
check.convergence(jags_out, threshold = 1.1, verbose = TRUE, ...)
jags_out |
mcmc.list object (from coda package) containing samples from MCMC chains. |
threshold |
Gelman-Rubin diagnostic parameter threshold. Default = 1.1 |
verbose |
If TRUE, print convergence result. Default = TRUE |
... |
Additional arguments to |
List length 3 containing the following: * convergence: Logical. Whether or not convergence was achieved. * diagnostics: Numerical value of Gelman-Rubin diagnostics for each parameter and multivariate diagnostic * error: Logical. Whether or not an error occured in the Gelman-Rubin calculation.
Default settings for PROSPECT inversion
default.settings.prospect
default.settings.prospect
An object of class list
of length 18.
Extract default parameter values from model.list
defparam(modname)
defparam(modname)
modname |
Model name. Must match |
Named vector of default parameter values
Calculates the log density of a univariate truncated normal variable
dtnorm(x, mu, sd, MIN)
dtnorm(x, mu, sd, MIN)
x |
A random variable |
mu |
The mean parameter of the distribution; NOTE this is not equal to the mean |
sd |
The standard deviation parameter of the distribution |
MIN |
Value at which the truncation takes place |
Numeric; log density of the distribution, or -1e15 if the x < MIN
Alexey Shiklomanov
This function provides a convenient way to call the ED radiative transfer module (EDR, which simulates full spectral return of an ED patch for a given point in time) directly from R.
EDR( img_path, ed2in_path, spectra_list, trait.values, soil_reflect_path = system.file("extdata", "soil_reflect_par.dat", package = "PEcAnRTM"), wood_reflect_path = system.file("extdata", "wood_reflect_par.dat", package = "PEcAnRTM"), par.wl = 400:2499, nir.wl = 2500, edr_exe_path = NULL, output.path = dirname(normalizePath(ed2in_path, mustWork = TRUE)), settings = list(model = list(revision = "git", config.header = NULL)), singularity_args = list(), clean = FALSE, stderr = TRUE, verbose_error = TRUE, ... )
EDR( img_path, ed2in_path, spectra_list, trait.values, soil_reflect_path = system.file("extdata", "soil_reflect_par.dat", package = "PEcAnRTM"), wood_reflect_path = system.file("extdata", "wood_reflect_par.dat", package = "PEcAnRTM"), par.wl = 400:2499, nir.wl = 2500, edr_exe_path = NULL, output.path = dirname(normalizePath(ed2in_path, mustWork = TRUE)), settings = list(model = list(revision = "git", config.header = NULL)), singularity_args = list(), clean = FALSE, stderr = TRUE, verbose_error = TRUE, ... )
img_path |
Path to Singularity container (usually a |
ed2in_path |
Path to ED2IN file. |
spectra_list |
List of spectral data matrices. Names must exactly match
the PFTs given in |
trait.values |
Named, hierarchical list of trait values for generating
config.xml file. Names must be PFTs, and must exactly match names of
|
soil_reflect_path |
Path to soil reflectance file (defaults to spectra
in package |
wood_reflect_path |
Path to wood reflectance file (defaults to spectra
in package |
par.wl |
Vector of wavelengths defining PAR region |
nir.wl |
Vector of wavelengths defining NIR region |
edr_exe_path |
If |
output.path |
Directory in which to execute the run. Defaults to
|
settings |
PEcAn settings list. Default is |
singularity_args |
Additional arguments to be passed to |
clean |
Logical. If |
stderr |
Logical. If |
verbose_error |
Logical. If |
... |
Additional arguments to |
Alexey Shiklomanov
Locate history file based on path and prefix, copy to specified output directory, and rename to correct time.
EDR.preprocess.history( history.path, output.path, datetime, history.prefix = "history" )
EDR.preprocess.history( history.path, output.path, datetime, history.prefix = "history" )
history.path |
Path to directory containing history file. |
output.path |
Directory in which to execute the run. |
datetime |
POSIX date and time for run |
history.prefix |
String describing the history file prefix in
|
Convert R list to a Fortran /data/
module block
fortran_data_module(dat, types, modname, fname = paste0(modname, ".f90"))
fortran_data_module(dat, types, modname, fname = paste0(modname, ".f90"))
dat |
Named list object for creating module. List names will be used to initialize FORTRAN variabiles. |
types |
Character vector of FORTAN types (e.g. real*8, integer) |
modname |
Name of the module. We suggest the format 'MOD_yourmodname'. |
fname |
Output file name. Defaults to 'yourmodname.f90' |
For models with large constants (e.g. absorption features in the PROSPECT model), it may be preferable to store these in FORTRAN90 modules. However, manually creating and formatting these files is tedious. This script allows you to automatically generate module files from R lists. It automatically interprets the object lengths as array dimensions (only vectors are supported right now – higher dimension arrays may be in the future) and splits long data into rows of 10. Currently, only numeric data are supported (i.e. no characters).
Alexey Shiklomanov
w <- 3.2 x <- 1:5 y <- 6:15 z <- seq(exp(1), pi, length.out=42) l <- list(x=x, y=y, z=z) ## NOTE that names must be explicitly declared l.types <- c('real','integer', 'real*4', 'real*8') fortran_data_module(l, l.types, 'testmod', file.path(tempdir(), "testmod.f90")) x <- runif(10) y <- rnorm(10) z <- rgamma(10, 3) d <- data.frame(x,y,z) ## NOTE that data.frames are just named lists d.types <- rep('real*8', ncol(d)) fortran_data_module(d, d.types, 'random', file.path(tempdir(), "random.f90"))
w <- 3.2 x <- 1:5 y <- 6:15 z <- seq(exp(1), pi, length.out=42) l <- list(x=x, y=y, z=z) ## NOTE that names must be explicitly declared l.types <- c('real','integer', 'real*4', 'real*8') fortran_data_module(l, l.types, 'testmod', file.path(tempdir(), "testmod.f90")) x <- runif(10) y <- rnorm(10) z <- rgamma(10, 3) d <- data.frame(x,y,z) ## NOTE that data.frames are just named lists d.types <- rep('real*8', ncol(d)) fortran_data_module(d, d.types, 'random', file.path(tempdir(), "random.f90"))
R wrapper for 4SAIL model
foursail(refl, tran, rsoil, param)
foursail(refl, tran, rsoil, param)
refl |
input leaf reflectance from 400-2500nm (can be measured or modeled) |
tran |
input leaf transmittance from 400-2500nm (can be measured or modeled) |
rsoil |
input soil reflectance spectra from 400-2500nm (can be measured or modeled) |
param |
Vector of SAIL parameter values: * LIDFa: Leaf angle distribution function - parameter a * LIDFb: Leaf angle distribution function - parameter b * TypeLIDF: Leaf angle distribution function type (1 or 2) * LAI: Leaf area index * q: Hot spot effect parameter * tts: Solar zenith angle * tto: Observer zenith angle * psi: Sun-sensor azimuth angle |
Spectra matrix (see spectra()
) (2101 x 4) of reflectance factors
for wavelengths 400 to 2500nm:
* bi-hemispherical reflectance (rddt)
* hemispherical directional (rsdt)
* directional hemispherical (rdot)
* bi-directional (rsot)
Shawn Serbin
Alexey Shiklomanov
This is the fundamental physical model underlying the PROSPECT family of leaf RTMs.
generalized_plate_model(k, refractive, N)
generalized_plate_model(k, refractive, N)
k |
Specific absorption coefficient (400 - 2500nm) |
refractive |
Refractive index (400 - 2500nm) |
N |
Effective number of mesophyll layers (see |
Generate autocorrelated spectral noise
generate.noise(n = 2101, sigma = 1e-04, fw = 201, fsd = 6)
generate.noise(n = 2101, sigma = 1e-04, fw = 201, fsd = 6)
n |
Length of output vector (default = 2101) |
sigma |
Gaussian noise standard deviation (default=1e-4) |
fw |
Filter width. Will be coerced to an odd number if even (default = 201). |
fsd |
Scaling factor for filter standard deviation (default = 6) |
Read EDR output
get.EDR.output(path = getwd())
get.EDR.output(path = getwd())
path |
Path to directory containing |
Use samplers from the BayesianTools package to fit models to data. Like
invert.auto
, this will continue to run until convergence is achieved
(based on Gelman diagnostic) and the result has enough samples (as
specified by the user; see Details).
invert_bt(observed, model, prior, custom_settings = list(), loglike = NULL)
invert_bt(observed, model, prior, custom_settings = list(), loglike = NULL)
observed |
Vector of observations. Ignored if |
model |
Function called by log-likelihood. Must be |
prior |
BayesianTools prior object. |
custom_settings |
Nested settings list. See Details. |
loglike |
Custom log likelihood function. If |
custom_settings
is a list of lists, containing the following:
common
– BayesianTools settings common to both the initial and subsequent samples.
init
– BayesianTools settings for just the first round of sampling.
This is most common for the initial number of iterations, which is the
minimum expected for convergence.
loop
– BayesianTools settings for iterations inside the convergence
checking while
loop. This is most commonly for setting a smaller
iteration count than in init
.
other
– Miscellaneous (non-BayesianTools) settings, including:
sampler
– String describing which sampler to use. Default is DEzs
use_mpsrf
– Use the multivariate PSRF to check convergence.
Default is FALSE
because it may be an excessively conservative
diagnostic.
min_samp
– Minimum number of samples after burnin before stopping.
Default is 5000.
max_iter
– Maximum total number of iterations. Default is 1e6.
lag.max
– Maximum lag to use for autocorrelation normalization.
Default is 10 * log10(n)
(same as stats::acf
function).
save_progress
– File name for saving samples between loop
iterations. If NULL
(default), do not save progress samples.
threshold
– Threshold for Gelman PSRF convergence diagnostic. Default is 1.1.
verbose_loglike
– Diagnostic messages in log likelihood output. Default is TRUE.
See the BayesianTools sampler documentation for what can go in the BayesianTools
settings lists.
Inversion with automatic convergence checking
invert.auto( observed, invert.options, return.samples = TRUE, save.samples = NULL, quiet = FALSE, parallel = TRUE, parallel.cores = NULL, parallel.output = "/dev/null" )
invert.auto( observed, invert.options, return.samples = TRUE, save.samples = NULL, quiet = FALSE, parallel = TRUE, parallel.cores = NULL, parallel.output = "/dev/null" )
observed |
Vector, matrix, or data frame (coerced to matrix) of
observed values. For spectral data, wavelengths are rows and spectra are
columns. Dimensions must align with the output of |
invert.options |
Parameters related to inversion. |
return.samples |
Include full samples list in output. Default = |
save.samples |
Save samples to file as the inversion proceeds (useful
for debugging). If |
quiet |
Suppress progress bar and status messages. Default=FALSE |
parallel |
Logical. Whether or not to run multiple chains in parallel
on multiple cores (default = |
parallel.cores |
Number of cores to use for parallelization. If
|
parallel.output |
Filename (or ” for stdout) for printing parallel
outputs. Use with caution. Default = |
Performs an inversion via the invert.custom
function with
multiple chains and automatic convergence checking. Convergence checks are
performed using the multivariate Gelman-Rubin diagnostic.
Parameters specific to invert.auto
are described here.
For the remaining parameters, see invert.custom()
.
model
– The model to be inverted. This should be an R function that
takes params
as input and returns one column of observed
(nrows should be the same). Constants should be implicitly included here.
nchains
– Number of independent chains.
inits.function
– Function for generating initial conditions.
ngibbs.max
– Maximum number of total iterations (per chain). DEFAULT = 5e6
ngibbs.min
– Minimum number of total iterations (per chain). DEFAULT = 5000.
ngibbs.step
– Number of iterations between convergence checks. Default = 1000.
run_first
– Function to run before running sampling. Takes parallel
inputs list containing runID, initial values, and resume (NULL) as an
argument.
calculate.burnin
– If TRUE
, use PEcAn.assim.batch::autoburin
function to calculate burnin. Otherwise, assume burnin is min(niter/2, iter_conv_check)
.
threshold
– Maximum value of the Gelman-Rubin diagnostic for
determining convergence. Default = 1.1
List including results
(summary statistics), samples
(mcmc.list
object, or NA
if return.samples=FALSE
), and other
information.
Performs an inversion of an arbitrary model using a modified Metropolis Hastings algorithm with block sampling. This may be slightly slower than the implementation in Fortran, but is much more customizable, as the model can be any R function.
invert.custom( observed, invert.options, quiet = FALSE, return.resume = FALSE, runID = NULL )
invert.custom( observed, invert.options, quiet = FALSE, return.resume = FALSE, runID = NULL )
observed |
Vector, matrix, or data frame (coerced to matrix) of
observed values. For spectral data, wavelengths are rows and spectra are
columns. Dimensions must align with the output of |
invert.options |
R list object containing inversion settings. See details. |
quiet |
Suppress progress bar and status messages. Default=FALSE |
return.resume |
If |
runID |
Run-unique ID. Useful for parallel runs. Default=NULL |
inversion.options
contains the following:
inits
– Vector of initial values of model parameters to be inverted.
ngibbs
– Number of MCMC iterations
prior.function
– Function for use as prior.
Should take a vector of parameters as input and return a single value – the
sum of their log-densities – as output.
param.mins
– Vector of minimum values for inversion parameters
param.maxs
– Vector of minimum values for inversion parameters
model
– The model to be inverted.
This should be an R function that takes params
and runID
as input and
returns one column of observed
(nrows should be the same).
Constants should be implicitly included here.
adapt
– Number of steps for adapting covariance matrix (i.e. adapt
every 'n' steps). Default=100
adj_min
– Minimum threshold for rescaling Jump standard deviation.
Default = 0.1.
target
– Target acceptance rate. Default=0.234, based on recommendation
for multivariate block sampling in Haario et al. 2001
do.lsq
– Perform least squares optimization first (see invert.lsq
),
and use outputs to initialize Metropolis Hastings.
This may improve mixing time, but risks getting caught in a local minimum.
Default=FALSE
catch_error
– If TRUE
(default), wrap model in tryCatch
to prevent sampling termination on model execution error.
Haario, Heikki; Saksman, Eero; Tamminen, Johanna. An adaptive Metropolis algorithm. Bernoulli 7 (2001), no. 2, 223–242. http://projecteuclid.org/euclid.bj/1080222083.
Performs a least-squares inversion of an arbitrary radiative transfer model
(passed as an R function). The inversion attempts to minimize the sum of
residual least squares between modeled and observed spectra via the
Levenberg-Marquardt algorithm (nls.lm
function from the minpack.lm
package).
invert.lsq(observed, inits, model, lower = NULL, upper = NULL)
invert.lsq(observed, inits, model, lower = NULL, upper = NULL)
observed |
Vector of observations (e.g. a reflectance spectrum). |
inits |
Vector of initial conditions for the parameters. |
model |
An R function that calls the RTM and returns the error to be minimized. Be sure to include constants here. |
lower |
Lower bounds on parameters (default=NULL, which means -Inf). |
upper |
Upper bounds on parameters (default=NULL, which means +Inf). |
Alexey Shiklomanov
Load object from an RData file
load.from.name(filename, filepath = ".")
load.from.name(filename, filepath = ".")
filename |
Full name (without path!) of RData file |
filepath |
Path of RData file (default='.') |
Functions for default priors Lognormal mean parameters
lognorm.mu(mean, sd)
lognorm.mu(mean, sd)
mean |
Sample mean |
sd |
Sample standard deviation |
Lognormal sigma parameter
lognorm.sigma(mean, sd)
lognorm.sigma(mean, sd)
mean |
Sample mean |
sd |
Sample standard deviation |
Matplot generic method
matplot(...)
matplot(...)
... |
arguments passed to methods |
Matplot default method
## Default S3 method: matplot(...)
## Default S3 method: matplot(...)
... |
arguments passed to |
Plot multiple spectra on same graph
## S3 method for class 'spectra' matplot(spectra, ...)
## S3 method for class 'spectra' matplot(spectra, ...)
spectra |
Vector ( |
... |
Additional arguments to |
Calculate effective sample size of vector based on its autocorrelation.
neff(x, ...)
neff(x, ...)
x |
A vector or time series |
... |
additional arguments passed to methods |
Shortcut lists for PROSPECT parameter names
params.prospect4
params.prospect4
An object of class character
of length 4.
PROSPECT 5 parameters
params.prospect5
params.prospect5
An object of class character
of length 5.
PROSPECT 5B parameters
params.prospect5b
params.prospect5b
An object of class character
of length 6.
PROSPECT D parameters
params.prospectd
params.prospectd
An object of class character
of length 7.
Creates a nested list whose components are suitable for passing to EDR.
params2edr(params, sep = ".", prospect = TRUE, version = 5)
params2edr(params, sep = ".", prospect = TRUE, version = 5)
params |
Named parameter vector |
sep |
Separator between PFT name and trait name. Must be a single character (default = "."). |
prospect |
Logical. If |
version |
PROSPECT version |
If prospect = TRUE
, parameters prefixed with prospect_
are passed to
prospect with the specified version
, and other parameters are passed to
trait.values
.
The regular expression defining the separation is greedy, i.e.
temperate.Early_Hardwood.SLA
will separate into temperate.Early_Hardwood
and SLA
(assuming the default sep = "."
). Therefore, it is crucial that
trait names not contain any sep
characters (neither ED nor PROSPECT
parameters should anyway). If this is a problem, use an alternate separator
(e.g. |
).
Note that using sep = "."
allows this function to directly invert the
default behavior of unlist
. That is, calling
unlist(params2edr(params, prospect = FALSE)$trait.values)
will return the input vector of
trait values. This makes unlist
a convenient way to go from a
trait.values
list to a properly formatted params
vector.
Because unused ED parameters in the config.xml
are ignored, the PROSPECT
parameters are saved in the trait.values
object as well, which may be
useful for debugging.
List containing spectra_list
and trait.values
, both objects needed by EDR.
Alexey Shiklomanov
Plot spectra vs. wavelength
## S3 method for class 'spectra' plot(spectra, type = "l", ...)
## S3 method for class 'spectra' plot(spectra, type = "l", ...)
spectra |
Vector ( |
type |
plot style, e.g. "l" for lines, "b" for lines and points |
... |
Additional arguments to |
Neatly print inversion results summary
print_results_summary(output)
print_results_summary(output)
output |
Output from |
Alexey Shiklomanov
Print method for spectra S3 class
## S3 method for class 'spectra' print(spectra, n = 10, ...)
## S3 method for class 'spectra' print(spectra, n = 10, ...)
spectra |
Object of class |
n |
Max number of rows to print (show first |
... |
Additional arguments to |
Default prior parameters for PROSPECT models
prior.defaultvals.prospect(sd.inflate = 3)
prior.defaultvals.prospect(sd.inflate = 3)
sd.inflate |
Standard deviation multiplier (default = 3) |
Default PROSPECT 5 prior function
priorfunc.prospect(pmu, psigma)
priorfunc.prospect(pmu, psigma)
pmu |
Lognormal mu parameter |
psigma |
Lognormal sigma parameter |
Assumes lognormal distribution for all parameters. NOTE that prior on N is shifted by 1.
Coupled PROSPECT-Two-stream model
pro2s(param, prospect.version = 5)
pro2s(param, prospect.version = 5)
param |
Model parameters, in the following order: N, Cab, (Car, Cbrown), Cw, Cm, solar zenith angle, LAI, soil_moisture |
prospect.version |
Version of PROSPECT to use (4, 5, or '5B'; default=5) |
Spectra matrix (see spectra()
) for wavelengths 400 to 2500nm containing the following columns:
alpha.c
– Direct ("collimated") reflectance
alpha.i
– Diffuse ("isotropic") reflectance
Tc
– Direct transmittance to background
Ti
– Diffuse reflectance to background
Ac
– Direct absorbance by canopy
Ai
– Diffuse absorbance by canopy
R wrapper for PRO4SAIL model
pro4sail(param)
pro4sail(param)
param |
Vector of PRO4SAIL parameter values: * N: Effective number of leaf layers (>1) * Cab: Leaf chlorophyll content (ug/cm2) (>0) * Car: Leaf carotenoid content (ug/cm2) (>0) * Cbrown: Leaf brown matter content (ug/cm2) (>0) * Cw: Leaf water content (cm) (>0) * Cm: Leaf dry matter content (ug/cm2) (>0) * LIDFa: Leaf angle distribution function - parameter a * LIDFb: Leaf angle distribution function - parameter b * TypeLIDF: Leaf angle distribution function type (1 or 2) * LAI: Leaf area index * q: Hot spot effect parameter * tts: Solar zenith angle * tto: Observer zenith angle * psi: Sun-sensor azimuth angle * psoil: Fraction of soil moisture |
Spectra matrix (see spectra()
) (2101 x 4) of reflectance factors
for wavelengths 400 to 2500nm:
* bi-hemispherical reflectance (rddt)
* hemispherical directional (rsdt)
* directional hemispherical (rdot)
* bi-directional (rsot)
Alexey Shiklomanov
R wrapper for PRO4SAILD model
pro4saild(param)
pro4saild(param)
param |
Vector of PRO4SAIL parameter values: * N: Effective number of leaf layers (>1) * Cab: Leaf chlorophyll content (ug/cm2) (>0) * Car: Leaf carotenoid content (ug/cm2) (>0) * Canth: Leaf anthocyanin content (ug/cm2) (>0) * Cbrown: Leaf brown matter content (ug/cm2) (>0) * Cw: Leaf water content (cm) (>0) * Cm: Leaf dry matter content (ug/cm2) (>0) * LIDFa: Leaf angle distribution function - parameter a * LIDFb: Leaf angle distribution function - parameter b * TypeLIDF: Leaf angle distribution function type (1 or 2) * LAI: Leaf area index * q: Hot spot effect parameter * tts: Solar zenith angle * tto: Observer zenith angle * psi: Sun-sensor azimuth angle * psoil: Fraction of soil moisture |
Spectra matrix (see spectra()
) (2101 x 4) of reflectance factors
for wavelengths 400 to 2500nm:
* bi-hemispherical reflectance (rddt)
* hemispherical directional (rsdt)
* directional hemispherical (rdot)
* bi-directional (rsot)
Alexey Shiklomanov, Shawn Serbin
R wrapper for PROSPECT models
prospect(param, version)
prospect(param, version)
param |
Vector of PROSPECT parameter values: * N: Effective number of leaf layers (>1) * Cab: Leaf chlorophyll content (ug/cm2) (>0) * (5) Car: Leaf carotenoid content (ug/cm2) (>0) * (D) Canth: Leaf anthocyanin content (ug/cm2) (>0) * (5B, D) Cbrown: Leaf brown matter content (ug/cm2) (>0) * Cw: Leaf water content (cm) (>0) * Cm: Leaf dry matter content (ug/cm2) (>0) |
version |
PROSPECT version: 4, 5, or '5B' |
Object of class spectra
(see spectra()
) with simulated
reflectance (column 1) and transmittance (column 2) from 400 to 2500 nm
Alexey Shiklomanov
Quick BayesianTools prior creator for PROSPECT model
prospect_bt_prior(version, custom_prior = list())
prospect_bt_prior(version, custom_prior = list())
version |
PROSPECT version: 4, 5, or '5B' |
custom_prior |
List containing |
Read and process RSR data from directory
read.rsr.folder(dir.path, type)
read.rsr.folder(dir.path, type)
dir.path |
Directory containing RSR data |
type |
Type of sensor. Options are: landsat, avhrr |
Convenient wrapper around base R's splinefun
and approxfun
. See
stats::splinefun()
and stats::approxfun()
documentation for information
on different spline methods and additional arguments.
resample(values, ...) ## Default S3 method: resample(values, from, to, method = "fmm", ...) ## S3 method for class 'matrix' resample(values, from, to, ...) ## S3 method for class 'spectra' resample(values, to, ...)
resample(values, ...) ## Default S3 method: resample(values, from, to, method = "fmm", ...) ## S3 method for class 'matrix' resample(values, from, to, ...) ## S3 method for class 'spectra' resample(values, to, ...)
values |
Vector or matrix of values, or object of class |
... |
Additional arguments to |
from |
X values for interpolation (for |
to |
Y values onto which to interpolate. For |
method |
One of the methods for |
Object of the same class as values
, resampled to the to
values.
Generate relative spectral response (RSR) matrix based on FWHM data
rsr.from.fwhm(wavelength, fwhm)
rsr.from.fwhm(wavelength, fwhm)
wavelength |
Vector of average band widths, as reported in FWHM data. |
fwhm |
Vector of full-width half maximum (FWHM) bandwidths, as reported in FWHM data. |
Generic log-likelihood generator for RTMs
rtm_loglike(nparams, model, observed, lag.max = NULL, verbose = TRUE, ...)
rtm_loglike(nparams, model, observed, lag.max = NULL, verbose = TRUE, ...)
nparams |
number of parameters in model |
model |
function to minimize |
observed |
vector of observations |
lag.max |
passed to |
verbose |
logical: print extra diagnostics during run? |
... |
additional arguments passed to model function |
Random sampling from one-sided truncated normal distribution
rtnorm(mu, sd, MIN)
rtnorm(mu, sd, MIN)
mu |
The mean parameter of the truncated normal. NOTE that this is NOT equal to the mean of the distribution |
sd |
The standard deviation parameter of the truncated normal. Again, NOTE tht this is not the SD of the distribution. |
MIN |
The minimum value, which defines the truncation. |
Draws a random number and, if it doesn't fall within the specified range, resample using an adjusted Normal CDF. This isn't performed immediately because CDF sampling calls three functions – qnorm, runif, and pnorm–and therefore is much less efficient than a simple random sample.
Numeric length one.
Alexey Shiklomanov
Sensor spectral response functions
sensor.list
sensor.list
An object of class character
of length 12.
Sensor list with proper names
sensor.proper
sensor.proper
An object of class character
of length 12.
Using an existing ED2IN file as a template, create a new ED2IN and history file configured for running EDR.
setup_edr( ed2in, output_dir, datetime = ISOdatetime(ed2in[["IYEARA"]], ed2in[["IMONTHA"]], ed2in[["IDATEA"]], 12, 0, 0, tz = "UTC"), ... )
setup_edr( ed2in, output_dir, datetime = ISOdatetime(ed2in[["IYEARA"]], ed2in[["IMONTHA"]], ed2in[["IDATEA"]], 12, 0, 0, tz = "UTC"), ... )
ed2in |
ED2IN list object (see PEcAn.ED2::read_ed2in). |
output_dir |
Directory in which run files will be stored |
datetime |
Date time object (or compliant string) at which to run EDR. Defaults to 12 noon on start date in ED2IN. |
... |
Additional arguments passed on to |
Path to EDR-configured ED2IN file.
Alexey Shiklomanov
Spectra S3 class
spectra(spectra, wavelengths = 400:2500) is_spectra(spectra)
spectra(spectra, wavelengths = 400:2500) is_spectra(spectra)
spectra |
Vector ( |
wavelengths |
Wavelengths of spectra. |
is_spectra()
: Test if object is class spectra
Convolution of spectra to sensor RSR
spectral.response(spec, sensor)
spectral.response(spec, sensor)
spec |
Full (1 nm) spectrum (vector) |
sensor |
Sensor name (string). See sensor.list |
spectra
objectStructure of spectra
object
## S3 method for class 'spectra' str(spectra, ...)
## S3 method for class 'spectra' str(spectra, ...)
spectra |
Object of class |
... |
additional arguments, currently ignored |
Fit multivariate normal to samples. Return means and covariance matrix as a long list (for easy construction of data.tables)
summary_mvnorm(samples)
summary_mvnorm(samples)
samples |
Matrix of MCMC samples. |
Calculate simple univariate summary statistics and return as named list
summary_simple(samples)
summary_simple(samples)
samples |
Matrix of MCMC samples |
Trim RSR matrix to wavelength limits
trim.rsr(rsr, wl.min = 400, wl.max = 2500)
trim.rsr(rsr, wl.min = 400, wl.max = 2500)
rsr |
RSR matrix |
wl.min |
Minimum wavelength (inclusive, default = 400) |
wl.max |
Maximum wavelength (inclusive, default = 2500) |
Retrieve wavelengths from spectra object
wavelengths(spectra)
wavelengths(spectra)
spectra |
Object of class |