SRI Class Reference

#include <SRI.hpp>

Inheritance diagram for SRI:

Inheritance graph
[legend]
Collaboration diagram for SRI:

Collaboration graph
[legend]
List of all members.

Detailed Description

class SRI encapsulates all the information associated with the solution of a set of simultaneous linear equations.

It is used in least squares estimation (linear and linearized) and is the basis of the preferred implementation of Kalman filtering. An SRI consists of just three things: (1) 'R', the 'information matrix', which is an upper triangular matrix of dimension N, equal to the inverse of the square root (or Cholesky decomposition) of the solution covariance matrix, (2) 'Z', the 'SRI state vector' of length N (parallels the components of R), (not to be confused with the regular state vector X), and (3) 'names', a Namelist used to label the elements of R and Z (parallels and labels rows and columns of R and elements of Z). A Namelist is part of class SRI because the manipulations of SRI (see functions below) requires a consistent way of manipulating the different individual elements of R and Z, in addition it allows the user to attach 'human-readable' labels to the elements of the state vector, which is useful in adding, dropping and bumping states, and it makes printed results more readable (see the LabelledMatrix class in Namelist.hpp).

The set of simultaneous equations represented by an SRI is R * X = Z, where X is the (unknown) state vector (the conventional solution vector) also of dimension N. The state X is solved for as X = inverse(R) * Z, and the covariance matrix of the state X is equal to transpose(inverse(R))*inverse(R).

Least squares estimation via SRI is very simple and efficient; it uses the Householder transformation to convert the problem to upper triangular form, and then uses very efficient algorithms to invert the information matrix to find the solution and its covariance. The usual matrix equation is H * X = D, where H is the 'design matrix' or the 'partials matrix', of dimension M x N, X is the (unknown) solution vector of length N, and D is the 'data' or 'measurement' vector of length M. In the least squares 'update' of the SRI, this set of information {H,D} is concatenated with the existing SRI {R,Z} to form an (N+M x N+1) matrix Q which has R in the upper left, Z upper right, H lower left and D lower right. This extended matrix is then subjected to a Householder transformation (see class Matrix), which will put (at least the first N columns of) Q into upper triangular form. The result is a new, updated SRI (R and Z) in the place of the old, while in place of D are residuals of fit corresponding to the measurements in D (the H part of Q is trashed). This result, in fact (see the reference), produces an updated SRI which gives precisely the usual least squares solution for the combined 'a priori SRI + new data' problem. This algorithm is called a 'measurement update' of the SRI.

It is most enlightening to think of the SRI and this process in terms of 'information'. The SRI contains all the 'information' which has come from updates that have been made to it using (H,D) pairs. Initially, the SRI is all zeros, which corresponds to 'no information'. This overcomes one serious problem with conventional least squares and the Kalman algorithm, namely that a 'zero information' starting value cannot be correctly expressed, because in that case the covariance matrix is singular and the state vector is indeterminate; in the SRI method this is perfectly consistent - the covariance matrix is singular because the information matrix (R) is zero, and thus the state is entirely indeterminate. As new 'information' (in the form of data D and partials matrix H pairs) is added to the SRI (via the Householder algorithm), the 'information' stored in R and Z is increased and they become non-zero. (By the way note that the number of rows in the {H,D} information is arbitrary - information can be added in 'batches' - M large - or one - M=1 - piece at a time.) When there is enough information, R becomes non-singular, and so can be inverted and the solution and covariance can be computed. As the amount of information becomes large, elements of R become large, and thus elements of the covariance (think of covariance as a measure of uncertainty - the opposite or inverse of information) become small.

The structure of the SRI method allows some powerful techniques to be used in manipulating, combining and separating state elements and the information associated with them in SRIs. For example, if the measurement updates have failed to increase the information about one particular state element, then that element, and its information, may be removed from the problem by deleting that element's row and column of R, and its element of Z (and then re-triangularizing the SRI). In general, any subset of an SRI may be separated, or the SRI split (see the routine of that name below - note the caveats) into two separate SRIs. For another example, SRI allows the information of a each state element to be selectively reduced or even zeroed, simply by multiplying the corresponding elements of R and Z by a factor; in Kalman filtering this is called a 'Q bump' of the element and is very important in some filtering applications. There are methods (see below) consistently to merge (operator+()), split, and permute elements of, SRIs.

