GNSS-SDR  0.0.19
An Open Source GNSS Software Defined Receiver
gnss_ephemeris.h
Go to the documentation of this file.
1 /*!
2  * \file gnss_ephemeris.h
3  * \brief Base class for GNSS Ephemeris
4  * \author Carles Fernandez, 2021. cfernandez(at)cttc.es
5  *
6  *
7  * -----------------------------------------------------------------------------
8  *
9  * GNSS-SDR is a Global Navigation Satellite System software-defined receiver.
10  * This file is part of GNSS-SDR.
11  *
12  * Copyright (C) 2010-2021 (see AUTHORS file for a list of contributors)
13  * SPDX-License-Identifier: GPL-3.0-or-later
14  *
15  * -----------------------------------------------------------------------------
16  */
17 
18 
19 #ifndef GNSS_SDR_GNSS_EPHEMERIS_H
20 #define GNSS_SDR_GNSS_EPHEMERIS_H
21 
22 #include <array>
23 #include <cstdint>
24 
25 /*!
26  * \brief Base class for GNSS ephemeris storage
27  */
29 {
30 public:
31  Gnss_Ephemeris() = default;
32 
33  /*!
34  * \brief Sets (\a satClkDrift) and (\a dtr), and returns the clock drift in
35  * seconds according to the User Algorithm for SV Clock Correction
36  * (IS-GPS-200M, 20.3.3.3.3.1, and Galileo OS SIS ICD, 5.1.4).
37  */
38  double sv_clock_drift(double transmitTime);
39 
40  /*!
41  * \brief Computes prediction of the Doppler shift for a given time and receiver's position and velocity.
42  * \f[
43  * f_{d} = - \mathbf{v} \frac{\mathbf{x}^{T}}{\left| \mathbf{x} \right| } \frac{f_{L}}{c}
44  * \f]
45  * where:
46  * \f[
47  * \mathbf{v} = \mathbf{v}_{sat} - \mathbf{v}_{rx}
48  * \f]
49  * \f[
50  * \mathbf{x} = \mathbf{x}_{sat} - \mathbf{x}_{rx}
51  * \f]
52  * \f[
53  * \left| \mathbf{x} \right| = \sqrt{\mathbf{x}\mathbf{x}^{T}}
54  * \f]
55  *
56  * @param[in] rx_time_s Time of Week in seconds
57  * @param[in] lat Receiver's latitude in degrees
58  * @param[in] lon Receiver's longitude in degrees
59  * @param[in] h Receiver's height in meters
60  * @param[in] ve Receiver's velocity in the East direction [m/s]
61  * @param[in] vn Receiver's velocity in the North direction [m/s]
62  * @param[in] vu Receiver's velocity in the Up direction [m/s]
63  * @param[in] band Signal band for which the Doppler will be computed
64  * (1: L1 C/A, E1B, BI1; 2: L2C, BI2; 3: BI3; 5: L5/E5a; 6: E6B; 7: E5b; 8: E5a+E5b)
65  */
66  double predicted_doppler(double rx_time_s, double lat, double lon, double h, double ve, double vn, double vu, int band) const;
67 
68  void satellitePosition(double transmitTime); //!< Computes the ECEF SV coordinates and ECEF velocity
69 
70  uint32_t PRN{}; //!< SV ID
71  double M_0{}; //!< Mean anomaly at reference time [rad]
72  double delta_n{}; //!< Mean motion difference from computed value [rad/sec]
73  double ecc{}; //!< Eccentricity
74  double sqrtA{}; //!< Square root of the semi-major axis [meters^1/2]
75  double OMEGA_0{}; //!< Longitude of ascending node of orbital plane at weekly epoch [rad]
76  double i_0{}; //!< Inclination angle at reference time [rad]
77  double omega{}; //!< Argument of perigee [rad]
78  double OMEGAdot{}; //!< Rate of right ascension [rad/sec]
79  double idot{}; //!< Rate of inclination angle [rad/sec]
80  double Cuc{}; //!< Amplitude of the cosine harmonic correction term to the argument of latitude [rad]
81  double Cus{}; //!< Amplitude of the sine harmonic correction term to the argument of latitude [rad]
82  double Crc{}; //!< Amplitude of the cosine harmonic correction term to the orbit radius [meters]
83  double Crs{}; //!< Amplitude of the sine harmonic correction term to the orbit radius [meters]
84  double Cic{}; //!< Amplitude of the cosine harmonic correction term to the angle of inclination [rad]
85  double Cis{}; //!< Amplitude of the sine harmonic correction term to the angle of inclination [rad]
86  int32_t toe{}; //!< Ephemeris reference time [s]
87 
88  // Clock correction parameters
89  int32_t toc{}; //!< Clock correction data reference Time of Week [sec]
90  double af0{}; //!< SV clock bias correction coefficient [s]
91  double af1{}; //!< SV clock drift correction coefficient [s/s]
92  double af2{}; //!< SV clock drift rate correction coefficient [s/s^2]
93 
94  double satClkDrift{}; //!< SV clock drift
95  double dtr{}; //!< Relativistic clock correction term
96 
97  // Time
98  int32_t WN{}; //!< Week number
99  int32_t tow{}; //!< Time of Week
100 
101  // satellite positions
102  double satpos_X{}; //!< Earth-fixed coordinate x of the satellite [m]. Intersection of the IERS Reference Meridian (IRM) and the plane passing through the origin and normal to the Z-axis.
103  double satpos_Y{}; //!< Earth-fixed coordinate y of the satellite [m]. Completes a right-handed, Earth-Centered, Earth-Fixed orthogonal coordinate system.
104  double satpos_Z{}; //!< Earth-fixed coordinate z of the satellite [m]. The direction of the IERS (International Earth Rotation and Reference Systems Service) Reference Pole (IRP).
105 
106  // Satellite velocity
107  double satvel_X{}; //!< Earth-fixed velocity coordinate x of the satellite [m]
108  double satvel_Y{}; //!< Earth-fixed velocity coordinate y of the satellite [m]
109  double satvel_Z{}; //!< Earth-fixed velocity coordinate z of the satellite [m]
110 
111 protected:
112  char System{}; //!< Character ID of the GNSS system. 'G': GPS. 'E': Galileo. 'B': BeiDou
113 
114 private:
115  void satellitePosVelComputation(double transmitTime, std::array<double, 7>& pos_vel_dtr) const;
116  double check_t(double time) const;
117  double sv_clock_relativistic_term(double transmitTime) const;
118 };
119 
120 #endif // GNSS_SDR_GNSS_EPHEMERIS_H
void satellitePosition(double transmitTime)
Computes the ECEF SV coordinates and ECEF velocity.
double sv_clock_drift(double transmitTime)
Sets (satClkDrift) and (dtr), and returns the clock drift in seconds according to the User Algorithm ...
double i_0
Inclination angle at reference time [rad].
double satvel_Y
Earth-fixed velocity coordinate y of the satellite [m].
double OMEGA_0
Longitude of ascending node of orbital plane at weekly epoch [rad].
double satpos_Z
Earth-fixed coordinate z of the satellite [m]. The direction of the IERS (International Earth Rotatio...
int32_t toc
Clock correction data reference Time of Week [sec].
double satvel_X
Earth-fixed velocity coordinate x of the satellite [m].
double Crc
Amplitude of the cosine harmonic correction term to the orbit radius [meters].
double af2
SV clock drift rate correction coefficient [s/s^2].
uint32_t PRN
SV ID.
double satpos_Y
Earth-fixed coordinate y of the satellite [m]. Completes a right-handed, Earth-Centered, Earth-Fixed orthogonal coordinate system.
double Cis
Amplitude of the sine harmonic correction term to the angle of inclination [rad]. ...
double af0
SV clock bias correction coefficient [s].
double ecc
Eccentricity.
double M_0
Mean anomaly at reference time [rad].
int32_t WN
Week number.
int32_t toe
Ephemeris reference time [s].
int32_t tow
Time of Week.
double Cus
Amplitude of the sine harmonic correction term to the argument of latitude [rad]. ...
double OMEGAdot
Rate of right ascension [rad/sec].
double delta_n
Mean motion difference from computed value [rad/sec].
double satvel_Z
Earth-fixed velocity coordinate z of the satellite [m].
double idot
Rate of inclination angle [rad/sec].
double satpos_X
Earth-fixed coordinate x of the satellite [m]. Intersection of the IERS Reference Meridian (IRM) and ...
Base class for GNSS ephemeris storage.
double Cic
Amplitude of the cosine harmonic correction term to the angle of inclination [rad].
double predicted_doppler(double rx_time_s, double lat, double lon, double h, double ve, double vn, double vu, int band) const
Computes prediction of the Doppler shift for a given time and receiver&#39;s position and velocity...
double Cuc
Amplitude of the cosine harmonic correction term to the argument of latitude [rad].
double omega
Argument of perigee [rad].
double af1
SV clock drift correction coefficient [s/s].
double Crs
Amplitude of the sine harmonic correction term to the orbit radius [meters].
char System
Character ID of the GNSS system. &#39;G&#39;: GPS. &#39;E&#39;: Galileo. &#39;B&#39;: BeiDou.
double satClkDrift
SV clock drift.
double dtr
Relativistic clock correction term.
double sqrtA
Square root of the semi-major axis [meters^1/2].