Base class for CSI models.

Details

Let Xt denote an N x d matrix of particle positions recorded at intervals of dt, where d = 1,2,3 is the number of recorded dimensions. A CSI model for Xt is of the form

Xt = R(phi) mu + Sigma^{1/2} Z,

where:

  • R(phi) is an N x p matrix of drift terms, possibly dependent on a parameter phi. In most cases though, R(phi) = t((1:N-1) * dt) is used to model linear drift.

  • mu is a p x d matrix of drift parameters. For linear drift, it represents the drift velocity in each of the d dimensions.

  • Sigma is a d x d variance matrix.

  • Z is an N x d matrix where each column is an iid realization of N values of a Gaussian continuous stationary increments (CSI) process. Such processes are completely determined by the mean square displacement (MSD) function

    MSD_Z(h) = E[ (Z_{i+h,j}, Z_{i,j})^2 ],

    where the MSD function also depends on the parameter phi.

The csi_model class is a base class for CSI models which requires the user to specify the drift and msd functions, based on which it provides generic methods for parameter inference, simulation, etc. For further details about the CSI model see vignette("subdiff").

Public fields

phi_names

Vector of kernel parameter names (on the original scale). Setting to NULL means there are no kernel parameters (n_phi = 0).

Active bindings

Xt

Particle trajectory. A matrix where row n is the position of the particle at time t = n * dt and each column is a measurement dimension. Only reallocates memory for the internal Toeplitz matrix if necessary.

dt

Interobservation time (scalar).

Methods


Method trans()

Transform kernel parameters from regular to computational basis.

Usage

csi_model$trans(phi)

Arguments

phi

Kernel parameters in the original basis.

Returns

Kernel parameters in the computational basis.


Method itrans()

Transform kernel parameters from computational to regular basis.

Usage

csi_model$itrans(psi)

Arguments

psi

Kernel parameters in the computational basis.

Returns

Kernel parameters in the original basis.


Method acf()

Increment autocorrelation function.

Usage

csi_model$acf(phi, dt, N)

Arguments

phi

Kernel parameters in the original basis.

dt

Interobservation time.

N

Number of trajectory increments.

Returns

A vector of N autocorrelations.


Method msd()

Position mean square displacement function.

Usage

csi_model$msd(phi, t)

Arguments

phi

Kernel parameters in the original basis.

t

Vector of time points at which to calculate the MSD.

Details

This method can be directly supplied by the derived class. Otherwise, it uses SuperGauss::acf2msd() to calculate the MSD from self$acf() at intervals of self$dt, and interpolates linearly between these timepoints at the desired values in t.

The self$msd() method is defined this way because the MSD of some CSI models (e.g., fsd and farma) is not defined in continuous time, and thus intrinsically depends on the interobservation time dt. Thus, the default self$msd() method throws an error if self$dt has not yet been set.

Returns

A vector of unscaled MSD values the same length as t. That is, for the drift-subtracted process X_sub(t) = X(t) - drift(t) * mu, the method returns eta(t) where

MSD_{X_sub}(t) = trace(Sigma) * eta(t).


Method drift()

Increment autocorrelation function.

Usage

csi_model$drift(phi, dt, N)

Arguments

phi

Kernel parameters in the original basis.

dt

Interobservation time.

N

Number of trajectory increments.

Returns

An N x n_drift matrix of drift basis functions.


Method get_omega()

Combine kernel parameters with conditional MLE of mu and Sigma.

Usage

csi_model$get_omega(psi)

Arguments

psi

Kernel parameters in the computational basis.

Returns

Named vector of full parameter set in the computational basis.


Method trans_full()

Convert from original to computational basis.

Usage

csi_model$trans_full(phi, mu, Sigma)

Arguments

phi

Kernel parameters in the original basis.

mu

Drift coefficients.

Sigma

Scale matrix.

Returns

Full parameter vector in the computational basis.


Method itrans_full()

Convert from computational to original basis.

Usage

