naginterfaces.library.tsa.multi_kalman_sqrt_invar¶
- naginterfaces.library.tsa.multi_kalman_sqrt_invar(transf, a, b, stq, q, c, r, s, tol=0.0)[source]¶
multi_kalman_sqrt_invar
performs a combined measurement and time update of one iteration of the time-invariant Kalman filter using a square root covariance filter.For full information please refer to the NAG Library document for g13eb
https://support.nag.com/numeric/nl/nagdoc_30.3/flhtml/g13/g13ebf.html
- Parameters
- transfstr, length 1
Indicates whether to transform the input matrix pair to lower observer Hessenberg form. The transformation will only be required on the first call to
multi_kalman_sqrt_invar
.The matrices in arrays and are transformed to lower observer Hessenberg form and the matrices in and are transformed as described in Notes.
The matrices in arrays , and should be as returned from a previous call to
multi_kalman_sqrt_invar
with .- afloat, array-like, shape
If , the state transition matrix, .
If , the transformed matrix as returned by a previous call to
multi_kalman_sqrt_invar
with .- bfloat, array-like, shape
If , the noise coefficient matrix .
If , the transformed matrix as returned by a previous call to
multi_kalman_sqrt_invar
with .- 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
If , the measurement coefficient matrix, .
If , the transformed matrix as returned by a previous call to
multi_kalman_sqrt_invar
with .- rfloat, array-like, shape
The lower triangular Cholesky factor of the measurement noise covariance matrix .
- sfloat, array-like, shape
If the lower triangular Cholesky factor of the state covariance matrix, .
If the lower triangular Cholesky factor of the covariance matrix of the transformed state vector as returned from a previous call to
multi_kalman_sqrt_invar
with .- 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
- afloat, ndarray, shape
If , the transformed matrix, , otherwise is unchanged.
- bfloat, ndarray, shape
If , the transformed matrix, , otherwise is unchanged.
- cfloat, ndarray, shape
If , the transformed matrix, , otherwise is unchanged.
- sfloat, ndarray, shape
The lower triangular Cholesky factor of the transformed state covariance matrix, .
- kfloat, ndarray, shape
The Kalman gain matrix for the transformed state vector premultiplied by the state transformed transition matrix, .
- hfloat, ndarray, shape
The lower triangular matrix .
- ufloat, ndarray, shape
If the transformation matrix , otherwise is not referenced.
- Raises
- NagValueError
- (errno )
On entry, .
Constraint: or .
- (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 matrices and are time invariant.
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
where is the Kalman gain matrix. The second stage is the time-update for , which is given by
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 triangularization is carried out via Householder transformations exploiting the zero pattern of the pre-array. The relationship between the Kalman gain matrix and is given by
In order to exploit the invariant parts of the model to simplify the computation of the results for the transformed state space are computed where is the transformation that reduces the matrix pair to lower observer Hessenberg form. That is, the matrix is computed such that the compound matrix
is a lower trapezoidal matrix. Further the matrix is transformed to . These transformations need only be computed once at the start of a series, and
multi_kalman_sqrt_invar
will, optionally, compute them.multi_kalman_sqrt_invar
returns transformed matrices , , and , the Cholesky factor of the updated transformed state covariance matrix (where ) and the matrix , valid for both transformed and original models, which is used in the computation of the likelihood for the model. Note that the covariance matrices and can be time-varying.
- 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