1
0
mirror of https://github.com/gnss-sdr/gnss-sdr synced 2025-01-28 09:54:51 +00:00

Commits from the Diff file created by Mara Branzanti on August 7th

git-svn-id: https://svn.code.sf.net/p/gnss-sdr/code/trunk@404 64b25241-fba3-4117-9849-534c7e92360d
This commit is contained in:
Javier Arribas 2013-08-16 12:45:30 +00:00
parent d070b4a389
commit 367105f666
8 changed files with 481 additions and 122 deletions

View File

@ -40,6 +40,12 @@
#include <utility> // std::pair
#include "MATH_CONSTANTS.h"
// Physical constants
const double GALILEO_PI = 3.1415926535898; //!< Pi as defined in GALILEO ICD
const double GALILEO_GM = 3.986004418e14; //!< Geocentric gravitational constant[m^3/s^2]
const double GALILEO_OMEGA_EARTH_DOT = 7.2921151467e-5; //!< Mean angular velocity of the Earth [rad/s]
const double GALILEO_C_m_s = 299792458.0; //!< The speed of light, [m/s]
const double GALILEO_F = -4.442807633e-10; //!< Constant, [s/(m)^(1/2)]
// carrier and code frequencies
const double Galileo_E1_FREQ_HZ = 1.57542e9; //!< E1 [Hz]
const double Galileo_E1_CODE_CHIP_RATE_HZ = 1.023e6; //!< Galileo E1 code rate [chips/s]
@ -151,7 +157,7 @@ const std::vector<std::pair<int,int>> Region5_5_bit({{47,1}}); //
const std::vector<std::pair<int,int>> BGD_E1E5a_5_bit({{48,10}}); //
const double BGD_E1E5a_5_LSB = TWO_N32;
const std::vector<std::pair<int,int>> BGD_E1E5b_5_bit({{58,10}}); //
const double BGD_E1E5b_5_LSB = TWO_N35;
const double BGD_E1E5b_5_LSB = TWO_N32;
const std::vector<std::pair<int,int>> E5b_HS_5_bit({{68,2}}); //
const std::vector<std::pair<int,int>> E1B_HS_5_bit({{70,2}}); //
const std::vector<std::pair<int,int>> E5b_DVS_5_bit({{72,1}}); //
@ -182,8 +188,8 @@ const std::vector<std::pair<int,int>>WN_a_7_bit({{11,2}});
const std::vector<std::pair<int,int>>t0a_7_bit({{13,10}});
const double t0a_7_LSB = 600;
const std::vector<std::pair<int,int>>SVID1_7_bit({{23,6}});
const std::vector<std::pair<int,int>>Delta_alpha_7_bit({{29,13}});
const double Delta_alpha_7_LSB = TWO_N9;
const std::vector<std::pair<int,int>>DELTA_A_7_bit({{29,13}});
const double DELTA_A_7_LSB = TWO_N9;
const std::vector<std::pair<int,int>>e_7_bit({{42,11}});
const double e_7_LSB = TWO_N16;
const std::vector<std::pair<int,int>>omega_7_bit({{53,16}});
@ -227,7 +233,7 @@ const std::vector<std::pair<int,int>>M0_9_bit({{23,16}});
const double M0_9_LSB = TWO_N15;
const std::vector<std::pair<int,int>>af0_9_bit({{39,16}});
const double af0_9_LSB = TWO_N19;
const std::vector<std::pair<int,int>>af1_9_bit({{54,13}});
const std::vector<std::pair<int,int>>af1_9_bit({{55,13}});
const double af1_9_LSB = TWO_N38;
const std::vector<std::pair<int,int>>E5b_HS_9_bit({{68,2}});
const std::vector<std::pair<int,int>>E1B_HS_9_bit({{70,2}});

View File

