File failed to load: https://isis.astrogeology.usgs.gov/9.0.0/Object/assets/jax/output/NativeMML/config.js
Isis Developer Reference
BundleAdjust.h
Go to the documentation of this file.
1#ifndef BundleAdjust_h
2#define BundleAdjust_h
3
9
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
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 const QString &cnetFile,
342 const QString &cubeList,
343 const QString &lidarDataFile,
344 bool printSummary = true);
345 BundleAdjust(BundleSettingsQsp bundleSettings,
346 QString &cnet,
347 SerialNumberList &snlist,
348 bool printSummary = true);
349 BundleAdjust(BundleSettingsQsp bundleSettings,
350 Control &cnet,
351 SerialNumberList &snlist,
352 bool bPrintSummary);
353 BundleAdjust(BundleSettingsQsp bundleSettings,
354 ControlNet &cnet,
355 SerialNumberList &snlist,
356 bool printSummary = true);
357 BundleAdjust(BundleSettingsQsp bundleSettings,
358 ControlNetQsp cnet,
359 const QString &cubeList,
360 bool printSummary = true);
361 BundleAdjust(BundleSettingsQsp bundleSettings,
362 Control &control,
363 QList<ImageList *> imgList,
364 bool printSummary);
367
369 bool isAborted();
370
371 public slots:
372 bool solveCholesky();
373 void abortBundle();
374 void outputBundleStatus(QString status);
375
376 // accessors
377
381 QString fileName(int index);
382 QString iterationSummaryGroup() const;
383 bool isConverged();
384 Table cMatrix(int index);
385 Table spVector(int index);
386 QString modelState(int index);
387 int numberOfImages() const;
388 double iteration() const;
389
390 signals:
391 void statusUpdate(QString);
392 void error(QString);
394 void pointUpdate(int);
395 void statusBarUpdate(QString);
396 void resultsReady(BundleSolutionInfo *bundleSolveInformation);
397 void finished();
398
399 private:
400 //TODO Should there be a resetBundle(BundleSettings bundleSettings) method
401 // that allows for rerunning with new settings? JWB
402 void init(Progress *progress = 0);
403 bool initializeNormalEquationsMatrix();
404 bool validateNetwork();
405 bool solveSystem();
406 void iterationSummary();
407 BundleSolutionInfo* bundleSolveInformation();
408 bool computeBundleStatistics();
409 void applyParameterCorrections();
410 bool errorPropagation();
411 void computeResiduals();
412 double computeVtpv();
413 bool computeRejectionLimit();
414 bool flagOutliers();
415
416 // normal equation matrices methods
417
418 bool formNormalEquations();
419 bool computePartials(LinearAlgebra::Matrix &coeffTarget,
420 LinearAlgebra::Matrix &coeffImage,
421 LinearAlgebra::Matrix &coeffPoint3D,
422 LinearAlgebra::Vector &coeffRHS,
423 BundleMeasure &measure,
424 BundleControlPoint &point);
425 bool formMeasureNormals(LinearAlgebra::MatrixUpperTriangular &N22,
429 LinearAlgebra::Matrix &coeffTarget,
430 LinearAlgebra::Matrix &coeffImage,
431 LinearAlgebra::Matrix &coeffPoint3D,
432 LinearAlgebra::Vector &coeffRHS,
433 int observationIndex);
434 int formPointNormals(LinearAlgebra::MatrixUpperTriangular &N22,
438 BundleControlPointQsp &point);
439 int formLidarPointNormals(LinearAlgebra::MatrixUpperTriangular &N22,
444 bool formWeightedNormals(LinearAlgebra::VectorCompressed &n1,
446
447 // dedicated matrix functions
448
449 void productAB(SparseBlockColumnMatrix &A,
451 void accumProductAlphaAB(double alpha,
455 bool invert3x3(LinearAlgebra::MatrixUpperTriangular &m);
456 bool productATransB(LinearAlgebra::MatrixUpperTriangular &N22,
459
460 // CHOLMOD library methods
461
462 bool initializeCHOLMODLibraryVariables();
463 bool freeCHOLMODLibraryVariables();
464 bool loadCholmodTriplet();
465
466 // member variables
467
468 BundleSettingsQsp m_bundleSettings;
469 BundleResults m_bundleResults;
471 ControlNetQsp m_controlNet;
472 QString m_cnetFileName;
473
474 QVector <BundleControlPointQsp> m_bundleControlPoints;
475 BundleLidarPointVector m_bundleLidarControlPoints;
476
477 QString m_lidarFileName;
478 LidarDataQsp m_lidarDataSet;
479 int m_numLidarConstraints;
480
481 BundleObservationVector m_bundleObservations;
484 SerialNumberList *m_serialNumberList;
485 BundleTargetBodyQsp m_bundleTargetBody;
487 bool m_abort;
488 QString m_iterationSummary;
490 bool m_printSummary;
493 bool m_cleanUp;
496 int m_rank;
497 int m_iteration;
498 double m_iterationTime;
499 int m_numberOfImagePartials;
500 QList<ImageList *> m_imageLists;
502
503 // ==========================================================================================
504 // === BEYOND THIS PLACE (THERE BE DRAGONS) all refers to the folded bundle solution. ===
505 // === Currently, everything uses the CHOLMOD library, ===
506 // === there is no dependence on the least-squares class. ===
507 // ==========================================================================================
508
510 boost::numeric::ublas::symmetric_matrix<
511 double,
512 boost::numeric::ublas::upper,
513 boost::numeric::ublas::column_major > m_normalInverse;
514 cholmod_common m_cholmodCommon;
519 SparseBlockMatrix m_sparseNormals;
523 cholmod_triplet *m_cholmodTriplet;
529 cholmod_sparse *m_cholmodNormal;
534 cholmod_factor *m_L;
538 LinearAlgebra::Vector m_imageSolution;
540
541 int m_previousNumberImagePartials;
544 };
545}
546
547#endif
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:838
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:819
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:809
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:90
Container class for BundleAdjustment results.
Definition BundleSolutionInfo.h:172
This represents an ISIS control net in a project-based GUI interface.
Definition Control.h:65
a control network
Definition ControlNet.h:258
Internalizes a list of images and allows for operations on the entire list.
Definition ImageList.h:53
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
QSharedPointer< LidarData > LidarDataQsp
Definition for a shared pointer to a LidarData object.
Definition LidarData.h:100
QSharedPointer< BundleLidarControlPoint > BundleLidarControlPointQsp
QSharedPointer to a BundleLidarControlPoint.
Definition BundleLidarControlPoint.h:74
QSharedPointer< BundleSettings > BundleSettingsQsp
Definition for a BundleSettingsQsp, a shared pointer to a BundleSettings object.
Definition BundleSettings.h:355
QSharedPointer< ControlNet > ControlNetQsp
Typedef for QSharedPointer to control network. This typedef is for future implementation of target bo...
Definition ControlNet.h:485
QSharedPointer< BundleTargetBody > BundleTargetBodyQsp
Definition for BundleTargetBodyQsp, a QSharedPointer to a BundleTargetBody.
Definition BundleTargetBody.h:187
QSharedPointer< BundleControlPoint > BundleControlPointQsp
Definition for BundleControlPointQSP, a shared pointer to a BundleControlPoint.
Definition BundleControlPoint.h:192