G13EAF (PDF version)
G13 Chapter Contents
G13 Chapter Introduction
NAG Library Manual

NAG Library Routine Document

G13EAF

Note:  before using this routine, please read the Users' Note for your implementation to check the interpretation of bold italicised terms and other implementation-dependent details.

+ Contents

    1  Purpose
    7  Accuracy

1  Purpose

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

2  Specification

SUBROUTINE G13EAF ( N, M, L, A, LDS, B, STQ, Q, LDQ, C, LDM, R, S, K, H, TOL, IWK, WK, IFAIL)
INTEGER  N, M, L, LDS, LDQ, LDM, IWK(M), IFAIL
REAL (KIND=nag_wp)  A(LDS,N), B(LDS,L), Q(LDQ,*), C(LDM,N), R(LDM,M), S(LDS,N), K(LDS,M), H(LDM,M), TOL, WK((N+M)*(N+M+L))
LOGICAL  STQ

3  Description

The Kalman filter arises from the state space model given by:
Xi+1=AiXi+BiWi, VarWi=Qi Yi=CiXi+Vi, VarVi=Ri
where Xi is the state vector of length n at time i, Yi is the observation vector of length m at time i, and Wi of length l and Vi of length m are the independent state noise and measurement noise respectively.
The estimate of Xi given observations Y1 to Yi-1 is denoted by X^ii-1 with state covariance matrix VarX^ii-1=Pii-1=SiSiT, while the estimate of Xi given observations Y1 to Yi is denoted by X^ii with covariance matrix VarX^ii=Pii. The update of the estimate, X^ii-1, from time i to time i+1, is computed in two stages. First, the measurement-update is given by
X^ii=X^ii-1+KiYi-CiX^ii-1 (1)
and
Pii=I-KiCiPii-1 (2)
where Ki=Pii-1 CiT CiPii-1 CiT +Ri -1 is the Kalman gain matrix. The second stage is the time-update for X which is given by
X^i+1i=AiX^ii+DiUi (3)
and
Pi+1i=AiPii AiT +BiQi BiT (4)
where DiUi 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
Ri1/2 CiSi 0 0 AiSi BiQi1/2 U= Hi1/2 0 0 Gi Si+1 0 (5)
where U is an orthogonal transformation triangularizing the left-hand pre-array to produce the right-hand post-array. The relationship between the Kalman gain matrix, Ki, and Gi is given by
AiKi=Gi Hi1/2 -1.
G13EAF requires the input of the lower triangular Cholesky factors of the noise covariance matrices Ri1/2 and, optionally, Qi1/2 and the lower triangular Cholesky factor of the current state covariance matrix, Si, and returns the product of the matrices Ai and Ki, AiKi, the Cholesky factor of the updated state covariance matrix Si+1 and the matrix Hi1/2 used in the computation of the likelihood for the model.

4  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

5  Parameters

1:     N – INTEGERInput
On entry: n, the size of the state vector.
Constraint: N1.
2:     M – INTEGERInput
On entry: m, the size of the observation vector.
Constraint: M1.
3:     L – INTEGERInput
On entry: l, the dimension of the state noise.
Constraint: L1.
4:     A(LDS,N) – REAL (KIND=nag_wp) arrayInput
On entry: the state transition matrix, Ai.
5:     LDS – INTEGERInput
On entry: the first dimension of the arrays A, B, S and K as declared in the (sub)program from which G13EAF is called.
Constraint: LDSN.
6:     B(LDS,L) – REAL (KIND=nag_wp) arrayInput
On entry: the noise coefficient matrix Bi.
7:     STQ – LOGICALInput
On entry: if STQ=.TRUE., the state noise covariance matrix Qi is assumed to be the identity matrix. Otherwise the lower triangular Cholesky factor, Qi1/2, must be provided in Q.
8:     Q(LDQ,*) – REAL (KIND=nag_wp) arrayInput
Note: the second dimension of the array Q must be at least L if STQ=.FALSE. and at least 1 if STQ=.TRUE..
On entry: if STQ=.FALSE., Q must contain the lower triangular Cholesky factor of the state noise covariance matrix, Qi1/2. Otherwise Q is not referenced.
9:     LDQ – INTEGERInput
On entry: the first dimension of the array Q as declared in the (sub)program from which G13EAF is called.
Constraints:
  • if STQ=.FALSE., LDQL;
  • otherwise LDQ1.