@ -2,7 +2,7 @@
* \file gps_almanac.h
* \brief Interface of a GPS ALMANAC storage
* \author Javier Arribas, 2013. jarribas(at)cttc.es
*
* \author Mara Branzanti 2013. mara.branzanti(at)gmail.com
* -------------------------------------------------------------------------
*
* Copyright (C) 2010-2013 (see AUTHORS file for a list of contributors)
@ -28,38 +28,71 @@
* -------------------------------------------------------------------------
*/
/*
* ToDo: rewrite the class Galileo Almanac (actually it is just a Gps Almanac copy!). Update also the in-line documentation!
*/
#ifndef GNSS_SDR_GALILEO_ALMANAC_H_
#define GNSS_SDR_GALILEO_ALMANAC_H_
/*!
* \brief This class is a storage for the GPS SV ALMANAC data as described in IS-GPS-200E
* \brief This class is a storage for the GALIELO ALMANAC data as described in GALILEO ICD
*
* See http://www.gps.gov/technical/icwg/IS-GPS-200E.pdf Appendix II
*/
/*
* ToDo: Rewrite the class to include all the parameters described in Galileo ICD (this is just a copy of GPS code!)
* See http:http://ec.europa.eu/enterprise/policies/satnav/galileo/files/galileo-os-sis-icd-issue1-revision1_en.pdf paragraph 5.1.10
*/
class Galileo_Almanac
{
public:
unsigned int i_satellite_PRN; //!< SV PRN NUMBER
double d_Delta_i;
double d_Toa; //!< Almanac data reference time of week (Ref. 20.3.3.4.3 IS-GPS-200E) [s]
double d_M_0; //!< Mean Anomaly at Reference Time [semi-circles]
double d_e_eccentricity; //!< Eccentricity [dimensionless]
double d_sqrt_A; //!< Square Root of the Semi-Major Axis [sqrt(m)]
double d_OMEGA0; //!< Longitude of Ascending Node of Orbit Plane at Weekly Epoch [semi-circles]
double d_OMEGA; //!< Argument of Perigee [semi-cicles]
double d_OMEGA_DOT; //!< Rate of Right Ascension [semi-circles/s]
int i_SV_health; // SV Health
double d_A_f0; //!< Coefficient 0 of code phase offset model [s]
double d_A_f1; //!< Coefficient 1 of code phase offset model [s/s]
/*Word type 7: Almanac for SVID1 (1/2), almanac reference time and almanac reference week number*/
int IOD_a_7;
double WN_a_7;
double t0a_7;
int SVID1_7;
double DELTA_A_7;
double e_7;
double omega_7;
double delta_i_7;
double Omega0_7;
double Omega_dot_7;
double M0_7;
/*Word type 8: Almanac for SVID1 (2/2) and SVID2 (1/2)*/
int IOD_a_8;
double af0_8;
double af1_8;
double E5b_HS_8;
double E1B_HS_8;
int SVID2_8;
double DELTA_A_8;
double e_8;
double omega_8;
double delta_i_8;
double Omega0_8;
double Omega_dot_8;
/*Word type 9: Almanac for SVID2 (2/2) and SVID3 (1/2)*/
int IOD_a_9;
double WN_a_9;
double t0a_9;
double M0_9;
double af0_9;
double af1_9;
double E5b_HS_9;
double E1B_HS_9;
int SVID3_9;
double DELTA_A_9;
double e_9;
double omega_9;
double delta_i_9;
/*Word type 10: Almanac for SVID3 (2/2)*/
int IOD_a_10;
double Omega0_10;
double Omega_dot_10;
double M0_10;
double af0_10;
double af1_10;
double E5b_HS_10;
double E1B_HS_10;
/*!
* Default constructor
*/

View File

@ -4,7 +4,7 @@
*
* See http://www.gps.gov/technical/icwg/IS-GPS-200E.pdf Appendix II
* \author Javier Arribas, 2013. jarribas(at)cttc.es
*
* \author Mara Branzanti 2013. mara.branzanti(at)gmail.com
* -------------------------------------------------------------------------
*
* Copyright (C) 2010-2013 (see AUTHORS file for a list of contributors)
@ -32,64 +32,237 @@
#include "galileo_ephemeris.h"
Galileo_Ephemeris::Galileo_Ephemeris()
{
i_satellite_PRN = 0;
d_TOW = 0;
d_Crs = 0;
d_Delta_n = 0;
d_M_0 = 0;
d_Cuc = 0;
d_e_eccentricity = 0;
d_Cus = 0;
d_sqrt_A = 0;
d_Toe = 0;
d_Toc = 0;
d_Cic = 0;
d_OMEGA0 = 0;
d_Cis = 0;
d_i_0 = 0;
d_Crc = 0;
d_OMEGA = 0;
d_OMEGA_DOT = 0;
d_IDOT = 0;
i_code_on_L2 = 0;
i_GPS_week = 0;
b_L2_P_data_flag = false;
i_SV_accuracy = 0;
i_SV_health = 0;
d_TGD = 0; //!< Estimated Group Delay Differential: L1-L2 correction term only for the benefit of "L1 P(Y)" or "L2 P(Y)" s users [s]
d_IODC = 0; //!< Issue of Data, Clock
i_AODO = 0; //!< Age of Data Offset (AODO) term for the navigation message correction table (NMCT) contained in subframe 4 (reference paragraph 20.3.3.5.1.9) [s]
double M0_1 = 0; // Mean anomaly at reference time [semi-circles]
double delta_n_3 = 0; // Mean motion difference from computed value [semi-circles/sec]
double e_1 = 0; // Eccentricity
double A_1 = 0; // Square root of the semi-major axis [metres^1/2]
double OMEGA_0_2 = 0; // Longitude of ascending node of orbital plane at weekly epoch [semi-circles]
double i_0_2 = 0; // Inclination angle at reference time [semi-circles]
double omega_2 = 0; // Argument of perigee [semi-circles]
double OMEGA_dot_3 = 0; // Rate of right ascension [semi-circles/sec]
double iDot_2 = 0; // Rate of inclination angle [semi-circles/sec]
double C_uc_3 = 0; // Amplitude of the cosine harmonic correction term to the argument of latitude [radians]
double C_us_3 = 0; // Amplitude of the sine harmonic correction term to the argument of latitude [radians]
double C_rc_3 = 0; // Amplitude of the cosine harmonic correction term to the orbit radius [meters]
double C_rs_3 = 0; // Amplitude of the sine harmonic correction term to the orbit radius [meters]
double C_ic_4 = 0; // Amplitude of the cosine harmonic correction term to the angle of inclination [radians]
double C_is_4 = 0; // Amplitude of the sine harmonic correction term to the angle of inclination [radians]
double t0e_1 = 0; // Ephemeris reference time [s]
b_fit_interval_flag = false;//!< indicates the curve-fit interval used by the CS (Block II/IIA/IIR/IIR-M/IIF) and SS (Block IIIA) in determining the ephemeris parameters, as follows: 0 = 4 hours, 1 = greater than 4 hours.
d_spare1 = 0;
d_spare2 = 0;
/*Clock correction parameters*/
double t0c_4 = 0; //Clock correction data reference Time of Week [sec]
double af0_4 = 0; //SV clock bias correction coefficient [s]
double af1_4 = 0; //SV clock drift correction coefficient [s/s]
double af2_4 = 0; //SV clock drift rate correction coefficient [s/s^2]
d_A_f0 = 0; //!< Coefficient 0 of code phase offset model [s]
d_A_f1 = 0; //!< Coefficient 1 of code phase offset model [s/s]
d_A_f2 = 0; //!< Coefficient 2 of code phase offset model [s/s^2]
b_integrity_status_flag = false;
b_alert_flag = false; //!< If true, indicates that the SV URA may be worse than indicated in d_SV_accuracy, use that SV at our own risk.
b_antispoofing_flag = false; //!< If true, the AntiSpoofing mode is ON in that SV
/*GST*/
double WN_5 = 0;
double TOW_5 = 0;
}
void Galileo_Ephemeris::satellitePosition(double transmitTime)
{
/*
* ToDo: Compute satellite position at transmit Time
double Galileo_Ephemeris::Galileo_System_Time(double WN, double TOW){
/* GALIELO SYSTEM TIME, ICD 5.1.2
* input parameter:
* WN: The Week Number is an integer counter that gives the sequential week number
from the origin of the Galileo time. It covers 4096 weeks (about 78 years).
Then the counter is reset to zero to cover additional period modulo 4096
TOW: The Time of Week is defined as the number of seconds that have occurred since
the transition from the previous week. The TOW covers an entire week from 0 to
604799 seconds and is reset to zero at the end of each week
WN and TOW are received in page 5
output:
t: it is the transmitted time in Galileo System Time (expressed in seconds)
The GST start epoch shall be 00:00 UT on Sunday 22nd August 1999 (midnight between 21st and 22nd August).
At the start epoch, GST shall be ahead of UTC by thirteen (13)
leap seconds. Since the next leap second was inserted at 01.01.2006, this implies that
as of 01.01.2006 GST is ahead of UTC by fourteen (14) leap seconds.
The epoch denoted in the navigation messages by TOW and WN
will be measured relative to the leading edge of the first chip of the
first code sequence of the first page symbol. The transmission timing of the navigation
message provided through the TOW is synchronised to each satellites version of Galileo System Time (GST).
*
*/
double t=0;
double sec_in_day = 86400;
double day_in_week = 7;
t = WN*sec_in_day*day_in_week+ TOW; // second from the origin of the Galileo time
return t;
}
double Galileo_Ephemeris::sv_clock_drift(double transmitTime){
/* Satellite Time Correction Algorithm, ICD 5.1.4
*
*/
double dt;
dt = transmitTime - t0c_4;
Galileo_satClkDrift = af0_4 + af1_4*dt + (af2_4 * dt)*(af2_4 * dt) + Galileo_dtr;
return Galileo_satClkDrift;
}
// compute the relativistic correction term
double Galileo_Ephemeris::sv_clock_relativistic_term(double transmitTime) //Satellite Time Correction Algorithm, ICD 5.1.4
{
double tk;
double a;
double n;
double n0;
double E;
double E_old;
double dE;
double M;
// Restore semi-major axis
a = A_1*A_1;
n0 = sqrt(GALILEO_GM / (a*a*a));
// Time from ephemeris reference epoch
//tk = check_t(transmitTime - d_Toe); this is tk for GPS; for Galileo it is different
//t = WN_5*86400*7 + TOW_5; //WN_5*86400*7 are the second from the origin of the Galileo time
tk = transmitTime - t0e_1;
// Corrected mean motion
n = n0 + delta_n_3;
// Mean anomaly
M = M0_1 + n * tk;
// Reduce mean anomaly to between 0 and 2pi
M = fmod((M + 2* GALILEO_PI), (2* GALILEO_PI));
// Initial guess of eccentric anomaly
E = M;
// --- Iteratively compute eccentric anomaly ----------------------------
for (int ii = 1; ii<20; ii++)
{
E_old = E;
E = M + e_1 * sin(E);
dE = fmod(E - E_old, 2*GALILEO_PI);
if (fabs(dE) < 1e-12)
{
//Necessary precision is reached, exit from the loop
break;
}
}
// Compute relativistic correction term
Galileo_dtr = GALILEO_F * e_1* A_1 * sin(E);
return Galileo_dtr;
}
void Galileo_Ephemeris::satellitePosition(double transmitTime) //when this function in used, the input must be the transmitted time (t) in second computed by Galileo_System_Time (above function)
{
double tk; // Time from ephemeris reference epoch
//double t; // Galileo System Time (ICD, paragraph 5.1.2)
double a; // Semi-major axis
double n; // Corrected mean motion
double n0; // Computed mean motion
double M; // Mean anomaly
double E; //Eccentric Anomaly (to be solved by iteration)
double E_old;
double dE;
double nu; //True anomaly
double phi; //argument of Latitude
double u; // Correct argument of latitude
double r; // Correct radius
double i;
double Omega;
// Find Galileo satellite's position ----------------------------------------------
// Restore semi-major axis
a = A_1*A_1;
// Computed mean motion
n0 = sqrt(GALILEO_GM / (a*a*a));
// Time from ephemeris reference epoch
//tk = check_t(transmitTime - d_Toe); this is tk for GPS; for Galileo it is different
//t = WN_5*86400*7 + TOW_5; //WN_5*86400*7 are the second from the origin of the Galileo time
tk = transmitTime - t0e_1;
// Corrected mean motion
n = n0 + delta_n_3;
// Mean anomaly
M = M0_1 + n * tk;
// Reduce mean anomaly to between 0 and 2pi
M = fmod((M + 2* GALILEO_PI), (2* GALILEO_PI));
// Initial guess of eccentric anomaly
E = M;
// --- Iteratively compute eccentric anomaly ----------------------------
for (int ii = 1; ii<20; ii++)
{
E_old = E;
E = M + e_1 * sin(E);
dE = fmod(E - E_old, 2*GALILEO_PI);
if (fabs(dE) < 1e-12)
{
//Necessary precision is reached, exit from the loop
break;
}
}
// Compute the true anomaly
double tmp_Y = sqrt(1.0 - e_1 * e_1) * sin(E);
double tmp_X = cos(E) - e_1;
nu = atan2(tmp_Y, tmp_X);
// Compute angle phi (argument of Latitude)
phi = nu + omega_2;
// Reduce phi to between 0 and 2*pi rad
phi = fmod((phi), (2*GALILEO_PI));
// Correct argument of latitude
u = phi + C_uc_3 * cos(2*phi) + C_us_3 * sin(2*phi);
// Correct radius
r = a * (1 - e_1*cos(E)) + C_rc_3 * cos(2*phi) + C_rs_3 * sin(2*phi);
// Correct inclination
i = i_0_2 + iDot_2 * tk + C_ic_4 * cos(2*phi) + C_is_4 * sin(2*phi);
// Compute the angle between the ascending node and the Greenwich meridian
Omega = OMEGA_0_2 + (OMEGA_dot_3 - GALILEO_OMEGA_EARTH_DOT)*tk - GALILEO_OMEGA_EARTH_DOT * t0e_1;
// Reduce to between 0 and 2*pi rad
Omega = fmod((Omega + 2*GALILEO_PI), (2*GALILEO_PI));
// --- Compute satellite coordinates in Earth-fixed coordinates
d_satpos_X = 0;
d_satpos_Y = 0;
d_satpos_Z = 0;
galileo_satpos_X = cos(u) * r * cos(Omega) - sin(u) * r * cos(i) * sin(Omega);
galileo_satpos_Y = cos(u) * r * sin(Omega) + sin(u) * r * cos(i) * cos(Omega); //***********************NOTE: in GALILEO ICD this expression is not correct because it has minus (- sin(u) * r * cos(i) * cos(Omega)) instead of plus
galileo_satpos_Z = sin(u) * r * sin(i);
// Satellite's velocity. Can be useful for Vector Tracking loops
double Omega_dot = OMEGA_dot_3 - GALILEO_OMEGA_EARTH_DOT;
galileo_satvel_X = - Omega_dot * (cos(u) * r + sin(u) * r * cos(i)) + galileo_satpos_X * cos(Omega) - galileo_satpos_Y * cos(i) * sin(Omega);
galileo_satvel_Y = Omega_dot * (cos(u) * r * cos(Omega) - sin(u) * r * cos(i) * sin(Omega)) + galileo_satpos_X * sin(Omega) + galileo_satpos_Y * cos(i) * cos(Omega);
galileo_satvel_Z = galileo_satpos_Y * sin(i);
// Satellite's velocity.
d_satvel_X =0;
d_satvel_Y = 0;
d_satvel_Z = 0;
}

View File

@ -1,8 +1,8 @@
/*!
* \file gps_navigation_message.h
* \file galileo_navigation_message.h
* \brief Interface of a GPS EPHEMERIS storage
* \author Javier Arribas, 2013. jarribas(at)cttc.es
*
* \author Mara Branzanti 2013. mara.branzanti(at)gmail.com
* -------------------------------------------------------------------------
*
* Copyright (C) 2010-2013 (see AUTHORS file for a list of contributors)
@ -36,22 +36,64 @@
#include <map>
#include "boost/assign.hpp"
#include <boost/serialization/nvp.hpp>
#include "Galileo_E1.h"
/*!
* \brief This class is a storage and orbital model functions for the Galileo SV ephemeris data as described in Galileo ICD
* \brief This class is a storage and orbital model functions for the Galileo SV ephemeris data as described in Galileo ICD paragraph 5.1.1
*
*/
/*
* ToDo: Rewrite the class to include all the parameters described in Galileo ICD (this is just a copy of GPS ephemeris!)
*/
class Galileo_Ephemeris
{
private:
public:
unsigned int i_satellite_PRN; // SV PRN NUMBER
/*Galileo ephemeris are 16 parameters and here are reported following the ICD order, paragraph 5.1.1.
The number in the name after underscore (_1, _2, _3 and so on) refers to the page were we can find that parameter */
double M0_1; // Mean anomaly at reference time [semi-circles]
double delta_n_3; // Mean motion difference from computed value [semi-circles/sec]
double e_1; // Eccentricity
double A_1; // Square root of the semi-major axis [metres^1/2]
double OMEGA_0_2; // Longitude of ascending node of orbital plane at weekly epoch [semi-circles]
double i_0_2; // Inclination angle at reference time [semi-circles]
double omega_2; // Argument of perigee [semi-circles]
double OMEGA_dot_3; // Rate of right ascension [semi-circles/sec]
double iDot_2; // Rate of inclination angle [semi-circles/sec]
double C_uc_3; // Amplitude of the cosine harmonic correction term to the argument of latitude [radians]
double C_us_3; // Amplitude of the sine harmonic correction term to the argument of latitude [radians]
double C_rc_3; // Amplitude of the cosine harmonic correction term to the orbit radius [meters]
double C_rs_3; // Amplitude of the sine harmonic correction term to the orbit radius [meters]
double C_ic_4; // Amplitude of the cosine harmonic correction term to the angle of inclination [radians]
double C_is_4; // Amplitude of the sine harmonic correction term to the angle of inclination [radians]
double t0e_1; // Ephemeris reference time [s]
/*Clock correction parameters*/
double t0c_4; //Clock correction data reference Time of Week [sec]
double af0_4; //SV clock bias correction coefficient [s]
double af1_4; //SV clock drift correction coefficient [s/s]
double af2_4; //SV clock drift rate correction coefficient [s/s^2]
/*GST*/
double WN_5; //Week number
double TOW_5; //Time of Week
double Galileo_satClkDrift;
double Galileo_dtr; // relativistic clock correction term
// satellite positions
double galileo_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.
double galileo_satpos_Y; //!< Earth-fixed coordinate y of the satellite [m]. Completes a right-handed, Earth-Centered, Earth-Fixed orthogonal coordinate system.
double galileo_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).
// Satellite velocity
double galileo_satvel_X; //!< Earth-fixed velocity coordinate x of the satellite [m]
double galileo_satvel_Y; //!< Earth-fixed velocity coordinate y of the satellite [m]
double galileo_satvel_Z; //!< Earth-fixed velocity coordinate z of the satellite [m]
/*The following parameters refers to GPS
unsigned int i_satellite_PRN; // SV PRN NUMBER
double d_TOW; //!< Time of GPS Week of the ephemeris set (taken from subframes TOW) [s]
double d_Crs; //!< Amplitude of the Sine Harmonic Correction Term to the Orbit Radius [m]
double d_Delta_n; //!< Mean Motion Difference From Computed Value [semi-circles/s]
@ -113,9 +155,9 @@ public:
template<class Archive>
/*
* \brief Serialize is a boost standard method to be called by the boost XML serialization. Here is used to save the ephemeris data on disk file.
*/
\\brief Serialize is a boost standard method to be called by the boost XML serialization. Here is used to save the ephemeris data on disk file.
void serialize(Archive& archive, const unsigned int version)
{
using boost::serialization::make_nvp;
@ -161,15 +203,19 @@ public:
archive & make_nvp("b_antispoofing_flag",b_antispoofing_flag); //!< If true, the AntiSpoofing mode is ON in that SV
}
/*!
* \brief Compute the ECEF SV coordinates and ECEF velocity
* [Insert here the reference in Galileo ICD]
\\brief Compute the ECEF SV coordinates and ECEF velocity
\\http://ec.europa.eu/enterprise/policies/satnav/galileo/open-service/
*/
void satellitePosition(double transmitTime);
/*!
* Default constructor
*/
double Galileo_System_Time(double WN, double TOW); // Galileo System Time (GST), ICD paragraph 5.1.2
double sv_clock_drift(double transmitTime); //Satellite Time Correction Algorithm, ICD 5.1.4
double sv_clock_relativistic_term(double transmitTime); //Satellite Time Correction Algorithm, ICD 5.1.4
//Default constructor
Galileo_Ephemeris();
};

