SRIleastSquares Class Reference

#include <SRIleastSquares.hpp>

Inheritance diagram for SRIleastSquares:

Inheritance graph
[legend]
Collaboration diagram for SRIleastSquares:

Collaboration graph
[legend]
List of all members.

Detailed Description

class SRIleastSquares inherits SRI and implements a general least squares algorithm using SRI, 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.

Definition at line 73 of file SRIleastSquares.hpp.

Public Member Functions

 SRIleastSquares (void) throw ()
 empty constructor
 SRIleastSquares (const unsigned int N) throw ()
 constructor given the dimension N.
 SRIleastSquares (const Namelist &NL) throw ()
 constructor given a Namelist; its dimension determines the SRI dimension.
 SRIleastSquares (const Matrix< double > &R, const Vector< double > &Z, const Namelist &NL) throw (MatrixException)
 explicit constructor - throw if the dimensions are inconsistent.
 SRIleastSquares (const SRIleastSquares &right) throw ()
 copy constructor
SRIleastSquaresoperator= (const SRIleastSquares &right) throw ()
 operator=
int dataUpdate (Vector< double > &D, Vector< double > &X, Matrix< double > &Cov, void(LSF)(Vector< double > &X, Vector< double > &f, Matrix< double > &P)) throw (MatrixException)
 A general least squares update, NOT the SRIF (Kalman) measurement update.
void zeroAll (void)
 remove all stored information by setting the SRI to zero (does not re-dimension).
bool isValid ()
 Return true if the solution is valid, i.e. if the problem is non-singular.
void Reset (const int N=0) throw (Exception)
 reset the computation, i.e.
Vector< double > Solution (void)
 Get the current solution vector.
int Iterations ()
 Get the number of iterations used in last call to leastSquaresEstimation().
double Convergence ()
 Get the convergence value found in last call to leastSquaresEstimation().
double ConditionNumber ()
 Get the condition number of the covariance matrix from last calls to leastSquaresEstimation() (Larger means 'closer to singular' except zero means condition number is infinite).

Public Attributes

int iterationsLimit
 limit on the number of iterations
double convergenceLimit
 limit on the RSS change in solution which produces success
double divergenceLimit
 upper limit on the RSS change in solution which produces an abort
bool doWeight
 if true, weight the equation using the inverse of covariance matrix on input - default is false
bool doRobust
 if true, weight the equation using robust statistical techniques - default is false
bool doSequential
 if true, save information for a sequential solution - default is false
bool doLinearize
 if true, equation F(X)=D is non-linear, the algorithm will be iterated, and LSF must return partials matrix and F(X).
bool doVerbose
 if true, output intermediate results

Friends

std::ostream & operator<< (std::ostream &s, const SRIleastSquares &srif)
 output operator


Constructor & Destructor Documentation

SRIleastSquares void   )  throw ()
 

empty constructor

Definition at line 61 of file SRIleastSquares.cpp.

SRIleastSquares const unsigned int  N  )  throw ()
 

constructor given the dimension N.

Parameters:
N dimension of the SRIleastSquares.

Definition at line 66 of file SRIleastSquares.cpp.

SRIleastSquares const Namelist NL  )  throw ()
 

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

Parameters:
NL Namelist for the SRIleastSquares.

Definition at line 77 of file SRIleastSquares.cpp.

SRIleastSquares const Matrix< double > &  R,
const Vector< double > &  Z,
const Namelist NL
throw (MatrixException)
 

explicit constructor - throw if the dimensions are inconsistent.

Parameters:
R Initial information matrix, an upper triangular matrix of dim N.
Z Initial information data vector, of length N.
NL Namelist for the SRIleastSquares, also of length N.
Exceptions:
MatrixException if dimensions are not consistent.

Definition at line 89 of file SRIleastSquares.cpp.

References GPSTK_THROW.

SRIleastSquares const SRIleastSquares right  )  throw () [inline]
 

copy constructor

Parameters:
right SRIleastSquares to be copied

Definition at line 100 of file SRIleastSquares.hpp.


Member Function Documentation

double ConditionNumber  )  [inline]
 

Get the condition number of the covariance matrix from last calls to leastSquaresEstimation() (Larger means 'closer to singular' except zero means condition number is infinite).

Definition at line 178 of file SRIleastSquares.hpp.