csi_model$itrans_full(omega)

Arguments

omega

Vector of parameters in the computational basis.

Returns

List with elements phi, mu, and Sigma.


Method nlp()

Evaluate the negative profile loglikelihood.

Usage

csi_model$nlp(psi)

Arguments

psi

Kernel parameters in the computational basis.

Returns

Value of the negative profile loglikelihood (scalar).


Method nu_hat()

Conditional MLE of nuisance parameters.

Usage

csi_model$nu_hat(phi)

Arguments

phi

Kernel parameters in the original basis.

Returns

A list with elements mu and Sigma.


Method loglik()

Evaluate the loglikelihood function.

Usage

csi_model$loglik(phi, mu, Sigma)

Arguments

phi

Kernel parameters in the original basis.

mu

Drift coefficients.

Sigma

Scale matrix.

Returns

Value of the loglikelihood (scalar).


Method fisher()

Calculate the observed Fisher information matrix.

Usage

csi_model$fisher(omega)

Arguments

omega

Vector of length n_omega of parameters in the computational basis.

Returns

The n_omega x n_omega observed Fisher information matrix.


Method get_vcov()

Convert a Fisher information matrix to a variance matrix.

Usage

csi_model$get_vcov(fi)

Arguments

fi

Fisher information matrix of size n_omega x n_omega with parameters in the computational basis.

Returns

Variance matrix of size n_omega x n_omega with parameters in the computational basis.


Method fit()

Calculate the maximum likelihood parameter values.

Usage

csi_model$fit(psi0, vcov = TRUE, ...)

Arguments

psi0

Vector of kernel parameter values (on the computational scale) to initialize the optimization if n_phi > 1. If n_phi == 1, a vector of length 2 giving the range in which to perform the optimum search.

vcov

Whether to also calculate the MLE variance estimate.

...

Additional arguments to stats::optim() or stats::optimize() for n_phi == 1.

Returns

A list with elements coef and optionally vcov containing the MLE in the computational basis and its variance estimate.


Method resid()

Calculate the model residuals.

Usage

csi_model$resid(phi, mu, Sigma)

Arguments

phi

Kernel parameters in the original basis.

mu

Drift coefficients.

Sigma

Scale matrix.

Returns

A matrix of residuals the same size as dX as calculated with csi_resid(), upon using the model's drift() and acf() specifications.


Method sim()

Simulate trajectories from the model.

Usage

csi_model$sim(phi, mu, Sigma, nsim = 1, fft = TRUE, nkeep, tol = 1e-06)

Arguments

phi

Kernel parameters in the original basis.

mu

Drift coefficients.

Sigma

Scale matrix.

nsim

Number of trajectories to simulate.

fft, nkeep, tol

Optional arguments to SuperGauss::rnormtz().

Returns

A matrix of size dim(Xt) or an array of size nrow(Xt) x ncol(dX) x nsim array when nsim > 1 of simulated trajectories, as calculated with csi_sim(), using the model's drift() and acf() specifications.


Method new()

Model object constructor.

Usage

csi_model$new(Xt, dt, drift = "linear", n_drift)

Arguments

Xt

Matrix of particle positions.

dt

Interobservation time.

drift

Drift specification. Either one of the strings "none", "linear", "quadratric", or a function with signature function(phi, dt, N).

n_drift

Integer number of drift terms. Ignored if drift is one of the default strings. Required otherwise.

Details

The value of n_phi is automatically determined from phi_names. But this means the latter must be set by derived$initialize() before super$initialize() is called. Otherwise an error is thrown.

The constructor can be called without Xt for accessing methods which don't require it, e.g.,

derived$new(dt = dt)$acf(phi, dt, N)

Development Notes

  • Should a validator be (optionally) called after the constructor?

  • Should psi necessarily contain the origin, to facilitate the validator?


Method clone()

The objects of this class are cloneable with this method.

Usage

csi_model$clone(deep = FALSE)

Arguments

deep

Whether to make a deep clone.