satellite-sgp4-mobility-model.cc
Go to the documentation of this file.
1 /* -*- Mode:C++; c-file-style:"gnu"; indent-tabs-mode:nil; -*- */
2 /*
3  * Copyright (c) 2016 INESC TEC
4  * Copyright (c) 2021 CNES
5  *
6  * This program is free software; you can redistribute it and/or modify
7  * it under the terms of the GNU General Public License version 2 as
8  * published by the Free Software Foundation;
9  *
10  * This program is distributed in the hope that it will be useful,
11  * but WITHOUT ANY WARRANTY; without even the implied warranty of
12  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13  * GNU General Public License for more details.
14  *
15  * You should have received a copy of the GNU General Public License
16  * along with this program; if not, write to the Free Software
17  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
18  *
19  * Code from https://gitlab.inesctec.pt/pmms/ns3-satellite
20  *
21  * Author: Pedro Silva <pmms@inesctec.pt>
22  * Author: Bastien Tauran <bastien.tauran@viveris.fr>
23  */
24 
26 
27 #include "vector-extensions.h"
28 
29 #include <ns3/boolean.h>
30 #include <ns3/log.h>
31 #include <ns3/simulator.h>
32 #include <ns3/string.h>
33 
34 NS_LOG_COMPONENT_DEFINE("sat-sgp4-mobility-model");
35 
36 namespace ns3
37 {
38 
39 NS_OBJECT_ENSURE_REGISTERED(SatSGP4MobilityModel);
40 
41 const gravconsttype SatSGP4MobilityModel::WGeoSys = wgs72; // recommended for SGP4/SDP4
42 const uint32_t SatSGP4MobilityModel::TleSatInfoWidth = 69;
43 
44 TypeId
46 {
47  static TypeId tid =
48  TypeId("ns3::SatSGP4MobilityModel")
49  .SetParent<SatMobilityModel>()
50  .AddConstructor<SatSGP4MobilityModel>()
51  .AddAttribute("StartDateStr",
52  "Absolute start time of simulation (UTC)",
53  StringValue("1992-01-01 00:00:00"),
54  MakeStringAccessor(&SatSGP4MobilityModel::m_startStr),
55  MakeStringChecker())
56  .AddAttribute("UpdatePositionEachRequest",
57  "Compute position each time a packet is transmitted",
58  BooleanValue(true),
60  MakeBooleanChecker())
61  .AddAttribute(
62  "UpdatePositionPeriod",
63  "Period of satellite position refresh, if UpdatePositionEachRequest is false",
64  TimeValue(Seconds(1)),
66  MakeTimeChecker());
67  return tid;
68 }
69 
70 TypeId
72 {
73  return GetTypeId();
74 }
75 
77  : m_startStr("1992-01-01 00:00:00"),
78  m_timeLastUpdate(Time::Min())
79 {
80  NS_LOG_FUNCTION(this);
81 
82  ObjectBase::ConstructSelf(AttributeConstructionList());
84 }
85 
87 {
88 }
89 
92 {
93  NS_LOG_FUNCTION(this);
94 
95  return m_start;
96 }
97 
98 void
100 {
101  NS_LOG_FUNCTION(this << t);
102 
103  m_start = t;
104 }
105 
106 Vector3D
108 {
109  NS_LOG_FUNCTION(this);
110 
111  JulianDate cur = m_start + Simulator::Now();
112 
113  double r[3], v[3];
114  double delta = (cur - GetTleEpoch()).GetMinutes();
115 
116  if (!IsInitialized())
117  return Vector3D();
118 
119  sgp4(WGeoSys, m_sgp4_record, delta, r, v);
120 
121  if (m_sgp4_record.error != 0)
122  return Vector3D();
123 
124  // velocity vector is in km/s so it needs to be converted to m/s
125  return 1000 * rvTemeTovItrf(Vector3D(r[0], r[1], r[2]), Vector3D(v[0], v[1], v[2]), cur);
126 }
127 
128 Vector
130 {
131  NS_LOG_FUNCTION(this);
132 
133  return DoGetGeoPosition().ToVector();
134 }
135 
136 void
137 SatSGP4MobilityModel::DoSetPosition(const Vector& position)
138 {
139  NS_LOG_FUNCTION(this << position);
140 
141  DoSetGeoPosition(GeoCoordinate(position));
142 }
143 
146 {
147  NS_LOG_FUNCTION(this);
148 
149  if ((m_updatePositionEachRequest == false) &&
150  (Simulator::Now() < m_timeLastUpdate + m_updatePositionPeriod))
151  {
152  return m_lastPosition;
153  }
154 
155  m_timeLastUpdate = Simulator::Now();
157 
158  double r[3], v[3];
159  double delta = (cur - GetTleEpoch()).GetMinutes();
160 
161  if (!IsInitialized())
162  return Vector3D();
163 
164  sgp4(WGeoSys, m_sgp4_record, delta, r, v);
165 
166  if (m_sgp4_record.error != 0)
167  return Vector3D();
168 
169  // vector r is in km so it needs to be converted to meters
170  m_lastPosition = rTemeTorItrf(Vector3D(r[0], r[1], r[2]), cur) * 1000;
171 
172  return m_lastPosition;
173 }
174 
175 void
177 {
178  NS_LOG_FUNCTION(this << position);
179 
180  m_lastPosition = position;
182 }
183 
184 bool
186 {
187  NS_LOG_FUNCTION(this);
188 
189  return ((m_sgp4_record.jdsatepoch > 0) && (m_tle1 != "") && (m_tle2 != ""));
190 }
191 
194 {
195  NS_LOG_FUNCTION(this);
196 
197  if (IsInitialized())
199 
200  return JulianDate();
201 }
202 
203 void
204 SatSGP4MobilityModel::SetTleInfo(const std::string& tle)
205 {
206  NS_LOG_FUNCTION(this << tle);
207 
208  uint32_t delimPos = tle.find("\n");
209  const std::string line1 = tle.substr(0, delimPos);
210  const std::string line2 = tle.substr(delimPos + 1, tle.size() - delimPos);
211 
212  double start, stop, delta;
213  char l1[TleSatInfoWidth + 1], l2[TleSatInfoWidth + 1];
214  double r[3], v[3];
215 
216  NS_ASSERT_MSG(line1.size() == TleSatInfoWidth && line2.size() == TleSatInfoWidth,
217  "Two-Line Element info lines must be of length" << TleSatInfoWidth << "!");
218 
219  m_tle1 = std::string(line1.c_str());
220  m_tle2 = std::string(line2.c_str());
221 
222  memcpy(l1, line1.c_str(), sizeof(l1));
223  memcpy(l2, line2.c_str(), sizeof(l2));
224 
225  // 'c' => catalog mode run
226  // 'e' => epoch time (relative to TLE lines)
227  // 'i' => improved mode of operation
228  twoline2rv(l1, l2, 'c', 'e', 'i', WGeoSys, start, stop, delta, m_sgp4_record);
229 
230  // call propagator to check if it has been properly initialized
231  sgp4(WGeoSys, m_sgp4_record, 0, r, v);
232 
233  if (m_start > GetTleEpoch() + Years(1))
234  {
235  NS_FATAL_ERROR("Simulation start date "
236  << (m_start - GetTleEpoch()).GetDays()
237  << " after TLE epoch. TLE are only valid for a few days");
238  }
239 
240  if (m_sgp4_record.error != 0)
241  {
242  NS_FATAL_ERROR("Error while loading TLE file");
243  }
244 
246 }
247 
248 Vector3D
249 SatSGP4MobilityModel::rTemeTorItrf(const Vector3D& rteme, const JulianDate& t)
250 {
251  Matrix pmt = PefToItrf(t); // PEF->ITRF matrix transposed
252  Matrix tmt = TemeToPef(t); // TEME->PEF matrix
253 
254  return pmt * (tmt * rteme);
255 }
256 
257 Vector3D
259  const Vector3D& vteme,
260  const JulianDate& t)
261 {
262  Matrix pmt = PefToItrf(t); // PEF->ITRF matrix transposed
263  Matrix tmt = TemeToPef(t); // TEME->PEF matrix
264  Vector3D w(0.0, 0.0, t.GetOmegaEarth());
265 
266  return pmt * ((tmt * vteme) - CrossProduct(w, tmt * rteme));
267 }
268 
271 {
272  std::pair<double, double> eop = t.GetPolarMotion();
273 
274  const double &xp = eop.first, &yp = eop.second;
275  const double cosxp = cos(xp), cosyp = cos(yp);
276  const double sinxp = sin(xp), sinyp = sin(yp);
277 
278  // [from AIAA-2006-6753 Report, Page 32, Appendix C - TEME Coordinate System]
279  //
280  // Matrix(ITRF<->PEF) = ROT1(yp)*ROT2(xp) [using c for cos, and s for sin]
281  //
282  // | 1 0 0 |*| c(xp) 0 -s(xp) |=| c(xp) 0 -s(xp) |
283  // | 0 c(yp) s(yp) | | 0 1 0 | | s(yp)*s(xp) c(yp) s(yp)*c(xp) |
284  // | 0 -s(yp) c(yp) | | s(xp) 0 c(xp) | | c(yp)*s(xp) -s(yp) c(yp)*c(xp) |
285  //
286 
287  // we return the transpose because it is what's needed
288  return Matrix(cosxp,
289  sinyp * sinxp,
290  cosyp * sinxp,
291  0,
292  cosyp,
293  -sinyp,
294  -sinxp,
295  sinyp * cosxp,
296  cosyp * cosxp);
297 }
298 
301 {
302  const double gmst = t.GetGmst();
303  const double cosg = cos(gmst), sing = sin(gmst);
304 
305  // [from AIAA-2006-6753 Report, Page 32, Appendix C - TEME Coordinate System]
306  //
307  // rPEF = ROT3(gmst)*rTEME
308  //
309  // | cos(gmst) sin(gmst) 0 |
310  // | -sin(gmst) cos(gmst) 0 |
311  // | 0 0 1 |
312  //
313 
314  return Matrix(cosg, sing, 0, -sing, cosg, 0, 0, 0, 1);
315 }
316 
318  double c01,
319  double c02,
320  double c10,
321  double c11,
322  double c12,
323  double c20,
324  double c21,
325  double c22)
326 {
327  m[0][0] = c00;
328  m[0][1] = c01;
329  m[0][2] = c02;
330  m[1][0] = c10;
331  m[1][1] = c11;
332  m[1][2] = c12;
333  m[2][0] = c20;
334  m[2][1] = c21;
335  m[2][2] = c22;
336 }
337 
338 Vector3D
340 {
341  return Vector3D(m[0][0] * v.x + m[0][1] * v.y + m[0][2] * v.z,
342  m[1][0] * v.x + m[1][1] * v.y + m[1][2] * v.z,
343  m[2][0] * v.x + m[2][1] * v.y + m[2][2] * v.z);
344 }
345 
346 } // namespace ns3
GeoCoordinate class is used to store and operate with geodetic coordinates.
Vector ToVector() const
Converts Geodetic coordinates to Cartesian coordinates.
Class to handle Julian days and provide respective Earth Orientation Parameters (EOP).
Definition: julian-date.h:78
double GetGmst(void) const
Retrieve the Greenwich Mean Sidereal Time.
Definition: julian-date.cc:183
double GetOmegaEarth(void) const
Retrieve Earth's angular velocity.
Definition: julian-date.cc:170
std::pair< double, double > GetPolarMotion(void) const
Retrieve the polar motion cofficients measured/predicted.
Definition: julian-date.cc:157
Keep track of the current position and velocity of an object in satellite network.
void NotifyGeoCourseChange(void) const
std::string m_startStr
Simulation absolute start time in string format.
virtual void DoSetGeoPosition(const GeoCoordinate &position)
static const gravconsttype WGeoSys
World Geodetic System (WGS) constants to be used by SGP4/SDP4 models.
JulianDate GetStartTime() const
Get the time instant considered as the simulation start.
virtual void DoSetPosition(const Vector &position)
static const uint32_t TleSatInfoWidth
Satellite's information line size defined by TLE data format.
GeoCoordinate m_lastPosition
Last saved satellite position.
void SetStartTime(const JulianDate &t)
Set the time instant considered as the simulation start.
static Matrix PefToItrf(const JulianDate &t)
Retrieve the matrix for converting from PEF to ITRF at a given time.
Time m_updatePositionPeriod
Period of satellite position refresh, if UpdatePositionEachRequest is false.
JulianDate GetTleEpoch(void) const
Retrieve the TLE epoch time.
bool m_updatePositionEachRequest
Compute position each time a packet is transmitted.
elsetrec m_sgp4_record
SGP4/SDP4 record.
static Vector3D rvTemeTovItrf(const Vector3D &rteme, const Vector3D &vteme, const JulianDate &t)
Retrieve the satellite's velocity vector in ITRF coordinates.
static Matrix TemeToPef(const JulianDate &t)
Retrieve the matrix for converting from TEME to PEF at a given time.
Time m_timeLastUpdate
Last position update time.
std::string m_tle2
satellite's TLE data.
static TypeId GetTypeId(void)
Get the type ID.
static Vector3D rTemeTorItrf(const Vector3D &rteme, const JulianDate &t)
Retrieve the satellite's position vector in ITRF coordinates.
virtual GeoCoordinate DoGetGeoPosition() const
bool IsInitialized(void) const
Check if the satellite has already been initialized.
void SetTleInfo(const std::string &tle)
Set satellite's TLE information required for its initialization.
JulianDate m_start
Simulation absolute start time.
SatArqSequenceNumber is handling the sequence numbers for the ARQ process.
Vector3D CrossProduct(const Vector3D &v1, const Vector3D &v2)
Cross product of two Vector3D objects.
void twoline2rv(char longstr1[130], char longstr2[130], char typerun, char typeinput, char opsmode, gravconsttype whichconst, double &startmfe, double &stopmfe, double &deltamin, elsetrec &satrec)
bool sgp4(gravconsttype whichconst, elsetrec &satrec, double tsince, double r[3], double v[3])
gravconsttype
@ wgs72
double jdsatepoch
Matrix data structure to make coordinate conversion code clearer and less verbose.
Vector3D operator*(const Vector3D &v) const