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