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

Last change on this file since 884 was 884, checked in by mervart, 16 years ago

* empty log message *

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