Position.cpp

Go to the documentation of this file.
00001 #pragma ident "$Id: Position.cpp 669 2007-07-05 16:52:08Z pben $"
00002 
00003 
00004 
00014 //============================================================================
00015 //
00016 //  This file is part of GPSTk, the GPS Toolkit.
00017 //
00018 //  The GPSTk is free software; you can redistribute it and/or modify
00019 //  it under the terms of the GNU Lesser General Public License as published
00020 //  by the Free Software Foundation; either version 2.1 of the License, or
00021 //  any later version.
00022 //
00023 //  The GPSTk is distributed in the hope that it will be useful,
00024 //  but WITHOUT ANY WARRANTY; without even the implied warranty of
00025 //  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00026 //  GNU Lesser General Public License for more details.
00027 //
00028 //  You should have received a copy of the GNU Lesser General Public
00029 //  License along with GPSTk; if not, write to the Free Software Foundation,
00030 //  Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00031 //  
00032 //  Copyright 2004, The University of Texas at Austin
00033 //
00034 //============================================================================
00035 
00036 #include "Position.hpp"
00037 #include "WGS84Geoid.hpp"
00038 #include "icd_200_constants.hpp"    // for TWO_PI, etc
00039 #include "geometry.hpp"             // for RAD_TO_DEG, etc
00040 #include "MiscMath.hpp"             // for RSS, SQRT
00041 
00042 namespace gpstk
00043 {
00044 
00045    using namespace std;
00046    using namespace StringUtils;
00047 
00048    // ----------- Part  1: coordinate systems --------------------------------
00049       // Labels for coordinate systems supported by Position
00050    static const char *SystemNames[] = {
00051       "Unknown",
00052       "Geodetic",
00053       "Geocentric",
00054       "Cartesian",
00055       "Spherical"};
00056 
00057       // return string giving name of coordinate system
00058    string Position::getSystemName()
00059       throw()
00060    { return SystemNames[system]; }
00061   
00062    // ----------- Part  2: tolerance -----------------------------------------
00063       // One millimeter tolerance.
00064    const double Position::ONE_MM_TOLERANCE = 0.001;
00065       // One centimeter tolerance.
00066    const double Position::ONE_CM_TOLERANCE = 0.01;
00067       // One micron tolerance.
00068    const double Position::ONE_UM_TOLERANCE = 0.000001;
00069    
00070       // Default tolerance for time equality in meters.
00071    double Position::POSITION_TOLERANCE = Position::ONE_MM_TOLERANCE;
00072 
00073       // Sets the tolerance for output and comparisons, for this object only.
00074       // See the constants in this file (e.g. ONE_MM_TOLERANCE)
00075       // for some easy to use tolerance values.
00076       // @param tol Tolerance in meters to be used by comparison operators.
00077    Position& Position::setTolerance(const double tol)
00078       throw()
00079    {
00080       tolerance = tol;
00081       return *this;
00082    }
00083 
00084    // ----------- Part  3: member functions: constructors --------------------
00085    //
00086       // Default constructor.
00087    Position::Position()
00088       throw()
00089    {
00090       WGS84Geoid WGS84;
00091       initialize(0.0,0.0,0.0,Unknown,&WGS84);
00092    }
00093 
00094    Position::Position(const double& a,
00095                       const double& b,
00096                       const double& c,
00097                       Position::CoordinateSystem s,
00098                       GeoidModel *geoid)
00099       throw(GeometryException)
00100    {
00101       try {
00102          initialize(a,b,c,s,geoid);
00103       }
00104       catch(GeometryException& ge) {
00105          GPSTK_RETHROW(ge);
00106       }
00107    }
00108 
00109    Position::Position(const double ABC[3],
00110                       CoordinateSystem s,
00111                       GeoidModel *geoid)
00112       throw(GeometryException)
00113    {
00114       double a=ABC[0];
00115       double b=ABC[1];
00116       double c=ABC[2];
00117       try {
00118          initialize(a,b,c,s,geoid);
00119       }
00120       catch(GeometryException& ge) {
00121          GPSTK_RETHROW(ge);
00122       }
00123    }
00124 
00125    Position::Position(const Triple& ABC,
00126                       CoordinateSystem s,
00127                       GeoidModel *geoid)
00128       throw(GeometryException)
00129    {
00130       double a=ABC[0];
00131       double b=ABC[1];
00132       double c=ABC[2];
00133       try {
00134          initialize(a,b,c,s,geoid);
00135       }
00136       catch(GeometryException& ge) {
00137          GPSTK_RETHROW(ge);
00138       }
00139    }
00140 
00141    Position::Position(const Xvt& xvt)
00142       throw()
00143    {
00144       double a=xvt.x[0];
00145       double b=xvt.x[1];
00146       double c=xvt.x[2];
00147       initialize(a,b,c,Cartesian);
00148    }
00149 
00150    // ----------- Part  4: member functions: arithmetic ----------------------
00151    //
00152    Position& Position::operator-=(const Position& right)
00153       throw()
00154    {
00155       Position r(right);
00156       CoordinateSystem savesys=system;    // save the original system
00157 
00158          // convert to cartestian and difference there
00159       transformTo(Cartesian);
00160       r.transformTo(Cartesian);
00161 
00162       for(int i=0; i<3; i++)
00163          theArray[i] -= r.theArray[i];
00164 
00165       transformTo(savesys);               // transform back to the original system
00166       return *this;
00167    }
00168 
00169    Position& Position::operator+=(const Position& right)
00170       throw()
00171    {
00172       Position r(right);
00173       CoordinateSystem savesys=system;    // save the original system
00174 
00175          // convert to cartestian and difference there
00176       transformTo(Cartesian);
00177       r.transformTo(Cartesian);
00178 
00179       for(int i=0; i<3; i++)
00180          theArray[i] += r.theArray[i];
00181 
00182       transformTo(savesys);               // transform back to the original system
00183       return *this;
00184    }
00185 
00186    Position operator-(const Position& left,
00187                             const Position& right)
00188       throw()
00189    {
00190       Position l(left),r(right);
00191          // convert both to Cartesian
00192       l.transformTo(Position::Cartesian);
00193       r.transformTo(Position::Cartesian);
00194          // difference
00195       l -= r;
00196 
00197       return l;
00198    }
00199 
00200    Position operator+(const Position& left,
00201                             const Position& right)
00202       throw()
00203    {
00204       Position l(left),r(right);
00205          // convert both to Cartesian
00206       l.transformTo(Position::Cartesian);
00207       r.transformTo(Position::Cartesian);
00208          // add
00209       l += r;
00210 
00211       return l;
00212    }
00213 
00214    // ----------- Part  5: member functions: comparisons ---------------------
00215    //
00216       // Equality operator. Returns false if geoid values differ.
00217    bool Position::operator==(const Position &right) const
00218       throw()
00219    {
00220       if(AEarth != right.AEarth || eccSquared != right.eccSquared)
00221          return false;
00222       if(range(*this,right) < tolerance)
00223          return true;
00224       else
00225          return false;
00226    }
00227 
00228       // Inequality operator. Returns true if geoid values differ.
00229    bool Position::operator!=(const Position &right) const
00230       throw()
00231    {
00232       return !(operator==(right));
00233    }
00234 
00235    // ----------- Part  6: member functions: coordinate transformations ------
00236    //
00237       // Transform coordinate system. Does nothing if sys already matches the
00238       // current value of member CoordinateSystem 'system'.
00239       // @param sys coordinate system into which *this is to be transformed.
00240       // @return *this
00241    Position Position::transformTo(CoordinateSystem sys)
00242       throw()
00243    {
00244       if(sys == Unknown || sys == system) return *this;
00245 
00246       // this copies geoid information and tolerance
00247       Position target(*this);
00248 
00249       // transform target.theArray and set target.system
00250       switch(system) {
00251          case Unknown:
00252             return *this;
00253          case Geodetic:
00254             // --------------- Geodetic to ... ------------------------
00255             switch(sys) {
00256                case Unknown: case Geodetic: return *this;
00257                case Geocentric:
00258                   convertGeodeticToGeocentric(*this,target,AEarth,eccSquared);
00259                   target.system = Geocentric;
00260                   break;
00261                case Cartesian:
00262                   convertGeodeticToCartesian(*this,target,AEarth,eccSquared);
00263                   target.system = Cartesian;
00264                   break;
00265                case Spherical:
00266                   convertGeodeticToGeocentric(*this,target,AEarth,eccSquared);
00267                   target.theArray[0] = 90 - target.theArray[0];   // geocen -> sph
00268                   target.system = Spherical;
00269                   break;
00270             }
00271             break;
00272          case Geocentric:
00273             // --------------- Geocentric to ... ----------------------
00274             switch(sys) {
00275                case Unknown: case Geocentric: return *this;
00276                case Geodetic:
00277                   convertGeocentricToGeodetic(*this,target,AEarth,eccSquared);
00278                   target.system = Geodetic;
00279                   break;
00280                case Cartesian:
00281                   convertGeocentricToCartesian(*this,target);
00282                   target.system = Cartesian;
00283                   break;
00284                case Spherical:
00285                   target.theArray[0] = 90 - target.theArray[0];   // geocen -> sph
00286                   target.system = Spherical;
00287                   break;
00288             }
00289             break;
00290          case Cartesian:
00291             // --------------- Cartesian to ... -----------------------
00292             switch(sys) {
00293                case Unknown: case Cartesian: return *this;
00294                case Geodetic:
00295                   convertCartesianToGeodetic(*this,target,AEarth,eccSquared);
00296                   target.system = Geodetic;
00297                   break;
00298                case Geocentric:
00299                   convertCartesianToGeocentric(*this,target);
00300                   target.system = Geocentric;
00301                   break;
00302                case Spherical:
00303                   convertCartesianToSpherical(*this,target);
00304                   target.system = Spherical;
00305                   break;
00306             }
00307             break;
00308          case Spherical:
00309             // --------------- Spherical to ... -----------------------
00310             switch(sys) {
00311                case Unknown: case Spherical: return *this;
00312                case Geodetic:
00313                   theArray[0] = 90 - theArray[0];   // sph -> geocen
00314                   convertGeocentricToGeodetic(*this,target,AEarth,eccSquared);
00315                   target.system = Geodetic;
00316                   break;
00317                case Geocentric:
00318                   target.theArray[0] = 90 - target.theArray[0];   // sph -> geocen
00319                   target.system = Geocentric;
00320                   break;
00321                case Cartesian:
00322                   convertSphericalToCartesian(*this,target);
00323                   target.system = Cartesian;
00324                   break;
00325             }
00326             break;
00327       }  // end switch(system)
00328 
00329       *this = target;
00330       return *this;
00331    }
00332   
00333    // ----------- Part  7: member functions: get -----------------------------
00334    // 
00335    // These routines retrieve coordinate values in all coordinate systems.
00336    // Note that calling these will transform the Position to another coordinate
00337    // system if that is required.
00338    //
00339       // Get X coordinate (meters)
00340    double Position::X() const
00341       throw()
00342    {
00343       if(system == Cartesian)
00344          return theArray[0];
00345       Position t(*this);
00346       t.transformTo(Cartesian);
00347       return t.theArray[0];
00348    }
00349 
00350       // Get Y coordinate (meters)
00351    double Position::Y() const
00352       throw()
00353    {
00354       if(system == Cartesian)
00355          return theArray[1];
00356       Position t(*this);
00357       t.transformTo(Cartesian);
00358       return t.theArray[1];
00359    }
00360 
00361       // Get Z coordinate (meters)
00362    double Position::Z() const
00363       throw()
00364    {
00365       if(system == Cartesian)
00366          return theArray[2];
00367       Position t(*this);
00368       t.transformTo(Cartesian);
00369       return t.theArray[2];
00370    }
00371 
00372       // Get geodetic latitude (degrees North).
00373    double Position::geodeticLatitude() const
00374       throw()
00375    {
00376       if(system == Geodetic)
00377          return theArray[0];
00378       Position t(*this);
00379       t.transformTo(Geodetic);
00380       return t.theArray[0];
00381    }
00382 
00383       // Get geocentric latitude (degrees North),
00384       // equal to 90 degress - theta in regular spherical coordinates.
00385    double Position::geocentricLatitude() const
00386       throw()
00387    {
00388       if(system == Geocentric)
00389          return theArray[0];
00390       Position t(*this);
00391       t.transformTo(Geocentric);
00392       return t.theArray[0];
00393    }
00394 
00395       // Get spherical coordinate theta in degrees
00396    double Position::theta() const
00397       throw()
00398    {
00399       if(system == Spherical)
00400          return theArray[0];
00401       Position t(*this);
00402       t.transformTo(Spherical);
00403       return t.theArray[0];
00404    }
00405 
00406       // Get spherical coordinate phi in degrees
00407    double Position::phi() const
00408       throw()
00409    {
00410       if(system == Spherical)
00411          return theArray[1];
00412       Position t(*this);
00413       t.transformTo(Spherical);
00414       return t.theArray[1];
00415    }
00416 
00417       // Get longitude (degrees East),
00418       // equal to phi in regular spherical coordinates.
00419    double Position::longitude() const
00420       throw()
00421    {
00422       if(system != Cartesian)
00423          return theArray[1];
00424       Position t(*this);
00425       t.transformTo(Spherical);
00426       return t.theArray[1];
00427    }
00428 
00429       // Get radius or distance from the center of Earth (meters),
00430       // Same as radius in spherical coordinates.
00431    double Position::radius() const
00432       throw()
00433    {
00434       if(system == Spherical || system == Geocentric)
00435          return theArray[2];
00436       Position t(*this);
00437       t.transformTo(Spherical);
00438       return t.theArray[2];
00439    }
00440 
00441       // Get height above ellipsoid (meters) (Geodetic).
00442    double Position::height() const
00443       throw()
00444    {
00445       if(system == Geodetic)
00446          return theArray[2];
00447       Position t(*this);
00448       t.transformTo(Geodetic);
00449       return t.theArray[2];
00450    }
00451 
00452    // ----------- Part  8: member functions: set -----------------------------
00453    //
00454       // Set the geoid values for this Position given a geoid.
00455       // @param geoid Pointer to the GeoidModel.
00456       // @throw GeometryException if input is NULL.
00457    void Position::setGeoidModel(const GeoidModel *geoid)
00458       throw(GeometryException)
00459    {
00460       if(!geoid)
00461       {
00462          GeometryException ge("Given GeoidModel pointer is NULL.");
00463          GPSTK_THROW(ge);
00464       }
00465       AEarth = geoid->a();
00466       eccSquared = geoid->eccSquared();
00467    }
00468 
00469       // Set the Position given geodetic coordinates, system is set to Geodetic.
00470       // @param lat geodetic latitude in degrees North
00471       // @param lon geodetic longitude in degrees East
00472       // @param ht height above the ellipsoid in meters
00473       // @return a reference to this object.
00474       // @throw GeometryException on invalid input
00475    Position& Position::setGeodetic(const double lat,
00476                                    const double lon,
00477                                    const double ht,
00478                                    const GeoidModel *geoid)
00479       throw(GeometryException)
00480    {
00481       if(lat > 90 || lat < -90)
00482       {
00483          GeometryException ge("Invalid latitude in setGeodetic: "
00484                                  + StringUtils::asString(lat));
00485          GPSTK_THROW(ge);
00486       }
00487       theArray[0] = lat;
00488 
00489       theArray[1] = lon;
00490       if(theArray[1] < 0)
00491          theArray[1] += 360*(1+(unsigned long)(theArray[1]/360));
00492       else if(theArray[1] >= 360)
00493          theArray[1] -= 360*(unsigned long)(theArray[1]/360);
00494 
00495       theArray[2] = ht;
00496 
00497       if(geoid) {
00498          AEarth = geoid->a();
00499          eccSquared = geoid->eccSquared();
00500       }
00501       system = Geodetic;
00502 
00503       return *this;
00504    }
00505 
00506       // Set the Position given geocentric coordinates, system is set to Geocentric
00507       // @param lat geocentric latitude in degrees North
00508       // @param lon geocentric longitude in degrees East
00509       // @param rad radius from the Earth's center in meters
00510       // @return a reference to this object.
00511       // @throw GeometryException on invalid input
00512    Position& Position::setGeocentric(const double lat,
00513                                      const double lon,
00514                                      const double rad)
00515       throw(GeometryException)
00516    {
00517       if(lat > 90 || lat < -90)
00518       {
00519          GeometryException ge("Invalid latitude in setGeocentric: "
00520                                  + StringUtils::asString(lat));
00521          GPSTK_THROW(ge);
00522       }
00523       if(rad < 0)
00524       {
00525          GeometryException ge("Invalid radius in setGeocentric: "
00526                                           + StringUtils::asString(rad));
00527          GPSTK_THROW(ge);
00528       }
00529       theArray[0] = lat;
00530       theArray[1] = lon;
00531       theArray[2] = rad;
00532 
00533       if(theArray[1] < 0)
00534          theArray[1] += 360*(1+(unsigned long)(theArray[1]/360));
00535       else if(theArray[1] >= 360)
00536          theArray[1] -= 360*(unsigned long)(theArray[1]/360);
00537       system = Geocentric;
00538 
00539       return *this;
00540    }
00541 
00542       // Set the Position given spherical coordinates, system is set to Spherical
00543       // @param theta angle from the Z-axis (degrees)
00544       // @param phi angle from the X-axis in the XY plane (degrees)
00545       // @param rad radius from the center in meters
00546       // @return a reference to this object.
00547       // @throw GeometryException on invalid input
00548    Position& Position::setSpherical(const double theta,
00549                                     const double phi,
00550                                     const double rad)
00551       throw(GeometryException)
00552    {
00553       if(theta < 0 || theta > 180)
00554       {
00555          GeometryException ge("Invalid theta in setSpherical: "
00556                                  + StringUtils::asString(theta));
00557          GPSTK_THROW(ge);
00558       }
00559       if(rad < 0)
00560       {
00561          GeometryException ge("Invalid radius in setSpherical: "
00562                                           + StringUtils::asString(rad));
00563          GPSTK_THROW(ge);
00564       }
00565 
00566       theArray[0] = theta;
00567       theArray[1] = phi;
00568       theArray[2] = rad;
00569 
00570       if(theArray[1] < 0)
00571          theArray[1] += 360*(1+(unsigned long)(theArray[1]/360));
00572       else if(theArray[1] >= 360)
00573          theArray[1] -= 360*(unsigned long)(theArray[1]/360);
00574       system = Spherical;
00575 
00576       return *this;
00577    }
00578 
00579       // Set the Position given ECEF coordinates, system is set to Cartesian.
00580       // @param X ECEF X coordinate in meters.
00581       // @param Y ECEF Y coordinate in meters.
00582       // @param Z ECEF Z coordinate in meters.
00583       // @return a reference to this object.
00584    Position& Position::setECEF(const double X,
00585                                const double Y,
00586                                const double Z)
00587       throw()
00588    {
00589       theArray[0] = X;
00590       theArray[1] = Y;
00591       theArray[2] = Z;
00592       system = Cartesian;
00593       return *this;
00594    }
00595 
00596    // ----------- Part 9: member functions: setToString, printf -------------
00597    //
00598       // setToString, similar to scanf, this function takes a string and a
00599       // format describing string in order to define Position
00600       // values.  The parameters it can take are listed below and
00601       // described above with the printf() function.
00602       //
00603       // The specification must be sufficient to define a Position.
00604       // The following table lists combinations that give valid
00605       // Positions. Anything more or other combinations will give
00606       // unknown (read as: "bad") results so don't try it.  Anything
00607       // less will throw an exception.
00608       //
00609       // @code
00610       //  %X %Y %Z  (cartesian or ECEF)
00611       //  %x %y %z  (cartesian or ECEF)
00612       //  %a %l %r  (geocentric)
00613       //  %A %L %h  (geodetic)
00614       //  %t %p %r  (spherical)
00615       // @endcode
00616       //
00617       // So
00618       // @code
00619       // pos.setToString("123.4342,9328.1982,-128987.399", "%X,%Y,%Z");
00620       // @endcode
00621       //
00622       // works but 
00623       //
00624       // @code
00625       // pos.setToString("123.4342,9328.1982", "%X,%Y");
00626       // @endcode
00627       // doesn't work (incomplete specification because it doesn't specify
00628       // a Position).
00629       //
00630       // Whitespace is unimportant here, the function will handle it.
00631       // The caller must ensure that that the extra characters in
00632       // the format string (ie '.' ',') are in the same relative
00633       // location as they are in the actual string, see the example above.
00634       //
00635       // @param str string from which to get the Position coordinates
00636       // @param fmt format to use to parse \c str.
00637       // @throw GeometryException if \c fmt is an incomplete or invalid specification
00638       // @throw FormatException if unable to scan \c str.
00639       // @throw StringException if an error occurs manipulating the
00640       // \c str or \c fmt strings.
00641       // @return a reference to this object.
00642    Position& Position::setToString(const std::string& str,
00643                                    const std::string& fmt)
00644       throw(GeometryException,
00645             DayTime::FormatException,
00646             StringUtils::StringException)
00647    {
00648       try {
00649             // make an object to return (so we don't fiddle with *this 
00650             // until it's necessary)
00651          Position toReturn;
00652          
00653             // flags indicated these defined
00654          bool hX=false, hY=false, hZ=false;
00655          bool hglat=false, hlon=false, hht=false;
00656          bool hclat=false, hrad=false;
00657          bool htheta=false, hphi=false;
00658             // store input values
00659          double x,y,z,glat,lon,ht,clat,rad,theta,phi;
00660             // copy format and input string to parse
00661          string f = fmt;
00662          string s = str;
00663          
00664             // parse strings...  As we process each part, it's removed from both
00665             // strings so when we reach 0, we're done
00666          while ( (s.size() > 0) && (f.size() > 0) )
00667          {
00668                // remove everything in f and s up to the first % in f
00669                // (these parts of the strings must be identical or this will break
00670                // after it tries to remove it!)
00671             while ( (s.length() != 0) && (f.length() != 0) && (f[0] != '%') )
00672             {
00673                   // remove that character now and other whitespace
00674                s.erase(0,1);
00675                f.erase(0,1);
00676                stripLeading(s);
00677                stripLeading(f);
00678             }
00679             
00680                // check just in case we hit the end of either string...
00681             if ( (s.length() == 0) || (f.length() == 0) )
00682                break;
00683             
00684                // lose the '%' in f...
00685             f.erase(0,1);
00686             
00687                // if the format string is like %03f, get '3' as the field
00688                // length.
00689             string::size_type fieldLength = string::npos;
00690             
00691             if (!isalpha(f[0]))
00692             {
00693                fieldLength = asInt(f);
00694                
00695                   // remove everything else up to the next character
00696                   // (in "%03f", that would be 'f')
00697                while ((!f.empty()) && (!isalpha(f[0])))
00698                   f.erase(0,1);
00699                if (f.empty())
00700                   break;
00701             }
00702             
00703                // finally, get the character that should end this field, if any
00704             char delimiter = 0;
00705             if (f.size() > 1)
00706             {
00707                if (f[1] != '%')
00708                {
00709                   delimiter = f[1];
00710                   
00711                   if (fieldLength == string::npos)
00712                      fieldLength = s.find(delimiter,0);
00713                }
00714                   // if the there is no delimiter character and the next field
00715                   // is another part of the time to parse, assume the length
00716                   // of this field is 1
00717                else if (fieldLength == string::npos)
00718                {
00719                   fieldLength = 1;
00720                }
00721             }
00722             
00723                // figure out the next string to be removed.  if there is a field
00724                // length, use that first
00725             string toBeRemoved = s.substr(0, fieldLength);
00726             
00727                // based on char at f[0], we know what to do...
00728             switch (f[0]) 
00729             {
00730           //%x   X() (meters)
00731           //%y   Y() (meters)
00732           //%z   Z() (meters)
00733           //%X   X()/1000 (kilometers)
00734           //%Y   Y()/1000 (kilometers)
00735           //%Z   Z()/1000 (kilometers)
00736                case 'X':
00737                   x = asDouble(toBeRemoved) * 1000;
00738                   hX = true;
00739                   break;
00740                case 'x':
00741                   x = asDouble(toBeRemoved);
00742                   hX = true;
00743                   break;
00744                case 'Y':
00745                   y = asDouble(toBeRemoved) * 1000;
00746                   hY = true;
00747                   break;
00748                case 'y':
00749                   y = asDouble(toBeRemoved);
00750                   hY = true;
00751                   break;
00752                case 'Z':
00753                   z = asDouble(toBeRemoved) * 1000;
00754                   hZ = true;
00755                   break;
00756                case 'z':
00757                   z = asDouble(toBeRemoved);
00758                   hZ = true;
00759                   break;
00760           //%A   geodeticLatitude() (degrees North)
00761           //%a   geocentricLatitude() (degrees North)
00762                case 'A':
00763                   glat = asDouble(toBeRemoved);
00764                   if(glat > 90. || glat < -90.) {
00765                      DayTime::FormatException f(
00766                            "Invalid geodetic latitude for setTostring: "
00767                            + toBeRemoved);
00768                      GPSTK_THROW(f);
00769                   }
00770                   hglat = true;
00771                   break;
00772                case 'a':
00773                   clat = asDouble(toBeRemoved);
00774                   if(clat > 90. || clat < -90.) {
00775                      DayTime::FormatException f(
00776                            "Invalid geocentric latitude for setTostring: "
00777                            + toBeRemoved);
00778                      GPSTK_THROW(f);
00779                   }
00780                   hclat = true;
00781                   break;
00782           //%L   longitude() (degrees East)
00783           //%l   longitude() (degrees East)
00784           //%w   longitude() (degrees West)
00785           //%W   longitude() (degrees West)
00786                case 'L':
00787                case 'l':
00788                   lon = asDouble(toBeRemoved);
00789                   if(lon < 0)
00790                      lon += 360*(1+(unsigned long)(lon/360));
00791                   else if(lon >= 360)
00792                      lon -= 360*(unsigned long)(lon/360);
00793                   hlon = true;
00794                   break;
00795                case 'w':
00796                case 'W':
00797                   lon = 360.0 - asDouble(toBeRemoved);
00798                   if(lon < 0)
00799                      lon += 360*(1+(unsigned long)(lon/360));
00800                   else if(lon >= 360)
00801                      lon -= 360*(unsigned long)(lon/360);
00802                   hlon = true;
00803                   break;
00804           //%t   theta() (degrees)
00805           //%T   theta() (radians)
00806                case 't':
00807                   theta = asDouble(toBeRemoved);
00808                   if(theta > 180. || theta < 0.) {
00809                      DayTime::FormatException f("Invalid theta for setTostring: "
00810                                                 + toBeRemoved);
00811                      GPSTK_THROW(f);
00812                   }
00813                   htheta = true;
00814                   break;
00815                case 'T':
00816                   theta = asDouble(toBeRemoved) * RAD_TO_DEG;
00817                   if(theta > 90. || theta < -90.) {
00818                      DayTime::FormatException f("Invalid theta for setTostring: "
00819                                                 + toBeRemoved);
00820                      GPSTK_THROW(f);
00821                   }
00822                   htheta = true;
00823                   break;
00824           //%p   phi() (degrees)
00825           //%P   phi() (radians)
00826                case 'p':
00827                   phi = asDouble(toBeRemoved);
00828                   if(phi < 0)
00829                      phi += 360*(1+(unsigned long)(phi/360));
00830                   else if(phi >= 360)
00831                      phi -= 360*(unsigned long)(phi/360);
00832                   hphi = true;
00833                   break;
00834                case 'P':
00835                   phi = asDouble(toBeRemoved) * RAD_TO_DEG;
00836                   if(phi < 0)
00837                      phi += 360*(1+(unsigned long)(phi/360));
00838                   else if(phi >= 360)
00839                      phi -= 360*(unsigned long)(phi/360);
00840                   hphi = true;
00841                   break;
00842           //%r   radius() meters
00843           //%R   radius()/1000 kilometers
00844           //%h   height() meters
00845           //%H   height()/1000 kilometers
00846                case 'r':
00847                   rad = asDouble(toBeRemoved);
00848                   if(rad < 0.0) {
00849                      DayTime::FormatException f("Invalid radius for setTostring: "
00850                                                 + toBeRemoved);
00851                      GPSTK_THROW(f);
00852                   }
00853                   hrad = true;
00854                   break;
00855                case 'R':
00856                   rad = asDouble(toBeRemoved) * 1000;
00857                   if(rad < 0.0) {
00858                      DayTime::FormatException f("Invalid radius for setTostring: "
00859                                                 + toBeRemoved);
00860                      GPSTK_THROW(f);
00861                   }
00862                   hrad = true;
00863                   break;
00864                case 'h':
00865                   ht = asDouble(toBeRemoved);
00866                   hht = true;
00867                   break;
00868                case 'H':
00869                   ht = asDouble(toBeRemoved) * 1000;
00870                   hht = true;
00871                   break;
00872                default: // do nothing
00873                   break;
00874             }
00875                // remove the part of s that we processed
00876             stripLeading(s,toBeRemoved,1);
00877             
00878                // remove the character we processed from f
00879             f.erase(0,1);    
00880             
00881                // check for whitespace again...
00882             stripLeading(f);
00883             stripLeading(s);
00884             
00885          }
00886          
00887          if ( s.length() != 0  ) 
00888          {
00889                // throw an error - something didn't get processed in the strings
00890             DayTime::FormatException fe(
00891                "Processing error - parts of strings left unread - " + s);
00892             GPSTK_THROW(fe);
00893          }
00894          
00895          if (f.length() != 0)
00896          {
00897                // throw an error - something didn't get processed in the strings
00898             DayTime::FormatException fe(
00899                "Processing error - parts of strings left unread - " + f);
00900             GPSTK_THROW(fe);
00901          }
00902          
00903             // throw if the specification is incomplete
00904          if ( !(hX && hY && hZ) && !(hglat && hlon && hht) &&
00905               !(hclat && hlon && hrad) && !(htheta && hphi && hrad)) {
00906             DayTime::FormatException fe("Incomplete specification for setToString");
00907             GPSTK_THROW(fe);
00908          }
00909 
00910             // define the Position toReturn
00911          if(hX && hY && hZ)
00912             toReturn.setECEF(x,y,z);
00913          else if(hglat && hlon && hht)
00914             toReturn.setGeodetic(glat,lon,ht);
00915          else if(hclat && hlon && hrad)
00916             toReturn.setGeocentric(clat,lon,rad);
00917          else if(htheta && hphi && hrad)
00918             toReturn.setSpherical(theta,phi,rad);
00919 
00920          *this = toReturn;
00921          return *this;
00922       }
00923       catch(gpstk::Exception& exc)
00924       {
00925          GeometryException ge(exc);
00926          ge.addText("Failed to convert string to Position");
00927          GPSTK_THROW(ge);
00928       }
00929       catch(std::exception& exc)
00930       {
00931          GeometryException ge(exc.what());
00932          ge.addText("Failed to convert string to Position");
00933          GPSTK_THROW(ge);
00934       }
00935    }
00936 
00937       // Format this Position into a string.
00938       //
00939       // Generate and return a string containing formatted
00940       // Position coordinates, formatted by the specification \c fmt.
00941       //
00942       // \li \%x   X() (meters)
00943       // \li \%y   Y() (meters)
00944       // \li \%z   Z() (meters)
00945       // \li \%X   X()/1000 (kilometers)
00946       // \li \%Y   Y()/1000 (kilometers)
00947       // \li \%Z   Z()/1000 (kilometers)
00948       // \li \%A   geodeticLatitude() (degrees North)
00949       // \li \%a   geocentricLatitude() (degrees North)
00950       // \li \%L   longitude() (degrees East)
00951       // \li \%l   longitude() (degrees East)
00952       // \li \%w   longitude() (degrees West)
00953       // \li \%W   longitude() (degrees West)
00954       // \li \%t   theta() (degrees)
00955       // \li \%T   theta() (radians)
00956       // \li \%p   phi() (degrees)
00957       // \li \%P   phi() (radians)
00958       // \li \%r   radius() meters
00959       // \li \%R   radius()/1000 kilometers
00960       // \li \%h   height() meters
00961       // \li \%H   height()/1000 kilometers
00962       //
00963       // @param fmt format to use for this time.
00964       // @return a string containing this Position in the
00965       // representation specified by \c fmt.
00966    std::string Position::printf(const char *fmt) const
00967       throw(StringUtils::StringException)
00968    {
00969       string rv = fmt;
00970       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?x"),
00971                           string("xf"), X());
00972       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?y"),
00973                           string("yf"), Y());
00974       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?z"),
00975                           string("zf"), Z());
00976       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?X"),
00977                           string("Xf"), X()/1000);
00978       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?Y"),
00979                           string("Yf"), Y()/1000);
00980       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?Z"),
00981                           string("Zf"), Z()/1000);
00982 
00983       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?A"),
00984                           string("Af"), geodeticLatitude());
00985       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?a"),
00986                           string("af"), geocentricLatitude());
00987       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?L"),
00988                           string("Lf"), longitude());
00989       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?l"),
00990                           string("lf"), longitude());
00991       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?w"),
00992                           string("wf"), 360-longitude());
00993       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?W"),
00994                           string("Wf"), 360-longitude());
00995 
00996       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?t"),
00997                           string("tf"), theta());
00998       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?T"),
00999                           string("Tf"), theta()*DEG_TO_RAD);
01000       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?p"),
01001                           string("pf"), phi());
01002       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?P"),
01003                           string("Pf"), phi()*DEG_TO_RAD);
01004       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?r"),
01005                           string("rf"), radius());
01006       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?R"),
01007                           string("Rf"), radius()/1000);
01008       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?h"),
01009                           string("hf"), height());
01010       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?H"),
01011                           string("Hf"), height()/1000);
01012       return rv;
01013    }
01014 
01015       // Returns the string that operator<<() would print.
01016    string Position::asString() const
01017       throw(StringUtils::StringException)
01018    {
01019       ostringstream o;
01020       o << *this;
01021       return o.str();
01022    }
01023 
01024    // ----------- Part 10: functions: fundamental conversions ---------------
01025    // 
01026       // Fundamental conversion from spherical to cartesian coordinates.
01027       // @param trp (input): theta, phi, radius
01028       // @param xyz (output): X,Y,Z in units of radius
01029       // Algorithm references: standard geometry.
01030    void Position::convertSphericalToCartesian(const Triple& tpr,
01031                                               Triple& xyz)
01032       throw()
01033    {
01034       double st=sin(tpr[0]*DEG_TO_RAD);
01035       xyz[0] = tpr[2]*st*cos(tpr[1]*DEG_TO_RAD);
01036       xyz[1] = tpr[2]*st*sin(tpr[1]*DEG_TO_RAD);
01037       xyz[2] = tpr[2]*cos(tpr[0]*DEG_TO_RAD);
01038    }
01039 
01040       // Fundamental routine to convert cartesian to spherical coordinates.
01041       // @param xyz (input): X,Y,Z
01042       // @param trp (output): theta, phi (deg), radius in units of input
01043       // Algorithm references: standard geometry.
01044    void Position::convertCartesianToSpherical(const Triple& xyz,
01045                                               Triple& tpr)
01046       throw()
01047    {
01048       tpr[2] = RSS(xyz[0],xyz[1],xyz[2]);
01049       if(tpr[2] <= Position::POSITION_TOLERANCE/5)
01050       {
01051             // zero-length Cartesian vector
01052          tpr[0] = 90;
01053          tpr[1] = 0;
01054          return;
01055       }
01056       tpr[0] = acos(xyz[2]/tpr[2]);
01057       tpr[0] *= RAD_TO_DEG;
01058       if(RSS(xyz[0],xyz[1]) < Position::POSITION_TOLERANCE/5) {       // pole
01059          tpr[1] = 0;
01060          return;
01061       }
01062       tpr[1] = atan2(xyz[1],xyz[0]);
01063       tpr[1] *= RAD_TO_DEG;
01064       if(tpr[1] < 0) tpr[1] += 360;
01065    }
01066 
01067       // Fundamental routine to convert cartesian (ECEF) to geodetic coordinates,
01068       // (Geoid specified by semi-major axis and eccentricity squared).
01069       // @param xyz (input): X,Y,Z in meters
01070       // @param llh (output): geodetic lat(deg N), lon(deg E),
01071       //                             height above ellipsoid (meters)
01072       // @param A (input) Earth semi-major axis
01073       // @param eccSq (input) square of Earth eccentricity
01074       // Algorithm references: 
01075    void Position::convertCartesianToGeodetic(const Triple& xyz,
01076                                              Triple& llh,
01077                                              const double A,
01078                                              const double eccSq)
01079       throw()
01080    {
01081       double p,slat,N,htold,latold;
01082       p = SQRT(xyz[0]*xyz[0]+xyz[1]*xyz[1]);
01083       if(p < Position::POSITION_TOLERANCE/5) {   // pole or origin
01084          llh[0] = llh[1] = 0;
01085          llh[2] = fabs(xyz[2]) - A;
01086          return;
01087       }
01088       llh[0] = atan2(xyz[2], p*(1.0-eccSq));
01089       llh[2] = 0;
01090       for(int i=0; i<5; i++) {
01091          slat = sin(llh[0]);
01092          N = A / SQRT(1.0 - eccSq*slat*slat);
01093          htold = llh[2];
01094          llh[2] = p/cos(llh[0]) - N;
01095          latold = llh[0];
01096          llh[0] = atan2(xyz[2], p*(1.0-eccSq*(N/(N+llh[2]))));
01097          if(fabs(llh[0]-latold) < 1.0e-9 && fabs(llh[2]-htold) < 1.0e-9 * A) break;
01098       }
01099       llh[1] = atan2(xyz[1],xyz[0]);
01100       if(llh[1] < 0.0) llh[1] += TWO_PI;
01101       llh[0] *= RAD_TO_DEG;
01102       llh[1] *= RAD_TO_DEG;
01103    }
01104 
01105       // Fundamental routine to convert geodetic to cartesian (ECEF) coordinates,
01106       // (Geoid specified by semi-major axis and eccentricity squared).
01107       // @param llh (input): geodetic lat(deg N), lon(deg E),
01108       //            height above ellipsoid (meters)
01109       // @param xyz (output): X,Y,Z in meters
01110       // @param A (input) Earth semi-major axis
01111       // @param eccSq (input) square of Earth eccentricity
01112       // Algorithm references: 
01113    void Position::convertGeodeticToCartesian(const Triple& llh,
01114                                              Triple& xyz,
01115                                              const double A,
01116                                              const double eccSq)
01117       throw()
01118    {
01119       double slat = sin(llh[0]*DEG_TO_RAD);
01120       double clat = cos(llh[0]*DEG_TO_RAD);
01121       double N = A/SQRT(1.0-eccSq*slat*slat);
01122       xyz[0] = (N+llh[2])*clat*cos(llh[1]*DEG_TO_RAD);
01123       xyz[1] = (N+llh[2])*clat*sin(llh[1]*DEG_TO_RAD);
01124       xyz[2] = (N*(1.0-eccSq)+llh[2])*slat;
01125    }
01126 
01127       // Fundamental routine to convert cartesian (ECEF) to geocentric coordinates.
01128       // @param xyz (input): X,Y,Z in meters
01129       // @param llr (output):
01130       //            geocentric lat(deg N),lon(deg E),radius (units of input)
01131    void Position::convertCartesianToGeocentric(const Triple& xyz,
01132                                                Triple& llr)
01133       throw()
01134    {
01135       convertCartesianToSpherical(xyz, llr);
01136       llr[0] = 90 - llr[0];         // convert theta to latitude
01137    }
01138 
01139       // Fundamental routine to convert geocentric to cartesian (ECEF) coordinates.
01140       // @param llr (input): geocentric lat(deg N),lon(deg E),radius
01141       // @param xyz (output): X,Y,Z (units of radius)
01142    void Position::convertGeocentricToCartesian(const Triple& llr,
01143                                                Triple& xyz)
01144       throw()
01145    {
01146       Triple llh(llr);
01147       llh[0] = 90 - llh[0];         // convert latitude to theta
01148       convertSphericalToCartesian(llh, xyz);
01149    }
01150 
01151       // Fundamental routine to convert geocentric to geodetic coordinates.
01152       // @param llr (input): geocentric Triple: lat(deg N),lon(deg E),radius (meters)
01153       // @param llh (output): geodetic latitude (deg N),
01154       //            longitude (deg E), and height above ellipsoid (meters)
01155       // @param A (input) Earth semi-major axis
01156       // @param eccSq (input) square of Earth eccentricity
01157    void Position::convertGeocentricToGeodetic(const Triple& llr,
01158                                                Triple& llh,
01159                                                const double A,
01160                                                const double eccSq)
01161       throw()
01162    {
01163       double cl,p,sl,slat,N,htold,latold;
01164       cl = sin((90-llr[0])*DEG_TO_RAD);
01165       sl = cos((90-llr[0])*DEG_TO_RAD);
01166       if(llr[2] <= Position::POSITION_TOLERANCE/5) {
01167          // radius is below tolerance, hence assign zero-length
01168          // arbitrarily set latitude = longitude = 0
01169          llh[0] = llh[1] = 0;
01170          llh[2] = -A;
01171          return;
01172       }
01173       else if(cl < 1.e-10) {
01174          // near pole ... note that 1mm/radius(Earth) = 1.5e-10
01175          if(llr[0] < 0) llh[0] = -90;
01176          else           llh[0] =  90;
01177          llh[1] = 0;
01178          llh[2] = llr[2] - A*SQRT(1-eccSq);
01179          return;
01180       }
01181       llh[0] = atan2(sl, cl*(1.0-eccSq));
01182       p = cl*llr[2];
01183       llh[2] = 0;
01184       for(int i=0; i<5; i++) {
01185          slat = sin(llh[0]);
01186          N = A / SQRT(1.0 - eccSq*slat*slat);
01187          htold = llh[2];
01188          llh[2] = p/cos(llh[0]) - N;
01189          latold = llh[0];
01190          llh[0] = atan2(sl, cl*(1.0-eccSq*(N/(N+llh[2]))));
01191          if(fabs(llh[0]-latold) < 1.0e-9 && fabs(llh[2]-htold) < 1.0e-9 * A) break;
01192       }
01193       llh[0] *= RAD_TO_DEG;
01194    }
01195 
01196       // Fundamental routine to convert geodetic to geocentric coordinates.
01197       // @param geodeticllh (input): geodetic latitude (deg N),
01198       //            longitude (deg E), and height above ellipsoid (meters)
01199       // @param llr (output): geocentric lat (deg N),lon (deg E),radius (meters)
01200       // @param A (input) Earth semi-major axis
01201       // @param eccSq (input) square of Earth eccentricity
01202    void Position::convertGeodeticToGeocentric(const Triple& llh,
01203                                               Triple& llr,
01204                                               const double A,
01205                                               const double eccSq)
01206       throw()
01207    {
01208       double slat = sin(llh[0]*DEG_TO_RAD);
01209       double N = A/SQRT(1.0-eccSq*slat*slat);
01210       // radius
01211       llr[2] = SQRT((N+llh[2])*(N+llh[2]) + N*eccSq*(N*eccSq-2*(N+llh[2]))*slat*slat);
01212       if(llr[2] <= Position::POSITION_TOLERANCE/5)
01213       {
01214             // radius is below tolerance, hence assign zero-length
01215             // arbitrarily set latitude = longitude = 0
01216          llr[0] = llr[1] = llr[2] = 0;
01217          return;
01218       }
01219       if(1-fabs(slat) < 1.e-10) {             // at the pole
01220          if(slat < 0) llr[0] = -90;
01221          else         llr[0] =  90;
01222          return;
01223       }
01224       // theta
01225       llr[0] = acos((N*(1-eccSq)+llh[2])*slat/llr[2]);
01226       llr[0] *= RAD_TO_DEG;
01227       llr[0] = 90 - llr[0];
01228    }
01229 
01230    // ----------- Part 11: operator<< and other useful functions -------------
01231    //
01232      // Stream output for Position objects.
01233      // @param s stream to append formatted Position to.
01234      // @param t Position to append to stream \c s.
01235      // @return reference to \c s.
01236    ostream& operator<<(ostream& s, const Position& p)
01237    {
01238       if(p.system == Position::Cartesian)
01239          s << p.printf("%.4x m %.4y m %.4z m");
01240       else if(p.system == Position::Geodetic)
01241          s << p.printf("%.8A degN %.8L degE %.4h m");
01242       else if(p.system == Position::Geocentric)
01243          s << p.printf("%.8a degN %.8L degE %.4r m");
01244       else if(p.system == Position::Spherical)
01245          s << p.printf("%.8t deg %.8p deg %.4r m");
01246       else
01247          s << " Unknown system! : " << p[0] << " " << p[1] << " " << p[2];
01248 
01249       return s;
01250    }
01251 
01252       // Compute the range in meters between this Position and
01253       // the Position passed as input.
01254       // @param right Position to which to find the range
01255       // @return the range (in meters)
01256       // @throw GeometryException if geoid values differ
01257    double range(const Position& A,
01258                 const Position& B)
01259       throw(GeometryException)
01260    {
01261       if(A.AEarth != B.AEarth || A.eccSquared != B.eccSquared)
01262       {
01263          GeometryException ge("Unequal geoids");
01264          GPSTK_THROW(ge);
01265       }
01266 
01267          Position L(A),R(B);
01268          L.transformTo(Position::Cartesian);
01269          R.transformTo(Position::Cartesian);
01270          double dif = RSS(L.X()-R.X(),L.Y()-R.Y(),L.Z()-R.Z());
01271          return dif;
01272       }
01273 
01274      // Compute the radius of the ellipsoidal Earth, given the geodetic latitude.
01275      // @param geolat geodetic latitude in degrees
01276      // @return the Earth radius (in meters)
01277    double Position::radiusEarth(const double geolat,
01278                                 const double A,
01279                                 const double eccSq)
01280       throw()
01281    {
01282       double slat=sin(DEG_TO_RAD*geolat);
01283       double e=(1.0-eccSq);
01284       double f=(1.0+(e*e-1.0)*slat*slat)/(1.0-eccSq*slat*slat);
01285       return (A * SQRT(f));
01286    }
01287 
01288       // A member function that computes the elevation of the input
01289       // (Target) position as seen from this Position.
01290       // @param Target the Position which is observed to have the
01291       //        computed elevation, as seen from this Position.
01292       // @return the elevation in degrees
01293    double Position::elevation(const Position& Target) const
01294       throw(GeometryException)
01295    {
01296       Position R(*this),S(Target);
01297       R.transformTo(Cartesian);
01298       S.transformTo(Cartesian);
01299       // use Triple:: functions in cartesian coordinates (only)
01300       double elevation;
01301       try {
01302          elevation = R.elvAngle(S);         
01303       }
01304       catch(GeometryException& ge)
01305       {
01306          GPSTK_RETHROW(ge);
01307       }
01308       return elevation;
01309    }
01310 
01311       // A member function that computes the elevation of the input
01312       // (Target) position as seen from this Position, using a Geodetic
01313       // (i.e. ellipsoidal) system.
01314       // @param Target the Position which is observed to have the
01315       //        computed elevation, as seen from this Position.
01316       // @return the elevation in degrees
01317    double Position::elevationGeodetic(const Position& Target) const
01318       throw(GeometryException)
01319    {
01320       Position R(*this),S(Target);
01321       double latGeodetic = R.getGeodeticLatitude()*DEG_TO_RAD;
01322       double longGeodetic = R.getLongitude()*DEG_TO_RAD;
01323       double localUp;
01324       double cosUp;
01325       R.transformTo(Cartesian);
01326       S.transformTo(Cartesian);
01327       Triple z;
01328       // Let's get the slant vector
01329       z = S.theArray - R.theArray;
01330 
01331       if (z.mag()<=1e-4) // if the positions are within .1 millimeter
01332       {
01333          GeometryException ge("Positions are within .1 millimeter");
01334          GPSTK_THROW(ge);
01335       }
01336 
01337       // Compute k vector in local North-East-Up (NEU) system
01338       Triple kVector(cos(latGeodetic)*cos(longGeodetic), cos(latGeodetic)*sin(longGeodetic), sin(latGeodetic));
01339       // Take advantage of dot method to get Up coordinate in local NEU system
01340       localUp = z.dot(kVector);
01341       // Let's get cos(z), being z the angle with respect to local vertical (Up);
01342       cosUp = localUp/z.mag();
01343 
01344       return 90.0 - ((::acos(cosUp))*RAD_TO_DEG);
01345    }
01346 
01347       // A member function that computes the azimuth of the input
01348       // (Target) position as seen from this Position.
01349       // @param Target the Position which is observed to have the
01350       //        computed azimuth, as seen from this Position.
01351       // @return the azimuth in degrees
01352    double Position::azimuth(const Position& Target) const
01353       throw(GeometryException)
01354    {
01355       Position R(*this),S(Target);
01356       R.transformTo(Cartesian);
01357       S.transformTo(Cartesian);
01358       // use Triple:: functions in cartesian coordinates (only)
01359       double az;
01360       try
01361       {
01362          az = R.azAngle(S);
01363          
01364       }
01365       catch(GeometryException& ge)
01366       {
01367          GPSTK_RETHROW(ge);
01368       }
01369       
01370       return az; 
01371    }
01372 
01373       // A member function that computes the azimuth of the input
01374       // (Target) position as seen from this Position, using a Geodetic
01375       // (i.e. ellipsoidal) system.
01376       // @param Target the Position which is observed to have the
01377       //        computed azimuth, as seen from this Position.
01378       // @return the azimuth in degrees
01379    double Position::azimuthGeodetic(const Position& Target) const
01380       throw(GeometryException)
01381    {
01382       Position R(*this),S(Target);
01383       double latGeodetic = R.getGeodeticLatitude()*DEG_TO_RAD;
01384       double longGeodetic = R.getLongitude()*DEG_TO_RAD;
01385       double localN, localE;
01386       R.transformTo(Cartesian);
01387       S.transformTo(Cartesian);
01388       Triple z;
01389       // Let's get the slant vector
01390       z = S.theArray - R.theArray;
01391 
01392       if (z.mag()<=1e-4) // if the positions are within .1 millimeter
01393       {
01394          GeometryException ge("Positions are within .1 millimeter");
01395          GPSTK_THROW(ge);
01396       }
01397       
01398       // Compute i vector in local North-East-Up (NEU) system
01399       Triple iVector(-sin(latGeodetic)*cos(longGeodetic), -sin(latGeodetic)*sin(longGeodetic), cos(latGeodetic));
01400       // Compute j vector in local North-East-Up (NEU) system
01401       Triple jVector(-sin(longGeodetic), cos(longGeodetic), 0);
01402 
01403       // Now, let's use dot product to get localN and localE unitary vectors
01404       localN = (z.dot(iVector))/z.mag();
01405       localE = (z.dot(jVector))/z.mag();
01406 
01407       // Let's test if computing azimuth has any sense
01408       double test = fabs(localN) + fabs(localE);
01409 
01410       // Warning: If elevation is very close to 90 degrees, we will return azimuth = 0.0
01411       if (test < 1.0e-16) return 0.0;
01412 
01413       double alpha = ((::atan2(localE, localN)) * RAD_TO_DEG);
01414       if (alpha < 0.0)
01415       {
01416          return alpha + 360.0;
01417       }
01418       else 
01419       {
01420          return alpha;
01421       }
01422    }
01423 
01424      // A member function that computes the point at which a signal, which
01425      // is received at *this Position and there is observed at the input
01426      // azimuth and elevation, crosses a model ionosphere that is taken to
01427      // be a uniform thin shell at the input height. This algorithm is done
01428      // in geocentric coordinates.
01429      // A member function that computes the point at which a signal, which
01430      // is received at *this Position and there is observed at the input
01431      // azimuth and elevation, crosses a model ionosphere that is taken to
01432      // be a uniform thin shell at the input height. This algorithm is done
01433      // in geocentric coordinates.
01434      // @param elev elevation angle of the signal at reception, in degrees
01435      // @param azim azimuth angle of the signal at reception, in degrees
01436      // @param ionoht height of the ionosphere, in meters
01437      // @return Position IPP the position of the ionospheric pierce point,
01438      //     in the same coordinate system as *this; *this is not modified.
01439    Position Position::getIonosphericPiercePoint(const double elev,
01440                                                 const double azim,
01441                                                 const double ionoht) const
01442       throw()
01443    {
01444       Position Rx(*this);
01445 
01446       // convert to Geocentric
01447       Rx.transformTo(Geocentric);
01448 
01449       // compute the geographic pierce point
01450       Position IPP(Rx);                   // copy system and geoid
01451       double el = elev * DEG_TO_RAD;
01452       // p is the angle subtended at Earth center by Rx and the IPP
01453       double p = PI/2.0 - el - asin(AEarth*cos(el)/(AEarth+ionoht));
01454       double lat = Rx.theArray[0] * DEG_TO_RAD;
01455       double az = azim * DEG_TO_RAD;
01456       IPP.theArray[0] = asin(sin(lat)*cos(p) + cos(lat)*sin(p)*cos(az));
01457       IPP.theArray[1] = Rx.theArray[1]*DEG_TO_RAD
01458          + asin(sin(p)*sin(az)/cos(IPP.theArray[0]));
01459 
01460       IPP.theArray[0] *= RAD_TO_DEG;
01461       IPP.theArray[1] *= RAD_TO_DEG;
01462       IPP.theArray[2] = AEarth + ionoht;
01463 
01464       // transform back
01465       IPP.transformTo(system);
01466 
01467       return IPP;
01468    }
01469 
01470    // ----------- Part 12: private functions and member data -----------------
01471    //
01472       // Initialization function, used by the constructors.
01473       // @param a coordinate [ X(m), or latitude (degrees N) ]
01474       // @param b coordinate [ Y(m), or longitude (degrees E) ]
01475       // @param c coordinate [ Z, height above ellipsoid or radius, in m ]
01476       // @param s CoordinateSystem, defaults to Cartesian
01477       // @param geiod pointer to a GeoidModel, default NULL (WGS84)
01478       // @throw GeometryException on invalid input.
01479    void Position::initialize(const double a,
01480                   const double b,
01481                   const double c,
01482                   Position::CoordinateSystem s,
01483                   GeoidModel *geoid)
01484       throw(GeometryException)
01485    {
01486       double bb(b);
01487       if(s == Geodetic || s==Geocentric)
01488       {
01489          if(a > 90 || a < -90)
01490          {
01491             GeometryException ge("Invalid latitude in constructor: "
01492                                     + StringUtils::asString(a));
01493             GPSTK_THROW(ge);
01494          }
01495          if(bb < 0)
01496             bb += 360*(1+(unsigned long)(bb/360));
01497          else if(bb >= 360)
01498             bb -= 360*(unsigned long)(bb/360);
01499       }
01500       if(s==Geocentric || s==Spherical)
01501       {
01502          if(c < 0)
01503          {
01504             GeometryException ge("Invalid radius in constructor: "
01505                                            + StringUtils::asString(c));
01506             GPSTK_THROW(ge);
01507          }
01508       }
01509       if(s==Spherical)
01510       {
01511          if(a < 0 || a > 180)
01512          {
01513             GeometryException ge("Invalid theta in constructor: "
01514                                     + StringUtils::asString(a));
01515             GPSTK_THROW(ge);
01516          }
01517          if(bb < 0)
01518             bb += 360*(1+(unsigned long)(bb/360));
01519          else if(bb >= 360)
01520             bb -= 360*(unsigned long)(bb/360);
01521       }
01522 
01523       theArray[0] = a;
01524       theArray[1] = bb;
01525       theArray[2] = c;
01526 
01527       if(geoid) {
01528          AEarth = geoid->a();
01529          eccSquared = geoid->eccSquared();
01530       }
01531       else {
01532          WGS84Geoid WGS84;
01533          AEarth = WGS84.a();
01534          eccSquared = WGS84.eccSquared();
01535       }
01536       system = s;
01537       tolerance = POSITION_TOLERANCE;
01538    }
01539 
01540 }  // namespace gpstk

Generated on Tue Jan 6 03:31:23 2009 for GPS ToolKit Software Library by  doxygen 1.3.9.1