calCov {magi} | R Documentation |
Calculate stationary Gaussian process kernel
Description
Covariance calculations for Gaussian process kernels. Currently supports matern, rbf, compact1, periodicMatern, generalMatern, and rationalQuadratic kernels. Can also return m_phi and other additional quantities useful for ODE inference.
Usage
calCov(
phi,
rInput,
signrInput,
bandsize = NULL,
complexity = 3,
kerneltype = "matern",
df,
noiseInjection = 1e-07
)
Arguments
phi |
the kernel hyper-parameters. See details for hyper-parameter specification for each |
rInput |
the distance matrix between all time points s and t, i.e., |s - t| |
signrInput |
the sign matrix of the time differences, i.e., sign(s - t) |
bandsize |
size for band matrix approximation. See details. |
complexity |
integer value for the complexity of the kernel calculations desired:
See details for their definitions. |
kerneltype |
must be one of |
df |
degrees of freedom, for |
noiseInjection |
a small value added to the diagonal elements of C and Kphi for numerical stability |
Details
The covariance formulae and the hyper-parameters phi
for the supported kernels are as follows. Stationary kernels have C(s,t) = C(r)
where r = |s-t|
is the distance between the two time points. Generally, the hyper-parameter phi[1]
controls the overall variance level while phi[2]
controls the bandwidth.
matern
This is the simplified Matern covariance with
df = 5/2
:C(r) = phi[1] * (1 + \sqrt 5 r/phi[2] + 5r^2/(3 phi[2]^2)) * \exp(-\sqrt 5 r/phi[2])
rbf
-
C(r) = phi[1] * \exp(-r^2/(2 phi[2]^2))
compact1
-
C(r) = phi[1] * \max(1-r/phi[2],0)^4 * (4r/phi[2]+1)
periodicMatern
-
Define
r' = | \sin(r \pi/phi[3])*2 |
. Then the covariance is given byC(r')
using the Matern formula. generalMatern
-
C(r) = phi[1] * 2^(1-df) / \Gamma(df) * ( \sqrt(2.0 * df) * r / phi[2] )^df * besselK( \sqrt(2.0 * df) * r / phi[2] , df)
where
besselK
is the modified Bessel function of the second kind. rationalQuadratic
-
C(r) = phi[1] * (1 + r^2/(2 df phi[2]^2))^(-df)
The kernel calculations available and their definitions are as follows:
- C
The covariance matrix corresponding to the distance matrix
rInput
.- Cprime
The cross-covariance matrix
d C(s,t) / ds
.- Cdoubleprime
The cross-covariance matrix
d^2 C(s,t) / ds dt
.- dCdphi
A list with the matrices
dC / dphi
for each element of phi.- Ceigen1over
The reciprocals of the eigenvalues of C.
- CeigenVec
Matrix of eigenvectors of C.
- Cinv
The inverse of C.
- mphi
The matrix
Cprime * Cinv
.- Kphi
The matrix
Cdoubleprime - Cprime * Kinv * t(Cprime)
.- Keigen1over
The reciprocals of the eigenvalues of Kphi.
- Kinv
The inverse of Kphi.
- mphiLeftHalf
The matrix
Cprime * CeigenVec
.- dCdphiCube
dC / dphi
as a 3-D array, with the third dimension corresponding to the elements of phi.
If bandsize
is a positive integer, additionally CinvBand, mphiBand, and KinvBand are provided in the return list, which are
band matrix approximations to Cinv, mphi, and Kinv with the specified bandsize
.
Value
A list containing the kernel calculations included by the value of complexity
.
Examples
foo <- outer(0:40, t(0:40), '-')[, 1, ]
r <- abs(foo)
signr <- -sign(foo)
calCov(c(0.2, 2), r, signr, bandsize = 20, kerneltype = "generalMatern", df = 2.01)