Position.cpp

Go to the documentation of this file.
00001 #pragma ident "$Id: Position.cpp 2033 2009-07-17 16:08:02Z afarris $"
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          stripLeading(s);
00664          stripLeading(f);
00665          
00666             // parse strings...  As we process each part, it's removed from both
00667             // strings so when we reach 0, we're done
00668          while ( (s.size() > 0) && (f.size() > 0) )
00669          {
00670                // remove everything in f and s up to the first % in f
00671                // (these parts of the strings must be identical or this will break
00672                // after it tries to remove it!)
00673             while ( (s.length() != 0) && (f.length() != 0) && (f[0] != '%') )
00674             {
00675                   // remove that character now and other whitespace
00676                s.erase(0,1);
00677                f.erase(0,1);
00678                stripLeading(s);
00679                stripLeading(f);
00680             }
00681             
00682                // check just in case we hit the end of either string...
00683             if ( (s.length() == 0) || (f.length() == 0) )
00684                break;
00685             
00686                // lose the '%' in f...
00687             f.erase(0,1);
00688             
00689                // if the format string is like %03f, get '3' as the field
00690                // length.
00691             string::size_type fieldLength = string::npos;
00692             
00693             if (!isalpha(f[0]))
00694             {
00695                fieldLength = asInt(f);
00696                
00697                   // remove everything else up to the next character
00698                   // (in "%03f", that would be 'f')
00699                while ((!f.empty()) && (!isalpha(f[0])))
00700                   f.erase(0,1);
00701                if (f.empty())
00702                   break;
00703             }
00704             
00705                // finally, get the character that should end this field, if any
00706             char delimiter = 0;
00707             if (f.size() > 1)
00708             {
00709                if (f[1] != '%')
00710                {
00711                   delimiter = f[1];
00712                   
00713                   if (fieldLength == string::npos)
00714                      fieldLength = s.find(delimiter,0);
00715                }
00716                   // if the there is no delimiter character and the next field
00717                   // is another part of the time to parse, assume the length
00718                   // of this field is 1
00719                else if (fieldLength == string::npos)
00720                {
00721                   fieldLength = 1;
00722                }
00723             }
00724             
00725                // figure out the next string to be removed.  if there is a field
00726                // length, use that first
00727             string toBeRemoved = s.substr(0, fieldLength);
00728             
00729                // based on char at f[0], we know what to do...
00730             switch (f[0]) 
00731             {
00732           //%x   X() (meters)
00733           //%y   Y() (meters)
00734           //%z   Z() (meters)
00735           //%X   X()/1000 (kilometers)
00736           //%Y   Y()/1000 (kilometers)
00737           //%Z   Z()/1000 (kilometers)
00738                case 'X':
00739                   x = asDouble(toBeRemoved) * 1000;
00740                   hX = true;
00741                   break;
00742                case 'x':
00743                   x = asDouble(toBeRemoved);
00744                   hX = true;
00745                   break;
00746                case 'Y':
00747                   y = asDouble(toBeRemoved) * 1000;
00748                   hY = true;
00749                   break;
00750                case 'y':
00751                   y = asDouble(toBeRemoved);
00752                   hY = true;
00753                   break;
00754                case 'Z':
00755                   z = asDouble(toBeRemoved) * 1000;
00756                   hZ = true;
00757                   break;
00758                case 'z':
00759                   z = asDouble(toBeRemoved);
00760                   hZ = true;
00761                   break;
00762           //%A   geodeticLatitude() (degrees North)
00763           //%a   geocentricLatitude() (degrees North)
00764                case 'A':
00765                   glat = asDouble(toBeRemoved);
00766                   if(glat > 90. || glat < -90.) {
00767                      DayTime::FormatException f(
00768                            "Invalid geodetic latitude for setTostring: "
00769                            + toBeRemoved);
00770                      GPSTK_THROW(f);
00771                   }
00772                   hglat = true;
00773                   break;
00774                case 'a':
00775                   clat = asDouble(toBeRemoved);
00776                   if(clat > 90. || clat < -90.) {
00777                      DayTime::FormatException f(
00778                            "Invalid geocentric latitude for setTostring: "
00779                            + toBeRemoved);
00780                      GPSTK_THROW(f);
00781                   }
00782                   hclat = true;
00783                   break;
00784           //%L   longitude() (degrees East)
00785           //%l   longitude() (degrees East)
00786           //%w   longitude() (degrees West)
00787           //%W   longitude() (degrees West)
00788                case 'L':
00789                case 'l':
00790                   lon = asDouble(toBeRemoved);
00791                   if(lon < 0)
00792                      lon += 360*(1+(unsigned long)(lon/360));
00793                   else if(lon >= 360)
00794                      lon -= 360*(unsigned long)(lon/360);
00795                   hlon = true;
00796                   break;
00797                case 'w':
00798                case 'W':
00799                   lon = 360.0 - asDouble(toBeRemoved);
00800                   if(lon < 0)
00801                      lon += 360*(1+(unsigned long)(lon/360));
00802                   else if(lon >= 360)
00803                      lon -= 360*(unsigned long)(lon/360);
00804                   hlon = true;
00805                   break;
00806           //%t   theta() (degrees)
00807           //%T   theta() (radians)
00808                case 't':
00809                   theta = asDouble(toBeRemoved);
00810                   if(theta > 180. || theta < 0.) {
00811                      DayTime::FormatException f("Invalid theta for setTostring: "
00812                                                 + toBeRemoved);
00813                      GPSTK_THROW(f);
00814                   }
00815                   htheta = true;
00816                   break;
00817                case 'T':
00818                   theta = asDouble(toBeRemoved) * RAD_TO_DEG;
00819                   if(theta > 90. || theta < -90.) {
00820                      DayTime::FormatException f("Invalid theta for setTostring: "
00821                                                 + toBeRemoved);
00822                      GPSTK_THROW(f);
00823                   }
00824                   htheta = true;
00825                   break;
00826           //%p   phi() (degrees)
00827           //%P   phi() (radians)
00828                case 'p':
00829                   phi = asDouble(toBeRemoved);
00830                   if(phi < 0)
00831                      phi += 360*(1+(unsigned long)(phi/360));
00832                   else if(phi >= 360)
00833                      phi -= 360*(unsigned long)(phi/360);
00834                   hphi = true;
00835                   break;
00836                case 'P':
00837                   phi = asDouble(toBeRemoved) * RAD_TO_DEG;
00838                   if(phi < 0)
00839                      phi += 360*(1+(unsigned long)(phi/360));
00840                   else if(phi >= 360)
00841                      phi -= 360*(unsigned long)(phi/360);
00842                   hphi = true;
00843                   break;
00844           //%r   radius() meters
00845           //%R   radius()/1000 kilometers
00846           //%h   height() meters
00847           //%H   height()/1000 kilometers
00848                case 'r':
00849                   rad = asDouble(toBeRemoved);
00850                   if(rad < 0.0) {
00851                      DayTime::FormatException f("Invalid radius for setTostring: "
00852                                                 + toBeRemoved);
00853                      GPSTK_THROW(f);
00854                   }
00855                   hrad = true;
00856                   break;
00857                case 'R':
00858                   rad = asDouble(toBeRemoved) * 1000;
00859                   if(rad < 0.0) {
00860                      DayTime::FormatException f("Invalid radius for setTostring: "
00861                                                 + toBeRemoved);
00862                      GPSTK_THROW(f);
00863                   }
00864                   hrad = true;
00865                   break;
00866                case 'h':
00867                   ht = asDouble(toBeRemoved);
00868                   hht = true;
00869                   break;
00870                case 'H':
00871                   ht = asDouble(toBeRemoved) * 1000;
00872                   hht = true;
00873                   break;
00874                default: // do nothing
00875                   break;
00876             }
00877                // remove the part of s that we processed
00878             stripLeading(s,toBeRemoved,1);
00879             
00880                // remove the character we processed from f
00881             f.erase(0,1);    
00882             
00883                // check for whitespace again...
00884             stripLeading(f);
00885             stripLeading(s);
00886             
00887          }
00888          
00889          if ( s.length() != 0  ) 
00890          {
00891                // throw an error - something didn't get processed in the strings
00892             DayTime::FormatException fe(
00893                "Processing error - parts of strings left unread - " + s);
00894             GPSTK_THROW(fe);
00895          }
00896          
00897          if (f.length() != 0)
00898          {
00899                // throw an error - something didn't get processed in the strings
00900             DayTime::FormatException fe(
00901                "Processing error - parts of strings left unread - " + f);
00902             GPSTK_THROW(fe);
00903          }
00904          
00905             // throw if the specification is incomplete
00906          if ( !(hX && hY && hZ) && !(hglat && hlon && hht) &&
00907               !(hclat && hlon && hrad) && !(htheta && hphi && hrad)) {
00908             DayTime::FormatException fe("Incomplete specification for setToString");
00909             GPSTK_THROW(fe);
00910          }
00911 
00912             // define the Position toReturn
00913          if(hX && hY && hZ)
00914             toReturn.setECEF(x,y,z);
00915          else if(hglat && hlon && hht)
00916             toReturn.setGeodetic(glat,lon,ht);
00917          else if(hclat && hlon && hrad)
00918             toReturn.setGeocentric(clat,lon,rad);
00919          else if(htheta && hphi && hrad)
00920             toReturn.setSpherical(theta,phi,rad);
00921 
00922          *this = toReturn;
00923          return *this;
00924       }
00925       catch(gpstk::Exception& exc)
00926       {
00927          GeometryException ge(exc);
00928          ge.addText("Failed to convert string to Position");
00929          GPSTK_THROW(ge);
00930       }
00931       catch(std::exception& exc)
00932       {
00933          GeometryException ge(exc.what());
00934          ge.addText("Failed to convert string to Position");
00935          GPSTK_THROW(ge);
00936       }
00937    }
00938 
00939       // Format this Position into a string.
00940       //
00941       // Generate and return a string containing formatted
00942       // Position coordinates, formatted by the specification \c fmt.
00943       //
00944       // \li \%x   X() (meters)
00945       // \li \%y   Y() (meters)
00946       // \li \%z   Z() (meters)
00947       // \li \%X   X()/1000 (kilometers)
00948       // \li \%Y   Y()/1000 (kilometers)
00949       // \li \%Z   Z()/1000 (kilometers)
00950       // \li \%A   geodeticLatitude() (degrees North)
00951       // \li \%a   geocentricLatitude() (degrees North)
00952       // \li \%L   longitude() (degrees East)
00953       // \li \%l   longitude() (degrees East)
00954       // \li \%w   longitude() (degrees West)
00955       // \li \%W   longitude() (degrees West)
00956       // \li \%t   theta() (degrees)
00957       // \li \%T   theta() (radians)
00958       // \li \%p   phi() (degrees)
00959       // \li \%P   phi() (radians)
00960       // \li \%r   radius() meters
00961       // \li \%R   radius()/1000 kilometers
00962       // \li \%h   height() meters
00963       // \li \%H   height()/1000 kilometers
00964       //
00965       // @param fmt format to use for this time.
00966       // @return a string containing this Position in the
00967       // representation specified by \c fmt.
00968    std::string Position::printf(const char *fmt) const
00969       throw(StringUtils::StringException)
00970    {
00971       string rv = fmt;
00972       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?x"),
00973                           string("xf"), X());
00974       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?y"),
00975                           string("yf"), Y());
00976       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?z"),
00977                           string("zf"), Z());
00978       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?X"),
00979                           string("Xf"), X()/1000);
00980       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?Y"),
00981                           string("Yf"), Y()/1000);
00982       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?Z"),
00983                           string("Zf"), Z()/1000);
00984 
00985       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?A"),
00986                           string("Af"), geodeticLatitude());
00987       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?a"),
00988                           string("af"), geocentricLatitude());
00989       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?L"),
00990                           string("Lf"), longitude());
00991       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?l"),
00992                           string("lf"), longitude());
00993       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?w"),
00994                           string("wf"), 360-longitude());
00995       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?W"),
00996                           string("Wf"), 360-longitude());
00997 
00998       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?t"),
00999                           string("tf"), theta());
01000       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?T"),
01001                           string("Tf"), theta()*DEG_TO_RAD);
01002       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?p"),
01003                           string("pf"), phi());
01004       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?P"),
01005                           string("Pf"), phi()*DEG_TO_RAD);
01006       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?r"),
01007                           string("rf"), radius());
01008       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?R"),
01009                           string("Rf"), radius()/1000);
01010       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?h"),
01011                           string("hf"), height());
01012       rv = formattedPrint(rv, string("%[ 0-]?[[:digit:]]*(\\.[[:digit:]]+)?H"),
01013                           string("Hf"), height()/1000);
01014       return rv;
01015    }
01016 
01017       // Returns the string that operator<<() would print.
01018    string Position::asString() const
01019       throw(StringUtils::StringException)
01020    {
01021       ostringstream o;
01022       o << *this;
01023       return o.str();
01024    }
01025 
01026    // ----------- Part 10: functions: fundamental conversions ---------------
01027    // 
01028       // Fundamental conversion from spherical to cartesian coordinates.
01029       // @param trp (input): theta, phi, radius
01030       // @param xyz (output): X,Y,Z in units of radius
01031       // Algorithm references: standard geometry.
01032    void Position::convertSphericalToCartesian(const Triple& tpr,
01033                                               Triple& xyz)
01034       throw()
01035    {
01036       double st=sin(tpr[0]*DEG_TO_RAD);
01037       xyz[0] = tpr[2]*st*cos(tpr[1]*DEG_TO_RAD);
01038       xyz[1] = tpr[2]*st*sin(tpr[1]*DEG_TO_RAD);
01039       xyz[2] = tpr[2]*cos(tpr[0]*DEG_TO_RAD);
01040    }
01041 
01042       // Fundamental routine to convert cartesian to spherical coordinates.
01043       // @param xyz (input): X,Y,Z
01044       // @param trp (output): theta, phi (deg), radius in units of input
01045       // Algorithm references: standard geometry.
01046    void Position::convertCartesianToSpherical(const Triple& xyz,
01047                                               Triple& tpr)
01048       throw()
01049    {
01050       tpr[2] = RSS(xyz[0],xyz[1],xyz[2]);
01051       if(tpr[2] <= Position::POSITION_TOLERANCE/5)
01052       {
01053             // zero-length Cartesian vector
01054          tpr[0] = 90;
01055          tpr[1] = 0;
01056          return;
01057       }
01058       tpr[0] = acos(xyz[2]/tpr[2]);
01059       tpr[0] *= RAD_TO_DEG;
01060       if(RSS(xyz[0],xyz[1]) < Position::POSITION_TOLERANCE/5) {       // pole
01061          tpr[1] = 0;
01062          return;
01063       }
01064       tpr[1] = atan2(xyz[1],xyz[0]);
01065       tpr[1] *= RAD_TO_DEG;
01066       if(tpr[1] < 0) tpr[1] += 360;
01067    }
01068 
01069       // Fundamental routine to convert cartesian (ECEF) to geodetic coordinates,
01070       // (Geoid specified by semi-major axis and eccentricity squared).
01071       // @param xyz (input): X,Y,Z in meters
01072       // @param llh (output): geodetic lat(deg N), lon(deg E),
01073       //                             height above ellipsoid (meters)
01074       // @param A (input) Earth semi-major axis
01075       // @param eccSq (input) square of Earth eccentricity
01076       // Algorithm references: 
01077    void Position::convertCartesianToGeodetic(const Triple& xyz,
01078                                              Triple& llh,
01079                                              const double A,
01080                                              const double eccSq)
01081       throw()
01082    {
01083       double p,slat,N,htold,latold;
01084       p = SQRT(xyz[0]*xyz[0]+xyz[1]*xyz[1]);
01085       if(p < Position::POSITION_TOLERANCE/5) {   // pole or origin
01086          llh[0] = llh[1] = 0;
01087          llh[2] = fabs(xyz[2]) - A;
01088          return;
01089       }
01090       llh[0] = atan2(xyz[2], p*(1.0-eccSq));
01091       llh[2] = 0;
01092       for(int i=0; i<5; i++) {
01093          slat = sin(llh[0]);
01094          N = A / SQRT(1.0 - eccSq*slat*slat);
01095          htold = llh[2];
01096          llh[2] = p/cos(llh[0]) - N;
01097          latold = llh[0];
01098          llh[0] = atan2(xyz[2], p*(1.0-eccSq*(N/(N+llh[2]))));
01099          if(fabs(llh[0]-latold) < 1.0e-9 && fabs(llh[2]-htold) < 1.0e-9 * A) break;
01100       }
01101       llh[1] = atan2(xyz[1],xyz[0]);
01102       if(llh[1] < 0.0) llh[1] += TWO_PI;
01103       llh[0] *= RAD_TO_DEG;
01104       llh[1] *= RAD_TO_DEG;
01105    }
01106 
01107       // Fundamental routine to convert geodetic to cartesian (ECEF) coordinates,
01108       // (Geoid specified by semi-major axis and eccentricity squared).
01109       // @param llh (input): geodetic lat(deg N), lon(deg E),
01110       //            height above ellipsoid (meters)
01111       // @param xyz (output): X,Y,Z in meters
01112       // @param A (input) Earth semi-major axis
01113       // @param eccSq (input) square of Earth eccentricity
01114       // Algorithm references: 
01115    void Position::convertGeodeticToCartesian(const Triple& llh,
01116                                              Triple& xyz,
01117                                              const double A,
01118                                              const double eccSq)
01119       throw()
01120    {
01121       double slat = sin(llh[0]*DEG_TO_RAD);
01122       double clat = cos(llh[0]*DEG_TO_RAD);
01123       double N = A/SQRT(1.0-eccSq*slat*slat);
01124       xyz[0] = (N+llh[2])*clat*cos(llh[1]*DEG_TO_RAD);
01125       xyz[1] = (N+llh[2])*clat*sin(llh[1]*DEG_TO_RAD);
01126       xyz[2] = (N*(1.0-eccSq)+llh[2])*slat;
01127    }
01128 
01129       // Fundamental routine to convert cartesian (ECEF) to geocentric coordinates.
01130       // @param xyz (input): X,Y,Z in meters
01131       // @param llr (output):
01132       //            geocentric lat(deg N),lon(deg E),radius (units of input)
01133    void Position::convertCartesianToGeocentric(const Triple& xyz,
01134                                                Triple& llr)
01135       throw()
01136    {
01137       convertCartesianToSpherical(xyz, llr);
01138       llr[0] = 90 - llr[0];         // convert theta to latitude
01139    }
01140 
01141       // Fundamental routine to convert geocentric to cartesian (ECEF) coordinates.
01142       // @param llr (input): geocentric lat(deg N),lon(deg E),radius
01143       // @param xyz (output): X,Y,Z (units of radius)
01144    void Position::convertGeocentricToCartesian(const Triple& llr,
01145                                                Triple& xyz)
01146       throw()
01147    {
01148       Triple llh(llr);
01149       llh[0] = 90 - llh[0];         // convert latitude to theta
01150       convertSphericalToCartesian(llh, xyz);
01151    }
01152 
01153       // Fundamental routine to convert geocentric to geodetic coordinates.
01154       // @param llr (input): geocentric Triple: lat(deg N),lon(deg E),radius (meters)
01155       // @param llh (output): geodetic latitude (deg N),
01156       //            longitude (deg E), and height above ellipsoid (meters)
01157       // @param A (input) Earth semi-major axis
01158       // @param eccSq (input) square of Earth eccentricity
01159    void Position::convertGeocentricToGeodetic(const Triple& llr,
01160                                                Triple& llh,
01161                                                const double A,
01162                                                const double eccSq)
01163       throw()
01164    {
01165       double cl,p,sl,slat,N,htold,latold;
01166       llh[1] = llr[1]; // longitude is the same in both systems
01167       cl = sin((90-llr[0])*DEG_TO_RAD);
01168       sl = cos((90-llr[0])*DEG_TO_RAD);
01169       if(llr[2] <= Position::POSITION_TOLERANCE/5) {
01170          // radius is below tolerance, hence assign zero-length
01171          // arbitrarily set latitude = longitude = 0
01172          llh[0] = llh[1] = 0;
01173          llh[2] = -A;
01174          return;
01175       }
01176       else if(cl < 1.e-10) {
01177          // near pole ... note that 1mm/radius(Earth) = 1.5e-10
01178          if(llr[0] < 0) llh[0] = -90;
01179          else           llh[0] =  90;
01180          llh[1] = 0;
01181          llh[2] = llr[2] - A*SQRT(1-eccSq);
01182          return;
01183       }
01184       llh[0] = atan2(sl, cl*(1.0-eccSq));
01185       p = cl*llr[2];
01186       llh[2] = 0;
01187       for(int i=0; i<5; i++) {
01188          slat = sin(llh[0]);
01189          N = A / SQRT(1.0 - eccSq*slat*slat);
01190          htold = llh[2];
01191          llh[2] = p/cos(llh[0]) - N;
01192          latold = llh[0];
01193          llh[0] = atan2(sl, cl*(1.0-eccSq*(N/(N+llh[2]))));
01194          if(fabs(llh[0]-latold) < 1.0e-9 && fabs(llh[2]-htold) < 1.0e-9 * A) break;
01195       }
01196       llh[0] *= RAD_TO_DEG;
01197    }
01198 
01199       // Fundamental routine to convert geodetic to geocentric coordinates.
01200       // @param geodeticllh (input): geodetic latitude (deg N),
01201       //            longitude (deg E), and height above ellipsoid (meters)
01202       // @param llr (output): geocentric lat (deg N),lon (deg E),radius (meters)
01203       // @param A (input) Earth semi-major axis
01204       // @param eccSq (input) square of Earth eccentricity
01205    void Position::convertGeodeticToGeocentric(const Triple& llh,
01206                                               Triple& llr,
01207                                               const double A,
01208                                               const double eccSq)
01209       throw()
01210    {
01211       double slat = sin(llh[0]*DEG_TO_RAD);
01212       double N = A/SQRT(1.0-eccSq*slat*slat);
01213       llr[1] = llh[1]; // longitude is the same in both systems
01214       // radius
01215       llr[2] = SQRT((N+llh[2])*(N+llh[2]) + N*eccSq*(N*eccSq-2*(N+llh[2]))*slat*slat);
01216       if(llr[2] <= Position::POSITION_TOLERANCE/5)
01217       {
01218             // radius is below tolerance, hence assign zero-length
01219             // arbitrarily set latitude = longitude = 0
01220          llr[0] = llr[1] = llr[2] = 0;
01221          return;
01222       }
01223       if(1-fabs(slat) < 1.e-10) {             // at the pole
01224          if(slat < 0) llr[0] = -90;
01225          else         llr[0] =  90;
01226          llr[1] = 0;
01227          return;
01228       }
01229       // theta
01230       llr[0] = acos((N*(1-eccSq)+llh[2])*slat/llr[2]);
01231       llr[0] *= RAD_TO_DEG;
01232       llr[0] = 90 - llr[0];
01233    }
01234 
01235    // ----------- Part 11: operator<< and other useful functions -------------
01236    //
01237      // Stream output for Position objects.
01238      // @param s stream to append formatted Position to.
01239      // @param t Position to append to stream \c s.
01240      // @return reference to \c s.
01241    ostream& operator<<(ostream& s, const Position& p)
01242    {
01243       if(p.system == Position::Cartesian)
01244          s << p.printf("%.4x m %.4y m %.4z m");
01245       else if(p.system == Position::Geodetic)
01246          s << p.printf("%.8A degN %.8L degE %.4h m");
01247       else if(p.system == Position::Geocentric)
01248          s << p.printf("%.8a degN %.8L degE %.4r m");
01249       else if(p.system == Position::Spherical)
01250          s << p.printf("%.8t deg %.8p deg %.4r m");
01251       else
01252          s << " Unknown system! : " << p[0] << " " << p[1] << " " << p[2];
01253 
01254       return s;
01255    }
01256 
01257       // Compute the range in meters between this Position and
01258       // the Position passed as input.
01259       // @param right Position to which to find the range
01260       // @return the range (in meters)
01261       // @throw GeometryException if geoid values differ
01262    double range(const Position& A,
01263                 const Position& B)
01264       throw(GeometryException)
01265    {
01266       if(A.AEarth != B.AEarth || A.eccSquared != B.eccSquared)
01267       {
01268          GeometryException ge("Unequal geoids");
01269          GPSTK_THROW(ge);
01270       }
01271 
01272          Position L(A),R(B);
01273          L.transformTo(Position::Cartesian);
01274          R.transformTo(Position::Cartesian);
01275          double dif = RSS(L.X()-R.X(),L.Y()-R.Y(),L.Z()-R.Z());
01276          return dif;
01277       }
01278 
01279      // Compute the radius of the ellipsoidal Earth, given the geodetic latitude.
01280      // @param geolat geodetic latitude in degrees
01281      // @return the Earth radius (in meters)
01282    double Position::radiusEarth(const double geolat,
01283                                 const double A,
01284                                 const double eccSq)
01285       throw()
01286    {
01287       double slat=sin(DEG_TO_RAD*geolat);
01288       double e=(1.0-eccSq);
01289       double f=(1.0+(e*e-1.0)*slat*slat)/(1.0-eccSq*slat*slat);
01290       return (A * SQRT(f));
01291    }
01292 
01293       // A member function that computes the elevation of the input
01294       // (Target) position as seen from this Position.
01295       // @param Target the Position which is observed to have the
01296       //        computed elevation, as seen from this Position.
01297       // @return the elevation in degrees
01298    double Position::elevation(const Position& Target) const
01299       throw(GeometryException)
01300    {
01301       Position R(*this),S(Target);
01302       R.transformTo(Cartesian);
01303       S.transformTo(Cartesian);
01304       // use Triple:: functions in cartesian coordinates (only)
01305       double elevation;
01306       try {
01307          elevation = R.elvAngle(S);         
01308       }
01309       catch(GeometryException& ge)
01310       {
01311          GPSTK_RETHROW(ge);
01312       }
01313       return elevation;
01314    }
01315 
01316       // A member function that computes the elevation of the input
01317       // (Target) position as seen from this Position, using a Geodetic
01318       // (i.e. ellipsoidal) system.
01319       // @param Target the Position which is observed to have the
01320       //        computed elevation, as seen from this Position.
01321       // @return the elevation in degrees
01322    double Position::elevationGeodetic(const Position& Target) const
01323       throw(GeometryException)
01324    {
01325       Position R(*this),S(Target);
01326       double latGeodetic = R.getGeodeticLatitude()*DEG_TO_RAD;
01327       double longGeodetic = R.getLongitude()*DEG_TO_RAD;
01328       double localUp;
01329       double cosUp;
01330       R.transformTo(Cartesian);
01331       S.transformTo(Cartesian);
01332       Triple z;
01333       // Let's get the slant vector
01334       z = S.theArray - R.theArray;
01335 
01336       if (z.mag()<=1e-4) // if the positions are within .1 millimeter
01337       {
01338          GeometryException ge("Positions are within .1 millimeter");
01339          GPSTK_THROW(ge);
01340       }
01341 
01342       // Compute k vector in local North-East-Up (NEU) system
01343       Triple kVector(cos(latGeodetic)*cos(longGeodetic), cos(latGeodetic)*sin(longGeodetic), sin(latGeodetic));
01344       // Take advantage of dot method to get Up coordinate in local NEU system
01345       localUp = z.dot(kVector);
01346       // Let's get cos(z), being z the angle with respect to local vertical (Up);
01347       cosUp = localUp/z.mag();
01348 
01349       return 90.0 - ((::acos(cosUp))*RAD_TO_DEG);
01350    }
01351 
01352       // A member function that computes the azimuth of the input
01353       // (Target) position as seen from this Position.
01354       // @param Target the Position which is observed to have the
01355       //        computed azimuth, as seen from this Position.
01356       // @return the azimuth in degrees
01357    double Position::azimuth(const Position& Target) const
01358       throw(GeometryException)
01359    {
01360       Position R(*this),S(Target);
01361       R.transformTo(Cartesian);
01362       S.transformTo(Cartesian);
01363       // use Triple:: functions in cartesian coordinates (only)
01364       double az;
01365       try
01366       {
01367          az = R.azAngle(S);
01368          
01369       }
01370       catch(GeometryException& ge)
01371       {
01372          GPSTK_RETHROW(ge);
01373       }
01374       
01375       return az; 
01376    }
01377 
01378       // A member function that computes the azimuth of the input
01379       // (Target) position as seen from this Position, using a Geodetic
01380       // (i.e. ellipsoidal) system.
01381       // @param Target the Position which is observed to have the
01382       //        computed azimuth, as seen from this Position.
01383       // @return the azimuth in degrees
01384    double Position::azimuthGeodetic(const Position& Target) const
01385       throw(GeometryException)
01386    {
01387       Position R(*this),S(Target);
01388       double latGeodetic = R.getGeodeticLatitude()*DEG_TO_RAD;
01389       double longGeodetic = R.getLongitude()*DEG_TO_RAD;
01390       double localN, localE;
01391       R.transformTo(Cartesian);
01392       S.transformTo(Cartesian);
01393       Triple z;
01394       // Let's get the slant vector
01395       z = S.theArray - R.theArray;
01396 
01397       if (z.mag()<=1e-4) // if the positions are within .1 millimeter
01398       {
01399          GeometryException ge("Positions are within .1 millimeter");
01400          GPSTK_THROW(ge);
01401       }
01402       
01403       // Compute i vector in local North-East-Up (NEU) system
01404       Triple iVector(-sin(latGeodetic)*cos(longGeodetic), -sin(latGeodetic)*sin(longGeodetic), cos(latGeodetic));
01405       // Compute j vector in local North-East-Up (NEU) system
01406       Triple jVector(-sin(longGeodetic), cos(longGeodetic), 0);
01407 
01408       // Now, let's use dot product to get localN and localE unitary vectors
01409       localN = (z.dot(iVector))/z.mag();
01410       localE = (z.dot(jVector))/z.mag();
01411 
01412       // Let's test if computing azimuth has any sense
01413       double test = fabs(localN) + fabs(localE);
01414 
01415       // Warning: If elevation is very close to 90 degrees, we will return azimuth = 0.0
01416       if (test < 1.0e-16) return 0.0;
01417 
01418       double alpha = ((::atan2(localE, localN)) * RAD_TO_DEG);
01419       if (alpha < 0.0)
01420       {
01421          return alpha + 360.0;
01422       }
01423       else 
01424       {
01425          return alpha;
01426       }
01427    }
01428 
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      // A member function that computes the point at which a signal, which
01435      // is received at *this Position and there is observed at the input
01436      // azimuth and elevation, crosses a model ionosphere that is taken to
01437      // be a uniform thin shell at the input height. This algorithm is done
01438      // in geocentric coordinates.
01439      // @param elev elevation angle of the signal at reception, in degrees
01440      // @param azim azimuth angle of the signal at reception, in degrees
01441      // @param ionoht height of the ionosphere, in meters
01442      // @return Position IPP the position of the ionospheric pierce point,
01443      //     in the same coordinate system as *this; *this is not modified.
01444    Position Position::getIonosphericPiercePoint(const double elev,
01445                                                 const double azim,
01446                                                 const double ionoht) const
01447       throw()
01448    {
01449       Position Rx(*this);
01450 
01451       // convert to Geocentric
01452       Rx.transformTo(Geocentric);
01453 
01454       // compute the geographic pierce point
01455       Position IPP(Rx);                   // copy system and geoid
01456       double el = elev * DEG_TO_RAD;
01457       // p is the angle subtended at Earth center by Rx and the IPP
01458       double p = PI/2.0 - el - asin(AEarth*cos(el)/(AEarth+ionoht));
01459       double lat = Rx.theArray[0] * DEG_TO_RAD;
01460       double az = azim * DEG_TO_RAD;
01461       IPP.theArray[0] = asin(sin(lat)*cos(p) + cos(lat)*sin(p)*cos(az));
01462       IPP.theArray[1] = Rx.theArray[1]*DEG_TO_RAD
01463          + asin(sin(p)*sin(az)/cos(IPP.theArray[0]));
01464 
01465       IPP.theArray[0] *= RAD_TO_DEG;
01466       IPP.theArray[1] *= RAD_TO_DEG;
01467       IPP.theArray[2] = AEarth + ionoht;
01468 
01469       // transform back
01470       IPP.transformTo(system);
01471 
01472       return IPP;
01473    }
01474 
01475    // ----------- Part 12: private functions and member data -----------------
01476    //
01477       // Initialization function, used by the constructors.
01478       // @param a coordinate [ X(m), or latitude (degrees N) ]
01479       // @param b coordinate [ Y(m), or longitude (degrees E) ]
01480       // @param c coordinate [ Z, height above ellipsoid or radius, in m ]
01481       // @param s CoordinateSystem, defaults to Cartesian
01482       // @param geiod pointer to a GeoidModel, default NULL (WGS84)
01483       // @throw GeometryException on invalid input.
01484    void Position::initialize(const double a,
01485                   const double b,
01486                   const double c,
01487                   Position::CoordinateSystem s,
01488                   GeoidModel *geoid)
01489       throw(GeometryException)
01490    {
01491       double bb(b);
01492       if(s == Geodetic || s==Geocentric)
01493       {
01494          if(a > 90 || a < -90)
01495          {
01496             GeometryException ge("Invalid latitude in constructor: "
01497                                     + StringUtils::asString(a));
01498             GPSTK_THROW(ge);
01499          }
01500          if(bb < 0)
01501             bb += 360*(1+(unsigned long)(bb/360));
01502          else if(bb >= 360)
01503             bb -= 360*(unsigned long)(bb/360);
01504       }
01505       if(s==Geocentric || s==Spherical)
01506       {
01507          if(c < 0)
01508          {
01509             GeometryException ge("Invalid radius in constructor: "
01510                                            + StringUtils::asString(c));
01511             GPSTK_THROW(ge);
01512          }
01513       }
01514       if(s==Spherical)
01515       {
01516          if(a < 0 || a > 180)
01517          {
01518             GeometryException ge("Invalid theta in constructor: "
01519                                     + StringUtils::asString(a));
01520             GPSTK_THROW(ge);
01521          }
01522          if(bb < 0)
01523             bb += 360*(1+(unsigned long)(bb/360));
01524          else if(bb >= 360)
01525             bb -= 360*(unsigned long)(bb/360);
01526       }
01527 
01528       theArray[0] = a;
01529       theArray[1] = bb;
01530       theArray[2] = c;
01531 
01532       if(geoid) {
01533          AEarth = geoid->a();
01534          eccSquared = geoid->eccSquared();
01535       }
01536       else {
01537          WGS84Geoid WGS84;
01538          AEarth = WGS84.a();
01539          eccSquared = WGS84.eccSquared();
01540       }
01541       system = s;
01542       tolerance = POSITION_TOLERANCE;
01543    }
01544 
01545 }  // namespace gpstk

Generated on Wed Mar 10 03:30:49 2010 for GPS ToolKit Software Library by  doxygen 1.3.9.1