Isis 3 Programmer Reference
Isis::LeastSquares Class Reference

Generic least square fitting class. More...

#include <LeastSquares.h>

Collaboration diagram for Isis::LeastSquares:
Collaboration graph

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 ()
 
std::vector< double > GetEpsilons () const
 
void SetParameterWeights (const std::vector< double > weights)
 
void SetNumberOfConstrainedParameters (int n)
 

Private Member Functions

void SolveSVD ()
 After all the data has been registered through AddKnown, invoke this method to solve the system of equations.
 
void SolveQRD ()
 After all the data has been registered through AddKnown, invoke this method to solve the system of equations with a QR decomposition of A = QR.
 
void SolveCholesky ()
 
int SolveSparse ()
 Solve using sparse class.
 
void FillSparseA (const std::vector< double > &data)
 Invoke this method for each set of knowns for sparse solutions.
 
bool ApplyParameterWeights ()
 

Private Attributes

arma::mat p_xSparse
 sparse solution matrix
 
std::vector< double > p_epsilonsSparse
 sparse vector of total parameter corrections
 
std::vector< double > p_parameterWeights
 vector of parameter weights
 
arma::SpMat< double > p_sparseA
 design matrix 'A'
 
arma::SpMat< double > p_normals
 normal equations matrix 'N'
 
arma::mat p_ATb
 right-hand side vector
 
arma::mat p_SLU_Factor
 decomposed normal equations matrix
 
bool p_jigsaw
 
bool p_sparse
 
bool p_solved
 Boolean value indicating solution is complete.
 
int p_currentFillRow
 
int p_sparseRows
 
int p_sparseCols
 
int p_constrainedParameters
 constrained parameters
 
int p_degreesOfFreedom
 degrees of freedom (redundancy)
 
double p_sigma0
 sigma nought - reference variance
 
std::vector< std::vector< double > > p_input
 A vector of the input variables to evaluate.
 
std::vector< double > p_expected
 A vector of the expected values when solved.
 
std::vector< double > p_sqrtWeight
 A vector of the square roots of the weights for each known value.
 
std::vector< double > p_residuals
 A vector of the residuals (or difference between expected and solved values).
 
Isis::BasisFunctionp_basis
 Pointer to the BasisFunction object.
 

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);
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);
Generic linear equation class.
Generic least square fitting class.
See also
PolynomialUnivariate
PolynomialBiVariate
BasisFunction
Author
2004-06-24 Jeff Anderson
History

2004-06-24 Jeff Anderson, Original version

2007-11-16 Tracie Sucharski and Debbie A. Cook Added SolveQRD method for a faster solve

2008-02-06 Tracie Sucharski, Renamed Solve to SolveSVD and make private. Add public Solve, with a new parameter, method which defaults to SVD. This calls the correct solution method. This was done for documentation purposes-clarifies the default solves by SVD.

2008-04-16 Debbie Cook / Tracie Sucharski, Added SolveSparse.

2008-06-09 Tracie Sucharski, Added conditional compilations for Solaris. We could not get SuperLu to build under Solaris due to a confilict with Complex definitions.

2009-04-06 Tracie Sucharski, Added return value to SolveSparse, which indicates a column containing all zeroes which causes superlu to bomb.

2009-12-21 Jeannie Walldren, Changed p_weight to p_sqrtweight. Modified Weight() and AddKnown() to take the square root of the given weight and add it to the p_sqrtweight vector.

2010-04-20 Debbie A. Cook, Replaced SparseReset with Reset to reset all solution methods

2010-11-22 Debbie A. Cook - Merged with Ken Edmundson version

2013-12-29 Jeannie Backer - Improved error messages. Fixes #962.

2019-09-05 Makayla Shepherd & Jesse Mapel - Changed sparse solution to use Armadillo library's SuperLU interface instead of GMM.

Definition at line 99 of file LeastSquares.h.

Member Enumeration Documentation

◆ SolveMethod

Enumerator
SVD 

Singular Value Decomposition.

QRD 

QR Decomposition.

SPARSE 

Sparse.

Definition at line 112 of file LeastSquares.h.

Constructor & Destructor Documentation

◆ LeastSquares()

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

Creates a LeastSquares Object.

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

Definition at line 23 of file LeastSquares.cpp.

References p_ATb, p_basis, p_epsilonsSparse, p_normals, p_parameterWeights, p_sigma0, p_solved, p_sparseA, p_xSparse, and Isis::IException::Programmer.

