USGS

Isis 3.0 Developer's Reference (API)

Home

Isis::LeastSquares Class Reference
[Math]

Generic least square fitting class. More...

#include <LeastSquares.h>

List of all members.

Public Types

enum  SolveMethod { SVD, QRD, SPARSE }

Public Member Functions

 LeastSquares (Isis::BasisFunction &basis, bool sparse=false, int sparseRows=0, int sparseCols=0, bool jigsaw=false)
 Creates a LeastSquares Object.
 ~LeastSquares ()
 Destroys the LeastSquares object.
void AddKnown (const std::vector< double > &input, double expected, double weight=1.0)
 Invoke this method for each set of knowns.
std::vector< double > GetInput (int row) const
 This method returns the data at the given row.
double GetExpected (int row) const
 This method returns the expected value at the given row.
int Rows () const
 This methods returns the number of rows in the matrix.
int Solve (Isis::LeastSquares::SolveMethod method=SVD)
 After all the data has been registered through AddKnown, invoke this method to solve the system of equations.
double Evaluate (const std::vector< double > &input)
 Invokes the BasisFunction Evaluate method.
std::vector< double > Residuals () const
 Returns a vector of residuals (errors).
double Residual (int i) const
 Returns the ith residual.
void Weight (int index, double weight)
 Reset the weight for the ith known.
int Knowns () const
 The number of knowns (or times AddKnown was invoked) linear combination of the variables.
double GetSigma0 ()
int GetDegreesOfFreedom ()
void Reset ()
void ResetSparse ()
bool SparseErrorPropagation ()
 Error propagation for sparse least-squares solution.
std::vector< double > GetEpsilons () const
void SetParameterWeights (const std::vector< double > weights)
void SetNumberOfConstrainedParameters (int n)
const gmm::row_matrix
< gmm::rsvector< double > > & 
GetCovarianceMatrix () const

Detailed Description

Generic least square fitting class.

This class can be used to solved systems of linear equations through least squares fitting. The solution is derived through singular value decomposition or QR decomposition. For example:

\[ x + y = 3 \]

\[ -2x + 3y = 1 \]

\[ 2x - y = 2 \]

Is a simple system of equations that can be solved using this class.

   Isis::BasisFunction b("Linear",2,2);
   Isis::LeastSquares lsq(b);

   std::vector<double> one;
   one.push_back(1.0);
   one.push_back(1.0);
   lsq.AddKnown(one,3.0);

   std::vector<double> two;
   two.push_back(-2.0);
   two.push_back(3.0);
   lsq.AddKnown(two,1.0);

   std::vector<double> tre;
   tre.push_back(2.0);
   tre.push_back(-1.0);
   lsq.AddKnown(tre,2.0);

   lsq.Solve();

   double out1 = lsq.Evaluate(one);
   double out2 = lsq.Evaluate(two);
   double out3 = lsq.Evaluate(tre);
See also:
PolynomialUnivariate
PolynomialBiVariate
BasisFunction
Author:
2004-06-24 Jeff Anderson

Member Enumeration Documentation

Enumerator:
SVD 

Singular Value Decomposition.

QRD 

QR Decomposition.

SPARSE 

Sparse.


Constructor & Destructor Documentation

Isis::LeastSquares::LeastSquares ( Isis::BasisFunction basis,
bool  sparse = false,
int  sparseRows = 0,
int  sparseCols = 0,
bool  jigsaw = false 
)

Creates a LeastSquares Object.

Parameters:
basis A BasisFunction. This parameter allows for the least squares fitting to be applied to arbitrary equations.

References _FILEINFO_, and Isis::IException::Programmer.

Isis::LeastSquares::~LeastSquares (  ) 

Destroys the LeastSquares object.


Member Function Documentation

void Isis::LeastSquares::AddKnown ( const std::vector< double > &  data,
double  result,
double  weight = 1.0 
)

Invoke this method for each set of knowns.

Given our example in the description, we have three knowns and expecteds. They are

\[ (1,1) = 3 \]

\[ (-2,3) = 1 \]

\[ (2,-1) = 2 \]

Parameters:
data A vector of knowns.
result The expected value for the knowns.
weight (Default = 1.0) How strongly to weight this known. Weight less than 1 increases residual for this known, while weight greater than 1 decreases the residual for this known.
Exceptions:
Isis::IException::Programmer - Number of elements in data does not match basis requirements

References _FILEINFO_, Isis::BasisFunction::Name(), Isis::IException::Programmer, and Isis::BasisFunction::Variables().

