Isis 3 Programmer Reference
BundleControlPoint.cpp
1 #include "BundleControlPoint.h"
2 
3 #include <QDebug>
4 
5 // boost lib
6 #include <boost/numeric/ublas/vector_proxy.hpp>
7 
8 // Isis lib
9 #include "ControlMeasure.h"
10 #include "Latitude.h"
11 #include "LinearAlgebra.h"
12 #include "Longitude.h"
13 #include "SparseBlockMatrix.h"
14 #include "SpecialPixel.h"
15 
16 namespace Isis {
17 
18 
26  ControlPoint *controlPoint) {
27  m_controlPoint = controlPoint;
28 
29  // setup vector of BundleMeasures for this control point
30  int numMeasures = controlPoint->GetNumMeasures();
31  for (int i = 0; i < numMeasures; i++) {
32  ControlMeasure *controlMeasure = controlPoint->GetMeasure(i);
33  if (controlMeasure->IsIgnored()) {
34  continue;
35  }
36 
37  addMeasure(controlMeasure);
38  }
39 
40  // initialize to 0.0
41  m_corrections.clear();
42  m_weights.clear();
43  m_nicVector.clear();
44 
45  // initialize to Null for consistency with other bundle classes...
46  m_aprioriSigmas.clear();
50  m_adjustedSigmas.clear();
54 
55  // initialize coordinate type from settings
56  // m_coordType = SurfacePoint::Latitudinal;
57  m_coordTypeReports = bundleSettings->controlPointCoordTypeReports();
58  m_coordTypeBundle = bundleSettings->controlPointCoordTypeBundle();
59  setWeights(bundleSettings);
60  }
61 
62 
70  copy(src);
71  }
72 
73 
78  }
79 
80 
87 
88  // sanity check
89  clear();
90 
91  m_controlPoint = src.m_controlPoint;
92 
93  int numMeasures = src.size();
94 
95  for (int i = 0; i < numMeasures; i++)
96  append(BundleMeasureQsp( new BundleMeasure(*(src.at(i))) ));
97 
101  m_weights = src.m_weights;
102  m_nicVector = src.m_nicVector;
104  m_coordTypeBundle = src.m_coordTypeBundle;
105  // *** TODO *** Why is the cholmodQMatrix not copied? Ask Ken
106  }
107 
108 
118 
119  BundleMeasureQsp bundleMeasure = BundleMeasureQsp( new BundleMeasure(controlMeasure, this) );
120 
121  append(bundleMeasure);
122 
123  return bundleMeasure;
124  }
125 
126 
133  m_controlPoint->ComputeResiduals();
134  }
135 
136 
143  m_controlPoint->SetAdjustedSurfacePoint(surfacePoint);
144  }
145 
146 
155  m_controlPoint->SetNumberOfRejectedMeasures(numRejected);
156  }
157 
158 
167  m_controlPoint->SetRejected(reject);
168  }
169 
170 
179 
180  double globalPointCoord1AprioriSigma = settings->globalPointCoord1AprioriSigma();
181  double globalPointCoord2AprioriSigma = settings->globalPointCoord2AprioriSigma();
182  double globalPointCoord3AprioriSigma = settings->globalPointCoord3AprioriSigma();
183 
184  if (m_controlPoint->GetType() == ControlPoint::Fixed) {
185  m_weights[0] = 1.0e+50;
186  m_weights[1] = 1.0e+50;
187  m_weights[2] = 1.0e+50;
188  // m_aprioriSigmas = Isis::Null by default
189  }
190 
191  if (m_controlPoint->GetType() == ControlPoint::Free) {
192  // Global sigmas are temporary and should not be stored in the control net. Currently
193  // global sigmas are always in meters. Make sure unit of weights is 1/(var unit squared),
194  // where var is a value being solved for in the adjustment. Latitude and longitude are in
195  // radians and everything else is in km.
196  if (!IsSpecial(globalPointCoord1AprioriSigma)) {
197  setSigmaWeightFromGlobals(globalPointCoord1AprioriSigma, 0);
198  } // else m_aprioriSigma = Isis::Null
199  // m_weights = 0.0
200  if (!IsSpecial(globalPointCoord2AprioriSigma)) {
201  m_aprioriSigmas[1] = globalPointCoord2AprioriSigma;
202  setSigmaWeightFromGlobals(globalPointCoord2AprioriSigma, 1);
203  }// else m_aprioriSigma = Isis::Null
204  // m_weights = 0.0
205  if (m_coordTypeBundle == SurfacePoint::Latitudinal && !settings->solveRadius()) {
206  m_weights[2] = 1.0e+50;
207  }
208  else {
209  if (!IsSpecial(globalPointCoord3AprioriSigma)) {
210  setSigmaWeightFromGlobals(globalPointCoord3AprioriSigma, 2);
211  }
212  }
213  }
214  if (m_controlPoint->GetType() == ControlPoint::Constrained) {
215  // *** Be careful...Is m_aprioriSigmas an output (for reports) or bundle variable? ***
216 
217  // Assuming the member variable sigmas are for output reports (internal use only) so use
218  // the report coordinate type to calculate. All point sigmas are in meters. Radius weights
219  // are in 1/km^2. Make x/y/z weights the same because BundleAdjust works in km.
220  // Weights and corrections go into the bundle so use bundle coordinate type.
221  if ( m_controlPoint->IsCoord1Constrained() ) {
222  m_aprioriSigmas[0] =
223  m_controlPoint->GetAprioriSurfacePoint().GetSigmaDistance
224  (m_coordTypeReports, SurfacePoint::One).meters();
225  m_weights[0] =
226  m_controlPoint->GetAprioriSurfacePoint().GetWeight(m_coordTypeBundle, SurfacePoint::One);
227  }
228  else if (!IsSpecial(globalPointCoord1AprioriSigma)) {
229  setSigmaWeightFromGlobals(globalPointCoord1AprioriSigma, 0);
230  } // else not constrained and global sigma is Null, then m_aprioriSigmas = Isis::Null
231  // m_weights = 0.0
232 
233  if ( m_controlPoint->IsCoord2Constrained() ) {
234  m_aprioriSigmas[1] =
235  m_controlPoint->GetAprioriSurfacePoint().GetSigmaDistance
236  (m_coordTypeReports, SurfacePoint::Two).meters();
237  m_weights[1] =
238  m_controlPoint->GetAprioriSurfacePoint().GetWeight(m_coordTypeBundle, SurfacePoint::Two);
239  }
240  else if (!IsSpecial(globalPointCoord2AprioriSigma)) {
241  setSigmaWeightFromGlobals(globalPointCoord2AprioriSigma, 1);
242  } // else not constrained and global sigma is Null, then m_aprioriSigmas = Isis::Null
243  // m_weights = 0.0
244 
245  if (m_coordTypeBundle == SurfacePoint::Latitudinal && !settings->solveRadius()) {
246  m_weights[2] = 1.0e+50;
247  }
248  else {
249  if ( m_controlPoint->IsCoord3Constrained() ) {
250  m_aprioriSigmas[2] =
251  m_controlPoint->GetAprioriSurfacePoint().GetSigmaDistance
252  (m_coordTypeReports, SurfacePoint::Three).meters();
253  m_weights[2] = m_controlPoint->GetAprioriSurfacePoint().GetWeight
254  (m_coordTypeBundle, SurfacePoint::Three);
255  }
256  else if (!IsSpecial(globalPointCoord3AprioriSigma)) {
257  setSigmaWeightFromGlobals(globalPointCoord3AprioriSigma, 2);
258  } // else not constrained and global sigma is Null, then m_aprioriSigmas = Isis::Null
259  // m_weights = 0.0
260  }
261  }
262  }
263 
264 
275  void BundleControlPoint::setSigmaWeightFromGlobals(double gSigma, int index) {
276  m_aprioriSigmas[index] = gSigma;
277 
278  switch (index) {
279  case 0:
280  if (m_coordTypeBundle == SurfacePoint::Latitudinal) {
281  m_aprioriSigmas[index] = gSigma;
282  double sigmaRadians =
283  m_controlPoint->GetAprioriSurfacePoint().MetersToLatitude(gSigma); // m to radians
284  if (sigmaRadians > DBL_EPSILON) {
285  m_weights[index] = 1. / (sigmaRadians*sigmaRadians); // 1/radians^2
286  }
287  }
288  else if (m_coordTypeBundle == SurfacePoint::Rectangular) {
289  double sigmakm = gSigma * .001; // km
290  if (sigmakm > DBL_EPSILON) {
291  m_weights[0] = 1./ (sigmakm*sigmakm); // 1/km^2
292  }
293  }
294  else {
295  IString msg ="Unknown surface point coordinate type enum ["
296  + toString(m_coordTypeBundle) + "]." ;
298  }
299  break;
300  case 1:
301  // if (!IsSpecial(globalPointCoord2AprioriSigma)) {
302  if (m_coordTypeBundle == SurfacePoint::Latitudinal) {
303  double sigmaRadians =
304  m_controlPoint->GetAprioriSurfacePoint().MetersToLongitude(gSigma); // m to radians
305  if (sigmaRadians > DBL_EPSILON) {
306  m_weights[1] = 1. / (sigmaRadians*sigmaRadians); // 1/radians^2
307  }
308  }
309  else if (m_coordTypeBundle == SurfacePoint::Rectangular) {
310  double sigmakm = gSigma * 0.001; // km
311  if (sigmakm > DBL_EPSILON) m_weights[1] = 1./ (sigmakm*sigmakm); // 1/km^2
312  }
313  else {
314  IString msg ="Unknown surface point coordinate type enum ["
315  + toString(m_coordTypeBundle) + "]." ;
317  }
318  break;
319 
320  case 2:
321  // Coordinate 2 is either latitudinal radius or rectangular z, both in m
322  {
323  double sigmakm = gSigma * .001; // km
324  m_weights[2] = 1./ (sigmakm*sigmakm); // 1/km^2
325  }
326  break;
327 
328  default:
329  IString msg ="Unknown coordinate index [" + toString(index) + "]." ;
331  }
332  }
333 
334 
341  m_controlPoint->ZeroNumberOfRejectedMeasures();
342  }
343 
344 
355  SparseBlockMatrix &sparseNormals,
356  LinearAlgebra::Vector &v1) {
357 
358  QMapIterator< int, LinearAlgebra::Matrix * > Qit(m_cholmodQMatrix);
359 
360  int subrangeStart, subrangeEnd;
361 
362  while ( Qit.hasNext() ) {
363  Qit.next();
364 
365  int columnIndex = Qit.key();
366 
367  subrangeStart = sparseNormals.at(columnIndex)->startColumn();
368  subrangeEnd = subrangeStart + Qit.value()->size2();
369 
370  m_nicVector += alpha * prod(*(Qit.value()),subrange(v1,subrangeStart,subrangeEnd));
371  }
372  }
373 
374 
383  SparseBlockMatrix &sparseNormals,
384  const BundleTargetBodyQsp target) {
385  if (!isRejected()) {
386 
387  // subtract product of Q and nj from NIC
388  productAlphaAV(-1.0, sparseNormals, imageSolution);
389 
390  // update adjusted surface point appropriately for the coordinate type
391  switch (m_coordTypeBundle) {
394  break;
397  break;
398  default:
399  IString msg ="Unknown surface point coordinate type enum ["
400  + toString(m_coordTypeBundle) + "]." ;
402  }
403 
404  // sum and save corrections
405  m_corrections(0) += m_nicVector[0];
406  m_corrections(1) += m_nicVector[1];
407  m_corrections(2) += m_nicVector[2];
408  }
409  }
410 
411 
418  return m_controlPoint;
419  }
420 
421 
428  return m_controlPoint->IsRejected();
429  }
430 
431 
438  return this->size();
439  }
440 
441 
450  return m_controlPoint->GetNumberOfRejectedMeasures();
451  }
452 
453 
462  return m_controlPoint->GetResidualRms();
463  }
464 
465 
472  return m_controlPoint->GetAdjustedSurfacePoint();
473  }
474 
475 
481  QString BundleControlPoint::id() const {
482  return m_controlPoint->GetId();
483  }
484 
485 
495  return m_controlPoint->GetType();
496  }
497 
498 
508  return m_coordTypeReports;
509  }
510 
511 
521  return m_coordTypeBundle;
522  }
523 
524 
525 // *** TODO *** Need to add bounded vectors to Linear Algebra class
533  boost::numeric::ublas::bounded_vector< double, 3 > &BundleControlPoint::corrections() {
534  return m_corrections;
535  }
536 
537 
544  boost::numeric::ublas::bounded_vector< double, 3 > &BundleControlPoint::aprioriSigmas() {
545  return m_aprioriSigmas;
546  // Be careful about units. Latitude and longitude sigmas are in radians. Radius and X/Y/Z are
547  // meters.
548 
549  }
550 
551 
558  boost::numeric::ublas::bounded_vector< double, 3 > &BundleControlPoint::adjustedSigmas() {
559  return m_adjustedSigmas;
560  }
561 
562 
569  boost::numeric::ublas::bounded_vector< double, 3 > &BundleControlPoint::weights() {
570  return m_weights;
571  }
572 
573 
579  boost::numeric::ublas::bounded_vector<double, 3> &BundleControlPoint::nicVector() {
580  return m_nicVector;
581  }
582 
583 
590  return m_cholmodQMatrix;
591  }
592 
593 
608  QString BundleControlPoint::formatBundleOutputSummaryString(bool errorPropagation) const {
609 
610  int numRays = numberOfMeasures(); // should this depend on the raw point, as written, or this->size()???
611  int numGoodRays = numRays - numberOfRejectedMeasures();
612  double ResidualRms = residualRms();
613 
614  // Return generic set of control point coordinates as floating point values
615  // Be careful because of the potential to confuse units.
616  SurfacePoint::CoordUnits units = SurfacePoint::Degrees;
617  if (m_coordTypeReports == SurfacePoint::Rectangular) units = SurfacePoint::Kilometers;
618  double cpc1 = m_controlPoint->GetAdjustedSurfacePoint().GetCoord(m_coordTypeReports,
619  SurfacePoint::One, units);
620  double cpc2 = m_controlPoint->GetAdjustedSurfacePoint().GetCoord(m_coordTypeReports,
621  SurfacePoint::Two, units);
622  double cpc3 = m_controlPoint->GetAdjustedSurfacePoint().GetCoord(m_coordTypeReports,
623  SurfacePoint::Three, SurfacePoint::Kilometers);
624 
625  QString pointType = ControlPoint::PointTypeToString(type()).toUpper();
626 
627  QString output = QString("%1%2%3 of %4%5%6%7%8%9%10%11\n")
628  .arg(id(), 16)
629  .arg(pointType, 15)
630  .arg(numGoodRays, 5)
631  .arg(numRays)
632  .arg(formatValue(ResidualRms, 6, 2))
633  .arg(formatValue(cpc1, 16, 8)) // deg km
634  .arg(formatValue(cpc2, 16, 8)) // deg OR km
635  .arg(formatValue(cpc3, 16, 8)) // km km
636  .arg(formatCoordAdjustedSigmaString(SurfacePoint::One, 16, 8, errorPropagation)) // m
637  .arg(formatCoordAdjustedSigmaString(SurfacePoint::Two, 16, 8, errorPropagation)) // m
638  .arg(formatCoordAdjustedSigmaString(SurfacePoint::Three, 16, 8, errorPropagation)); // m
639 
640  return output;
641  }
642 
643 
656  bool solveRadius) const {
657  QString output;
658 
659  switch (m_coordTypeReports) {
661  output = formatBundleLatitudinalOutputDetailString(errorPropagation, solveRadius);
662  break;
664  output = formatBundleRectangularOutputDetailString(errorPropagation);
665  break;
666  default:
667  IString msg ="Unknown surface point coordinate type enum [" + toString(m_coordTypeBundle) + "]." ;
669  }
670  return output;
671  }
672 
673 
695  bool errorPropagation,
696  bool solveRadius) const {
697 
698  int numRays = numberOfMeasures();
699  int numGoodRays = numRays - numberOfRejectedMeasures();
700  double lat = m_controlPoint->GetAdjustedSurfacePoint().GetLatitude().degrees();
701  double lon = m_controlPoint->GetAdjustedSurfacePoint().GetLongitude().degrees();
702  double rad = m_controlPoint->GetAdjustedSurfacePoint().GetLocalRadius().kilometers();
703 
704  // Use the local radius in meters, rad*1000., to convert radians to meters now instead of the
705  // target body equatorial radius DAC 09/17/2018
706  double rtm = rad * 1000.;
707 
708  // Coordinate corrections are currently set in BundleAdjust::applyParameterCorrections in the
709  // coordTypeBundle coordinates (radians for Latitudinal coordinates and km for Rectangular).
710  // point corrections and initial sigmas
711  double cor_lat_dd = 0.; // lat correction, decimal degs
712  double cor_lon_dd = 0.; // lon correction, decimal degs
713  double cor_rad_km = 0.; // radius correction, kilometers
714  double cor_lat_m = 0.; // lat correction, meters
715  double cor_lon_m = 0.; // lon correction, meters
716  double cor_rad_m = 0.; // radius correction, meters
717  double latInit = Isis::Null;
718  double lonInit = Isis::Null;
719  double radInit = Isis::Null;
720 
721  if (m_coordTypeBundle == SurfacePoint::Rectangular) {
722  double x = m_controlPoint->GetAdjustedSurfacePoint().GetX().kilometers();
723  double xCor = m_corrections(0); // km
724  double y = m_controlPoint->GetAdjustedSurfacePoint().GetY().kilometers();
725  double yCor = m_corrections(1); // km
726  double z = m_controlPoint->GetAdjustedSurfacePoint().GetZ().kilometers();
727  double zCor = m_corrections(2); // km
728 
729  if (!IsSpecial(x) && !IsSpecial(y) && !IsSpecial(z)) {
733  latInit = rectPoint.GetLatitude().degrees();
734  lonInit = rectPoint.GetLongitude().degrees();
735  radInit = rectPoint.GetLocalRadius().kilometers();
736  if (!IsSpecial(lat)) {
737  cor_lat_dd = (lat - latInit); // degrees
738  cor_lat_m = cor_lat_dd * DEG2RAD * rtm;
739  }
740  if (!IsSpecial(lon)) {
741  cor_lon_dd = (lon - lonInit); // degrees
742  cor_lon_m = cor_lon_dd * DEG2RAD * rtm * cos(lat*DEG2RAD); // lon corrections meters
743  }
744  if (!IsSpecial(rad)) {
745  cor_rad_km = rad - radInit;
746  cor_rad_m = cor_rad_km * 1000.;
747  }
748  }
749  }
750  else if (m_coordTypeBundle == SurfacePoint::Latitudinal) {
751  cor_lat_dd = m_corrections(0) * RAD2DEG; // lat correction, decimal degs
752  cor_lon_dd = m_corrections(1) * RAD2DEG; // lon correction, decimal degs
753  cor_rad_m = m_corrections(2) * 1000.0; // radius correction, meters
754 
755  cor_lat_m = m_controlPoint->GetAdjustedSurfacePoint().LatitudeToMeters(m_corrections(0));
756  cor_lon_m = m_controlPoint->GetAdjustedSurfacePoint().LongitudeToMeters(m_corrections(1));
757  cor_rad_km = m_corrections(2);
758 
759  if (!IsSpecial(lat)) {
760  latInit = lat - cor_lat_dd;
761  }
762 
763  if (!IsSpecial(lon)) {
764  lonInit = lon - cor_lon_dd;
765  }
766 
767  if (!IsSpecial(rad)) {
768  radInit = rad - m_corrections(2); // km
769  }
770  }
771  else {
772  IString msg ="Unknown surface point coordinate type enum [" +
773  toString(m_coordTypeBundle) + "]." ;
775  }
776 
777  QString pointType = ControlPoint::PointTypeToString(type()).toUpper();
778 
779  QString output;
780 
781  output = QString(" Label: %1\nStatus: %2\n Rays: %3 of %4\n")
782  .arg(id())
783  .arg(pointType)
784  .arg(numGoodRays)
785  .arg(numRays);
786 
787  QString labels = "\n Point Initial Total Total "
788  " Final Initial Final\n"
789  "Coordinate Value Correction Correction "
790  " Value Accuracy Accuracy\n"
791  " (dd/dd/km) (dd/dd/km) (Meters) "
792  " (dd/dd/km) (Meters) (Meters)\n";
793  output += labels;
794 
795  SurfacePoint::CoordIndex idx = SurfacePoint::One;
796  output += QString(" LATITUDE%1%2%3%4%5%6\n")
797  .arg(formatValue(latInit, 17, 8)) // deg
798  .arg(formatValue(cor_lat_dd, 21, 8)) // deg
799  .arg(formatValue(cor_lat_m, 20, 8)) // m
800  .arg(formatValue(lat, 20, 8)) // deg
801  .arg(formatCoordAprioriSigmaString(idx, 18,8, true)) // m
802  .arg(formatCoordAdjustedSigmaString(idx, 18,8,errorPropagation)); // m
803 
804  idx = SurfacePoint::Two;
805  output += QString(" LONGITUDE%1%2%3%4%5%6\n")
806  .arg(formatValue(lonInit, 17, 8)) // deg
807  .arg(formatValue(cor_lon_dd, 21, 8)) // deg
808  .arg(formatValue(cor_lon_m, 20, 8)) // m
809  .arg(formatValue(lon, 20, 8)) // deg
810  .arg(formatCoordAprioriSigmaString(idx, 18,8, true)) // m
811  .arg(formatCoordAdjustedSigmaString(idx, 18,8,errorPropagation)); // m
812 
813  idx = SurfacePoint::Three;
814  output += QString(" RADIUS%1%2%3%4%5%6\n\n")
815  .arg(formatValue(radInit, 17, 8)) // km
816  .arg(formatValue(cor_rad_km, 21, 8)) // km
817  .arg(formatValue(cor_rad_m, 20, 8)) // m
818  .arg(formatValue(rad, 20, 8)) // km
819  .arg(formatCoordAprioriSigmaString(idx, 18,8, solveRadius)) // m
820  .arg(formatCoordAdjustedSigmaString(idx, 18,8, errorPropagation)); // m
821 
822  return output;
823  }
824 
825 
838  (bool errorPropagation) const {
839  int numRays = numberOfMeasures();
840  int numGoodRays = numRays - numberOfRejectedMeasures();
841  double X = m_controlPoint->GetAdjustedSurfacePoint().GetX().kilometers();
842  double Y = m_controlPoint->GetAdjustedSurfacePoint().GetY().kilometers();
843  double Z = m_controlPoint->GetAdjustedSurfacePoint().GetZ().kilometers();
844 
845  // Coordinate corrections are currently set in applyParameterCorrections. Units depend on
846  // coordTypeBundle coordinates (radians for Latitudinal coordinates and km for Rectangular).
847  // point corrections and initial sigmas.....
848  double cor_X_m = 0.; // X correction, meters
849  double cor_Y_m = 0.; // Y correction, meters
850  double cor_Z_m = 0.; // Z correction, meters
851  double XInit = Isis::Null;
852  double YInit = Isis::Null;
853  double ZInit = Isis::Null;
854 
855  if (m_coordTypeBundle == SurfacePoint::Latitudinal) {
856  double lat = m_controlPoint->GetAdjustedSurfacePoint().GetLatitude().degrees();
857  double latcor = m_corrections(0) * RAD2DEG;
858  double lon = m_controlPoint->GetAdjustedSurfacePoint().GetLongitude().degrees();
859  double loncor = m_corrections(1) * RAD2DEG;
860  double rad = m_controlPoint->GetAdjustedSurfacePoint().GetLocalRadius().kilometers();
861  double radcor = m_corrections(2);
862 
863  if (!IsSpecial(lat) && !IsSpecial(lon) && !IsSpecial(rad)) {
864  SurfacePoint latPoint(Latitude(lat-latcor, Angle::Degrees),
865  Longitude(lon-loncor, Angle::Degrees),
866  Distance(rad-radcor, Distance::Kilometers));
867  XInit = latPoint.GetX().kilometers();
868  YInit = latPoint.GetY().kilometers();
869  ZInit = latPoint.GetZ().kilometers();
870  cor_X_m = Isis::Null;
871  if (!IsSpecial(X)) {
872  cor_X_m = (X - XInit); // kilometers
873  }
874  cor_Y_m = Isis::Null;
875  if (!IsSpecial(Y)) {
876  cor_Y_m = (Y - YInit); // kilometers
877  }
878  cor_Z_m = Isis::Null;
879  if (!IsSpecial(Z)) {
880  cor_Z_m = (Z - ZInit); // kilometers
881  }
882  }
883  }
884  else if (m_coordTypeBundle == SurfacePoint::Rectangular) {
885  cor_X_m = m_corrections(0); // X correction, kilometers
886  cor_Y_m = m_corrections(1); // Y correction, kilometers
887  cor_Z_m = m_corrections(2); // Z correction, kilometers
888  if (!IsSpecial(X)) {
889  XInit = X - m_corrections(0); // km
890  }
891  if (!IsSpecial(Y)) {
892  YInit = Y - m_corrections(1); // km
893  }
894  if (!IsSpecial(Z)) {
895  ZInit = Z - m_corrections(2); // km
896  }
897  }
898  else {
899  IString msg ="Unknown surface point coordinate type enum [" +
900  toString(m_coordTypeBundle) + "]." ;
902  }
903 
904  QString pointType = ControlPoint::PointTypeToString(type()).toUpper();
905 
906  QString output;
907 
908  output = QString(" Label: %1\nStatus: %2\n Rays: %3 of %4\n")
909  .arg(id())
910  .arg(pointType)
911  .arg(numGoodRays)
912  .arg(numRays);
913 
914  QString labels =
915  "\n Point Initial Total Final "
916  " Initial Final\n"
917  " Coordinate Value Correction Value "
918  " Accuracy Accuracy\n"
919  " (km/km/km) (km) (km/km/km) "
920  " (Meters) (Meters)\n";
921  output += labels;
922 
923  SurfacePoint::CoordIndex idx = SurfacePoint::One;
924  output += QString(" BODY-FIXED-X%1%2%3%4%5\n")
925  .arg(formatValue(XInit, 17, 8)) // km
926  .arg(formatValue(cor_X_m, 20, 8)) // km
927  .arg(formatValue(X, 20, 8)) // km
928  .arg(formatCoordAprioriSigmaString(idx, 18,8, true)) // m
929  .arg(formatCoordAdjustedSigmaString(idx, 18,8,errorPropagation)); // m
930 
931  idx = SurfacePoint::Two;
932  output += QString(" BODY-FIXED-Y%1%2%3%4%5\n")
933  .arg(formatValue(YInit, 17, 8)) // km
934  .arg(formatValue(cor_Y_m, 20, 8)) // km
935  .arg(formatValue(Y, 20, 8)) // km
936  .arg(formatCoordAprioriSigmaString(idx, 18, 8, true)) // m
937  .arg(formatCoordAdjustedSigmaString(idx, 18,8,errorPropagation)); // m
938 
939  idx = SurfacePoint::Three;
940  output += QString(" BODY-FIXED-Z%1%2%3%4%5\n\n")
941  .arg(formatValue(ZInit, 17, 8)) // km
942  .arg(formatValue(cor_Z_m, 20, 8)) // km
943  .arg(formatValue(Z, 20, 8)) // km
944  // Set solveRadius to true to avoid limiting output information for Z.
945  // Set radius does not really apply to rectangular coordinates.
946  .arg(formatCoordAprioriSigmaString(idx, 18,8, true)) // m
947  .arg(formatCoordAdjustedSigmaString(idx, 18,8,errorPropagation)); // m
948 
949  return output;
950  }
951 
952 
963  QString BundleControlPoint::formatValue(double value, int fieldWidth, int precision) const {
964  QString output;
965  IsSpecial(value) ?
966  output = QString("%1").arg("Null", fieldWidth) :
967  output = QString("%1").arg(value, fieldWidth, 'f', precision);
968  return output;
969  }
970 
971 
985  QString BundleControlPoint::formatAprioriSigmaString(SurfacePoint::CoordIndex index,
986  int fieldWidth, int precision,
987  bool solveRadius) const {
988  QString aprioriSigmaStr;
989  QString pointType = ControlPoint::PointTypeToString(type()).toUpper();
990  if (pointType == "CONSTRAINED" || !solveRadius) {
991  pointType = "N/A";
992  }
993  double sigma = m_aprioriSigmas[int(index)];
994  if (IsSpecial(sigma)) { // if globalAprioriSigma <= 0 (including Isis::Null), then m_aprioriSigmas = Null
995  aprioriSigmaStr = QString("%1").arg(pointType, fieldWidth);
996  }
997  else {
998  aprioriSigmaStr = QString("%1").arg(sigma, fieldWidth, 'f', precision);
999  }
1000  return aprioriSigmaStr;
1001  }
1002 
1003 
1014 QString BundleControlPoint::formatCoordAprioriSigmaString(SurfacePoint::CoordIndex index,
1015  int fieldWidth,
1016  int precision,
1017  bool solveRadius) const {
1018  return formatAprioriSigmaString(index, fieldWidth, precision, solveRadius);
1019  }
1020 
1021 
1035  QString BundleControlPoint::formatAdjustedSigmaString(SurfacePoint::CoordIndex index,
1036  int fieldWidth, int precision,
1037  bool errorPropagation) const {
1038  QString adjustedSigmaStr;
1039 
1040  if (!errorPropagation) {
1041  adjustedSigmaStr = QString("%1").arg("N/A",fieldWidth);
1042  }
1043  else {
1044  double sigma = Isis::Null;
1046 
1047  if (IsSpecial(sigma)) {
1048  adjustedSigmaStr = QString("%1").arg("N/A",fieldWidth);
1049  }
1050  else {
1051  adjustedSigmaStr = QString("%1").arg(sigma, fieldWidth, 'f', precision);
1052  }
1053  }
1054 
1055  return adjustedSigmaStr;
1056  }
1057 
1058 
1071 QString BundleControlPoint::formatCoordAdjustedSigmaString(SurfacePoint::CoordIndex index,
1072  int fieldWidth, int precision,
1073  bool errorPropagation) const {
1074  return formatAdjustedSigmaString(index, fieldWidth, precision, errorPropagation);
1075  }
1076 
1077 
1085  (const BundleTargetBodyQsp target) {
1086  SurfacePoint surfacepoint = adjustedSurfacePoint();
1087  double pointLat = surfacepoint.GetLatitude().degrees();
1088  double pointLon = surfacepoint.GetLongitude().degrees();
1089  double pointRad = surfacepoint.GetLocalRadius().meters();
1090 
1091  // get point parameter corrections
1092  double latCorrection = m_nicVector[0];
1093  double lonCorrection = m_nicVector[1];
1094  double radCorrection = m_nicVector[2];
1095 
1096  // convert to degrees and apply
1097  pointLat += RAD2DEG * latCorrection;
1098  pointLon += RAD2DEG * lonCorrection;
1099 
1100  // Make sure updated values are still in valid range(0 to 360 for lon and -90 to 90 for lat)
1101  if (pointLat < -90.0) {
1102  pointLat = -180.0 - pointLat;
1103  pointLon = pointLon + 180.0;
1104  }
1105  if (pointLat > 90.0) {
1106  pointLat = 180.0 - pointLat;
1107  pointLon = pointLon + 180.0;
1108  }
1109  while (pointLon > 360.0) {
1110  pointLon = pointLon - 360.0;
1111  }
1112  while (pointLon < 0.0) {
1113  pointLon = pointLon + 360.0;
1114  }
1115 
1116  // convert to meters and apply
1117  pointRad += 1000. * radCorrection;
1118 
1119  // ken testing - if solving for target body mean radius, set radius to current
1120  // mean radius value
1121  // Only allow radius options for Latitudinal coordinates
1122  if (target && (target->solveMeanRadius() || target->solveTriaxialRadii()) ) {
1123  if (target->solveMeanRadius()) {
1124  surfacepoint.SetSphericalCoordinates(Latitude(pointLat, Angle::Degrees),
1125  Longitude(pointLon, Angle::Degrees),
1126  target->meanRadius());
1127  }
1128  else if (target->solveTriaxialRadii()) {
1129  Distance localRadius = target->localRadius(Latitude(pointLat, Angle::Degrees),
1130  Longitude(pointLon, Angle::Degrees));
1131  surfacepoint.SetSphericalCoordinates(Latitude(pointLat, Angle::Degrees),
1132  Longitude(pointLon, Angle::Degrees),
1133  localRadius);
1134  }
1135  }
1136  else {
1137  surfacepoint.SetSphericalCoordinates(Latitude(pointLat, Angle::Degrees),
1138  Longitude(pointLon, Angle::Degrees),
1139  Distance(pointRad, Distance::Meters));
1140  }
1141  // Reset the point now that it has been updated
1142  setAdjustedSurfacePoint(surfacepoint);
1143  }
1144 
1145 
1152  SurfacePoint surfacepoint = adjustedSurfacePoint();
1153  double pointX = surfacepoint.GetX().meters();
1154  double pointY = surfacepoint.GetY().meters();
1155  double pointZ = surfacepoint.GetZ().meters();
1156 
1157  // get point parameter corrections
1158  double XCorrection = m_nicVector[0];
1159  double YCorrection = m_nicVector[1];
1160  double ZCorrection = m_nicVector[2];
1161 
1162  // Convert corrections to meters and apply
1163  pointX += 1000. * XCorrection;
1164  pointY += 1000. * YCorrection;
1165  pointZ += 1000. * ZCorrection;
1166 
1167  surfacepoint.SetRectangularCoordinates(
1171  // Reset the point now that it has been updated
1172  setAdjustedSurfacePoint(surfacepoint);
1173  }
1174 }
This class defines a body-fixed surface point.
Definition: SurfacePoint.h:148
double meters() const
Get the distance in meters.
Definition: Distance.cpp:97
void ZeroNumberOfRejectedMeasures()
Initialize the number of rejected measures to 0.
SparseBlockMatrix.
QString formatBundleRectangularOutputDetailString(bool errorPropagation) const
Formats a detailed output string table for this Rectangular BundleControlPoint.
ControlPoint::PointType type() const
Accesses BundleControlPoint&#39;s type.
const double Null
Value for an Isis Null pixel.
Definition: SpecialPixel.h:110
QString formatAdjustedSigmaString(SurfacePoint::CoordIndex, int fieldWidth, int precision, bool errorPropagation) const
Formats the adjusted sigma value indicated by the given type code.
The distance is being specified in kilometers.
Definition: Displacement.h:54
QString formatBundleOutputDetailString(bool errorPropagation, bool solveRadius=false) const
Formats a detailed output string table for this BundleControlPoint.
void SetNumberOfRejectedMeasures(int numRejected)
Set (update) the number of rejected measures for the control point.
boost::numeric::ublas::bounded_vector< double, 3 > m_adjustedSigmas
adjusted sigmas for point parameters
A Constrained point is a Control Point whose lat/lon/radius is somewhat established and should not be...
Definition: ControlPoint.h:391
double residualRms() const
Gets the root-mean-square (rms) of the BundleControlPoint&#39;s residuals.
A Fixed point is a Control Point whose lat/lon is well established and should not be changed...
Definition: ControlPoint.h:386
BundleControlPoint(BundleSettingsQsp bundleSettings, ControlPoint *point)
Constructs a BundleControlPoint object from a ControlPoint.
SurfacePoint::CoordinateType coordTypeBundle() const
Accesses BundleControlPoint&#39;s control point coordinate type for the bundle adjustment.
bool IsCoord3Constrained()
Return bool indicating if 3rd coordinate is Constrained or not.
QString formatAprioriSigmaString(SurfacePoint::CoordIndex index, int fieldWidth, int precision, bool solveRadius=false) const
Formats the apriori sigma value indicated by the given type code.
void setRejected(bool reject)
Sets this BundleControlPoint to rejected or not rejected.
int GetNumberOfRejectedMeasures() const
Get the number of rejected measures on the control point.
const ControlMeasure * GetMeasure(QString serialNumber) const
Get a control measure based on its cube&#39;s serial number.
BundleMeasureQsp addMeasure(ControlMeasure *controlMeasure)
Creates a BundleMeasure from the given ControlMeasure and appends it to this BundleControlPoint&#39;s mea...
void applyParameterCorrections(LinearAlgebra::Vector imageSolution, SparseBlockMatrix &sparseNormals, const BundleTargetBodyQsp target)
Apply the parameter corrections to the bundle control point.
double MetersToLatitude(double latLength)
This method returns an angular measure of a distance in the direction of and relative to the latitude...
void copy(const BundleControlPoint &src)
Copies given BundleControlPoint to this BundleControlPoint.
PointType
These are the valid &#39;types&#39; of point.
Definition: ControlPoint.h:379
double meters() const
Get the displacement in meters.
void SetRectangularCoordinates(const Displacement &x, const Displacement &y, const Displacement &z)
Set surface point in rectangular coordinates.
This class is designed to encapsulate the concept of a Latitude.
Definition: Latitude.h:63
void computeResiduals()
Computes the residuals for this BundleControlPoint.
boost::numeric::ublas::bounded_vector< double, 3 > & weights()
Accesses the 3 dimensional ordered vector of weight values associated with coordinate1, coordinate2, and coordinate3.
The distance is being specified in kilometers.
Definition: Distance.h:58
boost::numeric::ublas::bounded_vector< double, 3 > & aprioriSigmas()
Accesses the 3 dimensional ordered vector of apriori sigmas (apriori coordinate1, apriori coordinate2...
A container class for a ControlMeasure.
Definition: BundleMeasure.h:69
void updateAdjustedSurfacePointRectangularly()
pointer to the control point object this represents
QString formatValue(double value, int fieldWidth, int precision) const
Formats the given double precision value using the specified field width and precision.
QString toString(bool boolToConvert)
Global function to convert a boolean to a string.
Definition: IString.cpp:226
This error is for when a programmer made an API call that was illegal.
Definition: IException.h:162
double kilometers() const
Get the distance in kilometers.
Definition: Distance.cpp:118
The distance is being specified in meters.
Definition: Displacement.h:52
Distance measurement, usually in meters.
Definition: Distance.h:47
double GetResidualRms() const
Get rms of residuals.
Latitude GetLatitude() const
Return the body-fixed latitude for the surface point.
boost::numeric::ublas::bounded_vector< double, 3 > & adjustedSigmas()
Accesses the 3 dimensional ordered vector of adjusted sigmas (adjusted coordinate1, adjusted coordinate2, adjusted coordinate3).
double GetCoord(CoordinateType type, CoordIndex index, CoordUnits units)
This method returns a coordinate of a SurfacePoint.
double degrees() const
Get the angle in units of Degrees.
Definition: Angle.h:249
Unless noted otherwise, the portions of Isis written by the USGS are public domain.
Status SetAdjustedSurfacePoint(SurfacePoint newSurfacePoint)
Set or update the surface point relating to this control point.
boost::numeric::ublas::bounded_vector< double, 3 > & corrections()
Accesses the 3 dimensional ordered vector of correction values associated with coord1, coord2, and coord 3 (latitude, longitude, and radius or X, Y, and Z.
A Free point is a Control Point that identifies common measurements between two or more cubes...
Definition: ControlPoint.h:399
boost::numeric::ublas::bounded_vector< double, 3 > & nicVector()
Accesses the 3 dimensional ordered NIC vector.
Planetocentric latitudinal (lat/lon/rad) coordinates.
Definition: SurfacePoint.h:156
boost::numeric::ublas::bounded_vector< double, 3 > m_nicVector
array of NICs (see Brown, 1976)
This class is designed to encapsulate the concept of a Longitude.
Definition: Longitude.h:52
ControlPoint * rawControlPoint() const
Accessor for the raw ControlPoint object used for this BundleControlPoint.
boost::numeric::ublas::vector< double > Vector
Definition for an Isis::LinearAlgebra::Vector of doubles.
void updateAdjustedSurfacePointLatitudinally(const BundleTargetBodyQsp target)
Apply the parameter corrections to the bundle control point latitudinally.
Degrees are generally considered more human readable, 0-360 is one circle, however most math does not...
Definition: Angle.h:73
QString formatBundleLatitudinalOutputDetailString(bool errorPropagation, bool solveRadius=false) const
Formats a detailed output string table for this Latitudinal BundleControlPoint.
int numberOfRejectedMeasures() const
Accesses the number of rejected measures for this BundleControlPoint.
boost::numeric::ublas::bounded_vector< double, 3 > m_aprioriSigmas
apriori sigmas for point parameters
QString GetId() const
Return the Id of the control point.
double LongitudeToMeters(double longitude) const
This method returns a length in meters version of a delta longitude angle in radians relative to the ...
double MetersToLongitude(double lonLength)
This method returns an angular measure in radians of a distance in the direction of and relative to t...
bool IsCoord2Constrained()
Return bool indicating if 2nd coordinate is Constrained or not.
void zeroNumberOfRejectedMeasures()
Resets the number of rejected measures for this BundleControlPoint to zero.
void setAdjustedSurfacePoint(SurfacePoint surfacePoint)
Sets the adjusted surface point for this BundleControlPoint.
#define _FILEINFO_
Macro for the filename and line number.
Definition: IException.h:40
SurfacePoint::CoordinateType m_coordTypeReports
BundleControlPoint coordinate type.
boost::numeric::ublas::bounded_vector< double, 3 > m_weights
weights for point parameters
This class holds information about a control point that BundleAdjust needs to run correctly...
bool IsCoord1Constrained()
Return bool indicating if 1st coordinate is Constrained or not.
~BundleControlPoint()
Destructor for BundleControlPoint.
bool IsSpecial(const double d)
Returns if the input pixel is special.
Definition: SpecialPixel.h:212
A single control point.
Definition: ControlPoint.h:369
boost::numeric::ublas::bounded_vector< double, 3 > m_corrections
corrections to point parameters
void productAlphaAV(double alpha, SparseBlockMatrix &sparseNormals, LinearAlgebra::Vector &v1)
Perform the matrix multiplication v2 = alpha ( Q x v1 ).
SurfacePoint::CoordinateType coordTypeReports() const
Accesses BundleControlPoint&#39;s coordinate type for reports.
Body-fixed rectangular x/y/z coordinates.
Definition: SurfacePoint.h:157
CoordinateType
Defines the coordinate typ, units, and coordinate index for some of the output methods.
Definition: SurfacePoint.h:155
Distance GetSigmaDistance(CoordinateType type, CoordIndex index)
This method returns a sigma of a SurfacePoint coordinate as a Distance.
const double DEG2RAD
Multiplier for converting from degrees to radians.
Definition: Constants.h:59
QString formatCoordAdjustedSigmaString(SurfacePoint::CoordIndex, int fieldWidth, int precision, bool errorPropagation) const
Formats the adjusted coordinate sigma value.
Status SetRejected(bool rejected)
Set the jigsawRejected state.
double LatitudeToMeters(double latitude) const
This method returns a Displacement of an Angle relative to the current SurfacePoint latitude...
Unless noted otherwise, the portions of Isis written by the USGS are public domain.
SurfacePoint adjustedSurfacePoint() const
Accesses the adjusted SurfacePoint associated with this BundleControlPoint.
Longitude GetLongitude() const
Return the body-fixed longitude for the surface point.
int numberOfMeasures() const
Accesses number of measures associated with this BundleControlPoint.
void setSigmaWeightFromGlobals(double gSigma, int index)
Sets the member sigmas and weights from a global sigma.
Displacement is a signed length, usually in meters.
Definition: Displacement.h:43
QSharedPointer< BundleMeasure > BundleMeasureQsp
Definition for BundleMeasureQsp, a shared pointer to a BundleMeasure.
QString id() const
Accesses the Point ID associated with this BundleControlPoint.
Isis exception class.
Definition: IException.h:107
Adds specific functionality to C++ strings.
Definition: IString.h:181
QString formatCoordAprioriSigmaString(SurfacePoint::CoordIndex index, int fieldWidth, int precision, bool solveRadius=false) const
Formats the apriori coordinate 1 (latitude or X) sigma value.
Namespace for ISIS/Bullet specific routines.
Definition: Apollo.h:31
QString formatBundleOutputSummaryString(bool errorPropagation) const
Formats an output summary string for this BundleControlPoint.
SparseBlockRowMatrix m_cholmodQMatrix
The CholMod matrix associated with this point.
a control measurement
bool isRejected() const
Method used to determine whether this control point is rejected.
The distance is being specified in meters.
Definition: Distance.h:56
double GetWeight(CoordinateType type, CoordIndex index)
This method returns the weight of a SurfacePoint coordinate Note: At this time a units argument is no...
const double RAD2DEG
Multiplier for converting from radians to degrees.
Definition: Constants.h:60
double kilometers() const
Get the displacement in kilometers.
PointType GetType() const
void SetSphericalCoordinates(const Latitude &lat, const Longitude &lon, const Distance &radius)
Update spherical coordinates (lat/lon/radius)
Distance GetLocalRadius() const
Return the radius of the surface point.
void setWeights(const BundleSettingsQsp settings)
Sets the weights using the given BundleSettings QSharedPointer and a conversion value for meters to r...
void setNumberOfRejectedMeasures(int numRejected)
Sets the number of rejected measures for this BundleControlPoint.
SparseBlockRowMatrix & cholmodQMatrix()
Accesses the CholMod matrix associated with this BundleControlPoint.
SparseBlockRowMatrix.
static QString PointTypeToString(PointType type)
Obtain a string representation of a given PointType.
Status ComputeResiduals()
This method computes the BundleAdjust residuals for a point.
Unless noted otherwise, the portions of Isis written by the USGS are public domain.