GNSS-SDR 0.0.21
An Open Source GNSS Software Defined Receiver
Loading...
Searching...
No Matches
gnss_almanac.h
Go to the documentation of this file.
1/*!
2 * \file gnss_almanac.h
3 * \brief Base class for GNSS almanac storage
4 * \author Carles Fernandez, 2021. cfernandez(at)cttc.es
5 *
6 * -----------------------------------------------------------------------------
7 *
8 * GNSS-SDR is a Global Navigation Satellite System software-defined receiver.
9 * This file is part of GNSS-SDR.
10 *
11 * Copyright (C) 2010-2021 (see AUTHORS file for a list of contributors)
12 * SPDX-License-Identifier: GPL-3.0-or-later
13 *
14 * -----------------------------------------------------------------------------
15 */
16
17
18#ifndef GNSS_SDR_GNSS_ALMANAC_H
19#define GNSS_SDR_GNSS_ALMANAC_H
20
21#include <array>
22#include <cstdint>
23
24/** \addtogroup Core
25 * \{ */
26/** \addtogroup System_Parameters
27 * \{ */
28
29
30/*!
31 * \brief Base class for GNSS almanac storage
32 */
34{
35public:
36 /*!
37 * Default constructor
38 */
39 Gnss_Almanac() = default;
40
41 /*!
42 * \brief Computes prediction of the Doppler shift for a given time and receiver's position and velocity.
43 * \f[
44 * f_{d} = - \mathbf{v} \frac{\mathbf{x}^{T}}{\left| \mathbf{x} \right| } \frac{f_{L}}{c}
45 * \f]
46 * where:
47 * \f[
48 * \mathbf{v} = \mathbf{v}_{sat} - \mathbf{v}_{rx}
49 * \f]
50 * \f[
51 * \mathbf{x} = \mathbf{x}_{sat} - \mathbf{x}_{rx}
52 * \f]
53 * \f[
54 * \left| \mathbf{x} \right| = \sqrt{\mathbf{x}\mathbf{x}^{T}}
55 * \f]
56 *
57 * @param[in] rx_time_s Time of Week in seconds
58 * @param[in] lat Receiver's latitude in degrees
59 * @param[in] lon Receiver's longitude in degrees
60 * @param[in] h Receiver's height in meters
61 * @param[in] ve Receiver's velocity in the East direction [m/s]
62 * @param[in] vn Receiver's velocity in the North direction [m/s]
63 * @param[in] vu Receiver's velocity in the Up direction [m/s]
64 * @param[in] band Signal band for which the Doppler will be computed
65 * (1: L1 C/A, E1B, BI1; 2: L2C, BI2; 3: BI3; 5: L5/E5a; 6: E6B; 7: E5b; 8: E5a+E5b)
66 */
67 double predicted_doppler(double rx_time_s,
68 double lat,
69 double lon,
70 double h,
71 double ve,
72 double vn,
73 double vu,
74 int band) const;
75
76 /*!
77 * \brief Computes satellite Position and Velocity, in ECEF, for a given time (expressed in seconds of week)
78 */
79 void satellitePosVelComputation(double transmitTime, std::array<double, 7>& pos_vel_dtr) const;
80
81 uint32_t PRN{}; //!< SV PRN NUMBER
82 double delta_i{}; //!< Inclination Angle at Reference Time (relative to i_0 = 0.30 semi-circles)
83 int32_t toa{}; //!< Almanac data reference time of week [s]
84 int32_t WNa{}; //!< Almanac week number
85 double M_0{}; //!< Mean Anomaly at Reference Time [semi-circles]
86 double ecc{}; //!< Eccentricity [dimensionless]
87 double sqrtA{}; //!< Square Root of the Semi-Major Axis [sqrt(m)]
88 double OMEGA_0{}; //!< Longitude of Ascending Node of Orbit Plane at Weekly Epoch [semi-circles]
89 double omega{}; //!< Argument of Perigee [semi-cicles]
90 double OMEGAdot{}; //!< Rate of Right Ascension [semi-circles/s]
91 double af0{}; //!< Coefficient 0 of code phase offset model [s]
92 double af1{}; //!< Coefficient 1 of code phase offset model [s/s]
93
94protected:
95 char System{}; //!< Character ID of the GNSS system. 'G': GPS. 'E': Galileo. 'C': BeiDou
96private:
97 double check_t(double time) const;
98};
99
100
101/** \} */
102/** \} */
103#endif // GNSS_SDR_GNSS_ALMANAC_H
uint32_t PRN
SV PRN NUMBER.
double OMEGAdot
Rate of Right Ascension [semi-circles/s].
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's position and velocity.
double OMEGA_0
Longitude of Ascending Node of Orbit Plane at Weekly Epoch [semi-circles].
double M_0
Mean Anomaly at Reference Time [semi-circles].
double sqrtA
Square Root of the Semi-Major Axis [sqrt(m)].
char System
Character ID of the GNSS system. 'G': GPS. 'E': Galileo. 'C': BeiDou.
double af0
Coefficient 0 of code phase offset model [s].
int32_t WNa
Almanac week number.
Gnss_Almanac()=default
double omega
Argument of Perigee [semi-cicles].
int32_t toa
Almanac data reference time of week [s].
double af1
Coefficient 1 of code phase offset model [s/s].
double ecc
Eccentricity [dimensionless].
double delta_i
Inclination Angle at Reference Time (relative to i_0 = 0.30 semi-circles).
void satellitePosVelComputation(double transmitTime, std::array< double, 7 > &pos_vel_dtr) const
Computes satellite Position and Velocity, in ECEF, for a given time (expressed in seconds of week).