21 const std::array<double, 3>&
a,
24 const double&
Cis,
const double&
Cic,
const double&
Crs,
const double&
Crc,
53IRNSSEphemeris::IRNSSEphemeris(int32_t year, int32_t month, int32_t day, int32_t hour, int32_t minute,
double second,
double svClockBias,
double svClockDrift,
double svClockDriftRate,
54 double IODEC,
double Crs,
double delta_n,
double M_0,
55 double Cuc,
double e,
double Cus,
double sqrt_A,
56 double Toe,
double Cic,
double Omega_0,
double Cis,
57 double i_0,
double Crc,
double omega,
double Omega_dot,
58 double i_dot,
double ,
double IRNWeek,
double ,
59 double svAccuracy,
double svHealth,
double T_GD,
double ,
60 double ,
double ,
double ,
double )
64 IODEC(static_cast<size_t>(IODEC)),
65 a({ svClockBias, svClockDrift, svClockDriftRate }),
81 svAccuracy(svAccuracy),
82 svHealth(
static_cast<uint8_t
>(svHealth)),
96 LOG_DATA(
" toe {} (Time of ephemeris)",
toe.toGPSweekTow());
99 LOG_DATA(
" A {} [m] (Semi-major axis)", A);
100 auto n_0 = std::sqrt(mu / std::pow(A, 3));
101 LOG_DATA(
" n_0 {} [rad/s] (Computed mean motion)", n_0);
103 LOG_DATA(
" n {} [rad/s] (Corrected mean motion)", n);
108 InsTime transTime = transTime0;
109 LOG_DATA(
" Iterating Time at transmission");
111 double clkDrift = 0.0;
113 for (
size_t i = 0; i < 2; i++)
118 auto t_minus_toc =
static_cast<double>((transTime -
toc).count());
119 LOG_DATA(
" transTime - toc {} [s]", t_minus_toc);
122 double t_k =
static_cast<double>((transTime -
toe).count());
123 LOG_DATA(
" transTime - toe {} [s] (t_k = Time difference from ephemeris reference epoch)", t_k);
126 auto M_k =
M_0 + n * t_k;
127 LOG_DATA(
" M_k {} [s] (Mean anomaly)", M_k);
131 double E_k_old = 0.0;
133 for (
size_t i = 0; std::abs(E_k - E_k_old) > 1
e-13 && i < 10; i++)
136 E_k = M_k +
e * sin(E_k);
140 double dt_r = F *
e *
sqrt_A * std::sin(E_k);
141 LOG_DATA(
" dt_r {} [s] (Relativistic correction term)", dt_r);
144 dt_sv =
a[0] +
a[1] * t_minus_toc +
a[2] * std::pow(t_minus_toc, 2) + dt_r;
149 LOG_DATA(
" dt_sv {} [s] (SV PRN code phase time offset)", dt_sv);
152 clkDrift =
a[1] +
a[2] / 2.0 * t_minus_toc;
155 transTime = transTime0 - std::chrono::duration<double>(dt_sv);
159 return { .transmitTime = transTime, .bias = dt_sv, .drift = clkDrift };
164 Eigen::Vector3d e_pos = Eigen::Vector3d::Zero();
165 Eigen::Vector3d e_vel = Eigen::Vector3d::Zero();
166 Eigen::Vector3d e_accel = Eigen::Vector3d::Zero();
174 LOG_DATA(
" toe {} (Time of ephemeris)",
toe.toGPSweekTow());
177 LOG_DATA(
" A {} [m] (Semi-major axis)", A);
178 auto n_0 = std::sqrt(mu / std::pow(A, 3));
179 LOG_DATA(
" n_0 {} [rad/s] (Computed mean motion)", n_0);
181 LOG_DATA(
" n {} [rad/s] (Corrected mean motion)", n);
187 double t_k =
static_cast<double>((transTime -
toe).count());
188 LOG_DATA(
" t_k {} [s] (Time difference from ephemeris reference epoch)", t_k);
191 auto M_k =
M_0 + n * t_k;
192 LOG_DATA(
" M_k {} [s] (Mean anomaly)", M_k);
195 double E_k_old = 0.0;
197 LOG_DATA(
" E_k {} [rad] (Eccentric anomaly)", E_k);
198 for (
size_t i = 0; std::abs(E_k - E_k_old) > 1
e-13 && i < 10; i++)
201 E_k = E_k + (M_k - E_k +
e * std::sin(E_k)) / (1 -
e * std::cos(E_k));
202 LOG_DATA(
" E_k {} [rad] (Eccentric anomaly)", E_k);
207 auto v_k = std::atan2(std::sqrt(1 -
e *
e) * std::sin(E_k), (std::cos(E_k) -
e));
208 LOG_DATA(
" v_k {} [rad] (True Anomaly (unambiguous quadrant))", v_k);
209 auto Phi_k = v_k +
omega;
210 LOG_DATA(
" Phi_k {} [rad] (Argument of Latitude)", Phi_k);
213 auto delta_u_k =
Cus * std::sin(2 * Phi_k) +
Cuc * std::cos(2 * Phi_k);
214 LOG_DATA(
" delta_u_k {} [rad] (Argument of Latitude Correction)", delta_u_k);
215 auto delta_r_k =
Crs * std::sin(2 * Phi_k) +
Crc * std::cos(2 * Phi_k);
216 LOG_DATA(
" delta_r_k {} [m] (Radius Correction)", delta_r_k);
217 auto delta_i_k =
Cis * std::sin(2 * Phi_k) +
Cic * std::cos(2 * Phi_k);
218 LOG_DATA(
" delta_i_k {} [rad] (Inclination Correction)", delta_i_k);
220 auto u_k = Phi_k + delta_u_k;
221 LOG_DATA(
" u_k {} [rad] (Corrected Argument of Latitude)", u_k);
222 auto r_k = A * (1 -
e * std::cos(E_k)) + delta_r_k;
223 LOG_DATA(
" r_k {} [m] (Corrected Radius)", r_k);
224 auto i_k =
i_0 + delta_i_k +
i_dot * t_k;
225 LOG_DATA(
" i_k {} [rad] (Corrected Inclination)", i_k);
227 auto x_k_op = r_k * std::cos(u_k);
228 LOG_DATA(
" x_k_op {} [m] (Position in orbital plane)", x_k_op);
229 auto y_k_op = r_k * std::sin(u_k);
230 LOG_DATA(
" y_k_op {} [m] (Position in orbital plane)", y_k_op);
233 auto Omega_k =
Omega_0 + (
Omega_dot - Omega_e_dot) * t_k - Omega_e_dot *
static_cast<double>(
toe.toGPSweekTow(
IRNSST).tow);
234 LOG_DATA(
" Omega_k {} [rad] (Corrected longitude of ascending node)", Omega_k);
237 auto x_k = x_k_op * std::cos(Omega_k) - y_k_op * std::cos(i_k) * std::sin(Omega_k);
238 LOG_DATA(
" x_k {} [m] (Earth-fixed x coordinates)", x_k);
240 auto y_k = x_k_op * std::sin(Omega_k) + y_k_op * std::cos(i_k) * std::cos(Omega_k);
241 LOG_DATA(
" y_k {} [m] (Earth-fixed y coordinates)", y_k);
243 auto z_k = y_k_op * std::sin(i_k);
244 LOG_DATA(
" z_k {} [m] (Earth-fixed z coordinates)", z_k);
246 e_pos = Eigen::Vector3d{ x_k, y_k, z_k };
251 auto E_k_dot = n / (1 -
e * std::cos(E_k));
253 auto v_k_dot = E_k_dot * std::sqrt(1 -
e *
e) / (1 -
e * std::cos(E_k));
255 auto i_k_dot =
i_dot + 2 * v_k_dot * (
Cis * std::cos(2 * Phi_k) -
Cic * std::sin(2 * Phi_k));
257 auto u_k_dot = v_k_dot + 2 * v_k_dot * (
Cus * std::cos(2 * Phi_k) -
Cuc * std::sin(2 * Phi_k));
259 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));
261 auto Omega_k_dot =
Omega_dot - Omega_e_dot;
263 auto vx_k_op = r_k_dot * std::cos(u_k) - r_k * u_k_dot * std::sin(u_k);
265 auto vy_k_op = r_k_dot * std::sin(u_k) + r_k * u_k_dot * std::cos(u_k);
267 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)
268 - 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));
270 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)
271 - 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));
273 auto vz_k = vy_k_op * std::sin(i_k) + y_k_op * i_k_dot * std::cos(i_k);
277 e_vel = Eigen::Vector3d{ vx_k, vy_k, vz_k };
285 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))
286 + 2 * vy_k * Omega_e_dot + x_k * std::pow(Omega_e_dot, 2);
288 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))
289 + 2 * vx_k * Omega_e_dot + y_k * std::pow(Omega_e_dot, 2);
291 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));
293 e_accel = Eigen::Vector3d{ ax_k, ay_k, az_k };
297 return { .e_pos = e_pos,
299 .e_accel = e_accel };
IRNSS 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.
Corrections calcClockCorrections(const InsTime &recvTime, double dist, const Frequency &freq) const final
Calculates clock bias and drift of the satellite.
const std::bitset< 2 > svHealth
SV health.
const InsTime toc
Time of Clock.
double calcSatellitePositionVariance() const final
Calculates the Variance of the satellite position in [m^2].
bool isHealthy() const final
Checks whether the signal is healthy.
const double omega
Argument of perigee [rad].
const double i_dot
Rate of inclination angle [rad/s].
const double Crc
Amplitude of the cosine harmonic correction term to the orbit radius [m].
IRNSSEphemeris(const InsTime &toc, const InsTime &toe, const size_t &IODEC, 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, const double &T_GD)
Constructor.
const double e
Eccentricity [-].
const double Omega_dot
Rate of right ascension [rad/s].
const double T_GD
Total Group Delay.
PosVelAccel calcSatelliteData(const InsTime &transTime, Orbit::Calc calc) const final
Calculates position, velocity and acceleration of the satellite at transmission time.
const double delta_n
Mean motion difference from computed value [rad/s].
const double Cic
Amplitude of the cosine harmonic correction term to the angle of inclination [rad].
const double Cus
Amplitude of the sine harmonic correction term to the argument of latitude [rad].
const size_t IODEC
Issue of Data for Ephemeris and Clock.
const double Cuc
Amplitude of the cosine harmonic correction term to the argument of latitude [rad].
const double Cis
Amplitude of the sine harmonic correction term to the angle of inclination [rad].
const double svAccuracy
SV accuracy [m].
const double M_0
Mean anomaly at reference time [rad].
const double sqrt_A
Square root of the semi-major axis [m^1/2].
const std::array< double, 3 > a
const double i_0
Inclination angle at reference time [rad].
const double Omega_0
Longitude of Ascending Node of Orbit Plane at Weekly Epoch [rad].
const InsTime toe
Time of Ephemeris.
const double Crs
Amplitude of the sine harmonic correction term to the orbit radius [m].
static constexpr double R_E
Earth Equatorial Radius [m].
static constexpr double J2
Oblate Earth Gravity Coefficient [-].
static constexpr double MU
Gravitational constant GPS [m³/s²]. See is-gps-200m IS-GPS-200M p. 106.
static constexpr double omega_ie
Earth angular velocity GPS [rad/s]. See is-gps-200m IS-GPS-200M p. 106.
static constexpr double MU
Earth gravitational constant IRNSS [m³/s²].
static constexpr double F
Relativistic constant F for clock corrections [s/√m] (-2*√µ/c²)
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.
@ IRNSSEphemeris
IRNSS Broadcast Ephemeris.
@ IRNSST
Indian Regional Navigation Satellite System Time.
double ratioFreqSquared(Frequency f1, Frequency f2, int8_t num1, int8_t num2)
Calculates the ration of the frequencies squared γ
@ S01
SBAS L1 (1575.42 MHz).
double gpsUraIdx2Val(uint8_t idx)
Converts a GPS URA (user range accuracy) index to it's value.
uint8_t gpsUraVal2Idx(double val)
Converts a GPS URA (user range accuracy) value to it's index.
@ IRNSS
Indian Regional Navigation Satellite System.
Satellite clock corrections.
Satellite Position, Velocity and Acceleration.