Kalman filtering is an important application of SRI methods (actually it is called 'square root information filtering' or SRIF - 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 measurment update described above (which is actually just linear least squares) is half of the SRIF (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!). There are even SRI smoothing algorithms, corresponding to Kalman smoothers, which consist of a 'backwards' filter, implemented by applying a 'smoother update' to the SRI at each point in reverse order.

Ref: Bierman, G.J. "Factorization Methods for Discrete Sequential Estimation," Academic Press, 1977.

Definition at line 166 of file SRI.hpp.

Public Member Functions

 SRI (void) throw ()
 empty constructor
 throw ()
 constructor given the dimension N.
 throw ()
 constructor given a Namelist, its dimension determines the SRI dimension.
 SRI (const Matrix< double > &, const Vector< double > &, const Namelist &) throw (MatrixException)
 explicit constructor - throw if the dimensions are inconsistent.
 throw ()
 copy constructor
SRIoperator= (const SRI &right) throw ()
 operator=
void permute (const Namelist &) throw (MatrixException,VectorException)
 Permute the SRI elements to match the input Namelist, which may differ with the SRI Namelist by AT MOST A PERMUTATION; throw if this is not true.
void split (const Namelist &, SRI &) throw (MatrixException,VectorException)
 split an SRI into two others, this one matching the input Namelist, the other containing whatever is left.
SRIoperator+= (const Namelist &) throw (MatrixException,VectorException)
 extend this SRI to include the given Namelist, with no added information; names in the input namelist which are not unique are ignored.
void reshape (const Namelist &) throw (MatrixException,VectorException)
 reshape this SRI to match the input Namelist, by calling other member functions, including split(), operator+() and permute()
void merge (const SRI &S) throw (MatrixException,VectorException)
 merge an SRI into this one. NB names may be reordered in the result.
SRIoperator+= (const SRI &) throw (MatrixException,VectorException)
 merge this SRI with the given input SRI.
SRIappend (const SRI &S) throw (MatrixException,VectorException)
 append an SRI onto this SRI.
void zeroOne (const unsigned int n) throw ()
 Zero out the nth row of R and the nth element of Z, removing all information about that element.
void zeroAll (const int n=0) throw ()
 Zero out all the first n rows of R and elements of Z, removing all information about those elements.
void zeroState (void) throw ()
 Zero out (set all elements to zero) the state (Vector Z) only.
void shift (const Vector< double > &X0) throw (MatrixException)
 Shift the state vector by a constant vector X0; does not change information i.e.
void shiftZ (const Vector< double > &Z0) throw (MatrixException)
 Shift the SRI state vector (Z) by a constant vector Z0; does not change information.
void transform (const Matrix< double > &T, const Matrix< double > &invT=SRINullMatrix) throw (MatrixException,VectorException)
 Transform this SRI with the transformation matrix T; i.e.
void transformState (const Matrix< double > &invT) throw (MatrixException)
 Transform the state by the transformation matrix T; i.e.
void Qbump (const unsigned int &in, const double &q=0.0) throw (MatrixException,VectorException)
 Decrease the information in this SRI for, or 'Q bump', the element with the input index.
void stateFix (const unsigned int &index, const double &value) throw (MatrixException,VectorException)
 Fix the state element with the input index to the input value, and collapse the SRI by removing that element.
void stateFix (const Namelist &drops, const Vector< double > &values) throw (MatrixException,VectorException)
 Vector version of stateFix, with Namelist identifying the states.
void addAPriori (const Matrix< double > &Cov, const Vector< double > &X) throw (MatrixException)
 Add a priori or constraint information in the form of an ordinary state vector and covariance matrix.
void addAPrioriInformation (const Matrix< double > &ICov, const Vector< double > &X) throw (MatrixException)
 Add a priori or constraint information in the form of an information matrix (inverse covariance) and ordinary state.
void measurementUpdate (Matrix< double > &Partials, Vector< double > &Data) throw (MatrixException)
 SRIF (Kalman) measurement update, or least squares update Call the SRI measurement update for this SRI and the given input.
void getConditionNumber (double &small, double &big) throw (MatrixException)
 Compute the condition number, or rather the largest and smallest eigenvalues of the SRI matrix R (the condition number is the ratio of the largest and smallest eigenvalues).
void getState (Vector< double > &X, int *ptrSingularIndex=NULL) throw (MatrixException)
 Compute the state X without computing the covariance matrix C.
void getStateAndCovariance (Vector< double > &X, Matrix< double > &C, double *ptrSmall=NULL, double *ptrBig=NULL) throw (MatrixException,VectorException)
 Compute the state X and the covariance matrix C of the state, where C = transpose(inverse(R))*inverse(R) and X = inverse(R) * Z.
unsigned int size (void) const throw ()
 return the size of the SRI, which is the dimension of R(rows and columns), Z and names.
Namelist getNames (void) throw ()
 access the Namelist of the SRI
std::string getName (const unsigned int in) throw ()
 access the name of a specific state element, given its index.
bool setName (const unsigned int in, const std::string &name) throw ()
 assign the name of a specific state element, given its index; no effect, and return false, if the name is not unique; return true if successful.
unsigned int index (std::string &name) throw ()
 return the index of the name in the Namelist that matches the input, or -1 if not found.
Matrix< double > getR (void) throw ()
 access the R matrix
Vector< double > getZ (void) throw ()
 access the Z vector

Protected Attributes

Matrix< double > R
 Information matrix, an upper triangular (square) matrix.
Vector< double > Z
 SRI state vector, of length equal to the dimension (row and col) of R.
Namelist names
 Namelist parallel to R and Z, labelling the elements of the state vector.

Friends

SRI operator+ (const SRI &, const SRI &) throw (MatrixException,VectorException)
 merge two SRIs to produce a third. ? should this be operator&() ?
std::ostream & operator<< (std::ostream &s, const SRI &)
 output operator


Constructor & Destructor Documentation

SRI void   )  throw () [inline]
 