Referenced by Isis::SurfaceModel::AddTriplet(), Isis::ReseauDistortionMap::SetFocalPlane(), Isis::VimsGroundMap::SetGround(), Isis::SpicePosition::SetPolynomial(), Isis::VimsSkyMap::SetSky(), Isis::ReseauDistortionMap::SetUndistortedFocalPlane(), Isis::OverlapNormalization::Solve(), and Isis::Affine::Solve().

double Isis::LeastSquares::Evaluate ( const std::vector< double > &  data  ) 

Invokes the BasisFunction Evaluate method.

Parameters:
data The input variables to evaluate.
Returns:
The evaluation for the input variable.
Exceptions:
Isis::IException::Programmer - Unable to evaluate until a solution has been computed

References _FILEINFO_, Isis::BasisFunction::Evaluate(), and Isis::IException::Programmer.

Referenced by Isis::SurfaceModel::Evaluate(), Isis::ReseauDistortionMap::SetFocalPlane(), Isis::VimsGroundMap::SetGround(), Isis::VimsSkyMap::SetSky(), and Isis::ReseauDistortionMap::SetUndistortedFocalPlane().

const gmm::row_matrix<gmm::rsvector<double> >& Isis::LeastSquares::GetCovarianceMatrix (  )  const [inline]
int Isis::LeastSquares::GetDegreesOfFreedom (  )  [inline]
std::vector<double> Isis::LeastSquares::GetEpsilons (  )  const [inline]
double Isis::LeastSquares::GetExpected ( int  row  )  const

This method returns the expected value at the given row.

Parameters:
row 
Returns:
double

References _FILEINFO_, Isis::IException::Programmer, Rows(), and Isis::toString().

std::vector< double > Isis::LeastSquares::GetInput ( int  row  )  const

This method returns the data at the given row.

Parameters:
row 
Returns:
std::vector<double>

References _FILEINFO_, Isis::IException::Programmer, Rows(), and Isis::toString().

double Isis::LeastSquares::GetSigma0 (  )  [inline]
int Isis::LeastSquares::Knowns (  )  const [inline]

The number of knowns (or times AddKnown was invoked) linear combination of the variables.

Returns:
int The number of knowns

Referenced by Isis::VimsGroundMap::SetGround(), Isis::VimsSkyMap::SetSky(), and Isis::BundleAdjust::Solve().

void Isis::LeastSquares::Reset (  ) 
void Isis::LeastSquares::ResetSparse (  )  [inline]

References Reset().

double Isis::LeastSquares::Residual ( int  i  )  const

Returns the ith residual.

That is, the difference between the evaluation of a known with the solution against the expected value. There is one residual for each time AddKnown was invoked.

Parameters:
i The number of times AddKnown was invoked to be evaluated.
Returns:
The output value of the residual.
Exceptions:
Isis::IException::Programmer - Unable to return residuals until a solution has been computed

References _FILEINFO_, and Isis::IException::Programmer.

std::vector< double > Isis::LeastSquares::Residuals (  )  const

Returns a vector of residuals (errors).

That is, the difference between the evaluation of a known with the solution against the expected value.

Returns:
The vector of residuals.
Exceptions:
Isis::IException::Programmer - Unable to return residuals until a solution has been computed

References _FILEINFO_, and Isis::IException::Programmer.

Referenced by Isis::BundleAdjust::Solve().

int Isis::LeastSquares::Rows (  )  const

This methods returns the number of rows in the matrix.

Returns:
int

Referenced by GetExpected(), GetInput(), and Solve().

void Isis::LeastSquares::SetNumberOfConstrainedParameters ( int  n  )  [inline]
void Isis::LeastSquares::SetParameterWeights ( const std::vector< double >  weights  )  [inline]
int Isis::LeastSquares::Solve ( Isis::LeastSquares::SolveMethod  method = SVD  ) 
bool Isis::LeastSquares::SparseErrorPropagation (  ) 

Error propagation for sparse least-squares solution.

Computes the variance-covariance matrix of the parameters. This is the inverse of the normal equations matrix, scaled by the reference variance (also called variance of unit weight, etc).

Referenced by Isis::BundleAdjust::Solve().

void Isis::LeastSquares::Weight ( int  index,
double  weight 
)

Reset the weight for the ith known.

This weight will not be used unless the system is resolved using the Solve method.

Parameters:
index The position in the array to assign the given weight value
weight A weight factor to apply to the ith known. A weight less than one increase the residual for this known while a weight greater than one decrease the residual for this known.

The documentation for this class was generated from the following files: