Isis 3 Programmer Reference
BundleAdjust.cpp
1 
7 /* SPDX-License-Identifier: CC0-1.0 */
8 
9 #include "BundleAdjust.h"
10 
11 // std lib
12 #include <iomanip>
13 #include <iostream>
14 #include <sstream>
15 
16 // qt lib
17 #include <QCoreApplication>
18 #include <QDebug>
19 #include <QFile>
20 #include <QMutex>
21 
22 // boost lib
23 #include <boost/lexical_cast.hpp>
24 #include <boost/numeric/ublas/io.hpp>
25 #include <boost/numeric/ublas/matrix_sparse.hpp>
26 #include <boost/numeric/ublas/vector_proxy.hpp>
27 
28 // Isis lib
29 #include "Application.h"
30 #include "BundleObservation.h"
31 #include "IsisBundleObservation.h"
32 #include "BundleObservationSolveSettings.h"
33 #include "BundleResults.h"
34 #include "BundleSettings.h"
35 #include "BundleSolutionInfo.h"
36 #include "BundleTargetBody.h"
37 #include "Camera.h"
38 #include "CameraDetectorMap.h"
39 #include "CameraDistortionMap.h"
40 #include "CameraFocalPlaneMap.h"
41 #include "CameraGroundMap.h"
42 #include "CSMCamera.h"
43 #include "Control.h"
44 #include "ControlPoint.h"
45 #include "CorrelationMatrix.h"
46 #include "Distance.h"
47 #include "ImageList.h"
48 #include "iTime.h"
49 #include "Latitude.h"
50 #include "Longitude.h"
51 #include "MaximumLikelihoodWFunctions.h"
52 #include "SpecialPixel.h"
53 #include "StatCumProbDistDynCalc.h"
54 #include "SurfacePoint.h"
55 #include "Target.h"
56 
57 using namespace boost::numeric::ublas;
58 using namespace Isis;
59 
60 namespace Isis {
61 
62 
72  static void cholmodErrorHandler(int nStatus,
73  const char* file,
74  int nLineNo,
75  const char* message) {
76  QString errlog;
77 
78  errlog = "SPARSE: ";
79  errlog += message;
80 
81  PvlGroup gp(errlog);
82 
83  gp += PvlKeyword("File",file);
84  gp += PvlKeyword("Line_Number", toString(nLineNo));
85  gp += PvlKeyword("Status", toString(nStatus));
86 
87 // Application::Log(gp);
88 
89  errlog += ". (See print.prt for details)";
90 
91 // throw IException(IException::Unknown, errlog, _FILEINFO_);
92  }
93 
94 
104  BundleAdjust::BundleAdjust(BundleSettingsQsp bundleSettings,
105  const QString &cnetFile,
106  const QString &cubeList,
107  bool printSummary) {
108  m_abort = false;
109  Progress progress;
110  // initialize constructor dependent settings...
111  // m_printSummary, m_cleanUp, m_cnetFileName, m_controlNet,
112  // m_serialNumberList, m_bundleSettings
113  m_printSummary = printSummary;
114  m_cleanUp = true;
115  m_cnetFileName = cnetFile;
116  try {
117  m_controlNet = ControlNetQsp( new ControlNet(cnetFile, &progress,
118  bundleSettings->controlPointCoordTypeReports()) );
119  }
120  catch (IException &e) {
121  throw;
122  }
123  m_bundleResults.setOutputControlNet(m_controlNet);
124  m_serialNumberList = new SerialNumberList(cubeList);
125  m_bundleSettings = bundleSettings;
126  m_bundleTargetBody = bundleSettings->bundleTargetBody();
127 
128  init(&progress);
129  }
130 
131 
142  BundleAdjust::BundleAdjust(BundleSettingsQsp bundleSettings,
143  QString &cnetFile,
144  SerialNumberList &snlist,
145  bool printSummary) {
146  // initialize constructor dependent settings...
147  // m_printSummary, m_cleanUp, m_cnetFileName, m_controlNet,
148  // m_serialNumberList, m_bundleSettings
149  m_abort = false;
150  Progress progress;
151  m_printSummary = printSummary;
152  m_cleanUp = false;
153  m_cnetFileName = cnetFile;
154  try {
155  m_controlNet = ControlNetQsp( new ControlNet(cnetFile, &progress) );
156  }
157  catch (IException &e) {
158  throw;
159  }
160  m_bundleResults.setOutputControlNet(m_controlNet);
161  m_serialNumberList = &snlist;
162  m_bundleSettings = bundleSettings;
163  m_bundleTargetBody = bundleSettings->bundleTargetBody();
164 
165  init();
166  }
167 
168 
179  BundleAdjust::BundleAdjust(BundleSettingsQsp bundleSettings,
180  Control &cnet,
181  SerialNumberList &snlist,
182  bool printSummary) {
183  // initialize constructor dependent settings...
184  // m_printSummary, m_cleanUp, m_cnetFileName, m_controlNet,
185  // m_serialNumberList, m_bundleSettings
186  m_abort = false;
187  Progress progress;
188  m_printSummary = printSummary;
189  m_cleanUp = false;
190  m_cnetFileName = cnet.fileName();
191  try {
192  m_controlNet = ControlNetQsp( new ControlNet(cnet.fileName(), &progress) );
193  }
194  catch (IException &e) {
195  throw;
196  }
197  m_bundleResults.setOutputControlNet(m_controlNet);
198  m_serialNumberList = &snlist;
199  m_bundleSettings = bundleSettings;
200  m_bundleTargetBody = bundleSettings->bundleTargetBody();
201 
202  init();
203  }
204 
205 
216  BundleAdjust::BundleAdjust(BundleSettingsQsp bundleSettings,
217  ControlNet &cnet,
218  SerialNumberList &snlist,
219  bool printSummary) {
220  // initialize constructor dependent settings...
221  // m_printSummary, m_cleanUp, m_cnetFileName, m_controlNet,
222  // m_serialNumberList, m_bundleSettings
223  m_abort = false;
224  m_printSummary = printSummary;
225  m_cleanUp = false;
226  m_cnetFileName = "";
227  try {
228  m_controlNet = ControlNetQsp( new ControlNet(cnet) );
229  }
230  catch (IException &e) {
231  throw;
232  }
233  m_bundleResults.setOutputControlNet(m_controlNet);
234  m_serialNumberList = &snlist;
235  m_bundleSettings = bundleSettings;
236  m_bundleTargetBody = bundleSettings->bundleTargetBody();
237 
238  init();
239  }
240 
241 
250  BundleAdjust::BundleAdjust(BundleSettingsQsp bundleSettings,
251  ControlNetQsp cnet,
252  const QString &cubeList,
253  bool printSummary) {
254  m_abort = false;
255  m_printSummary = printSummary;
256  m_cleanUp = false;
257  m_cnetFileName = "";
258  try {
259  m_controlNet = cnet;
260  }
261  catch (IException &e) {
262  throw;
263  }
264  m_bundleResults.setOutputControlNet(m_controlNet);
265  m_serialNumberList = new SerialNumberList(cubeList);
266  m_bundleSettings = bundleSettings;
267  m_bundleTargetBody = bundleSettings->bundleTargetBody();
268 
269  init();
270  }
271 
272 
282  BundleAdjust::BundleAdjust(BundleSettingsQsp bundleSettings,
283  Control &control,
284  QList<ImageList *> imgLists,
285  bool printSummary) {
286  m_bundleSettings = bundleSettings;
287 
288  m_abort = false;
289  try {
290  m_controlNet = ControlNetQsp( new ControlNet(control.fileName()) );
291  }
292  catch (IException &e) {
293  throw;
294  }
295  m_bundleResults.setOutputControlNet(m_controlNet);
296 
297  m_imageLists = imgLists;
298 
299  // this is too slow and we need to get rid of the serial number list anyway
300  // should be unnecessary as Image class has serial number
301  // could hang on to image list until creating BundleObservations?
302  m_serialNumberList = new SerialNumberList;
303 
304  foreach (ImageList *imgList, imgLists) {
305  foreach (Image *image, *imgList) {
306  m_serialNumberList->add(image->fileName());
307 // m_serialNumberList->add(image->serialNumber(), image->fileName());
308  }
309  }
310 
311  m_bundleTargetBody = bundleSettings->bundleTargetBody();
312 
313  m_printSummary = printSummary;
314 
315  m_cleanUp = false;
316  m_cnetFileName = control.fileName();
317 
318  init();
319  }
320 
321 
330  BundleAdjust::~BundleAdjust() {
331  // If we have ownership
332  if (m_cleanUp) {
333  delete m_serialNumberList;
334  }
335 
336  freeCHOLMODLibraryVariables();
337 
338  }
339 
340 
372  void BundleAdjust::init(Progress *progress) {
373  emit(statusUpdate("Initialization"));
374  m_previousNumberImagePartials = 0;
375 
376  // initialize
377  //
378  // JWB
379  // - some of these not originally initialized.. better values???
380  m_iteration = 0;
381  m_rank = 0;
382  m_iterationSummary = "";
383 
384  // Get the cameras set up for all images
385  // NOTE - THIS IS NOT THE SAME AS "setImage" as called in BundleAdjust::computePartials
386  // this call only does initializations; sets measure's camera pointer, etc
387  // RENAME????????????
388  m_controlNet->SetImages(*m_serialNumberList, progress);
389 
390  // clear JigsawRejected flags
391  m_controlNet->ClearJigsawRejected();
392 
393  // initialize held variables
394  int numImages = m_serialNumberList->size();
395 
396  // matrix stuff
397  m_normalInverse.clear();
398  m_RHS.clear();
399  m_imageSolution.clear();
400 
401  // we don't want to call initializeCHOLMODLibraryVariables() here since mRank=0
402  // m_cholmodCommon, m_sparseNormals are not initialized
403  m_L = NULL;
404  m_cholmodNormal = NULL;
405  m_cholmodTriplet = NULL;
406 
407  // should we initialize objects m_xResiduals, m_yResiduals, m_xyResiduals
408 
409  // TESTING
410  // TODO: code below should go into a separate method???
411  // set up BundleObservations and assign solve settings for each from BundleSettings class
412  for (int i = 0; i < numImages; i++) {
413 
414  Camera *camera = m_controlNet->Camera(i);
415  QString observationNumber = m_serialNumberList->observationNumber(i);
416  QString instrumentId = m_serialNumberList->spacecraftInstrumentId(i);
417  QString serialNumber = m_serialNumberList->serialNumber(i);
418  QString fileName = m_serialNumberList->fileName(i);
419 
420  // create a new BundleImage and add to new (or existing if observation mode is on)
421  // BundleObservation
422  BundleImageQsp image = BundleImageQsp(new BundleImage(camera, serialNumber, fileName));
423 
424  if (!image) {
425  QString msg = "In BundleAdjust::init(): image " + fileName + "is null." + "\n";
426  throw IException(IException::Programmer, msg, _FILEINFO_);
427  }
428 
429  BundleObservationQsp observation =
430  m_bundleObservations.addNew(image, observationNumber, instrumentId, m_bundleSettings);
431 
432  if (!observation) {
433  QString msg = "In BundleAdjust::init(): observation "
434  + observationNumber + "is null." + "\n";
435  throw IException(IException::Programmer, msg, _FILEINFO_);
436  }
437  }
438 
439  // set up vector of BundleControlPoints
440  int numControlPoints = m_controlNet->GetNumPoints();
441  for (int i = 0; i < numControlPoints; i++) {
442  ControlPoint *point = m_controlNet->GetPoint(i);
443  if (point->IsIgnored()) {
444  continue;
445  }
446 
447  BundleControlPointQsp bundleControlPoint(new BundleControlPoint
448  (m_bundleSettings, point));
449  m_bundleControlPoints.append(bundleControlPoint);
450 
451  // set parent observation for each BundleMeasure
452 
453  int numMeasures = bundleControlPoint->size();
454  for (int j=0; j < numMeasures; j++) {
455  BundleMeasureQsp measure = bundleControlPoint->at(j);
456  QString cubeSerialNumber = measure->cubeSerialNumber();
457 
458  BundleObservationQsp observation =
459  m_bundleObservations.observationByCubeSerialNumber(cubeSerialNumber);
460  BundleImageQsp image = observation->imageByCubeSerialNumber(cubeSerialNumber);
461 
462  measure->setParentObservation(observation);
463  measure->setParentImage(image);
464  }
465  }
466 
467  //===========================================================================================//
468  //==== Use the bundle settings to initialize more member variables and set up solutions =====//
469  //===========================================================================================//
470 
471  // TODO: Need to have some validation code to make sure everything is
472  // on the up-and-up with the control network. Add checks for multiple
473  // networks, images without any points, and points on images removed from
474  // the control net (when we start adding software to remove points with high
475  // residuals) and ?. For "deltack" a single measure on a point is allowed
476  // so skip the test.
477  if (m_bundleSettings->validateNetwork()) {
478  validateNetwork();
479  }
480  m_bundleResults.maximumLikelihoodSetUp(m_bundleSettings->maximumLikelihoodEstimatorModels());
481 
482  //===========================================================================================//
483  //=============== End Bundle Settings =======================================================//
484  //===========================================================================================//
485 
486  //===========================================================================================//
487  //======================== initialize matrices and more parameters ==========================//
488  //===========================================================================================//
489 
490  // size of reduced normals matrix
491 
492  // TODO
493  // this should be determined from BundleSettings
494  // m_rank will be the sum of observation, target, and self-cal parameters
495  // TODO
496  m_rank = m_bundleObservations.numberParameters();
497 
498  if (m_bundleSettings->solveTargetBody()) {
499  m_rank += m_bundleSettings->numberTargetBodyParameters();
500 
501  if (m_bundleTargetBody->solveMeanRadius() || m_bundleTargetBody->solveTriaxialRadii()) {
502  outputBundleStatus("Warning: Solving for the target body radii (triaxial or mean) "
503  "is NOT possible and likely increases error in the solve.\n");
504  }
505 
506  if (m_bundleTargetBody->solveMeanRadius()){
507  // Check if MeanRadiusValue is abnormal compared to observation
508  bool isMeanRadiusValid = true;
509  double localRadius, aprioriRadius;
510 
511  // Arbitrary control point containing an observed localRadius
512  BundleControlPointQsp point = m_bundleControlPoints.at(0);
513  SurfacePoint surfacepoint = point->adjustedSurfacePoint();
514 
515  localRadius = surfacepoint.GetLocalRadius().meters();
516  aprioriRadius = m_bundleTargetBody->meanRadius().meters();
517 
518  // Ensure aprioriRadius is within some threshold
519  // Considered potentially inaccurate if it's off by atleast a factor of two
520  if (aprioriRadius >= 2 * localRadius || aprioriRadius <= localRadius / 2) {
521  isMeanRadiusValid = false;
522  }
523 
524  // Warn user for abnormal MeanRadiusValue
525  if (!isMeanRadiusValid) {
526  outputBundleStatus("Warning: User-entered MeanRadiusValue appears to be inaccurate. "
527  "This can cause a bundle failure.\n");
528  }
529  }
530  }
531 
532  int num3DPoints = m_bundleControlPoints.size();
533 
534  m_bundleResults.setNumberUnknownParameters(m_rank + 3 * num3DPoints);
535 
536  m_imageSolution.resize(m_rank);
537 
538  //~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
539  // initializations for cholmod
540  initializeCHOLMODLibraryVariables();
541 
542  // initialize normal equations matrix
543  initializeNormalEquationsMatrix();
544  }
545 
546 
566  bool BundleAdjust::validateNetwork() {
567 
568  outputBundleStatus("\nValidating network...");
569 
570  int imagesWithInsufficientMeasures = 0;
571  QString msg = "Images with one or less measures:\n";
572  int numObservations = m_bundleObservations.size();
573  for (int i = 0; i < numObservations; i++) {
574  int numImages = m_bundleObservations.at(i)->size();
575  for (int j = 0; j < numImages; j++) {
576  BundleImageQsp bundleImage = m_bundleObservations.at(i)->at(j);
577  int numMeasures = m_controlNet->
578  GetNumberOfValidMeasuresInImage(bundleImage->serialNumber());
579 
580  if (numMeasures > 1) {
581  continue;
582  }
583 
584  imagesWithInsufficientMeasures++;
585  msg += bundleImage->fileName() + ": " + toString(numMeasures) + "\n";
586  }
587  }
588 
589  if ( imagesWithInsufficientMeasures > 0 ) {
590  throw IException(IException::User, msg, _FILEINFO_);
591  }
592 
593  outputBundleStatus("\nValidation complete!...\n");
594 
595  return true;
596  }
597 
604  bool BundleAdjust::initializeCHOLMODLibraryVariables() {
605  if ( m_rank <= 0 ) {
606  return false;
607  }
608 
609  m_cholmodTriplet = NULL;
610 
611  cholmod_start(&m_cholmodCommon);
612 
613  // set user-defined cholmod error handler
614  m_cholmodCommon.error_handler = cholmodErrorHandler;
615 
616  // testing not using metis
617  m_cholmodCommon.nmethods = 1;
618  m_cholmodCommon.method[0].ordering = CHOLMOD_AMD;
619 
620  return true;
621  }
622 
623 
634  bool BundleAdjust::initializeNormalEquationsMatrix() {
635 
636  int nBlockColumns = m_bundleObservations.size();
637 
638  if (m_bundleSettings->solveTargetBody())
639  nBlockColumns += 1;
640 
641  m_sparseNormals.setNumberOfColumns(nBlockColumns);
642 
643  m_sparseNormals.at(0)->setStartColumn(0);
644 
645  int nParameters = 0;
646  if (m_bundleSettings->solveTargetBody()) {
647  nParameters += m_bundleSettings->numberTargetBodyParameters();
648  m_sparseNormals.at(1)->setStartColumn(nParameters);
649 
650  int observation = 0;
651  for (int i = 2; i < nBlockColumns; i++) {
652  nParameters += m_bundleObservations.at(observation)->numberParameters();
653  m_sparseNormals.at(i)->setStartColumn(nParameters);
654  observation++;
655  }
656  }
657  else {
658  for (int i = 0; i < nBlockColumns; i++) {
659  m_sparseNormals.at(i)->setStartColumn(nParameters);
660  nParameters += m_bundleObservations.at(i)->numberParameters();
661  }
662  }
663 
664  return true;
665  }
666 
667 
676  bool BundleAdjust::freeCHOLMODLibraryVariables() {
677 
678  cholmod_free_triplet(&m_cholmodTriplet, &m_cholmodCommon);
679  cholmod_free_sparse(&m_cholmodNormal, &m_cholmodCommon);
680  cholmod_free_factor(&m_L, &m_cholmodCommon);
681 
682  cholmod_finish(&m_cholmodCommon);
683 
684  return true;
685  }
686 
687 
697  BundleSolutionInfo* BundleAdjust::solveCholeskyBR() {
698  solveCholesky();
699  return bundleSolveInformation();
700  }
701 
702 
707  void BundleAdjust::abortBundle() {
708  m_abort = true;
709  }
710 
711 
726  bool BundleAdjust::solveCholesky() {
727  emit(statusBarUpdate("Solving"));
728  try {
729 
730  // throw error if a frame camera is included AND
731  // if m_bundleSettings->solveInstrumentPositionOverHermiteSpline()
732  // is set to true (can only use for line scan or radar)
733  // if (m_bundleSettings->solveInstrumentPositionOverHermiteSpline() == true) {
734  // int numImages = images();
735  // for (int i = 0; i < numImages; i++) {
736  // if (m_controlNet->Camera(i)->GetCameraType() == 0) {
737  // QString msg = "At least one sensor is a frame camera. "
738  // "Spacecraft Option OVERHERMITE is not valid for frame cameras\n";
739  // throw IException(IException::User, msg, _FILEINFO_);
740  // }
741  // }
742  // }
743 
744  // Compute the apriori coordinates for each nonheld point
745  m_controlNet->ComputeApriori(); // original location
746 
747  // ken testing - if solving for target mean radius, set point radius to current mean radius
748  // if solving for triaxial radii, set point radius to local radius
749  if (m_bundleTargetBody && m_bundleTargetBody->solveMeanRadius()) {
750  int numControlPoints = m_bundleControlPoints.size();
751  for (int i = 0; i < numControlPoints; i++) {
752  BundleControlPointQsp point = m_bundleControlPoints.at(i);
753  SurfacePoint surfacepoint = point->adjustedSurfacePoint();
754 
755  surfacepoint.ResetLocalRadius(m_bundleTargetBody->meanRadius());
756 
757  point->setAdjustedSurfacePoint(surfacepoint);
758  }
759  }
760 
761  // Only use target body solution options when using Latitudinal coordinates
762  if (m_bundleTargetBody && m_bundleTargetBody->solveTriaxialRadii()
763  && m_bundleSettings->controlPointCoordTypeBundle() == SurfacePoint::Latitudinal) {
764  int numControlPoints = m_bundleControlPoints.size();
765  for (int i = 0; i < numControlPoints; i++) {
766  BundleControlPointQsp point = m_bundleControlPoints.at(i);
767  SurfacePoint surfacepoint = point->adjustedSurfacePoint();
768 
769  Distance localRadius = m_bundleTargetBody->localRadius(surfacepoint.GetLatitude(),
770  surfacepoint.GetLongitude());
771  surfacepoint.ResetLocalRadius(localRadius);
772 
773  point->setAdjustedSurfacePoint(surfacepoint);
774  }
775  }
776 
777  // Beginning of iterations
778  m_iteration = 1;
779  double vtpv = 0.0;
780  double previousSigma0 = 0.0;
781 
782  // start the clock
783  clock_t solveStartClock = clock();
784 
785  for (;;) {
786 
787  emit iterationUpdate(m_iteration);
788 
789  // testing
790  if (m_abort) {
791  m_bundleResults.setConverged(false);
792  emit statusUpdate("\n aborting...");
793  emit finished();
794  return false;
795  }
796  // testing
797 
798  emit statusUpdate( QString("starting iteration %1\n").arg(m_iteration) );
799 
800  clock_t iterationStartClock = clock();
801 
802  // zero normals (after iteration 0)
803  if (m_iteration != 1) {
804  m_sparseNormals.zeroBlocks();
805  }
806 
807  // form normal equations -- computePartials is called in here.
808  if (!formNormalEquations()) {
809  m_bundleResults.setConverged(false);
810  break;
811  }
812 
813  // testing
814  if (m_abort) {
815  m_bundleResults.setConverged(false);
816  emit statusUpdate("\n aborting...");
817  emit finished();
818  return false;
819  }
820  // testing
821 
822  // solve the system
823  if (!solveSystem()) {
824  outputBundleStatus("\nsolve failed!");
825  m_bundleResults.setConverged(false);
826  break;
827  }
828 
829  // testing
830  if (m_abort) {
831  m_bundleResults.setConverged(false);
832  emit statusUpdate("\n aborting...");
833  emit finished();
834  return false;
835  }
836  // testing
837 
838  // apply parameter corrections
839  applyParameterCorrections();
840 
841  // testing
842  if (m_abort) {
843  m_bundleResults.setConverged(false);
844  emit statusUpdate("\n aborting...");
845  emit finished();
846  return false;
847  }
848  // testing
849 
850  // compute residuals
851  vtpv = computeResiduals();
852 
853  // flag outliers
854  if ( m_bundleSettings->outlierRejection() ) {
855  computeRejectionLimit();
856  flagOutliers();
857  }
858 
859  // testing
860  if (m_abort) {
861  m_bundleResults.setConverged(false);
862  emit statusUpdate("\n aborting...");
863  emit finished();
864  return false;
865  }
866  // testing
867 
868  // variance of unit weight (also reference variance, variance factor, etc.)
869  m_bundleResults.computeSigma0(vtpv, m_bundleSettings->convergenceCriteria());
870 
871  // Set up formatting for status updates with doubles (e.g. Sigma0, Elapsed Time)
872  int fieldWidth = 20;
873  char format = 'f';
874  int precision = 10;
875 
876  emit statusUpdate(QString("Iteration: %1 \n")
877  .arg(m_iteration));
878  emit statusUpdate(QString("Sigma0: %1 \n")
879  .arg(m_bundleResults.sigma0(),
880  fieldWidth,
881  format,
882  precision));
883  emit statusUpdate(QString("Observations: %1 \n")
884  .arg(m_bundleResults.numberObservations()));
885  emit statusUpdate(QString("Constrained Parameters:%1 \n")
886  .arg(m_bundleResults.numberConstrainedPointParameters()));
887  emit statusUpdate(QString("Unknowns: %1 \n")
888  .arg(m_bundleResults.numberUnknownParameters()));
889  emit statusUpdate(QString("Degrees of Freedom: %1 \n")
890  .arg(m_bundleResults.degreesOfFreedom()));
891  emit iterationUpdate(m_iteration);
892 
893  // check for convergence
894  if (m_bundleSettings->convergenceCriteria() == BundleSettings::Sigma0) {
895  if (fabs(previousSigma0 - m_bundleResults.sigma0())
896  <= m_bundleSettings->convergenceCriteriaThreshold()) {
897  // convergence detected
898 
899  // if maximum likelihood tiers are being processed,
900  // check to see if there's another tier to go.
901  if (m_bundleResults.maximumLikelihoodModelIndex()
902  < m_bundleResults.numberMaximumLikelihoodModels() - 1
903  && m_bundleResults.maximumLikelihoodModelIndex()
904  < 2) {
905  // TODO is this second condition redundant???
906  // should BundleResults require num models <= 3, so num models - 1 <= 2
907  if (m_bundleResults.numberMaximumLikelihoodModels()
908  > m_bundleResults.maximumLikelihoodModelIndex() + 1) {
909 
910  // If there is another tier left we will increment the index.
911  m_bundleResults.incrementMaximumLikelihoodModelIndex();
912  }
913  }
914  else { // otherwise iterations are complete
915  m_bundleResults.setConverged(true);
916  emit statusUpdate("Bundle has converged\n");
917  emit statusBarUpdate("Converged");
918  break;
919  }
920  }
921  }
922  else {
923  // bundleSettings.convergenceCriteria() == BundleSettings::ParameterCorrections
924  int numConvergedParams = 0;
925  int numImgParams = m_imageSolution.size();
926  for (int ij = 0; ij < numImgParams; ij++) {
927  if (fabs(m_imageSolution(ij)) > m_bundleSettings->convergenceCriteriaThreshold()) {
928  break;
929  }
930  else
931  numConvergedParams++;
932  }
933 
934  if ( numConvergedParams == numImgParams ) {
935  m_bundleResults.setConverged(true);
936  emit statusUpdate("Bundle has converged\n");
937  emit statusBarUpdate("Converged");
938  break;
939  }
940  }
941 
942  m_bundleResults.printMaximumLikelihoodTierInformation();
943  clock_t iterationStopClock = clock();
944  double iterationTime = (iterationStopClock - iterationStartClock)
945  / (double)CLOCKS_PER_SEC;
946  emit statusUpdate( QString("End of Iteration %1 \n").arg(m_iteration) );
947  emit statusUpdate( QString("Elapsed Time: %1 \n").arg(iterationTime,
948  fieldWidth,
949  format,
950  precision) );
951 
952  // check for maximum iterations
953  if (m_iteration >= m_bundleSettings->convergenceCriteriaMaximumIterations()) {
954  emit(statusBarUpdate("Max Iterations Reached"));
955  break;
956  }
957 
958  // restart the dynamic calculation of the cumulative probility distribution of residuals
959  // (in unweighted pixels) --so it will be up to date for the next iteration
960  if (!m_bundleResults.converged()) {
961  m_bundleResults.initializeResidualsProbabilityDistribution(101);
962  }
963  // TODO: is this necessary ???
964  // probably all ready initialized to 101 nodes in bundle settings constructor...
965 
966  // if we're using CHOLMOD and still going, release cholmod_factor
967  // (if we don't, memory leaks will occur), otherwise we need it for error propagation
968  if (!m_bundleResults.converged() || !m_bundleSettings->errorPropagation()) {
969  cholmod_free_factor(&m_L, &m_cholmodCommon);
970  }
971 
972  iterationSummary();
973 
974  m_iteration++;
975 
976  previousSigma0 = m_bundleResults.sigma0();
977  }
978 
979  if (m_bundleResults.converged() && m_bundleSettings->errorPropagation()) {
980  clock_t errorPropStartClock = clock();
981 
982  outputBundleStatus("\nStarting Error Propagation");
983 
984  errorPropagation();
985  emit statusUpdate("\n\nError Propagation Complete\n");
986  clock_t errorPropStopClock = clock();
987  m_bundleResults.setElapsedTimeErrorProp((errorPropStopClock - errorPropStartClock)
988  / (double)CLOCKS_PER_SEC);
989  }
990 
991  clock_t solveStopClock = clock();
992  m_bundleResults.setElapsedTime((solveStopClock - solveStartClock)
993  / (double)CLOCKS_PER_SEC);
994 
995  wrapUp();
996 
997  m_bundleResults.setIterations(m_iteration);
998  m_bundleResults.setObservations(m_bundleObservations);
999  m_bundleResults.setBundleControlPoints(m_bundleControlPoints);
1000 
1001  emit resultsReady(bundleSolveInformation());
1002 
1003  emit statusUpdate("\nBundle Complete\n");
1004 
1005  iterationSummary();
1006  }
1007  catch (IException &e) {
1008  m_bundleResults.setConverged(false);
1009  emit statusUpdate("\n aborting...");
1010  emit statusBarUpdate("Failed to Converge");
1011  emit finished();
1012  QString msg = "Could not solve bundle adjust.";
1013  throw IException(e, e.errorType(), msg, _FILEINFO_);
1014  }
1015 
1016  emit finished();
1017  return true;
1018  }
1019 
1020 
1029  BundleSolutionInfo *BundleAdjust::bundleSolveInformation() {
1030  BundleSolutionInfo *bundleSolutionInfo = new BundleSolutionInfo(m_bundleSettings,
1031  FileName(m_cnetFileName),
1032  m_bundleResults,
1033  imageLists());
1034  bundleSolutionInfo->setRunTime("");
1035  return bundleSolutionInfo;
1036  }
1037 
1038 
1050  bool BundleAdjust::formNormalEquations() {
1051  emit(statusBarUpdate("Forming Normal Equations"));
1052  bool status = false;
1053 
1054  m_bundleResults.setNumberObservations(0);// ???
1055  m_bundleResults.resetNumberConstrainedPointParameters();//???
1056 
1057  // Initialize auxilary matrices and vectors.
1058  static LinearAlgebra::Matrix coeffTarget;
1059  static LinearAlgebra::Matrix coeffImage;
1060  static LinearAlgebra::Matrix coeffPoint3D(2, 3);
1061  static LinearAlgebra::Vector coeffRHS(2);
1062  static boost::numeric::ublas::symmetric_matrix<double, upper> N22(3);
1064  static LinearAlgebra::Vector n2(3);
1065  boost::numeric::ublas::compressed_vector<double> n1(m_rank);
1066 
1067  m_RHS.resize(m_rank);
1068 
1069  // if solving for target body parameters, set size of coeffTarget
1070  // (note this size will not change through the adjustment).
1071  if (m_bundleSettings->solveTargetBody()) {
1072  int numTargetBodyParameters = m_bundleSettings->numberTargetBodyParameters();
1073  // TODO make sure numTargetBodyParameters is greater than 0
1074  coeffTarget.resize(2,numTargetBodyParameters);
1075  }
1076 
1077  // clear N12, n1, and nj
1078  N12.clear();
1079  n1.clear();
1080  m_RHS.clear();
1081 
1082  // clear static matrices
1083  coeffPoint3D.clear();
1084  coeffRHS.clear();
1085  N22.clear();
1086  n2.clear();
1087 
1088  // loop over 3D points
1089  int numGood3DPoints = 0;
1090  int numRejected3DPoints = 0;
1091  int pointIndex = 0;
1092  int num3DPoints = m_bundleControlPoints.size();
1093 
1094  outputBundleStatus("\n\n");
1095 
1096  for (int i = 0; i < num3DPoints; i++) {
1097  emit(pointUpdate(i+1));
1098  BundleControlPointQsp point = m_bundleControlPoints.at(i);
1099 
1100  if (point->isRejected()) {
1101  numRejected3DPoints++;
1102 
1103  pointIndex++;
1104  continue;
1105  }
1106 
1107  if ( i != 0 ) {
1108  N22.clear();
1109  N12.wipe();
1110  n2.clear();
1111  }
1112 
1113  // loop over measures for this point
1114  int numMeasures = point->size();
1115  for (int j = 0; j < numMeasures; j++) {
1116  BundleMeasureQsp measure = point->at(j);
1117 
1118  // flagged as "JigsawFail" implies this measure has been rejected
1119  // TODO IsRejected is obsolete -- replace code or add to ControlMeasure
1120  if (measure->isRejected()) {
1121  continue;
1122  }
1123 
1124  status = computePartials(coeffTarget, coeffImage, coeffPoint3D, coeffRHS, *measure,
1125  *point);
1126 
1127  if (!status) {
1128  // TODO should status be set back to true? JAM
1129  // TODO this measure should be flagged as rejected.
1130  continue;
1131  }
1132 
1133  // update number of observations
1134  int numObs = m_bundleResults.numberObservations();
1135  m_bundleResults.setNumberObservations(numObs + 2);
1136 
1137  formMeasureNormals(N22, N12, n1, n2, coeffTarget, coeffImage, coeffPoint3D, coeffRHS,
1138  measure->observationIndex());
1139 
1140  } // end loop over this points measures
1141 
1142  formPointNormals(N22, N12, n2, m_RHS, point);
1143 
1144  pointIndex++;
1145 
1146  numGood3DPoints++;
1147 
1148  } // end loop over 3D points
1149 
1150  // finally, form the reduced normal equations
1151  formWeightedNormals(n1, m_RHS);
1152 
1153  // update number of unknown parameters
1154  m_bundleResults.setNumberUnknownParameters(m_rank + 3 * numGood3DPoints);
1155 
1156  return status;
1157 }
1158 
1159 
1179  bool BundleAdjust::formMeasureNormals(symmetric_matrix<double, upper>&N22,
1181  compressed_vector<double> &n1,
1182  vector<double> &n2,
1183  matrix<double> &coeffTarget,
1184  matrix<double> &coeffImage,
1185  matrix<double> &coeffPoint3D,
1186  vector<double> &coeffRHS,
1187  int observationIndex) {
1188 
1189  static symmetric_matrix<double, upper> N11;
1190  static matrix<double> N11TargetImage;
1191 
1192  int blockIndex = observationIndex;
1193 
1194  // if we are solving for target body parameters
1195  int numTargetPartials = coeffTarget.size2();
1196  if (m_bundleSettings->solveTargetBody()) {
1197  blockIndex++;
1198 
1199  static vector<double> n1Target(numTargetPartials);
1200  n1Target.resize(numTargetPartials);
1201  n1Target.clear();
1202 
1203  // form N11 (normals for target body)
1204  N11.resize(numTargetPartials);
1205  N11.clear();
1206 
1207  N11 = prod(trans(coeffTarget), coeffTarget);
1208 
1209  // insert submatrix at column, row
1210  m_sparseNormals.insertMatrixBlock(0, 0, numTargetPartials, numTargetPartials);
1211 
1212  (*(*m_sparseNormals[0])[0]) += N11;
1213 
1214  // form portion of N11 between target and image
1215  N11TargetImage.resize(numTargetPartials, coeffImage.size2());
1216  N11TargetImage.clear();
1217  N11TargetImage = prod(trans(coeffTarget),coeffImage);
1218 
1219  m_sparseNormals.insertMatrixBlock(blockIndex, 0,
1220  numTargetPartials, coeffImage.size2());
1221  (*(*m_sparseNormals[blockIndex])[0]) += N11TargetImage;
1222 
1223  // form N12 target portion
1224  static matrix<double> N12Target(numTargetPartials, 3);
1225  N12Target.clear();
1226 
1227  N12Target = prod(trans(coeffTarget), coeffPoint3D);
1228 
1229  // insert N12Target into N12
1230  N12.insertMatrixBlock(0, numTargetPartials, 3);
1231  *N12[0] += N12Target;
1232 
1233  // form n1Target
1234  n1Target = prod(trans(coeffTarget), coeffRHS);
1235 
1236  // insert n1Target into n1
1237  for (int i = 0; i < numTargetPartials; i++) {
1238  n1(i) += n1Target(i);
1239  }
1240  }
1241 
1242 
1243 //~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ below is ok (2015-06-03)
1244 // TODO - if solving for target (and/or self-cal) have to use not observationIndex below but
1245 // observationIndex plus 1 or 2
1246 
1247  int numImagePartials = coeffImage.size2();
1248 
1249  static LinearAlgebra::Vector n1Image(numImagePartials);
1250  n1Image.resize(numImagePartials);
1251  n1Image.clear();
1252 
1253  // form N11 (normals for photo)
1254  N11.resize(numImagePartials);
1255  N11.clear();
1256 
1257  N11 = prod(trans(coeffImage), coeffImage);
1258 
1259  int t = m_sparseNormals.at(blockIndex)->startColumn();
1260 
1261  // insert submatrix at column, row
1262  m_sparseNormals.insertMatrixBlock(blockIndex, blockIndex,
1263  numImagePartials, numImagePartials);
1264 
1265  (*(*m_sparseNormals[blockIndex])[blockIndex]) += N11;
1266 
1267  // form N12Image
1268  static matrix<double> N12Image(numImagePartials, 3);
1269  N12Image.resize(numImagePartials, 3);
1270  N12Image.clear();
1271 
1272  N12Image = prod(trans(coeffImage), coeffPoint3D);
1273 
1274  // insert N12Image into N12
1275  N12.insertMatrixBlock(blockIndex, numImagePartials, 3);
1276  *N12[blockIndex] += N12Image;
1277 
1278  // form n1
1279  n1Image = prod(trans(coeffImage), coeffRHS);
1280 
1281  // insert n1Image into n1
1282  // TODO - MUST ACCOUNT FOR TARGET BODY PARAMETERS
1283  // WHEN INSERTING INTO n1 HERE!!!!!
1284  for (int i = 0; i < numImagePartials; i++) {
1285  n1(i + t) += n1Image(i);
1286  }
1287 
1288  // form N22
1289  N22 += prod(trans(coeffPoint3D), coeffPoint3D);
1290 
1291  // form n2
1292  n2 += prod(trans(coeffPoint3D), coeffRHS);
1293 
1294  return true;
1295  }
1296 
1297 
1315  bool BundleAdjust::formPointNormals(symmetric_matrix<double, upper>&N22,
1317  vector<double> &n2,
1318  vector<double> &nj,
1319  BundleControlPointQsp &bundleControlPoint) {
1320 
1321  boost::numeric::ublas::bounded_vector<double, 3> &NIC = bundleControlPoint->nicVector();
1322  SparseBlockRowMatrix &Q = bundleControlPoint->cholmodQMatrix();
1323 
1324  NIC.clear();
1325  Q.zeroBlocks();
1326 
1327  // weighting of 3D point parameters
1328  // Make sure weights are in the units corresponding to the bundle coordinate type
1329  boost::numeric::ublas::bounded_vector<double, 3> &weights
1330  = bundleControlPoint->weights();
1331  boost::numeric::ublas::bounded_vector<double, 3> &corrections
1332  = bundleControlPoint->corrections();
1333 
1334  if (weights(0) > 0.0) {
1335  N22(0,0) += weights(0);
1336  n2(0) += (-weights(0) * corrections(0));
1337  m_bundleResults.incrementNumberConstrainedPointParameters(1);
1338  }
1339 
1340  if (weights(1) > 0.0) {
1341  N22(1,1) += weights(1);
1342  n2(1) += (-weights(1) * corrections(1));
1343  m_bundleResults.incrementNumberConstrainedPointParameters(1);
1344  }
1345 
1346  if (weights(2) > 0.0) {
1347  N22(2,2) += weights(2);
1348  n2(2) += (-weights(2) * corrections(2));
1349  m_bundleResults.incrementNumberConstrainedPointParameters(1);
1350  }
1351 
1352  // invert N22
1353  invert3x3(N22);
1354 
1355  // save upper triangular covariance matrix for error propagation
1356  SurfacePoint SurfacePoint = bundleControlPoint->adjustedSurfacePoint();
1357  SurfacePoint.SetMatrix(m_bundleSettings->controlPointCoordTypeBundle(), N22);
1358  bundleControlPoint->setAdjustedSurfacePoint(SurfacePoint);
1359 
1360  // form Q (this is N22{-1} * N12{T})
1361  productATransB(N22, N12, Q);
1362 
1363  // form product of N22(inverse) and n2; store in NIC
1364  NIC = prod(N22, n2);
1365 
1366  // accumulate -R directly into reduced normal equations
1367  productAB(N12, Q);
1368 
1369  // accumulate -nj
1370  accumProductAlphaAB(-1.0, Q, n2, nj);
1371 
1372  return true;
1373  }
1374 
1375 
1387  bool BundleAdjust::formWeightedNormals(compressed_vector<double> &n1,
1388  vector<double> &nj) {
1389 
1390  m_bundleResults.resetNumberConstrainedImageParameters();
1391 
1392  int n = 0;
1393 
1394  for (int i = 0; i < m_sparseNormals.size(); i++) {
1395  LinearAlgebra::Matrix *diagonalBlock = m_sparseNormals.getBlock(i, i);
1396  if ( !diagonalBlock )
1397  continue;
1398 
1399  if (m_bundleSettings->solveTargetBody() && i == 0) {
1400  m_bundleResults.resetNumberConstrainedTargetParameters();
1401 
1402  // get parameter weights for target body
1403  vector<double> weights = m_bundleTargetBody->parameterWeights();
1404  vector<double> corrections = m_bundleTargetBody->parameterCorrections();
1405 
1406  int blockSize = diagonalBlock->size1();
1407  for (int j = 0; j < blockSize; j++) {
1408  if (weights[j] > 0.0) {
1409  (*diagonalBlock)(j,j) += weights[j];
1410  nj[n] -= weights[j] * corrections(j);
1411  m_bundleResults.incrementNumberConstrainedTargetParameters(1);
1412  }
1413  n++;
1414  }
1415  }
1416  else {
1417  BundleObservationQsp observation;
1418 
1419  // get parameter weights for this observation
1420  if (m_bundleSettings->solveTargetBody()) {
1421  observation = m_bundleObservations.at(i-1);
1422  }
1423  else {
1424  observation = m_bundleObservations.at(i);
1425  }
1426 
1427  LinearAlgebra::Vector weights = observation->parameterWeights();
1428  LinearAlgebra::Vector corrections = observation->parameterCorrections();
1429 
1430  int blockSize = diagonalBlock->size1();
1431  for (int j = 0; j < blockSize; j++) {
1432  if (weights(j) > 0.0) {
1433  (*diagonalBlock)(j,j) += weights(j);
1434  nj[n] -= weights(j) * corrections(j);
1435  m_bundleResults.incrementNumberConstrainedImageParameters(1);
1436  }
1437  n++;
1438  }
1439  }
1440  }
1441 
1442  // add n1 to nj
1443  nj += n1;
1444 
1445  return true;
1446  }
1447 
1448 
1458  bool BundleAdjust::productATransB(symmetric_matrix <double,upper> &N22,
1460  SparseBlockRowMatrix &Q) {
1461 
1462  QMapIterator< int, LinearAlgebra::Matrix * > N12it(N12);
1463 
1464  while ( N12it.hasNext() ) {
1465  N12it.next();
1466 
1467  int rowIndex = N12it.key();
1468 
1469  // insert submatrix in Q at block "rowIndex"
1470  Q.insertMatrixBlock(rowIndex, 3, N12it.value()->size1());
1471 
1472  *(Q[rowIndex]) = prod(N22,trans(*(N12it.value())));
1473  }
1474 
1475  return true;
1476  }
1477 
1478 
1488  void BundleAdjust::productAB(SparseBlockColumnMatrix &N12,
1489  SparseBlockRowMatrix &Q) {
1490  // iterators for N12 and Q
1491  QMapIterator<int, LinearAlgebra::Matrix*> N12it(N12);
1492  QMapIterator<int, LinearAlgebra::Matrix*> Qit(Q);
1493 
1494  // now multiply blocks and subtract from m_sparseNormals
1495  while ( N12it.hasNext() ) {
1496  N12it.next();
1497 
1498  int rowIndex = N12it.key();
1499  LinearAlgebra::Matrix *N12block = N12it.value();
1500 
1501  while ( Qit.hasNext() ) {
1502  Qit.next();
1503 
1504  int columnIndex = Qit.key();
1505 
1506  if ( rowIndex > columnIndex ) {
1507  continue;
1508  }
1509 
1510  LinearAlgebra::Matrix *Qblock = Qit.value();
1511 
1512  // insert submatrix at column, row
1513  m_sparseNormals.insertMatrixBlock(columnIndex, rowIndex,
1514  N12block->size1(), Qblock->size2());
1515 
1516  (*(*m_sparseNormals[columnIndex])[rowIndex]) -= prod(*N12block,*Qblock);
1517  }
1518  Qit.toFront();
1519  }
1520  }
1521 
1522 
1533  void BundleAdjust::accumProductAlphaAB(double alpha,
1535  vector<double> &n2,
1536  vector<double> &nj) {
1537 
1538  if (alpha == 0.0) {
1539  return;
1540  }
1541 
1542  int numParams;
1543 
1544  QMapIterator<int, LinearAlgebra::Matrix*> Qit(Q);
1545 
1546  while ( Qit.hasNext() ) {
1547  Qit.next();
1548 
1549  int columnIndex = Qit.key();
1550  LinearAlgebra::Matrix *Qblock = Qit.value();
1551 
1552  LinearAlgebra::Vector blockProduct = prod(trans(*Qblock),n2);
1553 
1554  numParams = m_sparseNormals.at(columnIndex)->startColumn();
1555 
1556  for (unsigned i = 0; i < blockProduct.size(); i++) {
1557  nj(numParams+i) += alpha*blockProduct(i);
1558  }
1559  }
1560  }
1561 
1562 
1572  bool BundleAdjust::solveSystem() {
1573 
1574  // load cholmod triplet
1575  if ( !loadCholmodTriplet() ) {
1576  QString msg = "CHOLMOD: Failed to load Triplet matrix";
1577  throw IException(IException::Programmer, msg, _FILEINFO_);
1578  }
1579 
1580  // convert triplet to sparse matrix
1581  m_cholmodNormal = cholmod_triplet_to_sparse(m_cholmodTriplet,
1582  m_cholmodTriplet->nnz,
1583  &m_cholmodCommon);
1584 
1585  // analyze matrix
1586  // TODO should we analyze just 1st iteration?
1587  m_L = cholmod_analyze(m_cholmodNormal, &m_cholmodCommon);
1588 
1589  // create cholmod cholesky factor
1590  // CHOLMOD will choose LLT or LDLT decomposition based on the characteristics of the matrix.
1591  cholmod_factorize(m_cholmodNormal, m_L, &m_cholmodCommon);
1592 
1593  // check for "matrix not positive definite" error
1594  if (m_cholmodCommon.status == CHOLMOD_NOT_POSDEF) {
1595  QString msg = "Matrix NOT positive-definite: failure at column " + toString((int) m_L->minor);
1596 // throw IException(IException::User, msg, _FILEINFO_);
1597  error(msg);
1598  emit(finished());
1599  return false;
1600  }
1601 
1602  // cholmod solution and right-hand side vectors
1603  cholmod_dense *x, *b;
1604 
1605  // initialize right-hand side vector
1606  b = cholmod_zeros(m_cholmodNormal->nrow, 1, m_cholmodNormal->xtype, &m_cholmodCommon);
1607 
1608  // copy right-hand side vector into b
1609  double *px = (double*)b->x;
1610  for (int i = 0; i < m_rank; i++) {
1611  px[i] = m_RHS[i];
1612  }
1613 
1614  // cholmod solve
1615  x = cholmod_solve(CHOLMOD_A, m_L, b, &m_cholmodCommon);
1616 
1617  // copy solution vector x out into m_imageSolution
1618  double *sx = (double*)x->x;
1619  for (int i = 0; i < m_rank; i++) {
1620  m_imageSolution[i] = sx[i];
1621  }
1622 
1623  // free cholmod structures
1624  cholmod_free_sparse(&m_cholmodNormal, &m_cholmodCommon);
1625  cholmod_free_dense(&b, &m_cholmodCommon);
1626  cholmod_free_dense(&x, &m_cholmodCommon);
1627 
1628  return true;
1629  }
1630 
1631 
1643  bool BundleAdjust::loadCholmodTriplet() {
1644 
1645  if ( m_iteration == 1 ) {
1646  int numElements = m_sparseNormals.numberOfElements();
1647  m_cholmodTriplet = cholmod_allocate_triplet(m_rank, m_rank, numElements,
1648  -1, CHOLMOD_REAL, &m_cholmodCommon);
1649 
1650  if ( !m_cholmodTriplet ) {
1651  outputBundleStatus("\nTriplet allocation failure\n");
1652  return false;
1653  }
1654  m_cholmodTriplet->nnz = 0;
1655  }
1656 
1657  int *tripletColumns = (int*) m_cholmodTriplet->i;
1658  int *tripletRows = (int*) m_cholmodTriplet->j;
1659  double *tripletValues = (double*)m_cholmodTriplet->x;
1660 
1661  double entryValue;
1662 
1663  int numEntries = 0;
1664 
1665  int numBlockcolumns = m_sparseNormals.size();
1666  for (int columnIndex = 0; columnIndex < numBlockcolumns; columnIndex++) {
1667 
1668  SparseBlockColumnMatrix *normalsColumn = m_sparseNormals[columnIndex];
1669 
1670  if ( !normalsColumn ) {
1671  QString status = "\nSparseBlockColumnMatrix retrieval failure at column " +
1672  QString::number(columnIndex);
1673  outputBundleStatus(status);
1674  return false;
1675  }
1676 
1677  int numLeadingColumns = normalsColumn->startColumn();
1678 
1679  QMapIterator< int, LinearAlgebra::Matrix * > it(*normalsColumn);
1680 
1681  while ( it.hasNext() ) {
1682  it.next();
1683 
1684  int rowIndex = it.key();
1685 
1686  // note: as the normal equations matrix is symmetric, the # of leading rows for a block is
1687  // equal to the # of leading columns for a block column at the "rowIndex" position
1688  int numLeadingRows = m_sparseNormals.at(rowIndex)->startColumn();
1689 
1690  LinearAlgebra::Matrix *normalsBlock = it.value();
1691  if ( !normalsBlock ) {
1692  QString status = "\nmatrix block retrieval failure at column ";
1693  status.append(columnIndex);
1694  status.append(", row ");
1695  status.append(rowIndex);
1696  outputBundleStatus(status);
1697  status = "Total # of block columns: " + QString::number(numBlockcolumns);
1698  outputBundleStatus(status);
1699  status = "Total # of blocks: " + QString::number(m_sparseNormals.numberOfBlocks());
1700  outputBundleStatus(status);
1701  return false;
1702  }
1703 
1704  if ( columnIndex == rowIndex ) { // diagonal block (upper-triangular)
1705  for (unsigned ii = 0; ii < normalsBlock->size1(); ii++) {
1706  for (unsigned jj = ii; jj < normalsBlock->size2(); jj++) {
1707  entryValue = normalsBlock->at_element(ii,jj);
1708  int entryColumnIndex = jj + numLeadingColumns;
1709  int entryRowIndex = ii + numLeadingRows;
1710 
1711  if ( m_iteration == 1 ) {
1712  tripletColumns[numEntries] = entryColumnIndex;
1713  tripletRows[numEntries] = entryRowIndex;
1714  m_cholmodTriplet->nnz++;
1715  }
1716 
1717  tripletValues[numEntries] = entryValue;
1718 
1719  numEntries++;
1720  }
1721  }
1722  }
1723  else { // off-diagonal block (square)
1724  for (unsigned ii = 0; ii < normalsBlock->size1(); ii++) {
1725  for (unsigned jj = 0; jj < normalsBlock->size2(); jj++) {
1726  entryValue = normalsBlock->at_element(ii,jj);
1727  int entryColumnIndex = jj + numLeadingColumns;
1728  int entryRowIndex = ii + numLeadingRows;
1729 
1730  if ( m_iteration ==1 ) {
1731  tripletColumns[numEntries] = entryRowIndex;
1732  tripletRows[numEntries] = entryColumnIndex;
1733  m_cholmodTriplet->nnz++;
1734  }
1735 
1736  tripletValues[numEntries] = entryValue;
1737 
1738  numEntries++;
1739  }
1740  }
1741  }
1742  }
1743  }
1744 
1745  return true;
1746  }
1747 
1748 
1757  bool BundleAdjust::cholmodInverse() {
1758  int i, j;
1759 
1760  // allocate matrix inverse
1761  m_normalInverse.resize(m_rank);
1762 
1763  cholmod_dense *x; // solution vector
1764  cholmod_dense *b; // right-hand side (column vectors of identity)
1765 
1766  b = cholmod_zeros ( m_rank, 1, CHOLMOD_REAL, &m_cholmodCommon ) ;
1767  double *pb = (double*)b->x;
1768 
1769  double *px = NULL;
1770 
1771  for (i = 0; i < m_rank; i++) {
1772  if ( i > 0 ) {
1773  pb[i-1] = 0.0;
1774  }
1775  pb[i] = 1.0;
1776 
1777  x = cholmod_solve ( CHOLMOD_A, m_L, b, &m_cholmodCommon ) ;
1778  px = (double*)x->x;
1779 
1780  // store solution in corresponding column of inverse
1781  for (j = 0; j <= i; j++) {
1782  m_normalInverse(j, i) = px[j];
1783  }
1784 
1785  cholmod_free_dense(&x,&m_cholmodCommon);
1786  }
1787 
1788  cholmod_free_dense(&b,&m_cholmodCommon);
1789 
1790  return true;
1791  }
1792 
1793 
1806  bool BundleAdjust::invert3x3(symmetric_matrix<double, upper> &m) {
1807  double det;
1808  double den;
1809 
1810  boost::numeric::ublas::symmetric_matrix< double, upper > c = m;
1811 
1812  den = m(0, 0) * (m(1, 1) * m(2, 2) - m(1, 2) * m(2, 1))
1813  - m(0, 1) * (m(1, 0) * m(2, 2) - m(1, 2) * m(2, 0))
1814  + m(0, 2) * (m(1, 0) * m(2, 1) - m(1, 1) * m(2, 0));
1815 
1816  // check for divide by zero
1817  if (fabs(den) < 1.0e-100) {
1818  return false;
1819  }
1820 
1821  det = 1.0 / den;
1822 
1823  m(0, 0) = (c(1, 1) * c(2, 2) - c(1, 2) * c(2, 1)) * det;
1824  m(0, 1) = (c(0, 2) * c(2, 1) - c(0, 1) * c(2, 2)) * det;
1825  m(0, 2) = (c(0, 1) * c(1, 2) - c(0, 2) * c(1, 1)) * det;
1826  m(1, 1) = (c(0, 0) * c(2, 2) - c(0, 2) * c(2, 0)) * det;
1827  m(1, 2) = (c(0, 2) * c(1, 0) - c(0, 0) * c(1, 2)) * det;
1828  m(2, 2) = (c(0, 0) * c(1, 1) - c(0, 1) * c(1, 0)) * det;
1829 
1830  return true;
1831  }
1832 
1833 
1853  bool BundleAdjust::computePartials(matrix<double> &coeffTarget,
1854  matrix<double> &coeffImage,
1855  matrix<double> &coeffPoint3D,
1856  vector<double> &coeffRHS,
1857  BundleMeasure &measure,
1858  BundleControlPoint &point) {
1859 
1860  Camera *measureCamera = measure.camera();
1861  BundleObservationQsp observation = measure.parentBundleObservation();
1862 
1863  int numImagePartials = observation->numberParameters();
1864 
1865  // we're saving the number of image partials in m_previousNumberImagePartials
1866  // to compare to the previous computePartials call to avoid unnecessary resizing of the
1867  // coeffImage matrix
1868  if (numImagePartials != m_previousNumberImagePartials) {
1869  coeffImage.resize(2,numImagePartials);
1870  m_previousNumberImagePartials = numImagePartials;
1871  }
1872 
1873  // No need to call SetImage for framing camera
1874  if (measureCamera->GetCameraType() != Camera::Framing) {
1875  // Set the Spice to the measured point. A framing camera exposes the entire image at one time.
1876  // It will have a single set of Spice for the entire image. Scanning cameras may populate a single
1877  // image with multiple exposures, each with a unique set of Spice. SetImage needs to be called
1878  // repeatedly for these images to point to the Spice for the current pixel.
1879  measureCamera->SetImage(measure.sample(), measure.line());
1880  }
1881 
1882  // CSM Cameras do not have a ground map
1883  if (measureCamera->GetCameraType() != Camera::Csm) {
1884  // Compute the look vector in instrument coordinates based on time of observation and apriori
1885  // lat/lon/radius. As of 05/15/2019, this call no longer does the back-of-planet test. An optional
1886  // bool argument was added CameraGroundMap::GetXY to turn off the test.
1887  double computedX, computedY;
1888  if (!(measureCamera->GroundMap()->GetXY(point.adjustedSurfacePoint(),
1889  &computedX, &computedY, false))) {
1890  QString msg = "Unable to map apriori surface point for measure ";
1891  msg += measure.cubeSerialNumber() + " on point " + point.id() + " into focal plane";
1892  throw IException(IException::User, msg, _FILEINFO_);
1893  }
1894  }
1895 
1896  if (m_bundleSettings->solveTargetBody()) {
1897  observation->computeTargetPartials(coeffTarget, measure, m_bundleSettings, m_bundleTargetBody);
1898  }
1899 
1900  observation->computeImagePartials(coeffImage, measure);
1901 
1902  // Complete partials calculations for 3D point (latitudinal or rectangular)
1903  // Retrieve the coordinate type (latitudinal or rectangular) and compute the partials for
1904  // the fixed point with respect to each coordinate in Body-Fixed
1905  SurfacePoint::CoordinateType coordType = m_bundleSettings->controlPointCoordTypeBundle();
1906  observation->computePoint3DPartials(coeffPoint3D, measure, coordType);
1907 
1908  // right-hand side (measured - computed)
1909  observation->computeRHSPartials(coeffRHS, measure);
1910 
1911  double deltaX = coeffRHS(0);
1912  double deltaY = coeffRHS(1);
1913 
1914  m_bundleResults.addResidualsProbabilityDistributionObservation(observation->computeObservationValue(measure, deltaX));
1915  m_bundleResults.addResidualsProbabilityDistributionObservation(observation->computeObservationValue(measure, deltaY));
1916 
1917  if (m_bundleResults.numberMaximumLikelihoodModels()
1918  > m_bundleResults.maximumLikelihoodModelIndex()) {
1919  // If maximum likelihood estimation is being used
1920  double residualR2ZScore = sqrt(deltaX * deltaX + deltaY * deltaY) / sqrt(2.0);
1921 
1922  // Dynamically build the cumulative probability distribution of the R^2 residual Z Scores
1923  m_bundleResults.addProbabilityDistributionObservation(residualR2ZScore);
1924 
1925  int currentModelIndex = m_bundleResults.maximumLikelihoodModelIndex();
1926  double observationWeight = m_bundleResults.maximumLikelihoodModelWFunc(currentModelIndex)
1927  .sqrtWeightScaler(residualR2ZScore);
1928  coeffImage *= observationWeight;
1929  coeffPoint3D *= observationWeight;
1930  coeffRHS *= observationWeight;
1931 
1932  if (m_bundleSettings->solveTargetBody()) {
1933  coeffTarget *= observationWeight;
1934  }
1935  }
1936 
1937  return true;
1938  }
1939 
1940 
1944  void BundleAdjust::applyParameterCorrections() {
1945  emit(statusBarUpdate("Updating Parameters"));
1946  int t = 0;
1947 
1948  // TODO - update target body parameters if in solution
1949  // note these come before BundleObservation parameters in normal equations matrix
1950  if (m_bundleSettings->solveTargetBody()) {
1951  int numTargetBodyParameters = m_bundleTargetBody->numberParameters();
1952 
1953  m_bundleTargetBody->applyParameterCorrections(subrange(m_imageSolution,0,
1954  numTargetBodyParameters));
1955 
1956  t += numTargetBodyParameters;
1957  }
1958 
1959  // Update spice for each BundleObservation
1960  int numObservations = m_bundleObservations.size();
1961  for (int i = 0; i < numObservations; i++) {
1962  BundleObservationQsp observation = m_bundleObservations.at(i);
1963 
1964  int numParameters = observation->numberParameters();
1965 
1966  observation->applyParameterCorrections(subrange(m_imageSolution,t,t+numParameters));
1967 
1968  if (m_bundleSettings->solveTargetBody()) {
1969  // TODO: needs to be updated for ISIS vs. CSM CSM has no updateBodyRotation]
1970  // TODO: this is no good.
1971  QSharedPointer<IsisBundleObservation> isisObservation = qSharedPointerDynamicCast<IsisBundleObservation>(observation);
1972  isisObservation->updateBodyRotation();
1973  }
1974 
1975  t += numParameters;
1976  }
1977  // TODO - When BundleXYZ gets merged into dev, go with Ken's version of merging the updating of
1978  // of the adjusted surface point into BundleControlPoint.
1979 
1980  int pointIndex = 0;
1981  int numControlPoints = m_bundleControlPoints.size();
1982 
1983  for (int i = 0; i < numControlPoints; i++) {
1984  BundleControlPointQsp point = m_bundleControlPoints.at(i);
1985 
1986  if (point->isRejected()) {
1987  pointIndex++;
1988  continue;
1989  }
1990 
1991  point->applyParameterCorrections(m_imageSolution, m_sparseNormals,
1992  m_bundleTargetBody);
1993  pointIndex++;
1994 
1995  } // end loop over point corrections
1996  }
1997 
1998 
2010  double BundleAdjust::computeResiduals() {
2011  emit(statusBarUpdate("Computing Residuals"));
2012  double vtpv = 0.0;
2013  double vtpvControl = 0.0;
2014  double vtpvImage = 0.0;
2015  double weight;
2016  double v, vx, vy;
2017 
2018  // clear residual stats vectors
2019  m_xResiduals.Reset();
2020  m_yResiduals.Reset();
2021  m_xyResiduals.Reset();
2022 
2023  // vtpv for image coordinates
2024  int numObjectPoints = m_bundleControlPoints.size();
2025 
2026  for (int i = 0; i < numObjectPoints; i++) {
2027 
2028  BundleControlPointQsp point = m_bundleControlPoints.at(i);
2029 
2030  point->computeResiduals();
2031 
2032  int numMeasures = point->numberOfMeasures();
2033  for (int j = 0; j < numMeasures; j++) {
2034  const BundleMeasureQsp measure = point->at(j);
2035 
2036  weight = 1.4 * (measure->camera())->PixelPitch();
2037  weight = 1.0 / weight;
2038  weight *= weight;
2039 
2040  vx = measure->focalPlaneMeasuredX() - measure->focalPlaneComputedX();
2041  vy = measure->focalPlaneMeasuredY() - measure->focalPlaneComputedY();
2042 
2043  // if rejected, don't include in statistics
2044  if (measure->isRejected()) {
2045  continue;
2046  }
2047 
2048  m_xResiduals.AddData(vx);
2049  m_yResiduals.AddData(vy);
2050  m_xyResiduals.AddData(vx);
2051  m_xyResiduals.AddData(vy);
2052 
2053  vtpv += vx * vx * weight + vy * vy * weight;
2054  }
2055  }
2056 
2057  // add vtpv from constrained 3D points
2058  int pointIndex = 0; // *** TODO *** This does not appear to be used. Delete? DAC 07-14-2017
2059  for (int i = 0; i < numObjectPoints; i++) {
2060  BundleControlPointQsp bundleControlPoint = m_bundleControlPoints.at(i);
2061 
2062  // get weight and correction vector for this point
2063  boost::numeric::ublas::bounded_vector<double, 3> weights = bundleControlPoint->weights();
2064  boost::numeric::ublas::bounded_vector<double, 3> corrections =
2065  bundleControlPoint->corrections();
2066 
2067  if ( weights(0) > 0.0 ) {
2068  vtpvControl += corrections(0) * corrections(0) * weights(0);
2069  }
2070  if ( weights(1) > 0.0 ) {
2071  vtpvControl += corrections(1) * corrections(1) * weights(1);
2072  }
2073  if ( weights(2) > 0.0 ) {
2074  vtpvControl += corrections(2) * corrections(2) * weights(2);
2075  }
2076 
2077  pointIndex++;
2078  }
2079 
2080  // add vtpv from constrained image parameters
2081  for (int i = 0; i < m_bundleObservations.size(); i++) {
2082  BundleObservationQsp observation = m_bundleObservations.at(i);
2083 
2084  // get weight and correction vector for this observation
2085  const LinearAlgebra::Vector &weights = observation->parameterWeights();
2086  const LinearAlgebra::Vector &corrections = observation->parameterCorrections();
2087 
2088  for (int j = 0; j < (int)corrections.size(); j++) {
2089  if (weights[j] > 0.0) {
2090  v = corrections[j];
2091  vtpvImage += v * v * weights[j];
2092  }
2093  }
2094  }
2095 
2096  // TODO - add vtpv from constrained target body parameters
2097  double vtpvTargetBody = 0.0;
2098  if ( m_bundleTargetBody) {
2099  vtpvTargetBody = m_bundleTargetBody->vtpv();
2100  }
2101 
2102  vtpv = vtpv + vtpvControl + vtpvImage + vtpvTargetBody;
2103 
2104  // Compute rms for all image coordinate residuals
2105  // separately for x, y, then x and y together
2106  m_bundleResults.setRmsXYResiduals(m_xResiduals.Rms(), m_yResiduals.Rms(), m_xyResiduals.Rms());
2107 
2108  return vtpv;
2109  }
2110 
2111 
2118  bool BundleAdjust::wrapUp() {
2119  // compute residuals in pixels
2120 
2121  // vtpv for image coordinates
2122  for (int i = 0; i < m_bundleControlPoints.size(); i++) {
2123  BundleControlPointQsp point = m_bundleControlPoints.at(i);
2124  point->computeResiduals();
2125  }
2126 
2127  computeBundleStatistics();
2128 
2129  return true;
2130  }
2131 
2132 
2147  bool BundleAdjust::computeRejectionLimit() {
2148  double vx, vy;
2149 
2150  int numResiduals = m_bundleResults.numberObservations() / 2;
2151 
2152  std::vector<double> residuals;
2153 
2154  residuals.resize(numResiduals);
2155 
2156  // load magnitude (squared) of residual vector
2157  int residualIndex = 0;
2158  int numObjectPoints = m_bundleControlPoints.size();
2159  for (int i = 0; i < numObjectPoints; i++) {
2160 
2161  const BundleControlPointQsp point = m_bundleControlPoints.at(i);
2162 
2163  if ( point->isRejected() ) {
2164  continue;
2165  }
2166 
2167  int numMeasures = point->numberOfMeasures();
2168  for (int j = 0; j < numMeasures; j++) {
2169 
2170  const BundleMeasureQsp measure = point->at(j);
2171 
2172  if ( measure->isRejected() ) {
2173  continue;
2174  }
2175 
2176  vx = measure->sampleResidual();
2177  vy = measure->lineResidual();
2178 
2179  residuals[residualIndex] = sqrt(vx*vx + vy*vy);
2180 
2181  residualIndex++;
2182  }
2183  }
2184 
2185  // sort vectors
2186  std::sort(residuals.begin(), residuals.end());
2187 
2188  double median;
2189  double medianDev;
2190  double mad;
2191 
2192  int midpointIndex = numResiduals / 2;
2193 
2194  if ( numResiduals % 2 == 0 ) {
2195  median = (residuals[midpointIndex - 1] + residuals[midpointIndex]) / 2;
2196  }
2197  else {
2198  median = residuals[midpointIndex];
2199  }
2200 
2201  // compute M.A.D.
2202  for (int i = 0; i < numResiduals; i++) {
2203  residuals[i] = fabs(residuals[i] - median);
2204  }
2205 
2206  std::sort(residuals.begin(), residuals.end());
2207 
2208  if ( numResiduals % 2 == 0 ) {
2209  medianDev = (residuals[midpointIndex - 1] + residuals[midpointIndex]) / 2;
2210  }
2211  else {
2212  medianDev = residuals[midpointIndex];
2213  }
2214 
2215  QString status = "\nmedian deviation: ";
2216  status.append(QString("%1").arg(medianDev));
2217  status.append("\n");
2218  outputBundleStatus(status);
2219 
2220  mad = 1.4826 * medianDev;
2221 
2222  status = "\nmad: ";
2223  status.append(QString("%1").arg(mad));
2224  status.append("\n");
2225  outputBundleStatus(status);
2226 
2227  m_bundleResults.setRejectionLimit(median
2228  + m_bundleSettings->outlierRejectionMultiplier() * mad);
2229 
2230  status = "\nRejection Limit: ";
2231  status.append(QString("%1").arg(m_bundleResults.rejectionLimit()));
2232  status.append("\n");
2233  outputBundleStatus(status);
2234 
2235  return true;
2236  }
2237 
2238 
2246  bool BundleAdjust::flagOutliers() {
2247  double vx, vy;
2248  int numRejected;
2249  int totalNumRejected = 0;
2250 
2251  int maxResidualIndex;
2252  double maxResidual;
2253  double sumSquares;
2254  double usedRejectionLimit = m_bundleResults.rejectionLimit();
2255 
2256  // TODO What to do if usedRejectionLimit is too low?
2257 
2258  int numComingBack = 0;
2259 
2260  int numObjectPoints = m_bundleControlPoints.size();
2261 
2262  outputBundleStatus("\n");
2263  for (int i = 0; i < numObjectPoints; i++) {
2264  BundleControlPointQsp point = m_bundleControlPoints.at(i);
2265 
2266  point->zeroNumberOfRejectedMeasures();
2267 
2268  numRejected = 0;
2269  maxResidualIndex = -1;
2270  maxResidual = -1.0;
2271 
2272  int numMeasures = point->numberOfMeasures();
2273  for (int j = 0; j < numMeasures; j++) {
2274 
2275  BundleMeasureQsp measure = point->at(j);
2276 
2277  vx = measure->sampleResidual();
2278  vy = measure->lineResidual();
2279 
2280  sumSquares = sqrt(vx*vx + vy*vy);
2281 
2282  // measure is good
2283  if ( sumSquares <= usedRejectionLimit ) {
2284 
2285  // was it previously rejected?
2286  if ( measure->isRejected() ) {
2287  QString status = "Coming back in: ";
2288  status.append(QString("%1").arg(point->id().toLatin1().data()));
2289  status.append("\r");
2290  outputBundleStatus(status);
2291  numComingBack++;
2292  m_controlNet->DecrementNumberOfRejectedMeasuresInImage(measure->cubeSerialNumber());
2293  }
2294 
2295  measure->setRejected(false);
2296  continue;
2297  }
2298 
2299  // if it's still rejected, skip it
2300  if ( measure->isRejected() ) {
2301  numRejected++;
2302  totalNumRejected++;
2303  continue;
2304  }
2305 
2306  if ( sumSquares > maxResidual ) {
2307  maxResidual = sumSquares;
2308  maxResidualIndex = j;
2309  }
2310  }
2311 
2312  // no observations above the current rejection limit for this 3D point
2313  if ( maxResidual == -1.0 || maxResidual <= usedRejectionLimit ) {
2314  point->setNumberOfRejectedMeasures(numRejected);
2315  continue;
2316  }
2317 
2318  // this is another kluge - if we only have two observations
2319  // we won't reject (for now)
2320  if ((numMeasures - (numRejected + 1)) < 2) {
2321  point->setNumberOfRejectedMeasures(numRejected);
2322  continue;
2323  }
2324 
2325  // otherwise, we have at least one observation
2326  // for this point whose residual is above the
2327  // current rejection limit - we'll flag the
2328  // worst of these as rejected
2329  BundleMeasureQsp rejected = point->at(maxResidualIndex);
2330  rejected->setRejected(true);
2331  numRejected++;
2332  point->setNumberOfRejectedMeasures(numRejected);
2333  m_controlNet->IncrementNumberOfRejectedMeasuresInImage(rejected->cubeSerialNumber());
2334  totalNumRejected++;
2335 
2336  // do we still have sufficient remaining observations for this 3D point?
2337  if ( ( numMeasures-numRejected ) < 2 ) {
2338  point->setRejected(true);
2339  QString status = "Rejecting Entire Point: ";
2340  status.append(QString("%1").arg(point->id().toLatin1().data()));
2341  status.append("\r");
2342  outputBundleStatus(status);
2343  }
2344  else
2345  point->setRejected(false);
2346 
2347  }
2348 
2349  int numberRejectedObservations = 2*totalNumRejected;
2350 
2351  QString status = "\nRejected Observations:";
2352  status.append(QString("%1").arg(numberRejectedObservations));
2353  status.append(" (Rejection Limit:");
2354  status.append(QString("%1").arg(usedRejectionLimit));
2355  status.append(")\n");
2356  outputBundleStatus(status);
2357 
2358  m_bundleResults.setNumberRejectedObservations(numberRejectedObservations);
2359 
2360  status = "\nMeasures that came back: ";
2361  status.append(QString("%1").arg(numComingBack));
2362  status.append("\n");
2363  outputBundleStatus(status);
2364 
2365  return true;
2366  }
2367 
2368 
2376  QList<ImageList *> BundleAdjust::imageLists() {
2377 
2378  if (m_imageLists.count() > 0) {
2379  return m_imageLists;
2380  }
2381  else if (m_serialNumberList->size() > 0) {
2382  ImageList *imgList = new ImageList;
2383  try {
2384  for (int i = 0; i < m_serialNumberList->size(); i++) {
2385  Image *image = new Image(m_serialNumberList->fileName(i));
2386  imgList->append(image);
2387  image->closeCube();
2388  }
2389  m_imageLists.append(imgList);
2390  }
2391  catch (IException &e) {
2392  QString msg = "Invalid image in serial number list\n";
2393  throw IException(IException::Programmer, msg, _FILEINFO_);
2394  }
2395  }
2396  else {
2397  QString msg = "No images used in bundle adjust\n";
2398  throw IException(IException::Programmer, msg, _FILEINFO_);
2399  }
2400 
2401  return m_imageLists;
2402  }
2403 
2404 
2425  bool BundleAdjust::errorPropagation() {
2426  emit(statusBarUpdate("Error Propagation"));
2427  // free unneeded memory
2428  cholmod_free_triplet(&m_cholmodTriplet, &m_cholmodCommon);
2429  cholmod_free_sparse(&m_cholmodNormal, &m_cholmodCommon);
2430 
2431  LinearAlgebra::Matrix T(3, 3);
2432  // *** TODO ***
2433  // Can any of the control point specific code be moved to BundleControlPoint?
2434 
2435  double sigma0Squared = m_bundleResults.sigma0() * m_bundleResults.sigma0();
2436 
2437  int numObjectPoints = m_bundleControlPoints.size();
2438 
2439  std::string currentTime = iTime::CurrentLocalTime().toLatin1().data();
2440 
2441  QString status = " Time: ";
2442  status.append(currentTime.c_str());
2443  status.append("\n\n");
2444  outputBundleStatus(status);
2445 
2446  // create and initialize array of 3x3 matrices for all object points
2447  std::vector< symmetric_matrix<double> > pointCovariances(numObjectPoints,
2448  symmetric_matrix<double>(3));
2449  for (int d = 0; d < numObjectPoints; d++) {
2450  pointCovariances[d].clear();
2451  }
2452 
2453  cholmod_dense *x; // solution vector
2454  cholmod_dense *b; // right-hand side (column vectors of identity)
2455 
2456  b = cholmod_zeros ( m_rank, 1, CHOLMOD_REAL, &m_cholmodCommon );
2457  double *pb = (double*)b->x;
2458 
2459  double *px = NULL;
2460 
2461  SparseBlockColumnMatrix inverseMatrix;
2462 
2463  // Create unique file name
2464  FileName matrixFile(m_bundleSettings->outputFilePrefix() + "inverseMatrix.dat");
2465  //???FileName matrixFile = FileName::createTempFile(m_bundleSettings.outputFilePrefix()
2466  //??? + "inverseMatrix.dat");
2467  // Create file handle
2468  QFile matrixOutput(matrixFile.expanded());
2469 
2470  // Check to see if creating the inverse correlation matrix is turned on
2471  if (m_bundleSettings->createInverseMatrix()) {
2472  // Open file to write to
2473  matrixOutput.open(QIODevice::WriteOnly);
2474  }
2475  QDataStream outStream(&matrixOutput);
2476 
2477  int i, j, k;
2478  int columnIndex = 0;
2479  int numColumns = 0;
2480  int numBlockColumns = m_sparseNormals.size();
2481  for (i = 0; i < numBlockColumns; i++) {
2482 
2483  // columns in this column block
2484  SparseBlockColumnMatrix *normalsColumn = m_sparseNormals.at(i);
2485  if (i == 0) {
2486  numColumns = normalsColumn->numberOfColumns();
2487  int numRows = normalsColumn->numberOfRows();
2488  inverseMatrix.insertMatrixBlock(i, numRows, numColumns);
2489  inverseMatrix.zeroBlocks();
2490  }
2491  else {
2492  if (normalsColumn->numberOfColumns() == numColumns) {
2493  int numRows = normalsColumn->numberOfRows();
2494  inverseMatrix.insertMatrixBlock(i, numRows, numColumns);
2495  inverseMatrix.zeroBlocks();
2496  }
2497  else {
2498  numColumns = normalsColumn->numberOfColumns();
2499 
2500  // reset inverseMatrix
2501  inverseMatrix.wipe();
2502 
2503  // insert blocks
2504  for (j = 0; j < (i+1); j++) {
2505  SparseBlockColumnMatrix *normalsRow = m_sparseNormals.at(j);
2506  int numRows = normalsRow->numberOfRows();
2507 
2508  inverseMatrix.insertMatrixBlock(j, numRows, numColumns);
2509  }
2510  }
2511  }
2512 
2513  int localCol = 0;
2514 
2515  // solve for inverse for nCols
2516  for (j = 0; j < numColumns; j++) {
2517  if ( columnIndex > 0 ) {
2518  pb[columnIndex - 1] = 0.0;
2519  }
2520  pb[columnIndex] = 1.0;
2521 
2522  x = cholmod_solve ( CHOLMOD_A, m_L, b, &m_cholmodCommon );
2523  px = (double*)x->x;
2524  int rp = 0;
2525 
2526  // store solution in corresponding column of inverse
2527  for (k = 0; k < inverseMatrix.size(); k++) {
2528  LinearAlgebra::Matrix *matrix = inverseMatrix.value(k);
2529 
2530  int sz1 = matrix->size1();
2531 
2532  for (int ii = 0; ii < sz1; ii++) {
2533  (*matrix)(ii,localCol) = px[ii + rp];
2534  }
2535  rp += matrix->size1();
2536  }
2537 
2538  columnIndex++;
2539  localCol++;
2540 
2541  cholmod_free_dense(&x,&m_cholmodCommon);
2542  }
2543 
2544  // save adjusted target body sigmas if solving for target
2545  if (m_bundleSettings->solveTargetBody() && i == 0) {
2546  vector< double > &adjustedSigmas = m_bundleTargetBody->adjustedSigmas();
2547  matrix< double > *targetCovMatrix = inverseMatrix.value(i);
2548 
2549  for (int z = 0; z < numColumns; z++)
2550  adjustedSigmas[z] = sqrt((*targetCovMatrix)(z,z))*m_bundleResults.sigma0();
2551  }
2552  // save adjusted image sigmas
2553  else {
2554  BundleObservationQsp observation;
2555  if (m_bundleSettings->solveTargetBody()) {
2556  observation = m_bundleObservations.at(i-1);
2557  }
2558  else {
2559  observation = m_bundleObservations.at(i);
2560  }
2561  vector< double > &adjustedSigmas = observation->adjustedSigmas();
2562  matrix< double > *imageCovMatrix = inverseMatrix.value(i);
2563  for ( int z = 0; z < numColumns; z++) {
2564  adjustedSigmas[z] = sqrt((*imageCovMatrix)(z,z))*m_bundleResults.sigma0();
2565  }
2566  }
2567 
2568  // Output the inverse matrix if requested
2569  if (m_bundleSettings->createInverseMatrix()) {
2570  outStream << inverseMatrix;
2571  }
2572 
2573  // now loop over all object points to sum contributions into 3x3 point covariance matrix
2574  int pointIndex = 0;
2575  for (j = 0; j < numObjectPoints; j++) {
2576  emit(pointUpdate(j+1));
2577  BundleControlPointQsp point = m_bundleControlPoints.at(pointIndex);
2578  if ( point->isRejected() ) {
2579  continue;
2580  }
2581 
2582  // only update point every 100 points
2583  if (j%100 == 0) {
2584  QString status = "\rError Propagation: Inverse Block ";
2585  status.append(QString::number(i+1));
2586  status.append(" of ");
2587  status.append(QString::number(numBlockColumns));
2588  status.append("; Point ");
2589  status.append(QString::number(j+1));
2590  status.append(" of ");
2591  status.append(QString::number(numObjectPoints));
2592  outputBundleStatus(status);
2593  }
2594 
2595  // get corresponding Q matrix
2596  // NOTE: we are getting a reference to the Q matrix stored
2597  // in the BundleControlPoint for speed (without the & it is dirt slow)
2598  SparseBlockRowMatrix &Q = point->cholmodQMatrix();
2599 
2600  T.clear();
2601 
2602  // get corresponding point covariance matrix
2603  boost::numeric::ublas::symmetric_matrix<double> &covariance = pointCovariances[pointIndex];
2604 
2605  // get firstQBlock - index i is the key into Q for firstQBlock
2606  LinearAlgebra::Matrix *firstQBlock = Q.value(i);
2607  if (!firstQBlock) {
2608  pointIndex++;
2609  continue;
2610  }
2611 
2612  // iterate over Q
2613  // secondQBlock is current map value
2614  QMapIterator< int, LinearAlgebra::Matrix * > it(Q);
2615  while ( it.hasNext() ) {
2616  it.next();
2617 
2618  int nKey = it.key();
2619 
2620  if (it.key() > i) {
2621  break;
2622  }
2623 
2624  LinearAlgebra::Matrix *secondQBlock = it.value();
2625 
2626  if ( !secondQBlock ) {// should never be NULL
2627  continue;
2628  }
2629 
2630  LinearAlgebra::Matrix *inverseBlock = inverseMatrix.value(it.key());
2631 
2632  if ( !inverseBlock ) {// should never be NULL
2633  continue;
2634  }
2635 
2636  T = prod(*inverseBlock, trans(*firstQBlock));
2637  T = prod(*secondQBlock,T);
2638 
2639  if (nKey != i) {
2640  T += trans(T);
2641  }
2642 
2643  try {
2644  covariance += T;
2645  }
2646 
2647  catch (std::exception &e) {
2648  outputBundleStatus("\n\n");
2649  QString msg = "Input data and settings are not sufficiently stable "
2650  "for error propagation.";
2651  throw IException(IException::User, msg, _FILEINFO_);
2652  }
2653  }
2654  pointIndex++;
2655  }
2656  }
2657 
2658  if (m_bundleSettings->createInverseMatrix()) {
2659  // Close the file.
2660  matrixOutput.close();
2661  // Save the location of the "covariance" matrix
2662  m_bundleResults.setCorrMatCovFileName(matrixFile);
2663  }
2664 
2665  // can free sparse normals now
2666  m_sparseNormals.wipe();
2667 
2668  // free b (right-hand side vector
2669  cholmod_free_dense(&b,&m_cholmodCommon);
2670 
2671  outputBundleStatus("\n\n");
2672 
2673  currentTime = Isis::iTime::CurrentLocalTime().toLatin1().data();
2674 
2675  status = "\rFilling point covariance matrices: Time ";
2676  status.append(currentTime.c_str());
2677  outputBundleStatus(status);
2678  outputBundleStatus("\n\n");
2679 
2680  // now loop over points again and set final covariance stuff
2681  // *** TODO *** Can this loop go into BundleControlPoint
2682  int pointIndex = 0;
2683  for (j = 0; j < numObjectPoints; j++) {
2684 
2685  BundleControlPointQsp point = m_bundleControlPoints.at(pointIndex);
2686 
2687  if ( point->isRejected() ) {
2688  continue;
2689  }
2690 
2691  if (j%100 == 0) {
2692  status = "\rError Propagation: Filling point covariance matrices ";
2693  status.append(QString("%1").arg(j+1));
2694  status.append(" of ");
2695  status.append(QString("%1").arg(numObjectPoints));
2696  status.append("\r");
2697  outputBundleStatus(status);
2698  }
2699 
2700  // get corresponding point covariance matrix
2701  boost::numeric::ublas::symmetric_matrix<double> &covariance = pointCovariances[pointIndex];
2702 
2703  // Update and reset the matrix
2704  // Get the Limiting Error Propagation uncertainties: sigmas for coordinate 1, 2, and 3 in meters
2705  //
2706  SurfacePoint SurfacePoint = point->adjustedSurfacePoint();
2707 
2708  // Get the TEP by adding the corresponding members of pCovar and covariance
2709  boost::numeric::ublas::symmetric_matrix <double,boost::numeric::ublas::upper> pCovar;
2710 
2711  if (m_bundleSettings->controlPointCoordTypeBundle() == SurfacePoint::Latitudinal) {
2712  pCovar = SurfacePoint.GetSphericalMatrix(SurfacePoint::Kilometers);
2713  }
2714  else {
2715  // Assume Rectangular coordinates
2716  pCovar = SurfacePoint.GetRectangularMatrix(SurfacePoint::Kilometers);
2717  }
2718  pCovar += covariance;
2719  pCovar *= sigma0Squared;
2720 
2721  // debug lines
2722  // if (j < 3) {
2723  // std::cout << " Adjusted surface point ..." << std::endl;
2724  // std:: cout << " sigmaLat (radians) = " << sqrt(pCovar(0,0)) << std::endl;
2725  // std:: cout << " sigmaLon (radians) = " << sqrt(pCovar(1,1)) << std::endl;
2726  // std:: cout << " sigmaRad (km) = " << sqrt(pCovar(2,2)) << std::endl;
2727  // std::cout << " Adjusted matrix = " << std::endl;
2728  // std::cout << " " << pCovar(0,0) << " " << pCovar(0,1) << " "
2729  // << pCovar(0,2) << std::endl;
2730  // std::cout << " " << pCovar(1,0) << " " << pCovar(1,1) << " "
2731  // << pCovar(1,2) << std::endl;
2732  // std::cout << " " << pCovar(2,0) << " " << pCovar(2,1) << " "
2733  // << pCovar(2,2) << std::endl;
2734  // }
2735  // end debug
2736 
2737  // Distance units are km**2
2738  SurfacePoint.SetMatrix(m_bundleSettings->controlPointCoordTypeBundle(),pCovar);
2739  point->setAdjustedSurfacePoint(SurfacePoint);
2740  // // debug lines
2741  // if (j < 3) {
2742  // boost::numeric::ublas::symmetric_matrix <double,boost::numeric::ublas::upper> recCovar;
2743  // recCovar = SurfacePoint.GetRectangularMatrix(SurfacePoint::Meters);
2744  // std:: cout << " sigmaLat (meters) = " <<
2745  // point->adjustedSurfacePoint().GetSigmaDistance(SurfacePoint::Latitudinal,
2746  // SurfacePoint::One).meters() << std::endl;
2747  // std:: cout << " sigmaLon (meters) = " <<
2748  // point->adjustedSurfacePoint().GetSigmaDistance(SurfacePoint::Latitudinal,
2749  // SurfacePoint::Two).meters() << std::endl;
2750  // std:: cout << " sigmaRad (km) = " << sqrt(pCovar(2,2)) << std::endl;
2751  // std::cout << "Rectangular matrix with radius in meters" << std::endl;
2752  // std::cout << " " << recCovar(0,0) << " " << recCovar(0,1) << " "
2753  // << recCovar(0,2) << std::endl;
2754  // std::cout << " " << recCovar(1,0) << " " << recCovar(1,1) << " "
2755  // << recCovar(1,2) << std::endl;
2756  // std::cout << " " << recCovar(2,0) << " " << recCovar(2,1) << " "
2757  // << recCovar(2,2) << std::endl;
2758  // }
2759  // // end debug
2760 
2761  pointIndex++;
2762  }
2763 
2764  return true;
2765  }
2766 
2767 
2773  ControlNetQsp BundleAdjust::controlNet() {
2774  return m_controlNet;
2775  }
2776 
2777 
2783  SerialNumberList *BundleAdjust::serialNumberList() {
2784  return m_serialNumberList;
2785  }
2786 
2787 
2793  int BundleAdjust::numberOfImages() const {
2794  return m_serialNumberList->size();
2795  }
2796 
2797 
2805  // TODO: probably don't need this, can get from BundleObservation
2806  QString BundleAdjust::fileName(int i) {
2807  return m_serialNumberList->fileName(i);
2808  }
2809 
2810 
2816  double BundleAdjust::iteration() const {
2817  return m_iteration;
2818  }
2819 
2820 
2831  Table BundleAdjust::cMatrix(int i) {
2832  return m_controlNet->Camera(i)->instrumentRotation()->Cache("InstrumentPointing");
2833  }
2834 
2835 
2846  Table BundleAdjust::spVector(int i) {
2847  return m_controlNet->Camera(i)->instrumentPosition()->Cache("InstrumentPosition");
2848  }
2849 
2850 
2859  QString BundleAdjust::modelState(int i) {
2860  Camera *imageCam = m_controlNet->Camera(i);
2861  if (imageCam->GetCameraType() != Camera::Csm) {
2862  QString msg = "Cannot get model state for image [" + toString(i) +
2863  "] because it is not a CSM camera model.";
2864  throw IException(IException::Programmer, msg, _FILEINFO_);
2865  }
2866 
2867  CSMCamera *csmCamera = dynamic_cast<CSMCamera*>(imageCam);
2868  return csmCamera->getModelState();
2869  }
2870 
2871 
2880  void BundleAdjust::iterationSummary() {
2881  QString iterationNumber;
2882 
2883  if ( m_bundleResults.converged() ) {
2884  iterationNumber = "Iteration" + toString(m_iteration) + ": Final";
2885  }
2886  else {
2887  iterationNumber = "Iteration" + toString(m_iteration);
2888  }
2889 
2890  PvlGroup summaryGroup(iterationNumber);
2891 
2892  summaryGroup += PvlKeyword("Sigma0",
2893  toString( m_bundleResults.sigma0() ) );
2894  summaryGroup += PvlKeyword("Observations",
2895  toString( m_bundleResults.numberObservations() ) );
2896  summaryGroup += PvlKeyword("Constrained_Point_Parameters",
2897  toString( m_bundleResults.numberConstrainedPointParameters() ) );
2898  summaryGroup += PvlKeyword("Constrained_Image_Parameters",
2899  toString( m_bundleResults.numberConstrainedImageParameters() ) );
2900  if (m_bundleSettings->bundleTargetBody()) {
2901  summaryGroup += PvlKeyword("Constrained_Target_Parameters",
2902  toString( m_bundleResults.numberConstrainedTargetParameters() ) );
2903  }
2904  summaryGroup += PvlKeyword("Unknown_Parameters",
2905  toString( m_bundleResults.numberUnknownParameters() ) );
2906  summaryGroup += PvlKeyword("Degrees_of_Freedom",
2907  toString( m_bundleResults.degreesOfFreedom() ) );
2908  summaryGroup += PvlKeyword( "Rejected_Measures",
2909  toString( m_bundleResults.numberRejectedObservations()/2) );
2910 
2911  if ( m_bundleResults.numberMaximumLikelihoodModels() >
2912  m_bundleResults.maximumLikelihoodModelIndex() ) {
2913  // if maximum likelihood estimation is being used
2914 
2915  summaryGroup += PvlKeyword("Maximum_Likelihood_Tier: ",
2916  toString( m_bundleResults.maximumLikelihoodModelIndex() ) );
2917  summaryGroup += PvlKeyword("Median_of_R^2_residuals: ",
2918  toString( m_bundleResults.maximumLikelihoodMedianR2Residuals() ) );
2919  }
2920 
2921  if ( m_bundleResults.converged() ) {
2922  summaryGroup += PvlKeyword("Converged", "TRUE");
2923  summaryGroup += PvlKeyword("TotalElapsedTime", toString( m_bundleResults.elapsedTime() ) );
2924 
2925  if (m_bundleSettings->errorPropagation()) {
2926  summaryGroup += PvlKeyword("ErrorPropagationElapsedTime",
2927  toString( m_bundleResults.elapsedTimeErrorProp() ) );
2928  }
2929  }
2930 
2931  std::ostringstream ostr;
2932  ostr << summaryGroup << std::endl;
2933  m_iterationSummary += QString::fromStdString( ostr.str() );
2934 
2935  if (m_printSummary && iApp != NULL) {
2936  Application::Log(summaryGroup);
2937  }
2938  else {
2939  std::cout << summaryGroup << std::endl;
2940  }
2941  }
2942 
2943 
2949  bool BundleAdjust::isConverged() {
2950  return m_bundleResults.converged();
2951  }
2952 
2953 
2959  bool BundleAdjust::isAborted() {
2960  return m_abort;
2961  }
2962 
2963 
2971  QString BundleAdjust::iterationSummaryGroup() const {
2972  return m_iterationSummary;
2973  }
2974 
2975 
2985  void BundleAdjust::outputBundleStatus(QString status) {
2986  if (iApp == NULL) { // in a function call
2987  printf("%s", status.toStdString().c_str());
2988  }
2989  else if (QCoreApplication::applicationName() != "ipce") {
2990  printf("%s", status.toStdString().c_str());
2991  }
2992  }
2993 
2994 
2995 
3031  bool BundleAdjust::computeBundleStatistics() {
3032 
3033  // use qvectors so that we can set the size.
3034  // this will be useful later when adding data.
3035  // data may added out of index order
3036  int numberImages = m_serialNumberList->size();
3037  QVector<Statistics> rmsImageSampleResiduals(numberImages);
3038  QVector<Statistics> rmsImageLineResiduals(numberImages);
3039  QVector<Statistics> rmsImageResiduals(numberImages);
3040 
3041  int numObjectPoints = m_bundleControlPoints.size();
3042  for (int i = 0; i < numObjectPoints; i++) {
3043 
3044  const BundleControlPointQsp point = m_bundleControlPoints.at(i);
3045 
3046  if (point->isRejected()) {
3047  continue;
3048  }
3049 
3050  int numMeasures = point->numberOfMeasures();
3051  for (int j = 0; j < numMeasures; j++) {
3052 
3053  const BundleMeasureQsp measure = point->at(j);
3054 
3055  if (measure->isRejected()) {
3056  continue;
3057  }
3058 
3059  double sampleResidual = fabs(measure->sampleResidual());
3060  double lineResidual = fabs(measure->lineResidual());
3061 
3062  // Determine the index for this measure's serial number
3063  int imageIndex = m_serialNumberList->serialNumberIndex(measure->cubeSerialNumber());
3064 
3065  // add residual data to the statistics object at the appropriate serial number index
3066  rmsImageSampleResiduals[imageIndex].AddData(sampleResidual);
3067  rmsImageLineResiduals[imageIndex].AddData(lineResidual);
3068  rmsImageResiduals[imageIndex].AddData(lineResidual);
3069  rmsImageResiduals[imageIndex].AddData(sampleResidual);
3070  }
3071  }
3072 
3073 
3074  if (m_bundleSettings->errorPropagation()) {
3075 
3076  // initialize body-fixed coordinate boundaries
3077 
3078  // Latitude or X
3079  Distance minSigmaCoord1Dist;
3080  QString minSigmaCoord1PointId = "";
3081 
3082  Distance maxSigmaCoord1Dist;
3083  QString maxSigmaCoord1PointId = "";
3084 
3085  // Longitude or Y
3086  Distance minSigmaCoord2Dist;
3087  QString minSigmaCoord2PointId = "";
3088 
3089  Distance maxSigmaCoord2Dist;
3090  QString maxSigmaCoord2PointId = "";
3091 
3092  // Radius or Z
3093  Distance minSigmaCoord3Dist;
3094  QString minSigmaCoord3PointId = "";
3095 
3096  Distance maxSigmaCoord3Dist;
3097  QString maxSigmaCoord3PointId = "";
3098 
3099  // compute stats for point sigmas
3100  Statistics sigmaCoord1Stats;
3101  Statistics sigmaCoord2Stats;
3102  Statistics sigmaCoord3Stats;
3103 
3104  Distance sigmaCoord1Dist, sigmaCoord2Dist, sigmaCoord3Dist;
3105  SurfacePoint::CoordinateType type = m_bundleSettings->controlPointCoordTypeReports();
3106 
3107  int numPoints = m_bundleControlPoints.size();
3108  // initialize max and min values to those from first valid point
3109  for (int i = 0; i < numPoints; i++) {
3110 
3111  const BundleControlPointQsp point = m_bundleControlPoints.at(i);
3112 
3113  maxSigmaCoord1Dist = point->adjustedSurfacePoint().GetSigmaDistance(type,
3114  SurfacePoint::One);
3115  minSigmaCoord1Dist = maxSigmaCoord1Dist;
3116 
3117  maxSigmaCoord2Dist = point->adjustedSurfacePoint().GetSigmaDistance(type,
3118  SurfacePoint::Two);
3119  minSigmaCoord2Dist = maxSigmaCoord2Dist;
3120 
3121  maxSigmaCoord1PointId = point->id();
3122  maxSigmaCoord2PointId = maxSigmaCoord1PointId;
3123  minSigmaCoord1PointId = maxSigmaCoord1PointId;
3124  minSigmaCoord2PointId = maxSigmaCoord1PointId;
3125 
3126  // Get stats for coordinate 3 if used
3127  if (m_bundleSettings->solveRadius() || type == SurfacePoint::Rectangular) {
3128  maxSigmaCoord3Dist = point->adjustedSurfacePoint().GetSigmaDistance(type,
3129  SurfacePoint::Three);
3130  minSigmaCoord3Dist = maxSigmaCoord3Dist;
3131 
3132  maxSigmaCoord3PointId = maxSigmaCoord1PointId;
3133  minSigmaCoord3PointId = maxSigmaCoord1PointId;
3134  }
3135  break;
3136  }
3137 
3138  for (int i = 0; i < numPoints; i++) {
3139 
3140  const BundleControlPointQsp point = m_bundleControlPoints.at(i);
3141 
3142  sigmaCoord1Dist = point->adjustedSurfacePoint().GetSigmaDistance(type,
3143  SurfacePoint::One);
3144  sigmaCoord2Dist = point->adjustedSurfacePoint().GetSigmaDistance(type,
3145  SurfacePoint::Two);
3146  sigmaCoord3Dist = point->adjustedSurfacePoint().GetSigmaDistance(type,
3147  SurfacePoint::Three);
3148 
3149  sigmaCoord1Stats.AddData(sigmaCoord1Dist.meters());
3150  sigmaCoord2Stats.AddData(sigmaCoord2Dist.meters());
3151  sigmaCoord3Stats.AddData(sigmaCoord3Dist.meters());
3152 
3153  if (sigmaCoord1Dist > maxSigmaCoord1Dist) {
3154  maxSigmaCoord1Dist = sigmaCoord1Dist;
3155  maxSigmaCoord1PointId = point->id();
3156  }
3157  if (sigmaCoord2Dist > maxSigmaCoord2Dist) {
3158  maxSigmaCoord2Dist = sigmaCoord2Dist;
3159  maxSigmaCoord2PointId = point->id();
3160  }
3161  if (m_bundleSettings->solveRadius() || type == SurfacePoint::Rectangular) {
3162  if (sigmaCoord3Dist > maxSigmaCoord3Dist) {
3163  maxSigmaCoord3Dist = sigmaCoord3Dist;
3164  maxSigmaCoord3PointId = point->id();
3165  }
3166  }
3167  if (sigmaCoord1Dist < minSigmaCoord1Dist) {
3168  minSigmaCoord1Dist = sigmaCoord1Dist;
3169  minSigmaCoord1PointId = point->id();
3170  }
3171  if (sigmaCoord2Dist < minSigmaCoord2Dist) {
3172  minSigmaCoord2Dist = sigmaCoord2Dist;
3173  minSigmaCoord2PointId = point->id();
3174  }
3175  if (m_bundleSettings->solveRadius() || type == SurfacePoint::Rectangular) {
3176  if (sigmaCoord3Dist < minSigmaCoord3Dist) {
3177  minSigmaCoord3Dist = sigmaCoord3Dist;
3178  minSigmaCoord3PointId = point->id();
3179  }
3180  }
3181  }
3182 
3183  // update bundle results
3184  m_bundleResults.resizeSigmaStatisticsVectors(numberImages);
3185 
3186  m_bundleResults.setSigmaCoord1Range(minSigmaCoord1Dist, maxSigmaCoord1Dist,
3187  minSigmaCoord1PointId, maxSigmaCoord1PointId);
3188 
3189  m_bundleResults.setSigmaCoord2Range(minSigmaCoord2Dist, maxSigmaCoord2Dist,
3190  minSigmaCoord2PointId, maxSigmaCoord2PointId);
3191 
3192  m_bundleResults.setSigmaCoord3Range(minSigmaCoord3Dist, maxSigmaCoord3Dist,
3193  minSigmaCoord3PointId, maxSigmaCoord3PointId);
3194 
3195  m_bundleResults.setRmsFromSigmaStatistics(sigmaCoord1Stats.Rms(),
3196  sigmaCoord2Stats.Rms(),
3197  sigmaCoord3Stats.Rms());
3198  }
3199  m_bundleResults.setRmsImageResidualLists(rmsImageLineResiduals.toList(),
3200  rmsImageSampleResiduals.toList(),
3201  rmsImageResiduals.toList());
3202 
3203  return true;
3204  }
3205 
3206 }
Isis::Statistics
This class is used to accumulate statistics on double arrays.
Definition: Statistics.h:94
Isis::Statistics::AddData
void AddData(const double *data, const unsigned int count)
Add an array of doubles to the accumulators and counters.
Definition: Statistics.cpp:141
Isis::BundleMeasure::cubeSerialNumber
QString cubeSerialNumber() const
Accesses the serial number of the cube containing this control measure.
Definition: BundleMeasure.cpp:251
Isis::BundleControlPoint::adjustedSurfacePoint
SurfacePoint adjustedSurfacePoint() const
Accesses the adjusted SurfacePoint associated with this BundleControlPoint.
Definition: BundleControlPoint.cpp:479
Isis::SparseBlockColumnMatrix
SparseBlockColumnMatrix.
Definition: SparseBlockMatrix.h:58
Isis::PvlKeyword
A single keyword-value pair.
Definition: PvlKeyword.h:82
Isis::CameraGroundMap::GetXY
virtual bool GetXY(const SurfacePoint &spoint, double *cudx, double *cudy, bool test=true)
Compute undistorted focal plane coordinate from ground position using current Spice from SetImage cal...
Definition: CameraGroundMap.cpp:152
QList
This is free and unencumbered software released into the public domain.
Definition: BoxcarCachingAlgorithm.h:13
Isis::SerialNumberList::size
int size() const
How many serial number / filename combos are in the list.
Definition: SerialNumberList.cpp:384
Isis::Camera::SetImage
virtual bool SetImage(const double sample, const double line)
Sets the sample/line values of the image to get the lat/lon values.
Definition: Camera.cpp:154
Isis::FileName
File name manipulation and expansion.
Definition: FileName.h:100
Isis::Camera::GetCameraType
virtual CameraType GetCameraType() const =0
Returns the type of camera that was created.
Isis::SparseBlockRowMatrix::zeroBlocks
void zeroBlocks()
Sets all elements of all matrix blocks to zero.
Definition: SparseBlockMatrix.cpp:611
Isis::ImageList::append
void append(Image *const &value)
Appends an image to the image list.
Definition: ImageList.cpp:153
Isis::BundleMeasure::camera
Camera * camera() const
Accesses the associated camera for this bundle measure.
Definition: BundleMeasure.cpp:128
Isis::iTime::CurrentLocalTime
static QString CurrentLocalTime()
Returns the current local time This time is taken directly from the system clock, so if the system cl...
Definition: iTime.cpp:513
Isis::Image::fileName
QString fileName() const
Get the file name of the cube that this image represents.
Definition: Image.cpp:340
Isis::SurfacePoint::GetLatitude
Latitude GetLatitude() const
Return the body-fixed latitude for the surface point.
Definition: SurfacePoint.cpp:1665
Isis::Camera::Camera
Camera(Cube &cube)
Constructs the Camera object.
Definition: Camera.cpp:54
Isis::BundleControlPoint
This class holds information about a control point that BundleAdjust needs to run correctly.
Definition: BundleControlPoint.h:91
Isis::Control
This represents an ISIS control net in a project-based GUI interface.
Definition: Control.h:66
QSharedPointer< BundleSettings >
Isis::SurfacePoint::ResetLocalRadius
void ResetLocalRadius(const Distance &radius)
This method resets the local radius of a SurfacePoint.
Definition: SurfacePoint.cpp:914
Isis::SparseBlockColumnMatrix::numberOfRows
int numberOfRows()
Returns total number of rows in map (this needs to be clarified and maybe rewritten).
Definition: SparseBlockMatrix.cpp:215
Isis::Camera
Definition: Camera.h:236
Isis::LinearAlgebra::Matrix
boost::numeric::ublas::matrix< double > Matrix
Definition for an Isis::LinearAlgebra::Matrix of doubles.
Definition: LinearAlgebra.h:100
Isis::BundleControlPoint::id
QString id() const
Accesses the Point ID associated with this BundleControlPoint.
Definition: BundleControlPoint.cpp:489
Isis::toString
QString toString(bool boolToConvert)
Global function to convert a boolean to a string.
Definition: IString.cpp:211
Isis::SerialNumberList
Serial Number list generator.
Definition: SerialNumberList.h:64
Isis::SparseBlockColumnMatrix::startColumn
int startColumn() const
Sets starting column for block in full matrix.
Definition: SparseBlockMatrix.cpp:157
Isis::Distance
Distance measurement, usually in meters.
Definition: Distance.h:34
Isis::SparseBlockColumnMatrix::insertMatrixBlock
bool insertMatrixBlock(int nColumnBlock, int nRows, int nCols)
Inserts a "newed" LinearAlgebra::Matrix pointer of size (nRows, nCols) into the map with the block co...
Definition: SparseBlockMatrix.cpp:121
Isis::ImageList
Internalizes a list of images and allows for operations on the entire list.
Definition: ImageList.h:55
Isis::ControlPoint
A single control point.
Definition: ControlPoint.h:354
Isis::SurfacePoint::SetMatrix
void SetMatrix(CoordinateType type, const boost::numeric::ublas::symmetric_matrix< double, boost::numeric::ublas::upper > &covar)
Set the covariance matrix.
Definition: SurfacePoint.cpp:749
Isis::LinearAlgebra::Vector
boost::numeric::ublas::vector< double > Vector
Definition for an Isis::LinearAlgebra::Vector of doubles.
Definition: LinearAlgebra.h:120
Isis::Statistics::Rms
double Rms() const
Computes and returns the rms.
Definition: Statistics.cpp:365
Isis::PvlGroup
Contains multiple PvlContainers.
Definition: PvlGroup.h:41
Isis::SurfacePoint::CoordinateType
CoordinateType
Defines the coordinate typ, units, and coordinate index for some of the output methods.
Definition: SurfacePoint.h:139
Isis::BundleSolutionInfo
Container class for BundleAdjustment results.
Definition: BundleSolutionInfo.h:159
Isis::BundleMeasure::sample
double sample() const
Accesses the current sample measurement for this control measure.
Definition: BundleMeasure.cpp:191
Isis::Table
Class for storing Table blobs information.
Definition: Table.h:61
Isis::SerialNumberList::add
void add(const QString &filename, bool def2filename=false)
Adds a new filename / serial number pair to the SerialNumberList.
Definition: SerialNumberList.cpp:121
Isis::BundleSolutionInfo::setRunTime
void setRunTime(QString runTime)
Sets the run time, and the name if a name is not already set.
Definition: BundleSolutionInfo.cpp:242
Isis::Image::closeCube
void closeCube()
Cleans up the Cube pointer.
Definition: Image.cpp:307
Isis::ControlNet
a control network
Definition: ControlNet.h:257
Isis::BundleMeasure
A container class for a ControlMeasure.
Definition: BundleMeasure.h:55
Isis::Image
This represents a cube in a project-based GUI interface.
Definition: Image.h:107
Isis::IException
Isis exception class.
Definition: IException.h:91
Isis::Progress
Program progress reporter.
Definition: Progress.h:42
Isis::IException::errorType
ErrorType errorType() const
Returns the source of the error for this exception.
Definition: IException.cpp:430
Isis::SurfacePoint::GetLongitude
Longitude GetLongitude() const
Return the body-fixed longitude for the surface point.
Definition: SurfacePoint.cpp:1685
Isis::CSMCamera::getModelState
QString getModelState() const
Get the CSM Model state string to re-create the CSM Model.
Definition: CSMCamera.cpp:969
Isis::SparseBlockRowMatrix::insertMatrixBlock
bool insertMatrixBlock(int nRowBlock, int nRows, int nCols)
Inserts a "newed" LinearAlgebra::Matrix pointer of size (nRows, nCols) into the map with the block ro...
Definition: SparseBlockMatrix.cpp:507
Isis::SurfacePoint::GetLocalRadius
Distance GetLocalRadius() const
Return the radius of the surface point.
Definition: SurfacePoint.cpp:1732
Isis::Camera::GroundMap
CameraGroundMap * GroundMap()
Returns a pointer to the CameraGroundMap object.
Definition: Camera.cpp:2856
Isis::Distance::meters
double meters() const
Get the distance in meters.
Definition: Distance.cpp:85
Isis::BundleImage
This class hold image information that BundleAdjust needs to run correctly.Definition for a BundleIma...
Definition: BundleImage.h:36
Isis::BundleMeasure::parentBundleObservation
QSharedPointer< BundleObservation > parentBundleObservation()
Accesses the parent BundleObservation for this bundle measure.
Definition: BundleMeasure.cpp:158
Isis::SparseBlockColumnMatrix::numberOfColumns
int numberOfColumns()
Returns total number of columns in map (NOTE: NOT the number of matrix blocks).
Definition: SparseBlockMatrix.cpp:190
Isis::CSMCamera
Definition: CSMCamera.h:25
Isis::SparseBlockRowMatrix
SparseBlockRowMatrix.
Definition: SparseBlockMatrix.h:125
QVector
This is free and unencumbered software released into the public domain.
Definition: Calculator.h:18
Isis::cholmodErrorHandler
static void cholmodErrorHandler(int nStatus, const char *file, int nLineNo, const char *message)
Custom error handler for CHOLMOD.
Definition: BundleAdjust.cpp:72
Isis::SparseBlockColumnMatrix::wipe
void wipe()
Deletes all pointer elements and removes them from the map.
Definition: SparseBlockMatrix.cpp:52
Isis::BundleMeasure::line
double line() const
Accesses the current line measurement for this control measure.
Definition: BundleMeasure.cpp:215
Isis::SurfacePoint
This class defines a body-fixed surface point.
Definition: SurfacePoint.h:132
Isis
This is free and unencumbered software released into the public domain.
Definition: Apollo.h:16
Isis::Control::fileName
QString fileName() const
Access the name of the control network file associated with this Control.
Definition: Control.cpp:272