0.4.1
Loading...
Searching...
No Matches
QZSSEphemeris.cpp
Go to the documentation of this file.
1// This file is part of INSTINCT, the INS Toolkit for Integrated
2// Navigation Concepts and Training by the Institute of Navigation of
3// the University of Stuttgart, Germany.
4//
5// This Source Code Form is subject to the terms of the Mozilla Public
6// License, v. 2.0. If a copy of the MPL was not distributed with this
7// file, You can obtain one at https://mozilla.org/MPL/2.0/.
8
9#include "QZSSEphemeris.hpp"
10
13
14#include "util/Logger.hpp"
15
16namespace NAV
17{
18
20 const size_t& IODE, const size_t& IODC,
21 const std::array<double, 3>& a,
22 const double& sqrt_A, const double& e, const double& i_0, const double& Omega_0, const double& omega, const double& M_0,
23 const double& delta_n, const double& Omega_dot, const double& i_dot, const double& Cus, const double& Cuc,
24 const double& Cis, const double& Cic, const double& Crs, const double& Crc,
25 const double& svAccuracy, uint8_t svHealth,
26 uint8_t L2ChannelCodes, bool L2DataFlagPCode,
27 const double& T_GD, bool fitIntervalFlag)
29 toc(toc),
30 toe(toe),
31 IODE(IODE),
32 IODC(IODC),
33 a(a),
35 e(e),
36 i_0(i_0),
38 omega(omega),
39 M_0(M_0),
42 i_dot(i_dot),
43 Cus(Cus),
44 Cuc(Cuc),
45 Cis(Cis),
46 Cic(Cic),
47 Crs(Crs),
48 Crc(Crc),
53 T_GD(T_GD),
55
56#ifdef TESTING
57
58QZSSEphemeris::QZSSEphemeris(int32_t year, int32_t month, int32_t day, int32_t hour, int32_t minute, double second, double svClockBias, double svClockDrift, double svClockDriftRate,
59 double IODE, double Crs, double delta_n, double M_0,
60 double Cuc, double e, double Cus, double sqrt_A,
61 double Toe, double Cic, double Omega_0, double Cis,
62 double i_0, double Crc, double omega, double Omega_dot,
63 double i_dot, double L2ChannelCodes, double GPSWeek, double L2DataFlagPCode,
64 double svAccuracy, double svHealth, double T_GD, double IODC,
65 double /* TransmissionTimeOfMessage */, double fitIntervalFlag, double /* spare1 */, double /* spare2 */)
66 : SatNavData(SatNavData::QZSSEphemeris, InsTime(year, month, day, hour, minute, second, SatelliteSystem(QZSS).getTimeSystem())),
67 toc(refTime),
68 toe(InsTime(0, static_cast<int32_t>(GPSWeek), Toe, SatelliteSystem(QZSS).getTimeSystem())),
69 IODE(static_cast<size_t>(IODE)),
70 IODC(static_cast<size_t>(IODC)),
71 a({ svClockBias, svClockDrift, svClockDriftRate }),
72 sqrt_A(sqrt_A),
73 e(e),
74 i_0(i_0),
75 Omega_0(Omega_0),
76 omega(omega),
77 M_0(M_0),
78 delta_n(delta_n),
79 Omega_dot(Omega_dot),
80 i_dot(i_dot),
81 Cus(Cus),
82 Cuc(Cuc),
83 Cis(Cis),
84 Cic(Cic),
85 Crs(Crs),
86 Crc(Crc),
87 svAccuracy(svAccuracy),
88 svHealth(static_cast<uint8_t>(svHealth)),
89 L2ChannelCodes(static_cast<uint8_t>(L2ChannelCodes)),
90 L2DataFlagPCode(static_cast<bool>(L2DataFlagPCode)),
91 T_GD(T_GD),
92 fitIntervalFlag(static_cast<bool>(fitIntervalFlag))
93{}
94
95#endif
96
97Clock::Corrections QZSSEphemeris::calcClockCorrections(const InsTime& recvTime, double dist, const Frequency& freq) const
98{
99 LOG_DATA("Calc Sat Clock corrections at receiver time {}", recvTime.toGPSweekTow());
100 // Earth gravitational constant [m³/s²] (WGS 84 value of the earth's gravitational constant for GPS user)
101 const auto mu = InsConst::QZSS::MU;
102 // Relativistic constant F for clock corrections [s/√m] (-2*√µ/c²)
103 const auto F = InsConst::QZSS::F;
104
105 LOG_DATA(" toe {} (Time of ephemeris)", toe.toGPSweekTow());
106
107 const auto A = sqrt_A * sqrt_A; // Semi-major axis [m]
108 LOG_DATA(" A {} [m] (Semi-major axis)", A);
109 auto n_0 = std::sqrt(mu / std::pow(A, 3)); // Computed mean motion [rad/s]
110 LOG_DATA(" n_0 {} [rad/s] (Computed mean motion)", n_0);
111 auto n = n_0 + delta_n; // Corrected mean motion [rad/s]
112 LOG_DATA(" n {} [rad/s] (Corrected mean motion)", n);
113
114 // Time at transmission
115 InsTime transTime0 = recvTime - std::chrono::duration<double>(dist / InsConst::C);
116
117 InsTime transTime = transTime0;
118 LOG_DATA(" Iterating Time at transmission");
119 double dt_sv = 0.0;
120 double clkDrift = 0.0;
121
122 for (size_t i = 0; i < 2; i++)
123 {
124 LOG_DATA(" transTime {} (Time at transmission)", transTime.toGPSweekTow());
125
126 // [s]
127 auto t_minus_toc = static_cast<double>((transTime - toc).count());
128 LOG_DATA(" transTime - toc {} [s]", t_minus_toc);
129
130 // Time difference from ephemeris reference epoch [s]
131 double t_k = static_cast<double>((transTime - toe).count());
132 LOG_DATA(" transTime - toe {} [s] (t_k = Time difference from ephemeris reference epoch)", t_k);
133
134 // Mean anomaly [rad]
135 auto M_k = M_0 + n * t_k;
136 LOG_DATA(" M_k {} [s] (Mean anomaly)", M_k);
137
138 // Eccentric anomaly [rad]
139 double E_k = M_k;
140 double E_k_old = 0.0;
141
142 for (size_t i = 0; std::abs(E_k - E_k_old) > 1e-13 && i < 10; i++)
143 {
144 E_k_old = E_k; // Kepler’s equation ( Mk = E_k − e sin E_k ) may be solved for Eccentric anomaly (E_k) by iteration:
145 E_k = M_k + e * sin(E_k);
146 }
147
148 // Relativistic correction term [s]
149 double dt_r = F * e * sqrt_A * std::sin(E_k);
150 LOG_DATA(" dt_r {} [s] (Relativistic correction term)", dt_r);
151
152 // SV PRN code phase time offset [s]
153 dt_sv = a[0] + a[1] * t_minus_toc + a[2] * std::pow(t_minus_toc, 2) + dt_r;
154
155 // See IS-GPS-200M GPS ICD, ch. 20.3.3.3.3.2, p.102
156 dt_sv -= ratioFreqSquared(J01, freq, -128, -128) * T_GD; // TODO: Check again
157
158 LOG_DATA(" dt_sv {} [s] (SV PRN code phase time offset)", dt_sv);
159
160 // Groves ch. 9.3.1, eq. 9.78, p. 391
161 clkDrift = a[1] + a[2] / 2.0 * t_minus_toc;
162
163 // Correct transmit time for the satellite clock bias
164 transTime = transTime0 - std::chrono::duration<double>(dt_sv);
165 }
166 LOG_DATA(" transTime {} (Time at transmission)", transTime.toGPSweekTow());
167
168 return { .transmitTime = transTime, .bias = dt_sv, .drift = clkDrift };
169}
170
172{
173 Eigen::Vector3d e_pos = Eigen::Vector3d::Zero();
174 Eigen::Vector3d e_vel = Eigen::Vector3d::Zero();
175 Eigen::Vector3d e_accel = Eigen::Vector3d::Zero();
176
177 LOG_DATA("Calc Sat Position at transmit time {}", transTime.toGPSweekTow());
178 // Earth gravitational constant [m³/s²] (WGS 84 value of the earth's gravitational constant for GPS user)
179 const auto mu = InsConst::QZSS::MU;
180 // Earth angular velocity [rad/s] (WGS 84 value of the earth's rotation rate)
181 const auto Omega_e_dot = InsConst::QZSS::omega_ie;
182
183 LOG_DATA(" toe {} (Time of ephemeris)", toe.toGPSweekTow());
184
185 const auto A = sqrt_A * sqrt_A; // Semi-major axis [m]
186 LOG_DATA(" A {} [m] (Semi-major axis)", A);
187 auto n_0 = std::sqrt(mu / std::pow(A, 3)); // Computed mean motion [rad/s]
188 LOG_DATA(" n_0 {} [rad/s] (Computed mean motion)", n_0);
189 auto n = n_0 + delta_n; // Corrected mean motion [rad/s]
190 LOG_DATA(" n {} [rad/s] (Corrected mean motion)", n);
191
192 // Eccentric anomaly [rad]
193 double E_k = 0.0;
194
195 // Time difference from ephemeris reference epoch [s]
196 double t_k = static_cast<double>((transTime - toe).count());
197 LOG_DATA(" t_k {} [s] (Time difference from ephemeris reference epoch)", t_k);
198
199 // Mean anomaly [rad]
200 auto M_k = M_0 + n * t_k;
201 LOG_DATA(" M_k {} [s] (Mean anomaly)", M_k);
202
203 E_k = M_k; // Initial Value [rad]
204 double E_k_old = 0.0;
205 LOG_DATA(" Iterating E_k");
206 LOG_DATA(" E_k {} [rad] (Eccentric anomaly)", E_k);
207 for (size_t i = 0; std::abs(E_k - E_k_old) > 1e-13 && i < 10; i++)
208 {
209 E_k_old = E_k; // Kepler’s equation ( Mk = E_k − e sin E_k ) may be solved for Eccentric anomaly (E_k) by iteration:
210 E_k = E_k + (M_k - E_k + e * std::sin(E_k)) / (1 - e * std::cos(E_k)); // - Refined Value, minimum of three iterations, (j=1,2,3)
211 LOG_DATA(" E_k {} [rad] (Eccentric anomaly)", E_k); // - Final Value (radians)
212 }
213
214 // auto v_k = 2.0 * std::atan(std::sqrt((1.0 + e) / (1.0 - e)) * std::tan(E_k / 2.0)); // True Anomaly (unambiguous quadrant) [rad] (GPS ICD algorithm)
215 // auto v_k = std::atan2(std::sqrt(1 - e * e) * std::sin(E_k) / (1 - e * std::cos(E_k)), (std::cos(E_k) - e) / (1 - e * std::cos(E_k))); // True Anomaly [rad] (GALILEO ICD algorithm)
216 auto v_k = std::atan2(std::sqrt(1 - e * e) * std::sin(E_k), (std::cos(E_k) - e)); // True Anomaly [rad] // simplified, since the denominators cancel out
217 LOG_DATA(" v_k {} [rad] (True Anomaly (unambiguous quadrant))", v_k);
218 auto Phi_k = v_k + omega; // Argument of Latitude [rad]
219 LOG_DATA(" Phi_k {} [rad] (Argument of Latitude)", Phi_k);
220
221 // Second Harmonic Perturbations
222 auto delta_u_k = Cus * std::sin(2 * Phi_k) + Cuc * std::cos(2 * Phi_k); // Argument of Latitude Correction [rad]
223 LOG_DATA(" delta_u_k {} [rad] (Argument of Latitude Correction)", delta_u_k);
224 auto delta_r_k = Crs * std::sin(2 * Phi_k) + Crc * std::cos(2 * Phi_k); // Radius Correction [m]
225 LOG_DATA(" delta_r_k {} [m] (Radius Correction)", delta_r_k);
226 auto delta_i_k = Cis * std::sin(2 * Phi_k) + Cic * std::cos(2 * Phi_k); // Inclination Correction [rad]
227 LOG_DATA(" delta_i_k {} [rad] (Inclination Correction)", delta_i_k);
228
229 auto u_k = Phi_k + delta_u_k; // Corrected Argument of Latitude [rad]
230 LOG_DATA(" u_k {} [rad] (Corrected Argument of Latitude)", u_k);
231 auto r_k = A * (1 - e * std::cos(E_k)) + delta_r_k; // Corrected Radius [m]
232 LOG_DATA(" r_k {} [m] (Corrected Radius)", r_k);
233 auto i_k = i_0 + delta_i_k + i_dot * t_k; // Corrected Inclination [rad]
234 LOG_DATA(" i_k {} [rad] (Corrected Inclination)", i_k);
235
236 auto x_k_op = r_k * std::cos(u_k); // Position in orbital plane [m]
237 LOG_DATA(" x_k_op {} [m] (Position in orbital plane)", x_k_op);
238 auto y_k_op = r_k * std::sin(u_k); // Position in orbital plane [m]
239 LOG_DATA(" y_k_op {} [m] (Position in orbital plane)", y_k_op);
240
241 // Corrected longitude of ascending node [rad]
242 auto Omega_k = Omega_0 + (Omega_dot - Omega_e_dot) * t_k - Omega_e_dot * static_cast<double>(toe.toGPSweekTow().tow);
243 LOG_DATA(" Omega_k {} [rad] (Corrected longitude of ascending node)", Omega_k);
244
245 // Earth-fixed x coordinates [m]
246 auto x_k = x_k_op * std::cos(Omega_k) - y_k_op * std::cos(i_k) * std::sin(Omega_k);
247 LOG_DATA(" x_k {} [m] (Earth-fixed x coordinates)", x_k);
248 // Earth-fixed y coordinates [m]
249 auto y_k = x_k_op * std::sin(Omega_k) + y_k_op * std::cos(i_k) * std::cos(Omega_k);
250 LOG_DATA(" y_k {} [m] (Earth-fixed y coordinates)", y_k);
251 // Earth-fixed z coordinates [m]
252 auto z_k = y_k_op * std::sin(i_k);
253 LOG_DATA(" z_k {} [m] (Earth-fixed z coordinates)", z_k);
254
255 e_pos = Eigen::Vector3d{ x_k, y_k, z_k };
256
257 if (calc & Calc_Velocity || calc & Calc_Acceleration)
258 {
259 // Eccentric Anomaly Rate [rad/s]
260 auto E_k_dot = n / (1 - e * std::cos(E_k));
261 // True Anomaly Rate [rad/s]
262 auto v_k_dot = E_k_dot * std::sqrt(1 - e * e) / (1 - e * std::cos(E_k));
263 // Corrected Inclination Angle Rate [rad/s]
264 auto i_k_dot = i_dot + 2 * v_k_dot * (Cis * std::cos(2 * Phi_k) - Cic * std::sin(2 * Phi_k));
265 // Corrected Argument of Latitude Rate [rad/s]
266 auto u_k_dot = v_k_dot + 2 * v_k_dot * (Cus * std::cos(2 * Phi_k) - Cuc * std::sin(2 * Phi_k));
267 // Corrected Radius Rate [m/s]
268 auto r_k_dot = e * A * E_k_dot * std::sin(E_k) + 2 * v_k_dot * (Crs * std::cos(2 * Phi_k) - Crc * std::sin(2 * Phi_k));
269 // Longitude of Ascending Node Rate [rad/s]
270 auto Omega_k_dot = Omega_dot - Omega_e_dot;
271 // In-plane x velocity [m/s]
272 auto vx_k_op = r_k_dot * std::cos(u_k) - r_k * u_k_dot * std::sin(u_k);
273 // In-plane y velocity [m/s]
274 auto vy_k_op = r_k_dot * std::sin(u_k) + r_k * u_k_dot * std::cos(u_k);
275 // Earth-Fixed x velocity [m/s]
276 auto vx_k = -x_k_op * Omega_k_dot * std::sin(Omega_k) + vx_k_op * std::cos(Omega_k) - vy_k_op * std::sin(Omega_k) * std::cos(i_k)
277 - y_k_op * (Omega_k_dot * std::cos(Omega_k) * std::cos(i_k) - i_k_dot * std::sin(Omega_k) * std::sin(i_k));
278 // Earth-Fixed y velocity [m/s]
279 auto vy_k = x_k_op * Omega_k_dot * std::cos(Omega_k) + vx_k_op * std::sin(Omega_k) + vy_k_op * std::cos(Omega_k) * std::cos(i_k)
280 - y_k_op * (Omega_k_dot * std::sin(Omega_k) * std::cos(i_k) + i_k_dot * std::cos(Omega_k) * std::sin(i_k));
281 // Earth-Fixed z velocity [m/s]
282 auto vz_k = vy_k_op * std::sin(i_k) + y_k_op * i_k_dot * std::cos(i_k);
283
284 if (calc & Calc_Velocity)
285 {
286 e_vel = Eigen::Vector3d{ vx_k, vy_k, vz_k };
287 }
288
289 if (calc & Calc_Acceleration)
290 {
291 // Oblate Earth acceleration Factor [m/s^2]
292 auto F = -(3.0 / 2.0) * InsConst::GPS::J2 * (mu / std::pow(r_k, 2)) * std::pow(InsConst::GPS::R_E / r_k, 2);
293 // Earth-Fixed x acceleration [m/s^2]
294 auto ax_k = -mu * (x_k / std::pow(r_k, 3)) + F * ((1.0 - 5.0 * std::pow(z_k / r_k, 2)) * (x_k / r_k))
295 + 2 * vy_k * Omega_e_dot + x_k * std::pow(Omega_e_dot, 2);
296 // Earth-Fixed y acceleration [m/s^2]
297 auto ay_k = -mu * (y_k / std::pow(r_k, 3)) + F * ((1.0 - 5.0 * std::pow(z_k / r_k, 2)) * (y_k / r_k))
298 + 2 * vx_k * Omega_e_dot + y_k * std::pow(Omega_e_dot, 2);
299 // Earth-Fixed z acceleration [m/s^2]
300 auto az_k = -mu * (z_k / std::pow(r_k, 3)) + F * ((3.0 - 5.0 * std::pow(z_k / r_k, 2)) * (z_k / r_k));
301
302 e_accel = Eigen::Vector3d{ ax_k, ay_k, az_k };
303 }
304 }
305
306 return { .e_pos = e_pos,
307 .e_vel = e_vel,
308 .e_accel = e_accel };
309}
310
311bool QZSSEphemeris::isHealthy() const // TODO Parse Signal Id as a parameter and differentiate depending on the bitset
312{
313 return svHealth.none();
314}
315
317{
318 // Getting the index and value again will discretize the URA values
319 return std::pow(gpsUraIdx2Val(gpsUraVal2Idx(svAccuracy)), 2);
320}
321
322} // namespace NAV
Holds all Constants.
GNSS helper functions.
Utility class for logging to console and file.
#define LOG_DATA
All output which occurs repeatedly every time observations are received.
Definition Logger.hpp:29
QZSS Ephemeris information.
Frequency definition for different satellite systems.
Definition Frequency.hpp:59
static constexpr double R_E
Earth Equatorial Radius [m].
static constexpr double J2
Oblate Earth Gravity Coefficient [-].
static constexpr double F
Relativistic constant F for clock corrections [s/√m] (-2*√µ/c²)
static constexpr double MU
Earth gravitational constant QZSS [m³/s²].
static constexpr double omega_ie
Earth angular velocity QZSS [rad/s].
static constexpr double C
Speed of light [m/s].
Definition Constants.hpp:34
The class is responsible for all time-related tasks.
Definition InsTime.hpp:710
constexpr InsTime_GPSweekTow toGPSweekTow(TimeSystem timesys=GPST) const
Converts this time object into a different format.
Definition InsTime.hpp:854
Calc
Calculation flags.
Definition Orbit.hpp:75
@ Calc_Velocity
Velocity calculation flag.
Definition Orbit.hpp:78
@ Calc_Acceleration
Acceleration calculation flag.
Definition Orbit.hpp:79
Broadcasted ephemeris message data.
const double Crc
Amplitude of the cosine harmonic correction term to the orbit radius [m].
const std::array< double, 3 > a
const double Cic
Amplitude of the cosine harmonic correction term to the angle of inclination [rad].
const size_t IODC
Issue of Data, Clock.
const double T_GD
Group delay between SV clock and L1C/A [s].
const double i_dot
Rate of inclination angle [rad/s].
const InsTime toc
Time of Clock.
const bool fitIntervalFlag
Fit Interval period of the ephemeris.
const double Omega_dot
Rate of right ascension [rad/s].
const std::bitset< 6 > svHealth
SV health.
const double Cuc
Amplitude of the cosine harmonic correction term to the argument of latitude [rad].
QZSSEphemeris(const InsTime &toc, const InsTime &toe, const size_t &IODE, const size_t &IODC, const std::array< double, 3 > &a, const double &sqrt_A, const double &e, const double &i_0, const double &Omega_0, const double &omega, const double &M_0, const double &delta_n, const double &Omega_dot, const double &i_dot, const double &Cus, const double &Cuc, const double &Cis, const double &Cic, const double &Crs, const double &Crc, const double &svAccuracy, uint8_t svHealth, uint8_t L2ChannelCodes, bool L2DataFlagPCode, const double &T_GD, bool fitIntervalFlag)
Constructor.
bool isHealthy() const final
Checks whether the signal is healthy.
double calcSatellitePositionVariance() const final
Calculates the Variance of the satellite position in [m^2].
const double M_0
Mean anomaly at reference time [rad].
Corrections calcClockCorrections(const InsTime &recvTime, double dist, const Frequency &freq) const final
Calculates clock bias and drift of the satellite.
const double Crs
Amplitude of the sine harmonic correction term to the orbit radius [m].
const double sqrt_A
Square root of the semi-major axis [m^1/2].
const bool L2DataFlagPCode
Data Flag for L2 P-Code.
const size_t IODE
Issue of Data, Ephemeris.
const double Cis
Amplitude of the sine harmonic correction term to the angle of inclination [rad].
const double svAccuracy
SV accuracy [m].
PosVelAccel calcSatelliteData(const InsTime &transTime, Orbit::Calc calc) const final
Calculates position, velocity and acceleration of the satellite at transmission time.
const double Omega_0
Longitude of the ascending node at reference time [rad].
const double omega
Argument of perigee [rad].
const InsTime toe
Time of Ephemeris.
const double delta_n
Mean motion difference from computed value [rad/s].
const uint8_t L2ChannelCodes
Indicate which code(s) is (are) commanded ON for the in-phase component of the L2 channel.
const double Cus
Amplitude of the sine harmonic correction term to the argument of latitude [rad].
const double e
Eccentricity [-].
const double i_0
Inclination angle at reference time [rad].
Satellite Navigation data (to calculate SatNavData and clock)
SatNavData(Type type, const InsTime &refTime)
Constructor.
@ QZSSEphemeris
QZSS Broadcast Ephemeris.
double ratioFreqSquared(Frequency f1, Frequency f2, int8_t num1, int8_t num2)
Calculates the ration of the frequencies squared γ
Definition Functions.cpp:24
@ J01
QZSS L1 (1575.42 MHz).
Definition Frequency.hpp:47
double gpsUraIdx2Val(uint8_t idx)
Converts a GPS URA (user range accuracy) index to it's value.
Definition Functions.cpp:54
uint8_t gpsUraVal2Idx(double val)
Converts a GPS URA (user range accuracy) value to it's index.
Definition Functions.cpp:47
@ QZSS
Quasi-Zenith Satellite System.
Satellite clock corrections.
Definition Clock.hpp:28
Satellite Position, Velocity and Acceleration.
Definition Orbit.hpp:40
Satellite System type.