0.4.1
Loading...
Searching...
No Matches
GalileoEphemeris.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
10
13
14#include "util/Logger.hpp"
15
16namespace NAV
17{
18
21
23 const size_t& IODnav,
24 const std::array<double, 3>& a,
25 const double& sqrt_A, const double& e, const double& i_0, const double& Omega_0, const double& omega, const double& M_0,
26 const double& delta_n, const double& Omega_dot, const double& i_dot, const double& Cus, const double& Cuc,
27 const double& Cis, const double& Cic, const double& Crs, const double& Crc,
28 const std::bitset<10>& dataSource, const double& signalAccuracy, const SvHealth& svHealth,
29 const double& BGD_E1_E5a, const double& BGD_E1_E5b)
31 toc(toc),
32 toe(toe),
34 a(a),
36 e(e),
37 i_0(i_0),
39 omega(omega),
40 M_0(M_0),
43 i_dot(i_dot),
44 Cus(Cus),
45 Cuc(Cuc),
46 Cis(Cis),
47 Cic(Cic),
48 Crs(Crs),
49 Crc(Crc),
55
56#ifdef TESTING
57
58GalileoEphemeris::GalileoEphemeris(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 IODnav, 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 dataSource, double GALWeek, double /* spare1 */,
64 double signalAccuracy, double svHealth, double BGD_E1_E5a, double BGD_E1_E5b,
65 double /* TransmissionTimeOfMessage */, double /* spare2 */, double /* spare3 */, double /* spare4 */)
66 : SatNavData(SatNavData::GalileoEphemeris, InsTime(year, month, day, hour, minute, second, GST)),
67 toc(refTime),
68 toe(InsTime(0, static_cast<int32_t>(GALWeek), Toe, GST)),
69 IODnav(static_cast<size_t>(IODnav)),
70 a({ svClockBias, svClockDrift, svClockDriftRate }),
71 sqrt_A(sqrt_A),
72 e(e),
73 i_0(i_0),
74 Omega_0(Omega_0),
75 omega(omega),
76 M_0(M_0),
77 delta_n(delta_n),
78 Omega_dot(Omega_dot),
79 i_dot(i_dot),
80 Cus(Cus),
81 Cuc(Cuc),
82 Cis(Cis),
83 Cic(Cic),
84 Crs(Crs),
85 Crc(Crc),
86 dataSource(static_cast<uint16_t>(dataSource)),
87 signalAccuracy(signalAccuracy),
88 svHealth({ .E5a_DataValidityStatus = static_cast<GalileoEphemeris::SvHealth::DataValidityStatus>((static_cast<uint16_t>(svHealth) & 0b000001000) >> 3),
89 .E5b_DataValidityStatus = static_cast<GalileoEphemeris::SvHealth::DataValidityStatus>((static_cast<uint16_t>(svHealth) & 0b001000000) >> 6),
90 .E1B_DataValidityStatus = static_cast<GalileoEphemeris::SvHealth::DataValidityStatus>((static_cast<uint16_t>(svHealth) & 0b000000001) >> 0),
91 .E5a_SignalHealthStatus = static_cast<GalileoEphemeris::SvHealth::SignalHealthStatus>((static_cast<uint16_t>(svHealth) & 0b000110000) >> 4),
92 .E5b_SignalHealthStatus = static_cast<GalileoEphemeris::SvHealth::SignalHealthStatus>((static_cast<uint16_t>(svHealth) & 0b110000000) >> 7),
93 .E1B_SignalHealthStatus = static_cast<GalileoEphemeris::SvHealth::SignalHealthStatus>((static_cast<uint16_t>(svHealth) & 0b000000110) >> 1) }),
94 BGD_E1_E5a(BGD_E1_E5a),
95 BGD_E1_E5b(BGD_E1_E5b)
96{}
97
98#endif
99
100Clock::Corrections GalileoEphemeris::calcClockCorrections(const InsTime& recvTime, double dist, const Frequency& freq) const
101{
102 LOG_DATA("Calc Sat Clock corrections at receiver time {}", recvTime.toGPSweekTow());
103 // Earth gravitational constant [m³/s²] (WGS 84 value of the earth's gravitational constant for GPS user)
104 const auto mu = InsConst::GAL::MU;
105 // Relativistic constant F for clock corrections [s/√m] (-2*√µ/c²)
106 const auto F = InsConst::GAL::F;
107
108 LOG_DATA(" toe {} (Time of ephemeris)", toe.toGPSweekTow());
109
110 const auto A = sqrt_A * sqrt_A; // Semi-major axis [m]
111 LOG_DATA(" A {} [m] (Semi-major axis)", A);
112 auto n_0 = std::sqrt(mu / std::pow(A, 3)); // Computed mean motion [rad/s]
113 LOG_DATA(" n_0 {} [rad/s] (Computed mean motion)", n_0);
114 auto n = n_0 + delta_n; // Corrected mean motion [rad/s]
115 LOG_DATA(" n {} [rad/s] (Corrected mean motion)", n);
116
117 // Time at transmission
118 InsTime transTime0 = recvTime - std::chrono::duration<double>(dist / InsConst::C);
119
120 InsTime transTime = transTime0;
121 LOG_DATA(" Iterating Time at transmission");
122 double dt_sv = 0.0;
123 double clkDrift = 0.0;
124
125 for (size_t i = 0; i < 2; i++)
126 {
127 LOG_DATA(" transTime {} (Time at transmission)", transTime.toGPSweekTow());
128
129 // [s]
130 auto t_minus_toc = static_cast<double>((transTime - toc).count());
131 LOG_DATA(" transTime - toc {} [s]", t_minus_toc);
132
133 // Time difference from ephemeris reference epoch [s]
134 double t_k = static_cast<double>((transTime - toe).count());
135 LOG_DATA(" transTime - toe {} [s] (t_k = Time difference from ephemeris reference epoch)", t_k);
136
137 // Mean anomaly [rad]
138 auto M_k = M_0 + n * t_k;
139 LOG_DATA(" M_k {} [s] (Mean anomaly)", M_k);
140
141 // Eccentric anomaly [rad]
142 double E_k = M_k;
143 double E_k_old = 0.0;
144
145 for (size_t i = 0; std::abs(E_k - E_k_old) > 1e-13 && i < 10; i++)
146 {
147 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:
148 E_k = M_k + e * sin(E_k);
149 }
150
151 // Relativistic correction term [s]
152 double dt_r = F * e * sqrt_A * std::sin(E_k);
153 LOG_DATA(" dt_r {} [s] (Relativistic correction term)", dt_r);
154
155 // SV PRN code phase time offset [s]
156 dt_sv = a[0] + a[1] * t_minus_toc + a[2] * std::pow(t_minus_toc, 2) + dt_r;
157
158 // See GAL-ICD-2.0 GAL ICD, ch. 5.1.5, p.47
159 dt_sv -= ratioFreqSquared(E01, freq, -128, -128) * (freq == E05 ? BGD_E1_E5a : BGD_E1_E5b); // TODO: Check again
160
161 LOG_DATA(" dt_sv {} [s] (SV PRN code phase time offset)", dt_sv);
162
163 // Groves ch. 9.3.1, eq. 9.78, p. 391
164 clkDrift = a[1] + a[2] / 2.0 * t_minus_toc;
165
166 // Correct transmit time for the satellite clock bias
167 transTime = transTime0 - std::chrono::duration<double>(dt_sv);
168 }
169 LOG_DATA(" transTime {} (Time at transmission)", transTime.toGPSweekTow());
170
171 return { .transmitTime = transTime, .bias = dt_sv, .drift = clkDrift };
172}
173
175{
176 Eigen::Vector3d e_pos = Eigen::Vector3d::Zero();
177 Eigen::Vector3d e_vel = Eigen::Vector3d::Zero();
178 Eigen::Vector3d e_accel = Eigen::Vector3d::Zero();
179
180 LOG_DATA("Calc Sat Position at transmit time {}", transTime.toGPSweekTow());
181 // Earth gravitational constant [m³/s²] (WGS 84 value of the earth's gravitational constant for GPS user)
182 const auto mu = InsConst::GAL::MU;
183 // Earth angular velocity [rad/s] (WGS 84 value of the earth's rotation rate)
184 const auto Omega_e_dot = InsConst::GAL::omega_ie;
185
186 LOG_DATA(" toe {} (Time of ephemeris)", toe.toGPSweekTow());
187
188 const auto A = sqrt_A * sqrt_A; // Semi-major axis [m]
189 LOG_DATA(" A {} [m] (Semi-major axis)", A);
190 auto n_0 = std::sqrt(mu / std::pow(A, 3)); // Computed mean motion [rad/s]
191 LOG_DATA(" n_0 {} [rad/s] (Computed mean motion)", n_0);
192 auto n = n_0 + delta_n; // Corrected mean motion [rad/s]
193 LOG_DATA(" n {} [rad/s] (Corrected mean motion)", n);
194
195 // Eccentric anomaly [rad]
196 double E_k = 0.0;
197
198 // Time difference from ephemeris reference epoch [s]
199 double t_k = static_cast<double>((transTime - toe).count());
200 LOG_DATA(" t_k {} [s] (Time difference from ephemeris reference epoch)", t_k);
201
202 // Mean anomaly [rad]
203 auto M_k = M_0 + n * t_k;
204 LOG_DATA(" M_k {} [s] (Mean anomaly)", M_k);
205
206 E_k = M_k; // Initial Value [rad]
207 double E_k_old = 0.0;
208 LOG_DATA(" Iterating E_k");
209 LOG_DATA(" E_k {} [rad] (Eccentric anomaly)", E_k);
210 for (size_t i = 0; std::abs(E_k - E_k_old) > 1e-13 && i < 10; i++)
211 {
212 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:
213 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)
214 LOG_DATA(" E_k {} [rad] (Eccentric anomaly)", E_k); // - Final Value (radians)
215 }
216
217 // 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)
218 // 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)
219 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
220 LOG_DATA(" v_k {} [rad] (True Anomaly (unambiguous quadrant))", v_k);
221 auto Phi_k = v_k + omega; // Argument of Latitude [rad]
222 LOG_DATA(" Phi_k {} [rad] (Argument of Latitude)", Phi_k);
223
224 // Second Harmonic Perturbations
225 auto delta_u_k = Cus * std::sin(2 * Phi_k) + Cuc * std::cos(2 * Phi_k); // Argument of Latitude Correction [rad]
226 LOG_DATA(" delta_u_k {} [rad] (Argument of Latitude Correction)", delta_u_k);
227 auto delta_r_k = Crs * std::sin(2 * Phi_k) + Crc * std::cos(2 * Phi_k); // Radius Correction [m]
228 LOG_DATA(" delta_r_k {} [m] (Radius Correction)", delta_r_k);
229 auto delta_i_k = Cis * std::sin(2 * Phi_k) + Cic * std::cos(2 * Phi_k); // Inclination Correction [rad]
230 LOG_DATA(" delta_i_k {} [rad] (Inclination Correction)", delta_i_k);
231
232 auto u_k = Phi_k + delta_u_k; // Corrected Argument of Latitude [rad]
233 LOG_DATA(" u_k {} [rad] (Corrected Argument of Latitude)", u_k);
234 auto r_k = A * (1 - e * std::cos(E_k)) + delta_r_k; // Corrected Radius [m]
235 LOG_DATA(" r_k {} [m] (Corrected Radius)", r_k);
236 auto i_k = i_0 + delta_i_k + i_dot * t_k; // Corrected Inclination [rad]
237 LOG_DATA(" i_k {} [rad] (Corrected Inclination)", i_k);
238
239 auto x_k_op = r_k * std::cos(u_k); // Position in orbital plane [m]
240 LOG_DATA(" x_k_op {} [m] (Position in orbital plane)", x_k_op);
241 auto y_k_op = r_k * std::sin(u_k); // Position in orbital plane [m]
242 LOG_DATA(" y_k_op {} [m] (Position in orbital plane)", y_k_op);
243
244 // Corrected longitude of ascending node [rad]
245 auto Omega_k = Omega_0 + (Omega_dot - Omega_e_dot) * t_k - Omega_e_dot * static_cast<double>(toe.toGPSweekTow().tow);
246 LOG_DATA(" Omega_k {} [rad] (Corrected longitude of ascending node)", Omega_k);
247
248 // Earth-fixed x coordinates [m]
249 auto x_k = x_k_op * std::cos(Omega_k) - y_k_op * std::cos(i_k) * std::sin(Omega_k);
250 LOG_DATA(" x_k {} [m] (Earth-fixed x coordinates)", x_k);
251 // Earth-fixed y coordinates [m]
252 auto y_k = x_k_op * std::sin(Omega_k) + y_k_op * std::cos(i_k) * std::cos(Omega_k);
253 LOG_DATA(" y_k {} [m] (Earth-fixed y coordinates)", y_k);
254 // Earth-fixed z coordinates [m]
255 auto z_k = y_k_op * std::sin(i_k);
256 LOG_DATA(" z_k {} [m] (Earth-fixed z coordinates)", z_k);
257
258 e_pos = Eigen::Vector3d{ x_k, y_k, z_k };
259
260 if (calc & Calc_Velocity || calc & Calc_Acceleration)
261 {
262 // Eccentric Anomaly Rate [rad/s]
263 auto E_k_dot = n / (1 - e * std::cos(E_k));
264 // True Anomaly Rate [rad/s]
265 auto v_k_dot = E_k_dot * std::sqrt(1 - e * e) / (1 - e * std::cos(E_k));
266 // Corrected Inclination Angle Rate [rad/s]
267 auto i_k_dot = i_dot + 2 * v_k_dot * (Cis * std::cos(2 * Phi_k) - Cic * std::sin(2 * Phi_k));
268 // Corrected Argument of Latitude Rate [rad/s]
269 auto u_k_dot = v_k_dot + 2 * v_k_dot * (Cus * std::cos(2 * Phi_k) - Cuc * std::sin(2 * Phi_k));
270 // Corrected Radius Rate [m/s]
271 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));
272 // Longitude of Ascending Node Rate [rad/s]
273 auto Omega_k_dot = Omega_dot - Omega_e_dot;
274 // In-plane x velocity [m/s]
275 auto vx_k_op = r_k_dot * std::cos(u_k) - r_k * u_k_dot * std::sin(u_k);
276 // In-plane y velocity [m/s]
277 auto vy_k_op = r_k_dot * std::sin(u_k) + r_k * u_k_dot * std::cos(u_k);
278 // Earth-Fixed x velocity [m/s]
279 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)
280 - 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));
281 // Earth-Fixed y velocity [m/s]
282 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)
283 - 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));
284 // Earth-Fixed z velocity [m/s]
285 auto vz_k = vy_k_op * std::sin(i_k) + y_k_op * i_k_dot * std::cos(i_k);
286
287 if (calc & Calc_Velocity)
288 {
289 e_vel = Eigen::Vector3d{ vx_k, vy_k, vz_k };
290 }
291
292 if (calc & Calc_Acceleration)
293 {
294 // Oblate Earth acceleration Factor [m/s^2]
295 auto F = -(3.0 / 2.0) * InsConst::GPS::J2 * (mu / std::pow(r_k, 2)) * std::pow(InsConst::GPS::R_E / r_k, 2);
296 // Earth-Fixed x acceleration [m/s^2]
297 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))
298 + 2 * vy_k * Omega_e_dot + x_k * std::pow(Omega_e_dot, 2);
299 // Earth-Fixed y acceleration [m/s^2]
300 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))
301 + 2 * vx_k * Omega_e_dot + y_k * std::pow(Omega_e_dot, 2);
302 // Earth-Fixed z acceleration [m/s^2]
303 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));
304
305 e_accel = Eigen::Vector3d{ ax_k, ay_k, az_k };
306 }
307 }
308
309 return { .e_pos = e_pos,
310 .e_vel = e_vel,
311 .e_accel = e_accel };
312}
313
315{
316 return svHealth.E5a_DataValidityStatus == SvHealth::DataValidityStatus::NavigationDataValid
319 && svHealth.E5a_SignalHealthStatus == SvHealth::SignalHealthStatus::SignalOK
320 && svHealth.E5b_SignalHealthStatus == SvHealth::SignalHealthStatus::SignalOK
321 && svHealth.E1B_SignalHealthStatus == SvHealth::SignalHealthStatus::SignalOK;
322}
323
325{
326 return std::pow(signalAccuracy, 2);
327}
328
329} // namespace NAV
Holds all Constants.
GNSS helper functions.
Galileo Ephemeris information.
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
Frequency definition for different satellite systems.
Definition Frequency.hpp:59
Broadcasted ephemeris message data.
InsTime toe
Time of Ephemeris.
double calcSatellitePositionVariance() const final
Calculates the Variance of the satellite position in [m^2].
double Omega_0
Longitude of the ascending node at reference time [rad].
double BGD_E1_E5b
E1-E5b Broadcast Group Delay [s].
bool isHealthy() const final
Checks whether the signal is healthy.
double Omega_dot
Rate of change of right ascension [rad/s].
PosVelAccel calcSatelliteData(const InsTime &transTime, Orbit::Calc calc) const final
Calculates position, velocity and acceleration of the satellite at transmission time.
double e
Eccentricity [-].
std::bitset< 10 > dataSource
Data sources.
double Cis
Amplitude of the sine harmonic correction term to the angle of inclination [rad].
double omega
Argument of perigee [rad].
double M_0
Mean anomaly at reference time [rad].
double i_dot
Rate of change of inclination [rad/s].
double Cuc
Amplitude of the cosine harmonic correction term to the argument of latitude [rad].
double Cus
Amplitude of the sine harmonic correction term to the argument of latitude [rad].
Corrections calcClockCorrections(const InsTime &recvTime, double dist, const Frequency &freq) const final
Calculates clock bias and drift of the satellite.
double sqrt_A
Square root of the semi-major axis [m^1/2].
GalileoEphemeris(const InsTime &toc)
Default Constructor.
double Cic
Amplitude of the cosine harmonic correction term to the angle of inclination [rad].
size_t IODnav
Issue of Data of the nav batch.
double Crc
Amplitude of the cosine harmonic correction term to the orbit radius [m].
std::array< double, 3 > a
double Crs
Amplitude of the sine harmonic correction term to the orbit radius [m].
double delta_n
Mean motion difference from computed value [rad/s].
InsTime toc
Time of Clock.
double BGD_E1_E5a
E1-E5a Broadcast Group Delay [s].
SvHealth svHealth
Signal Health.
double i_0
Inclination angle at reference time [rad].
static constexpr double MU
Earth gravitational constant GALILEO [m³/s²].
static constexpr double omega_ie
Earth angular velocity GALILEO [rad/s].
static constexpr double F
Relativistic constant F for clock corrections [s/√m] (-2*√µ/c²)
static constexpr double R_E
Earth Equatorial Radius [m].
static constexpr double J2
Oblate Earth Gravity Coefficient [-].
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
Satellite Navigation data (to calculate SatNavData and clock)
SatNavData(Type type, const InsTime &refTime)
Constructor.
@ GalileoEphemeris
Galileo Broadcast Ephemeris.
@ GST
Galileo System Time.
double ratioFreqSquared(Frequency f1, Frequency f2, int8_t num1, int8_t num2)
Calculates the ration of the frequencies squared γ
Definition Functions.cpp:24
@ E01
Galileo, "E1" (1575.42 MHz).
Definition Frequency.hpp:31
@ E05
Galileo E5a (1176.45 MHz).
Definition Frequency.hpp:32
Satellite clock corrections.
Definition Clock.hpp:28
Navigation Data Validity and Signal Health Status.
DataValidityStatus
Navigation Data Validity.
@ NavigationDataValid
Navigation data valid.
SignalHealthStatus
Signal Health Status.
Satellite Position, Velocity and Acceleration.
Definition Orbit.hpp:40