empty constructor

Definition at line 169 of file SRI.hpp.

SRI const Matrix< double > &  ,
const Vector< double > &  ,
const Namelist
throw (MatrixException)
 

explicit constructor - throw if the dimensions are inconsistent.

Definition at line 95 of file SRI.cpp.

References GPSTK_THROW.


Member Function Documentation

void addAPriori const Matrix< double > &  Cov,
const Vector< double > &  X
throw (MatrixException)
 

Add a priori or constraint information in the form of an ordinary state vector and covariance matrix.

The matrix must be non-singular.

Parameters:
Cov Covariance matrix of same dimension as this SRIFilter
X State vector of same dimension as this SRIFilter
Exceptions:
if input is invalid: dimensions are wrong or Cov is singular.

Definition at line 831 of file SRI.cpp.

References GPSTK_THROW, and gpstk::inverse().

void addAPrioriInformation const Matrix< double > &  ICov,
const Vector< double > &  X
throw (MatrixException)
 

Add a priori or constraint information in the form of an information matrix (inverse covariance) and ordinary state.

ICov must be non-singular.

Parameters:
ICov Inverse covariance matrix of same dimension as this SRIFilter
X State vector of same dimension as this SRIFilter
Exceptions:
if input is invalid: dimensions are wrong.

Definition at line 854 of file SRI.cpp.

References GPSTK_THROW, Cholesky::L, gpstk::SrifMU(), and gpstk::transpose().

SRI & append const SRI S  )  throw (MatrixException,VectorException)
 

append an SRI onto this SRI.

Similar to opertor+= but simpler; input SRI is simply appended, first using operator+=(Namelist), then filling the new portions of R and Z, all without final Householder transformation of result. Do not allow a name that is already present to be added: throw.

Definition at line 358 of file SRI.cpp.

References GPSTK_RETHROW, and GPSTK_THROW.

void getConditionNumber double &  small,
double &  big
throw (MatrixException)
 

Compute the condition number, or rather the largest and smallest eigenvalues of the SRI matrix R (the condition number is the ratio of the largest and smallest eigenvalues).

Note that the condition number of the covariance matrix would be the square of the condition number of R.

Definition at line 882 of file SRI.cpp.

References GPSTK_RETHROW, SVD::S, and SVD::sort().

std::string getName const unsigned int  in  )  throw () [inline]
 

access the name of a specific state element, given its index.

returns 'out-of-range' if the index is out of range.

Definition at line 394 of file SRI.hpp.

References Namelist::getName().

Namelist getNames void   )  throw () [inline]
 

access the Namelist of the SRI

Definition at line 388 of file SRI.hpp.

Matrix<double> getR void   )  throw () [inline]
 

access the R matrix

Definition at line 413 of file SRI.hpp.

void getState Vector< double > &  X,
int *  ptrSingularIndex = NULL
throw (MatrixException)
 

Compute the state X without computing the covariance matrix C.

R*X=Z so X=inverse(R)*Z; in this routine the state is computed explicitly, without forming the inverse of R, using the fact that R is upper triangular. NB. The matrix is singular if and only if one or more of the diagonal elements is zero; in the case the routine will still have valid entries in the state vector for index greater than the largest index with zero diagonal

Parameters:
X State vector (output)
ptrSingularIndex if ptr is non-null, on output *ptr will be the largest index of singularity
Exceptions:
SingularMatrixException if R is singular.

Definition at line 905 of file SRI.cpp.

References GPSTK_THROW, Vector::size(), and gpstk::sum().

