Isis 3 Programmer Reference
LeastSquares.h
1 #ifndef LeastSquares_h
2 #define LeastSquares_h
3 
8 /* SPDX-License-Identifier: CC0-1.0 */
9 #include <vector>
10 
11 #include "tnt/tnt_array2d.h"
12 
13 #include <armadillo>
14 
15 #include "BasisFunction.h"
16 
17 namespace Isis {
99  class LeastSquares {
100  public:
101 
102  LeastSquares(Isis::BasisFunction &basis, bool sparse = false,
103  int sparseRows = 0, int sparseCols = 0, bool jigsaw = false);
104  ~LeastSquares();
105  void AddKnown(const std::vector<double> &input, double expected,
106  double weight = 1.0);
107 
108  std::vector<double> GetInput(int row) const;
109  double GetExpected(int row) const;
110  int Rows() const;
111 
112  enum SolveMethod { SVD,
113  QRD,
114  SPARSE
115  };
116 
118  double Evaluate(const std::vector<double> &input);
119  std::vector<double> Residuals() const;
120  double Residual(int i) const;
121  void Weight(int index, double weight);
122 
129  int Knowns() const {
130  return p_expected.size();
131  };
132 
133  double GetSigma0() { return p_sigma0; }
134  int GetDegreesOfFreedom() { return p_degreesOfFreedom; }
135  void Reset ();
136 
137  void ResetSparse() { Reset(); }
138  std::vector<double> GetEpsilons () const { return p_epsilonsSparse; }
139  void SetParameterWeights(const std::vector<double> weights) { p_parameterWeights = weights; }
140  void SetNumberOfConstrainedParameters(int n) { p_constrainedParameters = n; }
141 
142  private:
143  void SolveSVD();
144  void SolveQRD();
145  void SolveCholesky () {}
146 
147  int SolveSparse();
148  void FillSparseA(const std::vector<double> &data);
149  bool ApplyParameterWeights();
150 
151  arma::mat p_xSparse;
152  std::vector<double> p_epsilonsSparse;
153  std::vector<double> p_parameterWeights;
155  arma::SpMat<double> p_sparseA;
156  arma::SpMat<double> p_normals;
157  arma::mat p_ATb;
158  arma::mat p_SLU_Factor;
160  bool p_jigsaw;
161  bool p_sparse;
162  bool p_solved;
164  int p_currentFillRow;
165  int p_sparseRows;
166  int p_sparseCols;
170  double p_sigma0;
172  std::vector<std::vector<double> > p_input;
174  std::vector<double> p_expected;
176  std::vector<double> p_sqrtWeight;
179  std::vector<double> p_residuals;
184 
185 
186 
187  };
188 };
189 
190 
191 #endif
Isis::LeastSquares::LeastSquares
LeastSquares(Isis::BasisFunction &basis, bool sparse=false, int sparseRows=0, int sparseCols=0, bool jigsaw=false)
Creates a LeastSquares Object.
Definition: LeastSquares.cpp:23
Isis::LeastSquares::p_epsilonsSparse
std::vector< double > p_epsilonsSparse
sparse vector of total parameter corrections
Definition: LeastSquares.h:152
Isis::LeastSquares::Solve
int Solve(Isis::LeastSquares::SolveMethod method=SVD)
After all the data has been registered through AddKnown, invoke this method to solve the system of eq...
Definition: LeastSquares.cpp:205
Isis::LeastSquares::p_xSparse
arma::mat p_xSparse
sparse solution matrix
Definition: LeastSquares.h:151
Isis::LeastSquares::SolveSVD
void SolveSVD()
After all the data has been registered through AddKnown, invoke this method to solve the system of eq...
Definition: LeastSquares.cpp:236
Isis::LeastSquares::p_constrainedParameters
int p_constrainedParameters
constrained parameters
Definition: LeastSquares.h:167
Isis::LeastSquares::GetInput
std::vector< double > GetInput(int row) const
This method returns the data at the given row.
Definition: LeastSquares.cpp:158
Isis::LeastSquares::Knowns
int Knowns() const
The number of knowns (or times AddKnown was invoked) linear combination of the variables.
Definition: LeastSquares.h:129
Isis::BasisFunction
Generic linear equation class.
Definition: BasisFunction.h:48
Isis::LeastSquares::FillSparseA
void FillSparseA(const std::vector< double > &data)
Invoke this method for each set of knowns for sparse solutions.
Definition: LeastSquares.cpp:137
Isis::LeastSquares::p_parameterWeights
std::vector< double > p_parameterWeights
vector of parameter weights
Definition: LeastSquares.h:153
Isis::LeastSquares::p_expected
std::vector< double > p_expected
A vector of the expected values when solved.
Definition: LeastSquares.h:174
Isis::LeastSquares::SPARSE
@ SPARSE
Sparse.
Definition: LeastSquares.h:114
Isis::LeastSquares::Evaluate
double Evaluate(const std::vector< double > &input)
Invokes the BasisFunction Evaluate method.
Definition: LeastSquares.cpp:553
Isis::LeastSquares::SolveMethod
SolveMethod
Definition: LeastSquares.h:112
Isis::LeastSquares::Rows
int Rows() const
This methods returns the number of rows in the matrix.
Definition: LeastSquares.cpp:187
Isis::LeastSquares::Residual
double Residual(int i) const
Returns the ith residual.
Definition: LeastSquares.cpp:590
Isis::LeastSquares::p_input
std::vector< std::vector< double > > p_input
A vector of the input variables to evaluate.
Definition: LeastSquares.h:172
Isis::LeastSquares::~LeastSquares
~LeastSquares()
Destroys the LeastSquares object.
Definition: LeastSquares.cpp:62
Isis::LeastSquares
Generic least square fitting class.
Definition: LeastSquares.h:99
Isis::LeastSquares::p_sparseA
arma::SpMat< double > p_sparseA
design matrix 'A'
Definition: LeastSquares.h:155
Isis::LeastSquares::p_normals
arma::SpMat< double > p_normals
normal equations matrix 'N'
Definition: LeastSquares.h:156
Isis::LeastSquares::p_basis
Isis::BasisFunction * p_basis
Pointer to the BasisFunction object.
Definition: LeastSquares.h:183
Isis::LeastSquares::GetExpected
double GetExpected(int row) const
This method returns the expected value at the given row.
Definition: LeastSquares.cpp:173
Isis::LeastSquares::p_residuals
std::vector< double > p_residuals
A vector of the residuals (or difference between expected and solved values).
Definition: LeastSquares.h:179
Isis::LeastSquares::p_sqrtWeight
std::vector< double > p_sqrtWeight
A vector of the square roots of the weights for each known value.
Definition: LeastSquares.h:176
Isis::LeastSquares::p_ATb
arma::mat p_ATb
right-hand side vector
Definition: LeastSquares.h:157
Isis::LeastSquares::QRD
@ QRD
QR Decomposition.
Definition: LeastSquares.h:113
Isis::LeastSquares::p_SLU_Factor
arma::mat p_SLU_Factor
decomposed normal equations matrix
Definition: LeastSquares.h:158
Isis::LeastSquares::p_sigma0
double p_sigma0
sigma nought - reference variance
Definition: LeastSquares.h:170
Isis::LeastSquares::p_solved
bool p_solved
Boolean value indicating solution is complete.
Definition: LeastSquares.h:162
Isis::LeastSquares::SolveSparse
int SolveSparse()
Solve using sparse class.
Definition: LeastSquares.cpp:428
Isis::LeastSquares::AddKnown
void AddKnown(const std::vector< double > &input, double expected, double weight=1.0)
Invoke this method for each set of knowns.
Definition: LeastSquares.cpp:96
Isis::LeastSquares::SVD
@ SVD
Singular Value Decomposition.
Definition: LeastSquares.h:112
Isis::LeastSquares::p_degreesOfFreedom
int p_degreesOfFreedom
degrees of freedom (redundancy)
Definition: LeastSquares.h:168
Isis
This is free and unencumbered software released into the public domain.
Definition: Apollo.h:16
Isis::LeastSquares::Weight
void Weight(int index, double weight)
Reset the weight for the ith known.
Definition: LeastSquares.cpp:614
Isis::LeastSquares::SolveQRD
void SolveQRD()
After all the data has been registered through AddKnown, invoke this method to solve the system of eq...
Definition: LeastSquares.cpp:339
Isis::LeastSquares::Residuals
std::vector< double > Residuals() const
Returns a vector of residuals (errors).
Definition: LeastSquares.cpp:570