Isis 3 Programmer Reference
|
Generic least square fitting class. More...
#include <LeastSquares.h>
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. More... | |
~LeastSquares () | |
Destroys the LeastSquares object. More... | |
void | AddKnown (const std::vector< double > &input, double expected, double weight=1.0) |
Invoke this method for each set of knowns. More... | |
std::vector< double > | GetInput (int row) const |
This method returns the data at the given row. More... | |
double | GetExpected (int row) const |
This method returns the expected value at the given row. More... | |
int | Rows () const |
This methods returns the number of rows in the matrix. More... | |
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. More... | |
double | Evaluate (const std::vector< double > &input) |
Invokes the BasisFunction Evaluate method. More... | |
std::vector< double > | Residuals () const |
Returns a vector of residuals (errors). More... | |
double | Residual (int i) const |
Returns the ith residual. More... | |
void | Weight (int index, double weight) |
Reset the weight for the ith known. More... | |
int | Knowns () const |
The number of knowns (or times AddKnown was invoked) linear combination of the variables. More... | |
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. More... | |
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. More... | |
void | SolveCholesky () |
int | SolveSparse () |
Solve using sparse class. More... | |
void | FillSparseA (const std::vector< double > &data) |
Invoke this method for each set of knowns for sparse solutions. More... | |
bool | ApplyParameterWeights () |
Private Attributes | |
arma::mat | p_xSparse |
sparse solution matrix More... | |
std::vector< double > | p_epsilonsSparse |
sparse vector of total parameter corrections More... | |
std::vector< double > | p_parameterWeights |
vector of parameter weights More... | |
arma::SpMat< double > | p_sparseA |
design matrix 'A' More... | |
arma::SpMat< double > | p_normals |
normal equations matrix 'N' More... | |
arma::mat | p_ATb |
right-hand side vector More... | |
arma::mat | p_SLU_Factor |
decomposed normal equations matrix More... | |
bool | p_jigsaw |
bool | p_sparse |
bool | p_solved |
Boolean value indicating solution is complete. More... | |
int | p_currentFillRow |
int | p_sparseRows |
int | p_sparseCols |
int | p_constrainedParameters |
constrained parameters More... | |
int | p_degreesOfFreedom |
degrees of freedom (redundancy) More... | |
double | p_sigma0 |
sigma nought - reference variance More... | |
std::vector< std::vector< double > > | p_input |
A vector of the input variables to evaluate. More... | |
std::vector< double > | p_expected |
A vector of the expected values when solved. More... | |
std::vector< double > | p_sqrtWeight |
A vector of the square roots of the weights for each known value. More... | |
std::vector< double > | p_residuals |
A vector of the residuals (or difference between expected and solved values). More... | |
Isis::BasisFunction * | p_basis |
Pointer to the BasisFunction object. More... | |
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.
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 115 of file LeastSquares.h.
Enumerator | |
---|---|
SVD | Singular Value Decomposition. |
QRD | QR Decomposition. |
SPARSE | Sparse. |
Definition at line 128 of file LeastSquares.h.
Isis::LeastSquares::LeastSquares | ( | Isis::BasisFunction & | basis, |
bool | sparse = false , |
||
int | sparseRows = 0 , |
||
int | sparseCols = 0 , |
||
bool | jigsaw = false |
||
) |
Creates a LeastSquares Object.
basis | A BasisFunction. This parameter allows for the least squares fitting to be applied to arbitrary equations. |
Definition at line 39 of file LeastSquares.cpp.
References _FILEINFO_, p_ATb, p_basis, p_epsilonsSparse, p_normals, p_parameterWeights, p_sigma0, p_solved, p_sparseA, p_xSparse, and Isis::IException::Programmer.
Isis::LeastSquares::~LeastSquares | ( | ) |
Destroys the LeastSquares object.
Definition at line 78 of file LeastSquares.cpp.
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 \]
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. |
Isis::IException::Programmer | - Number of elements in data does not match basis requirements |
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 112 of file LeastSquares.cpp.
References _FILEINFO_, FillSparseA(), Isis::BasisFunction::Name(), p_basis, p_expected, p_input, p_sqrtWeight, Isis::IException::Programmer, and Isis::BasisFunction::Variables().
Referenced by Isis::SurfaceModel::AddTriplet(), Isis::ReseauDistortionMap::SetFocalPlane(), Isis::VimsGroundMap::SetGround(), Isis::SpicePosition::SetPolynomial(), Isis::SpiceRotation::SetPolynomial(), Isis::VimsSkyMap::SetSky(), Isis::ReseauDistortionMap::SetUndistortedFocalPlane(), and Isis::Affine::Solve().
double Isis::LeastSquares::Evaluate | ( | const std::vector< double > & | data | ) |
Invokes the BasisFunction Evaluate method.
data | The input variables to evaluate. |
Isis::IException::Programmer | - Unable to evaluate until a solution has been computed |
Definition at line 569 of file LeastSquares.cpp.
References _FILEINFO_, Isis::BasisFunction::Evaluate(), p_basis, p_solved, and Isis::IException::Programmer.
Referenced by Isis::SurfaceModel::Evaluate(), Isis::ReseauDistortionMap::SetFocalPlane(), Isis::VimsGroundMap::SetGround(), Isis::VimsSkyMap::SetSky(), and Isis::ReseauDistortionMap::SetUndistortedFocalPlane().
|
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.
data | A vector of knowns. |
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 153 of file LeastSquares.cpp.
References Isis::BasisFunction::Expand(), p_basis, p_sparseA, p_sqrtWeight, and Isis::BasisFunction::Term().
Referenced by AddKnown().
double Isis::LeastSquares::GetExpected | ( | int | row | ) | const |
This method returns the expected value at the given row.
row |
Definition at line 189 of file LeastSquares.cpp.
References _FILEINFO_, p_expected, 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.
row |
Definition at line 174 of file LeastSquares.cpp.
References _FILEINFO_, p_input, Isis::IException::Programmer, Rows(), and Isis::toString().
|
inline |
The number of knowns (or times AddKnown was invoked) linear combination of the variables.
Definition at line 145 of file LeastSquares.h.
References p_expected.
Referenced by Isis::VimsGroundMap::SetGround(), and Isis::VimsSkyMap::SetSky().
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.
i | The number of times AddKnown was invoked to be evaluated. |
Isis::IException::Programmer | - Unable to return residuals until a solution has been computed |
Definition at line 606 of file LeastSquares.cpp.
References _FILEINFO_, p_residuals, p_solved, 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.
Isis::IException::Programmer | - Unable to return residuals until a solution has been computed |
Definition at line 586 of file LeastSquares.cpp.
References _FILEINFO_, p_residuals, p_solved, and Isis::IException::Programmer.
int Isis::LeastSquares::Rows | ( | ) | const |
This methods returns the number of rows in the matrix.
Definition at line 203 of file LeastSquares.cpp.
References p_input.
Referenced by GetExpected(), GetInput(), and 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.
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 221 of file LeastSquares.cpp.
References _FILEINFO_, p_solved, QRD, Rows(), SolveQRD(), SolveSparse(), SolveSVD(), SPARSE, SVD, and Isis::IException::Unknown.
Referenced by Isis::ReseauDistortionMap::SetFocalPlane(), Isis::VimsGroundMap::SetGround(), Isis::SpicePosition::SetPolynomial(), Isis::SpiceRotation::SetPolynomial(), Isis::VimsSkyMap::SetSky(), Isis::ReseauDistortionMap::SetUndistortedFocalPlane(), Isis::SurfaceModel::Solve(), and Isis::Affine::Solve().
|
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.
Definition at line 355 of file LeastSquares.cpp.
References _FILEINFO_, 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().
|
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.
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 444 of file LeastSquares.cpp.
References _FILEINFO_, 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().
|
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.
Definition at line 252 of file LeastSquares.cpp.
References _FILEINFO_, 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().
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.
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. |
Definition at line 630 of file LeastSquares.cpp.
References p_sqrtWeight.
|
private |
right-hand side vector
Definition at line 173 of file LeastSquares.h.
Referenced by LeastSquares(), and SolveSparse().
|
private |
Pointer to the BasisFunction object.
Definition at line 199 of file LeastSquares.h.
Referenced by AddKnown(), Evaluate(), FillSparseA(), LeastSquares(), SolveQRD(), SolveSparse(), and SolveSVD().
|
private |
|
private |
degrees of freedom (redundancy)
Definition at line 184 of file LeastSquares.h.
Referenced by SolveSparse(), and SolveSVD().
|
private |
sparse vector of total parameter corrections
Definition at line 168 of file LeastSquares.h.
Referenced by LeastSquares(), and SolveSparse().
|
private |
A vector of the expected values when solved.
Definition at line 190 of file LeastSquares.h.
Referenced by AddKnown(), GetExpected(), Knowns(), SolveQRD(), SolveSparse(), and SolveSVD().
|
private |
A vector of the input variables to evaluate.
Definition at line 188 of file LeastSquares.h.
Referenced by AddKnown(), GetInput(), Rows(), SolveQRD(), and SolveSVD().
|
private |
normal equations matrix 'N'
Definition at line 172 of file LeastSquares.h.
Referenced by LeastSquares(), and SolveSparse().
|
private |
vector of parameter weights
Definition at line 169 of file LeastSquares.h.
Referenced by LeastSquares(), and SolveSparse().
|
private |
A vector of the residuals (or difference between expected and solved values).
Definition at line 195 of file LeastSquares.h.
Referenced by Residual(), Residuals(), SolveQRD(), SolveSparse(), and SolveSVD().
|
private |
sigma nought - reference variance
Definition at line 186 of file LeastSquares.h.
Referenced by LeastSquares(), SolveSparse(), and SolveSVD().
|
private |
decomposed normal equations matrix
Definition at line 174 of file LeastSquares.h.
|
private |
Boolean value indicating solution is complete.
Definition at line 178 of file LeastSquares.h.
Referenced by Evaluate(), LeastSquares(), Residual(), Residuals(), Solve(), SolveQRD(), SolveSparse(), and SolveSVD().
|
private |
design matrix 'A'
Definition at line 171 of file LeastSquares.h.
Referenced by FillSparseA(), LeastSquares(), and SolveSparse().
|
private |
A vector of the square roots of the weights for each known value.
Definition at line 192 of file LeastSquares.h.
Referenced by AddKnown(), FillSparseA(), SolveQRD(), SolveSparse(), SolveSVD(), and Weight().
|
private |
sparse solution matrix
Definition at line 167 of file LeastSquares.h.
Referenced by LeastSquares(), and SolveSparse().