naginterfaces.library.tsa.multi_​kalman_​sqrt_​var

naginterfaces.library.tsa.multi_kalman_sqrt_var(a, b, stq, q, c, r, s, tol=0.0)[source]

multi_kalman_sqrt_var performs a combined measurement and time update of one iteration of the time-varying Kalman filter using a square root covariance filter.

For full information please refer to the NAG Library document for g13ea

https://support.nag.com/numeric/nl/nagdoc_30.2/flhtml/g13/g13eaf.html

Parameters
afloat, array-like, shape

The state transition matrix, .

bfloat, array-like, shape

The noise coefficient matrix .

stqbool

If , the state noise covariance matrix is assumed to be the identity matrix. Otherwise the lower triangular Cholesky factor, , must be provided in .

qfloat, array-like, shape

Note: the required extent for this argument in dimension 1 is determined as follows: if : ; otherwise: .

Note: the required extent for this argument in dimension 2 is determined as follows: if : ; otherwise: .

If , must contain the lower triangular Cholesky factor of the state noise covariance matrix, . Otherwise is not referenced.

cfloat, array-like, shape

The measurement coefficient matrix, .

rfloat, array-like, shape

The lower triangular Cholesky factor of the measurement noise covariance matrix .

sfloat, array-like, shape

The lower triangular Cholesky factor of the state covariance matrix, .

tolfloat, optional

The tolerance used to test for the singularity of . If , then is used instead. The inverse of the condition number of is estimated by a call to lapacklin.dtrcon. If this estimate is less than then is assumed to be singular.

Returns
sfloat, ndarray, shape

The lower triangular Cholesky factor of the state covariance matrix, .

kfloat, ndarray, shape

The Kalman gain matrix, , premultiplied by the state transition matrix, , .

hfloat, ndarray, shape

The lower triangular matrix .

Raises
NagValueError
(errno )

On entry, .

Constraint: .

(errno )

On entry, .

Constraint: .

(errno )

On entry, .

Constraint: .

(errno )

On entry, .

Constraint: .

(errno )

The matrix is singular.

Notes

In the NAG Library the traditional C interface for this routine uses a different algorithmic base. Please contact NAG if you have any questions about compatibility.

The Kalman filter arises from the state space model given by:

where is the state vector of length at time , is the observation vector of length at time , and of length and of length are the independent state noise and measurement noise respectively.

The estimate of given observations to is denoted by with state covariance matrix , while the estimate of given observations to is denoted by with covariance matrix . The update of the estimate, , from time to time , is computed in two stages. First, the measurement-update is given by

and

where is the Kalman gain matrix. The second stage is the time-update for which is given by

and

where represents any deterministic control used.

The square root covariance filter algorithm provides a stable method for computing the Kalman gain matrix and the state covariance matrix. The algorithm can be summarised as

where is an orthogonal transformation triangularizing the left-hand pre-array to produce the right-hand post-array. The relationship between the Kalman gain matrix, , and is given by

multi_kalman_sqrt_var requires the input of the lower triangular Cholesky factors of the noise covariance matrices and, optionally, and the lower triangular Cholesky factor of the current state covariance matrix, , and returns the product of the matrices and , , the Cholesky factor of the updated state covariance matrix and the matrix used in the computation of the likelihood for the model.

References

Vanbegin, M, van Dooren, P and Verhaegen, M H G, 1989, Algorithm 675: FORTRAN subroutines for computing the square root covariance filter and square root information filter in dense or Hessenberg forms, ACM Trans. Math. Software (15), 243–256

Verhaegen, M H G and van Dooren, P, 1986, Numerical aspects of different Kalman filter implementations, IEEE Trans. Auto. Contr. (AC-31), 907–917