1   /* Copyright 2002-2024 CS GROUP
2    * Licensed to CS GROUP (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  
20  import java.util.HashMap;
21  import java.util.List;
22  import java.util.Map;
23  
24  import org.hipparchus.CalculusFieldElement;
25  import org.hipparchus.geometry.euclidean.threed.FieldVector3D;
26  import org.hipparchus.util.FastMath;
27  import org.hipparchus.util.MathUtils;
28  import org.orekit.annotation.DefaultDataContext;
29  import org.orekit.attitudes.AttitudeProvider;
30  import org.orekit.attitudes.FieldAttitude;
31  import org.orekit.attitudes.FrameAlignedProvider;
32  import org.orekit.data.DataContext;
33  import org.orekit.errors.OrekitException;
34  import org.orekit.errors.OrekitMessages;
35  import org.orekit.frames.Frame;
36  import org.orekit.orbits.FieldCartesianOrbit;
37  import org.orekit.orbits.FieldOrbit;
38  import org.orekit.propagation.FieldSpacecraftState;
39  import org.orekit.propagation.analytical.FieldAbstractAnalyticalPropagator;
40  import org.orekit.propagation.analytical.tle.generation.TleGenerationAlgorithm;
41  import org.orekit.time.FieldAbsoluteDate;
42  import org.orekit.time.TimeScale;
43  import org.orekit.utils.FieldPVCoordinates;
44  import org.orekit.utils.PVCoordinates;
45  import org.orekit.utils.ParameterDriver;
46  import org.orekit.utils.TimeSpanMap;
47  
48  
49  /** This class provides elements to propagate TLE's.
50   * <p>
51   * The models used are SGP4 and SDP4, initially proposed by NORAD as the unique convenient
52   * propagator for TLE's. Inputs and outputs of this propagator are only suited for
53   * NORAD two lines elements sets, since it uses estimations and mean values appropriate
54   * for TLE's only.
55   * </p>
56   * <p>
57   * Deep- or near- space propagator is selected internally according to NORAD recommendations
58   * so that the user has not to worry about the used computation methods. One instance is created
59   * for each TLE (this instance can only be get using {@link #selectExtrapolator(FieldTLE, CalculusFieldElement[])} method,
60   * and can compute {@link PVCoordinates position and velocity coordinates} at any
61   * time. Maximum accuracy is guaranteed in a 24h range period before and after the provided
62   * TLE epoch (of course this accuracy is not really measurable nor predictable: according to
63   * <a href="https://www.celestrak.com/">CelesTrak</a>, the precision is close to one kilometer
64   * and error won't probably rise above 2 km).
65   * </p>
66   * <p>This implementation is largely inspired from the paper and source code <a
67   * href="https://www.celestrak.com/publications/AIAA/2006-6753/">Revisiting Spacetrack
68   * Report #3</a> and is fully compliant with its results and tests cases.</p>
69   * @author Felix R. Hoots, Ronald L. Roehrich, December 1980 (original fortran)
70   * @author David A. Vallado, Paul Crawford, Richard Hujsak, T.S. Kelso (C++ translation and improvements)
71   * @author Fabien Maussion (java translation)
72   * @author Thomas Paulet (field translation)
73   * @since 11.0
74   * @see FieldTLE
75   * @param <T> type of the field elements
76   */
77  public abstract class FieldTLEPropagator<T extends CalculusFieldElement<T>> extends FieldAbstractAnalyticalPropagator<T> {
78  
79      // CHECKSTYLE: stop VisibilityModifier check
80  
81      /** Initial state. */
82      protected FieldTLE<T> tle;
83  
84      /** UTC time scale. */
85      protected final TimeScale utc;
86  
87      /** final RAAN. */
88      protected T xnode;
89  
90      /** final semi major axis. */
91      protected T a;
92  
93      /** final eccentricity. */
94      protected T e;
95  
96      /** final inclination. */
97      protected T i;
98  
99      /** final perigee argument. */
100     protected T omega;
101 
102     /** L from SPTRCK #3. */
103     protected T xl;
104 
105     /** original recovered semi major axis. */
106     protected T a0dp;
107 
108     /** original recovered mean motion. */
109     protected T xn0dp;
110 
111     /** cosinus original inclination. */
112     protected T cosi0;
113 
114     /** cos io squared. */
115     protected T theta2;
116 
117     /** sinus original inclination. */
118     protected T sini0;
119 
120     /** common parameter for mean anomaly (M) computation. */
121     protected T xmdot;
122 
123     /** common parameter for perigee argument (omega) computation. */
124     protected T omgdot;
125 
126     /** common parameter for raan (OMEGA) computation. */
127     protected T xnodot;
128 
129     /** original eccentricity squared. */
130     protected T e0sq;
131     /** 1 - e2. */
132     protected T beta02;
133 
134     /** sqrt (1 - e2). */
135     protected T beta0;
136 
137     /** perigee, expressed in KM and ALTITUDE. */
138     protected T perige;
139 
140     /** eta squared. */
141     protected T etasq;
142 
143     /** original eccentricity * eta. */
144     protected T eeta;
145 
146     /** s* new value for the contant s. */
147     protected T s4;
148 
149     /** tsi from SPTRCK #3. */
150     protected T tsi;
151 
152     /** eta from SPTRCK #3. */
153     protected T eta;
154 
155     /** coef for SGP C3 computation. */
156     protected T coef;
157 
158     /** coef for SGP C5 computation. */
159     protected T coef1;
160 
161     /** C1 from SPTRCK #3. */
162     protected T c1;
163 
164     /** C2 from SPTRCK #3. */
165     protected T c2;
166 
167     /** C4 from SPTRCK #3. */
168     protected T c4;
169 
170     /** common parameter for raan (OMEGA) computation. */
171     protected T xnodcf;
172 
173     /** 3/2 * C1. */
174     protected T t2cof;
175 
176     // CHECKSTYLE: resume VisibilityModifier check
177 
178     /** TLE frame. */
179     private final Frame teme;
180 
181     /** Spacecraft masses (kg) mapped to TLEs. */
182     private Map<FieldTLE<T>, T> masses;
183 
184     /** All TLEs. */
185     private TimeSpanMap<FieldTLE<T>> tles;
186 
187     /** Protected constructor for derived classes.
188      *
189      * <p>This constructor uses the {@link DataContext#getDefault() default data context}.
190      *
191      * @param initialTLE the unique TLE to propagate
192      * @param attitudeProvider provider for attitude computation
193      * @param mass spacecraft mass (kg)
194      * @param parameters SGP4 and SDP4 model parameters
195      * @see #FieldTLEPropagator(FieldTLE, AttitudeProvider, CalculusFieldElement, Frame, CalculusFieldElement[])
196      */
197     @DefaultDataContext
198     protected FieldTLEPropagator(final FieldTLE<T> initialTLE, final AttitudeProvider attitudeProvider, final T mass,
199                                  final T[] parameters) {
200         this(initialTLE, attitudeProvider, mass, DataContext.getDefault().getFrames().getTEME(), parameters);
201     }
202 
203     /** Protected constructor for derived classes.
204      * @param initialTLE the unique TLE to propagate
205      * @param attitudeProvider provider for attitude computation
206      * @param mass spacecraft mass (kg)
207      * @param teme the TEME frame to use for propagation.
208      * @param parameters SGP4 and SDP4 model parameters
209      */
210     protected FieldTLEPropagator(final FieldTLE<T> initialTLE, final AttitudeProvider attitudeProvider, final T mass,
211                                  final Frame teme, final T[] parameters) {
212         super(initialTLE.getE().getField(), attitudeProvider);
213         setStartDate(initialTLE.getDate());
214         this.utc       = initialTLE.getUtc();
215         initializeTle(initialTLE);
216         this.teme      = teme;
217         this.tles      = new TimeSpanMap<>(tle);
218         this.masses    = new HashMap<>();
219         this.masses.put(tle, mass);
220 
221         initializeCommons(parameters);
222         sxpInitialize(parameters);
223         // set the initial state
224         final FieldOrbit<T> orbit = propagateOrbit(initialTLE.getDate(), parameters);
225         final FieldAttitude<T> attitude = attitudeProvider.getAttitude(orbit, orbit.getDate(), orbit.getFrame());
226         super.resetInitialState(new FieldSpacecraftState<>(orbit, attitude, mass));
227     }
228 
229     /** Selects the extrapolator to use with the selected TLE.
230      *
231      * <p>This method uses the {@link DataContext#getDefault() default data context}.
232      *
233      * @param tle the TLE to propagate.
234      * @param parameters SGP4 and SDP4 model parameters
235      * @return the correct propagator.
236      * @param <T> elements type
237      * @see #selectExtrapolator(FieldTLE, Frame, CalculusFieldElement[])
238      */
239     @DefaultDataContext
240     public static <T extends CalculusFieldElement<T>> FieldTLEPropagator<T> selectExtrapolator(final FieldTLE<T> tle, final T[] parameters) {
241         return selectExtrapolator(tle, DataContext.getDefault().getFrames().getTEME(), parameters);
242     }
243 
244     /** Selects the extrapolator to use with the selected TLE.
245      *
246      *<p>This method uses the {@link DataContext#getDefault() default data context}.
247      *
248      * @param tle the TLE to propagate.
249      * @param teme TEME frame.
250      * @param parameters SGP4 and SDP4 model parameters
251      * @return the correct propagator.
252      * @param <T> elements type
253      */
254     public static <T extends CalculusFieldElement<T>> FieldTLEPropagator<T> selectExtrapolator(final FieldTLE<T> tle, final Frame teme, final T[] parameters) {
255         return selectExtrapolator(
256                 tle,
257                 FrameAlignedProvider.of(teme),
258                 tle.getE().getField().getZero().newInstance(DEFAULT_MASS),
259                 teme,
260                 parameters);
261     }
262 
263     /** Selects the extrapolator to use with the selected TLE.
264      *
265      * <p>This method uses the {@link DataContext#getDefault() default data context}.
266      *
267      * @param tle the TLE to propagate.
268      * @param attitudeProvider provider for attitude computation
269      * @param mass spacecraft mass (kg)
270      * @param parameters SGP4 and SDP4 model parameters
271      * @return the correct propagator.
272      * @param <T> elements type
273      * @see #selectExtrapolator(FieldTLE, AttitudeProvider, CalculusFieldElement, Frame, CalculusFieldElement[])
274      */
275     @DefaultDataContext
276     public static <T extends CalculusFieldElement<T>> FieldTLEPropagator<T> selectExtrapolator(final FieldTLE<T> tle,
277                                                    final AttitudeProvider attitudeProvider,
278                                                    final T mass,
279                                                    final T[] parameters) {
280         return selectExtrapolator(tle, attitudeProvider, mass,
281                 DataContext.getDefault().getFrames().getTEME(), parameters);
282     }
283 
284     /** Selects the extrapolator to use with the selected TLE.
285      *
286      * @param tle the TLE to propagate.
287      * @param attitudeProvider provider for attitude computation
288      * @param mass spacecraft mass (kg)
289      * @param teme the TEME frame to use for propagation.
290      * @param parameters SGP4 and SDP4 model parameters
291      * @return the correct propagator.
292      * @param <T> elements type
293      */
294     public static <T extends CalculusFieldElement<T>> FieldTLEPropagator<T> selectExtrapolator(final FieldTLE<T> tle,
295                                                    final AttitudeProvider attitudeProvider,
296                                                    final T mass,
297                                                    final Frame teme,
298                                                    final T[] parameters) {
299 
300         final T a1 = tle.getMeanMotion().multiply(60.0).reciprocal().multiply(TLEConstants.XKE).pow(TLEConstants.TWO_THIRD);
301         final T cosi0 = FastMath.cos(tle.getI());
302         final T temp1 = cosi0.multiply(cosi0.multiply(3.0)).subtract(1.0).multiply(1.5 * TLEConstants.CK2);
303         final T temp2 = tle.getE().multiply(tle.getE()).negate().add(1.0).pow(-1.5);
304         final T temp = temp1.multiply(temp2);
305         final T delta1 = temp.divide(a1.multiply(a1));
306         final T a0 = a1.multiply(delta1.multiply(delta1.multiply(
307                         delta1.multiply(134.0 / 81.0).add(1.0)).add(TLEConstants.ONE_THIRD)).negate().add(1.0));
308         final T delta0 = temp.divide(a0.multiply(a0));
309 
310         // recover original mean motion :
311         final T xn0dp = tle.getMeanMotion().multiply(60.0).divide(delta0.add(1.0));
312 
313         // Period >= 225 minutes is deep space
314         if (MathUtils.TWO_PI / (xn0dp.multiply(TLEConstants.MINUTES_PER_DAY).getReal()) >= (1.0 / 6.4)) {
315             return new FieldDeepSDP4<>(tle, attitudeProvider, mass, teme, parameters);
316         } else {
317             return new FieldSGP4<>(tle, attitudeProvider, mass, teme, parameters);
318         }
319     }
320 
321     /** Get the Earth gravity coefficient used for TLE propagation.
322      * @return the Earth gravity coefficient.
323      */
324     public static double getMU() {
325         return TLEConstants.MU;
326     }
327 
328     /** Get the extrapolated position and velocity from an initial TLE.
329      * @param date the final date
330      * @param parameters values of the model
331      * @return the final PVCoordinates
332      */
333     public FieldPVCoordinates<T> getPVCoordinates(final FieldAbsoluteDate<T> date, final T[] parameters) {
334 
335         sxpPropagate(date.durationFrom(tle.getDate()).divide(60.0), parameters);
336 
337         // Compute PV with previous calculated parameters
338         return computePVCoordinates();
339     }
340 
341     /** Computation of the first commons parameters.
342      * @param parameters SGP4 and SDP4 model parameters
343      */
344     private void initializeCommons(final T[] parameters) {
345 
346         final T zero = tle.getDate().getField().getZero();
347         final T bStar = parameters[0];
348         final T a1 = tle.getMeanMotion().multiply(60.0).reciprocal().multiply(TLEConstants.XKE).pow(TLEConstants.TWO_THIRD);
349         cosi0 = FastMath.cos(tle.getI());
350         theta2 = cosi0.multiply(cosi0);
351         final T x3thm1 = theta2.multiply(3.0).subtract(1.0);
352         e0sq = tle.getE().square();
353         beta02 = e0sq.negate().add(1.0);
354         beta0 = FastMath.sqrt(beta02);
355         final T tval = x3thm1.multiply(1.5 * TLEConstants.CK2).divide(beta0.multiply(beta02));
356         final T delta1 = tval.divide(a1.multiply(a1));
357         final T a0 = a1.multiply(delta1.multiply(
358                      delta1.multiply(134.0 / 81.0).add(1.0).multiply(delta1).add(TLEConstants.ONE_THIRD)).negate().add(1.0));
359         final T delta0 = tval.divide(a0.multiply(a0));
360 
361         // recover original mean motion and semi-major axis :
362         xn0dp = tle.getMeanMotion().multiply(60.0).divide(delta0.add(1.0));
363         a0dp = a0.divide(delta0.negate().add(1.0));
364 
365         // Values of s and qms2t :
366         s4 = zero.newInstance(TLEConstants.S);  // unmodified value for s
367         T q0ms24 = zero.newInstance(TLEConstants.QOMS2T); // unmodified value for q0ms2T
368 
369         perige = a0dp.multiply(tle.getE().negate().add(1.0)).subtract(TLEConstants.NORMALIZED_EQUATORIAL_RADIUS).multiply(
370                                                                                                 TLEConstants.EARTH_RADIUS); // perige
371 
372         //  For perigee below 156 km, the values of s and qoms2t are changed :
373         if (perige.getReal() < 156.0) {
374             if (perige.getReal() <= 98.0) {
375                 s4 = zero.newInstance(20.0);
376             } else {
377                 s4 = perige.subtract(78.0);
378             }
379             final T temp_val = s4.negate().add(120.0).multiply(TLEConstants.NORMALIZED_EQUATORIAL_RADIUS / TLEConstants.EARTH_RADIUS);
380             final T temp_val_squared = temp_val.multiply(temp_val);
381             q0ms24 = temp_val_squared.square();
382             s4 = s4.divide(TLEConstants.EARTH_RADIUS).add(TLEConstants.NORMALIZED_EQUATORIAL_RADIUS); // new value for q0ms2T and s
383         }
384 
385         final T pinv = a0dp.multiply(beta02).reciprocal();
386         final T pinvsq = pinv.square();
387         tsi = a0dp.subtract(s4).reciprocal();
388         eta = a0dp.multiply(tle.getE()).multiply(tsi);
389         etasq = eta.square();
390         eeta = tle.getE().multiply(eta);
391 
392         final T psisq = etasq.negate().add(1.0).abs(); // abs because pow 3.5 needs positive value
393         final T tsi_squared = tsi.multiply(tsi);
394         coef = q0ms24.multiply(tsi_squared.square());
395         coef1 = coef.divide(psisq.pow(3.5));
396 
397         // C2 and C1 coefficients computation :
398         c2 = coef1.multiply(xn0dp).multiply(a0dp.multiply(
399                            etasq.multiply(1.5).add(eeta.multiply(etasq.add(4.0))).add(1.0)).add(
400                            tsi.divide(psisq).multiply(x3thm1).multiply(0.75 * TLEConstants.CK2).multiply(
401                            etasq.multiply(etasq.add(8.0)).multiply(3.0).add(8.0))));
402         c1 = bStar.multiply(c2);
403         sini0 = FastMath.sin(tle.getI());
404 
405         final T x1mth2 = theta2.negate().add(1.0);
406 
407         // C4 coefficient computation :
408         c4 = xn0dp.multiply(coef1).multiply(a0dp).multiply(2.0).multiply(beta02).multiply(
409                            eta.multiply(etasq.multiply(0.5).add(2.0)).add(tle.getE().multiply(etasq.multiply(2.0).add(0.5))).subtract(
410                            tsi.divide(a0dp.multiply(psisq)).multiply(2 * TLEConstants.CK2).multiply(
411                            x3thm1.multiply(-3).multiply(etasq.multiply(eeta.multiply(-0.5).add(1.5)).add(eeta.multiply(-2.0)).add(1.0)).add(
412                            x1mth2.multiply(0.75).multiply(etasq.multiply(2.0).subtract(eeta.multiply(etasq.add(1.0)))).multiply(FastMath.cos(tle.getPerigeeArgument().multiply(2.0)))))));
413 
414         final T theta4 = theta2.multiply(theta2);
415         final T temp1  = pinvsq.multiply(xn0dp).multiply(3 * TLEConstants.CK2);
416         final T temp2  = temp1.multiply(pinvsq).multiply(TLEConstants.CK2);
417         final T temp3  = pinvsq.multiply(pinvsq).multiply(xn0dp).multiply(1.25 * TLEConstants.CK4);
418 
419         // atmospheric and gravitation coefs :(Mdf and OMEGAdf)
420         xmdot = xn0dp.add(
421                 temp1.multiply(0.5).multiply(beta0).multiply(x3thm1)).add(
422                 temp2.multiply(0.0625).multiply(beta0).multiply(
423                 theta2.multiply(78.0).negate().add(13.0).add(theta4.multiply(137.0))));
424 
425         final T x1m5th = theta2.multiply(5.0).negate().add(1.0);
426 
427         omgdot = temp1.multiply(-0.5).multiply(x1m5th).add(
428                  temp2.multiply(0.0625).multiply(theta2.multiply(114.0).negate().add(
429                  theta4.multiply(395.0)).add(7.0))).add(
430                  temp3.multiply(theta2.multiply(36.0).negate().add(theta4.multiply(49.0)).add(3.0)));
431 
432         final T xhdot1 = temp1.negate().multiply(cosi0);
433 
434         xnodot = xhdot1.add(temp2.multiply(0.5).multiply(theta2.multiply(19.0).negate().add(4.0)).add(
435                  temp3.multiply(2.0).multiply(theta2.multiply(7.0).negate().add(3.0))).multiply(cosi0));
436         xnodcf = beta02.multiply(xhdot1).multiply(c1).multiply(3.5);
437         t2cof = c1.multiply(1.5);
438 
439     }
440 
441     /** Retrieves the position and velocity.
442      * @return the computed PVCoordinates.
443      */
444     private FieldPVCoordinates<T> computePVCoordinates() {
445 
446         final T zero = tle.getDate().getField().getZero();
447         // Long period periodics
448         final T axn = e.multiply(FastMath.cos(omega));
449         T temp = a.multiply(e.multiply(e).negate().add(1.0)).reciprocal();
450         final T xlcof = sini0.multiply(0.125 * TLEConstants.A3OVK2).multiply(
451                              cosi0.multiply(5.0).add(3.0).divide(cosi0.add(1.0)));
452         final T aycof = sini0.multiply(0.25 * TLEConstants.A3OVK2);
453         final T xll   = temp.multiply(xlcof).multiply(axn);
454         final T aynl  = temp.multiply(aycof);
455         final T xlt   = xl.add(xll);
456         final T ayn   = e.multiply(FastMath.sin(omega)).add(aynl);
457         final T elsq  = axn.square().add(ayn.square());
458         final T capu  = MathUtils.normalizeAngle(xlt.subtract(xnode), zero.getPi());
459         T epw    = capu;
460         T ecosE  = zero;
461         T esinE  = zero;
462         T sinEPW = zero;
463         T cosEPW = zero;
464 
465         // Dundee changes:  items dependent on cosio get recomputed:
466         final T cosi0Sq = cosi0.square();
467         final T x3thm1  = cosi0Sq.multiply(3.0).subtract(1.0);
468         final T x1mth2  = cosi0Sq.negate().add(1.0);
469         final T x7thm1  = cosi0Sq.multiply(7.0).subtract(1.0);
470 
471         if (e.getReal() > (1 - 1e-6)) {
472             throw new OrekitException(OrekitMessages.TOO_LARGE_ECCENTRICITY_FOR_PROPAGATION_MODEL, e);
473         }
474 
475         // Solve Kepler's' Equation.
476         final double newtonRaphsonEpsilon = 1e-12;
477         for (int j = 0; j < 10; j++) {
478 
479             boolean doSecondOrderNewtonRaphson = true;
480 
481             sinEPW = FastMath.sin( epw);
482             cosEPW = FastMath.cos( epw);
483             ecosE  = axn.multiply(cosEPW).add(ayn.multiply(sinEPW));
484             esinE  = axn.multiply(sinEPW).subtract(ayn.multiply(cosEPW));
485             final T f = capu.subtract(epw).add(esinE);
486             if (FastMath.abs(f.getReal()) < newtonRaphsonEpsilon) {
487                 break;
488             }
489             final T fdot = ecosE.negate().add(1.0);
490             T delta_epw = f.divide(fdot);
491             if (j == 0) {
492                 final T maxNewtonRaphson = e.abs().multiply(1.25);
493                 doSecondOrderNewtonRaphson = false;
494                 if (delta_epw.getReal() > maxNewtonRaphson.getReal()) {
495                     delta_epw = maxNewtonRaphson;
496                 } else if (delta_epw.getReal() < -maxNewtonRaphson.getReal()) {
497                     delta_epw = maxNewtonRaphson.negate();
498                 } else {
499                     doSecondOrderNewtonRaphson = true;
500                 }
501             }
502             if (doSecondOrderNewtonRaphson) {
503                 delta_epw = f.divide(fdot.add(esinE.multiply(0.5).multiply(delta_epw)));
504             }
505             epw = epw.add(delta_epw);
506         }
507 
508         // Short period preliminary quantities
509         temp = elsq.negate().add(1.0);
510         final T pl = a.multiply(temp);
511         final T r  = a.multiply(ecosE.negate().add(1.0));
512         T temp2 = a.divide(r);
513         final T betal = FastMath.sqrt(temp);
514         temp = esinE.divide(betal.add(1.0));
515         final T cosu  = temp2.multiply(cosEPW.subtract(axn).add(ayn.multiply(temp)));
516         final T sinu  = temp2.multiply(sinEPW.subtract(ayn).subtract(axn.multiply(temp)));
517         final T u     = FastMath.atan2(sinu, cosu);
518         final T sin2u = sinu.multiply(cosu).multiply(2.0);
519         final T cos2u = cosu.multiply(cosu).multiply(2.0).subtract(1.0);
520         final T temp1 = pl.reciprocal().multiply(TLEConstants.CK2);
521         temp2         = temp1.divide(pl);
522 
523         // Update for short periodics
524         final T rk = r.multiply(temp2.multiply(betal).multiply(x3thm1).multiply(-1.5).add(1.0)).add(
525                      temp1.multiply(x1mth2).multiply(cos2u).multiply(0.5));
526         final T uk = u.subtract(temp2.multiply(x7thm1).multiply(sin2u).multiply(0.25));
527         final T xnodek = xnode.add(temp2.multiply(cosi0).multiply(sin2u).multiply(1.5));
528         final T xinck = i.add(temp2.multiply(cosi0).multiply(sini0).multiply(cos2u).multiply(1.5));
529 
530         // Orientation vectors
531         final T sinuk  = FastMath.sin(uk);
532         final T cosuk  = FastMath.cos(uk);
533         final T sinik  = FastMath.sin(xinck);
534         final T cosik  = FastMath.cos(xinck);
535         final T sinnok = FastMath.sin(xnodek);
536         final T cosnok = FastMath.cos(xnodek);
537         final T xmx    = sinnok.negate().multiply(cosik);
538         final T xmy    = cosnok.multiply(cosik);
539         final T ux     = xmx.multiply(sinuk).add(cosnok.multiply(cosuk));
540         final T uy     = xmy.multiply(sinuk).add(sinnok.multiply(cosuk));
541         final T uz     = sinik.multiply(sinuk);
542 
543         // Position and velocity
544         final T cr = rk.multiply(1000 * TLEConstants.EARTH_RADIUS);
545         final FieldVector3D<T> pos = new FieldVector3D<>(cr.multiply(ux), cr.multiply(uy), cr.multiply(uz));
546 
547         final T sqrtA  = FastMath.sqrt(a);
548         final T rdot   = sqrtA.multiply(esinE.divide(r)).multiply(TLEConstants.XKE);
549         final T rfdot  = FastMath.sqrt(pl).divide(r).multiply(TLEConstants.XKE);
550         final T xn     = a.multiply(sqrtA).reciprocal().multiply(TLEConstants.XKE);
551         final T rdotk  = rdot.subtract(xn.multiply(temp1).multiply(x1mth2).multiply(sin2u));
552         final T rfdotk = rfdot.add(xn.multiply(temp1).multiply(x1mth2.multiply(cos2u).add(x3thm1.multiply(1.5))));
553         final T vx     = xmx.multiply(cosuk).subtract(cosnok.multiply(sinuk));
554         final T vy     = xmy.multiply(cosuk).subtract(sinnok.multiply(sinuk));
555         final T vz     = sinik.multiply(cosuk);
556 
557         final double cv = 1000.0 * TLEConstants.EARTH_RADIUS / 60.0;
558         final FieldVector3D<T> vel = new FieldVector3D<>(rdotk.multiply(ux).add(rfdotk.multiply(vx)).multiply(cv),
559                                                           rdotk.multiply(uy).add(rfdotk.multiply(vy)).multiply(cv),
560                                                           rdotk.multiply(uz).add(rfdotk.multiply(vz)).multiply(cv));
561         return new FieldPVCoordinates<>(pos, vel);
562 
563     }
564 
565     /** {@inheritDoc} */
566     @Override
567     public List<ParameterDriver> getParametersDrivers() {
568         return tle.getParametersDrivers();
569     }
570 
571     /** Initialization proper to each propagator (SGP or SDP).
572      * @param parameters model parameters
573      */
574     protected abstract void sxpInitialize(T[] parameters);
575 
576     /** Propagation proper to each propagator (SGP or SDP).
577      * @param t the offset from initial epoch (min)
578      * @param parameters model parameters
579      */
580     protected abstract void sxpPropagate(T t, T[] parameters);
581 
582     /** {@inheritDoc}
583      * <p>
584      * For TLE propagator, calling this method is only recommended
585      * for covariance propagation when the new <code>state</code>
586      * differs from the previous one by only adding the additional
587      * state containing the derivatives.
588      * </p>
589      */
590     public void resetInitialState(final FieldSpacecraftState<T> state) {
591         super.resetInitialState(state);
592         resetTle(state);
593         masses = new HashMap<>();
594         masses.put(tle, state.getMass());
595         tles = new TimeSpanMap<>(tle);
596     }
597 
598     /** {@inheritDoc} */
599     protected void resetIntermediateState(final FieldSpacecraftState<T> state, final boolean forward) {
600         resetTle(state);
601         if (forward) {
602             tles.addValidAfter(tle, state.getDate().toAbsoluteDate(), false);
603         } else {
604             tles.addValidBefore(tle, state.getDate().toAbsoluteDate(), false);
605         }
606         stateChanged(state);
607         masses.put(tle, state.getMass());
608     }
609 
610     /** Reset internal TLE from a SpacecraftState.
611      * @param state spacecraft state on which to base new TLE
612      */
613     private void resetTle(final FieldSpacecraftState<T> state) {
614         final TleGenerationAlgorithm algorithm = TLEPropagator.getDefaultTleGenerationAlgorithm(utc, teme);
615         final FieldTLE<T> newTle = algorithm.generate(state, tle);
616         initializeTle(newTle);
617     }
618 
619     /** Initialize internal TLE.
620      * @param newTle tle to replace current one
621      */
622     private void initializeTle(final FieldTLE<T> newTle) {
623         tle = newTle;
624         final T[] parameters = tle.getParameters(tle.getDate().getField());
625         initializeCommons(parameters);
626         sxpInitialize(parameters);
627     }
628 
629     /** {@inheritDoc} */
630     protected T getMass(final FieldAbsoluteDate<T> date) {
631         return masses.get(tles.get(date.toAbsoluteDate()));
632     }
633 
634     /** {@inheritDoc} */
635     public FieldOrbit<T> propagateOrbit(final FieldAbsoluteDate<T> date, final T[] parameters) {
636         final FieldTLE<T> closestTle = tles.get(date.toAbsoluteDate());
637         if (!tle.equals(closestTle)) {
638             initializeTle(closestTle);
639         }
640         final T mu = date.getField().getZero().newInstance(TLEConstants.MU);
641         return new FieldCartesianOrbit<>(getPVCoordinates(date, parameters), teme, date, mu);
642     }
643 
644     /** Get the underlying TLE.
645      * If there has been calls to #resetInitialState or #resetIntermediateState,
646      * it will not be the same as given to the constructor.
647      * @return underlying TLE
648      */
649     public FieldTLE<T> getTLE() {
650         return tle;
651     }
652 
653     /** {@inheritDoc} */
654     public Frame getFrame() {
655         return teme;
656     }
657 
658 }