10:   C(LDM,N) – REAL (KIND=nag_wp) arrayInput
On entry: the measurement coefficient matrix, Ci.
11:   LDM – INTEGERInput
On entry: the first dimension of the arrays C, R and H as declared in the (sub)program from which G13EAF is called.
Constraint: LDMM.
12:   R(LDM,M) – REAL (KIND=nag_wp) arrayInput
On entry: the lower triangular Cholesky factor of the measurement noise covariance matrix Ri1/2.
13:   S(LDS,N) – REAL (KIND=nag_wp) arrayInput/Output
On entry: the lower triangular Cholesky factor of the state covariance matrix, Si.
On exit: the lower triangular Cholesky factor of the state covariance matrix, Si+1.
14:   K(LDS,M) – REAL (KIND=nag_wp) arrayOutput
On exit: the Kalman gain matrix, Ki, premultiplied by the state transition matrix, Ai, AiKi.
15:   H(LDM,M) – REAL (KIND=nag_wp) arrayOutput
On exit: the lower triangular matrix Hi1/2.
16:   TOL – REAL (KIND=nag_wp)Input
On entry: the tolerance used to test for the singularity of Hi1/2. If 0.0TOL<m2×machine precision, then m2×machine precision is used instead. The inverse of the condition number of H1/2 is estimated by a call to F07TGF (DTRCON). If this estimate is less than TOL then H1/2 is assumed to be singular.
Suggested value: TOL=0.0.
Constraint: TOL0.0.
17:   IWK(M) – INTEGER arrayWorkspace
18:   WK(N+M×N+M+L) – REAL (KIND=nag_wp) arrayWorkspace
19:   IFAIL – INTEGERInput/Output
On entry: IFAIL must be set to 0, -1​ or ​1. If you are unfamiliar with this parameter you should refer to Section 3.3 in the Essential Introduction for details.
For environments where it might be inappropriate to halt program execution when an error is detected, the value -1​ or ​1 is recommended. If the output of error messages is undesirable, then the value 1 is recommended. Otherwise, if you are not familiar with this parameter, the recommended value is 0. When the value -1​ or ​1 is used it is essential to test the value of IFAIL on exit.
On exit: IFAIL=0 unless the routine detects an error or a warning has been flagged (see Section 6).

6  Error Indicators and Warnings

If on entry IFAIL=0 or -1, explanatory error messages are output on the current error message unit (as defined by X04AAF).
Errors or warnings detected by the routine:
IFAIL=1
On entry,N<1,
orM<1,
orL<1,
orLDS<N,
orLDM<M,
orSTQ=.TRUE. and LDQ<1,
orSTQ=.FALSE. and LDQ<L,
orTOL<0.0.
IFAIL=2
The matrix Hi1/2 is singular.

7  Accuracy

The use of the square root algorithm improves the stability of the computations as compared with the direct coding of the Kalman filter. The accuracy will depend on the model.

8  Further Comments

For models with time-invariant A,B and C, G13EBF can be used.
The estimate of the state vector X^i+1i can be computed from X^ii-1 by
X^i+1i=AiX^ii-1+AKiri
where
ri=Yi-CiX^ii- 1
are the independent one step prediction residuals. The required matrix-vector multiplications can be performed by F06PAF (DGEMV).
If Wi and Vi are independent multivariate Normal variates then the log-likelihood for observations i=1,2,,t is given by
lθ = κ - 12 i=1 t l n detHi - 12 i=1 t Yi - Ci X ii-1 T H i -1 Yi - Ci X ii-1
where κ is a constant.
The Cholesky factors of the covariance matrices can be computed using F07FDF (DPOTRF).
Note that the model
Xi+1=AiXi+Wi, VarWi=Qi Yi=CiXi+Vi, VarVi=Ri
can be specified either with B set to the identity matrix and STQ=.FALSE. and the matrix Q1/2 input in Q or with STQ=.TRUE. and B set to Q1/2.
The algorithm requires 76n3+n2 52m+l+n12l2+m2 operations and is backward stable (see Verhaegen and van Dooren (1986)).

9  Example

This example first inputs the number of updates to be computed and the problem sizes. The initial state vector and state covariance matrix are input followed by the model matrices Ai,Bi,Ci,Ri and optionally Qi. The Cholesky factors of the covariance matrices can be computed if required. The model matrices can be input at each update or only once at the first step. At each update the observed values are input and the residuals are computed and printed and the estimate of the state vector, X^ii-1, and the deviance are updated. The deviance is -2×log-likelihood ignoring the constant. After the final update the state covariance matrix is computed from S and printed along with final estimate of the state vector and the value of the deviance.
The data is for a two-dimensional time series to which a VARMA1,1 has been fitted. For the specification of a VARMA model as a state space model see the G13 Chapter Introduction. The initial value of P, P0, is the solution to
P0=A1P0 A1T +B1Q1 B1T .
For convenience, the mean of each series is input before the first update and subtracted from the observations before the measurement update is computed.

9.1  Program Text

Program Text (g13eafe.f90)

9.2  Program Data

Program Data (g13eafe.d)

9.3  Program Results

Program Results (g13eafe.r)


G13EAF (PDF version)
G13 Chapter Contents
G13 Chapter Introduction
NAG Library Manual

© The Numerical Algorithms Group Ltd, Oxford, UK. 2012