infer.timedeppar {timedeppar} | R Documentation |
Jointly infer constant and time-dependent parameters of a dynamic model given time-series data
Description
This function draws a Markov Chain from the posterior of constant and time-dependent parameters
(following Ornstein-Uhlenbeck processes) of a dynamic model.
The dynamic model is specified by a function that calculates the log likelihood for given time-series
data.
The Ornstein-Uhlenbeck processes of time-dependent processes are characterized by there mean (mean
),
standard deviation (sd
), and a rate parameter (gamma
) that quantifies temporal correlation.
Usage
infer.timedeppar(
loglikeli = NULL,
loglikeli.keepstate = FALSE,
param.ini = list(),
param.range = list(),
param.log = logical(0),
param.logprior = NULL,
param.ou.ini = numeric(0),
param.ou.fixed = numeric(0),
param.ou.logprior = NULL,
task = c("start", "continue", "restart"),
n.iter = NA,
cov.prop.const.ini = NA,
cov.prop.ou.ini = NA,
scale.prop.const.ini = NA,
scale.prop.ou.ini = NA,
control = list(),
res.infer.timedeppar = list(),
verbose = 0,
file.save = "",
...
)
Arguments
loglikeli |
function that calculates the log likelihood of the model for given constant or
time-dependent parameters and given observational data. |
loglikeli.keepstate |
boolean to indicate which kind of interface to the likelihood function is used.
See argument |
param.ini |
named list of initial vallues of parameters to be estimated.
scalar initial values for constant parameters,
two-column matrices for time and parameter values for time-dependent parameters
(values of time-dependent parameters may be NA and are then drawn from the
Ornstein-Uhlenbeck process).
The list |
param.range |
named list of ranges (2 element vectors with minimum and maximum) of parameters that are constrained (non-logarithmic for all parameters) |
param.log |
named vector of logicals indicating if inference should be done on the log scale (parameters are still given and returned on non-log scales). For time-dependent parameters, selecting this option implies the use of a lognormal marginal for the Ornstein-Uhlenbeck process. This means that the parameter is modelled as exp(Ornstein-Uhlenbeck), but mean and standard devistion of the process are still on non-log scales. |
param.logprior |
function to calculate the (joint) log prior of all estimaged constant parameters of the model. The function gets as its argument a named vector of the values of the estimaged constant parameters to allow the function to identify for which parameters a joint prior is required in the current setting). |
param.ou.ini |
named vector of initial values of parameters of the Ornstein-Uhlenbeck processes of
time-dependent parameters; see description of argument |
param.ou.fixed |
named vector of values of parameters of the Ornstein-Uhlenbeck processes of
time-dependent parameters that are kept fixed rather than being extimated.
If all process parameters are kept fixed, these names are
|
param.ou.logprior |
function to calculate the (joint) log prior of all estimated parameters of
the Ornstein-Uhlenbeck processes of a single time-dependent parameter.
The function gets as its argument a named vector of the values of the process
parameters to estimated.
These names are a subset of
|
task |
Which task to perform (default value: "start"): |
n.iter |
number of iterations of the Markov chain to be performed (default value: 10000). |
cov.prop.const.ini |
scaled covariance matrix of the proposal distribution for the Metropolis step of
constant parameters.
The proposal distribution of the Metropolis step is a normal distribution centered
at the last point of the chain with a covariance matrix equal to
|
cov.prop.ou.ini |
list of scaled covariance matrices of the proposal distributions for the Metropolis
step of the parameters of the Ornstein-Uhlenbeck processes of time-dependent
parameters.
The proposal distribution of the Metropolis step for the process parameters of the
time-dependent parameter i is a normal distribution centered
at the last point of the chain with a covariance matrix equal to
|
scale.prop.const.ini |
scale factor for the covariance matrix of the Metropolis step for constant parameters
with a proposal distribution equal to a normal distribution centered at the previous
point of the chain and a covariance matrix equal to
|
scale.prop.ou.ini |
vector of scale factors for the covariance matrices of the Metropolis step for parameters
of Ornstein-Uhlenbeck processes with a proposal distribution equal to a normal
distribution centered at the previous point of the chain and a covariance matrix
for the time dependent parameter i equal to |
control |
list of control parameters of the algorithm:
|
res.infer.timedeppar |
results of a previous call to this function.
These results are ignored if the argument |
verbose |
integer parameter indicating the level of progress reporting: |
file.save |
if non-empty string, the intermediate results are saved to this file as variable |
... |
additional parameters passed to the function |
Value
class of type timedeppar
with the following elements:
-
package
: package timedeppar: version and date, -
func
: function called (infer.timedeppar), -
date
: date of call, -
dot.args
: arguments passed to the likelihood function (included for reproducibility of results), -
task
: task that was performed (start
,restart
orcontinue
), -
file
: name of file to which output was written, -
param.ini
: initial values of likelihood parameters (constant and time-dependent), -
param.ou.ini
: initial values of Ornstein-Uhlenbeck process parameters that are estimated, -
param.ou.fixed
: values of Ornstein-Uhlenbeck process parameters that are not estimated, -
loglikeli
: function that was passed to calculate the log likelihood of the observations, -
loglikeli.keepstate
: boolean indicating whether or not the state from the previous run should be kept (this allows only partial time evaluation when only part of the input was replaced), -
param.logprior
: function that was passed to calculate the joint log prior of the constant likelihood parameters, -
param.ou.logprior
: function that was passed to calculate the joint log prior of the estimated Ornstein-Uhlenbeck process parameters (in case of multiple Ornstein-Uhlenbeck processes the function has to return the prior for the correct process; this can be identified by the names of the argument), -
param.range
: parameter ranges, -
param.log
: named logical vector of indicators for log inference, -
control
: named list of control parameters as passed to the call (or read from a previous call), -
n.iter
: number of iterations peformed (note that the size of the sample will be n.iter/control$thin), -
sample.diag
: list of samples of proposals, log acceptance ratios, and interval lengths of time-dependent parameters (only available if the control variablesave.diag
is set toTRUE
), -
sample.param.timedep
: list of samples of time dependent parameters (first row contains time points), -
sample.param.ou
: sample of Ornstein-Uhlenbeck process parameters, -
sample.param.const
: sample of constant parameters, -
sample.logpdf
: sample of prior, Ornstein-Uhlenbeck and posterior pdf, -
acceptfreq.constpar
: acceptance frequency of constant parameters after adaptation phase, -
acceptfreq.oupar
: acceptance frequencies of Ornstein-Uhlenbeck process parameters after adaptation phase, -
acceptfreq.timedeppar
: acceptance frequencies of time-depenent parameters, -
param.maxpost
: parameters at the maximum posterior (constant and time-dependent parameters), -
param.ou.maxpost
: Ornstein-Uhlenbeck process parameters at the maximum posterior, -
cov.prop.const
: final covariance matrix used for proposal distribution of constant parameters, -
cov.prop.ou
: list of final covariance matrices used for proposal distribution of Ornstein-Uhlenbeck process paramemters, -
scale.prop.const
: final scale of proposal distribution of constant parameters, -
scale.prop.ou
: final scale of proposal distribution of Ornstein-Uhlenbeck process parameters, -
sys.time
: run time used for the previous inference job.
References
Reichert, P.
timedeppar: An R package for inferring stochastic, time-dependent model parameters
in preparation, 2020.
Reichert, P., Ammann, L. and Fenicia, F.
Potential and challenges of investigating intrinsic uncertainty of hydrological models with time-dependent, stochastic parameters.
Water Resources Research 57(8), e2020WR028311, 2021.
doi:10.1029/2020WR028311
Reichert, P. and Mieleitner, J.
Analyzing input and structural uncertainty of nonlinear dynamic models with stochastic, time-dependent parameters.
Water Resources Research, 45, W10402, 2009.
doi:10.1029/2009WR007814
Tomassini, L., Reichert, P., Kuensch, H.-R. Buser, C., Knutti, R. and Borsuk, M.E.
A smoothing algorithm for estimating stochastic, continuous-time model parameters
and its application to a simple climate model.
Journal of the Royal Statistical Society: Series C (Applied Statistics) 58, 679-704, 2009.
doi:10.1111/j.1467-9876.2009.00678.x
See Also
plot.timedeppar
for visualizing results.
calc.acceptfreq
for calculating (apparent) acceptance frequencies.
calc.logpdf
for calculating log pdf values (prior, internal, posterior) from the results.
get.param
for extracting individual parameters from the Markov chain.
get.parsamp
for extracting subsamples of the Markov chain.
readres.timedeppar
for reading saved results from a previous run.
randOU
for sampling from an Ornstein-Uhlenbeck process.
logpdfOU
for calculating the probability density of a sample from an Ornstein-Uhlenbeck process.
Examples
# Simple example for re-inferring parameters of an Ornstein-Uhlenbeck process
# with observational noise from synthetically generated data
# ---------------------------------------------------------------------------
# load package:
if ( !require("timedeppar") ) { install.packages("timedeppar"); library(timedeppar) }
# choose model parameters:
y_mean <- 0
y_sd <- 1
y_gamma <- 10
obs_sd <- 0.2
# choose control parameters of numerical algorithm:
n.iter <- 100 # this is just to demonstrate how it works and is compatible
# with the computation time requirements for examples in CRAN
# n.iter <- 50000 # please go for a sample size like this
# # for getting a reasonable sample
n.interval <- 25 # increase if rejection frequency of stoch. par. too high
fract.adapt <- 0.4
n.adapt <- floor(fract.adapt*n.iter)
# synthetically generate data:
set.seed(123)
data <- randOU(mean=y_mean,sd=y_sd,gamma=y_gamma,t=seq(from=0,to=2,length.out=101))
data$yobs <- data$y + rnorm(nrow(data),mean=0,sd=obs_sd)
# define observational likelihood:
loglikeli <- function(param,data)
{
# get parameter y at time points of observations::
y <- param$y
if ( is.matrix(y) | is.data.frame(y) ) y <- approx(x=y[,1],y=y[,2],xout=data[,1])$y
# calculate likelihood:
loglikeli <- sum(dnorm(data[,"yobs"],mean=y,sd=param$obs_sd,log=TRUE))
# return result:
return(loglikeli)
}
# sample from the posterior of y, mu_y, sd_y and sd_obs assuming a uniform prior:
res <- infer.timedeppar(loglikeli = loglikeli,
param.ini = list(y=randOU(mean=y_mean,sd=y_sd,gamma=y_gamma,
t=seq(from=0,to=2,length.out=501)),
obs_sd=obs_sd),
param.log = c(y=FALSE,obs_sd=TRUE),
param.ou.ini = c(y_mean=0,y_sd=1),
param.ou.fixed = c(y_gamma=10),
n.iter = n.iter,
control = list(n.interval = n.interval,
n.adapt = n.adapt),
data = data)
# plot results using pre-defined options:
# pdf(paste0("infer_OU_",n.iter,"_",n.adapt,"_original.pdf"),width=8,height=12)
plot(res,
labels=expression(y = y,
y_mean = mu[y],
y_sd = sigma[y],
y_gamma = gamma[y]))
# dev.off()
# plot time series and data:
# pdf(paste0("infer_OU_",n.iter,"_",n.adapt,"_comparison.pdf"),width=8,height=6)
t <- res$sample.param.timedep$y[1,]
sample <- res$sample.param.timedep$y[(1+n.adapt+1):nrow(res$sample.param.timedep$y),]
q <- apply(sample,2,quantile,probs=c(0.025,0.5,0.975))
plot(numeric(0),numeric(0),type="n",xaxs="i",yaxs="i",
xlim=range(t),ylim=2.5*c(-1,1),xlab="t",ylab="y")
polygon(c(t,rev(t),t[1]),c(q[1,],rev(q[3,]),q[1,1]),
col="grey80",border=NA)
lines(t,q[2,])
lines(data$t,data$y,col="red")
points(data$t,data$yobs,pch=19,cex=0.8)
legend("bottomright",
legend=c("original process","noisy data","inferred median","inferred 95% range"),
lwd=c(1,NA,1,5),lty=c("solid",NA,"solid","solid"),col=c("red","black","black","grey80"),
pch=c(NA,19,NA,NA),cex=0.8)
# dev.off()