Isis 3 Programmer Reference
BundleAdjust.h
Go to the documentation of this file.
1 #ifndef BundleAdjust_h
2 #define BundleAdjust_h
3 
26 // Qt lib
27 #include <QObject> // parent class
28 
29 // std lib
30 #include <vector>
31 #include <fstream>
32 
33 // cholmod lib
34 #include <cholmod.h>
35 
36 // Isis lib
37 #include "BundleControlPoint.h"
40 #include "BundleResults.h"
41 #include "BundleSettings.h"
42 #include "BundleSolutionInfo.h"
43 #include "BundleTargetBody.h"
44 #include "Camera.h"
45 #include "CameraGroundMap.h"
46 #include "ControlMeasure.h"
47 #include "ControlNet.h"
48 #include "LinearAlgebra.h"
49 #include "MaximumLikelihoodWFunctions.h" // why not just forward declare???
50 #include "ObservationNumberList.h"
51 #include "Progress.h"
52 #include "SerialNumberList.h"
53 #include "SparseBlockMatrix.h"
54 #include "Statistics.h"
55 
56 template< typename T > class QList;
57 template< typename A, typename B > class QMap;
58 
59 namespace Isis {
60  class Control;
61  class ImageList;
62 
333  class BundleAdjust : public QObject {
334  Q_OBJECT
335  public:
336  BundleAdjust(BundleSettingsQsp bundleSettings,
337  const QString &cnetFile,
338  const QString &cubeList,
339  bool printSummary = true);
340  BundleAdjust(BundleSettingsQsp bundleSettings,
341  QString &cnet,
342  SerialNumberList &snlist,
343  bool printSummary = true);
344  BundleAdjust(BundleSettingsQsp bundleSettings,
345  Control &cnet,
346  SerialNumberList &snlist,
347  bool bPrintSummary);
348  BundleAdjust(BundleSettingsQsp bundleSettings,
349  ControlNet &cnet,
350  SerialNumberList &snlist,
351  bool printSummary = true);
352  BundleAdjust(BundleSettingsQsp bundleSettings,
353  ControlNetQsp cnet,
354  const QString &cubeList,
355  bool printSummary = true);
356  BundleAdjust(BundleSettingsQsp bundleSettings,
357  Control &control,
358  QList<ImageList *> imgList,
359  bool printSummary);
360  ~BundleAdjust();
362 
364  bool isAborted();
365 
366  public slots:
367  bool solveCholesky();
368  void abortBundle();
369  void outputBundleStatus(QString status);
370 
371  // accessors
372 
375  QString fileName(int index);
376  QString iterationSummaryGroup() const;
377  bool isConverged();
378  Table cMatrix(int index);
379  Table spVector(int index);
380  int numberOfImages() const;
381  double iteration() const;
382 
383  signals:
384  void statusUpdate(QString);
385  void error(QString);
386  void iterationUpdate(int);
387  void pointUpdate(int);
388  void statusBarUpdate(QString);
389  void resultsReady(BundleSolutionInfo *bundleSolveInformation);
390  void finished();
391 
392  private:
393  //TODO Should there be a resetBundle(BundleSettings bundleSettings) method
394  // that allows for rerunning with new settings? JWB
395  void init(Progress *progress = 0);
397  bool validateNetwork();
398  bool solveSystem();
399  void iterationSummary();
403  bool errorPropagation();
404  double computeResiduals();
405  bool computeRejectionLimit();
406  bool flagOutliers();
407 
408  // normal equation matrices methods
409 
410  bool formNormalEquations();
411  bool computePartials(LinearAlgebra::Matrix &coeffTarget,
412  LinearAlgebra::Matrix &coeffImage,
413  LinearAlgebra::Matrix &coeffPoint3D,
414  LinearAlgebra::Vector &coeffRHS,
415  BundleMeasure &measure,
416  BundleControlPoint &point);
417  bool formMeasureNormals(boost::numeric::ublas::symmetric_matrix<
418  double, boost::numeric::ublas::upper > &N22,
420  boost::numeric::ublas::compressed_vector< double > &n1,
422  LinearAlgebra::Matrix &coeffTarget,
423  LinearAlgebra::Matrix &coeffImage,
424  LinearAlgebra::Matrix &coeffPoint3D,
425  LinearAlgebra::Vector &coeffRHS,
426  int observationIndex);
427  bool formPointNormals(boost::numeric::ublas::symmetric_matrix<
428  double, boost::numeric::ublas::upper > &N22,
432  BundleControlPointQsp &point);
433  bool formWeightedNormals(boost::numeric::ublas::compressed_vector< double > &n1,
435 
436  // dedicated matrix functions
437 
440  void accumProductAlphaAB(double alpha,
444  bool invert3x3(boost::numeric::ublas::symmetric_matrix<
445  double, boost::numeric::ublas::upper > &m);
446  bool productATransB(boost::numeric::ublas::symmetric_matrix<
447  double, boost::numeric::ublas::upper > &N22,
450  void productAlphaAV(double alpha,
451  boost::numeric::ublas::bounded_vector< double, 3 > &v2,
454 
455  // CHOLMOD library methods
456 
459  bool cholmodInverse();
460  bool loadCholmodTriplet();
461  bool wrapUp();
462 
463  // member variables
464 
472  QString m_cnetFileName;
473  QVector <BundleControlPointQsp> m_bundleControlPoints;
481  BundleTargetBodyQsp m_bundleTargetBody;
483  bool m_abort;
484  QString m_iterationSummary;
489  bool m_cleanUp;
492  int m_rank;
495  QList<ImageList *> m_imageLists;
498  // ==========================================================================================
499  // === BEYOND THIS PLACE (THERE BE DRAGONS) all refers to the folded bundle solution. ===
500  // === Currently, everything uses the CHOLMOD library, ===
501  // === there is no dependence on the least-squares class. ===
502  // ==========================================================================================
503 
505  boost::numeric::ublas::symmetric_matrix<
506  double,
507  boost::numeric::ublas::upper,
508  boost::numeric::ublas::column_major > m_normalInverse;
509  cholmod_common m_cholmodCommon;
518  cholmod_triplet *m_cholmodTriplet;
524  cholmod_sparse *m_cholmodNormal;
529  cholmod_factor *m_L;
539  };
540 }
541 
542 #endif
This represents an ISIS control net in a project-based GUI interface.
Definition: Control.h:79
boost::numeric::ublas::symmetric_matrix< double, boost::numeric::ublas::upper, boost::numeric::ublas::column_major > m_normalInverse
!< The lists of images used in the bundle.
Definition: BundleAdjust.h:508
This class is a container class for BundleObservations.
bool initializeCHOLMODLibraryVariables()
Initializations for CHOLMOD sparse matrix package.
Unless noted otherwise, the portions of Isis written by the USGS are public domain.
SparseBlockMatrix.
QString fileName(int index)
Return the ith filename in the cube list file given to constructor.
Internalizes a list of images and allows for operations on the entire list.
Definition: ImageList.h:55
bool flagOutliers()
Flags outlier measures and control points.
bool m_printSummary
!< Summary of the most recently completed iteration.
Definition: BundleAdjust.h:486
int m_previousNumberImagePartials
!< The image parameter solution vector.
Definition: BundleAdjust.h:536
SparseBlockColumnMatrix.
void init(Progress *progress=0)
Initialize all solution parameters.
bool isConverged()
Returns if the BundleAdjust converged.
Container class for BundleAdjustment results.
SerialNumberList * serialNumberList()
Returns a pointer to the serial number list.
int m_numberOfImagePartials
number of image-related partials.
Definition: BundleAdjust.h:494
bool cholmodInverse()
Compute inverse of normal equations matrix for CHOLMOD.
bool formWeightedNormals(boost::numeric::ublas::compressed_vector< double > &n1, LinearAlgebra::Vector &nj)
Apply weighting for spacecraft position, velocity, acceleration and camera angles, angular velocities, angular accelerations if so stipulated (legalese).
BundleSolutionInfo * bundleSolveInformation()
Create a BundleSolutionInfo containing the settings and results from the bundle adjustment.
bool wrapUp()
Compute the residuals for each adjusted point and store adjustment results in m_bundleResults.
bool freeCHOLMODLibraryVariables()
Free CHOLMOD library variables.
SparseBlockMatrix m_sparseNormals
!< The right hand side of the normal equations.
Definition: BundleAdjust.h:514
Table spVector(int index)
Return a table spacecraft vector for the ith cube in the cube list given to the constructor.
Statistics m_xyResiduals
xy residual statistics.
Definition: BundleAdjust.h:470
Unless noted otherwise, the portions of Isis written by the USGS are public domain.
A container class for a ControlMeasure.
Definition: BundleMeasure.h:69
boost::numeric::ublas::matrix< double > Matrix
Definition for an Isis::LinearAlgebra::Matrix of doubles.
Unless noted otherwise, the portions of Isis written by the USGS are public domain.
bool initializeNormalEquationsMatrix()
Initialize Normal Equations matrix (m_sparseNormals).
SerialNumberList * m_serialNumberList
!< Vector of observations.
Definition: BundleAdjust.h:480
bool validateNetwork()
control network validation - on the very real chance that the net has not been checked before running...
BundleSettingsQsp m_bundleSettings
Contains the solve settings.
Definition: BundleAdjust.h:465
void abortBundle()
Flag to abort when bundle is threaded.
This class is used to accumulate statistics on double arrays.
Definition: Statistics.h:107
LinearAlgebra::Vector m_imageSolution
!< The lower triangular L matrix from Cholesky decomposition.
Definition: BundleAdjust.h:533
bool solveCholesky()
Compute the least squares bundle adjustment solution using Cholesky decomposition.
boost::numeric::ublas::vector< double > Vector
Definition for an Isis::LinearAlgebra::Vector of doubles.
bool formMeasureNormals(boost::numeric::ublas::symmetric_matrix< double, boost::numeric::ublas::upper > &N22, SparseBlockColumnMatrix &N12, boost::numeric::ublas::compressed_vector< double > &n1, LinearAlgebra::Vector &n2, LinearAlgebra::Matrix &coeffTarget, LinearAlgebra::Matrix &coeffImage, LinearAlgebra::Matrix &coeffPoint3D, LinearAlgebra::Vector &coeffRHS, int observationIndex)
Form the auxilary normal equation matrices for a measure.
Program progress reporter.
Definition: Progress.h:58
~BundleAdjust()
Destroys BundleAdjust object, deallocates pointers (if we have ownership), and frees variables from c...
a control network
Definition: ControlNet.h:271
bool loadCholmodTriplet()
Load sparse normal equations matrix into CHOLMOD triplet.
Table cMatrix(int index)
Return a table cmatrix for the ith cube in the cube list given to the constructor.
BundleAdjust(BundleSettingsQsp bundleSettings, const QString &cnetFile, const QString &cubeList, bool printSummary=true)
Construct a BundleAdjust object from the given settings, control network file, and cube list...
int m_rank
!< If the serial number lists need to be deleted by the destructor.
Definition: BundleAdjust.h:492
ControlNetQsp controlNet()
Returns a pointer to the output control network.
cholmod_triplet * m_cholmodTriplet
!< The sparse block normal equations matrix.
Definition: BundleAdjust.h:518
bool solveSystem()
Compute the solution to the normal equations using the CHOLMOD library.
bool productATransB(boost::numeric::ublas::symmetric_matrix< double, boost::numeric::ublas::upper > &N22, SparseBlockColumnMatrix &N12, SparseBlockRowMatrix &Q)
Perform the matrix multiplication Q = N22 x N12(transpose)
void accumProductAlphaAB(double alpha, SparseBlockRowMatrix &A, LinearAlgebra::Vector &B, LinearAlgebra::Vector &C)
Performs the matrix multiplication nj = nj + alpha (Q x n2).
void applyParameterCorrections()
apply parameter corrections for solution.
Statistics m_yResiduals
y residual statistics.
Definition: BundleAdjust.h:469
This class holds information about a control point that BundleAdjust needs to run correctly...
bool m_cleanUp
!< If the iteration summaries should be output to the log file.
Definition: BundleAdjust.h:489
BundleSolutionInfo * solveCholeskyBR()
Compute the least squares bundle adjustment solution using Cholesky decomposition.
BundleResults m_bundleResults
Stores the results of the bundle adjust.
Definition: BundleAdjust.h:466
LinearAlgebra::Vector m_RHS
!< Contains object parameters, statistics, and workspace used by the CHOLMOD library.
Definition: BundleAdjust.h:512
bool isAborted()
Returns if the BundleAdjust has been aborted.
int m_iteration
The current iteration.
Definition: BundleAdjust.h:493
QString m_cnetFileName
The control net filename.
Definition: BundleAdjust.h:472
Statistics m_xResiduals
x residual statistics.
Definition: BundleAdjust.h:468
void outputBundleStatus(QString status)
Slot for deltack and jigsaw to output the bundle status.
QList< ImageList * > imageLists()
This method returns the image list used in the bundle adjust.
bool formNormalEquations()
Form the least-squares normal equations matrix via cholmod.
bool formPointNormals(boost::numeric::ublas::symmetric_matrix< double, boost::numeric::ublas::upper > &N22, SparseBlockColumnMatrix &N12, LinearAlgebra::Vector &n2, LinearAlgebra::Vector &nj, BundleControlPointQsp &point)
Compute the Q matrix and NIC vector for a control point.
ControlNetQsp m_controlNet
Output control net.
Definition: BundleAdjust.h:471
int numberOfImages() const
Returns the number of images.
BundleObservationVector m_bundleObservations
!< Vector of control points.
Definition: BundleAdjust.h:477
Class for storing Table blobs information.
Definition: Table.h:77
void productAB(SparseBlockColumnMatrix &A, SparseBlockRowMatrix &B)
Perform the matrix multiplication C = N12 x Q.
bool computeRejectionLimit()
Compute rejection limit.
Namespace for ISIS/Bullet specific routines.
Definition: Apollo.h:31
A container class for statistical results from a BundleAdjust solution.
Definition: BundleResults.h:96
void iterationSummary()
Creates an iteration summary and an iteration group for the solution summary.
QString iterationSummaryGroup() const
Returns the iteration summary string.
cholmod_sparse * m_cholmodNormal
!< The CHOLMOD triplet representation of the sparse normal equations matrix.
Definition: BundleAdjust.h:524
double iteration() const
Returns what iteration the BundleAdjust is currently on.
An image bundle adjustment object.
Definition: BundleAdjust.h:333
bool errorPropagation()
Error propagation for solution.
double computeResiduals()
This method computes the focal plane residuals for the measures.
Serial Number list generator.
SparseBlockRowMatrix.
cholmod_factor * m_L
!< The CHOLMOD sparse normal equations matrix used by cholmod_factorize to solve the system...
Definition: BundleAdjust.h:529
bool computePartials(LinearAlgebra::Matrix &coeffTarget, LinearAlgebra::Matrix &coeffImage, LinearAlgebra::Matrix &coeffPoint3D, LinearAlgebra::Vector &coeffRHS, BundleMeasure &measure, BundleControlPoint &point)
Compute partial derivatives and weighted residuals for a measure.
bool m_abort
!< Contains information about the target body.
Definition: BundleAdjust.h:483
bool computeBundleStatistics()
Compute Bundle statistics and store them in m_bundleResults.
bool invert3x3(boost::numeric::ublas::symmetric_matrix< double, boost::numeric::ublas::upper > &m)
Dedicated quick inverse of 3x3 matrix.