Isis 3 Programmer Reference
BundleControlPoint.cpp
1 
7 /* SPDX-License-Identifier: CC0-1.0 */
8 
9 #include "BundleControlPoint.h"
10 
11 #include <QDebug>
12 
13 // boost lib
14 #include <boost/numeric/ublas/vector_proxy.hpp>
15 
16 // Isis lib
17 #include "ControlMeasure.h"
18 #include "Latitude.h"
19 #include "LinearAlgebra.h"
20 #include "Longitude.h"
21 #include "SparseBlockMatrix.h"
22 #include "SpecialPixel.h"
23 
24 namespace Isis {
25 
26 
34  ControlPoint *controlPoint) {
35  m_controlPoint = controlPoint;
36 
37  // setup vector of BundleMeasures for this control point
38  int numMeasures = controlPoint->GetNumMeasures();
39  for (int i = 0; i < numMeasures; i++) {
40  ControlMeasure *controlMeasure = controlPoint->GetMeasure(i);
41  if (controlMeasure->IsIgnored()) {
42  continue;
43  }
44 
45  addMeasure(controlMeasure);
46  }
47 
48  // initialize to 0.0
49  m_corrections.clear();
50  m_weights.clear();
51  m_nicVector.clear();
52 
53  // initialize to Null for consistency with other bundle classes...
54  m_aprioriSigmas.clear();
58  m_adjustedSigmas.clear();
62 
63  // initialize coordinate type from settings
64  // m_coordType = SurfacePoint::Latitudinal;
65  m_coordTypeReports = bundleSettings->controlPointCoordTypeReports();
66  m_coordTypeBundle = bundleSettings->controlPointCoordTypeBundle();
67  setWeights(bundleSettings);
68  }
69 
70 
78  copy(src);
79  }
80 
81 
86  }
87 
88 
95 
96  // sanity check
97  clear();
98 
99  m_controlPoint = src.m_controlPoint;
100 
101  int numMeasures = src.size();
102 
103  for (int i = 0; i < numMeasures; i++)
104  append(BundleMeasureQsp( new BundleMeasure(*(src.at(i))) ));
105 
109  m_weights = src.m_weights;
110  m_nicVector = src.m_nicVector;
112  m_coordTypeBundle = src.m_coordTypeBundle;
113  // *** TODO *** Why is the cholmodQMatrix not copied? Ask Ken
114  }
115 
116 
126 
127  BundleMeasureQsp bundleMeasure = BundleMeasureQsp( new BundleMeasure(controlMeasure, this) );
128 
129  append(bundleMeasure);
130 
131  return bundleMeasure;
132  }
133 
134 
141  m_controlPoint->ComputeResiduals();
142  }
143 
144 
151  m_controlPoint->SetAdjustedSurfacePoint(surfacePoint);
152  }
153 
154 
163  m_controlPoint->SetNumberOfRejectedMeasures(numRejected);
164  }
165 
166 
175  m_controlPoint->SetRejected(reject);
176  }
177 
178 
187 
188  double globalPointCoord1AprioriSigma = settings->globalPointCoord1AprioriSigma();
189  double globalPointCoord2AprioriSigma = settings->globalPointCoord2AprioriSigma();
190  double globalPointCoord3AprioriSigma = settings->globalPointCoord3AprioriSigma();
191 
192  if (m_controlPoint->GetType() == ControlPoint::Fixed) {
193  m_weights[0] = 1.0e+50;
194  m_weights[1] = 1.0e+50;
195  m_weights[2] = 1.0e+50;
196  // m_aprioriSigmas = Isis::Null by default
197  }
198 
199  if (m_controlPoint->GetType() == ControlPoint::Free) {
200  // Global sigmas are temporary and should not be stored in the control net. Currently
201  // global sigmas are always in meters. Make sure unit of weights is 1/(var unit squared),
202  // where var is a value being solved for in the adjustment. Latitude and longitude are in
203  // radians and everything else is in km.
204  if (!IsSpecial(globalPointCoord1AprioriSigma)) {
205  setSigmaWeightFromGlobals(globalPointCoord1AprioriSigma, 0);
206  } // else m_aprioriSigma = Isis::Null
207  // m_weights = 0.0
208  if (!IsSpecial(globalPointCoord2AprioriSigma)) {
209  m_aprioriSigmas[1] = globalPointCoord2AprioriSigma;
210  setSigmaWeightFromGlobals(globalPointCoord2AprioriSigma, 1);
211  }// else m_aprioriSigma = Isis::Null
212  // m_weights = 0.0
213  if (m_coordTypeBundle == SurfacePoint::Latitudinal && !settings->solveRadius()) {
214  m_weights[2] = 1.0e+50;
215  }
216  else {
217  if (!IsSpecial(globalPointCoord3AprioriSigma)) {
218  setSigmaWeightFromGlobals(globalPointCoord3AprioriSigma, 2);
219  }
220  }
221  }
222  if (m_controlPoint->GetType() == ControlPoint::Constrained) {
223  // *** Be careful...Is m_aprioriSigmas an output (for reports) or bundle variable? ***
224 
225  // Assuming the member variable sigmas are for output reports (internal use only) so use
226  // the report coordinate type to calculate. All point sigmas are in meters. Radius weights
227  // are in 1/km^2. Make x/y/z weights the same because BundleAdjust works in km.
228  // Weights and corrections go into the bundle so use bundle coordinate type.
229  if ( m_controlPoint->IsCoord1Constrained() ) {
230  m_aprioriSigmas[0] =
231  m_controlPoint->GetAprioriSurfacePoint().GetSigmaDistance
232  (m_coordTypeReports, SurfacePoint::One).meters();
233  m_weights[0] =
234  m_controlPoint->GetAprioriSurfacePoint().GetWeight(m_coordTypeBundle, SurfacePoint::One);
235  }
236  else if (!IsSpecial(globalPointCoord1AprioriSigma)) {
237  setSigmaWeightFromGlobals(globalPointCoord1AprioriSigma, 0);
238  } // else not constrained and global sigma is Null, then m_aprioriSigmas = Isis::Null
239  // m_weights = 0.0
240 
241  if ( m_controlPoint->IsCoord2Constrained() ) {
242  m_aprioriSigmas[1] =
243  m_controlPoint->GetAprioriSurfacePoint().GetSigmaDistance
244  (m_coordTypeReports, SurfacePoint::Two).meters();
245  m_weights[1] =
246  m_controlPoint->GetAprioriSurfacePoint().GetWeight(m_coordTypeBundle, SurfacePoint::Two);
247  }
248  else if (!IsSpecial(globalPointCoord2AprioriSigma)) {
249  setSigmaWeightFromGlobals(globalPointCoord2AprioriSigma, 1);
250  } // else not constrained and global sigma is Null, then m_aprioriSigmas = Isis::Null
251  // m_weights = 0.0
252 
253  if (m_coordTypeBundle == SurfacePoint::Latitudinal && !settings->solveRadius()) {
254  m_weights[2] = 1.0e+50;
255  }
256  else {
257  if ( m_controlPoint->IsCoord3Constrained() ) {
258  m_aprioriSigmas[2] =
259  m_controlPoint->GetAprioriSurfacePoint().GetSigmaDistance
260  (m_coordTypeReports, SurfacePoint::Three).meters();
261  m_weights[2] = m_controlPoint->GetAprioriSurfacePoint().GetWeight
262  (m_coordTypeBundle, SurfacePoint::Three);
263  }
264  else if (!IsSpecial(globalPointCoord3AprioriSigma)) {
265  setSigmaWeightFromGlobals(globalPointCoord3AprioriSigma, 2);
266  } // else not constrained and global sigma is Null, then m_aprioriSigmas = Isis::Null
267  // m_weights = 0.0
268  }
269  }
270  }
271 
272 
283  void BundleControlPoint::setSigmaWeightFromGlobals(double gSigma, int index) {
284  m_aprioriSigmas[index] = gSigma;
285 
286  switch (index) {
287  case 0:
288  if (m_coordTypeBundle == SurfacePoint::Latitudinal) {
289  m_aprioriSigmas[index] = gSigma;
290  double sigmaRadians =
291  m_controlPoint->GetAprioriSurfacePoint().MetersToLatitude(gSigma); // m to radians
292  if (sigmaRadians > DBL_EPSILON) {
293  m_weights[index] = 1. / (sigmaRadians*sigmaRadians); // 1/radians^2
294  }
295  }
296  else if (m_coordTypeBundle == SurfacePoint::Rectangular) {
297  double sigmakm = gSigma * .001; // km
298  if (sigmakm > DBL_EPSILON) {
299  m_weights[0] = 1./ (sigmakm*sigmakm); // 1/km^2
300  }
301  }
302  else {
303  IString msg ="Unknown surface point coordinate type enum ["
304  + toString(m_coordTypeBundle) + "]." ;
305  throw IException(IException::Programmer, msg, _FILEINFO_);
306  }
307  break;
308  case 1:
309  // if (!IsSpecial(globalPointCoord2AprioriSigma)) {
310  if (m_coordTypeBundle == SurfacePoint::Latitudinal) {
311  double sigmaRadians =
312  m_controlPoint->GetAprioriSurfacePoint().MetersToLongitude(gSigma); // m to radians
313  if (sigmaRadians > DBL_EPSILON) {
314  m_weights[1] = 1. / (sigmaRadians*sigmaRadians); // 1/radians^2
315  }
316  }
317  else if (m_coordTypeBundle == SurfacePoint::Rectangular) {
318  double sigmakm = gSigma * 0.001; // km
319  if (sigmakm > DBL_EPSILON) m_weights[1] = 1./ (sigmakm*sigmakm); // 1/km^2
320  }
321  else {
322  IString msg ="Unknown surface point coordinate type enum ["
323  + toString(m_coordTypeBundle) + "]." ;
324  throw IException(IException::Programmer, msg, _FILEINFO_);
325  }
326  break;
327 
328  case 2:
329  // Coordinate 2 is either latitudinal radius or rectangular z, both in m
330  {
331  double sigmakm = gSigma * .001; // km
332  m_weights[2] = 1./ (sigmakm*sigmakm); // 1/km^2
333  }
334  break;
335 
336  default:
337  IString msg ="Unknown coordinate index [" + toString(index) + "]." ;
338  throw IException(IException::Programmer, msg, _FILEINFO_);
339  }
340  }
341 
342 
349  m_controlPoint->ZeroNumberOfRejectedMeasures();
350  }
351 
352 
363  SparseBlockMatrix &sparseNormals,
364  LinearAlgebra::Vector &v1) {
365 
366  QMapIterator< int, LinearAlgebra::Matrix * > Qit(m_cholmodQMatrix);
367 
368  int subrangeStart, subrangeEnd;
369 
370  while ( Qit.hasNext() ) {
371  Qit.next();
372 
373  int columnIndex = Qit.key();
374 
375  subrangeStart = sparseNormals.at(columnIndex)->startColumn();
376  subrangeEnd = subrangeStart + Qit.value()->size2();
377 
378  m_nicVector += alpha * prod(*(Qit.value()),subrange(v1,subrangeStart,subrangeEnd));
379  }
380  }
381 
382 
391  SparseBlockMatrix &sparseNormals,
392  const BundleTargetBodyQsp target) {
393  if (!isRejected()) {
394 
395  // subtract product of Q and nj from NIC
396  productAlphaAV(-1.0, sparseNormals, imageSolution);
397 
398  // update adjusted surface point appropriately for the coordinate type
399  switch (m_coordTypeBundle) {
402  break;
405  break;
406  default:
407  IString msg ="Unknown surface point coordinate type enum ["
408  + toString(m_coordTypeBundle) + "]." ;
409  throw IException(IException::Programmer, msg, _FILEINFO_);
410  }
411 
412  // sum and save corrections
413  m_corrections(0) += m_nicVector[0];
414  m_corrections(1) += m_nicVector[1];
415  m_corrections(2) += m_nicVector[2];
416  }
417  }
418 
419 
426  return m_controlPoint;
427  }
428 
429 
436  return m_controlPoint->IsRejected();
437  }
438 
439 
446  return this->size();
447  }
448 
449 
458  return m_controlPoint->GetNumberOfRejectedMeasures();
459  }
460 
461 
470  return m_controlPoint->GetResidualRms();
471  }
472 
473 
480  return m_controlPoint->GetAdjustedSurfacePoint();
481  }
482 
483 
489  QString BundleControlPoint::id() const {
490  return m_controlPoint->GetId();
491  }
492 
493 
503  return m_controlPoint->GetType();
504  }
505 
506 
516  return m_coordTypeReports;
517  }
518 
519 
529  return m_coordTypeBundle;
530  }
531 
532 
533 // *** TODO *** Need to add bounded vectors to Linear Algebra class
541  boost::numeric::ublas::bounded_vector< double, 3 > &BundleControlPoint::corrections() {
542  return m_corrections;
543  }
544 
545 
552  boost::numeric::ublas::bounded_vector< double, 3 > &BundleControlPoint::aprioriSigmas() {
553  return m_aprioriSigmas;
554  // Be careful about units. Latitude and longitude sigmas are in radians. Radius and X/Y/Z are
555  // meters.
556 
557  }
558 
559 
566  boost::numeric::ublas::bounded_vector< double, 3 > &BundleControlPoint::adjustedSigmas() {
567  return m_adjustedSigmas;
568  }
569 
570 
577  boost::numeric::ublas::bounded_vector< double, 3 > &BundleControlPoint::weights() {
578  return m_weights;
579  }
580 
581 
587  boost::numeric::ublas::bounded_vector<double, 3> &BundleControlPoint::nicVector() {
588  return m_nicVector;
589  }
590 
591 
598  return m_cholmodQMatrix;
599  }
600 
601 
616  QString BundleControlPoint::formatBundleOutputSummaryString(bool errorPropagation) const {
617 
618  int numRays = numberOfMeasures(); // should this depend on the raw point, as written, or this->size()???
619  int numGoodRays = numRays - numberOfRejectedMeasures();
620  double ResidualRms = residualRms();
621 
622  // Return generic set of control point coordinates as floating point values
623  // Be careful because of the potential to confuse units.
624  SurfacePoint::CoordUnits units = SurfacePoint::Degrees;
625  if (m_coordTypeReports == SurfacePoint::Rectangular) units = SurfacePoint::Kilometers;
626  double cpc1 = m_controlPoint->GetAdjustedSurfacePoint().GetCoord(m_coordTypeReports,
627  SurfacePoint::One, units);
628  double cpc2 = m_controlPoint->GetAdjustedSurfacePoint().GetCoord(m_coordTypeReports,
629  SurfacePoint::Two, units);
630  double cpc3 = m_controlPoint->GetAdjustedSurfacePoint().GetCoord(m_coordTypeReports,
631  SurfacePoint::Three, SurfacePoint::Kilometers);
632 
633  QString pointType = ControlPoint::PointTypeToString(type()).toUpper();
634 
635  QString output = QString("%1%2%3 of %4%5%6%7%8%9%10%11\n")
636  .arg(id(), 16)
637  .arg(pointType, 15)
638  .arg(numGoodRays, 5)
639  .arg(numRays)
640  .arg(formatValue(ResidualRms, 6, 2))
641  .arg(formatValue(cpc1, 16, 8)) // deg km
642  .arg(formatValue(cpc2, 16, 8)) // deg OR km
643  .arg(formatValue(cpc3, 16, 8)) // km km
644  .arg(formatCoordAdjustedSigmaString(SurfacePoint::One, 16, 8, errorPropagation)) // m
645  .arg(formatCoordAdjustedSigmaString(SurfacePoint::Two, 16, 8, errorPropagation)) // m
646  .arg(formatCoordAdjustedSigmaString(SurfacePoint::Three, 16, 8, errorPropagation)); // m
647 
648  return output;
649  }
650 
651 
664  bool solveRadius) const {
665  QString output;
666 
667  switch (m_coordTypeReports) {
669  output = formatBundleLatitudinalOutputDetailString(errorPropagation, solveRadius);
670  break;
672  output = formatBundleRectangularOutputDetailString(errorPropagation);
673  break;
674  default:
675  IString msg ="Unknown surface point coordinate type enum [" + toString(m_coordTypeBundle) + "]." ;
676  throw IException(IException::Programmer, msg, _FILEINFO_);
677  }
678  return output;
679  }
680 
681 
703  bool errorPropagation,
704  bool solveRadius) const {
705 
706  int numRays = numberOfMeasures();
707  int numGoodRays = numRays - numberOfRejectedMeasures();
708  double lat = m_controlPoint->GetAdjustedSurfacePoint().GetLatitude().degrees();
709  double lon = m_controlPoint->GetAdjustedSurfacePoint().GetLongitude().degrees();
710  double rad = m_controlPoint->GetAdjustedSurfacePoint().GetLocalRadius().kilometers();
711 
712  // Use the local radius in meters, rad*1000., to convert radians to meters now instead of the
713  // target body equatorial radius DAC 09/17/2018
714  double rtm = rad * 1000.;
715 
716  // Coordinate corrections are currently set in BundleAdjust::applyParameterCorrections in the
717  // coordTypeBundle coordinates (radians for Latitudinal coordinates and km for Rectangular).
718  // point corrections and initial sigmas
719  double cor_lat_dd = 0.; // lat correction, decimal degs
720  double cor_lon_dd = 0.; // lon correction, decimal degs
721  double cor_rad_km = 0.; // radius correction, kilometers
722  double cor_lat_m = 0.; // lat correction, meters
723  double cor_lon_m = 0.; // lon correction, meters
724  double cor_rad_m = 0.; // radius correction, meters
725  double latInit = Isis::Null;
726  double lonInit = Isis::Null;
727  double radInit = Isis::Null;
728 
729  if (m_coordTypeBundle == SurfacePoint::Rectangular) {
730  double x = m_controlPoint->GetAdjustedSurfacePoint().GetX().kilometers();
731  double xCor = m_corrections(0); // km
732  double y = m_controlPoint->GetAdjustedSurfacePoint().GetY().kilometers();
733  double yCor = m_corrections(1); // km
734  double z = m_controlPoint->GetAdjustedSurfacePoint().GetZ().kilometers();
735  double zCor = m_corrections(2); // km
736 
737  if (!IsSpecial(x) && !IsSpecial(y) && !IsSpecial(z)) {
741  latInit = rectPoint.GetLatitude().degrees();
742  lonInit = rectPoint.GetLongitude().degrees();
743  radInit = rectPoint.GetLocalRadius().kilometers();
744  if (!IsSpecial(lat)) {
745  cor_lat_dd = (lat - latInit); // degrees
746  cor_lat_m = cor_lat_dd * DEG2RAD * rtm;
747  }
748  if (!IsSpecial(lon)) {
749  cor_lon_dd = (lon - lonInit); // degrees
750  cor_lon_m = cor_lon_dd * DEG2RAD * rtm * cos(lat*DEG2RAD); // lon corrections meters
751  }
752  if (!IsSpecial(rad)) {
753  cor_rad_km = rad - radInit;
754  cor_rad_m = cor_rad_km * 1000.;
755  }
756  }
757  }
758  else if (m_coordTypeBundle == SurfacePoint::Latitudinal) {
759  cor_lat_dd = m_corrections(0) * RAD2DEG; // lat correction, decimal degs
760  cor_lon_dd = m_corrections(1) * RAD2DEG; // lon correction, decimal degs
761  cor_rad_m = m_corrections(2) * 1000.0; // radius correction, meters
762 
763  cor_lat_m = m_controlPoint->GetAdjustedSurfacePoint().LatitudeToMeters(m_corrections(0));
764  cor_lon_m = m_controlPoint->GetAdjustedSurfacePoint().LongitudeToMeters(m_corrections(1));
765  cor_rad_km = m_corrections(2);
766 
767  if (!IsSpecial(lat)) {
768  latInit = lat - cor_lat_dd;
769  }
770 
771  if (!IsSpecial(lon)) {
772  lonInit = lon - cor_lon_dd;
773  }
774 
775  if (!IsSpecial(rad)) {
776  radInit = rad - m_corrections(2); // km
777  }
778  }
779  else {
780  IString msg ="Unknown surface point coordinate type enum [" +
781  toString(m_coordTypeBundle) + "]." ;
782  throw IException(IException::Programmer, msg, _FILEINFO_);
783  }
784 
785  QString pointType = ControlPoint::PointTypeToString(type()).toUpper();
786 
787  QString output;
788 
789  output = QString(" Label: %1\nStatus: %2\n Rays: %3 of %4\n")
790  .arg(id())
791  .arg(pointType)
792  .arg(numGoodRays)
793  .arg(numRays);
794 
795  QString labels = "\n Point Initial Total Total "
796  " Final Initial Final\n"
797  "Coordinate Value Correction Correction "
798  " Value Accuracy Accuracy\n"
799  " (dd/dd/km) (dd/dd/km) (Meters) "
800  " (dd/dd/km) (Meters) (Meters)\n";
801  output += labels;
802 
803  SurfacePoint::CoordIndex idx = SurfacePoint::One;
804  output += QString(" LATITUDE%1%2%3%4%5%6\n")
805  .arg(formatValue(latInit, 17, 8)) // deg
806  .arg(formatValue(cor_lat_dd, 21, 8)) // deg
807  .arg(formatValue(cor_lat_m, 20, 8)) // m
808  .arg(formatValue(lat, 20, 8)) // deg
809  .arg(formatCoordAprioriSigmaString(idx, 18,8, true)) // m
810  .arg(formatCoordAdjustedSigmaString(idx, 18,8,errorPropagation)); // m
811 
812  idx = SurfacePoint::Two;
813  output += QString(" LONGITUDE%1%2%3%4%5%6\n")
814  .arg(formatValue(lonInit, 17, 8)) // deg
815  .arg(formatValue(cor_lon_dd, 21, 8)) // deg
816  .arg(formatValue(cor_lon_m, 20, 8)) // m
817  .arg(formatValue(lon, 20, 8)) // deg
818  .arg(formatCoordAprioriSigmaString(idx, 18,8, true)) // m
819  .arg(formatCoordAdjustedSigmaString(idx, 18,8,errorPropagation)); // m
820 
821  idx = SurfacePoint::Three;
822  output += QString(" RADIUS%1%2%3%4%5%6\n\n")
823  .arg(formatValue(radInit, 17, 8)) // km
824  .arg(formatValue(cor_rad_km, 21, 8)) // km
825  .arg(formatValue(cor_rad_m, 20, 8)) // m
826  .arg(formatValue(rad, 20, 8)) // km
827  .arg(formatCoordAprioriSigmaString(idx, 18,8, solveRadius)) // m
828  .arg(formatCoordAdjustedSigmaString(idx, 18,8, errorPropagation)); // m
829 
830  return output;
831  }
832 
833 
846  (bool errorPropagation) const {
847  int numRays = numberOfMeasures();
848  int numGoodRays = numRays - numberOfRejectedMeasures();
849  double X = m_controlPoint->GetAdjustedSurfacePoint().GetX().kilometers();
850  double Y = m_controlPoint->GetAdjustedSurfacePoint().GetY().kilometers();
851  double Z = m_controlPoint->GetAdjustedSurfacePoint().GetZ().kilometers();
852 
853  // Coordinate corrections are currently set in applyParameterCorrections. Units depend on
854  // coordTypeBundle coordinates (radians for Latitudinal coordinates and km for Rectangular).
855  // point corrections and initial sigmas.....
856  double cor_X_m = 0.; // X correction, meters
857  double cor_Y_m = 0.; // Y correction, meters
858  double cor_Z_m = 0.; // Z correction, meters
859  double XInit = Isis::Null;
860  double YInit = Isis::Null;
861  double ZInit = Isis::Null;
862 
863  if (m_coordTypeBundle == SurfacePoint::Latitudinal) {
864  double lat = m_controlPoint->GetAdjustedSurfacePoint().GetLatitude().degrees();
865  double latcor = m_corrections(0) * RAD2DEG;
866  double lon = m_controlPoint->GetAdjustedSurfacePoint().GetLongitude().degrees();
867  double loncor = m_corrections(1) * RAD2DEG;
868  double rad = m_controlPoint->GetAdjustedSurfacePoint().GetLocalRadius().kilometers();
869  double radcor = m_corrections(2);
870 
871  if (!IsSpecial(lat) && !IsSpecial(lon) && !IsSpecial(rad)) {
872  SurfacePoint latPoint(Latitude(lat-latcor, Angle::Degrees),
873  Longitude(lon-loncor, Angle::Degrees),
874  Distance(rad-radcor, Distance::Kilometers));
875  XInit = latPoint.GetX().kilometers();
876  YInit = latPoint.GetY().kilometers();
877  ZInit = latPoint.GetZ().kilometers();
878  cor_X_m = Isis::Null;
879  if (!IsSpecial(X)) {
880  cor_X_m = (X - XInit); // kilometers
881  }
882  cor_Y_m = Isis::Null;
883  if (!IsSpecial(Y)) {
884  cor_Y_m = (Y - YInit); // kilometers
885  }
886  cor_Z_m = Isis::Null;
887  if (!IsSpecial(Z)) {
888  cor_Z_m = (Z - ZInit); // kilometers
889  }
890  }
891  }
892  else if (m_coordTypeBundle == SurfacePoint::Rectangular) {
893  cor_X_m = m_corrections(0); // X correction, kilometers
894  cor_Y_m = m_corrections(1); // Y correction, kilometers
895  cor_Z_m = m_corrections(2); // Z correction, kilometers
896  if (!IsSpecial(X)) {
897  XInit = X - m_corrections(0); // km
898  }
899  if (!IsSpecial(Y)) {
900  YInit = Y - m_corrections(1); // km
901  }
902  if (!IsSpecial(Z)) {
903  ZInit = Z - m_corrections(2); // km
904  }
905  }
906  else {
907  IString msg ="Unknown surface point coordinate type enum [" +
908  toString(m_coordTypeBundle) + "]." ;
909  throw IException(IException::Programmer, msg, _FILEINFO_);
910  }
911 
912  QString pointType = ControlPoint::PointTypeToString(type()).toUpper();
913 
914  QString output;
915 
916  output = QString(" Label: %1\nStatus: %2\n Rays: %3 of %4\n")
917  .arg(id())
918  .arg(pointType)
919  .arg(numGoodRays)
920  .arg(numRays);
921 
922  QString labels =
923  "\n Point Initial Total Final "
924  " Initial Final\n"
925  " Coordinate Value Correction Value "
926  " Accuracy Accuracy\n"
927  " (km/km/km) (km) (km/km/km) "
928  " (Meters) (Meters)\n";
929  output += labels;
930 
931  SurfacePoint::CoordIndex idx = SurfacePoint::One;
932  output += QString(" BODY-FIXED-X%1%2%3%4%5\n")
933  .arg(formatValue(XInit, 17, 8)) // km
934  .arg(formatValue(cor_X_m, 20, 8)) // km
935  .arg(formatValue(X, 20, 8)) // km
936  .arg(formatCoordAprioriSigmaString(idx, 18,8, true)) // m
937  .arg(formatCoordAdjustedSigmaString(idx, 18,8,errorPropagation)); // m
938 
939  idx = SurfacePoint::Two;
940  output += QString(" BODY-FIXED-Y%1%2%3%4%5\n")
941  .arg(formatValue(YInit, 17, 8)) // km
942  .arg(formatValue(cor_Y_m, 20, 8)) // km
943  .arg(formatValue(Y, 20, 8)) // km
944  .arg(formatCoordAprioriSigmaString(idx, 18, 8, true)) // m
945  .arg(formatCoordAdjustedSigmaString(idx, 18,8,errorPropagation)); // m
946 
947  idx = SurfacePoint::Three;
948  output += QString(" BODY-FIXED-Z%1%2%3%4%5\n\n")
949  .arg(formatValue(ZInit, 17, 8)) // km
950  .arg(formatValue(cor_Z_m, 20, 8)) // km
951  .arg(formatValue(Z, 20, 8)) // km
952  // Set solveRadius to true to avoid limiting output information for Z.
953  // Set radius does not really apply to rectangular coordinates.
954  .arg(formatCoordAprioriSigmaString(idx, 18,8, true)) // m
955  .arg(formatCoordAdjustedSigmaString(idx, 18,8,errorPropagation)); // m
956 
957  return output;
958  }
959 
960 
971  QString BundleControlPoint::formatValue(double value, int fieldWidth, int precision) const {
972  QString output;
973  IsSpecial(value) ?
974  output = QString("%1").arg("Null", fieldWidth) :
975  output = QString("%1").arg(value, fieldWidth, 'f', precision);
976  return output;
977  }
978 
979 
993  QString BundleControlPoint::formatAprioriSigmaString(SurfacePoint::CoordIndex index,
994  int fieldWidth, int precision,
995  bool solveRadius) const {
996  QString aprioriSigmaStr;
997  QString pointType = ControlPoint::PointTypeToString(type()).toUpper();
998  if (pointType == "CONSTRAINED" || !solveRadius) {
999  pointType = "N/A";
1000  }
1001  double sigma = m_aprioriSigmas[int(index)];
1002  if (IsSpecial(sigma)) { // if globalAprioriSigma <= 0 (including Isis::Null), then m_aprioriSigmas = Null
1003  aprioriSigmaStr = QString("%1").arg(pointType, fieldWidth);
1004  }
1005  else {
1006  aprioriSigmaStr = QString("%1").arg(sigma, fieldWidth, 'f', precision);
1007  }
1008  return aprioriSigmaStr;
1009  }
1010 
1011 
1022 QString BundleControlPoint::formatCoordAprioriSigmaString(SurfacePoint::CoordIndex index,
1023  int fieldWidth,
1024  int precision,
1025  bool solveRadius) const {
1026  return formatAprioriSigmaString(index, fieldWidth, precision, solveRadius);
1027  }
1028 
1029 
1043  QString BundleControlPoint::formatAdjustedSigmaString(SurfacePoint::CoordIndex index,
1044  int fieldWidth, int precision,
1045  bool errorPropagation) const {
1046  QString adjustedSigmaStr;
1047 
1048  if (!errorPropagation) {
1049  adjustedSigmaStr = QString("%1").arg("N/A",fieldWidth);
1050  }
1051  else {
1052  double sigma = Isis::Null;
1054 
1055  if (IsSpecial(sigma)) {
1056  adjustedSigmaStr = QString("%1").arg("N/A",fieldWidth);
1057  }
1058  else {
1059  adjustedSigmaStr = QString("%1").arg(sigma, fieldWidth, 'f', precision);
1060  }
1061  }
1062 
1063  return adjustedSigmaStr;
1064  }
1065 
1066 
1079 QString BundleControlPoint::formatCoordAdjustedSigmaString(SurfacePoint::CoordIndex index,
1080  int fieldWidth, int precision,
1081  bool errorPropagation) const {
1082  return formatAdjustedSigmaString(index, fieldWidth, precision, errorPropagation);
1083  }
1084 
1085 
1093  (const BundleTargetBodyQsp target) {
1094  SurfacePoint surfacepoint = adjustedSurfacePoint();
1095  double pointLat = surfacepoint.GetLatitude().degrees();
1096  double pointLon = surfacepoint.GetLongitude().degrees();
1097  double pointRad = surfacepoint.GetLocalRadius().meters();
1098 
1099  // get point parameter corrections
1100  double latCorrection = m_nicVector[0];
1101  double lonCorrection = m_nicVector[1];
1102  double radCorrection = m_nicVector[2];
1103 
1104  // convert to degrees and apply
1105  pointLat += RAD2DEG * latCorrection;
1106  pointLon += RAD2DEG * lonCorrection;
1107 
1108  // Make sure updated values are still in valid range(0 to 360 for lon and -90 to 90 for lat)
1109  if (pointLat < -90.0) {
1110  pointLat = -180.0 - pointLat;
1111  pointLon = pointLon + 180.0;
1112  }
1113  if (pointLat > 90.0) {
1114  pointLat = 180.0 - pointLat;
1115  pointLon = pointLon + 180.0;
1116  }
1117  while (pointLon > 360.0) {
1118  pointLon = pointLon - 360.0;
1119  }
1120  while (pointLon < 0.0) {
1121  pointLon = pointLon + 360.0;
1122  }
1123 
1124  // convert to meters and apply
1125  pointRad += 1000. * radCorrection;
1126 
1127  // ken testing - if solving for target body mean radius, set radius to current
1128  // mean radius value
1129  // Only allow radius options for Latitudinal coordinates
1130  if (target && (target->solveMeanRadius() || target->solveTriaxialRadii()) ) {
1131  if (target->solveMeanRadius()) {
1132  surfacepoint.SetSphericalCoordinates(Latitude(pointLat, Angle::Degrees),
1133  Longitude(pointLon, Angle::Degrees),
1134  target->meanRadius());
1135  }
1136  else if (target->solveTriaxialRadii()) {
1137  Distance localRadius = target->localRadius(Latitude(pointLat, Angle::Degrees),
1138  Longitude(pointLon, Angle::Degrees));
1139  surfacepoint.SetSphericalCoordinates(Latitude(pointLat, Angle::Degrees),
1140  Longitude(pointLon, Angle::Degrees),
1141  localRadius);
1142  }
1143  }
1144  else {
1145  surfacepoint.SetSphericalCoordinates(Latitude(pointLat, Angle::Degrees),
1146  Longitude(pointLon, Angle::Degrees),
1147  Distance(pointRad, Distance::Meters));
1148  }
1149  // Reset the point now that it has been updated
1150  setAdjustedSurfacePoint(surfacepoint);
1151  }
1152 
1153 
1160  SurfacePoint surfacepoint = adjustedSurfacePoint();
1161  double pointX = surfacepoint.GetX().meters();
1162  double pointY = surfacepoint.GetY().meters();
1163  double pointZ = surfacepoint.GetZ().meters();
1164 
1165  // get point parameter corrections
1166  double XCorrection = m_nicVector[0];
1167  double YCorrection = m_nicVector[1];
1168  double ZCorrection = m_nicVector[2];
1169 
1170  // Convert corrections to meters and apply
1171  pointX += 1000. * XCorrection;
1172  pointY += 1000. * YCorrection;
1173  pointZ += 1000. * ZCorrection;
1174 
1175  surfacepoint.SetRectangularCoordinates(
1179  // Reset the point now that it has been updated
1180  setAdjustedSurfacePoint(surfacepoint);
1181  }
1182 }
Isis::BundleControlPoint::setAdjustedSurfacePoint
void setAdjustedSurfacePoint(SurfacePoint surfacePoint)
Sets the adjusted surface point for this BundleControlPoint.
Definition: BundleControlPoint.cpp:150
Isis::BundleControlPoint::~BundleControlPoint
~BundleControlPoint()
Destructor for BundleControlPoint.
Definition: BundleControlPoint.cpp:85
Isis::Distance::kilometers
double kilometers() const
Get the distance in kilometers.
Definition: Distance.cpp:106
Isis::BundleControlPoint::applyParameterCorrections
void applyParameterCorrections(LinearAlgebra::Vector imageSolution, SparseBlockMatrix &sparseNormals, const BundleTargetBodyQsp target)
Apply the parameter corrections to the bundle control point.
Definition: BundleControlPoint.cpp:390
Isis::Angle::Degrees
@ Degrees
Degrees are generally considered more human readable, 0-360 is one circle, however most math does not...
Definition: Angle.h:56
Isis::BundleControlPoint::formatAdjustedSigmaString
QString formatAdjustedSigmaString(SurfacePoint::CoordIndex, int fieldWidth, int precision, bool errorPropagation) const
Formats the adjusted sigma value indicated by the given type code.
Definition: BundleControlPoint.cpp:1043
Isis::SurfacePoint::MetersToLatitude
double MetersToLatitude(double latLength)
This method returns an angular measure of a distance in the direction of and relative to the latitude...
Definition: SurfacePoint.cpp:1248
Isis::BundleControlPoint::numberOfRejectedMeasures
int numberOfRejectedMeasures() const
Accesses the number of rejected measures for this BundleControlPoint.
Definition: BundleControlPoint.cpp:457
Isis::BundleControlPoint::updateAdjustedSurfacePointRectangularly
void updateAdjustedSurfacePointRectangularly()
pointer to the control point object this represents
Definition: BundleControlPoint.cpp:1159
Isis::ControlPoint::PointTypeToString
static QString PointTypeToString(PointType type)
Obtain a string representation of a given PointType.
Definition: ControlPoint.cpp:1333
Isis::Displacement::Meters
@ Meters
The distance is being specified in meters.
Definition: Displacement.h:40
Isis::BundleControlPoint::adjustedSurfacePoint
SurfacePoint adjustedSurfacePoint() const
Accesses the adjusted SurfacePoint associated with this BundleControlPoint.
Definition: BundleControlPoint.cpp:479
Isis::ControlPoint::ZeroNumberOfRejectedMeasures
void ZeroNumberOfRejectedMeasures()
Initialize the number of rejected measures to 0.
Definition: ControlPoint.cpp:2065
Isis::ControlPoint::IsCoord3Constrained
bool IsCoord3Constrained()
Return bool indicating if 3rd coordinate is Constrained or not.
Definition: ControlPoint.cpp:1651
Isis::BundleControlPoint::type
ControlPoint::PointType type() const
Accesses BundleControlPoint's type.
Definition: BundleControlPoint.cpp:502
Isis::ControlPoint::GetMeasure
const ControlMeasure * GetMeasure(QString serialNumber) const
Get a control measure based on its cube's serial number.
Definition: ControlPoint.cpp:416
Isis::Latitude
This class is designed to encapsulate the concept of a Latitude.
Definition: Latitude.h:51
Isis::BundleControlPoint::m_coordTypeReports
SurfacePoint::CoordinateType m_coordTypeReports
BundleControlPoint coordinate type.
Definition: BundleControlPoint.h:182
Isis::BundleControlPoint::weights
boost::numeric::ublas::bounded_vector< double, 3 > & weights()
Accesses the 3 dimensional ordered vector of weight values associated with coordinate1,...
Definition: BundleControlPoint.cpp:577
Isis::BundleControlPoint::m_cholmodQMatrix
SparseBlockRowMatrix m_cholmodQMatrix
The CholMod matrix associated with this point.
Definition: BundleControlPoint.h:180
Isis::DEG2RAD
const double DEG2RAD
Multiplier for converting from degrees to radians.
Definition: Constants.h:43
Isis::ControlPoint::GetResidualRms
double GetResidualRms() const
Get rms of residuals.
Definition: ControlPoint.cpp:2161
Isis::BundleControlPoint::setRejected
void setRejected(bool reject)
Sets this BundleControlPoint to rejected or not rejected.
Definition: BundleControlPoint.cpp:174
Isis::ControlPoint::SetAdjustedSurfacePoint
Status SetAdjustedSurfacePoint(SurfacePoint newSurfacePoint)
Set or update the surface point relating to this control point.
Definition: ControlPoint.cpp:692
Isis::SurfacePoint::SetRectangularCoordinates
void SetRectangularCoordinates(const Displacement &x, const Displacement &y, const Displacement &z)
Set surface point in rectangular coordinates.
Definition: SurfacePoint.cpp:327
Isis::ControlPoint::SetNumberOfRejectedMeasures
void SetNumberOfRejectedMeasures(int numRejected)
Set (update) the number of rejected measures for the control point.
Definition: ControlPoint.cpp:2078
Isis::BundleControlPoint::setWeights
void setWeights(const BundleSettingsQsp settings)
Sets the weights using the given BundleSettings QSharedPointer and a conversion value for meters to r...
Definition: BundleControlPoint.cpp:186
Isis::BundleControlPoint::m_weights
boost::numeric::ublas::bounded_vector< double, 3 > m_weights
weights for point parameters
Definition: BundleControlPoint.h:176
Isis::BundleControlPoint::isRejected
bool isRejected() const
Method used to determine whether this control point is rejected.
Definition: BundleControlPoint.cpp:435
Isis::ControlPoint::PointType
PointType
These are the valid 'types' of point.
Definition: ControlPoint.h:364
Isis::SurfacePoint::GetLatitude
Latitude GetLatitude() const
Return the body-fixed latitude for the surface point.
Definition: SurfacePoint.cpp:1665
Isis::BundleControlPoint
This class holds information about a control point that BundleAdjust needs to run correctly.
Definition: BundleControlPoint.h:91
Isis::ControlPoint::GetNumberOfRejectedMeasures
int GetNumberOfRejectedMeasures() const
Get the number of rejected measures on the control point.
Definition: ControlPoint.cpp:2089
Isis::ControlPoint::GetId
QString GetId() const
Return the Id of the control point.
Definition: ControlPoint.cpp:1306
Isis::BundleControlPoint::corrections
boost::numeric::ublas::bounded_vector< double, 3 > & corrections()
Accesses the 3 dimensional ordered vector of correction values associated with coord1,...
Definition: BundleControlPoint.cpp:541
QSharedPointer< BundleSettings >
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::BundleControlPoint::formatBundleLatitudinalOutputDetailString
QString formatBundleLatitudinalOutputDetailString(bool errorPropagation, bool solveRadius=false) const
Formats a detailed output string table for this Latitudinal BundleControlPoint.
Definition: BundleControlPoint.cpp:702
Isis::IsSpecial
bool IsSpecial(const double d)
Returns if the input pixel is special.
Definition: SpecialPixel.h:197
Isis::BundleControlPoint::residualRms
double residualRms() const
Gets the root-mean-square (rms) of the BundleControlPoint's residuals.
Definition: BundleControlPoint.cpp:469
Isis::BundleControlPoint::aprioriSigmas
boost::numeric::ublas::bounded_vector< double, 3 > & aprioriSigmas()
Accesses the 3 dimensional ordered vector of apriori sigmas (apriori coordinate1, apriori coordinate2...
Definition: BundleControlPoint.cpp:552
Isis::SurfacePoint::GetCoord
double GetCoord(CoordinateType type, CoordIndex index, CoordUnits units)
This method returns a coordinate of a SurfacePoint.
Definition: SurfacePoint.cpp:962
Isis::Distance
Distance measurement, usually in meters.
Definition: Distance.h:34
Isis::Longitude
This class is designed to encapsulate the concept of a Longitude.
Definition: Longitude.h:40
Isis::BundleControlPoint::formatBundleRectangularOutputDetailString
QString formatBundleRectangularOutputDetailString(bool errorPropagation) const
Formats a detailed output string table for this Rectangular BundleControlPoint.
Definition: BundleControlPoint.cpp:846
Isis::BundleMeasureQsp
QSharedPointer< BundleMeasure > BundleMeasureQsp
Definition for BundleMeasureQsp, a shared pointer to a BundleMeasure.
Definition: BundleMeasure.h:99
Isis::SurfacePoint::LongitudeToMeters
double LongitudeToMeters(double longitude) const
This method returns a length in meters version of a delta longitude angle in radians relative to the ...
Definition: SurfacePoint.cpp:1336
Isis::BundleControlPoint::setSigmaWeightFromGlobals
void setSigmaWeightFromGlobals(double gSigma, int index)
Sets the member sigmas and weights from a global sigma.
Definition: BundleControlPoint.cpp:283
Isis::Distance::Kilometers
@ Kilometers
The distance is being specified in kilometers.
Definition: Distance.h:45
Isis::Displacement
Displacement is a signed length, usually in meters.
Definition: Displacement.h:31
Isis::BundleControlPoint::BundleControlPoint
BundleControlPoint(BundleSettingsQsp bundleSettings, ControlPoint *point)
Constructs a BundleControlPoint object from a ControlPoint.
Definition: BundleControlPoint.cpp:33
Isis::ControlPoint
A single control point.
Definition: ControlPoint.h:354
Isis::BundleControlPoint::coordTypeReports
SurfacePoint::CoordinateType coordTypeReports() const
Accesses BundleControlPoint's coordinate type for reports.
Definition: BundleControlPoint.cpp:515
Isis::ControlPoint::Fixed
@ Fixed
A Fixed point is a Control Point whose lat/lon is well established and should not be changed.
Definition: ControlPoint.h:371
Isis::Distance::Meters
@ Meters
The distance is being specified in meters.
Definition: Distance.h:43
Isis::ControlPoint::IsCoord2Constrained
bool IsCoord2Constrained()
Return bool indicating if 2nd coordinate is Constrained or not.
Definition: ControlPoint.cpp:1641
Isis::Displacement::meters
double meters() const
Get the displacement in meters.
Definition: Displacement.cpp:73
Isis::SurfacePoint::Rectangular
@ Rectangular
Body-fixed rectangular x/y/z coordinates.
Definition: SurfacePoint.h:141
Isis::BundleControlPoint::nicVector
boost::numeric::ublas::bounded_vector< double, 3 > & nicVector()
Accesses the 3 dimensional ordered NIC vector.
Definition: BundleControlPoint.cpp:587
Isis::LinearAlgebra::Vector
boost::numeric::ublas::vector< double > Vector
Definition for an Isis::LinearAlgebra::Vector of doubles.
Definition: LinearAlgebra.h:120
Isis::SurfacePoint::CoordinateType
CoordinateType
Defines the coordinate typ, units, and coordinate index for some of the output methods.
Definition: SurfacePoint.h:139
Isis::BundleControlPoint::formatCoordAprioriSigmaString
QString formatCoordAprioriSigmaString(SurfacePoint::CoordIndex index, int fieldWidth, int precision, bool solveRadius=false) const
Formats the apriori coordinate 1 (latitude or X) sigma value.
Definition: BundleControlPoint.cpp:1022
Isis::BundleControlPoint::computeResiduals
void computeResiduals()
Computes the residuals for this BundleControlPoint.
Definition: BundleControlPoint.cpp:140
Isis::BundleControlPoint::m_adjustedSigmas
boost::numeric::ublas::bounded_vector< double, 3 > m_adjustedSigmas
adjusted sigmas for point parameters
Definition: BundleControlPoint.h:174
Isis::SurfacePoint::LatitudeToMeters
double LatitudeToMeters(double latitude) const
This method returns a Displacement of an Angle relative to the current SurfacePoint latitude.
Definition: SurfacePoint.cpp:1310
Isis::ControlPoint::Constrained
@ Constrained
A Constrained point is a Control Point whose lat/lon/radius is somewhat established and should not be...
Definition: ControlPoint.h:376
Isis::BundleControlPoint::m_aprioriSigmas
boost::numeric::ublas::bounded_vector< double, 3 > m_aprioriSigmas
apriori sigmas for point parameters
Definition: BundleControlPoint.h:172
Isis::SurfacePoint::MetersToLongitude
double MetersToLongitude(double lonLength)
This method returns an angular measure in radians of a distance in the direction of and relative to t...
Definition: SurfacePoint.cpp:1276
Isis::BundleControlPoint::formatValue
QString formatValue(double value, int fieldWidth, int precision) const
Formats the given double precision value using the specified field width and precision.
Definition: BundleControlPoint.cpp:971
Isis::BundleControlPoint::zeroNumberOfRejectedMeasures
void zeroNumberOfRejectedMeasures()
Resets the number of rejected measures for this BundleControlPoint to zero.
Definition: BundleControlPoint.cpp:348
Isis::BundleControlPoint::formatBundleOutputSummaryString
QString formatBundleOutputSummaryString(bool errorPropagation) const
Formats an output summary string for this BundleControlPoint.
Definition: BundleControlPoint.cpp:616
Isis::BundleMeasure
A container class for a ControlMeasure.
Definition: BundleMeasure.h:55
Isis::IException
Isis exception class.
Definition: IException.h:91
Isis::BundleControlPoint::rawControlPoint
ControlPoint * rawControlPoint() const
Accessor for the raw ControlPoint object used for this BundleControlPoint.
Definition: BundleControlPoint.cpp:425
Isis::SurfacePoint::GetLongitude
Longitude GetLongitude() const
Return the body-fixed longitude for the surface point.
Definition: SurfacePoint.cpp:1685
Isis::Null
const double Null
Value for an Isis Null pixel.
Definition: SpecialPixel.h:95
Isis::Displacement::kilometers
double kilometers() const
Get the displacement in kilometers.
Definition: Displacement.cpp:94
Isis::BundleControlPoint::addMeasure
BundleMeasureQsp addMeasure(ControlMeasure *controlMeasure)
Creates a BundleMeasure from the given ControlMeasure and appends it to this BundleControlPoint's mea...
Definition: BundleControlPoint.cpp:125
Isis::BundleControlPoint::setNumberOfRejectedMeasures
void setNumberOfRejectedMeasures(int numRejected)
Sets the number of rejected measures for this BundleControlPoint.
Definition: BundleControlPoint.cpp:162
Isis::BundleControlPoint::numberOfMeasures
int numberOfMeasures() const
Accesses number of measures associated with this BundleControlPoint.
Definition: BundleControlPoint.cpp:445
Isis::BundleControlPoint::updateAdjustedSurfacePointLatitudinally
void updateAdjustedSurfacePointLatitudinally(const BundleTargetBodyQsp target)
Apply the parameter corrections to the bundle control point latitudinally.
Definition: BundleControlPoint.cpp:1093
Isis::SurfacePoint::GetLocalRadius
Distance GetLocalRadius() const
Return the radius of the surface point.
Definition: SurfacePoint.cpp:1732
Isis::IException::Programmer
@ Programmer
This error is for when a programmer made an API call that was illegal.
Definition: IException.h:146
Isis::ControlPoint::IsCoord1Constrained
bool IsCoord1Constrained()
Return bool indicating if 1st coordinate is Constrained or not.
Definition: ControlPoint.cpp:1630
Isis::SurfacePoint::SetSphericalCoordinates
void SetSphericalCoordinates(const Latitude &lat, const Longitude &lon, const Distance &radius)
Update spherical coordinates (lat/lon/radius)
Definition: SurfacePoint.cpp:553
Isis::BundleControlPoint::adjustedSigmas
boost::numeric::ublas::bounded_vector< double, 3 > & adjustedSigmas()
Accesses the 3 dimensional ordered vector of adjusted sigmas (adjusted coordinate1,...
Definition: BundleControlPoint.cpp:566
Isis::BundleControlPoint::formatCoordAdjustedSigmaString
QString formatCoordAdjustedSigmaString(SurfacePoint::CoordIndex, int fieldWidth, int precision, bool errorPropagation) const
Formats the adjusted coordinate sigma value.
Definition: BundleControlPoint.cpp:1079
Isis::BundleControlPoint::formatAprioriSigmaString
QString formatAprioriSigmaString(SurfacePoint::CoordIndex index, int fieldWidth, int precision, bool solveRadius=false) const
Formats the apriori sigma value indicated by the given type code.
Definition: BundleControlPoint.cpp:993
Isis::BundleControlPoint::productAlphaAV
void productAlphaAV(double alpha, SparseBlockMatrix &sparseNormals, LinearAlgebra::Vector &v1)
Perform the matrix multiplication v2 = alpha ( Q x v1 ).
Definition: BundleControlPoint.cpp:362
Isis::ControlPoint::Free
@ Free
A Free point is a Control Point that identifies common measurements between two or more cubes.
Definition: ControlPoint.h:384
Isis::Distance::meters
double meters() const
Get the distance in meters.
Definition: Distance.cpp:85
Isis::Angle::degrees
double degrees() const
Get the angle in units of Degrees.
Definition: Angle.h:232
Isis::ControlPoint::GetType
PointType GetType() const
Definition: ControlPoint.cpp:1401
Isis::ControlPoint::SetRejected
Status SetRejected(bool rejected)
Set the jigsawRejected state.
Definition: ControlPoint.cpp:539
Isis::SurfacePoint::GetSigmaDistance
Distance GetSigmaDistance(CoordinateType type, CoordIndex index)
This method returns a sigma of a SurfacePoint coordinate as a Distance.
Definition: SurfacePoint.cpp:1094
Isis::BundleControlPoint::m_corrections
boost::numeric::ublas::bounded_vector< double, 3 > m_corrections
corrections to point parameters
Definition: BundleControlPoint.h:170
Isis::BundleControlPoint::formatBundleOutputDetailString
QString formatBundleOutputDetailString(bool errorPropagation, bool solveRadius=false) const
Formats a detailed output string table for this BundleControlPoint.
Definition: BundleControlPoint.cpp:663
Isis::SurfacePoint::Latitudinal
@ Latitudinal
Planetocentric latitudinal (lat/lon/rad) coordinates.
Definition: SurfacePoint.h:140
Isis::BundleControlPoint::cholmodQMatrix
SparseBlockRowMatrix & cholmodQMatrix()
Accesses the CholMod matrix associated with this BundleControlPoint.
Definition: BundleControlPoint.cpp:597
Isis::IString
Adds specific functionality to C++ strings.
Definition: IString.h:165
Isis::SparseBlockRowMatrix
SparseBlockRowMatrix.
Definition: SparseBlockMatrix.h:125
Isis::BundleControlPoint::coordTypeBundle
SurfacePoint::CoordinateType coordTypeBundle() const
Accesses BundleControlPoint's control point coordinate type for the bundle adjustment.
Definition: BundleControlPoint.cpp:528
Isis::BundleControlPoint::m_nicVector
boost::numeric::ublas::bounded_vector< double, 3 > m_nicVector
array of NICs (see Brown, 1976)
Definition: BundleControlPoint.h:178
Isis::SurfacePoint::GetWeight
double GetWeight(CoordinateType type, CoordIndex index)
This method returns the weight of a SurfacePoint coordinate Note: At this time a units argument is no...
Definition: SurfacePoint.cpp:1491
Isis::SparseBlockMatrix
SparseBlockMatrix.
Definition: SparseBlockMatrix.h:186
Isis::RAD2DEG
const double RAD2DEG
Multiplier for converting from radians to degrees.
Definition: Constants.h:44
Isis::Displacement::Kilometers
@ Kilometers
The distance is being specified in kilometers.
Definition: Displacement.h:42
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::ControlPoint::ComputeResiduals
Status ComputeResiduals()
This method computes the BundleAdjust residuals for a point.
Definition: ControlPoint.cpp:1028
Isis::BundleControlPoint::copy
void copy(const BundleControlPoint &src)
Copies given BundleControlPoint to this BundleControlPoint.
Definition: BundleControlPoint.cpp:94
Isis::ControlMeasure
a control measurement
Definition: ControlMeasure.h:175