source: ntrip/trunk/BNC/RTCM3/ephemeris.cpp@ 2257

Last change on this file since 2257 was 2257, checked in by mervart, 14 years ago

* empty log message *

File size: 8.0 KB
Line 
1#include <math.h>
2#include <sstream>
3#include <iostream>
4#include <iomanip>
5#include <cstring>
6
7#include <newmatio.h>
8
9#include "ephemeris.h"
10#include "bncutils.h"
11#include "timeutils.h"
12
13using namespace std;
14
15//
16////////////////////////////////////////////////////////////////////////////
17bool t_eph::isNewerThan(const t_eph* eph) const {
18 if (_GPSweek > eph->_GPSweek ||
19 (_GPSweek == eph->_GPSweek && _GPSweeks > eph->_GPSweeks)) {
20 return true;
21 }
22 else {
23 return false;
24 }
25}
26
27// Set GPS Satellite Position
28////////////////////////////////////////////////////////////////////////////
29void t_ephGPS::set(const gpsephemeris* ee) {
30 ostringstream prn;
31 prn << 'G' << setfill('0') << setw(2) << ee->satellite;
32
33 _prn = prn.str();
34
35 // TODO: check if following two lines are correct
36 _GPSweek = ee->GPSweek;
37 _GPSweeks = ee->TOE;
38
39 _TOW = ee->TOW;
40 _TOC = ee->TOC;
41 _TOE = ee->TOE;
42 _IODE = ee->IODE;
43 _IODC = ee->IODC;
44
45 _clock_bias = ee->clock_bias ;
46 _clock_drift = ee->clock_drift ;
47 _clock_driftrate = ee->clock_driftrate;
48
49 _Crs = ee->Crs;
50 _Delta_n = ee->Delta_n;
51 _M0 = ee->M0;
52 _Cuc = ee->Cuc;
53 _e = ee->e;
54 _Cus = ee->Cus;
55 _sqrt_A = ee->sqrt_A;
56 _Cic = ee->Cic;
57 _OMEGA0 = ee->OMEGA0;
58 _Cis = ee->Cis;
59 _i0 = ee->i0;
60 _Crc = ee->Crc;
61 _omega = ee->omega;
62 _OMEGADOT = ee->OMEGADOT;
63 _IDOT = ee->IDOT;
64
65 _TGD = ee->TGD;
66}
67
68// Compute GPS Satellite Position (virtual)
69////////////////////////////////////////////////////////////////////////////
70void t_ephGPS::position(int GPSweek, double GPSweeks,
71 double* xc,
72 double* vv) const {
73
74 static const double secPerWeek = 7 * 86400.0;
75 static const double omegaEarth = 7292115.1467e-11;
76 static const double gmWGS = 398.6005e12;
77
78 memset(xc, 0, 4*sizeof(double));
79 memset(vv, 0, 3*sizeof(double));
80
81 double a0 = _sqrt_A * _sqrt_A;
82 if (a0 == 0) {
83 return;
84 }
85
86 double n0 = sqrt(gmWGS/(a0*a0*a0));
87 double tk = GPSweeks - _TOE;
88 if (GPSweek != _GPSweek) {
89 tk += (GPSweek - _GPSweek) * secPerWeek;
90 }
91 double n = n0 + _Delta_n;
92 double M = _M0 + n*tk;
93 double E = M;
94 double E_last;
95 do {
96 E_last = E;
97 E = M + _e*sin(E);
98 } while ( fabs(E-E_last)*a0 > 0.001 );
99 double v = 2.0*atan( sqrt( (1.0 + _e)/(1.0 - _e) )*tan( E/2 ) );
100 double u0 = v + _omega;
101 double sin2u0 = sin(2*u0);
102 double cos2u0 = cos(2*u0);
103 double r = a0*(1 - _e*cos(E)) + _Crc*cos2u0 + _Crs*sin2u0;
104 double i = _i0 + _IDOT*tk + _Cic*cos2u0 + _Cis*sin2u0;
105 double u = u0 + _Cuc*cos2u0 + _Cus*sin2u0;
106 double xp = r*cos(u);
107 double yp = r*sin(u);
108 double OM = _OMEGA0 + (_OMEGADOT - omegaEarth)*tk -
109 omegaEarth*_TOE;
110
111 double sinom = sin(OM);
112 double cosom = cos(OM);
113 double sini = sin(i);
114 double cosi = cos(i);
115 xc[0] = xp*cosom - yp*cosi*sinom;
116 xc[1] = xp*sinom + yp*cosi*cosom;
117 xc[2] = yp*sini;
118
119 double tc = GPSweeks - _TOC;
120 if (GPSweek != _GPSweek) {
121 tc += (GPSweek - _GPSweek) * secPerWeek;
122 }
123 xc[3] = _clock_bias + _clock_drift*tc + _clock_driftrate*tc*tc
124 - 4.442807633e-10 * _e * sqrt(a0) *sin(E);
125
126 // Velocity
127 // --------
128 double tanv2 = tan(v/2);
129 double dEdM = 1 / (1 - _e*cos(E));
130 double dotv = sqrt((1.0 + _e)/(1.0 - _e)) / cos(E/2)/cos(E/2) / (1 + tanv2*tanv2)
131 * dEdM * n;
132 double dotu = dotv + (-_Cuc*sin2u0 + _Cus*cos2u0)*2*dotv;
133 double dotom = _OMEGADOT - omegaEarth;
134 double doti = _IDOT + (-_Cic*sin2u0 + _Cis*cos2u0)*2*dotv;
135 double dotr = a0 * _e*sin(E) * dEdM * n
136 + (-_Crc*sin2u0 + _Crs*cos2u0)*2*dotv;
137 double dotx = dotr*cos(u) - r*sin(u)*dotu;
138 double doty = dotr*sin(u) + r*cos(u)*dotu;
139
140 vv[0] = cosom *dotx - cosi*sinom *doty // dX / dr
141 - xp*sinom*dotom - yp*cosi*cosom*dotom // dX / dOMEGA
142 + yp*sini*sinom*doti; // dX / di
143
144 vv[1] = sinom *dotx + cosi*cosom *doty
145 + xp*cosom*dotom - yp*cosi*sinom*dotom
146 - yp*sini*cosom*doti;
147
148 vv[2] = sini *doty + yp*cosi *doti;
149}
150
151
152// Derivative of the state vector using a simple force model (static)
153////////////////////////////////////////////////////////////////////////////
154ColumnVector t_ephGlo::glo_deriv(double /* tt */, const ColumnVector& xv) {
155
156 // State vector components
157 // -----------------------
158 ColumnVector rr = xv.rows(1,3);
159 ColumnVector vv = xv.rows(4,6);
160
161 // Acceleration
162 // ------------
163 static const double GM = 398.60044e12;
164 static const double AE = 6378136.0;
165 static const double OMEGA = 7292115.e-11;
166 static const double C20 = -1082.63e-6;
167
168 double rho = rr.norm_Frobenius();
169 double t1 = -GM/(rho*rho*rho);
170 double t2 = 3.0/2.0 * C20 * (GM*AE*AE) / (rho*rho*rho*rho*rho);
171 double t3 = OMEGA * OMEGA;
172 double t4 = 2.0 * OMEGA;
173 double z2 = rr(3) * rr(3);
174
175 // Vector of derivatives
176 // ---------------------
177 ColumnVector va(6);
178 va(1) = vv(1);
179 va(2) = vv(2);
180 va(3) = vv(3);
181 va(4) = (t1 + t2*(1.0-5.0*z2/(rho*rho)) + t3) * rr(1) + t4*vv(2);
182 va(5) = (t1 + t2*(1.0-5.0*z2/(rho*rho)) + t3) * rr(2) - t4*vv(1);
183 va(6) = (t1 + t2*(3.0-5.0*z2/(rho*rho)) ) * rr(3);
184
185 return va;
186}
187
188// Compute Glonass Satellite Position (virtual)
189////////////////////////////////////////////////////////////////////////////
190void t_ephGlo::position(int GPSweek, double GPSweeks,
191 double* xc, double* vv) const {
192
193 static const double secPerWeek = 7 * 86400.0;
194 static const double nominalStep = 10.0;
195
196 memset(xc, 0, 4*sizeof(double));
197 memset(vv, 0, 3*sizeof(double));
198
199 double dtPos = GPSweeks - _tt;
200 if (GPSweek != _GPSweek) {
201 dtPos += (GPSweek - _GPSweek) * secPerWeek;
202 }
203
204 int nSteps = int(fabs(dtPos) / nominalStep) + 1;
205 double step = dtPos / nSteps;
206
207 for (int ii = 1; ii <= nSteps; ii++) {
208 _xv = rungeKutta4(_tt, _xv, step, glo_deriv);
209 _tt += step;
210 }
211
212 // Position and Velocity
213 // ---------------------
214 xc[0] = _xv(1);
215 xc[1] = _xv(2);
216 xc[2] = _xv(3);
217
218 vv[0] = _xv(4);
219 vv[1] = _xv(5);
220 vv[2] = _xv(6);
221
222 // Clock Correction
223 // ----------------
224 double dtClk = GPSweeks - _GPSweeks;
225 if (GPSweek != _GPSweek) {
226 dtClk += (GPSweek - _GPSweek) * secPerWeek;
227 }
228 xc[3] = -_tau + _gamma * dtClk;
229}
230
231// IOD of Glonass Ephemeris (virtual)
232////////////////////////////////////////////////////////////////////////////
233int t_ephGlo::IOD() const {
234
235 bool old = false;
236
237 if (old) { // 5 LSBs of iod are equal to 5 LSBs of tb
238 unsigned int tb = int(fmod(_GPSweeks,86400.0)); //sec of day
239 const int shift = sizeof(tb) * 8 - 5;
240 unsigned int iod = tb << shift;
241 return (iod >> shift);
242 }
243 else {
244 return int(fmod(_tki, 3600)) / 30;
245 }
246}
247
248// Set Glonass Ephemeris
249////////////////////////////////////////////////////////////////////////////
250void t_ephGlo::set(const glonassephemeris* ee) {
251
252 ostringstream prn;
253 prn << 'R' << setfill('0') << setw(2) << ee->almanac_number;
254 _prn = prn.str();
255
256 int ww = ee->GPSWeek;
257 int tow = ee->GPSTOW;
258 updatetime(&ww, &tow, ee->tb*1000, 0); // Moscow -> GPS
259
260 _GPSweek = ww;
261 _GPSweeks = tow;
262 _E = ee->E;
263 _tau = ee->tau;
264 _gamma = ee->gamma;
265 _x_pos = ee->x_pos;
266 _x_velocity = ee->x_velocity;
267 _x_acceleration = ee->x_acceleration;
268 _y_pos = ee->y_pos;
269 _y_velocity = ee->y_velocity;
270 _y_acceleration = ee->y_acceleration;
271 _z_pos = ee->z_pos;
272 _z_velocity = ee->z_velocity;
273 _z_acceleration = ee->z_acceleration;
274 _health = 0;
275 _frequency_number = ee->frequency_number;
276 _tki = ee->tk-3*60*60; if (_tki < 0) _tki += 86400;
277
278 // Initialize status vector
279 // ------------------------
280 _tt = _GPSweeks;
281
282 _xv(1) = _x_pos * 1.e3;
283 _xv(2) = _y_pos * 1.e3;
284 _xv(3) = _z_pos * 1.e3;
285 _xv(4) = _x_velocity * 1.e3;
286 _xv(5) = _y_velocity * 1.e3;
287 _xv(6) = _z_velocity * 1.e3;
288}
Note: See TracBrowser for help on using the repository browser.