1   /* Copyright 2002-2019 CS Systèmes d'Information
2    * Licensed to CS Systèmes d'Information (CS) under one or more
3    * contributor license agreements.  See the NOTICE file distributed with
4    * this work for additional information regarding copyright ownership.
5    * CS licenses this file to You under the Apache License, Version 2.0
6    * (the "License"); you may not use this file except in compliance with
7    * the License.  You may obtain a copy of the License at
8    *
9    *   http://www.apache.org/licenses/LICENSE-2.0
10   *
11   * Unless required by applicable law or agreed to in writing, software
12   * distributed under the License is distributed on an "AS IS" BASIS,
13   * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14   * See the License for the specific language governing permissions and
15   * limitations under the License.
16   */
17  package org.orekit.propagation.analytical.tle;
18  
19  import org.hipparchus.geometry.euclidean.threed.Vector3D;
20  import org.hipparchus.util.FastMath;
21  import org.hipparchus.util.MathUtils;
22  import org.orekit.attitudes.Attitude;
23  import org.orekit.attitudes.AttitudeProvider;
24  import org.orekit.errors.OrekitException;
25  import org.orekit.errors.OrekitMessages;
26  import org.orekit.frames.Frame;
27  import org.orekit.frames.FramesFactory;
28  import org.orekit.orbits.CartesianOrbit;
29  import org.orekit.orbits.Orbit;
30  import org.orekit.propagation.SpacecraftState;
31  import org.orekit.propagation.analytical.AbstractAnalyticalPropagator;
32  import org.orekit.time.AbsoluteDate;
33  import org.orekit.utils.PVCoordinates;
34  
35  
36  /** This class provides elements to propagate TLE's.
37   * <p>
38   * The models used are SGP4 and SDP4, initially proposed by NORAD as the unique convenient
39   * propagator for TLE's. Inputs and outputs of this propagator are only suited for
40   * NORAD two lines elements sets, since it uses estimations and mean values appropriate
41   * for TLE's only.
42   * </p>
43   * <p>
44   * Deep- or near- space propagator is selected internally according to NORAD recommendations
45   * so that the user has not to worry about the used computation methods. One instance is created
46   * for each TLE (this instance can only be get using {@link #selectExtrapolator(TLE)} method,
47   * and can compute {@link PVCoordinates position and velocity coordinates} at any
48   * time. Maximum accuracy is guaranteed in a 24h range period before and after the provided
49   * TLE epoch (of course this accuracy is not really measurable nor predictable: according to
50   * <a href="http://www.celestrak.com/">CelesTrak</a>, the precision is close to one kilometer
51   * and error won't probably rise above 2 km).
52   * </p>
53   * <p>This implementation is largely inspired from the paper and source code <a
54   * href="http://www.celestrak.com/publications/AIAA/2006-6753/">Revisiting Spacetrack
55   * Report #3</a> and is fully compliant with its results and tests cases.</p>
56   * @author Felix R. Hoots, Ronald L. Roehrich, December 1980 (original fortran)
57   * @author David A. Vallado, Paul Crawford, Richard Hujsak, T.S. Kelso (C++ translation and improvements)
58   * @author Fabien Maussion (java translation)
59   * @see TLE
60   */
61  public abstract class TLEPropagator extends AbstractAnalyticalPropagator {
62  
63      // CHECKSTYLE: stop VisibilityModifier check
64  
65      /** Initial state. */
66      protected final TLE tle;
67  
68      /** final RAAN. */
69      protected double xnode;
70  
71      /** final semi major axis. */
72      protected double a;
73  
74      /** final eccentricity. */
75      protected double e;
76  
77      /** final inclination. */
78      protected double i;
79  
80      /** final perigee argument. */
81      protected double omega;
82  
83      /** L from SPTRCK #3. */
84      protected double xl;
85  
86      /** original recovered semi major axis. */
87      protected double a0dp;
88  
89      /** original recovered mean motion. */
90      protected double xn0dp;
91  
92      /** cosinus original inclination. */
93      protected double cosi0;
94  
95      /** cos io squared. */
96      protected double theta2;
97  
98      /** sinus original inclination. */
99      protected double sini0;
100 
101     /** common parameter for mean anomaly (M) computation. */
102     protected double xmdot;
103 
104     /** common parameter for perigee argument (omega) computation. */
105     protected double omgdot;
106 
107     /** common parameter for raan (OMEGA) computation. */
108     protected double xnodot;
109 
110     /** original eccentricity squared. */
111     protected double e0sq;
112     /** 1 - e2. */
113     protected double beta02;
114 
115     /** sqrt (1 - e2). */
116     protected double beta0;
117 
118     /** perigee, expressed in KM and ALTITUDE. */
119     protected double perige;
120 
121     /** eta squared. */
122     protected double etasq;
123 
124     /** original eccentricity * eta. */
125     protected double eeta;
126 
127     /** s* new value for the contant s. */
128     protected double s4;
129 
130     /** tsi from SPTRCK #3. */
131     protected double tsi;
132 
133     /** eta from SPTRCK #3. */
134     protected double eta;
135 
136     /** coef for SGP C3 computation. */
137     protected double coef;
138 
139     /** coef for SGP C5 computation. */
140     protected double coef1;
141 
142     /** C1 from SPTRCK #3. */
143     protected double c1;
144 
145     /** C2 from SPTRCK #3. */
146     protected double c2;
147 
148     /** C4 from SPTRCK #3. */
149     protected double c4;
150 
151     /** common parameter for raan (OMEGA) computation. */
152     protected double xnodcf;
153 
154     /** 3/2 * C1. */
155     protected double t2cof;
156 
157     // CHECKSTYLE: resume VisibilityModifier check
158 
159     /** TLE frame. */
160     private final Frame teme;
161 
162     /** Spacecraft mass (kg). */
163     private final double mass;
164 
165     /** Protected constructor for derived classes.
166      * @param initialTLE the unique TLE to propagate
167      * @param attitudeProvider provider for attitude computation
168      * @param mass spacecraft mass (kg)
169      */
170     protected TLEPropagator(final TLE initialTLE, final AttitudeProvider attitudeProvider,
171                             final double mass) {
172         super(attitudeProvider);
173         setStartDate(initialTLE.getDate());
174         this.tle  = initialTLE;
175         this.teme = FramesFactory.getTEME();
176         this.mass = mass;
177         initializeCommons();
178         sxpInitialize();
179         // set the initial state
180         final Orbit orbit = propagateOrbit(initialTLE.getDate());
181         final Attitude attitude = attitudeProvider.getAttitude(orbit, orbit.getDate(), orbit.getFrame());
182         super.resetInitialState(new SpacecraftState(orbit, attitude, mass));
183     }
184 
185     /** Selects the extrapolator to use with the selected TLE.
186      * @param tle the TLE to propagate.
187      * @return the correct propagator.
188      */
189     public static TLEPropagator selectExtrapolator(final TLE tle) {
190         return selectExtrapolator(tle, DEFAULT_LAW, DEFAULT_MASS);
191     }
192 
193     /** Selects the extrapolator to use with the selected TLE.
194      * @param tle the TLE to propagate.
195      * @param attitudeProvider provider for attitude computation
196      * @param mass spacecraft mass (kg)
197      * @return the correct propagator.
198      */
199     public static TLEPropagator selectExtrapolator(final TLE tle, final AttitudeProvider attitudeProvider,
200                                                    final double mass) {
201 
202         final double a1 = FastMath.pow( TLEConstants.XKE / (tle.getMeanMotion() * 60.0), TLEConstants.TWO_THIRD);
203         final double cosi0 = FastMath.cos(tle.getI());
204         final double temp = TLEConstants.CK2 * 1.5 * (3 * cosi0 * cosi0 - 1.0) *
205                             FastMath.pow(1.0 - tle.getE() * tle.getE(), -1.5);
206         final double delta1 = temp / (a1 * a1);
207         final double a0 = a1 * (1.0 - delta1 * (TLEConstants.ONE_THIRD + delta1 * (delta1 * 134.0 / 81.0 + 1.0)));
208         final double delta0 = temp / (a0 * a0);
209 
210         // recover original mean motion :
211         final double xn0dp = tle.getMeanMotion() * 60.0 / (delta0 + 1.0);
212 
213         // Period >= 225 minutes is deep space
214         if (MathUtils.TWO_PI / (xn0dp * TLEConstants.MINUTES_PER_DAY) >= (1.0 / 6.4)) {
215             return new DeepSDP4(tle, attitudeProvider, mass);
216         } else {
217             return new SGP4(tle, attitudeProvider, mass);
218         }
219     }
220 
221     /** Get the Earth gravity coefficient used for TLE propagation.
222      * @return the Earth gravity coefficient.
223      */
224     public static double getMU() {
225         return TLEConstants.MU;
226     }
227 
228     /** Get the extrapolated position and velocity from an initial TLE.
229      * @param date the final date
230      * @return the final PVCoordinates
231      */
232     public PVCoordinates getPVCoordinates(final AbsoluteDate date) {
233 
234         sxpPropagate(date.durationFrom(tle.getDate()) / 60.0);
235 
236         // Compute PV with previous calculated parameters
237         return computePVCoordinates();
238     }
239 
240     /** Computation of the first commons parameters.
241      */
242     private void initializeCommons() {
243 
244         final double a1 = FastMath.pow(TLEConstants.XKE / (tle.getMeanMotion() * 60.0), TLEConstants.TWO_THIRD);
245         cosi0 = FastMath.cos(tle.getI());
246         theta2 = cosi0 * cosi0;
247         final double x3thm1 = 3.0 * theta2 - 1.0;
248         e0sq = tle.getE() * tle.getE();
249         beta02 = 1.0 - e0sq;
250         beta0 = FastMath.sqrt(beta02);
251         final double tval = TLEConstants.CK2 * 1.5 * x3thm1 / (beta0 * beta02);
252         final double delta1 = tval / (a1 * a1);
253         final double a0 = a1 * (1.0 - delta1 * (TLEConstants.ONE_THIRD + delta1 * (1.0 + 134.0 / 81.0 * delta1)));
254         final double delta0 = tval / (a0 * a0);
255 
256         // recover original mean motion and semi-major axis :
257         xn0dp = tle.getMeanMotion() * 60.0 / (delta0 + 1.0);
258         a0dp = a0 / (1.0 - delta0);
259 
260         // Values of s and qms2t :
261         s4 = TLEConstants.S;  // unmodified value for s
262         double q0ms24 = TLEConstants.QOMS2T; // unmodified value for q0ms2T
263 
264         perige = (a0dp * (1 - tle.getE()) - TLEConstants.NORMALIZED_EQUATORIAL_RADIUS) * TLEConstants.EARTH_RADIUS; // perige
265 
266         //  For perigee below 156 km, the values of s and qoms2t are changed :
267         if (perige < 156.0) {
268             if (perige <= 98.0) {
269                 s4 = 20.0;
270             } else {
271                 s4 = perige - 78.0;
272             }
273             final double temp_val = (120.0 - s4) * TLEConstants.NORMALIZED_EQUATORIAL_RADIUS / TLEConstants.EARTH_RADIUS;
274             final double temp_val_squared = temp_val * temp_val;
275             q0ms24 = temp_val_squared * temp_val_squared;
276             s4 = s4 / TLEConstants.EARTH_RADIUS + TLEConstants.NORMALIZED_EQUATORIAL_RADIUS; // new value for q0ms2T and s
277         }
278 
279         final double pinv = 1.0 / (a0dp * beta02);
280         final double pinvsq = pinv * pinv;
281         tsi = 1.0 / (a0dp - s4);
282         eta = a0dp * tle.getE() * tsi;
283         etasq = eta * eta;
284         eeta = tle.getE() * eta;
285 
286         final double psisq = FastMath.abs(1.0 - etasq); // abs because pow 3.5 needs positive value
287         final double tsi_squared = tsi * tsi;
288         coef = q0ms24 * tsi_squared * tsi_squared;
289         coef1 = coef / FastMath.pow(psisq, 3.5);
290 
291         // C2 and C1 coefficients computation :
292         c2 = coef1 * xn0dp * (a0dp * (1.0 + 1.5 * etasq + eeta * (4.0 + etasq)) +
293              0.75 * TLEConstants.CK2 * tsi / psisq * x3thm1 * (8.0 + 3.0 * etasq * (8.0 + etasq)));
294         c1 = tle.getBStar() * c2;
295         sini0 = FastMath.sin(tle.getI());
296 
297         final double x1mth2 = 1.0 - theta2;
298 
299         // C4 coefficient computation :
300         c4 = 2.0 * xn0dp * coef1 * a0dp * beta02 * (eta * (2.0 + 0.5 * etasq) +
301              tle.getE() * (0.5 + 2.0 * etasq) -
302              2 * TLEConstants.CK2 * tsi / (a0dp * psisq) *
303              (-3.0 * x3thm1 * (1.0 - 2.0 * eeta + etasq * (1.5 - 0.5 * eeta)) +
304               0.75 * x1mth2 * (2.0 * etasq - eeta * (1.0 + etasq)) * FastMath.cos(2.0 * tle.getPerigeeArgument())));
305 
306         final double theta4 = theta2 * theta2;
307         final double temp1 = 3 * TLEConstants.CK2 * pinvsq * xn0dp;
308         final double temp2 = temp1 * TLEConstants.CK2 * pinvsq;
309         final double temp3 = 1.25 * TLEConstants.CK4 * pinvsq * pinvsq * xn0dp;
310 
311         // atmospheric and gravitation coefs :(Mdf and OMEGAdf)
312         xmdot = xn0dp +
313                 0.5 * temp1 * beta0 * x3thm1 +
314                 0.0625 * temp2 * beta0 * (13.0 - 78.0 * theta2 + 137.0 * theta4);
315 
316         final double x1m5th = 1.0 - 5.0 * theta2;
317 
318         omgdot = -0.5 * temp1 * x1m5th +
319                  0.0625 * temp2 * (7.0 - 114.0 * theta2 + 395.0 * theta4) +
320                  temp3 * (3.0 - 36.0 * theta2 + 49.0 * theta4);
321 
322         final double xhdot1 = -temp1 * cosi0;
323 
324         xnodot = xhdot1 + (0.5 * temp2 * (4.0 - 19.0 * theta2) + 2.0 * temp3 * (3.0 - 7.0 * theta2)) * cosi0;
325         xnodcf = 3.5 * beta02 * xhdot1 * c1;
326         t2cof = 1.5 * c1;
327 
328     }
329 
330     /** Retrieves the position and velocity.
331      * @return the computed PVCoordinates.
332      */
333     private PVCoordinates computePVCoordinates() {
334 
335         // Long period periodics
336         final double axn = e * FastMath.cos(omega);
337         double temp = 1.0 / (a * (1.0 - e * e));
338         final double xlcof = 0.125 * TLEConstants.A3OVK2 * sini0 * (3.0 + 5.0 * cosi0) / (1.0 + cosi0);
339         final double aycof = 0.25 * TLEConstants.A3OVK2 * sini0;
340         final double xll = temp * xlcof * axn;
341         final double aynl = temp * aycof;
342         final double xlt = xl + xll;
343         final double ayn = e * FastMath.sin(omega) + aynl;
344         final double elsq = axn * axn + ayn * ayn;
345         final double capu = MathUtils.normalizeAngle(xlt - xnode, FastMath.PI);
346         double epw = capu;
347         double ecosE = 0;
348         double esinE = 0;
349         double sinEPW = 0;
350         double cosEPW = 0;
351 
352         // Dundee changes:  items dependent on cosio get recomputed:
353         final double cosi0Sq = cosi0 * cosi0;
354         final double x3thm1 = 3.0 * cosi0Sq - 1.0;
355         final double x1mth2 = 1.0 - cosi0Sq;
356         final double x7thm1 = 7.0 * cosi0Sq - 1.0;
357 
358         if (e > (1 - 1e-6)) {
359             throw new OrekitException(OrekitMessages.TOO_LARGE_ECCENTRICITY_FOR_PROPAGATION_MODEL, e);
360         }
361 
362         // Solve Kepler's' Equation.
363         final double newtonRaphsonEpsilon = 1e-12;
364         for (int j = 0; j < 10; j++) {
365 
366             boolean doSecondOrderNewtonRaphson = true;
367 
368             sinEPW = FastMath.sin( epw);
369             cosEPW = FastMath.cos( epw);
370             ecosE = axn * cosEPW + ayn * sinEPW;
371             esinE = axn * sinEPW - ayn * cosEPW;
372             final double f = capu - epw + esinE;
373             if (FastMath.abs(f) < newtonRaphsonEpsilon) {
374                 break;
375             }
376             final double fdot = 1.0 - ecosE;
377             double delta_epw = f / fdot;
378             if (j == 0) {
379                 final double maxNewtonRaphson = 1.25 * FastMath.abs(e);
380                 doSecondOrderNewtonRaphson = false;
381                 if (delta_epw > maxNewtonRaphson) {
382                     delta_epw = maxNewtonRaphson;
383                 } else if (delta_epw < -maxNewtonRaphson) {
384                     delta_epw = -maxNewtonRaphson;
385                 } else {
386                     doSecondOrderNewtonRaphson = true;
387                 }
388             }
389             if (doSecondOrderNewtonRaphson) {
390                 delta_epw = f / (fdot + 0.5 * esinE * delta_epw);
391             }
392             epw += delta_epw;
393         }
394 
395         // Short period preliminary quantities
396         temp = 1.0 - elsq;
397         final double pl = a * temp;
398         final double r = a * (1.0 - ecosE);
399         double temp2 = a / r;
400         final double betal = FastMath.sqrt(temp);
401         temp = esinE / (1.0 + betal);
402         final double cosu = temp2 * (cosEPW - axn + ayn * temp);
403         final double sinu = temp2 * (sinEPW - ayn - axn * temp);
404         final double u = FastMath.atan2(sinu, cosu);
405         final double sin2u = 2.0 * sinu * cosu;
406         final double cos2u = 2.0 * cosu * cosu - 1.0;
407         final double temp1 = TLEConstants.CK2 / pl;
408         temp2 = temp1 / pl;
409 
410         // Update for short periodics
411         final double rk = r * (1.0 - 1.5 * temp2 * betal * x3thm1) + 0.5 * temp1 * x1mth2 * cos2u;
412         final double uk = u - 0.25 * temp2 * x7thm1 * sin2u;
413         final double xnodek = xnode + 1.5 * temp2 * cosi0 * sin2u;
414         final double xinck = i + 1.5 * temp2 * cosi0 * sini0 * cos2u;
415 
416         // Orientation vectors
417         final double sinuk = FastMath.sin(uk);
418         final double cosuk = FastMath.cos(uk);
419         final double sinik = FastMath.sin(xinck);
420         final double cosik = FastMath.cos(xinck);
421         final double sinnok = FastMath.sin(xnodek);
422         final double cosnok = FastMath.cos(xnodek);
423         final double xmx = -sinnok * cosik;
424         final double xmy = cosnok * cosik;
425         final double ux = xmx * sinuk + cosnok * cosuk;
426         final double uy = xmy * sinuk + sinnok * cosuk;
427         final double uz = sinik * sinuk;
428 
429         // Position and velocity
430         final double cr = 1000 * rk * TLEConstants.EARTH_RADIUS;
431         final Vector3D pos = new Vector3D(cr * ux, cr * uy, cr * uz);
432 
433         final double rdot   = TLEConstants.XKE * FastMath.sqrt(a) * esinE / r;
434         final double rfdot  = TLEConstants.XKE * FastMath.sqrt(pl) / r;
435         final double xn     = TLEConstants.XKE / (a * FastMath.sqrt(a));
436         final double rdotk  = rdot - xn * temp1 * x1mth2 * sin2u;
437         final double rfdotk = rfdot + xn * temp1 * (x1mth2 * cos2u + 1.5 * x3thm1);
438         final double vx     = xmx * cosuk - cosnok * sinuk;
439         final double vy     = xmy * cosuk - sinnok * sinuk;
440         final double vz     = sinik * cosuk;
441 
442         final double cv = 1000.0 * TLEConstants.EARTH_RADIUS / 60.0;
443         final Vector3D vel = new Vector3D(cv * (rdotk * ux + rfdotk * vx),
444                                           cv * (rdotk * uy + rfdotk * vy),
445                                           cv * (rdotk * uz + rfdotk * vz));
446 
447         return new PVCoordinates(pos, vel);
448 
449     }
450 
451     /** Initialization proper to each propagator (SGP or SDP).
452      */
453     protected abstract void sxpInitialize();
454 
455     /** Propagation proper to each propagator (SGP or SDP).
456      * @param t the offset from initial epoch (min)
457      */
458     protected abstract void sxpPropagate(double t);
459 
460     /** {@inheritDoc} */
461     public void resetInitialState(final SpacecraftState state) {
462         throw new OrekitException(OrekitMessages.NON_RESETABLE_STATE);
463     }
464 
465     /** {@inheritDoc} */
466     protected void resetIntermediateState(final SpacecraftState state, final boolean forward) {
467         throw new OrekitException(OrekitMessages.NON_RESETABLE_STATE);
468     }
469 
470     /** {@inheritDoc} */
471     protected double getMass(final AbsoluteDate date) {
472         return mass;
473     }
474 
475     /** {@inheritDoc} */
476     protected Orbit propagateOrbit(final AbsoluteDate date) {
477         return new CartesianOrbit(getPVCoordinates(date), teme, date, TLEConstants.MU);
478     }
479 
480     /** Get the underlying TLE.
481      * @return underlying TLE
482      */
483     public TLE getTLE() {
484         return tle;
485     }
486 
487     /** {@inheritDoc} */
488     public Frame getFrame() {
489         return teme;
490     }
491 
492 }