setup {CollocInfer}R Documentation

Setup Functions for proc and lik objects

Description

These functions set up lik and proc objects of squared error and multinormal processes.

Usage

LS.setup(pars,coefs=NULL,fn,basisvals=NULL,lambda,fd.obj=NULL,
        more=NULL,data=NULL,weights=NULL,times=NULL,quadrature=NULL,
        likfn = make.id(), likmore = NULL,eps=1e-6,
        posproc=FALSE,poslik=FALSE,discrete=FALSE,names=NULL,sparse=FALSE)

multinorm.setup(pars,coefs=NULL,fn,basisvals=NULL,var=c(1,0.01),fd.obj=NULL,
        more=NULL,data=NULL,times=NULL,quadrature=NULL,eps=1e-6,posproc=FALSE,
        poslik=FALSE,discrete=FALSE,names=NULL,sparse=FALSE)

Arguments

pars

Initial values of parameters to be estimated processes.

coefs

Vector giving the current estimate of the coefficients in the spline.

fn

A function giving the right hand side of a differential/difference equation. The function should have arguments

  • times The times at which the RHS is being evaluated.

  • x The state values at those times.

  • p Parameters to be entered in the system.

  • more An object containing additional inputs to fn

It should return a matrix of the same dimension of x giving the right hand side values.

If fn is given as a single function, its derivatives are estimated by finite-differencing with stepsize eps. Alternatively, a list can be supplied with elements:

  • fn Function to calculate the right hand side should accept a matrix of state values at .

  • dfdx Function to calculate the derivative with respect to x

  • dfdp Function to calculate the derivative with respect to p

  • d2fdx2 Function to calculate the second derivative with respect to x

  • d2fdxdp Function to calculate the second derivative with respect to x and p

These functions take the same arguments as fn and should output multidimensional arrays with the dimensions ordered according to time, state, deriv1, deriv2; here derivatives with respect to x always precede derivatives with respect to p.

fn can also be given as a pomp object (see the pomp package), in which case it is interfaced to CollocInfer through pomp.skeleton using a finite differencing.

basisvals

Values of the collocation basis to be used. This can either be a basis object from the fda package, or a list elements:

  • bvals.obs A matrix giving the values of the basis at the observation times

  • bvals A matrix giving the values of the basis at the quadrature times

  • dbvals A matrix giving the derivative of the basis at the quadrature times

For discrete systems, it may also be specified as a matrix, in which case bvals$bvals is obtained by deleting the last row and bvals$dbvals is obtained by deleting the first/

If left as NULL, it is taken from fd.obj for discrete=FALSE and defaults to an identity matrix of the same dimension as the number of observations for discrete=TRUE systems.

lambda

(LS.setup only) Penalty value trading off fidelity to data with fidelity to differential equations.

var

(profile.Cproc or profile.Dproc) A vector of length 2, giving

fd.obj

(Optional) A functional data object; if this is non-null, coefs and basisvals is extracted from here.

more

An object specifying additional arguments to fn.

data

The data to be used, this can be a matrix, or a three-dimensional array. If the latter, the middle dimension is taken to be replicates. The data are returned, if replicated they are returned in a concatenated form.

weights

(LS.setup only)

times

Vector observation times for the data. If the data are replicated, times are returned in a concatenated form.

quadrature

Quadrature points, should contain two elements (if not NULL)

  • qpts Quadrature points; defaults to midpoints between knots

  • qwts Quadrature weights; defaults to normalizing by the length of qpts.

eps

Finite differencing step size, if needed.

posproc

Should the state vector be constrained to be positive? If this is the case, the state is represented by an exponentiated basis expansion in the proc object.

poslik

Should the state be exponentiated before being compared to the data? When the state is represented on the log scale TRUE, this is an alternative to taking the log of the data.

discrete

Is this a discrete or continuous-time system?

names

The names of the state variables if not given by the column names of coefs.

sparse

Should sparse matrices be used for basis values? This option can save memory when ProfileGausNewt and SplineEstNewtRaph are called. Otherwise sparse matrices will be converted to full matrices and this can slow the code down.

likfn

Defines a map from the trajectory to the observations. This should be in the same form as fn. If a function is given, derivatives are estimated by finite differencing, otherwise a list is expected to provide the same derivatives as fn. If poslik=TRUE, the states are exponentiated before the likfn is evaluated and the derivatives are updated to account for this. Defaults to the identity transform.

likmore

A list containing additional inputs to likfn if needed, otherwise set to NULL

Details

These functions provide basic setup utilities for the collocation inference methods. They define lik and proc objects for sum of squared errors and multivariate normal log likelihoods with nonlinear transfer functions describing the evolution of the state vector.

Value

A list with elements

coefs

Starting values for coefs

lik

The lik object generated

proc

The proc item generated

data

The data matrix, concatenated if from a 3d array.

times

The vector of observation times, concatenated if data is a 3d array.

See Also

inneropt, outeropt, Profile.LS, Profile.multinorm, Smooth.LS, Smooth.multinorm

Examples


# FitzHugh-Nagumo

t = seq(0,20,0.05)            # Observation times

pars = c(0.2,0.2,3)           # Parameter vector
names(pars) = c('a','b','c')

knots = seq(0,20,0.2)         # Create a basis
norder = 3
nbasis = length(knots) + norder - 2
range = c(0,20)

bbasis = create.bspline.basis(range=range,nbasis=nbasis,norder=norder,breaks=knots)

lambda = 10000               # Penalty value

coefs = matrix(0,nbasis,2)   # Coefficient matrix

profile.obj = LS.setup(pars=pars,coefs=coefs,fn=make.fhn(),basisvals=bbasis,
                       lambda=lambda,times=t)


# Using multinorm

var = c(1,0.01)

profile.obj = multinorm.setup(pars=pars,coefs=coefs,fn=make.fhn(),
                                        basisvals=bbasis,var=var,times=t)


# Henon - discrete

hpars = c(1.4,0.3)
t = 1:200

coefs = matrix(0,200,2)
lambda = 10000

profile.obj = LS.setup(pars=hpars,coefs=coefs,fn=make.Henon(),basisvals=NULL,
                             lambda=lambda,times=t,discrete=TRUE)

[Package CollocInfer version 1.0.4 Index]