Isis 3 Programmer Reference
LeastSquares.cpp
Go to the documentation of this file.
1 
23 #include "jama/jama_svd.h"
24 #include "jama/jama_qr.h"
25 
26 #include <armadillo>
27 
28 #include "LeastSquares.h"
29 #include "IException.h"
30 #include "IString.h"
31 
32 namespace Isis {
40  int sparseRows, int sparseCols, bool jigsaw) {
41  p_jigsaw = jigsaw;
42  p_basis = &basis;
43  p_solved = false;
44  p_sparse = sparse;
45  p_sigma0 = 0.;
46 
47  p_sparseRows = sparseRows;
48  p_sparseCols = sparseCols;
49 
50 
51  if (p_sparse) {
52 
53  // make sure sparse nrows/ncols have been set
54  if (sparseRows == 0 || sparseCols == 0) {
55  QString msg = "If solving using sparse matrices, you must enter the "
56  "number of rows/columns";
58  }
59 
60 
61  p_sparseA.set_size(sparseRows, sparseCols);
62  p_normals.set_size(sparseCols, sparseCols);
63  p_ATb.resize(sparseCols, 1);
64  p_xSparse.resize(sparseCols);
65 
66  if( p_jigsaw ) {
67  p_epsilonsSparse.resize(sparseCols);
68  std::fill_n(p_epsilonsSparse.begin(), sparseCols, 0.0);
69 
70  p_parameterWeights.resize(sparseCols);
71  }
72 
73  }
74  p_currentFillRow = -1;
75  }
76 
79  }
80 
112  void LeastSquares::AddKnown(const std::vector<double> &data, double result,
113  double weight) {
114  if((int) data.size() != p_basis->Variables()) {
115  QString msg = "Number of elements in data does not match basis [" +
116  p_basis->Name() + "] requirements";
118  }
119 
120  p_expected.push_back(result);
121 
122  if (weight == 1) {
123  p_sqrtWeight.push_back(weight);
124  }
125  else {
126  p_sqrtWeight.push_back(sqrt(weight));
127  }
128 
129  if(p_sparse) {
130  FillSparseA(data);
131  }
132  else {
133  p_input.push_back(data);
134  }
135  }
136 
137 
138 
153  void LeastSquares::FillSparseA(const std::vector<double> &data) {
154 
155  p_basis->Expand(data);
156 
157  p_currentFillRow++;
158 
159  int ncolumns = (int)data.size();
160 
161  for(int c = 0; c < ncolumns; c++) {
162  p_sparseA(p_currentFillRow, c) = p_basis->Term(c) * p_sqrtWeight[p_currentFillRow];
163  }
164  }
165 
166 
174  std::vector<double> LeastSquares::GetInput(int row) const {
175  if((row >= Rows()) || (row < 0)) {
176  QString msg = "Index out of bounds [Given = " + toString(row) + "]";
178  }
179  return p_input[row];
180  }
181 
189  double LeastSquares::GetExpected(int row) const {
190  if((row >= Rows()) || (row < 0)) {
191  QString msg = "Index out of bounds [Given = " + toString(row) + "]";
193  }
194  return p_expected[row];
195  }
196 
203  int LeastSquares::Rows() const {
204  return (int)p_input.size();
205  }
206 
222 
223  if((method == SPARSE && p_sparseRows == 0) ||
224  (method != SPARSE && Rows() == 0 )) {
225  p_solved = false;
226  QString msg = "No solution available because no input data was provided";
228  }
229 
230  if(method == SVD) {
231  SolveSVD();
232  }
233  else if(method == QRD) {
234  SolveQRD();
235  }
236  else if(method == SPARSE) {
237  int column = SolveSparse();
238  return column;
239  }
240  return 0;
241  }
242 
253 
254  // We are solving Ax=b ... start by creating A
255  TNT::Array2D<double> A(p_input.size(), p_basis->Coefficients());
256  for(int r = 0; r < A.dim1(); r++) {
257  p_basis->Expand(p_input[r]);
258  for(int c = 0; c < A.dim2(); c++) {
259  A[r][c] = p_basis->Term(c) * p_sqrtWeight[r];
260  }
261  }
262 
263  // Ok use singular value decomposition to solve for the coefficients
264  // A = [U][S][V'] where [U] is MxN, [S] is NxN, [V'] is NxN transpose
265  // of [V]. We are solving for [A]x=b and need inverse of [A] such
266  // that x = [invA]b. Since inverse may not exist we use the
267  // pseudo-inverse [A+] from SVD which is [A+] = [V][invS][U']
268  // Our coefficents are then x = [A+]b where b is p_b.
269  JAMA::SVD<double> svd(A);
270 
271  TNT::Array2D<double> V;
272  svd.getV(V);
273 
274  // The inverse of S is the 1 over each diagonal element of S
275  TNT::Array2D<double> invS;
276  svd.getS(invS);
277 
278  for(int i = 0; i < invS.dim1(); i++) {
279  if(invS[i][i] != 0.0) invS[i][i] = 1.0 / invS[i][i];
280  }
281 
282  // Transpose U
283  TNT::Array2D<double> U;
284  svd.getU(U);
285  TNT::Array2D<double> transU(U.dim2(), U.dim1());
286 
287  for(int r = 0; r < U.dim1(); r++) {
288  for(int c = 0; c < U.dim2(); c++) {
289  transU[c][r] = U[r][c];
290  }
291  }
292 
293  // Now multiply everything together to get [A+]
294  TNT::Array2D<double> VinvS = TNT::matmult(V, invS);
295  TNT::Array2D<double> Aplus = TNT::matmult(VinvS, transU);
296 
297  // Using Aplus and our b we can solve for the coefficients
298  TNT::Array2D<double> b(p_expected.size(), 1);
299 
300  for(int r = 0; r < (int)p_expected.size(); r++) {
301  b[r][0] = p_expected[r] * p_sqrtWeight[r];
302  }
303 
304  TNT::Array2D<double> coefs = TNT::matmult(Aplus, b);
305 
306  // If the rank of the matrix is not large enough we don't
307  // have enough coefficients for the solution
308  if (coefs.dim1() < p_basis->Coefficients()) {
309  QString msg = "Unable to solve least-squares using SVD method. No "
310  "solution available. Not enough knowns or knowns are "
311  "co-linear ... [Unknowns = "
312  + toString(p_basis->Coefficients()) + "] [Knowns = "
313  + toString(coefs.dim1()) + "]";
315  }
316 
317  // Set the coefficients in our basis equation
318  std::vector<double> bcoefs;
319  for (int i = 0; i < coefs.dim1(); i++) bcoefs.push_back(coefs[i][0]);
320 
321  p_basis->SetCoefficients(bcoefs);
322 
323  // Compute the errors
324  for(int i = 0; i < (int)p_input.size(); i++) {
325  double value = p_basis->Evaluate(p_input[i]);
326  p_residuals.push_back(value - p_expected[i]);
328  }
329  // calculate degrees of freedom (or redundancy)
330  // DOF = # observations + # constrained parameters - # unknown parameters
331  p_degreesOfFreedom = p_basis->Coefficients() - coefs.dim1();
332 
333  if( p_degreesOfFreedom > 0.0 ) {
335  }
336 
337  // check for p_sigma0 < 0
338  p_sigma0 = sqrt(p_sigma0);
339 
340  // All done
341  p_solved = true;
342  }
343 
356 
357  // We are solving Ax=b ... start by creating an MxN matrix, A
358  TNT::Array2D<double> A(p_input.size(), p_basis->Coefficients());
359  for(int r = 0; r < A.dim1(); r++) {
360  p_basis->Expand(p_input[r]);
361  for(int c = 0; c < A.dim2(); c++) {
362  A[r][c] = p_basis->Term(c) * p_sqrtWeight[r];
363  }
364  }
365 
366  // Ok use to solve for the coefficients
367  // [A] = [Q][R] where [Q] is MxN and orthogonal and [R] is an NxN,
368  // upper triangular matrix. TNT provides the solve method that inverts
369  // [Q] and backsolves [R] to get the coefficients in the vector x.
370  // That is, we solve the system Rx = Q^T b
371  JAMA::QR<double> qr(A);
372 
373  // Using A and our b we can solve for the coefficients
374  TNT::Array1D<double> b(p_expected.size());
375  for(int r = 0; r < (int)p_expected.size(); r++) {
376  b[r] = p_expected[r] * p_sqrtWeight[r];
377  }// by construction, we know the size of b is equal to M, so b is conformant
378 
379  // Check to make sure the matrix is full rank before solving
380  // -- rectangular matrices must be full rank in order for the solve method
381  // to be successful
382  int full = qr.isFullRank();
383  if(full == 0) {
384  QString msg = "Unable to solve-least squares using QR Decomposition. "
385  "The upper triangular R matrix is not full rank";
387  }
388 
389  TNT::Array1D<double> coefs = qr.solve(b);
390 
391  // Set the coefficients in our basis equation
392  std::vector<double> bcoefs;
393  for(int i = 0; i < coefs.dim1(); i++) {
394  bcoefs.push_back(coefs[i]);
395  }
396  p_basis->SetCoefficients(bcoefs);
397 
398  // Compute the errors
399  for(int i = 0; i < (int)p_input.size(); i++) {
400  double value = p_basis->Evaluate(p_input[i]);
401  p_residuals.push_back(value - p_expected[i]);
402  }
403 
404  // All done
405  p_solved = true;
406  }
407 
408 
409 
410 
445 
446  // form "normal equations" matrix by multiplying ATA
448 
449  // Create the right-hand-side column vector 'b'
450  arma::mat b(p_sparseRows, 1);
451 
452  // multiply each element of 'b' by it's associated weight
453  for ( int r = 0; r < p_sparseRows; r++ )
454  b(r,0) = p_expected[r] * p_sqrtWeight[r];
455 
456  // form ATb
457  p_ATb = p_sparseA.t()*b;
458 
459  // apply parameter weighting if Jigsaw (bundle adjustment)
460  if ( p_jigsaw ) {
461  for( int i = 0; i < p_sparseCols; i++) {
462  double weight = p_parameterWeights[i];
463 
464  if( weight <= 0.0 )
465  continue;
466 
467  p_normals(i, i) += weight;
468  p_ATb(i, 0) -= p_epsilonsSparse[i]*weight;
469  }
470  }
471 
472  bool status = spsolve(p_xSparse, p_normals, p_ATb, "superlu");
473 
474  if (status == false) {
475  QString msg = "Could not solve sparse least squares problem.";
477  }
478 
479  // Set the coefficients in our basis equation
480  p_basis->SetCoefficients(arma::conv_to< std::vector<double> >::from(p_xSparse));
481 
482  // if Jigsaw (bundle adjustment)
483  // add corrections into epsilon vector (keeping track of total corrections)
484  if ( p_jigsaw ) {
485  for( int i = 0; i < p_sparseCols; i++ )
486  p_epsilonsSparse[i] += p_xSparse[i];
487  }
488 
489  // Compute the image coordinate residuals and sum into Sigma0
490  // (note this is exactly what was being done before, but with less overhead - I think)
491  // ultimately, we should not be using the A matrix but forming the normals
492  // directly. Then we'll have to compute the residuals by back projection
493 
494  p_residuals.resize(p_sparseRows);
495  p_residuals = arma::conv_to< std::vector<double> >::from(p_sparseA*p_xSparse);
496  p_sigma0 = 0.0;
497 
498  for ( int i = 0; i < p_sparseRows; i++ ) {
500  p_residuals[i] -= p_expected[i];
502  }
503 
504  // if Jigsaw (bundle adjustment)
505  // add contibution to Sigma0 from constrained parameters
506  if ( p_jigsaw ) {
507  double constrained_vTPv = 0.0;
508 
509  for ( int i = 0; i < p_sparseCols; i++ ) {
510  double weight = p_parameterWeights[i];
511 
512  if ( weight <= 0.0 )
513  continue;
514 
515  constrained_vTPv += p_epsilonsSparse[i]*p_epsilonsSparse[i]*weight;
516  }
517  p_sigma0 += constrained_vTPv;
518  }
519  // calculate degrees of freedom (or redundancy)
520  // DOF = # observations + # constrained parameters - # unknown parameters
521  p_degreesOfFreedom = p_sparseRows + p_constrainedParameters - p_sparseCols;
522 
523  if( p_degreesOfFreedom <= 0.0 ) {
524  p_sigma0 = 1.0;
525  }
526  else {
528  }
529 
530  // check for p_sigma0 < 0
531  p_sigma0 = sqrt(p_sigma0);
532 
533  // All done
534  p_solved = true;
535  return 0;
536  }
537 
538 
539  void LeastSquares::Reset ()
540  {
541  if ( p_sparse ) {
542  p_sparseA.zeros();
543  p_ATb.zeros();
544  p_normals.zeros();
545  p_currentFillRow = -1;
546  }
547  else {
548  p_input.clear();
549  }
550  p_sigma0 = 0.;
551  p_residuals.clear();
552  p_expected.clear();
553  p_sqrtWeight.clear();
554  p_solved = false;
555  }
556 
557 
558 
569  double LeastSquares::Evaluate(const std::vector<double> &data) {
570  if(!p_solved) {
571  QString msg = "Unable to evaluate until a solution has been computed";
573  }
574  return p_basis->Evaluate(data);
575  }
576 
586  std::vector<double> LeastSquares::Residuals() const {
587  if(!p_solved) {
588  QString msg = "Unable to return residuals until a solution has been computed";
590  }
591  return p_residuals;
592  }
593 
606  double LeastSquares::Residual(int i) const {
607  if(!p_solved) {
608  QString msg = "Unable to return residuals until a solution has been computed";
610  }
611  return p_residuals[i];
612  }
613 
630  void LeastSquares::Weight(int index, double weight) {
631  if(weight == 1) {
632  p_sqrtWeight[index] = weight;
633  }
634  else {
635  p_sqrtWeight[index] = sqrt(weight);
636  }
637  }
638 
639 } // end namespace isis
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
void SetCoefficients(const std::vector< double > &coefs)
Set the coefficients for the equation.
QString Name() const
Returns the name of the equation.
Definition: BasisFunction.h:96
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 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.
QString toString(bool boolToConvert)
Global function to convert a boolean to a string.
Definition: IString.cpp:226
double Term(int c) const
Returns the cth term.
This error is for when a programmer made an API call that was illegal.
Definition: IException.h:162
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.
virtual void Expand(const std::vector< double > &vars)
This is the function you should replace depending on your needs.
std::vector< double > p_residuals
A vector of the residuals (or difference between expected and solved values).
Definition: LeastSquares.h:195
#define _FILEINFO_
Macro for the filename and line number.
Definition: IException.h:40
A type of error that cannot be classified as any of the other error types.
Definition: IException.h:134
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
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
double Evaluate(const std::vector< double > &vars)
Compute the equation using the input variables.
Generic linear equation class.
Definition: BasisFunction.h:64
Isis::BasisFunction * p_basis
Pointer to the BasisFunction object.
Definition: LeastSquares.h:199
Isis exception class.
Definition: IException.h:107
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
int Variables() const
Returns the number of variables in the equation.
Definition: BasisFunction.h:88
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 Coefficients() const
Returns the number of coefficients for the equation.
Definition: BasisFunction.h:80
int p_constrainedParameters
constrained parameters
Definition: LeastSquares.h:183