View File

@ -2,7 +2,7 @@
* \file gps_iono.h
* \brief Interface of a GPS IONOSPHERIC MODEL storage
* \author Javier Arribas, 2013. jarribas(at)cttc.es
*
* \author Mara Branzanti 2013. mara.branzanti(at)gmail.com
* -------------------------------------------------------------------------
*
* Copyright (C) 2010-2013 (see AUTHORS file for a list of contributors)
@ -34,30 +34,41 @@
/*!
* \brief This class is a storage for the GALILEO IONOSPHERIC data as described in Galileo ICD
* \brief This class is a storage for the GALILEO IONOSPHERIC data as described in Galileo ICD paragraph 5.1.6
*
* See [Update ref]
*/
/*
* ToDo: Rewrite the class to include all the parameters described in Galileo ICD (this is just a copy of GPS code!)
* See http://ec.europa.eu/enterprise/policies/satnav/galileo/files/galileo-os-sis-icd-issue1-revision1_en.pdf
*/
class Galileo_Iono
{
private:
public:
// valid flag
// valid flag
bool valid;
// Ionospheric parameters
double d_alpha0; //!< Coefficient 0 of a cubic equation representing the amplitude of the vertical delay [s]
/*Ionospheric correction*/
/*Az*/
double ai0_5; //Effective Ionisation Level 1st order parameter [sfu]
double ai1_5; //Effective Ionisation Level 2st order parameter [sfu/degree]
double ai2_5; //Effective Ionisation Level 3st order parameter [sfu/degree]
/*Ionospheric disturbance flag*/
bool Region1_flag_5; // Ionospheric Disturbance Flag for region 1
bool Region2_flag_5; // Ionospheric Disturbance Flag for region 2
bool Region3_flag_5; // Ionospheric Disturbance Flag for region 3
bool Region4_flag_5; // Ionospheric Disturbance Flag for region 4
bool Region5_flag_5; // Ionospheric Disturbance Flag for region 5
// Ionospheric parameters GPS
/*double d_alpha0; //!< Coefficient 0 of a cubic equation representing the amplitude of the vertical delay [s]
double d_alpha1; //!< Coefficient 1 of a cubic equation representing the amplitude of the vertical delay [s/semi-circle]
double d_alpha2; //!< Coefficient 2 of a cubic equation representing the amplitude of the vertical delay [s(semi-circle)^2]
double d_alpha3; //!< Coefficient 3 of a cubic equation representing the amplitude of the vertical delay [s(semi-circle)^3]
double d_beta0; //!< Coefficient 0 of a cubic equation representing the period of the model [s]
double d_beta1; //!< Coefficient 1 of a cubic equation representing the period of the model [s/semi-circle]
double d_beta2; //!< Coefficient 2 of a cubic equation representing the period of the model [s(semi-circle)^2]
double d_beta3; //!< Coefficient 3 of a cubic equation representing the period of the model [s(semi-circle)^3]
double d_beta3; //!< Coefficient 3 of a cubic equation representing the period of the model [s(semi-circle)^3]*/
/*!
* Default constructor
*/

View File

@ -769,8 +769,8 @@ int Galileo_Navigation_Message::page_jk_decoder(char *data_jk)
SVID1_7= (double)read_navigation_unsigned(data_jk_bits, SVID1_7_bit);
std::cout << "SVID1_7= " << SVID1_7 << std::endl;
Delta_alpha_7= (double)read_navigation_unsigned(data_jk_bits, Delta_alpha_7_bit);
Delta_alpha_7= Delta_alpha_7 * Delta_alpha_7_LSB;
Delta_alpha_7= (double)read_navigation_unsigned(data_jk_bits, DELTA_A_7_bit);
Delta_alpha_7= Delta_alpha_7 * DELTA_A_7_LSB;
std::cout << "Delta_alpha_7= " << Delta_alpha_7 << std::endl;
e_7= (double)read_navigation_unsigned(data_jk_bits, e_7_bit);

View File

@ -32,15 +32,84 @@
Galileo_Utc_Model::Galileo_Utc_Model()
{
valid = false;
d_A1 = 0;
d_A0 = 0;
d_t_OT = 0;
i_WN_T = 0;
d_DeltaT_LS = 0;
i_WN_LSF = 0;
i_DN = 0;
d_DeltaT_LSF = 0;
//valid = false;
/*Word type 6: GST-UTC conversion parameters*/
A0_6 = 0;
A1_6 = 0;
Delta_tLS_6 = 0;
t0t_6 = 0;
WNot_6 = 0;
WN_LSF_6 = 0;
DN_6 = 0;
Delta_tLSF_6 = 0;
}
double Galileo_Utc_Model::GST_to_UTC_time(double t_e, int WN)
{
double t_Utc;
double t_Utc_daytime;
double Delta_t_Utc = Delta_tLS_6 + A0_6 + A1_6 * (t_e - t0t_6 + 604800 * (double)(WN - WNot_6));
// Determine if the effectivity time of the leap second event is in the past
int weeksToLeapSecondEvent = WN_LSF_6 - WN;
if ((weeksToLeapSecondEvent) >= 0) // is not in the past
{
//Detect if the effectivity time and user's time is within six hours = 6 * 60 *60 = 21600 s
int secondOfLeapSecondEvent = DN_6 * 24 * 60 * 60;
if (weeksToLeapSecondEvent > 0)
{
t_Utc_daytime = fmod(t_e - Delta_t_Utc, 86400);
}
else //we are in the same week than the leap second event
{
if (abs(t_e - secondOfLeapSecondEvent) > 21600)
{
/* 5.1.7a
* Whenever the leap second adjusted time indicated by the WN_LSF and the DN values
* is not in the past (relative to the user's present time), and the user's
* present time does not fall in the time span which starts at six hours prior
* to the effective time and ends at six hours after the effective time,
* the GST/Utc relationship is given by
*/
t_Utc_daytime = fmod(t_e - Delta_t_Utc, 86400);
}
else
{
/* 5.1.7b
* Whenever the user's current time falls within the time span of six hours
* prior to the leap second adjustment to six hours after the adjustment time, ,
* the effective time is computed according to the following equations:
*/
int W = fmod(t_e - Delta_t_Utc - 43200, 86400) + 43200;
t_Utc_daytime = fmod(W, 86400 + Delta_tLSF_6 - Delta_tLS_6);
//implement something to handle a leap second event!
}
if ( (t_e - secondOfLeapSecondEvent) > 21600)
{
Delta_t_Utc = Delta_tLSF_6 + A0_6 + A1_6 * (t_e - t0t_6 + 604800*(double)(WN - WNot_6));
t_Utc_daytime = fmod(t_e - Delta_t_Utc, 86400);
}
}
}
else // the effectivity time is in the past
{
/* 5.1.7c
* Whenever the leap second adjustment time, as indicated by the WN_LSF and DN values,
* is in the past (relative to the users current time) and the users present time does not
* fall in the time span which starts six hours prior to the leap second adjustment time and
* ends six hours after the adjustment time, the effective time is computed according to
* the following equation:
*/
Delta_t_Utc = Delta_tLSF_6 + A0_6 + A1_6 * (t_e - t0t_6 + 604800 * (double)(WN - WNot_6));
t_Utc_daytime = fmod(t_e - Delta_t_Utc, 86400);
}
double secondsOfWeekBeforeToday = 43200 * floor(t_e / 43200);
t_Utc = secondsOfWeekBeforeToday + t_Utc_daytime;
return t_Utc;
}

View File

@ -3,6 +3,7 @@
* \brief Interface of a GPS UTC MODEL storage
* \author Javier Arribas, 2013. jarribas(at)cttc.es
*
*
* -------------------------------------------------------------------------
*
* Copyright (C) 2010-2013 (see AUTHORS file for a list of contributors)
@ -32,21 +33,37 @@
#ifndef GNSS_SDR_GALILEO_UTC_MODEL_H_
#define GNSS_SDR_GALILEO_UTC_MODEL_H_
#include "Galileo_E1.h"
/*!
* \brief This class is a storage for the GALILEO UTC MODEL data as described in Galileo ICD
*
* See [Update ref]
*/
/*
* ToDo: Rewrite the class to include all the parameters described in Galileo ICD (this is just a copy of GPS code!)
* http://ec.europa.eu/enterprise/policies/satnav/galileo/files/galileo-os-sis-icd-issue1-revision1_en.pdf
* paragraph 5.1.7
*/
class Galileo_Utc_Model
{
public:
bool valid;
//bool valid;
/*Word type 6: GST-UTC conversion parameters*/
double A0_6;
double A1_6;
double Delta_tLS_6;
double t0t_6;
double WNot_6;
double WN_LSF_6;
double DN_6;
double Delta_tLSF_6;
//double TOW_6;
/*GST*/
//double WN_5; //Week number
//double TOW_5; //Time of Week
//this is from gps
/*bool valid;
// UTC parameters
double d_A1; //!< 1st order term of a model that relates GPS and UTC time (ref. 20.3.3.5.2.4 IS-GPS-200E) [s/s]
double d_A0; //!< Constant of a model that relates GPS and UTC time (ref. 20.3.3.5.2.4 IS-GPS-200E) [s]
@ -56,10 +73,14 @@ public:
int i_WN_LSF; //!< Week number at the end of which the leap second becomes effective [weeks]
int i_DN; //!< Day number (DN) at the end of which the leap second becomes effective [days]
double d_DeltaT_LSF; //!< Scheduled future or recent past (relative to NAV message upload) value of the delta time due to leap seconds [s]
*/
/*!
double GST_to_UTC_time(double t_e, int WN);
/*!
* Default constructor
*/
Galileo_Utc_Model();
};