Isis 3 Programmer Reference
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 #include <armadillo>
30 
31 #include "BasisFunction.h"
32 
33 namespace Isis {
115  class LeastSquares {
116  public:
117 
118  LeastSquares(Isis::BasisFunction &basis, bool sparse = false,
119  int sparseRows = 0, int sparseCols = 0, bool jigsaw = false);
120  ~LeastSquares();
121  void AddKnown(const std::vector<double> &input, double expected,
122  double weight = 1.0);
123 
124  std::vector<double> GetInput(int row) const;
125  double GetExpected(int row) const;
126  int Rows() const;
127 
128  enum SolveMethod { SVD,
129  QRD,
131  };
132 
134  double Evaluate(const std::vector<double> &input);
135  std::vector<double> Residuals() const;
136  double Residual(int i) const;
137  void Weight(int index, double weight);
138 
145  int Knowns() const {
146  return p_expected.size();
147  };
148 
149  double GetSigma0() { return p_sigma0; }
150  int GetDegreesOfFreedom() { return p_degreesOfFreedom; }
151  void Reset ();
152 
153  void ResetSparse() { Reset(); }
154  std::vector<double> GetEpsilons () const { return p_epsilonsSparse; }
155  void SetParameterWeights(const std::vector<double> weights) { p_parameterWeights = weights; }
156  void SetNumberOfConstrainedParameters(int n) { p_constrainedParameters = n; }
157 
158  private:
159  void SolveSVD();
160  void SolveQRD();
161  void SolveCholesky () {}
162 
163  int SolveSparse();
164  void FillSparseA(const std::vector<double> &data);
165  bool ApplyParameterWeights();
166 
167  arma::mat p_xSparse;
168  std::vector<double> p_epsilonsSparse;
169  std::vector<double> p_parameterWeights;
171  arma::SpMat<double> p_sparseA;
172  arma::SpMat<double> p_normals;
173  arma::mat p_ATb;
174  arma::mat p_SLU_Factor;
176  bool p_jigsaw;
177  bool p_sparse;
178  bool p_solved;
180  int p_currentFillRow;
181  int p_sparseRows;
182  int p_sparseCols;
186  double p_sigma0;
188  std::vector<std::vector<double> > p_input;
190  std::vector<double> p_expected;
192  std::vector<double> p_sqrtWeight;
195  std::vector<double> p_residuals;
200 
201 
202 
203  };
204 };
205 
206 
207 #endif
arma::mat p_xSparse
sparse solution matrix
Definition: LeastSquares.h:167
void Weight(int index, double weight)
Reset the weight for the ith known.
QR Decomposition.
Definition: LeastSquares.h:129
std::vector< double > Residuals() const
Returns a vector of residuals (errors).
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:168
double Residual(int i) const
Returns the ith residual.
int Knowns() const
The number of knowns (or times AddKnown was invoked) linear combination of the variables.
Definition: LeastSquares.h:145
int Rows() const
This methods returns the number of rows in the matrix.
std::vector< double > GetInput(int row) const
This method returns the data at the given row.
std::vector< double > p_expected
A vector of the expected values when solved.
Definition: LeastSquares.h:190
~LeastSquares()
Destroys the LeastSquares object.
std::vector< std::vector< double > > p_input
A vector of the input variables to evaluate.
Definition: LeastSquares.h:188
void FillSparseA(const std::vector< double > &data)
Invoke this method for each set of knowns for sparse solutions.
std::vector< double > p_parameterWeights
vector of parameter weights
Definition: LeastSquares.h:169
Singular Value Decomposition.
Definition: LeastSquares.h:128
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:195
double GetExpected(int row) const
This method returns the expected value at the given row.
arma::SpMat< double > p_sparseA
design matrix &#39;A&#39;
Definition: LeastSquares.h:171
Generic least square fitting class.
Definition: LeastSquares.h:115
arma::SpMat< double > p_normals
normal equations matrix &#39;N&#39;
Definition: LeastSquares.h:172
int SolveSparse()
Solve using sparse class.
double p_sigma0
sigma nought - reference variance
Definition: LeastSquares.h:186
void AddKnown(const std::vector< double > &input, double expected, double weight=1.0)
Invoke this method for each set of knowns.
int p_degreesOfFreedom
degrees of freedom (redundancy)
Definition: LeastSquares.h:184
arma::mat p_ATb
right-hand side vector
Definition: LeastSquares.h:173
Generic linear equation class.
Definition: BasisFunction.h:64
Isis::BasisFunction * p_basis
Pointer to the BasisFunction object.
Definition: LeastSquares.h:199
Namespace for ISIS/Bullet specific routines.
Definition: Apollo.h:31
std::vector< double > p_sqrtWeight
A vector of the square roots of the weights for each known value.
Definition: LeastSquares.h:192
arma::mat p_SLU_Factor
decomposed normal equations matrix
Definition: LeastSquares.h:174
bool p_solved
Boolean value indicating solution is complete.
Definition: LeastSquares.h:178
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 p_constrainedParameters
constrained parameters
Definition: LeastSquares.h:183