#include <SRIFilter.hpp>
Inheritance diagram for SRIFilter:


SRIFilter may be used for Kalman filtering, smoothing, or for simple least squares, including weighted, linear or linearized, robust and/or sequential algorithms.
At any point the state X and covariance P are related to the SRI by X = inverse(R) * z , P = inverse(R) * inverse(transpose(R)), or R = upper triangular square root (Cholesky decomposition) of the inverse of P, and z = R * X.
The SRIFilter implements Kalman filter algorithm, which includes sequential least squares (measurement update), dynamic propagation (time update), and smoothing (technically the term 'Kalman filter algorithm' is reserved for the classical algorithm just as Kalman presented it, in terms of a state vector and its covariance matrix).
The SRIFilter measurment update (which is actually just linear least squares) is half of the SRIFilter (Kalman filter) - there is a 'time update' that propagates the SRI (and thus the state and covariance) forward in time using the dynamical model of the filter. These are algebraically equivalent to the classical Kalman algorithm, but are more efficient and numerically stable (actually the Kalman algorithm has been shown to be numerically unstable!). The SRIFilter smoothing algorithms consists of a 'backwards' filter, implemented by applying a 'smoother update' to the SRIFilter at each point in reverse order.
Ref: Bierman, G.J. "Factorization Methods for Discrete Sequential Estimation," Academic Press, 1977.
Definition at line 93 of file SRIFilter.hpp.
Public Member Functions | |
| SRIFilter (void) throw () | |
| empty constructor | |
| SRIFilter (const unsigned int N) throw () | |
| constructor given the dimension N. | |
| SRIFilter (const Namelist &NL) throw () | |
| constructor given a Namelist; its dimension determines the SRI dimension. | |
| SRIFilter (const Matrix< double > &R, const Vector< double > &Z, const Namelist &NL) throw (MatrixException) | |
| explicit constructor - throw if the dimensions are inconsistent. | |
| SRIFilter (const SRIFilter &right) throw () | |
| copy constructor | |
| SRIFilter & | operator= (const SRIFilter &right) throw () |
| operator= | |
| void | measurementUpdate (const Matrix< double > &H, Vector< double > &D, const Matrix< double > &CM=SRINullMatrix) throw (MatrixException,VectorException) |
| SRIF (Kalman) simple linear measurement update with optional weight matrix. | |
| void | timeUpdate (Matrix< double > &PhiInv, Matrix< double > &Rw, Matrix< double > &G, Vector< double > &zw, Matrix< double > &Rwx) throw (MatrixException) |
| SRIF (Kalman) time update This routine uses the Householder transformation to propagate the SRIFilter state and covariance through a time step. | |
| void | smootherUpdate (Matrix< double > &Phi, Matrix< double > &Rw, Matrix< double > &G, Vector< double > &zw, Matrix< double > &Rwx) throw (MatrixException) |
| SRIF (Kalman) smoother update This routine uses the Householder transformation to propagate the SRIF state and covariance through a smoother (backward filter) step. | |
| void | zeroAll (void) |
| remove all stored information by setting the SRI to zero (does not re-dimension). | |
| void | Reset (const int N=0) |
| reset the computation, i.e. | |
Static Public Member Functions | |
| void | DMsmootherUpdate (Matrix< double > &P, Vector< double > &X, Matrix< double > &Phinv, Matrix< double > &Rw, Matrix< double > &G, Vector< double > &Zw, Matrix< double > &Rwx) throw (MatrixException) |
| Covariance/State version of the Kalman smoother update (Dyer-McReynolds). | |
| void | DMsmootherUpdateWithControl (Matrix< double > &P, Vector< double > &X, Matrix< double > &Phinv, Matrix< double > &Rw, Matrix< double > &G, Vector< double > &Zw, Matrix< double > &Rwx, Vector< double > &U) throw (MatrixException) |
| Modification for case with control vector: Xj+1 = Phi*Xj + Gwj + u. | |
|
|
empty constructor
Definition at line 65 of file SRIFilter.cpp. |
|
|
constructor given the dimension N.
Definition at line 73 of file SRIFilter.cpp. |
|
|
constructor given a Namelist; its dimension determines the SRI dimension.
Definition at line 84 of file SRIFilter.cpp. |
|
||||||||||||||||
|
explicit constructor - throw if the dimensions are inconsistent.
Definition at line 96 of file SRIFilter.cpp. References GPSTK_THROW. |
|
|
copy constructor
Definition at line 120 of file SRIFilter.hpp. |
|
||||||||||||||||||||||||||||||||
|
Covariance/State version of the Kalman smoother update (Dyer-McReynolds). This routine implements the Dyer-McReynolds form of the state and covariance recursions which constitute the backward filter of the Square Root Information Smoother; it is equivalent to the SRI form implemented in SRIFilter::smootherUpdate(). NB. This routine does NOT use the SRIFilter object; it is implemented as a member function to be consistent with other updates.
Ref: Bierman, G.J. "Factorization Methods for Discrete Sequential Estimation," Academic Press, 1977. Definition at line 206 of file SRIFilter.cpp. References GPSTK_RETHROW. |
|
||||||||||||||||||||||||||||||||||||
|
Modification for case with control vector: Xj+1 = Phi*Xj + Gwj + u.
|
|
||||||||||||||||
|
SRIF (Kalman) simple linear measurement update with optional weight matrix.
Definition at line 133 of file SRIFilter.cpp. References GPSTK_RETHROW, GPSTK_THROW, gpstk::inverse(), Cholesky::L, gpstk::SrifMU(), and gpstk::SRINullMatrix. |
|
|
operator=
Definition at line 120 of file SRIFilter.cpp. |
|
|
reset the computation, i.e. remove all stored information, and optionally change the dimension. If N is not input, the dimension is not changed.
Definition at line 231 of file SRIFilter.cpp. References Vector::resize(), Matrix::resize(), and Matrix::rows(). |
|
||||||||||||||||||||||||
|
SRIF (Kalman) smoother update This routine uses the Householder transformation to propagate the SRIF state and covariance through a smoother (backward filter) step. If the existing SRI state is of dimension N, and the number of noise parameter is Ns, then the inputs must be dimensioned as indicated.
In this implementation of the backward filter, the Householder transformation is applied to the following matrix [dimensions are shown in ()]: _ (Ns) (N) (1) _ _ _ (Ns) | Rw+Rwx*G Rwx*Phi zw | ==> | Rw Rwx zw | (N) | R*G R*Phi z | ==> | 0 R z | .
Ref: Bierman, G.J. "Factorization Methods for Discrete Sequential Estimation," Academic Press, 1977, pg 216. Definition at line 194 of file SRIFilter.cpp. References GPSTK_RETHROW. |
|
||||||||||||||||||||||||
|
SRIF (Kalman) time update This routine uses the Householder transformation to propagate the SRIFilter state and covariance through a time step. If the existing SRI state is of dimension n, and the number of noise parameter is ns, then the inputs must be dimensioned as indicated.
The matrix Rwx is related to the sensitivity of the state estimate to the unmodeled parameters in Zw. The sensitivity matrix is Sen = -inverse(Rw)*Rwx, where perturbation in model X = Sen * diagonal(a priori sigmas of parameter uncertainties). The quantities Rw, Rwx and Zw on output are to be saved and used in the sqrt information fixed interval smoother (SRIS), during the backward filter process. ------------------------------------------------------------------- Ref: Bierman, G.J. "Factorization Methods for Discrete Sequential Estimation," Academic Press, 1977, pg 121. Definition at line 181 of file SRIFilter.cpp. References GPSTK_RETHROW. |
|
|
remove all stored information by setting the SRI to zero (does not re-dimension).
Definition at line 221 of file SRIFilter.cpp. |
1.3.9.1