void getStateAndCovariance Vector< double > &  X,
Matrix< double > &  C,
double *  ptrSmall = NULL,
double *  ptrBig = NULL
throw (MatrixException,VectorException)
 

Compute the state X and the covariance matrix C of the state, where C = transpose(inverse(R))*inverse(R) and X = inverse(R) * Z.

Optional pointer arguments will return smallest and largest eigenvalues of the R matrix, which is a measure of singularity.

Parameters:
X State vector (output)
C Covariance of the state vector (output)
ptrSmall Pointer to double, on output *ptrSmall set to smallest eigenvalue of R
ptrBig Pointer to double, on output *ptrBig set to largest eigenvalue of R
Exceptions:
SingularMatrixException if R is singular. NB this is the most efficient way to invert the SRI equation.

Definition at line 931 of file SRI.cpp.

References GPSTK_RETHROW, gpstk::inverseUT(), and gpstk::UTtimesTranspose().

Vector<double> getZ void   )  throw () [inline]
 

access the Z vector

Definition at line 418 of file SRI.hpp.

unsigned int index std::string &  name  )  throw () [inline]
 

return the index of the name in the Namelist that matches the input, or -1 if not found.

Definition at line 408 of file SRI.hpp.

References Namelist::index().

void measurementUpdate Matrix< double > &  Partials,
Vector< double > &  Data
throw (MatrixException) [inline]
 

SRIF (Kalman) measurement update, or least squares update Call the SRI measurement update for this SRI and the given input.

See doc. for SrifMU().

Definition at line 337 of file SRI.hpp.

References gpstk::SrifMU().

void merge const SRI S  )  throw (MatrixException,VectorException) [inline]
 

merge an SRI into this one. NB names may be reordered in the result.

Definition at line 221 of file SRI.hpp.

SRI & operator+= const SRI  )  throw (MatrixException,VectorException)
 

merge this SRI with the given input SRI.

? should this be operator&=() ? NB may reorder the names in the resulting Namelist.

Definition at line 398 of file SRI.cpp.

References Householder::A, ConstMatrixBase::colCopy(), GPSTK_RETHROW, GPSTK_THROW, Namelist::index(), Namelist::labels, and Namelist::size().

SRI & operator+= const Namelist  )  throw (MatrixException,VectorException)
 

extend this SRI to include the given Namelist, with no added information; names in the input namelist which are not unique are ignored.

Definition at line 293 of file SRI.cpp.

References GPSTK_RETHROW, SRI::R, and SRI::Z.

SRI & operator= const SRI right  )  throw ()
 

operator=

Definition at line 127 of file SRI.cpp.

void permute const Namelist  )  throw (MatrixException,VectorException)
 

Permute the SRI elements to match the input Namelist, which may differ with the SRI Namelist by AT MOST A PERMUTATION; throw if this is not true.

Definition at line 141 of file SRI.cpp.

References GPSTK_RETHROW, GPSTK_THROW, gpstk::identical(), gpstk::SrifMU(), and gpstk::transpose().

void Qbump const unsigned int &  in,
const double &  q = 0.0
throw (MatrixException,VectorException)
 

Decrease the information in this SRI for, or 'Q bump', the element with the input index.

This means that the uncertainty and the state element given by the index are divided by the input factor q; the default input is zero, which means zero out the information (q = infinite). A Q bump by factor q is equivalent to 'de-weighting' the element by q. No effect if in is out of range.

Definition at line 644 of file SRI.cpp.

References Householder::A, ConstMatrixBase::colCopy(), and GPSTK_RETHROW.

void reshape const Namelist  )  throw (MatrixException,VectorException)
 

reshape this SRI to match the input Namelist, by calling other member functions, including split(), operator+() and permute()

Definition at line 327 of file SRI.cpp.

References GPSTK_RETHROW, and gpstk::StringUtils::split().

bool setName const unsigned int  in,
const std::string &  name
throw () [inline]
 

assign the name of a specific state element, given its index; no effect, and return false, if the name is not unique; return true if successful.

Definition at line 401 of file SRI.hpp.

References Namelist::setName().

void shift const Vector< double > &  X0  )  throw (MatrixException)
 

Shift the state vector by a constant vector X0; does not change information i.e.

let R * X = Z => R' * (X-X0) = Z' throw on invalid input dimension

Definition at line 523 of file SRI.cpp.

References GPSTK_THROW.

void shiftZ const Vector< double > &  Z0  )  throw (MatrixException)
 

Shift the SRI state vector (Z) by a constant vector Z0; does not change information.

i.e. let Z => Z-Z0 throw on invalid input dimension

Definition at line 540 of file SRI.cpp.

References GPSTK_THROW.

