1   /* Copyright 2002-2018 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      * @exception OrekitException if some specific error occurs
170      */
171     protected TLEPropagator(final TLE initialTLE, final AttitudeProvider attitudeProvider,
172                             final double mass)
173         throws OrekitException {
174         super(attitudeProvider);
175         setStartDate(initialTLE.getDate());
176         this.tle  = initialTLE;
177         this.teme = FramesFactory.getTEME();
178         this.mass = mass;
179         initializeCommons();
180         sxpInitialize();
181         // set the initial state
182         final Orbit orbit = propagateOrbit(initialTLE.getDate());
183         final Attitude attitude = attitudeProvider.getAttitude(orbit, orbit.getDate(), orbit.getFrame());
184         super.resetInitialState(new SpacecraftState(orbit, attitude, mass));
185     }
186 
187     /** Selects the extrapolator to use with the selected TLE.
188      * @param tle the TLE to propagate.
189      * @return the correct propagator.
190      * @exception OrekitException if the underlying model cannot be initialized
191      */
192     public static TLEPropagator selectExtrapolator(final TLE tle) throws OrekitException {
193         return selectExtrapolator(tle, DEFAULT_LAW, DEFAULT_MASS);
194     }
195 
196     /** Selects the extrapolator to use with the selected TLE.
197      * @param tle the TLE to propagate.
198      * @param attitudeProvider provider for attitude computation
199      * @param mass spacecraft mass (kg)
200      * @return the correct propagator.
201      * @exception OrekitException if the underlying model cannot be initialized
202      */
203     public static TLEPropagator selectExtrapolator(final TLE tle, final AttitudeProvider attitudeProvider,
204                                                    final double mass) throws OrekitException {
205 
206         final double a1 = FastMath.pow( TLEConstants.XKE / (tle.getMeanMotion() * 60.0), TLEConstants.TWO_THIRD);
207         final double cosi0 = FastMath.cos(tle.getI());
208         final double temp = TLEConstants.CK2 * 1.5 * (3 * cosi0 * cosi0 - 1.0) *
209                             FastMath.pow(1.0 - tle.getE() * tle.getE(), -1.5);
210         final double delta1 = temp / (a1 * a1);
211         final double a0 = a1 * (1.0 - delta1 * (TLEConstants.ONE_THIRD + delta1 * (delta1 * 134.0 / 81.0 + 1.0)));
212         final double delta0 = temp / (a0 * a0);
213 
214         // recover original mean motion :
215         final double xn0dp = tle.getMeanMotion() * 60.0 / (delta0 + 1.0);
216 
217         // Period >= 225 minutes is deep space
218         if (MathUtils.TWO_PI / (xn0dp * TLEConstants.MINUTES_PER_DAY) >= (1.0 / 6.4)) {
219             return new DeepSDP4(tle, attitudeProvider, mass);
220         } else {
221             return new SGP4(tle, attitudeProvider, mass);
222         }
223     }
224 
225     /** Get the Earth gravity coefficient used for TLE propagation.
226      * @return the Earth gravity coefficient.
227      */
228     public static double getMU() {
229         return TLEConstants.MU;
230     }
231 
232     /** Get the extrapolated position and velocity from an initial TLE.
233      * @param date the final date
234      * @return the final PVCoordinates
235      * @exception OrekitException if propagation cannot be performed at given date
236      */
237     public PVCoordinates getPVCoordinates(final AbsoluteDate date)
238         throws OrekitException {
239 
240         sxpPropagate(date.durationFrom(tle.getDate()) / 60.0);
241 
242         // Compute PV with previous calculated parameters
243         return computePVCoordinates();
244     }
245 
246     /** Computation of the first commons parameters.
247      */
248     private void initializeCommons() {
249 
250         final double a1 = FastMath.pow(TLEConstants.XKE / (tle.getMeanMotion() * 60.0), TLEConstants.TWO_THIRD);
251         cosi0 = FastMath.cos(tle.getI());
252         theta2 = cosi0 * cosi0;
253         final double x3thm1 = 3.0 * theta2 - 1.0;
254         e0sq = tle.getE() * tle.getE();
255         beta02 = 1.0 - e0sq;
256         beta0 = FastMath.sqrt(beta02);
257         final double tval = TLEConstants.CK2 * 1.5 * x3thm1 / (beta0 * beta02);
258         final double delta1 = tval / (a1 * a1);
259         final double a0 = a1 * (1.0 - delta1 * (TLEConstants.ONE_THIRD + delta1 * (1.0 + 134.0 / 81.0 * delta1)));
260         final double delta0 = tval / (a0 * a0);
261 
262         // recover original mean motion and semi-major axis :
263         xn0dp = tle.getMeanMotion() * 60.0 / (delta0 + 1.0);
264         a0dp = a0 / (1.0 - delta0);
265 
266         // Values of s and qms2t :
267         s4 = TLEConstants.S;  // unmodified value for s
268         double q0ms24 = TLEConstants.QOMS2T; // unmodified value for q0ms2T
269 
270         perige = (a0dp * (1 - tle.getE()) - TLEConstants.NORMALIZED_EQUATORIAL_RADIUS) * TLEConstants.EARTH_RADIUS; // perige
271 
272         //  For perigee below 156 km, the values of s and qoms2t are changed :
273         if (perige < 156.0) {
274             if (perige <= 98.0) {
275                 s4 = 20.0;
276             } else {
277                 s4 = perige - 78.0;
278             }
279             final double temp_val = (120.0 - s4) * TLEConstants.NORMALIZED_EQUATORIAL_RADIUS / TLEConstants.EARTH_RADIUS;
280             final double temp_val_squared = temp_val * temp_val;
281             q0ms24 = temp_val_squared * temp_val_squared;
282             s4 = s4 / TLEConstants.EARTH_RADIUS + TLEConstants.NORMALIZED_EQUATORIAL_RADIUS; // new value for q0ms2T and s
283         }
284 
285         final double pinv = 1.0 / (a0dp * beta02);
286         final double pinvsq = pinv * pinv;
287         tsi = 1.0 / (a0dp - s4);
288         eta = a0dp * tle.getE() * tsi;
289         etasq = eta * eta;
290         eeta = tle.getE() * eta;
291 
292         final double psisq = FastMath.abs(1.0 - etasq); // abs because pow 3.5 needs positive value
293         final double tsi_squared = tsi * tsi;
294         coef = q0ms24 * tsi_squared * tsi_squared;
295         coef1 = coef / FastMath.pow(psisq, 3.5);
296 
297         // C2 and C1 coefficients computation :
298         c2 = coef1 * xn0dp * (a0dp * (1.0 + 1.5 * etasq + eeta * (4.0 + etasq)) +
299              0.75 * TLEConstants.CK2 * tsi / psisq * x3thm1 * (8.0 + 3.0 * etasq * (8.0 + etasq)));
300         c1 = tle.getBStar() * c2;
301         sini0 = FastMath.sin(tle.getI());
302 
303         final double x1mth2 = 1.0 - theta2;
304 
305         // C4 coefficient computation :
306         c4 = 2.0 * xn0dp * coef1 * a0dp * beta02 * (eta * (2.0 + 0.5 * etasq) +
307              tle.getE() * (0.5 + 2.0 * etasq) -
308              2 * TLEConstants.CK2 * tsi / (a0dp * psisq) *
309              (-3.0 * x3thm1 * (1.0 - 2.0 * eeta + etasq * (1.5 - 0.5 * eeta)) +
310               0.75 * x1mth2 * (2.0 * etasq - eeta * (1.0 + etasq)) * FastMath.cos(2.0 * tle.getPerigeeArgument())));
311 
312         final double theta4 = theta2 * theta2;
313         final double temp1 = 3 * TLEConstants.CK2 * pinvsq * xn0dp;
314         final double temp2 = temp1 * TLEConstants.CK2 * pinvsq;
315         final double temp3 = 1.25 * TLEConstants.CK4 * pinvsq * pinvsq * xn0dp;
316 
317         // atmospheric and gravitation coefs :(Mdf and OMEGAdf)
318         xmdot = xn0dp +
319                 0.5 * temp1 * beta0 * x3thm1 +
320                 0.0625 * temp2 * beta0 * (13.0 - 78.0 * theta2 + 137.0 * theta4);
321 
322         final double x1m5th = 1.0 - 5.0 * theta2;
323 
324         omgdot = -0.5 * temp1 * x1m5th +
325                  0.0625 * temp2 * (7.0 - 114.0 * theta2 + 395.0 * theta4) +
326                  temp3 * (3.0 - 36.0 * theta2 + 49.0 * theta4);
327 
328         final double xhdot1 = -temp1 * cosi0;
329 
330         xnodot = xhdot1 + (0.5 * temp2 * (4.0 - 19.0 * theta2) + 2.0 * temp3 * (3.0 - 7.0 * theta2)) * cosi0;
331         xnodcf = 3.5 * beta02 * xhdot1 * c1;
332         t2cof = 1.5 * c1;
333 
334     }
335 
336     /** Retrieves the position and velocity.
337      * @return the computed PVCoordinates.
338      * @exception OrekitException if current orbit is out of supported range
339      * (too large eccentricity, too low perigee ...)
340      */
341     private PVCoordinates computePVCoordinates() throws OrekitException {
342 
343         // Long period periodics
344         final double axn = e * FastMath.cos(omega);
345         double temp = 1.0 / (a * (1.0 - e * e));
346         final double xlcof = 0.125 * TLEConstants.A3OVK2 * sini0 * (3.0 + 5.0 * cosi0) / (1.0 + cosi0);
347         final double aycof = 0.25 * TLEConstants.A3OVK2 * sini0;
348         final double xll = temp * xlcof * axn;
349         final double aynl = temp * aycof;
350         final double xlt = xl + xll;
351         final double ayn = e * FastMath.sin(omega) + aynl;
352         final double elsq = axn * axn + ayn * ayn;
353         final double capu = MathUtils.normalizeAngle(xlt - xnode, FastMath.PI);
354         double epw = capu;
355         double ecosE = 0;
356         double esinE = 0;
357         double sinEPW = 0;
358         double cosEPW = 0;
359 
360         // Dundee changes:  items dependent on cosio get recomputed:
361         final double cosi0Sq = cosi0 * cosi0;
362         final double x3thm1 = 3.0 * cosi0Sq - 1.0;
363         final double x1mth2 = 1.0 - cosi0Sq;
364         final double x7thm1 = 7.0 * cosi0Sq - 1.0;
365 
366         if (e > (1 - 1e-6)) {
367             throw new OrekitException(OrekitMessages.TOO_LARGE_ECCENTRICITY_FOR_PROPAGATION_MODEL, e);
368         }
369 
370         // Solve Kepler's' Equation.
371         final double newtonRaphsonEpsilon = 1e-12;
372         for (int j = 0; j < 10; j++) {
373 
374             boolean doSecondOrderNewtonRaphson = true;
375 
376             sinEPW = FastMath.sin( epw);
377             cosEPW = FastMath.cos( epw);
378             ecosE = axn * cosEPW + ayn * sinEPW;
379             esinE = axn * sinEPW - ayn * cosEPW;
380             final double f = capu - epw + esinE;
381             if (FastMath.abs(f) < newtonRaphsonEpsilon) {
382                 break;
383             }
384             final double fdot = 1.0 - ecosE;
385             double delta_epw = f / fdot;
386             if (j == 0) {
387                 final double maxNewtonRaphson = 1.25 * FastMath.abs(e);
388                 doSecondOrderNewtonRaphson = false;
389                 if (delta_epw > maxNewtonRaphson) {
390                     delta_epw = maxNewtonRaphson;
391                 } else if (delta_epw < -maxNewtonRaphson) {
392                     delta_epw = -maxNewtonRaphson;
393                 } else {
394                     doSecondOrderNewtonRaphson = true;
395                 }
396             }
397             if (doSecondOrderNewtonRaphson) {
398                 delta_epw = f / (fdot + 0.5 * esinE * delta_epw);
399             }
400             epw += delta_epw;
401         }
402 
403         // Short period preliminary quantities
404         temp = 1.0 - elsq;
405         final double pl = a * temp;
406         final double r = a * (1.0 - ecosE);
407         double temp2 = a / r;
408         final double betal = FastMath.sqrt(temp);
409         temp = esinE / (1.0 + betal);
410         final double cosu = temp2 * (cosEPW - axn + ayn * temp);
411         final double sinu = temp2 * (sinEPW - ayn - axn * temp);
412         final double u = FastMath.atan2(sinu, cosu);
413         final double sin2u = 2.0 * sinu * cosu;
414         final double cos2u = 2.0 * cosu * cosu - 1.0;
415         final double temp1 = TLEConstants.CK2 / pl;
416         temp2 = temp1 / pl;
417 
418         // Update for short periodics
419         final double rk = r * (1.0 - 1.5 * temp2 * betal * x3thm1) + 0.5 * temp1 * x1mth2 * cos2u;
420         final double uk = u - 0.25 * temp2 * x7thm1 * sin2u;
421         final double xnodek = xnode + 1.5 * temp2 * cosi0 * sin2u;
422         final double xinck = i + 1.5 * temp2 * cosi0 * sini0 * cos2u;
423 
424         // Orientation vectors
425         final double sinuk = FastMath.sin(uk);
426         final double cosuk = FastMath.cos(uk);
427         final double sinik = FastMath.sin(xinck);
428         final double cosik = FastMath.cos(xinck);
429         final double sinnok = FastMath.sin(xnodek);
430         final double cosnok = FastMath.cos(xnodek);
431         final double xmx = -sinnok * cosik;
432         final double xmy = cosnok * cosik;
433         final double ux = xmx * sinuk + cosnok * cosuk;
434         final double uy = xmy * sinuk + sinnok * cosuk;
435         final double uz = sinik * sinuk;
436 
437         // Position and velocity
438         final double cr = 1000 * rk * TLEConstants.EARTH_RADIUS;
439         final Vector3D pos = new Vector3D(cr * ux, cr * uy, cr * uz);
440 
441         final double rdot   = TLEConstants.XKE * FastMath.sqrt(a) * esinE / r;
442         final double rfdot  = TLEConstants.XKE * FastMath.sqrt(pl) / r;
443         final double xn     = TLEConstants.XKE / (a * FastMath.sqrt(a));
444         final double rdotk  = rdot - xn * temp1 * x1mth2 * sin2u;
445         final double rfdotk = rfdot + xn * temp1 * (x1mth2 * cos2u + 1.5 * x3thm1);
446         final double vx     = xmx * cosuk - cosnok * sinuk;
447         final double vy     = xmy * cosuk - sinnok * sinuk;
448         final double vz     = sinik * cosuk;
449 
450         final double cv = 1000.0 * TLEConstants.EARTH_RADIUS / 60.0;
451         final Vector3D vel = new Vector3D(cv * (rdotk * ux + rfdotk * vx),
452                                           cv * (rdotk * uy + rfdotk * vy),
453                                           cv * (rdotk * uz + rfdotk * vz));
454 
455         return new PVCoordinates(pos, vel);
456 
457     }
458 
459     /** Initialization proper to each propagator (SGP or SDP).
460      * @exception OrekitException if some specific error occurs
461      */
462     protected abstract void sxpInitialize() throws OrekitException;
463 
464     /** Propagation proper to each propagator (SGP or SDP).
465      * @param t the offset from initial epoch (min)
466      * @exception OrekitException if current state cannot be propagated
467      */
468     protected abstract void sxpPropagate(double t) throws OrekitException;
469 
470     /** {@inheritDoc} */
471     public void resetInitialState(final SpacecraftState state)
472         throws OrekitException {
473         throw new OrekitException(OrekitMessages.NON_RESETABLE_STATE);
474     }
475 
476     /** {@inheritDoc} */
477     protected void resetIntermediateState(final SpacecraftState state, final boolean forward)
478         throws OrekitException {
479         throw new OrekitException(OrekitMessages.NON_RESETABLE_STATE);
480     }
481 
482     /** {@inheritDoc} */
483     protected double getMass(final AbsoluteDate date) {
484         return mass;
485     }
486 
487     /** {@inheritDoc} */
488     protected Orbit propagateOrbit(final AbsoluteDate date) throws OrekitException {
489         return new CartesianOrbit(getPVCoordinates(date), teme, date, TLEConstants.MU);
490     }
491 
492     /** Get the underlying TLE.
493      * @return underlying TLE
494      */
495     public TLE getTLE() {
496         return tle;
497     }
498 
499     /** {@inheritDoc} */
500     public Frame getFrame() {
501         return teme;
502     }
503 
504 }