Isis Developer Reference
BundleAdjust.h
Go to the documentation of this file.
1#ifndef BundleAdjust_h
2#define BundleAdjust_h
3
10/* SPDX-License-Identifier: CC0-1.0 */
11
12#include <QObject> // parent class
13
14// std lib
15#include <vector>
16#include <fstream>
17
18// cholmod lib
19#include <cholmod.h>
20
21// Isis lib
22#include "BundleControlPoint.h"
26#include "BundleResults.h"
27#include "BundleSettings.h"
28#include "BundleSolutionInfo.h"
29#include "BundleTargetBody.h"
30#include "Camera.h"
31#include "CameraGroundMap.h"
32#include "ControlMeasure.h"
33#include "ControlNet.h"
34#include "LidarData.h"
35#include "LinearAlgebra.h"
38#include "Progress.h"
39#include "SerialNumberList.h"
40#include "SparseBlockMatrix.h"
41#include "Statistics.h"
42
43template< typename T > class QList;
44template< typename A, typename B > class QMap;
45
46namespace Isis {
47 class Control;
48 class ImageList;
49
329 class BundleAdjust : public QObject {
330 Q_OBJECT
331 public:
332 BundleAdjust(BundleSettingsQsp bundleSettings,
333 const QString &cnetFile,
334 const QString &cubeList,
335 bool printSummary = true);
336 BundleAdjust(BundleSettingsQsp bundleSettings,
337 const QString &cnetFile,
338 const QString &cubeList,
339 const QString &lidarDataFile,
340 bool printSummary = true);
341 BundleAdjust(BundleSettingsQsp bundleSettings,
342 QString &cnet,
343 SerialNumberList &snlist,
344 bool printSummary = true);
345 BundleAdjust(BundleSettingsQsp bundleSettings,
346 Control &cnet,
347 SerialNumberList &snlist,
348 bool bPrintSummary);
349 BundleAdjust(BundleSettingsQsp bundleSettings,
350 ControlNet &cnet,
351 SerialNumberList &snlist,
352 bool printSummary = true);
353 BundleAdjust(BundleSettingsQsp bundleSettings,
354 ControlNetQsp cnet,
355 const QString &cubeList,
356 bool printSummary = true);
357 BundleAdjust(BundleSettingsQsp bundleSettings,
358 Control &control,
359 QList<ImageList *> imgList,
360 bool printSummary);
363
364 QList<ImageList *> imageLists();
365 bool isAborted();
366
367 public slots:
368 bool solveCholesky();
369 void abortBundle();
370 void outputBundleStatus(QString status);
371
372 // accessors
373
377 QString fileName(int index);
378 QString iterationSummaryGroup() const;
379 bool isConverged();
380 Table cMatrix(int index);
381 Table spVector(int index);
382 QString modelState(int index);
383 int numberOfImages() const;
384 double iteration() const;
385
386 signals:
387 void statusUpdate(QString);
388 void error(QString);
390 void pointUpdate(int);
391 void statusBarUpdate(QString);
392 void resultsReady(BundleSolutionInfo *bundleSolveInformation);
393 void finished();
394
395 private:
396 //TODO Should there be a resetBundle(BundleSettings bundleSettings) method
397 // that allows for rerunning with new settings? JWB
398 void init(Progress *progress = 0);
399 bool initializeNormalEquationsMatrix();
400 bool validateNetwork();
401 bool solveSystem();
402 void iterationSummary();
403 BundleSolutionInfo* bundleSolveInformation();
404 bool computeBundleStatistics();
405 void applyParameterCorrections();
406 bool errorPropagation();
407 void computeResiduals();
408 double computeVtpv();
409 bool computeRejectionLimit();
410 bool flagOutliers();
411
412 // normal equation matrices methods
413
414 bool formNormalEquations();
415 bool computePartials(LinearAlgebra::Matrix &coeffTarget,
416 LinearAlgebra::Matrix &coeffImage,
417 LinearAlgebra::Matrix &coeffPoint3D,
418 LinearAlgebra::Vector &coeffRHS,
419 BundleMeasure &measure,
420 BundleControlPoint &point);
421 bool formMeasureNormals(LinearAlgebra::MatrixUpperTriangular &N22,
425 LinearAlgebra::Matrix &coeffTarget,
426 LinearAlgebra::Matrix &coeffImage,
427 LinearAlgebra::Matrix &coeffPoint3D,
428 LinearAlgebra::Vector &coeffRHS,
429 int observationIndex);
430 int formPointNormals(LinearAlgebra::MatrixUpperTriangular &N22,
434 BundleControlPointQsp &point);
435 int formLidarPointNormals(LinearAlgebra::MatrixUpperTriangular &N22,
440 bool formWeightedNormals(LinearAlgebra::VectorCompressed &n1,
442
443 // dedicated matrix functions
444
445 void productAB(SparseBlockColumnMatrix &A,
447 void accumProductAlphaAB(double alpha,
451 bool invert3x3(LinearAlgebra::MatrixUpperTriangular &m);
452 bool productATransB(LinearAlgebra::MatrixUpperTriangular &N22,
455
456 // CHOLMOD library methods
457
458 bool initializeCHOLMODLibraryVariables();
459 bool freeCHOLMODLibraryVariables();
460 bool loadCholmodTriplet();
461
462 // member variables
463
464 BundleSettingsQsp m_bundleSettings;
465 BundleResults m_bundleResults;
467 ControlNetQsp m_controlNet;
468 QString m_cnetFileName;
469
470 QVector <BundleControlPointQsp> m_bundleControlPoints;
471 BundleLidarPointVector m_bundleLidarControlPoints;
472
473 QString m_lidarFileName;
474 LidarDataQsp m_lidarDataSet;
475 int m_numLidarConstraints;
476
477 BundleObservationVector m_bundleObservations;
480 SerialNumberList *m_serialNumberList;
481 BundleTargetBodyQsp m_bundleTargetBody;
483 bool m_abort;
484 QString m_iterationSummary;
486 bool m_printSummary;
489 bool m_cleanUp;
492 int m_rank;
493 int m_iteration;
494 double m_iterationTime;
495 int m_numberOfImagePartials;
496 QList<ImageList *> m_imageLists;
499 // ==========================================================================================
500 // === BEYOND THIS PLACE (THERE BE DRAGONS) all refers to the folded bundle solution. ===
501 // === Currently, everything uses the CHOLMOD library, ===
502 // === there is no dependence on the least-squares class. ===
503 // ==========================================================================================
504
506 boost::numeric::ublas::symmetric_matrix<
507 double,
508 boost::numeric::ublas::upper,
509 boost::numeric::ublas::column_major > m_normalInverse;
510 cholmod_common m_cholmodCommon;
515 SparseBlockMatrix m_sparseNormals;
519 cholmod_triplet *m_cholmodTriplet;
525 cholmod_sparse *m_cholmodNormal;
530 cholmod_factor *m_L;
534 LinearAlgebra::Vector m_imageSolution;
537 int m_previousNumberImagePartials;
540 };
541}
542
543#endif
An image bundle adjustment object.
Definition BundleAdjust.h:329
QString modelState(int index)
Return the updated model state for the ith cube in the cube list given to the constructor.
Definition BundleAdjust.cpp:2974
int numberOfImages() const
Returns the number of images.
Definition BundleAdjust.cpp:2907
void statusBarUpdate(QString)
bool solveCholesky()
Compute the least squares bundle adjustment solution using Cholesky decomposition.
Definition BundleAdjust.cpp:831
void statusUpdate(QString)
SerialNumberList * serialNumberList()
Returns a pointer to the serial number list.
Definition BundleAdjust.cpp:2897
bool isConverged()
Returns if the BundleAdjust converged.
Definition BundleAdjust.cpp:3068
QString fileName(int index)
Return the ith filename in the cube list file given to constructor.
Definition BundleAdjust.cpp:2921
void abortBundle()
Flag to abort when bundle is threaded.
Definition BundleAdjust.cpp:812
void iterationUpdate(int)
double iteration() const
Returns what iteration the BundleAdjust is currently on.
Definition BundleAdjust.cpp:2931
void outputBundleStatus(QString status)
Slot for deltack and jigsaw to output the bundle status.
Definition BundleAdjust.cpp:3104
void pointUpdate(int)
LidarDataQsp lidarData()
Returns a pointer to the output lidar data file.
Definition BundleAdjust.cpp:2887
bool isAborted()
Returns if the BundleAdjust has been aborted.
Definition BundleAdjust.cpp:3078
Table cMatrix(int index)
Return the updated instrument pointing table for the ith cube in the cube list given to the construct...
Definition BundleAdjust.cpp:2946
QList< ImageList * > imageLists()
This method returns the image list used in the bundle adjust.
Definition BundleAdjust.cpp:2480
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.
Definition BundleAdjust.cpp:107
QString iterationSummaryGroup() const
Returns the iteration summary string.
Definition BundleAdjust.cpp:3090
void resultsReady(BundleSolutionInfo *bundleSolveInformation)
void error(QString)
BundleSolutionInfo * solveCholeskyBR()
Compute the least squares bundle adjustment solution using Cholesky decomposition.
Definition BundleAdjust.cpp:802
ControlNetQsp controlNet()
Returns a pointer to the output control network.
Definition BundleAdjust.cpp:2877
~BundleAdjust()
Destroys BundleAdjust object, deallocates pointers (if we have ownership), and frees variables from c...
Definition BundleAdjust.cpp:383
Table spVector(int index)
Return the updated instrument position table for the ith cube in the cube list given to the construct...
Definition BundleAdjust.cpp:2961
This class holds information about a control point that BundleAdjust needs to run correctly.
Definition BundleControlPoint.h:91
This class is a container class for BundleLidarControlPoints.
Definition BundleLidarPointVector.h:31
A container class for a ControlMeasure.
Definition BundleMeasure.h:55
This class is a container class for BundleObservations.
Definition BundleObservationVector.h:57
A container class for statistical results from a BundleAdjust solution.
Definition BundleResults.h:91
Container class for BundleAdjustment results.
Definition BundleSolutionInfo.h:162
This represents an ISIS control net in a project-based GUI interface.
Definition Control.h:66
a control network
Definition ControlNet.h:258
Internalizes a list of images and allows for operations on the entire list.
Definition ImageList.h:55
boost::numeric::ublas::compressed_vector< double > VectorCompressed
Definition for an Isis::LinearAlgebra::VectorCompressed of doubles.
Definition LinearAlgebra.h:142
boost::numeric::ublas::vector< double > Vector
Definition for an Isis::LinearAlgebra::Vector of doubles.
Definition LinearAlgebra.h:132
boost::numeric::ublas::matrix< double > Matrix
Definition for an Isis::LinearAlgebra::Matrix of doubles.
Definition LinearAlgebra.h:102
boost::numeric::ublas::symmetric_matrix< double, boost::numeric::ublas::upper > MatrixUpperTriangular
Definition for an Isis::LinearAlgebra::MatrixUpperTriangular of doubles with an upper configuration.
Definition LinearAlgebra.h:122
Program progress reporter.
Definition Progress.h:42
Serial Number list generator.
Definition SerialNumberList.h:64
SparseBlockColumnMatrix.
Definition SparseBlockMatrix.h:58
SparseBlockMatrix.
Definition SparseBlockMatrix.h:186
SparseBlockRowMatrix.
Definition SparseBlockMatrix.h:125
Class for storing Table blobs information.
Definition Table.h:61
This is free and unencumbered software released into the public domain.
Definition BoxcarCachingAlgorithm.h:13
This is free and unencumbered software released into the public domain.
Definition CubeIoHandler.h:22
This is free and unencumbered software released into the public domain.
Definition Apollo.h:16