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

Last change on this file since 929 was 929, checked in by weber, 16 years ago

* empty log message *

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