double Convergence  )  [inline]
 

Get the convergence value found in last call to leastSquaresEstimation().

Returns:
the convergence value

Definition at line 173 of file SRIleastSquares.hpp.

int dataUpdate Vector< double > &  D,
Vector< double > &  X,
Matrix< double > &  Cov,
void(LSF)(Vector< double > &X, Vector< double > &f, Matrix< double > &P) 
throw (MatrixException)
 

A general least squares update, NOT the SRIF (Kalman) measurement update.

Given data and measurement covariance, compute a solution and covariance using the appropriate least squares algorithm.

Parameters:
D Data vector, length M Input: raw data Output: post-fit residuals
X Solution vector, length N Input: nominal solution X0 (zero when doLinearized is false) Output: final solution
Cov Covariance matrix, dimension (N,N) Input: (If doWeight is true) inverse measurement covariance or weight matrix(M,M) Output: Solution covariance matrix (N,N)
LSF Pointer to a function which is used to define the equation to be solved. LSF arguments are: X Nominal solution (input) f Values of the equation f(X) (length M) (output) P Partials matrix df/dX evaluated at X (dimension M,N) (output) When doLinearize is false, LSF ignores X and returns the (constant) partials matrix P and zero for f(X).
Exceptions:
MatrixException if the input is inconsistent Return values: 0 ok -1 Problem is underdetermined (M<N) // TD -- naturalized sol? -2 Problem is singular -3 Algorithm failed to converge -4 Algorithm diverged

Definition at line 271 of file SRIleastSquares.cpp.

References gpstk::StringUtils::asString(), GPSTK_RETHROW, GPSTK_THROW, gpstk::inverse(), Cholesky::L, gpstk::Robust::MedianAbsoluteDeviation(), LabelledVector::message(), gpstk::RMS(), RobustTuningT, Vector::size(), and gpstk::SrifMU().

bool isValid void   )  [inline]
 

Return true if the solution is valid, i.e. if the problem is non-singular.

Definition at line 154 of file SRIleastSquares.hpp.

int Iterations  )  [inline]
 

Get the number of iterations used in last call to leastSquaresEstimation().

Returns:
the number of iterations

Definition at line 169 of file SRIleastSquares.hpp.

SRIleastSquares & operator= const SRIleastSquares right  )  throw ()
 

operator=

Parameters:
right SRIleastSquares to be copied

Definition at line 113 of file SRIleastSquares.cpp.

void Reset const int  N = 0  )  throw (Exception)
 

reset the computation, i.e.

remove all stored information, and optionally change the dimension. If N is not input, the dimension is not changed.

Parameters:
N new SRIleastSquares dimension (optional).

Definition at line 542 of file SRIleastSquares.cpp.

References GPSTK_RETHROW, and Vector::resize().

Vector<double> Solution void   )  [inline]
 

Get the current solution vector.

Returns:
current solution vector

Definition at line 165 of file SRIleastSquares.hpp.

void zeroAll void   ) 
 

remove all stored information by setting the SRI to zero (does not re-dimension).

Definition at line 530 of file SRIleastSquares.cpp.


Friends And Related Function Documentation

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

output operator


Member Data Documentation

double convergenceLimit
 

limit on the RSS change in solution which produces success

Definition at line 185 of file SRIleastSquares.hpp.

double divergenceLimit
 

upper limit on the RSS change in solution which produces an abort

Definition at line 188 of file SRIleastSquares.hpp.

bool doLinearize
 

if true, equation F(X)=D is non-linear, the algorithm will be iterated, and LSF must return partials matrix and F(X).

  • default is false

Definition at line 203 of file SRIleastSquares.hpp.

bool doRobust
 

if true, weight the equation using robust statistical techniques - default is false

Definition at line 196 of file SRIleastSquares.hpp.

bool doSequential
 

if true, save information for a sequential solution - default is false

Definition at line 199 of file SRIleastSquares.hpp.

bool doVerbose
 

if true, output intermediate results

Definition at line 206 of file SRIleastSquares.hpp.

bool doWeight
 

if true, weight the equation using the inverse of covariance matrix on input - default is false

Definition at line 192 of file SRIleastSquares.hpp.

int iterationsLimit
 

limit on the number of iterations

Definition at line 182 of file SRIleastSquares.hpp.


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