24 const std::array<double, 3>&
a,
27 const double&
Cis,
const double&
Cic,
const double&
Crs,
const double&
Crc,
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 ,
64 double signalAccuracy,
double svHealth,
double BGD_E1_E5a,
double BGD_E1_E5b,
65 double ,
double ,
double ,
double )
68 toe(
InsTime(0, static_cast<int32_t>(GALWeek), Toe,
GST)),
69 IODnav(static_cast<size_t>(IODnav)),
70 a({ svClockBias, svClockDrift, svClockDriftRate }),
86 dataSource(
static_cast<uint16_t
>(dataSource)),
87 signalAccuracy(signalAccuracy),
94 BGD_E1_E5a(BGD_E1_E5a),
95 BGD_E1_E5b(BGD_E1_E5b)
108 LOG_DATA(
" toe {} (Time of ephemeris)",
toe.toGPSweekTow());
111 LOG_DATA(
" A {} [m] (Semi-major axis)", A);
112 auto n_0 = std::sqrt(mu / std::pow(A, 3));
113 LOG_DATA(
" n_0 {} [rad/s] (Computed mean motion)", n_0);
115 LOG_DATA(
" n {} [rad/s] (Corrected mean motion)", n);
120 InsTime transTime = transTime0;
121 LOG_DATA(
" Iterating Time at transmission");
123 double clkDrift = 0.0;
125 for (
size_t i = 0; i < 2; i++)
130 auto t_minus_toc =
static_cast<double>((transTime -
toc).count());
131 LOG_DATA(
" transTime - toc {} [s]", t_minus_toc);
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);
138 auto M_k =
M_0 + n * t_k;
139 LOG_DATA(
" M_k {} [s] (Mean anomaly)", M_k);
143 double E_k_old = 0.0;
145 for (
size_t i = 0; std::abs(E_k - E_k_old) > 1
e-13 && i < 10; i++)
148 E_k = M_k +
e * sin(E_k);
152 double dt_r = F *
e *
sqrt_A * std::sin(E_k);
153 LOG_DATA(
" dt_r {} [s] (Relativistic correction term)", dt_r);
156 dt_sv =
a[0] +
a[1] * t_minus_toc +
a[2] * std::pow(t_minus_toc, 2) + dt_r;
161 LOG_DATA(
" dt_sv {} [s] (SV PRN code phase time offset)", dt_sv);
164 clkDrift =
a[1] +
a[2] / 2.0 * t_minus_toc;
167 transTime = transTime0 - std::chrono::duration<double>(dt_sv);
171 return { .transmitTime = transTime, .bias = dt_sv, .drift = clkDrift };
176 Eigen::Vector3d e_pos = Eigen::Vector3d::Zero();
177 Eigen::Vector3d e_vel = Eigen::Vector3d::Zero();
178 Eigen::Vector3d e_accel = Eigen::Vector3d::Zero();
186 LOG_DATA(
" toe {} (Time of ephemeris)",
toe.toGPSweekTow());
189 LOG_DATA(
" A {} [m] (Semi-major axis)", A);
190 auto n_0 = std::sqrt(mu / std::pow(A, 3));
191 LOG_DATA(
" n_0 {} [rad/s] (Computed mean motion)", n_0);
193 LOG_DATA(
" n {} [rad/s] (Corrected mean motion)", n);
199 double t_k =
static_cast<double>((transTime -
toe).count());
200 LOG_DATA(
" t_k {} [s] (Time difference from ephemeris reference epoch)", t_k);
203 auto M_k =
M_0 + n * t_k;
204 LOG_DATA(
" M_k {} [s] (Mean anomaly)", M_k);
207 double E_k_old = 0.0;
209 LOG_DATA(
" E_k {} [rad] (Eccentric anomaly)", E_k);
210 for (
size_t i = 0; std::abs(E_k - E_k_old) > 1
e-13 && i < 10; i++)
213 E_k = E_k + (M_k - E_k +
e * std::sin(E_k)) / (1 -
e * std::cos(E_k));
214 LOG_DATA(
" E_k {} [rad] (Eccentric anomaly)", E_k);
219 auto v_k = std::atan2(std::sqrt(1 -
e *
e) * std::sin(E_k), (std::cos(E_k) -
e));
220 LOG_DATA(
" v_k {} [rad] (True Anomaly (unambiguous quadrant))", v_k);
221 auto Phi_k = v_k +
omega;
222 LOG_DATA(
" Phi_k {} [rad] (Argument of Latitude)", Phi_k);
225 auto delta_u_k =
Cus * std::sin(2 * Phi_k) +
Cuc * std::cos(2 * Phi_k);
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);
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);
230 LOG_DATA(
" delta_i_k {} [rad] (Inclination Correction)", delta_i_k);
232 auto u_k = Phi_k + delta_u_k;
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;
235 LOG_DATA(
" r_k {} [m] (Corrected Radius)", r_k);
236 auto i_k =
i_0 + delta_i_k +
i_dot * t_k;
237 LOG_DATA(
" i_k {} [rad] (Corrected Inclination)", i_k);
239 auto x_k_op = r_k * std::cos(u_k);
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);
242 LOG_DATA(
" y_k_op {} [m] (Position in orbital plane)", y_k_op);
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);
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);
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);
255 auto z_k = y_k_op * std::sin(i_k);
256 LOG_DATA(
" z_k {} [m] (Earth-fixed z coordinates)", z_k);
258 e_pos = Eigen::Vector3d{ x_k, y_k, z_k };
263 auto E_k_dot = n / (1 -
e * std::cos(E_k));
265 auto v_k_dot = E_k_dot * std::sqrt(1 -
e *
e) / (1 -
e * std::cos(E_k));
267 auto i_k_dot =
i_dot + 2 * v_k_dot * (
Cis * std::cos(2 * Phi_k) -
Cic * std::sin(2 * Phi_k));
269 auto u_k_dot = v_k_dot + 2 * v_k_dot * (
Cus * std::cos(2 * Phi_k) -
Cuc * std::sin(2 * Phi_k));
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));
273 auto Omega_k_dot =
Omega_dot - Omega_e_dot;
275 auto vx_k_op = r_k_dot * std::cos(u_k) - r_k * u_k_dot * std::sin(u_k);
277 auto vy_k_op = r_k_dot * std::sin(u_k) + r_k * u_k_dot * std::cos(u_k);
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));
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));
285 auto vz_k = vy_k_op * std::sin(i_k) + y_k_op * i_k_dot * std::cos(i_k);
289 e_vel = Eigen::Vector3d{ vx_k, vy_k, vz_k };
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);
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);
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));
305 e_accel = Eigen::Vector3d{ ax_k, ay_k, az_k };
309 return { .e_pos = e_pos,
311 .e_accel = e_accel };
Galileo Ephemeris information.
Utility class for logging to console and file.
#define LOG_DATA
All output which occurs repeatedly every time observations are received.
Frequency definition for different satellite systems.
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].
The class is responsible for all time-related tasks.
constexpr InsTime_GPSweekTow toGPSweekTow(TimeSystem timesys=GPST) const
Converts this time object into a different format.
@ Calc_Velocity
Velocity calculation flag.
@ Calc_Acceleration
Acceleration calculation flag.
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 γ
@ E01
Galileo, "E1" (1575.42 MHz).
@ E05
Galileo E5a (1176.45 MHz).
Satellite clock corrections.
Navigation Data Validity and Signal Health Status.
DataValidityStatus
Navigation Data Validity.
@ NavigationDataValid
Navigation data valid.
SignalHealthStatus
Signal Health Status.
Satellite Position, Velocity and Acceleration.