Isis 3.0 Programmer Reference
Back | Home
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 
271  class BundleAdjust : public QObject {
272  Q_OBJECT
273  public:
274  BundleAdjust(BundleSettingsQsp bundleSettings,
275  const QString &cnetFile,
276  const QString &cubeList,
277  bool printSummary = true);
278  BundleAdjust(BundleSettingsQsp bundleSettings,
279  QString &cnet,
280  SerialNumberList &snlist,
281  bool printSummary = true);
282  BundleAdjust(BundleSettingsQsp bundleSettings,
283  Control &cnet,
284  SerialNumberList &snlist,
285  bool bPrintSummary);
286  BundleAdjust(BundleSettingsQsp bundleSettings,
287  ControlNet &cnet,
288  SerialNumberList &snlist,
289  bool printSummary = true);
290  BundleAdjust(BundleSettingsQsp bundleSettings,
291  ControlNetQsp cnet,
292  const QString &cubeList,
293  bool printSummary = true);
294  BundleAdjust(BundleSettingsQsp bundleSettings,
295  Control &control,
296  QList<ImageList *> imgList,
297  bool printSummary);
298  ~BundleAdjust();
300 
301  public slots:
302  bool solveCholesky();
303  void abortBundle();
304  void outputBundleStatus(QString status);
305 
306  // accessors
307 
310  QString fileName(int index);
311  QString iterationSummaryGroup() const;
312  bool isConverged();
313  Table cMatrix(int index);
314  Table spVector(int index);
315  int numberOfImages() const;
316  double iteration() const;
317 
318  signals:
319  void statusUpdate(QString);
320  void error(QString);
321  void iterationUpdate(int, double);
322  void resultsReady(BundleSolutionInfo *bundleSolveInformation);
323  void finished();
324 
325  private:
326  //TODO Should there be a resetBundle(BundleSettings bundleSettings) method
327  // that allows for rerunning with new settings? JWB
328  void init(Progress *progress = 0);
329  bool validateNetwork();
330  bool solveSystem();
331  void iterationSummary();
335  bool errorPropagation();
336  double computeResiduals();
337  bool computeRejectionLimit();
338  bool flagOutliers();
339 
340  // normal equation matrices methods
341 
342  bool formNormalEquations();
343  bool computePartials(LinearAlgebra::Matrix &coeffTarget,
344  LinearAlgebra::Matrix &coeffImage,
345  LinearAlgebra::Matrix &coeffPoint3D,
346  LinearAlgebra::Vector &coeffRHS,
347  BundleMeasure &measure,
348  BundleControlPoint &point);
349  bool formMeasureNormals(boost::numeric::ublas::symmetric_matrix<
350  double, boost::numeric::ublas::upper > &N22,
352  boost::numeric::ublas::compressed_vector< double > &n1,
354  LinearAlgebra::Matrix &coeffTarget,
355  LinearAlgebra::Matrix &coeffImage,
356  LinearAlgebra::Matrix &coeffPoint3D,
357  LinearAlgebra::Vector &coeffRHS,
358  int observationIndex);
359  bool formPointNormals(boost::numeric::ublas::symmetric_matrix<
360  double, boost::numeric::ublas::upper > &N22,
364  BundleControlPointQsp &point);
365  bool formWeightedNormals(boost::numeric::ublas::compressed_vector< double > &n1,
367 
368  // dedicated matrix functions
369 
372  void accumProductAlphaAB(double alpha,
376  bool invert3x3(boost::numeric::ublas::symmetric_matrix<
377  double, boost::numeric::ublas::upper > &m);
378  bool productATransB(boost::numeric::ublas::symmetric_matrix<
379  double, boost::numeric::ublas::upper > &N22,
382  void productAlphaAV(double alpha,
383  boost::numeric::ublas::bounded_vector< double, 3 > &v2,
386 
387  // CHOLMOD library methods
388 
391  bool cholmodInverse();
392  bool loadCholmodTriplet();
393  bool wrapUp();
394 
395  // member variables
396 
404  QString m_cnetFileName;
405  QVector <BundleControlPointQsp> m_bundleControlPoints;
413  BundleTargetBodyQsp m_bundleTargetBody;
416  bool m_abort;
417  QString m_iterationSummary;
422  bool m_cleanUp;
425  int m_rank;
428  double m_radiansToMeters;
433  // ==========================================================================================
434  // === BEYOND THIS PLACE (THERE BE DRAGONS) all refers to the folded bundle solution. ===
435  // === Currently, everything uses the CHOLMOD library, ===
436  // === there is no dependence on the least-squares class. ===
437  // ==========================================================================================
438 
440  boost::numeric::ublas::symmetric_matrix<
441  double,
442  boost::numeric::ublas::upper,
443  boost::numeric::ublas::column_major > m_normalInverse;
444  cholmod_common m_cholmodCommon;
453  cholmod_triplet *m_cholmodTriplet;
459  cholmod_sparse *m_cholmodNormal;
464  cholmod_factor *m_L;
470  };
471 }
472 
473 #endif
474 
This represents an ISIS control net in a project-based GUI interface.
Definition: Control.h:57
boost::numeric::ublas::symmetric_matrix< double, boost::numeric::ublas::upper, boost::numeric::ublas::column_major > m_normalInverse
!&lt; The body specific meters to radians conversion factor.
Definition: BundleAdjust.h:443
This class is a container class for BundleObservations.
BundleSolutionInfo bundleSolveInformation()
Create a BundleSolutionInfo containing the settings and results from the bundle adjustment.
bool initializeCHOLMODLibraryVariables()
Initializations for CHOLMOD sparse matrix package.
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:44
QSharedPointer< BundleSettings > BundleSettingsQsp
Definition for a BundleSettingsQsp, a shared pointer to a BundleSettings object.
bool flagOutliers()
Flags outlier measures and control points.
bool m_printSummary
!&lt; Summary of the most recently completed iteration.
Definition: BundleAdjust.h:419
SparseBlockColumnMatrix.
void init(Progress *progress=0)
Initialize all solution parameters.
QSharedPointer< ControlNet > ControlNetQsp
This typedef is for future implementation of target body.
Definition: ControlNet.h:446
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:427
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).
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
!&lt; The right hand side of the normal equations.
Definition: BundleAdjust.h:449
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:402
QString iterationSummaryGroup() const
Returns the iteration summary string.
QSharedPointer< BundleTargetBody > BundleTargetBodyQsp
Definition for BundleTargetBodyQsp, a QSharedPointer to a BundleTargetBody.
double m_metersToRadians
!&lt; The body specific radians to meters conversion factor.
Definition: BundleAdjust.h:430
double iteration() const
Returns what iteration the BundleAdjust is currently on.
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.
Distance measurement, usually in meters.
Definition: Distance.h:47
Unless noted otherwise, the portions of Isis written by the USGS are public domain.
SerialNumberList * m_serialNumberList
!&lt; Vector of observations.
Definition: BundleAdjust.h:412
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:397
void abortBundle()
Flag to abort when bundle is threaded.
This class is used to accumulate statistics on double arrays.
Definition: Statistics.h:109
LinearAlgebra::Vector m_imageSolution
!&lt; The lower triangular L matrix from Cholesky decomposition.
Definition: BundleAdjust.h:468
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:207
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
!&lt; If the serial number lists need to be deleted by the destructor.
Definition: BundleAdjust.h:425
ControlNetQsp controlNet()
Returns a pointer to the output control network.
cholmod_triplet * m_cholmodTriplet
!&lt; The sparse block normal equations matrix.
Definition: BundleAdjust.h:453
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:401
This class holds information about a control point that BundleAdjust needs to run corretly...
bool m_cleanUp
!&lt; If the iteration summaries should be output to the log file.
Definition: BundleAdjust.h:422
BundleResults m_bundleResults
Stores the results of the &lt; bundle adjust.
Definition: BundleAdjust.h:398
LinearAlgebra::Vector m_RHS
!&lt; Contains object parameters, statistics, and workspace used by the CHOLMOD library.
Definition: BundleAdjust.h:447
int m_iteration
The current iteration.
Definition: BundleAdjust.h:426
QString m_cnetFileName
The control net filename.
Definition: BundleAdjust.h:404
int numberOfImages() const
Returns the number of images.
Statistics m_xResiduals
x residual statistics.
Definition: BundleAdjust.h:400
void outputBundleStatus(QString status)
Slot for deltack and jigsaw to output the bundle status.
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.
BundleSolutionInfo solveCholeskyBR()
Compute the least squares bundle adjustment solution using Cholesky decomposition.
ControlNetQsp m_controlNet
Output control net.
Definition: BundleAdjust.h:403
BundleObservationVector m_bundleObservations
!&lt; Vector of control points.
Definition: BundleAdjust.h:409
Class for storing Table blobs information.
Definition: Table.h:74
void productAB(SparseBlockColumnMatrix &A, SparseBlockRowMatrix &B)
Perform the matrix multiplication C = N12 x Q.
bool computeRejectionLimit()
Compute rejection limit.
A container class for statistical results from a BundleAdjust solution.
Definition: BundleResults.h:90
void iterationSummary()
Creates an iteration summary and an iteration group for the solution summary.
cholmod_sparse * m_cholmodNormal
!&lt; The CHOLMOD triplet representation of the sparse normal equations matrix.
Definition: BundleAdjust.h:459
An image bundle adjustment object.
Definition: BundleAdjust.h:271
bool errorPropagation()
Error propagation for solution.
double computeResiduals()
This method computes the focal plane residuals for the measures.
QSharedPointer< BundleControlPoint > BundleControlPointQsp
Definition for BundleControlPointQSP, a shared pointer to a BundleControlPoint.
Serial Number list generator.
SparseBlockRowMatrix.
cholmod_factor * m_L
!&lt; The CHOLMOD sparse normal equations matrix used by cholmod_factorize to solve the system...
Definition: BundleAdjust.h:464
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.
void productAlphaAV(double alpha, boost::numeric::ublas::bounded_vector< double, 3 > &v2, SparseBlockRowMatrix &Q, LinearAlgebra::Vector &v1)
Perform the matrix multiplication v2 = alpha ( Q x v1 ).
bool m_abort
If the bundle should abort.
Definition: BundleAdjust.h:416
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.
Distance m_bodyRadii[3]
!&lt; Contains information about the target body.
Definition: BundleAdjust.h:415

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:14:59