Isis 3.0 Programmer Reference
Back | Home
LeastSquares.h
Go to the documentation of this file.
1 #ifndef LeastSquares_h
2 #define LeastSquares_h
3 
25 #include <vector>
26 
27 #include "tnt/tnt_array2d.h"
28 
29 
30 #ifndef __sun__
31 // used to ignore warnings generated by gmm.h when building on clang
32 #include "gmm/gmm.h"
33 #endif
34 
35 #include "BasisFunction.h"
36 
37 namespace Isis {
117  class LeastSquares {
118  public:
119 
120  LeastSquares(Isis::BasisFunction &basis, bool sparse = false,
121  int sparseRows = 0, int sparseCols = 0, bool jigsaw = false);
122  ~LeastSquares();
123  void AddKnown(const std::vector<double> &input, double expected,
124  double weight = 1.0);
125 
126  std::vector<double> GetInput(int row) const;
127  double GetExpected(int row) const;
128  int Rows() const;
129 
130  enum SolveMethod { SVD,
131  QRD,
133  };
134 
136  double Evaluate(const std::vector<double> &input);
137  std::vector<double> Residuals() const;
138  double Residual(int i) const;
139  void Weight(int index, double weight);
140 
147  int Knowns() const {
148  return p_expected.size();
149  };
150 
151  double GetSigma0() { return p_sigma0; }
152  int GetDegreesOfFreedom() { return p_degreesOfFreedom; }
153  void Reset ();
154 
155 #ifndef __sun__
156  void ResetSparse() { Reset(); }
157  bool SparseErrorPropagation();
158  std::vector<double> GetEpsilons () const { return p_epsilonsSparse; }
159  void SetParameterWeights(const std::vector<double> weights) { p_parameterWeights = weights; }
160  void SetNumberOfConstrainedParameters(int n) { p_constrainedParameters = n; }
161  const gmm::row_matrix<gmm::rsvector<double> >& GetCovarianceMatrix () const { return p_normals; }
162 #endif
163 
164  private:
165  void SolveSVD();
166  void SolveQRD();
167  void SolveCholesky () {}
168 
169 #ifndef __sun__
170  int SolveSparse();
171  void FillSparseA(const std::vector<double> &data);
172  bool ApplyParameterWeights();
173 
174  std::vector<double> p_xSparse;
175  std::vector<double> p_epsilonsSparse;
176  std::vector<double> p_parameterWeights;
178  gmm::row_matrix<gmm::rsvector<double> > p_sparseA;
179  private:
180  gmm::row_matrix<gmm::rsvector<double> > p_normals;
181  private:
182  gmm::dense_matrix<double> p_ATb;
183  gmm::SuperLU_factor<double> p_SLU_Factor;
184 #endif
185  bool p_jigsaw;
186  bool p_sparse;
187  bool p_solved;
189  int p_currentFillRow;
190  int p_sparseRows;
191  int p_sparseCols;
195  double p_sigma0;
197  std::vector<std::vector<double> > p_input;
199  std::vector<double> p_expected;
201  std::vector<double> p_sqrtWeight;
204  std::vector<double> p_residuals;
209 
210 
211 
212  };
213 };
214 
215 
216 #endif
void Weight(int index, double weight)
Reset the weight for the ith known.
QR Decomposition.
Definition: LeastSquares.h:131
gmm::row_matrix< gmm::rsvector< double > > p_sparseA
design matrix &#39;A&#39;
Definition: LeastSquares.h:178
void SolveSVD()
After all the data has been registered through AddKnown, invoke this method to solve the system of eq...
LeastSquares(Isis::BasisFunction &basis, bool sparse=false, int sparseRows=0, int sparseCols=0, bool jigsaw=false)
Creates a LeastSquares Object.
std::vector< double > p_epsilonsSparse
sparse vector of total parameter corrections
Definition: LeastSquares.h:175
std::vector< double > p_expected
A vector of the expected values when solved.
Definition: LeastSquares.h:199
std::vector< double > Residuals() const
Returns a vector of residuals (errors).
~LeastSquares()
Destroys the LeastSquares object.
std::vector< std::vector< double > > p_input
A vector of the input variables to evaluate.
Definition: LeastSquares.h:197
void FillSparseA(const std::vector< double > &data)
Invoke this method for each set of knowns for sparse solutions.
bool SparseErrorPropagation()
Error propagation for sparse least-squares solution.
std::vector< double > p_xSparse
sparse solution vector
Definition: LeastSquares.h:174
std::vector< double > p_parameterWeights
vector of parameter weights
Definition: LeastSquares.h:176
Singular Value Decomposition.
Definition: LeastSquares.h:130
double Evaluate(const std::vector< double > &input)
Invokes the BasisFunction Evaluate method.
std::vector< double > p_residuals
A vector of the residuals (or difference between expected and solved values).
Definition: LeastSquares.h:204
gmm::row_matrix< gmm::rsvector< double > > p_normals
normal equations matrix &#39;N&#39;
Definition: LeastSquares.h:180
gmm::dense_matrix< double > p_ATb
right-hand side vector
Definition: LeastSquares.h:182
Generic least square fitting class.
Definition: LeastSquares.h:117
gmm::SuperLU_factor< double > p_SLU_Factor
decomposed normal equations matrix
Definition: LeastSquares.h:183
int SolveSparse()
Solve using sparse class.
double p_sigma0
sigma nought - reference variance
Definition: LeastSquares.h:195
void AddKnown(const std::vector< double > &input, double expected, double weight=1.0)
Invoke this method for each set of knowns.
double Residual(int i) const
Returns the ith residual.
int p_degreesOfFreedom
degrees of freedom (redundancy)
Definition: LeastSquares.h:193
std::vector< double > GetInput(int row) const
This method returns the data at the given row.
Generic linear equation class.
Definition: BasisFunction.h:64
double GetExpected(int row) const
This method returns the expected value at the given row.
Isis::BasisFunction * p_basis
Pointer to the BasisFunction object.
Definition: LeastSquares.h:208
std::vector< double > p_sqrtWeight
A vector of the square roots of the weights for each known value.
Definition: LeastSquares.h:201
int Knowns() const
The number of knowns (or times AddKnown was invoked) linear combination of the variables.
Definition: LeastSquares.h:147
bool p_solved
Boolean value indicating solution is complete.
Definition: LeastSquares.h:187
void SolveQRD()
After all the data has been registered through AddKnown, invoke this method to solve the system of eq...
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...
int Rows() const
This methods returns the number of rows in the matrix.
int p_constrainedParameters
constrained parameters
Definition: LeastSquares.h:192

U.S. Department of the Interior | U.S. Geological Survey
ISIS | Privacy & Disclaimers | Astrogeology Research Program
To contact us, please post comments and questions on the ISIS Support Center
File Modified: 07/12/2023 23:22:00