Isis 3 Programmer Reference
SurfacePoint.cpp
1 #include "SurfacePoint.h"
2 
3 #include <SpiceUsr.h>
4 
5 #include "IException.h"
6 #include "IString.h"
7 #include "Latitude.h"
8 #include "Longitude.h"
9 #include "SpecialPixel.h"
10 
11 using namespace boost::numeric::ublas;
12 using namespace std;
13 
14 namespace Isis {
15 
20  SurfacePoint::SurfacePoint() {
21  InitCovariance();
22  InitPoint();
23  }
24 
29  SurfacePoint::SurfacePoint(const SurfacePoint &other) {
30  p_localRadius = other.p_localRadius;
31 
32  if(other.p_x) {
33  p_x = new Displacement(*other.p_x);
34  }
35  else {
36  p_x = NULL;
37  }
38 
39  if(other.p_y) {
40  p_y = new Displacement(*other.p_y);
41  }
42  else {
43  p_y = NULL;
44  }
45 
46  if(other.p_z) {
47  p_z = new Displacement(*other.p_z);
48  }
49  else {
50  p_z = NULL;
51  }
52 
53  if(other.p_rectCovar) {
54  p_rectCovar = new symmetric_matrix<double, upper>(*other.p_rectCovar);
55  }
56  else {
57  p_rectCovar = NULL;
58  }
59 
60  if(other.p_sphereCovar) {
61  p_sphereCovar = new symmetric_matrix<double, upper>(*other.p_sphereCovar);
62  }
63  else {
64  p_sphereCovar = NULL;
65  }
66  }
67 
68 
76  SurfacePoint::SurfacePoint(const Latitude &lat, const Longitude &lon,
77  const Distance &radius) {
78  InitCovariance();
79  InitPoint();
80  SetSphericalPoint(lat, lon, radius);
81  }
82 
83 
84 
99  SurfacePoint::SurfacePoint(const Latitude &lat, const Longitude &lon,
100  const Distance &radius, const Angle &latSigma, const Angle &lonSigma,
101  const Distance &radiusSigma) {
102  InitCovariance();
103  InitPoint();
104  SetSpherical(lat, lon, radius, latSigma, lonSigma, radiusSigma);
105  }
106 
107 
113  SurfacePoint::SurfacePoint(const Latitude &lat, const Longitude &lon,
114  const Distance &radius, const symmetric_matrix<double, upper> &covar) {
115  InitCovariance();
116  InitPoint();
117  SetSpherical(lat, lon, radius, covar);
118  }
119 
120 
128  SurfacePoint::SurfacePoint(const Displacement &x, const Displacement &y,
129  const Displacement &z) {
130  InitCovariance();
131  InitPoint();
132  SetRectangular(x, y, z);
133  }
134 
135 
150  SurfacePoint::SurfacePoint(const Displacement &x, const Displacement &y,
151  const Displacement &z, const Distance &xSigma, const Distance &ySigma,
152  const Distance &zSigma) {
153  InitCovariance();
154  InitPoint();
155  SetRectangular(x, y, z, xSigma, ySigma, zSigma);
156  }
157 
158 
168  SurfacePoint::SurfacePoint(const Displacement &x, const Displacement &y,
169  const Displacement &z, const symmetric_matrix<double, upper> &covar) {
170  InitCovariance();
171  InitPoint();
172  SetRectangular(x, y, z, covar);
173  }
174 
175 
180  SurfacePoint::~SurfacePoint() {
181  FreeAllocatedMemory();
182  }
183 
184 
189  void SurfacePoint::InitCovariance() {
190  p_rectCovar = NULL;
191  p_sphereCovar = NULL;
192  }
193 
194 
199  void SurfacePoint::InitPoint() {
200  p_x = NULL;
201  p_y = NULL;
202  p_z = NULL;
203  p_localRadius = Distance();
204  }
205 
206 
219  void SurfacePoint::SetRectangularPoint(const Displacement &x,
220  const Displacement &y, const Displacement &z) {
221 
222  if (!x.isValid() || !y.isValid() || !z.isValid()) {
223  IString msg = "x, y, and z must be set to valid displacements. One or "
224  "more coordinates have been set to an invalid displacement.";
225  throw IException(IException::User, msg, _FILEINFO_);
226  }
227 
228  if(p_x) {
229  *p_x = x;
230  }
231  else {
232  p_x = new Displacement(x);
233  }
234 
235  if(p_y) {
236  *p_y = y;
237  }
238  else {
239  p_y = new Displacement(y);
240  }
241 
242  if(p_z) {
243  *p_z = z;
244  }
245  else {
246  p_z = new Displacement(z);
247  }
248 
249  // Added 07-30-2018 to avoid trying to set an invalid SurfacePoint. This breaks the ringspt and cnetnewradii tests. Also the pole...mean test.
250  // if (!(this->Valid())) {
251  // IString msg = "All coodinates of the point are zero. At least one coordinate"
252  // " must be nonzero to be a valid SurfacePoint.";
253  // throw IException(IException::User, msg, _FILEINFO_);
254  // }
255 
256  if (!p_localRadius.isValid()) {
257  ComputeLocalRadius();
258  p_localRadius = GetLocalRadius();
259  }
260  }
261 
262 
277  void SurfacePoint::SetRectangular(const Displacement &x,
278  const Displacement &y, const Displacement &z, const Distance &xSigma,
279  const Distance &ySigma, const Distance &zSigma) {
280 
281  // Wipe out current local radius to ensure it will be recalculated in SetRectangularPoint
282  p_localRadius = Distance();
283 
284  SetRectangularPoint(x, y, z);
285 
286  if (xSigma.isValid() && ySigma.isValid() && zSigma.isValid())
287  SetRectangularSigmas(xSigma, ySigma, zSigma);
288  }
289 
290 
302  void SurfacePoint::SetRectangular(Displacement x, Displacement y, Displacement z,
303  const symmetric_matrix<double,upper>& covar) {
304  // Wipe out current local radius to ensure it will be recalulated in SetRectangularPoint
305  p_localRadius = Distance();
306 
307  SetRectangularPoint(x, y, z);
308  SetRectangularMatrix(covar);
309  }
310 
311 
321  void SurfacePoint::SetRectangularCoordinates(const Displacement &x,
322  const Displacement &y,
323  const Displacement &z) {
324  // Wipe out current local radius to ensure it will be recalculated in SetRectangularPoint
325  p_localRadius = Distance();
326 
327  SetRectangularPoint(x, y, z);
328  }
329 
330 
344  void SurfacePoint::SetRectangularSigmas(const Distance &xSigma,
345  const Distance &ySigma,
346  const Distance &zSigma) {
347  if (!xSigma.isValid() || !ySigma.isValid() || !zSigma.isValid()) {
348  IString msg = "x sigma, y sigma , and z sigma must be set to valid "
349  "distances. One or more sigmas have been set to an invalid distance.";
350  throw IException(IException::User, msg, _FILEINFO_);
351  }
352 
353  symmetric_matrix<double,upper> covar(3);
354  covar.clear();
355  covar(0,0) = xSigma.meters() * xSigma.meters();
356  covar(1,1) = ySigma.meters() * ySigma.meters();
357  covar(2,2) = zSigma.meters() * zSigma.meters();
358  SetRectangularMatrix(covar, Meters);
359  }
360 
361 
370  void SurfacePoint::SetRectangularMatrix(const symmetric_matrix<double, upper> &covar,
371  SurfacePoint::CoordUnits units) {
372 
373  // Make sure the point is set first
374  if (!Valid()) {
375  IString msg = "A point must be set before a variance/covariance matrix "
376  "can be set.";
377  throw IException(IException::Programmer, msg, _FILEINFO_);
378  }
379 
380  // Make sure the units are valid
381  if (units != Kilometers && units != Meters) {
382  IString msg = "Units must be Kilometers or Meters";
383  throw IException(IException::Programmer, msg, _FILEINFO_);
384  }
385 
386  // covar units are km**2
387  if (p_rectCovar) {
388  *p_rectCovar = covar;
389  }
390  else {
391  p_rectCovar = new symmetric_matrix<double, upper>(covar);
392  }
393 
394  if (units == Meters) {
395  // Convert input matrix to km to hold in memory
396  (*p_rectCovar)(0,0) = covar(0,0)/1.0e6;
397  (*p_rectCovar)(0,1) = covar(0,1)/1.0e6;
398  (*p_rectCovar)(0,2) = covar(0,2)/1.0e6;
399  (*p_rectCovar)(1,1) = covar(1,1)/1.0e6;
400  (*p_rectCovar)(1,2) = covar(1,2)/1.0e6;
401  (*p_rectCovar)(2,2) = covar(2,2)/1.0e6;
402  }
403 
404  SpiceDouble rectMat[3][3];
405 
406  // Compute the local radius of the surface point in meters. We will convert to km before saving the matrix.
407  double x2 = p_x->meters() * p_x->meters();
408  double y2 = p_y->meters() * p_y->meters();
409  double z = p_z->meters();
410  double radius = sqrt(x2 + y2 + z*z);
411 
412  // *** TODO *** Replace this section with LinearAlgebra multiply calls and avoid having to create a Spice matrix
413  rectMat[0][0] = covar(0,0);
414  rectMat[0][1] = rectMat[1][0] = covar(0,1);
415  rectMat[0][2] = rectMat[2][0] = covar(0,2);
416  rectMat[1][1] = covar(1,1);
417  rectMat[1][2] = rectMat[2][1] = covar(1,2);
418  rectMat[2][2] = covar(2,2);
419 
420  // Compute the Jacobian in meters. Don't deal with unit mismatch yet to preserve precision.
421  SpiceDouble J[3][3];
422  double zOverR = p_z->meters() / radius;
423  double r2 = radius*radius;
424  double denom = r2*radius*sqrt(1.0 - (zOverR*zOverR));
425  J[0][0] = -p_x->meters() * p_z->meters() / denom;
426  J[0][1] = -p_y->meters() * p_z->meters() / denom;
427  J[0][2] = (r2 - p_z->meters() * p_z->meters()) / denom;
428  J[1][0] = -p_y->meters() / (x2 + y2);
429  J[1][1] = p_x->meters() / (x2 + y2);
430  J[1][2] = 0.0;
431  J[2][0] = p_x->meters() / radius; // This row is unitless
432  J[2][1] = p_y->meters() / radius;
433  J[2][2] = p_z->meters() / radius;
434 
435  if(!p_sphereCovar)
436  p_sphereCovar = new symmetric_matrix<double, upper>(3);
437 
438  SpiceDouble mat[3][3];
439  mxm_c (J, rectMat, mat);
440  mxmt_c (mat, J, mat);
441  if (units == Kilometers) {
442  // Now take care of unit mismatch between rect matrix in km and Jacobian in m
443  (*p_sphereCovar)(0,0) = mat[0][0] * 1.0e6;
444  (*p_sphereCovar)(0,1) = mat[0][1] * 1.0e6;
445  (*p_sphereCovar)(0,2) = mat[0][2] * 1000.0;
446  (*p_sphereCovar)(1,1) = mat[1][1] * 1.0e6;
447 
448  (*p_sphereCovar)(1,2) = mat[1][2] * 1000.0;
449  (*p_sphereCovar)(2,2) = mat[2][2];
450  }
451  else { // (units == Meters)
452  // Convert matrix lengths from m to km
453  (*p_sphereCovar)(0,0) = mat[0][0];
454  (*p_sphereCovar)(0,1) = mat[0][1];
455  (*p_sphereCovar)(0,2) = mat[0][2] / 1000.0;
456  (*p_sphereCovar)(1,1) = mat[1][1];
457  (*p_sphereCovar)(1,2) = mat[1][2] / 1000.0;
458  (*p_sphereCovar)(2,2) = mat[2][2] / 1.0e6;
459  }
460  }
461 
462 
474  void SurfacePoint::SetSphericalPoint(const Latitude &lat,
475  const Longitude &lon,
476  const Distance &radius) {
477 // Is error checking necessary or does Latitude, Longitude, and Distance handle it?????
478  if (!lat.isValid() || !lon.isValid() || !radius.isValid()) {
479  IString msg = "Latitude, longitude, or radius is an invalid value.";
480  throw IException(IException::User, msg, _FILEINFO_);
481  }
482 
483  SpiceDouble dlat = (double) lat.radians();
484  SpiceDouble dlon = (double) lon.radians();
485  SpiceDouble dradius = radius.kilometers();
486 
487  SpiceDouble rect[3];
488  latrec_c ( dradius, dlon, dlat, rect);
489 
490  // Set local radius now since we have it and avoid calculating it later
491  p_localRadius = radius;
492 
493  SetRectangularPoint(Displacement(rect[0], Displacement::Kilometers),
494  Displacement(rect[1], Displacement::Kilometers),
495  Displacement(rect[2], Displacement::Kilometers));
496  }
497 
498 
512  void SurfacePoint::SetSpherical(const Latitude &lat, const Longitude &lon,
513  const Distance &radius, const Angle &latSigma, const Angle &lonSigma,
514  const Distance &radiusSigma) {
515 
516  SetSphericalPoint(lat, lon, radius);
517 
518  if (latSigma.isValid() && lonSigma.isValid() && radiusSigma.isValid())
519  SetSphericalSigmas(latSigma, lonSigma, radiusSigma);
520  }
521 
522 
532  void SurfacePoint::SetSpherical(const Latitude &lat, const Longitude &lon,
533  const Distance &radius, const symmetric_matrix<double, upper> &covar) {
534  SetSphericalPoint(lat, lon, radius);
535  SetSphericalMatrix(covar);
536  }
537 
538 
547  void SurfacePoint::SetSphericalCoordinates(const Latitude &lat,
548  const Longitude &lon, const Distance &radius) {
549  SetSphericalPoint(lat, lon, radius);
550  }
551 
552 
567  void SurfacePoint::SetSphericalSigmas(const Angle &latSigma,
568  const Angle &lonSigma,
569  const Distance &radiusSigma) {
570  if (latSigma.isValid() && lonSigma.isValid() && radiusSigma.isValid()) {
571  symmetric_matrix<double,upper> covar(3);
572  covar.clear();
573 
574  double sphericalCoordinate;
575  sphericalCoordinate = (double) latSigma.radians();
576  covar(0,0) = sphericalCoordinate*sphericalCoordinate;
577  sphericalCoordinate = (double) lonSigma.radians();
578  covar(1,1) = sphericalCoordinate*sphericalCoordinate;
579  sphericalCoordinate = (double) radiusSigma.meters();
580  covar(2,2) = sphericalCoordinate*sphericalCoordinate;
581 
582  // Use default set units for radius of meters*meters
583  SetSphericalMatrix(covar);
584  }
585  else {
586  delete p_sphereCovar;
587  p_sphereCovar = NULL;
588 
589  delete p_rectCovar;
590  p_rectCovar = NULL;
591  }
592  }
593 
594 
611  void SurfacePoint::SetSphericalSigmasDistance(const Distance &latSigma,
612  const Distance &lonSigma, const Distance &radiusSigma) {
613 
614  if (!Valid()) {
615  IString msg = "Cannot set spherical sigmas on an invalid surface point";
616  throw IException(IException::Programmer, msg, _FILEINFO_);
617  }
618 
619  SetSphericalSigmas(Angle(MetersToLatitude(latSigma.meters()),Angle::Radians),
620  Angle(MetersToLongitude(lonSigma.meters()), Angle::Radians), radiusSigma);
621  }
622 
623 
631  void SurfacePoint::SetSphericalMatrix(
632  const symmetric_matrix<double, upper> & covar,
633  SurfacePoint::CoordUnits units) {
634 
635  // Make sure the point is set first
636  if (!Valid()) {
637  IString msg = "A point must be set before a variance/covariance matrix "
638  "can be set.";
639  throw IException(IException::Programmer, msg, _FILEINFO_);
640  }
641 
642  // Make sure the units are valid
643  if (units != Kilometers && units != Meters) {
644  IString msg = "Units must be Kilometers or Meters";
645  throw IException(IException::Programmer, msg, _FILEINFO_);
646  }
647 
648  // Get the radius of the point in the same units as the input matrix when saving the input matrix
649  double radius = (double) GetLocalRadius().kilometers();
650 
651  // Save the spherical matrix in km and km**2
652  if (p_sphereCovar) {
653  *p_sphereCovar = covar;
654  }
655  else {
656  p_sphereCovar = new symmetric_matrix<double, upper>(covar);
657  }
658 
659  if (units == Meters) {
660  // Convert input matrix to km to store
661  (*p_sphereCovar)(0,0) = covar(0,0);
662  (*p_sphereCovar)(0,1) = covar(0,1);
663  (*p_sphereCovar)(0,2) = covar(0,2) / 1000.;
664  (*p_sphereCovar)(1,1) = covar(1,1);
665  (*p_sphereCovar)(1,2) = covar(1,2) / 1000.;
666  (*p_sphereCovar)(2,2) = covar(2,2) / 1.0e6;
667  radius = (double) GetLocalRadius().meters();
668  }
669 
670  // ***TODO*** Consider using LinearAlgebra matrix multiply and avoid creating SpiceDouble matrix.
671  SpiceDouble sphereMat[3][3];
672 
673  sphereMat[0][0] = covar(0,0);
674  sphereMat[0][1] = sphereMat[1][0] = covar(0,1);
675  sphereMat[0][2] = sphereMat[2][0] = covar(0,2);
676  sphereMat[1][1] = covar(1,1);
677  sphereMat[1][2] = sphereMat[2][1] = covar(1,2);
678  sphereMat[2][2] = covar(2,2);
679 
680 // std::cout<<"Ocovar = "<<sphereMat[0][0]<<" "<<sphereMat[0][1]<<" "<<sphereMat[0][2]<<std::endl
681 // <<" "<<sphereMat[1][0]<<" "<<sphereMat[1][1]<<" "<<sphereMat[1][2]<<std::endl
682 // <<" "<<sphereMat[2][0]<<" "<<sphereMat[2][1]<<" "<<sphereMat[2][2]<<std::endl;
683 
684  // Get the lat/lon/radius of the point
685  double lat = (double) GetLatitude().radians();
686  double lon = (double) GetLongitude().radians();
687 
688  // Compute the Jacobian in same units as input matrix.
689  SpiceDouble J[3][3];
690  double cosPhi = cos(lat);
691  double sinPhi = sin(lat);
692  double cosLamda = cos(lon);
693  double sinLamda = sin(lon);
694  double rcosPhi = radius*cosPhi;
695  double rsinPhi = radius*sinPhi;
696  J[0][0] = -rsinPhi * cosLamda;
697  J[0][1] = -rcosPhi * sinLamda;
698  J[0][2] = cosPhi * cosLamda;
699  J[1][0] = -rsinPhi * sinLamda;
700  J[1][1] = rcosPhi * cosLamda;
701  J[1][2] = cosPhi * sinLamda;
702  J[2][0] = rcosPhi;
703  J[2][1] = 0.0;
704  J[2][2] = sinPhi;
705 
706  if(!p_rectCovar)
707  p_rectCovar = new symmetric_matrix<double, upper>(3);
708 
709  SpiceDouble mat[3][3];
710  mxm_c (J, sphereMat, mat);
711  mxmt_c (mat, J, mat);
712 
713  if (units == Kilometers) {
714  (*p_rectCovar)(0,0) = mat[0][0];
715  (*p_rectCovar)(0,1) = mat[0][1];
716  (*p_rectCovar)(0,2) = mat[0][2];
717  (*p_rectCovar)(1,1) = mat[1][1];
718  (*p_rectCovar)(1,2) = mat[1][2];
719  (*p_rectCovar)(2,2) = mat[2][2];
720  }
721  else { // (units == Meters)
722  // Convert to km
723  (*p_rectCovar)(0,0) = mat[0][0]/1.0e6;
724  (*p_rectCovar)(0,1) = mat[0][1]/1.0e6;
725  (*p_rectCovar)(0,2) = mat[0][2]/1.0e6;
726  (*p_rectCovar)(1,1) = mat[1][1]/1.0e6;
727  (*p_rectCovar)(1,2) = mat[1][2]/1.0e6;
728  (*p_rectCovar)(2,2) = mat[2][2]/1.0e6;
729  }
730 // std::cout<<"Rcovar = "<<p_rectCovar(0,0)<<" "<<p_rectCovar(0,1)<<" "<<p_rectCovar(0,2)<<std::endl
731 // <<" "<<p_rectCovar(1,0)<<" "<<p_rectCovar(1,1)<<" "<<p_rectCovar(1,2)<<std::endl
732 // <<" "<<p_rectCovar(2,0)<<" "<<p_rectCovar(2,1)<<" "<<p_rectCovar(2,2)<<std::endl;
733  }
734 
735 
743  void SurfacePoint::SetMatrix(CoordinateType type, const symmetric_matrix<double, upper>& covar) {
744 
745  switch (type) {
746  case Latitudinal:
747  SetSphericalMatrix(covar, SurfacePoint::Kilometers);
748  break;
749  case Rectangular:
750  SetRectangularMatrix(covar, SurfacePoint::Kilometers);
751  break;
752  default:
753  IString msg ="Unknown surface point coordinate type enum [" + toString(type) + "]." ;
754  throw IException(IException::Programmer, msg, _FILEINFO_);
755  }
756  }
757 
758 
767  std::vector<double> SurfacePoint::Partial(CoordinateType type, CoordIndex index) {
768  std::vector<double> derivative(3);
769  switch (type) {
770  case Latitudinal:
771  derivative = LatitudinalDerivative(index);
772  break;
773  case Rectangular:
774  derivative = RectangularDerivative(index);
775  break;
776  default:
777  IString msg ="Unknown surface point coordinate type enum [" + toString(type) + "]." ;
778  throw IException(IException::Programmer, msg, _FILEINFO_);
779  }
780  return derivative;
781  }
782 
783 
792  std::vector<double> SurfacePoint::LatitudinalDerivative(CoordIndex index) {
793  std::vector<double> derivative(3);
794  double rlat = GetLatitude().radians();
795  double rlon = GetLongitude().radians();
796  double sinLon = sin(rlon);
797  double cosLon = cos(rlon);
798  double sinLat = sin(rlat);
799  double cosLat = cos(rlat);
800  double radkm = GetLocalRadius().kilometers();
801 
802  switch (index) {
803  case One:
804  derivative[0] = -radkm * sinLat * cosLon;
805  derivative[1] = -radkm * sinLon * sinLat;
806  derivative[2] = radkm * cosLat;
807  break;
808  case Two:
809  derivative[0] = -radkm * cosLat * sinLon;
810  derivative[1] = radkm * cosLat * cosLon;
811  derivative[2] = 0.0;
812  break;
813  case Three:
814  derivative[0] = cosLon * cosLat;
815  derivative[1] = sinLon * cosLat;
816  derivative[2] = sinLat;
817  break;
818  default:
819  IString msg = "Invalid coordinate index (must be less than 3).";
820  throw IException(IException::Programmer, msg, _FILEINFO_);
821  }
822  return derivative;
823  }
824 
825 
834  std::vector<double> SurfacePoint::RectangularDerivative(CoordIndex index) {
835  std::vector<double> derivative(3,0.0);
836 
837  switch (index) {
838  case One:
839  derivative[0] = 1.;
840  break;
841  case Two:
842  derivative[1] = 1.;
843  break;
844  case Three:
845  derivative[2] = 1.;
846  break;
847  default:
848  IString msg = "Invalid coordinate index (must be less than 3).";
849  throw IException(IException::Programmer, msg, _FILEINFO_);
850  }
851  return derivative;
852  }
853 
854 
864  void SurfacePoint::ToNaifArray(double naifOutput[3]) const {
865  if(Valid()) {
866  naifOutput[0] = p_x->kilometers();
867  naifOutput[1] = p_y->kilometers();
868  naifOutput[2] = p_z->kilometers();
869  }
870  else {
871  IString msg = "Cannot convert an invalid surface point to a naif array";
872  throw IException(IException::Programmer, msg, _FILEINFO_);
873  }
874  }
875 
885  void SurfacePoint::FromNaifArray(const double naifValues[3]) {
886  if(p_x && p_y && p_z) {
887  p_x->setKilometers(naifValues[0]);
888  p_y->setKilometers(naifValues[1]);
889  p_z->setKilometers(naifValues[2]);
890  }
891  else {
892  p_x = new Displacement(naifValues[0], Displacement::Kilometers);
893  p_y = new Displacement(naifValues[1], Displacement::Kilometers);
894  p_z = new Displacement(naifValues[2], Displacement::Kilometers);
895  }
896 
897  ComputeLocalRadius();
898  p_localRadius = GetLocalRadius();
899  }
900 
901 
908  void SurfacePoint::ResetLocalRadius(const Distance &radius) {
909 
910  if (!radius.isValid()) {
911  IString msg = "Radius value must be a valid Distance.";
912  throw IException(IException::Programmer, msg, _FILEINFO_);
913  }
914 
915  if (!p_x || !p_y || !p_z || !p_x->isValid() || !p_y->isValid() ||
916  !p_z->isValid()) {
917  IString msg = "In order to reset the local radius, a Surface Point must "
918  "already be set.";
919  throw IException(IException::Programmer, msg, _FILEINFO_);
920  }
921 
922  // Set latitudinal coordinates
923  SpiceDouble lat = (double) GetLatitude().radians();
924  SpiceDouble lon = (double) GetLongitude().radians();
925 
926  p_localRadius = radius;
927 
928  // Set rectangular coordinates
929  SpiceDouble rect[3];
930  latrec_c ((SpiceDouble) radius.kilometers(), lon, lat, rect);
931  p_x->setKilometers(rect[0]);
932  p_y->setKilometers(rect[1]);
933  p_z->setKilometers(rect[2]);
934 
935  // TODO What should be done to the variance/covariance matrix when the
936  // radius is reset??? With Bundle updates will this functionality be
937  // obsolete??? Ask Ken
938  }
939 
940 
941  bool SurfacePoint::Valid() const {
942  static const Displacement zero(0, Displacement::Meters);
943  return p_x && p_y && p_z && p_x->isValid() && p_y->isValid() && p_z->isValid() &&
944  (*p_x != zero || *p_y != zero || *p_z != zero);
945  }
946 
947 
956  double SurfacePoint::GetCoord(CoordinateType type, CoordIndex index, CoordUnits units) {
957  // TODO *** Is there a better way to satisfy the compiler that a value will be initialized?
958  // Don't the enums take care of preventing any other possibilities? no
959  double value = 0.;
960 
961  switch (type) {
962 
963  case Latitudinal:
964 
965  switch (index) {
966 
967  case One: // Latitude
968  value = LatToDouble(GetLatitude(), units);
969  break;
970 
971  case Two: // Longitude
972  value = LonToDouble(GetLongitude(), units);
973  break;
974 
975  case Three: // LocalRadius
976  value = DistanceToDouble(GetLocalRadius(), units);
977  break;
978 
979  default:
980  IString msg = "Invalid coordinate index (must be less than 3).";
981  throw IException(IException::Programmer, msg, _FILEINFO_);
982  }
983  break;
984 
985  case Rectangular:
986 
987  switch (index) {
988 
989  case One: // X
990  value = DisplacementToDouble(GetX(), units);
991  break;
992 
993  case Two: // Y
994  value = DisplacementToDouble(GetY(), units);
995  break;
996 
997  case Three: // Z
998  value = DisplacementToDouble(GetZ(), units);
999  break;
1000 
1001  default:
1002  IString msg = "Invalid coordinate index enum (must be less than 3).";
1003  throw IException(IException::Programmer, msg, _FILEINFO_);
1004  }
1005  break;
1006 
1007  default:
1008  IString msg ="Unknown surface point coordinate type enum [" + toString(type) + "]." ;
1009  throw IException(IException::Programmer, msg, _FILEINFO_);
1010  }
1011  return value;
1012  }
1013 
1014 
1023  double SurfacePoint::GetSigma(CoordinateType type, CoordIndex index, CoordUnits units) {
1024  double value = 0; // See first TODO in GetCoord
1025 
1026  switch (type) {
1027 
1028  case Latitudinal:
1029 
1030  switch (index) {
1031 
1032  case One:
1033  value = DistanceToDouble(GetLatSigmaDistance(), units);
1034  break;
1035 
1036  case Two:
1037  value = DistanceToDouble(GetLonSigmaDistance(), units);
1038  break;
1039 
1040  case Three:
1041  value = DistanceToDouble(GetLocalRadiusSigma(), units);
1042  break;
1043 
1044  default:
1045  IString msg = "Invalid coordinate index (must be less than 3).";
1046  throw IException(IException::Programmer, msg, _FILEINFO_);
1047  }
1048  break;
1049 
1050  case Rectangular:
1051 
1052  switch (index) {
1053 
1054  case One:
1055  value = DistanceToDouble(GetXSigma(), units);
1056  break;
1057 
1058  case Two:
1059  value = DistanceToDouble(GetYSigma(), units);
1060  break;
1061 
1062  case Three:
1063  value = DistanceToDouble(GetZSigma(), units);
1064  break;
1065 
1066  default:
1067  IString msg = "Invalid coordinate index (must be less than 3).";
1068  throw IException(IException::Programmer, msg, _FILEINFO_);
1069  }
1070  break;
1071 
1072  default:
1073  IString msg ="Unknown surface point coordinate type enum [" + toString(type) + "]." ;
1074  throw IException(IException::Programmer, msg, _FILEINFO_);
1075  }
1076 
1077  return value;
1078  }
1079 
1080 
1088  Distance SurfacePoint::GetSigmaDistance(CoordinateType type,
1089  CoordIndex index) {
1090  Distance dist = Distance(); // See first TODO in GetCoord
1091 
1092  switch (type) {
1093 
1094  case Latitudinal:
1095 
1096  switch (index) {
1097 
1098  case One:
1099  dist = GetLatSigmaDistance();
1100  break;
1101 
1102  case Two:
1103  dist = GetLonSigmaDistance();
1104  break;
1105 
1106  case Three:
1107  dist = GetLocalRadiusSigma();
1108  break;
1109 
1110  default:
1111  IString msg = "Invalid coordinate index (must be less than 3).";
1112  throw IException(IException::Programmer, msg, _FILEINFO_);
1113  }
1114  break;
1115 
1116  case Rectangular:
1117 
1118  switch (index) {
1119 
1120  case One:
1121  dist = GetXSigma();
1122  break;
1123 
1124  case Two:
1125  dist = GetYSigma();
1126  break;
1127 
1128  case Three:
1129  dist = GetZSigma();
1130  break;
1131 
1132  default:
1133  IString msg = "Invalid coordinate index (must be less than 3).";
1134  throw IException(IException::Programmer, msg, _FILEINFO_);
1135  }
1136  break;
1137 
1138  default:
1139  IString msg ="Unknown surface point coordinate type enum [" + toString(type) + "]." ;
1140  throw IException(IException::Programmer, msg, _FILEINFO_);
1141  }
1142 
1143  return dist;
1144  }
1145 
1146 
1154  double SurfacePoint::DisplacementToDouble(Displacement disp, CoordUnits units) {
1155  double value;
1156 
1157  switch (units) {
1158 
1159  case Meters:
1160  value = disp.meters();
1161  break;
1162 
1163  case Kilometers:
1164  value = disp.kilometers();
1165  break;
1166 
1167  default:
1168  IString msg = "Unrecognized unit for a Displacement (must be meters or kilometers).";
1169  throw IException(IException::Programmer, msg, _FILEINFO_);
1170  }
1171  return value;
1172  }
1173 
1174 
1182  double SurfacePoint::DistanceToDouble(Distance dist, CoordUnits units) {
1183  double value;
1184 
1185  switch (units) {
1186 
1187  case Meters:
1188  value = dist.meters();
1189  break;
1190 
1191  case Kilometers:
1192  value = dist.kilometers();
1193  break;
1194 
1195  default:
1196  IString msg = "Unrecognized unit for a Distance (not meters or kilometers).";
1197  throw IException(IException::Programmer, msg, _FILEINFO_);
1198  }
1199  return value;
1200  }
1201 
1202 
1210  double SurfacePoint::LatToDouble(Latitude lat, CoordUnits units) {
1211  double value;
1212 
1213  switch (units) {
1214 
1215  case Radians:
1216  value = GetLatitude().radians();
1217  break;
1218 
1219  case Degrees:
1220  value = GetLatitude().degrees();
1221  break;
1222 
1223  default:
1224  IString msg = "Invalid unit for latitude (must be Radians or Degrees).";
1225  throw IException(IException::Programmer, msg, _FILEINFO_);
1226  }
1227  return value;
1228  }
1229 
1230 
1242  double SurfacePoint::MetersToLatitude(double latLength) {
1243  if (Valid() && !IsSpecial(latLength)) {
1244  // Convert angle measure in meters to radians relative to latitude of SurfacePoint.
1245  return latLength / GetLocalRadius().meters();
1246  }
1247  else {
1248  // Return Null to be consistent with the bundle classes
1249  return Isis::Null;
1250  }
1251  }
1252 
1253 
1270  double SurfacePoint::MetersToLongitude(double deltaLonMeters) {
1271 
1272  if (Valid() && !IsSpecial(deltaLonMeters)) {
1273  double convFactor = cos((double)GetLatitude().radians());
1274  double deltaLonRadians;
1275 
1276  // Convert angle displacement to radians relative to longitude of SurfacePoint.
1277  if (convFactor > 1.0e-50) {
1278  deltaLonRadians = deltaLonMeters / (convFactor*GetLocalRadius().meters());
1279  }
1280  else {
1281  // Brent Archinal suggested setting sigma to pi in the case of a point near the pole
1282  deltaLonRadians = PI;
1283  }
1284  return deltaLonRadians;
1285  }
1286  else {
1287  // Return Null to be consistent with the bundle classes
1288  return Isis::Null;
1289  }
1290  }
1291 
1292 
1304  double SurfacePoint::LatitudeToMeters(double relativeLat) const {
1305  // Returns are consistent with the bundle classes
1306  if (relativeLat == 0.) {
1307  return 0.;
1308  }
1309  else if (Valid() && !IsSpecial(relativeLat) && GetLocalRadius().isValid() ) {
1310  return relativeLat * GetLocalRadius().meters();
1311  }
1312  else {
1313  return Isis::Null;
1314  }
1315  }
1316 
1317 
1330  double SurfacePoint::LongitudeToMeters(double deltaLonRadians) const{
1331  // Returns are consistent with the bundle classes
1332  double deltaLonMeters = Isis::Null;
1333 
1334  if (deltaLonRadians == 0.) {
1335  deltaLonMeters = 0.;
1336  }
1337  else if (Valid() && !IsSpecial(deltaLonRadians) && GetLocalRadius().isValid() ) {
1338  // Convert from radians to meters and return
1339  double scalingRadius = cos(GetLatitude().radians()) * GetLocalRadius().meters();
1340  deltaLonMeters = deltaLonRadians * scalingRadius;
1341  }
1342 
1343  return deltaLonMeters;
1344  }
1345 
1346 
1357  SurfacePoint::stringToCoordinateType(QString type) {
1358  if (type.compare("LATITUDINAL", Qt::CaseInsensitive) == 0) {
1359  return SurfacePoint::Latitudinal;
1360  }
1361  else if (type.compare("RECTANGULAR", Qt::CaseInsensitive) == 0) {
1362  return SurfacePoint::Rectangular;
1363  }
1364  else {
1365  throw IException(IException::Programmer,
1366  "Unknown coordinate type for a SurfacePoint [" + type + "].",
1367  _FILEINFO_);
1368  }
1369  }
1370 
1371 
1383  QString SurfacePoint::coordinateTypeToString(
1385  switch (type) {
1386 
1387  case Latitudinal:
1388  return "Latitudinal";
1389  break;
1390 
1391  case Rectangular:
1392  return "Rectangular";
1393  break;
1394 
1395  default:
1396  IString msg ="Unknown surface point coordinate type enum [" + toString(type) + "]." ;
1397  throw IException(IException::Programmer, msg, _FILEINFO_);
1398  }
1399  }
1400 
1401 
1409  double SurfacePoint::LonToDouble(Longitude lon, CoordUnits units) {
1410  double value;
1411 
1412  switch (units) {
1413 
1414  case Radians:
1415  value = GetLongitude().radians();
1416  break;
1417 
1418  case Degrees:
1419  value = GetLongitude().degrees();
1420  break;
1421 
1422  default:
1423  IString msg = "Invalid unit for longitude (must be Radians of Degrees).";
1424  throw IException(IException::Programmer, msg, _FILEINFO_);
1425  }
1426  return value;
1427  }
1428 
1429 
1430  Displacement SurfacePoint::GetX() const {
1431  if(!p_x) return Displacement();
1432 
1433  return *p_x;
1434  }
1435 
1436 
1437  Displacement SurfacePoint::GetY() const {
1438  if(!p_y) return Displacement();
1439 
1440  return *p_y;
1441  }
1442 
1443 
1444  Displacement SurfacePoint::GetZ() const {
1445  if(!p_z) return Displacement();
1446 
1447  return *p_z;
1448  }
1449 
1450 
1451  Distance SurfacePoint::GetXSigma() const {
1452  if(!p_rectCovar) return Distance();
1453 
1454  return Distance(sqrt((*p_rectCovar)(0, 0)), Distance::Kilometers);
1455  }
1456 
1457 
1458  Distance SurfacePoint::GetYSigma() const {
1459  if(!p_rectCovar) return Distance();
1460 
1461  return Distance(sqrt((*p_rectCovar)(1, 1)), Distance::Kilometers);
1462  }
1463 
1464 
1465  Distance SurfacePoint::GetZSigma() const {
1466  if(!p_rectCovar) return Distance();
1467 
1468  return Distance(sqrt((*p_rectCovar)(2, 2)), Distance::Kilometers);
1469  }
1470 
1471 
1485  double SurfacePoint::GetWeight(CoordinateType type, CoordIndex index) {
1486  double value = 0; // See first TODO in GetCoord
1487 
1488  switch (type) {
1489 
1490  case Latitudinal:
1491 
1492  switch (index) {
1493 
1494  case One:
1495  value = GetLatWeight(); // 1/radians**2
1496  break;
1497 
1498  case Two:
1499  value = GetLonWeight(); // 1/radians**2
1500  break;
1501 
1502  case Three:
1503  value = GetLocalRadiusWeight(); // 1/km**2
1504  break;
1505 
1506  default:
1507  IString msg = "Invalid coordinate index (must be less than 3).";
1508  throw IException(IException::Programmer, msg, _FILEINFO_);
1509  }
1510  break;
1511 
1512  case Rectangular:
1513 
1514  switch (index) {
1515 
1516  case One:
1517  value = GetXWeight(); // 1/km**2
1518  break;
1519 
1520  case Two:
1521  value = GetYWeight(); // 1/km**2
1522  break;
1523 
1524  case Three:
1525  value = GetZWeight(); // 1/km**2
1526  break;
1527 
1528  default:
1529  IString msg = "Invalid coordinate index (must be less than 3).";
1530  throw IException(IException::Programmer, msg, _FILEINFO_);
1531  }
1532  break;
1533 
1534  default:
1535  IString msg ="Unknown surface point coordinate type enum [" + toString(type) + "]." ;
1536  throw IException(IException::Programmer, msg, _FILEINFO_);
1537  }
1538  return value;
1539  }
1540 
1541 
1547  double SurfacePoint::GetXWeight() const {
1548 
1549  double dXSigma = GetXSigma().kilometers();
1550 
1551  if (dXSigma <= 0.0 ) {
1552  IString msg = "SurfacePoint::GetXWeight(): Sigma <= 0.0";
1553  throw IException(IException::Programmer, msg, _FILEINFO_);
1554  }
1555 
1556  return 1.0/(dXSigma*dXSigma);
1557  }
1558 
1559 
1565  double SurfacePoint::GetYWeight() const {
1566 
1567  double dYSigma = GetYSigma().kilometers();
1568 
1569  if (dYSigma <= 0.0 ) {
1570  IString msg = "SurfacePoint::GetYWeight(): Sigma <= 0.0";
1571  throw IException(IException::Programmer, msg, _FILEINFO_);
1572  }
1573 
1574  return 1.0/(dYSigma*dYSigma);
1575  }
1576 
1577 
1583  double SurfacePoint::GetZWeight() const {
1584 
1585  double dZSigma = GetZSigma().kilometers();
1586 
1587  if (dZSigma <= 0.0 ) {
1588  IString msg = "SurfacePoint::GetZWeight(): Sigma <= 0.0";
1589  throw IException(IException::Programmer, msg, _FILEINFO_);
1590  }
1591 
1592  return 1.0/(dZSigma*dZSigma);
1593  }
1594 
1595 
1596  symmetric_matrix<double, upper> SurfacePoint::GetRectangularMatrix
1597  (SurfacePoint::CoordUnits units) const {
1598  if(!p_rectCovar) {
1599  symmetric_matrix<double, upper> tmp(3);
1600  tmp.clear();
1601  return tmp;
1602  }
1603 
1604  // Make sure the units are valid
1605  if (units != Kilometers && units != Meters) {
1606  IString msg = "Units must be Kilometers or Meters";
1607  throw IException(IException::Programmer, msg, _FILEINFO_);
1608  }
1609 
1610  symmetric_matrix<double, upper> covar(3);
1611  covar.clear();
1612 
1613  switch (units) {
1614 
1615  case Meters:
1616  // Convert member matrix to Meters to return
1617  covar(0,0) = (*p_rectCovar)(0,0)*1.0e6;
1618  covar(0,1) = (*p_rectCovar)(0,1)*1.0e6;
1619  covar(0,2) = (*p_rectCovar)(0,2)*1.0e6;
1620  covar(1,1) = (*p_rectCovar)(1,1)*1.0e6;
1621  covar(1,2) = (*p_rectCovar)(1,2)*1.0e6;
1622  covar(2,2) = (*p_rectCovar)(2,2)*1.0e6;
1623  return covar;
1624  break;
1625 
1626  case Kilometers:
1627  return *p_rectCovar;
1628  break;
1629 
1630  default:
1631  IString msg = "Unrecognized unit for a Displacement (must be meters or kilometers).";
1632  throw IException(IException::Programmer, msg, _FILEINFO_);
1633  }
1634 
1635  return *p_rectCovar;
1636  }
1637 
1638 
1639  Angle SurfacePoint::GetLatSigma() const {
1640  if(!p_sphereCovar)
1641  return Angle();
1642 
1643  return Angle(sqrt((*p_sphereCovar)(0, 0)), Angle::Radians);
1644  }
1645 
1646 
1647  Angle SurfacePoint::GetLonSigma() const {
1648  if(!p_sphereCovar)
1649  return Angle();
1650 
1651  return Angle(sqrt((*p_sphereCovar)(1, 1)), Angle::Radians);
1652  }
1653 
1654 
1659  Latitude SurfacePoint::GetLatitude() const {
1660  if (!Valid())
1661  return Latitude();
1662 
1663  // TODO Scale for accuracy with coordinate of largest magnitude
1664  double x = p_x->meters();
1665  double y = p_y->meters();
1666  double z = p_z->meters();
1667 
1668  if (x != 0. || y != 0. || z != 0.)
1669  return Latitude(atan2(z, sqrt(x*x + y*y) ), Angle::Radians);
1670  else
1671  return Latitude();
1672  }
1673 
1674 
1679  Longitude SurfacePoint::GetLongitude() const {
1680  if (!Valid())
1681  return Longitude();
1682 
1683  double x = p_x->meters();
1684  double y = p_y->meters();
1685 
1686  if(x == 0.0 && y == 0.0) {
1687  return Longitude(0, Angle::Radians);
1688  }
1689 
1690  double lon = atan2(y, x);
1691  if (lon < 0) {
1692  lon += 2 * PI;
1693  }
1694 
1695  return Longitude(lon, Angle::Radians);
1696  }
1697 
1698 
1703  void SurfacePoint::ComputeLocalRadius() {
1704  static const Displacement zero(0, Displacement::Meters);
1705  if (Valid()) {
1706  double x = p_x->meters();
1707  double y = p_y->meters();
1708  double z = p_z->meters();
1709 
1710  p_localRadius = Distance(sqrt(x*x + y*y + z*z), Distance::Meters);
1711  }
1712  else if (*p_x == zero && *p_y == zero && *p_z == zero) { // for backwards compatability
1713  p_localRadius = Distance(0., Distance::Meters);
1714  }
1715  else { // Invalid point
1716  IString msg = "SurfacePoint::Can't compute local radius on invalid point";
1717  throw IException(IException::Programmer, msg, _FILEINFO_);
1718  }
1719  }
1720 
1721 
1726  Distance SurfacePoint::GetLocalRadius() const {
1727  return p_localRadius;
1728  }
1729 
1730 
1741  Distance SurfacePoint::GetLatSigmaDistance() const {
1742  double d = LatitudeToMeters(GetLatSigma().radians());
1743 
1744  if (d > 1.0e-50) {
1745  return Distance(d, Distance::Meters);
1746  }
1747  else {
1748  return Distance();
1749  }
1750  }
1751 
1752 
1757  Distance SurfacePoint::GetLonSigmaDistance() const{
1758  // return lonSigmaDistance;
1759  double d = LongitudeToMeters(GetLonSigma().radians());
1760 // TODO What do we do when the scaling radius is 0 (at the pole)?
1761  // if (d > DBL_EPSILON) {
1762  return Distance(d, Distance::Meters);
1763  // }
1764  // else { // Too close to the pole
1765  // return Distance();
1766  // }
1767  }
1768 
1769 
1770  Distance SurfacePoint::GetLocalRadiusSigma() const {
1771  if(!p_sphereCovar)
1772  return Distance();
1773 
1774  return Distance(sqrt((*p_sphereCovar)(2, 2)), Distance::Kilometers);
1775  }
1776 
1777 
1778  symmetric_matrix<double, upper> SurfacePoint::GetSphericalMatrix
1779  (SurfacePoint::CoordUnits units) const {
1780  if(!p_sphereCovar) {
1781  symmetric_matrix<double, upper> tmp(3);
1782  tmp.clear();
1783  return tmp;
1784  }
1785 
1786  // Make sure the units are valid
1787  if (units != Kilometers && units != Meters) {
1788  IString msg = "Units must be Kilometers or Meters";
1789  throw IException(IException::Programmer, msg, _FILEINFO_);
1790  }
1791 
1792  symmetric_matrix<double, upper> covar(3);
1793  covar.clear();
1794 
1795  switch (units) {
1796 
1797  case Meters:
1798  // Convert member matrix to Meters to return
1799  covar(0,0) = (*p_sphereCovar)(0,0);
1800  covar(0,1) = (*p_sphereCovar)(0,1);
1801  covar(0,2) = (*p_sphereCovar)(0,2)*1000.;
1802  covar(1,1) = (*p_sphereCovar)(1,1);
1803  covar(1,2) = (*p_sphereCovar)(1,2)*1000.;
1804  covar(2,2) = (*p_sphereCovar)(2,2)*1.0e6;
1805  return covar;
1806  break;
1807 
1808  case Kilometers:
1809  return *p_sphereCovar;
1810  break;
1811 
1812  default:
1813  IString msg = "Unrecognized unit for a Displacement (must be meters or kilometers).";
1814  throw IException(IException::Programmer, msg, _FILEINFO_);
1815  }
1816 
1817  return *p_sphereCovar;
1818  }
1819 
1820 
1826  double SurfacePoint::GetLatWeight() const {
1827  double dlatSigma = GetLatSigma().radians();
1828 
1829  if( dlatSigma <= 0.0 ) {
1830  IString msg = "SurfacePoint::GetLatWeight(): Sigma <= 0.0";
1831  throw IException(IException::Programmer, msg, _FILEINFO_);
1832  }
1833 
1834  return 1.0/(dlatSigma*dlatSigma);
1835  }
1836 
1842  double SurfacePoint::GetLonWeight() const {
1843  double dlonSigma = GetLonSigma().radians();
1844 
1845  if( dlonSigma <= 0.0 ) {
1846  IString msg = "SurfacePoint::GetLonWeight(): Sigma <= 0.0";
1847  throw IException(IException::Programmer, msg, _FILEINFO_);
1848  }
1849 
1850  return 1.0/(dlonSigma*dlonSigma);
1851  }
1852 
1858  double SurfacePoint::GetLocalRadiusWeight() const {
1859 
1860  double dlocalRadiusSigma = GetLocalRadiusSigma().kilometers();
1861 
1862  if (dlocalRadiusSigma <= 0.0 ) {
1863  IString msg = "SurfacePoint::GetRadWeight(): Sigma <= 0.0";
1864  throw IException(IException::Programmer, msg, _FILEINFO_);
1865  }
1866 
1867  return 1.0/(dlocalRadiusSigma*dlocalRadiusSigma);
1868  }
1869 
1874  Distance SurfacePoint::GetDistanceToPoint(const SurfacePoint &other) const {
1875 
1876  if(!Valid() || !other.Valid())
1877  return Distance();
1878 
1879  return GetDistanceToPoint(other,
1880  ((GetLocalRadius() + other.GetLocalRadius()) / 2.0));
1881  }
1882 
1883 
1891  Distance SurfacePoint::GetDistanceToPoint(const SurfacePoint &other,
1892  const Distance &sphereRadius) const {
1893  if(!Valid() || !other.Valid())
1894  return Distance();
1895 
1896  // Convert lat/lon values to radians
1897  const Angle &latitude = GetLatitude();
1898  const Angle &longitude = GetLongitude();
1899  const Angle &otherLatitude = other.GetLatitude();
1900  const Angle &otherLongitude = other.GetLongitude();
1901 
1902  // The harvestine method:
1903  // http://en.wikipedia.org/wiki/Haversine_formula
1904  Angle deltaLat = latitude - otherLatitude;
1905  Angle deltaLon = longitude - otherLongitude;
1906 
1907  double haversinLat = sin(deltaLat.radians() / 2.0);
1908  haversinLat *= haversinLat;
1909 
1910  double haversinLon = sin(deltaLon.radians() / 2.0);
1911  haversinLon *= haversinLon;
1912 
1913  double a = haversinLat + cos(latitude.radians()) *
1914  cos(otherLatitude.radians()) *
1915  haversinLon;
1916 
1917  double c = 2 * atan(sqrt(a) / sqrt(1 - a));
1918 
1919  return sphereRadius * c;
1920  }
1921 
1922 
1923  bool SurfacePoint::operator==(const SurfacePoint &other) const {
1924  bool equal = true;
1925 
1926  if(equal && p_x && p_y && p_z && other.p_x && other.p_y && other.p_z) {
1927  equal = equal && *p_x == *other.p_x;
1928  equal = equal && *p_y == *other.p_y;
1929  equal = equal && *p_z == *other.p_z;
1930  }
1931  else if(equal) {
1932  equal = equal && p_x == NULL && other.p_x == NULL;
1933  equal = equal && p_y == NULL && other.p_y == NULL;
1934  equal = equal && p_z == NULL && other.p_z == NULL;
1935  }
1936 
1937  if(equal && p_rectCovar) {
1938  equal = equal && (*p_rectCovar)(0, 0) == (*other.p_rectCovar)(0, 0);
1939  equal = equal && (*p_rectCovar)(0, 1) == (*other.p_rectCovar)(0, 1);
1940  equal = equal && (*p_rectCovar)(0, 2) == (*other.p_rectCovar)(0, 2);
1941  equal = equal && (*p_rectCovar)(1, 1) == (*other.p_rectCovar)(1, 1);
1942  equal = equal && (*p_rectCovar)(1, 2) == (*other.p_rectCovar)(1, 2);
1943  equal = equal && (*p_rectCovar)(2, 2) == (*other.p_rectCovar)(2, 2);
1944  }
1945  else if(equal) {
1946  equal = equal && p_rectCovar == NULL && other.p_rectCovar == NULL;
1947  }
1948 
1949  if(equal && p_sphereCovar) {
1950  equal = equal && (*p_sphereCovar)(0, 0) == (*other.p_sphereCovar)(0, 0);
1951  equal = equal && (*p_sphereCovar)(0, 1) == (*other.p_sphereCovar)(0, 1);
1952  equal = equal && (*p_sphereCovar)(0, 2) == (*other.p_sphereCovar)(0, 2);
1953  equal = equal && (*p_sphereCovar)(1, 1) == (*other.p_sphereCovar)(1, 1);
1954  equal = equal && (*p_sphereCovar)(1, 2) == (*other.p_sphereCovar)(1, 2);
1955  equal = equal && (*p_sphereCovar)(2, 2) == (*other.p_sphereCovar)(2, 2);
1956  }
1957  else if(equal) {
1958  equal = equal && p_sphereCovar == NULL && other.p_sphereCovar == NULL;
1959  }
1960 
1961  return equal;
1962  }
1963 
1964  SurfacePoint &SurfacePoint::operator=(const SurfacePoint &other) {
1965  // The lazy way of doing this (free all memory and copy) is too expensive
1966  // in the default case!
1967  if(p_x && other.p_x &&
1968  p_y && other.p_y &&
1969  p_z && other.p_z &&
1970  !p_rectCovar && !other.p_rectCovar &&
1971  !p_sphereCovar && !other.p_sphereCovar) {
1972  *p_x = *other.p_x;
1973  *p_y = *other.p_y;
1974  *p_z = *other.p_z;
1975  }
1976  else {
1977  FreeAllocatedMemory();
1978 
1979  if(other.p_x) {
1980  p_x = new Displacement(*other.p_x);
1981  }
1982 
1983  if(other.p_y) {
1984  p_y = new Displacement(*other.p_y);
1985  }
1986 
1987  if(other.p_z) {
1988  p_z = new Displacement(*other.p_z);
1989  }
1990 
1991  if(other.p_rectCovar) {
1992  p_rectCovar = new symmetric_matrix<double, upper>(*other.p_rectCovar);
1993  }
1994 
1995  if(other.p_sphereCovar) {
1996  p_sphereCovar = new symmetric_matrix<double, upper>(*other.p_sphereCovar);
1997  }
1998  }
1999 
2000  // Finally initialize local radius to avoid using a previous value
2001  p_localRadius = other.GetLocalRadius();
2002 
2003  return *this;
2004  }
2005 
2006  void SurfacePoint::FreeAllocatedMemory() {
2007  if(p_x) {
2008  delete p_x;
2009  p_x = NULL;
2010  }
2011 
2012  if(p_y) {
2013  delete p_y;
2014  p_y = NULL;
2015  }
2016 
2017  if(p_z) {
2018  delete p_z;
2019  p_z = NULL;
2020  }
2021 
2022  if(p_rectCovar) {
2023  delete p_rectCovar;
2024  p_rectCovar = NULL;
2025  }
2026 
2027  if(p_sphereCovar) {
2028  delete p_sphereCovar;
2029  p_sphereCovar = NULL;
2030  }
2031  }
2032 }
This class defines a body-fixed surface point.
Definition: SurfacePoint.h:148
double meters() const
Get the distance in meters.
Definition: Distance.cpp:97
boost::numeric::ublas::symmetric_matrix< double, boost::numeric::ublas::upper > * p_sphereCovar
3x3 upper triangular covariance matrix ocentric coordinates
Definition: SurfacePoint.h:331
const double Null
Value for an Isis Null pixel.
Definition: SpecialPixel.h:110
boost::numeric::ublas::symmetric_matrix< double, boost::numeric::ublas::upper > * p_rectCovar
3x3 upper triangular covariance matrix rectangular coordinates
Definition: SurfacePoint.h:328
const double PI
The mathematical constant PI.
Definition: Constants.h:56
double radians() const
Convert an angle to a double.
Definition: Angle.h:243
double meters() const
Get the displacement in meters.
bool isValid() const
Test if this displacement has been initialized or not.
Namespace for the standard library.
This class is designed to encapsulate the concept of a Latitude.
Definition: Latitude.h:63
QString toString(bool boolToConvert)
Global function to convert a boolean to a string.
Definition: IString.cpp:226
double kilometers() const
Get the distance in kilometers.
Definition: Distance.cpp:118
Distance measurement, usually in meters.
Definition: Distance.h:47
Latitude GetLatitude() const
Return the body-fixed latitude for the surface point.
This class is designed to encapsulate the concept of a Longitude.
Definition: Longitude.h:52
#define _FILEINFO_
Macro for the filename and line number.
Definition: IException.h:40
void setKilometers(double distanceInKilometers)
Set the distance in kilometers.
Definition: Distance.cpp:130
bool IsSpecial(const double d)
Returns if the input pixel is special.
Definition: SpecialPixel.h:212
bool isValid() const
Test if this distance has been initialized or not.
Definition: Distance.cpp:204
CoordinateType
Defines the coordinate typ, units, and coordinate index for some of the output methods.
Definition: SurfacePoint.h:155
Defines an angle and provides unit conversions.
Definition: Angle.h:62
Unless noted otherwise, the portions of Isis written by the USGS are public domain.
Longitude GetLongitude() const
Return the body-fixed longitude for the surface point.
Displacement is a signed length, usually in meters.
Definition: Displacement.h:43
Isis exception class.
Definition: IException.h:107
Adds specific functionality to C++ strings.
Definition: IString.h:181
Namespace for ISIS/Bullet specific routines.
Definition: Apollo.h:31
bool isValid() const
This indicates whether we have a legitimate angle stored or are in an unset, or invalid, state.
Definition: Angle.cpp:110
double kilometers() const
Get the displacement in kilometers.
Distance GetLocalRadius() const
Return the radius of the surface point.
Unless noted otherwise, the portions of Isis written by the USGS are public domain.