Changeset 2221 in ntrip for trunk/BNC/RTCM3
- Timestamp:
- Jan 12, 2010, 8:36:29 AM (15 years ago)
- Location:
- trunk/BNC/RTCM3
- Files:
-
- 2 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/BNC/RTCM3/ephemeris.cpp
r1239 r2221 5 5 6 6 #include "ephemeris.h" 7 #include "bncutils.h" 7 8 #include "timeutils.h" 8 9 … … 270 271 out << endl; 271 272 } 273 274 // Derivative of the state vector using a simple force model (static) 275 //////////////////////////////////////////////////////////////////////////// 276 ColumnVector t_ephGlo::glo_deriv(double /* tt */, const ColumnVector& xv) { 277 278 // State vector components 279 // ----------------------- 280 ColumnVector rr = xv.rows(1,3); 281 ColumnVector vv = xv.rows(4,6); 282 283 // Acceleration 284 // ------------ 285 static const double GM = 398.60044e12; 286 static const double AE = 6378136.0; 287 static const double OMEGA = 7292115.e-11; 288 static const double C20 = -1082.63e-6; 289 290 double rho = rr.norm_Frobenius(); 291 double t1 = -GM/(rho*rho*rho); 292 double t2 = 3.0/2.0 * C20 * (GM*AE*AE) / (rho*rho*rho*rho*rho); 293 double t3 = OMEGA * OMEGA; 294 double t4 = 2.0 * OMEGA; 295 double z2 = rr(3) * rr(3); 296 297 // Vector of derivatives 298 // --------------------- 299 ColumnVector va(6); 300 va(1) = vv(1); 301 va(2) = vv(2); 302 va(3) = vv(3); 303 va(4) = (t1 + t2*(1.0-5.0*z2/(rho*rho)) + t3) * rr(1) + t4*vv(2); 304 va(5) = (t1 + t2*(1.0-5.0*z2/(rho*rho)) + t3) * rr(2) - t4*vv(1); 305 va(6) = (t1 + t2*(3.0-5.0*z2/(rho*rho)) ) * rr(3); 306 307 return va; 308 } 309 310 // Compute Glonass Satellite Position (virtual) 311 //////////////////////////////////////////////////////////////////////////// 312 void t_ephGlo::position(int GPSweek, double GPSweeks, 313 double* xc, double* vv) const { 314 315 static const double secPerWeek = 7 * 86400.0; 316 static const double nominalStep = 10.0; 317 318 memset(xc, 0, 4*sizeof(double)); 319 memset(vv, 0, 3*sizeof(double)); 320 321 double dtPos = GPSweeks - _tt; 322 if (GPSweek != _GPSweek) { 323 dtPos += (GPSweek - _GPSweek) * secPerWeek; 324 } 325 326 int nSteps = int(fabs(dtPos) / nominalStep) + 1; 327 double step = dtPos / nSteps; 328 329 for (int ii = 1; ii <= nSteps; ii++) { 330 _xv = rungeKutta4(_tt, _xv, step, glo_deriv); 331 _tt += step; 332 } 333 334 // Position and Velocity 335 // --------------------- 336 xc[0] = _xv(1); 337 xc[1] = _xv(2); 338 xc[2] = _xv(3); 339 340 vv[0] = _xv(4); 341 vv[1] = _xv(5); 342 vv[2] = _xv(6); 343 344 // Clock Correction 345 // ---------------- 346 double dtClk = GPSweeks - _GPSweeks; 347 if (GPSweek != _GPSweek) { 348 dtClk += (GPSweek - _GPSweek) * secPerWeek; 349 } 350 xc[3] = -_tau + _gamma * dtClk; 351 } 352 353 // IOD of Glonass Ephemeris (virtual) 354 //////////////////////////////////////////////////////////////////////////// 355 int t_ephGlo::IOD() const { 356 357 bool old = false; 358 359 if (old) { // 5 LSBs of iod are equal to 5 LSBs of tb 360 unsigned int tb = int(fmod(_GPSweeks,86400.0)); //sec of day 361 const int shift = sizeof(tb) * 8 - 5; 362 unsigned int iod = tb << shift; 363 return (iod >> shift); 364 } 365 else { 366 return int(fmod(_tki, 3600)) / 30; 367 } 368 } 369 370 // Print Glonass Ephemeris (virtual) 371 //////////////////////////////////////////////////////////////////////////// 372 void t_ephGlo::print(std::ostream& out) const { 373 374 } 375 376 // Set Glonass Ephemeris 377 //////////////////////////////////////////////////////////////////////////// 378 void t_ephGlo::set(const glonassephemeris* ee) { 379 380 } -
trunk/BNC/RTCM3/ephemeris.h
r2040 r2221 1 1 #ifndef EPHEMERIS_H 2 2 #define EPHEMERIS_H 3 4 #include <newmat.h> 3 5 4 6 #include <stdio.h> … … 19 21 20 22 virtual void position(int GPSweek, double GPSweeks, 21 22 23 double* xc, 24 double* vv) const = 0; 23 25 24 26 void position(int GPSweek, double GPSweeks, 25 27 double& xx, double& yy, double& zz, double& cc) const { 26 28 double tmp_xx[4]; 27 29 double tmp_vv[4]; … … 49 51 public: 50 52 t_ephGPS() { } 51 ~t_ephGPS() {}53 virtual ~t_ephGPS() {} 52 54 double TOC() const {return _TOC;} 53 55 … … 55 57 56 58 void set(int prn, 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 59 int GPSWeek, 60 double toc, double toe, double tot, 61 double IODE, double IODC, 62 double clock_bias, double clock_drift, double clock_driftrate, 63 double OMEGA0, double OMEGADOT, 64 double i0, double IDOT, 65 double omega, 66 double M0, double Delta_n, 67 double sqrt_A, 68 double e, 69 double Crc, double Crs, 70 double Cic, double Cis, 71 double Cuc, double Cus, 72 double TGD, 73 int health); 72 74 73 v oid position(int GPSweek, double GPSweeks,74 75 75 virtual void position(int GPSweek, double GPSweeks, 76 double* xc, 77 double* vv) const; 76 78 77 int IOD() const { return static_cast<int>(_IODC); }79 virtual int IOD() const { return static_cast<int>(_IODC); } 78 80 79 v oid print(std::ostream& out) const;81 virtual void print(std::ostream& out) const; 80 82 81 83 private: … … 109 111 }; 110 112 113 class t_ephGlo : public t_eph { 114 public: 115 t_ephGlo() { _gps_utc = 0.0; _xv.ReSize(6); } 116 117 virtual ~t_ephGlo() {} 118 119 virtual void position(int GPSweek, double GPSweeks, 120 double* xc, 121 double* vv) const; 122 123 virtual int IOD() const; 124 125 virtual void print(std::ostream& out) const; 126 127 void set(const glonassephemeris* ee); 128 129 private: 130 static ColumnVector glo_deriv(double /* tt */, const ColumnVector& xv); 131 132 mutable double _tt; // time in seconds of GPSweek 133 mutable ColumnVector _xv; // status vector (position, velocity) at time _tt 134 135 double _gps_utc; // GPS - UTC in seconds 136 double _E; // [days] 137 double _tau; // [s] 138 double _gamma; // 139 double _x_pos; // [km] 140 double _x_velocity; // [km/s] 141 double _x_acceleration; // [km/s^2] 142 double _y_pos; // [km] 143 double _y_velocity; // [km/s] 144 double _y_acceleration; // [km/s^2] 145 double _z_pos; // [km] 146 double _z_velocity; // [km/s] 147 double _z_acceleration; // [km/s^2] 148 double _health; // 0 = O.K. 149 double _frequency_number; // ICD-GLONASS data position 150 double _tki; // message frame time 151 }; 152 111 153 #endif
Note:
See TracChangeset
for help on using the changeset viewer.