◆ ~LeastSquares()

Isis::LeastSquares::~LeastSquares ( )

Destroys the LeastSquares object.

Definition at line 62 of file LeastSquares.cpp.

Member Function Documentation

◆ AddKnown()

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
dataA vector of knowns.
resultThe 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
History

2008-04-22 Tracie Sucharski, Fill sparse matrix.

2009-12-21 Jeannie Walldren, Modified code to add the square root of the weight to the vector p_sqrtweight.

Definition at line 96 of file LeastSquares.cpp.

References FillSparseA(), Isis::BasisFunction::Name(), p_basis, p_expected, p_input, p_sqrtWeight, Isis::IException::Programmer, and Isis::BasisFunction::Variables().

Referenced by Isis::SurfaceModel::AddTriplet(), and Isis::OverlapNormalization::Solve().

◆ Evaluate()

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

Invokes the BasisFunction Evaluate method.

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

Definition at line 553 of file LeastSquares.cpp.

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

Referenced by Isis::SurfaceModel::Evaluate().

◆ FillSparseA()

void Isis::LeastSquares::FillSparseA ( const std::vector< double > & data)
private

Invoke this method for each set of knowns for sparse solutions.

The A sparse matrix must be filled as we go or we will quickly run out of memory for large solutions. So, expand the basis function, apply weights (which is done in the Solve method for the non-sparse case.

Parameters
dataA vector of knowns.
History

2008-04-22 Tracie Sucharski - New method for sparse solutions.

2009-12-21 Jeannie Walldren - Changed variable name from p_weight to p_sqrtweight.

Definition at line 137 of file LeastSquares.cpp.

References Isis::BasisFunction::Expand(), p_basis, p_sparseA, p_sqrtWeight, and Isis::BasisFunction::Term().

Referenced by AddKnown().

◆ GetDegreesOfFreedom()

int Isis::LeastSquares::GetDegreesOfFreedom ( )
inline

Definition at line 134 of file LeastSquares.h.

◆ GetEpsilons()

std::vector< double > Isis::LeastSquares::GetEpsilons ( ) const
inline

Definition at line 138 of file LeastSquares.h.

◆ GetExpected()

double Isis::LeastSquares::GetExpected ( int row) const

This method returns the expected value at the given row.

Parameters
row
Returns
double

Definition at line 173 of file LeastSquares.cpp.

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

◆ GetInput()

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

This method returns the data at the given row.

Parameters
row
Returns
std::vector<double>

Definition at line 158 of file LeastSquares.cpp.

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

◆ GetSigma0()

double Isis::LeastSquares::GetSigma0 ( )
inline

Definition at line 133 of file LeastSquares.h.

◆ Knowns()

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

Definition at line 129 of file LeastSquares.h.

References p_expected.

◆ Reset()

void Isis::LeastSquares::Reset ( )

Definition at line 523 of file LeastSquares.cpp.

◆ ResetSparse()

void Isis::LeastSquares::ResetSparse ( )
inline

Definition at line 137 of file LeastSquares.h.

◆ Residual()

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
iThe 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

Definition at line 590 of file LeastSquares.cpp.

References p_residuals, p_solved, and Isis::IException::Programmer.

◆ Residuals()

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

Definition at line 570 of file LeastSquares.cpp.

References p_residuals, p_solved, and Isis::IException::Programmer.

◆ Rows()

int Isis::LeastSquares::Rows ( ) const

This methods returns the number of rows in the matrix.

Returns
int

Definition at line 187 of file LeastSquares.cpp.

References p_input.

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

◆ SetNumberOfConstrainedParameters()

void Isis::LeastSquares::SetNumberOfConstrainedParameters ( int n)
inline

Definition at line 140 of file LeastSquares.h.

◆ SetParameterWeights()

void Isis::LeastSquares::SetParameterWeights ( const std::vector< double > weights)
inline

Definition at line 139 of file LeastSquares.h.

◆ Solve()

int Isis::LeastSquares::Solve ( Isis::LeastSquares::SolveMethod method = SVD)

After all the data has been registered through AddKnown, invoke this method to solve the system of equations.

You can then use the Evaluate and Residual methods freely.

History

2008-04-16 Debbie Cook / Tracie Sucharski, Added SolveSparse.

2009-04-08 Tracie Sucharski - Added return value which will pass on what is returned from SolveSparse which is a column number of a column that contained all zeros.

2010-12-12 Debbie A. Cook Fixed "no data" test for SPARSE case

Definition at line 205 of file LeastSquares.cpp.

References p_solved, QRD, Rows(), SolveQRD(), SolveSparse(), SolveSVD(), SPARSE, SVD, and Isis::IException::Unknown.

Referenced by Isis::SurfaceModel::Solve(), and Isis::OverlapNormalization::Solve().

◆ SolveCholesky()

void Isis::LeastSquares::SolveCholesky ( )
inlineprivate

Definition at line 145 of file LeastSquares.h.

◆ SolveQRD()

void Isis::LeastSquares::SolveQRD ( )
private

After all the data has been registered through AddKnown, invoke this method to solve the system of equations with a QR decomposition of A = QR.

You can then use the Evaluate and Residual methods freely. The QR decomposition is only slightly less reliable than the SVD, but much faster.

History
2009-12-21 Jeannie Walldren - Changed variable name from p_weight to p_sqrtweight.

Definition at line 339 of file LeastSquares.cpp.

References Isis::BasisFunction::Coefficients(), Isis::BasisFunction::Evaluate(), Isis::BasisFunction::Expand(), p_basis, p_expected, p_input, p_residuals, p_solved, p_sqrtWeight, Isis::BasisFunction::SetCoefficients(), Isis::BasisFunction::Term(), and Isis::IException::Unknown.

Referenced by Solve().

◆ SolveSparse()

int Isis::LeastSquares::SolveSparse ( )
private

Solve using sparse class.

After all the data has been registered through AddKnown, invoke this method to solve the system of equations Nx = b, where

N = ATPA b = ATPl, and

N is the "normal equations matrix; A is the so-called "design" matrix; P is the observation weight matrix (typically diagonal); l is the "computed - observed" column vector;

The solution is achieved using a sparse matrix formulation of the LU decomposition of the normal equations.

You can then use the Evaluate and Residual methods freely.

History

2008-04-16 Debbie Cook / Tracie Sucharski, New method

2008-04-23 Tracie Sucharski, Fill sparse matrix as we go in AddKnown method rather than in the solve method, otherwise we run out of memory very quickly.

2009-04-08 Tracie Sucharski - Added return value which is a column number of a column that contained all zeros.

2009-12-21 Jeannie Walldren - Changed variable name from p_weight to p_sqrtweight.

2010 Ken Edmundson

2010-11-20 Debbie A. Cook Merged Ken Edmundson verion with system version

2011-03-17 Ken Edmundson Corrected computation of residuals

Definition at line 428 of file LeastSquares.cpp.

References p_ATb, p_basis, p_constrainedParameters, p_degreesOfFreedom, p_epsilonsSparse, p_expected, p_normals, p_parameterWeights, p_residuals, p_sigma0, p_solved, p_sparseA, p_sqrtWeight, p_xSparse, Isis::BasisFunction::SetCoefficients(), and Isis::IException::Unknown.

Referenced by Solve().

◆ SolveSVD()

void Isis::LeastSquares::SolveSVD ( )
private

After all the data has been registered through AddKnown, invoke this method to solve the system of equations.

You can then use the Evaluate and Residual methods freely.

History
2009-12-21 Jeannie Walldren - Changed variable name from p_weight to p_sqrtweight.

Definition at line 236 of file LeastSquares.cpp.

References Isis::BasisFunction::Coefficients(), Isis::BasisFunction::Evaluate(), Isis::BasisFunction::Expand(), p_basis, p_degreesOfFreedom, p_expected, p_input, p_residuals, p_sigma0, p_solved, p_sqrtWeight, Isis::BasisFunction::SetCoefficients(), Isis::BasisFunction::Term(), Isis::toString(), and Isis::IException::Unknown.

Referenced by Solve().

◆ Weight()

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
indexThe position in the array to assign the given weight value
weightA 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.
History
2009-12-21 Jeannie Walldren, Modified code to add the square root of the weight to the vector p_sqrtweight.

Definition at line 614 of file LeastSquares.cpp.

References p_sqrtWeight.

Member Data Documentation

◆ p_ATb

arma::mat Isis::LeastSquares::p_ATb
private

right-hand side vector

Definition at line 157 of file LeastSquares.h.

Referenced by LeastSquares(), and SolveSparse().

◆ p_basis

Isis::BasisFunction* Isis::LeastSquares::p_basis
private

Pointer to the BasisFunction object.

Definition at line 183 of file LeastSquares.h.

Referenced by AddKnown(), Evaluate(), FillSparseA(), LeastSquares(), SolveQRD(), SolveSparse(), and SolveSVD().

◆ p_constrainedParameters

int Isis::LeastSquares::p_constrainedParameters
private

constrained parameters

Definition at line 167 of file LeastSquares.h.

Referenced by SolveSparse().

◆ p_currentFillRow

int Isis::LeastSquares::p_currentFillRow
private

Definition at line 164 of file LeastSquares.h.

◆ p_degreesOfFreedom

int Isis::LeastSquares::p_degreesOfFreedom
private

degrees of freedom (redundancy)

Definition at line 168 of file LeastSquares.h.

Referenced by SolveSparse(), and SolveSVD().

◆ p_epsilonsSparse

std::vector<double> Isis::LeastSquares::p_epsilonsSparse
private

sparse vector of total parameter corrections

Definition at line 152 of file LeastSquares.h.

Referenced by LeastSquares(), and SolveSparse().

◆ p_expected

std::vector<double> Isis::LeastSquares::p_expected
private

A vector of the expected values when solved.

Definition at line 174 of file LeastSquares.h.

Referenced by AddKnown(), GetExpected(), Knowns(), SolveQRD(), SolveSparse(), and SolveSVD().

◆ p_input

std::vector<std::vector<double> > Isis::LeastSquares::p_input
private

A vector of the input variables to evaluate.

Definition at line 172 of file LeastSquares.h.

Referenced by AddKnown(), GetInput(), Rows(), SolveQRD(), and SolveSVD().

◆ p_jigsaw

bool Isis::LeastSquares::p_jigsaw
private

Definition at line 160 of file LeastSquares.h.

◆ p_normals

arma::SpMat<double> Isis::LeastSquares::p_normals
private

normal equations matrix 'N'

Definition at line 156 of file LeastSquares.h.

Referenced by LeastSquares(), and SolveSparse().

◆ p_parameterWeights

std::vector<double> Isis::LeastSquares::p_parameterWeights
private

vector of parameter weights

Definition at line 153 of file LeastSquares.h.

Referenced by LeastSquares(), and SolveSparse().

◆ p_residuals

std::vector<double> Isis::LeastSquares::p_residuals
private

A vector of the residuals (or difference between expected and solved values).

Definition at line 179 of file LeastSquares.h.

Referenced by Residual(), Residuals(), SolveQRD(), SolveSparse(), and SolveSVD().

◆ p_sigma0

double Isis::LeastSquares::p_sigma0
private

sigma nought - reference variance

Definition at line 170 of file LeastSquares.h.

Referenced by LeastSquares(), SolveSparse(), and SolveSVD().

◆ p_SLU_Factor

arma::mat Isis::LeastSquares::p_SLU_Factor
private

decomposed normal equations matrix

Definition at line 158 of file LeastSquares.h.

◆ p_solved

bool Isis::LeastSquares::p_solved
private

Boolean value indicating solution is complete.

Definition at line 162 of file LeastSquares.h.

Referenced by Evaluate(), LeastSquares(), Residual(), Residuals(), Solve(), SolveQRD(), SolveSparse(), and SolveSVD().

◆ p_sparse

bool Isis::LeastSquares::p_sparse
private

Definition at line 161 of file LeastSquares.h.

◆ p_sparseA

arma::SpMat<double> Isis::LeastSquares::p_sparseA
private

design matrix 'A'

Definition at line 155 of file LeastSquares.h.

Referenced by FillSparseA(), LeastSquares(), and SolveSparse().

◆ p_sparseCols

int Isis::LeastSquares::p_sparseCols
private

Definition at line 166 of file LeastSquares.h.

◆ p_sparseRows

int Isis::LeastSquares::p_sparseRows
private

Definition at line 165 of file LeastSquares.h.

◆ p_sqrtWeight

std::vector<double> Isis::LeastSquares::p_sqrtWeight
private

A vector of the square roots of the weights for each known value.

Definition at line 176 of file LeastSquares.h.

Referenced by AddKnown(), FillSparseA(), SolveQRD(), SolveSparse(), SolveSVD(), and Weight().

◆ p_xSparse

arma::mat Isis::LeastSquares::p_xSparse
private

sparse solution matrix

Definition at line 151 of file LeastSquares.h.

Referenced by LeastSquares(), and SolveSparse().


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