File failed to load: https://isis.astrogeology.usgs.gov/9.0.0/Object/assets/jax/output/NativeMML/config.js
Isis 3 Programmer Reference
BundleAdjust.h
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"
23#include "BundleLidarPointVector.h"
24#include "BundleObservationSolveSettings.h"
25#include "BundleObservationVector.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"
36#include "MaximumLikelihoodWFunctions.h"
37#include "ObservationNumberList.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);
393 void iterationUpdate(int);
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);
404 bool validateNetwork();
405 bool solveSystem();
406 void iterationSummary();
410 bool errorPropagation();
411 void computeResiduals();
412 double computeVtpv();
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);
429 LinearAlgebra::Matrix &coeffTarget,
430 LinearAlgebra::Matrix &coeffImage,
431 LinearAlgebra::Matrix &coeffPoint3D,
432 LinearAlgebra::Vector &coeffRHS,
433 int observationIndex);
438 BundleControlPointQsp &point);
446
447 // dedicated matrix functions
448
451 void accumProductAlphaAB(double alpha,
459
460 // CHOLMOD library methods
461
464 bool loadCholmodTriplet();
465
466 // member variables
467
473
474 QVector <BundleControlPointQsp> m_bundleControlPoints;
476
480
481 BundleObservationVector m_bundleObservations;
485 BundleTargetBodyQsp m_bundleTargetBody;
487 bool m_abort;
488 QString m_iterationSummary;
496 int m_rank;
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;
523 cholmod_triplet *m_cholmodTriplet;
529 cholmod_sparse *m_cholmodNormal;
534 cholmod_factor *m_L;
540
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.
bool m_cleanUp
!< If the iteration summaries should be output to the log file.
int m_rank
!< If the serial number lists need to be deleted by the destructor.
SparseBlockMatrix m_sparseNormals
!< The right hand side of the normal equations.
int numberOfImages() const
Returns the number of images.
void productAB(SparseBlockColumnMatrix &A, SparseBlockRowMatrix &B)
Perform the matrix multiplication C = N12 x Q.
bool solveCholesky()
Compute the least squares bundle adjustment solution using Cholesky decomposition.
bool computeBundleStatistics()
Compute Bundle statistics and store them in m_bundleResults.
double m_iterationTime
Time for last iteration.
int formLidarPointNormals(LinearAlgebra::MatrixUpperTriangular &N22, SparseBlockColumnMatrix &N12, LinearAlgebra::Vector &n2, LinearAlgebra::Vector &nj, BundleLidarControlPointQsp &point)
Compute the Q matrix and NIC vector for a control point.
bool formNormalEquations()
Form the least-squares normal equations matrix via cholmod.
SerialNumberList * serialNumberList()
Returns a pointer to the serial number list.
void init(Progress *progress=0)
Initialize all solution parameters.
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 isConverged()
Returns if the BundleAdjust converged.
QString fileName(int index)
Return the ith filename in the cube list file given to constructor.
void abortBundle()
Flag to abort when bundle is threaded.
QString m_lidarFileName
Input lidar point filename.
int m_previousNumberImagePartials
!< The image parameter solution vector.
double iteration() const
Returns what iteration the BundleAdjust is currently on.
void outputBundleStatus(QString status)
Slot for deltack and jigsaw to output the bundle status.
bool flagOutliers()
Flags outlier measures and control points.
QVector< BundleControlPointQsp > m_bundleControlPoints
Vector of control points.
LidarDataQsp lidarData()
Returns a pointer to the output lidar data file.
bool isAborted()
Returns if the BundleAdjust has been aborted.
BundleSettingsQsp m_bundleSettings
Contains the solve settings.
bool invert3x3(LinearAlgebra::MatrixUpperTriangular &m)
Dedicated quick inverse of 3x3 matrix.
int m_numLidarConstraints
TODO: temp.
Table cMatrix(int index)
Return the updated instrument pointing table for the ith cube in the cube list given to the construct...
double computeVtpv()
Computes vtpv, the weighted sum of squares of residuals.
QList< ImageList * > imageLists()
This method returns the image list used in the bundle adjust.
void applyParameterCorrections()
Apply parameter corrections for current iteration.
bool initializeCHOLMODLibraryVariables()
Initializations for CHOLMOD sparse matrix package.
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.
int m_numberOfImagePartials
number of image-related partials.
LinearAlgebra::Vector m_imageSolution
!< The lower triangular L matrix from Cholesky decomposition.
cholmod_triplet * m_cholmodTriplet
!< The sparse block normal equations matrix.
ControlNetQsp m_controlNet
Output control net.
bool validateNetwork()
control network validation - on the very real chance that the net has not been checked before running...
SerialNumberList * m_serialNumberList
!< Vector of observations.
bool loadCholmodTriplet()
Load sparse normal equations matrix into CHOLMOD triplet.
LidarDataQsp m_lidarDataSet
Output lidar data.
cholmod_factor * m_L
!< The CHOLMOD sparse normal equations matrix used by cholmod_factorize to solve the system.
bool solveSystem()
Compute the solution to the normal equations using the CHOLMOD library.
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.
void accumProductAlphaAB(double alpha, SparseBlockRowMatrix &A, LinearAlgebra::Vector &B, LinearAlgebra::Vector &C)
Performs the matrix multiplication nj = nj + alpha (Q x n2).
bool initializeNormalEquationsMatrix()
Initialize Normal Equations matrix (m_sparseNormals).
QString iterationSummaryGroup() const
Returns the iteration summary string.
LinearAlgebra::Vector m_RHS
!< Contains object parameters, statistics, and workspace used by the CHOLMOD library.
bool errorPropagation()
Error propagation for solution.
BundleResults m_bundleResults
Stores the results of the bundle adjust.
QString m_cnetFileName
The control net filename.
bool formMeasureNormals(LinearAlgebra::MatrixUpperTriangular &N22, SparseBlockColumnMatrix &N12, LinearAlgebra::VectorCompressed &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.
BundleSolutionInfo * bundleSolveInformation()
Create a BundleSolutionInfo containing the settings and results from the bundle adjustment.
int formPointNormals(LinearAlgebra::MatrixUpperTriangular &N22, SparseBlockColumnMatrix &N12, LinearAlgebra::Vector &n2, LinearAlgebra::Vector &nj, BundleControlPointQsp &point)
Compute the Q matrix and NIC vector for a control point.
bool m_abort
!< Contains information about the target body.
BundleLidarPointVector m_bundleLidarControlPoints
Vector of lidar points.
void computeResiduals()
Compute image measure residuals.
bool m_printSummary
!< Summary of the most recently completed iteration.
bool computeRejectionLimit()
Compute rejection limit.
bool productATransB(LinearAlgebra::MatrixUpperTriangular &N22, SparseBlockColumnMatrix &N12, SparseBlockRowMatrix &Q)
Perform the matrix multiplication Q = N22 x N12(transpose)
void iterationSummary()
Creates an iteration summary and an iteration group for the solution summary.
int m_iteration
The current iteration.
bool freeCHOLMODLibraryVariables()
Free CHOLMOD library variables.
bool formWeightedNormals(LinearAlgebra::VectorCompressed &n1, LinearAlgebra::Vector &nj)
Apply weighting for spacecraft position, velocity, acceleration and camera angles,...
BundleSolutionInfo * solveCholeskyBR()
Compute the least squares bundle adjustment solution using Cholesky decomposition.
ControlNetQsp controlNet()
Returns a pointer to the output control network.
cholmod_sparse * m_cholmodNormal
!< The CHOLMOD triplet representation of the sparse normal equations matrix.
~BundleAdjust()
Destroys BundleAdjust object, deallocates pointers (if we have ownership), and frees variables from c...
Table spVector(int index)
Return the updated instrument position table for the ith cube in the cube list given to the construct...
This class holds information about a control point that BundleAdjust needs to run correctly.
This class is a container class for BundleLidarControlPoints.
A container class for a ControlMeasure.
This class is a container class for BundleObservations.
A container class for statistical results from a BundleAdjust solution.
Container class for BundleAdjustment results.
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.
boost::numeric::ublas::vector< double > Vector
Definition for an Isis::LinearAlgebra::Vector of doubles.
boost::numeric::ublas::matrix< double > Matrix
Definition for an Isis::LinearAlgebra::Matrix of doubles.
boost::numeric::ublas::symmetric_matrix< double, boost::numeric::ublas::upper > MatrixUpperTriangular
Definition for an Isis::LinearAlgebra::MatrixUpperTriangular of doubles with an upper configuration.
Program progress reporter.
Definition Progress.h:42
Serial Number list generator.
SparseBlockColumnMatrix.
Class for storing Table blobs information.
Definition Table.h:61
This is free and unencumbered software released into the public domain.
This is free and unencumbered software released into the public domain.
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.
QSharedPointer< BundleSettings > BundleSettingsQsp
Definition for a BundleSettingsQsp, a shared pointer to a BundleSettings object.
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.
QSharedPointer< BundleControlPoint > BundleControlPointQsp
Definition for BundleControlPointQSP, a shared pointer to a BundleControlPoint.