unsigned int size void   )  const throw () [inline]
 

return the size of the SRI, which is the dimension of R(rows and columns), Z and names.

Definition at line 383 of file SRI.hpp.

void split const Namelist ,
SRI
throw (MatrixException,VectorException)
 

split an SRI into two others, this one matching the input Namelist, the other containing whatever is left.

The input Namelist must be a non-trivial subset of this->names; throw MatrixException if it is not. NB. Interpreting the results of a split() and merge (operator+()) operations should be done very carefully; remember that the SRI contains both solution and noise, and that the results of these operations are not always as expected, particularly note that split() and operator+() are usually NOT reversible.

Definition at line 239 of file SRI.cpp.

References GPSTK_RETHROW, GPSTK_THROW, Namelist::labels, SRI::R, Vector::resize(), Vector::size(), Namelist::swap(), and SRI::Z.

void stateFix const Namelist drops,
const Vector< double > &  values
throw (MatrixException,VectorException)
 

Vector version of stateFix, with Namelist identifying the states.

Fix the given state elements to the input value, and collapse the SRI by removing those elements. No effect if name is not found.

Definition at line 722 of file SRI.cpp.

References GPSTK_RETHROW, GPSTK_THROW, and Vector::size().

void stateFix const unsigned int &  index,
const double &  value
throw (MatrixException,VectorException)
 

Fix the state element with the input index to the input value, and collapse the SRI by removing that element.

No effect if index is out of range.

Definition at line 685 of file SRI.cpp.

References GPSTK_RETHROW.

throw  ) 
 

copy constructor

throw  ) 
 

constructor given a Namelist, its dimension determines the SRI dimension.

throw  ) 
 

constructor given the dimension N.

void transform const Matrix< double > &  T,
const Matrix< double > &  invT = SRINullMatrix
throw (MatrixException,VectorException)
 

Transform this SRI with the transformation matrix T; i.e.

R -> T * R * inverse(T) and Z -> T * Z. The matrix inverse(T) may optionally be supplied as input, otherwise it is computed from T. NB names in this SRI are most likely changed; but this routine does not change the Namelist. Throw MatrixException if the input has the wrong dimension or cannot be inverted.

Definition at line 560 of file SRI.cpp.

References GPSTK_RETHROW, GPSTK_THROW, gpstk::inverseSVD(), gpstk::SrifMU(), and gpstk::SRINullMatrix.

void transformState const Matrix< double > &  invT  )  throw (MatrixException)
 

Transform the state by the transformation matrix T; i.e.

X -> T*X, without transforming the SRI; this is done by right multiplying R by inverse(T), which is the input. Thus R -> R*inverse(T), so R*inverse(T)*T*X = Z. Input is the _inverse_ of the transformation. throw MatrixException if input dimensions are wrong.

Definition at line 609 of file SRI.cpp.

References Householder::A, and GPSTK_THROW.

void zeroAll const int  n = 0  )  throw ()
 

Zero out all the first n rows of R and elements of Z, removing all information about those elements.

Default value of the input is 0, meaning zero out the entire SRI.

Definition at line 500 of file SRI.cpp.

void zeroOne const unsigned int  n  )  throw ()
 

Zero out the nth row of R and the nth element of Z, removing all information about that element.

Definition at line 483 of file SRI.cpp.

void zeroState void   )  throw () [inline]
 

Zero out (set all elements to zero) the state (Vector Z) only.

Definition at line 254 of file SRI.hpp.


Friends And Related Function Documentation

SRI operator+ const SRI Sleft,
const SRI Sright
throw (MatrixException,VectorException) [friend]
 

merge two SRIs to produce a third. ? should this be operator&() ?

Definition at line 463 of file SRI.cpp.

std::ostream& operator<< std::ostream &  s,
const SRI
[friend]
 

output operator


Member Data Documentation

Namelist names [protected]
 

Namelist parallel to R and Z, labelling the elements of the state vector.

Definition at line 435 of file SRI.hpp.

Referenced by gpstk::operator<<().

Matrix<double> R [protected]
 

Information matrix, an upper triangular (square) matrix.

Definition at line 429 of file SRI.hpp.

Referenced by SRI::operator+=(), gpstk::operator<<(), and SRI::split().

Vector<double> Z [protected]
 

SRI state vector, of length equal to the dimension (row and col) of R.

Definition at line 432 of file SRI.hpp.

Referenced by SRI::operator+=(), gpstk::operator<<(), and SRI::split().


The documentation for this class was generated from the following files:
Generated on Wed Feb 8 03:31:28 2012 for GPS ToolKit Software Library by  doxygen 1.3.9.1