source: ntrip/trunk/BNS/bnseph.cpp@ 955

Last change on this file since 955 was 932, checked in by weber, 17 years ago

* empty log message *

File size: 10.5 KB
Line 
1/* -------------------------------------------------------------------------
2 * BKG NTRIP Server
3 * -------------------------------------------------------------------------
4 *
5 * Class: bnseph
6 *
7 * Purpose: Retrieve broadcast ephemeris from BNC
8 *
9 * Author: L. Mervart
10 *
11 * Created: 29-Mar-2008
12 *
13 * Changes:
14 *
15 * -----------------------------------------------------------------------*/
16
17#include <iostream>
18
19#include "bnseph.h"
20#include "bnsutils.h"
21
22using namespace std;
23
24// Constructor
25////////////////////////////////////////////////////////////////////////////
26t_bnseph::t_bnseph(QObject* parent) : QThread(parent) {
27
28 this->setTerminationEnabled(true);
29
30 _socket = 0;
31}
32
33// Destructor
34////////////////////////////////////////////////////////////////////////////
35t_bnseph::~t_bnseph() {
36 delete _socket;
37}
38
39// Connect the Socket
40////////////////////////////////////////////////////////////////////////////
41void t_bnseph::reconnect() {
42
43 delete _socket;
44
45 QSettings settings;
46 QString host = "localhost";
47 int port = settings.value("ephPort").toInt();
48
49 _socket = new QTcpSocket();
50 _socket->connectToHost(host, port);
51
52 const int timeOut = 10*1000; // 10 seconds
53 if (!_socket->waitForConnected(timeOut)) {
54 emit(newMessage("bnseph::run Connect Timeout"));
55 }
56}
57
58// Start
59////////////////////////////////////////////////////////////////////////////
60void t_bnseph::run() {
61
62 emit(newMessage("bnseph::run Start"));
63
64 while (true) {
65 if (_socket == 0 || _socket->state() != QAbstractSocket::ConnectedState) {
66 reconnect();
67 }
68 if (_socket && _socket->state() == QAbstractSocket::ConnectedState) {
69 readEph();
70 }
71 else {
72 msleep(10);
73 }
74 }
75}
76
77// Read One Ephemeris
78////////////////////////////////////////////////////////////////////////////
79void t_bnseph::readEph() {
80
81 t_eph* eph = 0;
82 QByteArray line = waitForLine(_socket);
83
84 QTextStream in(line);
85 QString prn;
86
87 in >> prn;
88
89 int numlines = 0;
90 if (prn.indexOf('R') != -1) {
91 eph = new t_ephGlo();
92 numlines = 4;
93 }
94 else {
95 eph = new t_ephGPS();
96 numlines = 8;
97 }
98
99 QStringList lines;
100 lines << line;
101
102 for (int ii = 2; ii <= numlines; ii++) {
103 QByteArray line = waitForLine(_socket);
104 lines << line;
105 }
106
107 eph->read(lines);
108
109 emit(newEph(eph));
110}
111
112// Compare Time
113////////////////////////////////////////////////////////////////////////////
114bool t_eph::isNewerThan(const t_eph* eph) const {
115 if (_GPSweek > eph->_GPSweek ||
116 (_GPSweek == eph->_GPSweek && _GPSweeks > eph->_GPSweeks)) {
117 return true;
118 }
119 else {
120 return false;
121 }
122}
123
124// Read GPS Ephemeris
125////////////////////////////////////////////////////////////////////////////
126void t_ephGPS::read(const QStringList& lines) {
127
128 for (int ii = 1; ii <= lines.size(); ii++) {
129 QTextStream in(lines.at(ii-1).toAscii());
130
131 if (ii == 1) {
132 double year, month, day, hour, minute, second;
133 in >> _prn >> year >> month >> day >> hour >> minute >> second
134 >> _clock_bias >> _clock_drift >> _clock_driftrate;
135
136 if (year < 100) year += 2000;
137
138 QDateTime dateTime(QDate(int(year), int(month), int(day)),
139 QTime(int(hour), int(minute), int(second)), Qt::UTC);
140
141 GPSweekFromDateAndTime(dateTime, _GPSweek, _GPSweeks);
142 _TOC = _GPSweeks;
143 }
144 else if (ii == 2) {
145 in >> _IODE >> _Crs >> _Delta_n >> _M0;
146 }
147 else if (ii == 3) {
148 in >> _Cuc >> _e >> _Cus >> _sqrt_A;
149 }
150 else if (ii == 4) {
151 in >> _TOE >> _Cic >> _OMEGA0 >> _Cis;
152 }
153 else if (ii == 5) {
154 in >> _i0 >> _Crc >> _omega >> _OMEGADOT;
155 }
156 else if (ii == 6) {
157 in >> _IDOT;
158 }
159 else if (ii == 7) {
160 double hlp, health;
161 in >> hlp >> health >> _TGD >> _IODC;
162 }
163 else if (ii == 8) {
164 in >> _TOW;
165 }
166 }
167}
168
169// Compute GPS Satellite Position
170////////////////////////////////////////////////////////////////////////////
171void t_ephGPS::position(int GPSweek, double GPSweeks, ColumnVector& xc,
172 ColumnVector& vv) const {
173
174 const static double secPerWeek = 7 * 86400.0;
175 const static double omegaEarth = 7292115.1467e-11;
176 const static double gmWGS = 398.6005e12;
177
178 if (xc.Nrows() < 4) {
179 xc.ReSize(4);
180 }
181 xc = 0.0;
182
183 if (vv.Nrows() < 3) {
184 vv.ReSize(3);
185 }
186 vv = 0.0;
187
188 double a0 = _sqrt_A * _sqrt_A;
189 if (a0 == 0) {
190 return;
191 }
192
193 double n0 = sqrt(gmWGS/(a0*a0*a0));
194 double tk = GPSweeks - _TOE;
195 if (GPSweek != _GPSweek) {
196 tk += (GPSweek - _GPSweek) * secPerWeek;
197 }
198 double n = n0 + _Delta_n;
199 double M = _M0 + n*tk;
200 double E = M;
201 double E_last;
202 do {
203 E_last = E;
204 E = M + _e*sin(E);
205 } while ( fabs(E-E_last)*a0 > 0.001 );
206 double v = 2.0*atan( sqrt( (1.0 + _e)/(1.0 - _e) )*tan( E/2 ) );
207 double u0 = v + _omega;
208 double sin2u0 = sin(2*u0);
209 double cos2u0 = cos(2*u0);
210 double r = a0*(1 - _e*cos(E)) + _Crc*cos2u0 + _Crs*sin2u0;
211 double i = _i0 + _IDOT*tk + _Cic*cos2u0 + _Cis*sin2u0;
212 double u = u0 + _Cuc*cos2u0 + _Cus*sin2u0;
213 double xp = r*cos(u);
214 double yp = r*sin(u);
215 double OM = _OMEGA0 + (_OMEGADOT - omegaEarth)*tk -
216 omegaEarth*_TOE;
217
218 double sinom = sin(OM);
219 double cosom = cos(OM);
220 double sini = sin(i);
221 double cosi = cos(i);
222 xc(1) = xp*cosom - yp*cosi*sinom;
223 xc(2) = xp*sinom + yp*cosi*cosom;
224 xc(3) = yp*sini;
225
226 double tc = GPSweeks - _TOC;
227 if (GPSweek != _GPSweek) {
228 tc += (GPSweek - _GPSweek) * secPerWeek;
229 }
230 xc(4) = _clock_bias + _clock_drift*tc + _clock_driftrate*tc*tc
231 - 4.442807633e-10 * _e * sqrt(a0) *sin(E);
232
233 // Velocity
234 // --------
235 double tanv2 = tan(v/2);
236 double dEdM = 1 / (1 - _e*cos(E));
237 double dotv = sqrt((1.0 + _e)/(1.0 - _e)) / cos(E/2)/cos(E/2) / (1 + tanv2*tanv2)
238 * dEdM * n;
239 double dotu = dotv + (-_Cuc*sin2u0 + _Cus*cos2u0)*2*dotv;
240 double dotom = _OMEGADOT - omegaEarth;
241 double doti = _IDOT + (-_Cic*sin2u0 + _Cis*cos2u0)*2*dotv;
242 double dotr = a0 * _e*sin(E) * dEdM * n
243 + (-_Crc*sin2u0 + _Crs*cos2u0)*2*dotv;
244 double dotx = dotr*cos(u) - r*sin(u)*dotu;
245 double doty = dotr*sin(u) + r*cos(u)*dotu;
246
247 vv(1) = cosom *dotx - cosi*sinom *doty // dX / dr
248 - xp*sinom*dotom - yp*cosi*cosom*dotom // dX / dOMEGA
249 + yp*sini*sinom*doti; // dX / di
250
251 vv(2) = sinom *dotx + cosi*cosom *doty
252 + xp*cosom*dotom - yp*cosi*sinom*dotom
253 - yp*sini*cosom*doti;
254
255 vv(3) = sini *doty + yp*cosi *doti;
256}
257
258// Read Glonass Ephemeris
259////////////////////////////////////////////////////////////////////////////
260void t_ephGlo::read(const QStringList& lines) {
261
262 for (int ii = 1; ii <= lines.size(); ii++) {
263 QTextStream in(lines.at(ii-1).toAscii());
264
265 if (ii == 1) {
266 double year, month, day, hour, minute, second;
267 in >> _prn >> year >> month >> day >> hour >> minute >> second
268 >> _tau >> _gamma;
269
270 _tau = -_tau;
271
272 if (year < 100) year += 2000;
273
274 QDateTime dateTime(QDate(int(year), int(month), int(day)),
275 QTime(int(hour), int(minute), int(second)), Qt::UTC);
276
277 GPSweekFromDateAndTime(dateTime, _GPSweek, _GPSweeks);
278
279 //// beg test
280 //// _gps_utc = 14.0;
281 //// end test
282
283 _GPSweeks += _gps_utc;
284 }
285 else if (ii == 2) {
286 in >>_x_pos >> _x_velocity >> _x_acceleration >> _health;
287 }
288 else if (ii == 3) {
289 in >>_y_pos >> _y_velocity >> _y_acceleration >> _frequency_number;
290 }
291 else if (ii == 4) {
292 in >>_z_pos >> _z_velocity >> _z_acceleration >> _E;
293 }
294 }
295
296 // Initialize status vector
297 // ------------------------
298 _tt = _GPSweeks;
299
300 _xv(1) = _x_pos * 1.e3;
301 _xv(2) = _y_pos * 1.e3;
302 _xv(3) = _z_pos * 1.e3;
303 _xv(4) = _x_velocity * 1.e3;
304 _xv(5) = _y_velocity * 1.e3;
305 _xv(6) = _z_velocity * 1.e3;
306}
307
308// Derivative of the state vector using a simple force model (static)
309////////////////////////////////////////////////////////////////////////////
310ColumnVector t_ephGlo::glo_deriv(double /* tt */, const ColumnVector& xv) {
311
312 // State vector components
313 // -----------------------
314 ColumnVector rr = xv.rows(1,3);
315 ColumnVector vv = xv.rows(4,6);
316
317 // Acceleration
318 // ------------
319 const static double GM = 398.60044e12;
320 const static double AE = 6378136.0;
321 const static double OMEGA = 7292115.e-11;
322 const static double C20 = -1082.63e-6;
323
324 double rho = rr.norm_Frobenius();
325 double t1 = -GM/(rho*rho*rho);
326 double t2 = 3.0/2.0 * C20 * (GM*AE*AE) / (rho*rho*rho*rho*rho);
327 double t3 = OMEGA * OMEGA;
328 double t4 = 2.0 * OMEGA;
329 double z2 = rr(3) * rr(3);
330
331 // Vector of derivatives
332 // ---------------------
333 ColumnVector va(6);
334 va(1) = vv(1);
335 va(2) = vv(2);
336 va(3) = vv(3);
337 va(4) = (t1 + t2*(1.0-5.0*z2/(rho*rho)) + t3) * rr(1) + t4*vv(2);
338 va(5) = (t1 + t2*(1.0-5.0*z2/(rho*rho)) + t3) * rr(2) - t4*vv(1);
339 va(6) = (t1 + t2*(3.0-5.0*z2/(rho*rho)) ) * rr(3);
340
341 return va;
342}
343
344// Compute Glonass Satellite Position
345////////////////////////////////////////////////////////////////////////////
346void t_ephGlo::position(int GPSweek, double GPSweeks, ColumnVector& xc,
347 ColumnVector& vv) const {
348
349 const static double secPerWeek = 7 * 86400.0;
350 const static double nominalStep = 10.0;
351
352 double dtPos = GPSweeks - _tt;
353 if (GPSweek != _GPSweek) {
354 dtPos += (GPSweek - _GPSweek) * secPerWeek;
355 }
356
357 int nSteps = int(fabs(dtPos) / nominalStep) + 1;
358 double step = dtPos / nSteps;
359
360 for (int ii = 1; ii <= nSteps; ii++) {
361 _xv = rungeKutta4(_tt, _xv, step, glo_deriv);
362 _tt += step;
363 }
364
365 // Position and Velocity
366 // ---------------------
367 xc(1) = _xv(1);
368 xc(2) = _xv(2);
369 xc(3) = _xv(3);
370
371 vv(1) = _xv(4);
372 vv(2) = _xv(5);
373 vv(3) = _xv(6);
374
375 // Clock Correction
376 // ----------------
377 double dtClk = GPSweeks - _GPSweeks;
378 if (GPSweek != _GPSweek) {
379 dtClk += (GPSweek - _GPSweek) * secPerWeek;
380 }
381 xc(4) = -_tau + _gamma * dtClk;
382}
383
384// Glonass IOD
385////////////////////////////////////////////////////////////////////////////
386int t_ephGlo::IOD() const {
387 //// return int(fmod(_GPSweeks,86400.0)) / 600;
388
389 //// unsigned int tb = int(fmod(_GPSweeks,86400.0)) * 1000; // msec of day
390 unsigned int tb = int(fmod(_GPSweeks,86400.0)); //sec of day
391
392 // 5 LSBs of iod are equal to 5 LSBs of tb, remaining bits are zero
393 // ----------------------------------------------------------------
394 const int shift = sizeof(tb) * 8 - 5;
395 unsigned int iod = tb << shift;
396 return (iod >> shift);
397}
Note: See TracBrowser for help on using the repository browser.