TLEPropagator.java

  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. import java.util.ArrayList;
  19. import java.util.Collections;
  20. import java.util.HashMap;
  21. import java.util.List;
  22. import java.util.Map;

  23. import org.hipparchus.geometry.euclidean.threed.Vector3D;
  24. import org.hipparchus.linear.RealMatrix;
  25. import org.hipparchus.util.FastMath;
  26. import org.hipparchus.util.MathUtils;
  27. import org.hipparchus.util.SinCos;
  28. import org.orekit.annotation.DefaultDataContext;
  29. import org.orekit.attitudes.Attitude;
  30. import org.orekit.attitudes.AttitudeProvider;
  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.CartesianOrbit;
  37. import org.orekit.orbits.Orbit;
  38. import org.orekit.propagation.AbstractMatricesHarvester;
  39. import org.orekit.propagation.MatricesHarvester;
  40. import org.orekit.propagation.SpacecraftState;
  41. import org.orekit.propagation.analytical.AbstractAnalyticalPropagator;
  42. import org.orekit.propagation.analytical.tle.generation.FixedPointTleGenerationAlgorithm;
  43. import org.orekit.propagation.analytical.tle.generation.TleGenerationAlgorithm;
  44. import org.orekit.time.AbsoluteDate;
  45. import org.orekit.time.TimeScale;
  46. import org.orekit.utils.DoubleArrayDictionary;
  47. import org.orekit.utils.PVCoordinates;
  48. import org.orekit.utils.ParameterDriver;
  49. import org.orekit.utils.TimeSpanMap;
  50. import org.orekit.utils.TimeSpanMap.Span;

  51. /** This class provides elements to propagate TLE's.
  52.  * <p>
  53.  * The models used are SGP4 and SDP4, initially proposed by NORAD as the unique convenient
  54.  * propagator for TLE's. Inputs and outputs of this propagator are only suited for
  55.  * NORAD two lines elements sets, since it uses estimations and mean values appropriate
  56.  * for TLE's only.
  57.  * </p>
  58.  * <p>
  59.  * Deep- or near- space propagator is selected internally according to NORAD recommendations
  60.  * so that the user has not to worry about the used computation methods. One instance is created
  61.  * for each TLE (this instance can only be get using {@link #selectExtrapolator(TLE)} method,
  62.  * and can compute {@link PVCoordinates position and velocity coordinates} at any
  63.  * time. Maximum accuracy is guaranteed in a 24h range period before and after the provided
  64.  * TLE epoch (of course this accuracy is not really measurable nor predictable: according to
  65.  * <a href="https://www.celestrak.com/">CelesTrak</a>, the precision is close to one kilometer
  66.  * and error won't probably rise above 2 km).
  67.  * </p>
  68.  * <p>This implementation is largely inspired from the paper and source code <a
  69.  * href="https://www.celestrak.com/publications/AIAA/2006-6753/">Revisiting Spacetrack
  70.  * Report #3</a> and is fully compliant with its results and tests cases.</p>
  71.  * @author Felix R. Hoots, Ronald L. Roehrich, December 1980 (original fortran)
  72.  * @author David A. Vallado, Paul Crawford, Richard Hujsak, T.S. Kelso (C++ translation and improvements)
  73.  * @author Fabien Maussion (java translation)
  74.  * @see TLE
  75.  */
  76. public abstract class TLEPropagator extends AbstractAnalyticalPropagator {

  77.     // CHECKSTYLE: stop VisibilityModifier check

  78.     /** Initial state. */
  79.     protected TLE tle;

  80.     /** UTC time scale. */
  81.     protected final TimeScale utc;

  82.     /** final RAAN. */
  83.     protected double xnode;

  84.     /** final semi major axis. */
  85.     protected double a;

  86.     /** final eccentricity. */
  87.     protected double e;

  88.     /** final inclination. */
  89.     protected double i;

  90.     /** final perigee argument. */
  91.     protected double omega;

  92.     /** L from SPTRCK #3. */
  93.     protected double xl;

  94.     /** original recovered semi major axis. */
  95.     protected double a0dp;

  96.     /** original recovered mean motion. */
  97.     protected double xn0dp;

  98.     /** cosinus original inclination. */
  99.     protected double cosi0;

  100.     /** cos io squared. */
  101.     protected double theta2;

  102.     /** sinus original inclination. */
  103.     protected double sini0;

  104.     /** common parameter for mean anomaly (M) computation. */
  105.     protected double xmdot;

  106.     /** common parameter for perigee argument (omega) computation. */
  107.     protected double omgdot;

  108.     /** common parameter for raan (OMEGA) computation. */
  109.     protected double xnodot;

  110.     /** original eccentricity squared. */
  111.     protected double e0sq;
  112.     /** 1 - e2. */
  113.     protected double beta02;

  114.     /** sqrt (1 - e2). */
  115.     protected double beta0;

  116.     /** perigee, expressed in KM and ALTITUDE. */
  117.     protected double perige;

  118.     /** eta squared. */
  119.     protected double etasq;

  120.     /** original eccentricity * eta. */
  121.     protected double eeta;

  122.     /** s* new value for the contant s. */
  123.     protected double s4;

  124.     /** tsi from SPTRCK #3. */
  125.     protected double tsi;

  126.     /** eta from SPTRCK #3. */
  127.     protected double eta;

  128.     /** coef for SGP C3 computation. */
  129.     protected double coef;

  130.     /** coef for SGP C5 computation. */
  131.     protected double coef1;

  132.     /** C1 from SPTRCK #3. */
  133.     protected double c1;

  134.     /** C2 from SPTRCK #3. */
  135.     protected double c2;

  136.     /** C4 from SPTRCK #3. */
  137.     protected double c4;

  138.     /** common parameter for raan (OMEGA) computation. */
  139.     protected double xnodcf;

  140.     /** 3/2 * C1. */
  141.     protected double t2cof;

  142.     // CHECKSTYLE: resume VisibilityModifier check

  143.     /** TLE frame. */
  144.     private final Frame teme;

  145.     /** Spacecraft masses (kg) mapped to TLEs. */
  146.     private Map<TLE, Double> masses;

  147.     /** All TLEs. */
  148.     private TimeSpanMap<TLE> tles;

  149.     /** Protected constructor for derived classes.
  150.      *
  151.      * <p>This constructor uses the {@link DataContext#getDefault() default data context}.
  152.      *
  153.      * @param initialTLE the unique TLE to propagate
  154.      * @param attitudeProvider provider for attitude computation
  155.      * @param mass spacecraft mass (kg)
  156.      * @see #TLEPropagator(TLE, AttitudeProvider, double, Frame)
  157.      */
  158.     @DefaultDataContext
  159.     protected TLEPropagator(final TLE initialTLE, final AttitudeProvider attitudeProvider,
  160.                             final double mass) {
  161.         this(initialTLE, attitudeProvider, mass,
  162.                 DataContext.getDefault().getFrames().getTEME());
  163.     }

  164.     /** Protected constructor for derived classes.
  165.      * @param initialTLE the unique TLE to propagate
  166.      * @param attitudeProvider provider for attitude computation
  167.      * @param mass spacecraft mass (kg)
  168.      * @param teme the TEME frame to use for propagation.
  169.      * @since 10.1
  170.      */
  171.     protected TLEPropagator(final TLE initialTLE,
  172.                             final AttitudeProvider attitudeProvider,
  173.                             final double mass,
  174.                             final Frame teme) {
  175.         super(attitudeProvider);
  176.         setStartDate(initialTLE.getDate());
  177.         this.utc       = initialTLE.getUtc();
  178.         initializeTle(initialTLE);
  179.         this.teme      = teme;
  180.         this.tles      = new TimeSpanMap<>(tle);
  181.         this.masses    = new HashMap<>();
  182.         this.masses.put(tle, mass);

  183.         // set the initial state
  184.         final Orbit orbit = propagateOrbit(initialTLE.getDate());
  185.         final Attitude attitude = attitudeProvider.getAttitude(orbit, orbit.getDate(), orbit.getFrame());
  186.         super.resetInitialState(new SpacecraftState(orbit, attitude, mass));
  187.     }

  188.     /** Selects the extrapolator to use with the selected TLE.
  189.      *
  190.      * <p>This method uses the {@link DataContext#getDefault() default data context}.
  191.      *
  192.      * @param tle the TLE to propagate.
  193.      * @return the correct propagator.
  194.      * @see #selectExtrapolator(TLE, Frame)
  195.      */
  196.     @DefaultDataContext
  197.     public static TLEPropagator selectExtrapolator(final TLE tle) {
  198.         return selectExtrapolator(tle, DataContext.getDefault().getFrames().getTEME());
  199.     }

  200.     /** Selects the extrapolator to use with the selected TLE.
  201.      * @param tle the TLE to propagate.
  202.      * @param teme TEME frame.
  203.      * @return the correct propagator.
  204.      * @since 10.1
  205.      * @see #selectExtrapolator(TLE, Frame, AttitudeProvider)
  206.      */
  207.     public static TLEPropagator selectExtrapolator(final TLE tle, final Frame teme) {
  208.         return selectExtrapolator(tle, teme, FrameAlignedProvider.of(teme));
  209.     }

  210.     /** Selects the extrapolator to use with the selected TLE.
  211.      * @param tle the TLE to propagate.
  212.      * @param teme TEME frame.
  213.      * @param attitudeProvider provider for attitude computation
  214.      * @return the correct propagator.
  215.      * @since 12.2
  216.      */
  217.     public static TLEPropagator selectExtrapolator(final TLE tle, final Frame teme, final AttitudeProvider attitudeProvider) {
  218.         return selectExtrapolator(tle, attitudeProvider, DEFAULT_MASS, teme);
  219.     }

  220.     /** Selects the extrapolator to use with the selected TLE.
  221.      *
  222.      * <p>This method uses the {@link DataContext#getDefault() default data context}.
  223.      *
  224.      * @param tle the TLE to propagate.
  225.      * @param attitudeProvider provider for attitude computation
  226.      * @param mass spacecraft mass (kg)
  227.      * @return the correct propagator.
  228.      * @see #selectExtrapolator(TLE, AttitudeProvider, double, Frame)
  229.      */
  230.     @DefaultDataContext
  231.     public static TLEPropagator selectExtrapolator(final TLE tle, final AttitudeProvider attitudeProvider,
  232.                                                    final double mass) {
  233.         return selectExtrapolator(tle, attitudeProvider, mass,
  234.                                   DataContext.getDefault().getFrames().getTEME());
  235.     }

  236.     /** Selects the extrapolator to use with the selected TLE.
  237.      * @param tle the TLE to propagate.
  238.      * @param attitudeProvider provider for attitude computation
  239.      * @param mass spacecraft mass (kg)
  240.      * @param teme the TEME frame to use for propagation.
  241.      * @return the correct propagator.
  242.      * @since 10.1
  243.      */
  244.     public static TLEPropagator selectExtrapolator(final TLE tle,
  245.                                                    final AttitudeProvider attitudeProvider,
  246.                                                    final double mass,
  247.                                                    final Frame teme) {

  248.         final double a1 = FastMath.pow( TLEConstants.XKE / (tle.getMeanMotion() * 60.0), TLEConstants.TWO_THIRD);
  249.         final double cosi0 = FastMath.cos(tle.getI());
  250.         final double temp = TLEConstants.CK2 * 1.5 * (3 * cosi0 * cosi0 - 1.0) *
  251.                             FastMath.pow(1.0 - tle.getE() * tle.getE(), -1.5);
  252.         final double delta1 = temp / (a1 * a1);
  253.         final double a0 = a1 * (1.0 - delta1 * (TLEConstants.ONE_THIRD + delta1 * (delta1 * 134.0 / 81.0 + 1.0)));
  254.         final double delta0 = temp / (a0 * a0);

  255.         // recover original mean motion :
  256.         final double xn0dp = tle.getMeanMotion() * 60.0 / (delta0 + 1.0);

  257.         // Period >= 225 minutes is deep space
  258.         if (MathUtils.TWO_PI / (xn0dp * TLEConstants.MINUTES_PER_DAY) >= (1.0 / 6.4)) {
  259.             return new DeepSDP4(tle, attitudeProvider, mass, teme);
  260.         } else {
  261.             return new SGP4(tle, attitudeProvider, mass, teme);
  262.         }
  263.     }

  264.     /** Get the Earth gravity coefficient used for TLE propagation.
  265.      * @return the Earth gravity coefficient.
  266.      */
  267.     public static double getMU() {
  268.         return TLEConstants.MU;
  269.     }

  270.     /** Get the extrapolated position and velocity from an initial TLE.
  271.      * @param date the final date
  272.      * @return the final PVCoordinates
  273.      */
  274.     public PVCoordinates getPVCoordinates(final AbsoluteDate date) {

  275.         sxpPropagate(date.durationFrom(tle.getDate()) / 60.0);

  276.         // Compute PV with previous calculated parameters
  277.         return computePVCoordinates();
  278.     }

  279.     /** Computation of the first commons parameters.
  280.      */
  281.     private void initializeCommons() {

  282.         // Sine and cosine of inclination
  283.         final SinCos scI0 = FastMath.sinCos(tle.getI());

  284.         final double a1 = FastMath.pow(TLEConstants.XKE / (tle.getMeanMotion() * 60.0), TLEConstants.TWO_THIRD);
  285.         cosi0 = scI0.cos();
  286.         theta2 = cosi0 * cosi0;
  287.         final double x3thm1 = 3.0 * theta2 - 1.0;
  288.         e0sq = tle.getE() * tle.getE();
  289.         beta02 = 1.0 - e0sq;
  290.         beta0 = FastMath.sqrt(beta02);
  291.         final double tval = TLEConstants.CK2 * 1.5 * x3thm1 / (beta0 * beta02);
  292.         final double delta1 = tval / (a1 * a1);
  293.         final double a0 = a1 * (1.0 - delta1 * (TLEConstants.ONE_THIRD + delta1 * (1.0 + 134.0 / 81.0 * delta1)));
  294.         final double delta0 = tval / (a0 * a0);

  295.         // recover original mean motion and semi-major axis :
  296.         xn0dp = tle.getMeanMotion() * 60.0 / (delta0 + 1.0);
  297.         a0dp = a0 / (1.0 - delta0);

  298.         // Values of s and qms2t :
  299.         s4 = TLEConstants.S;  // unmodified value for s
  300.         double q0ms24 = TLEConstants.QOMS2T; // unmodified value for q0ms2T

  301.         perige = (a0dp * (1 - tle.getE()) - TLEConstants.NORMALIZED_EQUATORIAL_RADIUS) * TLEConstants.EARTH_RADIUS; // perige

  302.         //  For perigee below 156 km, the values of s and qoms2t are changed :
  303.         if (perige < 156.0) {
  304.             if (perige <= 98.0) {
  305.                 s4 = 20.0;
  306.             } else {
  307.                 s4 = perige - 78.0;
  308.             }
  309.             final double temp_val = (120.0 - s4) * TLEConstants.NORMALIZED_EQUATORIAL_RADIUS / TLEConstants.EARTH_RADIUS;
  310.             final double temp_val_squared = temp_val * temp_val;
  311.             q0ms24 = temp_val_squared * temp_val_squared;
  312.             s4 = s4 / TLEConstants.EARTH_RADIUS + TLEConstants.NORMALIZED_EQUATORIAL_RADIUS; // new value for q0ms2T and s
  313.         }

  314.         final double pinv = 1.0 / (a0dp * beta02);
  315.         final double pinvsq = pinv * pinv;
  316.         tsi = 1.0 / (a0dp - s4);
  317.         eta = a0dp * tle.getE() * tsi;
  318.         etasq = eta * eta;
  319.         eeta = tle.getE() * eta;

  320.         final double psisq = FastMath.abs(1.0 - etasq); // abs because pow 3.5 needs positive value
  321.         final double tsi_squared = tsi * tsi;
  322.         coef = q0ms24 * tsi_squared * tsi_squared;
  323.         coef1 = coef / FastMath.pow(psisq, 3.5);

  324.         // C2 and C1 coefficients computation :
  325.         c2 = coef1 * xn0dp * (a0dp * (1.0 + 1.5 * etasq + eeta * (4.0 + etasq)) +
  326.              0.75 * TLEConstants.CK2 * tsi / psisq * x3thm1 * (8.0 + 3.0 * etasq * (8.0 + etasq)));
  327.         c1 = tle.getBStar() * c2;
  328.         sini0 = scI0.sin();

  329.         final double x1mth2 = 1.0 - theta2;

  330.         // C4 coefficient computation :
  331.         c4 = 2.0 * xn0dp * coef1 * a0dp * beta02 * (eta * (2.0 + 0.5 * etasq) +
  332.              tle.getE() * (0.5 + 2.0 * etasq) -
  333.              2 * TLEConstants.CK2 * tsi / (a0dp * psisq) *
  334.              (-3.0 * x3thm1 * (1.0 - 2.0 * eeta + etasq * (1.5 - 0.5 * eeta)) +
  335.               0.75 * x1mth2 * (2.0 * etasq - eeta * (1.0 + etasq)) * FastMath.cos(2.0 * tle.getPerigeeArgument())));

  336.         final double theta4 = theta2 * theta2;
  337.         final double temp1 = 3 * TLEConstants.CK2 * pinvsq * xn0dp;
  338.         final double temp2 = temp1 * TLEConstants.CK2 * pinvsq;
  339.         final double temp3 = 1.25 * TLEConstants.CK4 * pinvsq * pinvsq * xn0dp;

  340.         // atmospheric and gravitation coefs :(Mdf and OMEGAdf)
  341.         xmdot = xn0dp +
  342.                 0.5 * temp1 * beta0 * x3thm1 +
  343.                 0.0625 * temp2 * beta0 * (13.0 - 78.0 * theta2 + 137.0 * theta4);

  344.         final double x1m5th = 1.0 - 5.0 * theta2;

  345.         omgdot = -0.5 * temp1 * x1m5th +
  346.                  0.0625 * temp2 * (7.0 - 114.0 * theta2 + 395.0 * theta4) +
  347.                  temp3 * (3.0 - 36.0 * theta2 + 49.0 * theta4);

  348.         final double xhdot1 = -temp1 * cosi0;

  349.         xnodot = xhdot1 + (0.5 * temp2 * (4.0 - 19.0 * theta2) + 2.0 * temp3 * (3.0 - 7.0 * theta2)) * cosi0;
  350.         xnodcf = 3.5 * beta02 * xhdot1 * c1;
  351.         t2cof = 1.5 * c1;

  352.     }

  353.     /** Retrieves the position and velocity.
  354.      * @return the computed PVCoordinates.
  355.      */
  356.     private PVCoordinates computePVCoordinates() {

  357.         // Sine and cosine of final perigee argument
  358.         final SinCos scOmega = FastMath.sinCos(omega);

  359.         // Long period periodics
  360.         final double axn = e * scOmega.cos();
  361.         double temp = 1.0 / (a * (1.0 - e * e));
  362.         final double xlcof = 0.125 * TLEConstants.A3OVK2 * sini0 * (3.0 + 5.0 * cosi0) / (1.0 + cosi0);
  363.         final double aycof = 0.25 * TLEConstants.A3OVK2 * sini0;
  364.         final double xll = temp * xlcof * axn;
  365.         final double aynl = temp * aycof;
  366.         final double xlt = xl + xll;
  367.         final double ayn = e * scOmega.sin() + aynl;
  368.         final double elsq = axn * axn + ayn * ayn;
  369.         final double capu = MathUtils.normalizeAngle(xlt - xnode, FastMath.PI);
  370.         double epw = capu;
  371.         double ecosE = 0;
  372.         double esinE = 0;
  373.         double sinEPW = 0;
  374.         double cosEPW = 0;

  375.         // Dundee changes:  items dependent on cosio get recomputed:
  376.         final double cosi0Sq = cosi0 * cosi0;
  377.         final double x3thm1 = 3.0 * cosi0Sq - 1.0;
  378.         final double x1mth2 = 1.0 - cosi0Sq;
  379.         final double x7thm1 = 7.0 * cosi0Sq - 1.0;

  380.         if (e > (1 - 1e-6)) {
  381.             throw new OrekitException(OrekitMessages.TOO_LARGE_ECCENTRICITY_FOR_PROPAGATION_MODEL, e);
  382.         }

  383.         // Solve Kepler's' Equation.
  384.         final double newtonRaphsonEpsilon = 1e-12;
  385.         for (int j = 0; j < 10; j++) {

  386.             boolean doSecondOrderNewtonRaphson = true;

  387.             final SinCos scEPW = FastMath.sinCos(epw);
  388.             sinEPW = scEPW.sin();
  389.             cosEPW = scEPW.cos();
  390.             ecosE = axn * cosEPW + ayn * sinEPW;
  391.             esinE = axn * sinEPW - ayn * cosEPW;
  392.             final double f = capu - epw + esinE;
  393.             if (FastMath.abs(f) < newtonRaphsonEpsilon) {
  394.                 break;
  395.             }
  396.             final double fdot = 1.0 - ecosE;
  397.             double delta_epw = f / fdot;
  398.             if (j == 0) {
  399.                 final double maxNewtonRaphson = 1.25 * FastMath.abs(e);
  400.                 doSecondOrderNewtonRaphson = false;
  401.                 if (delta_epw > maxNewtonRaphson) {
  402.                     delta_epw = maxNewtonRaphson;
  403.                 } else if (delta_epw < -maxNewtonRaphson) {
  404.                     delta_epw = -maxNewtonRaphson;
  405.                 } else {
  406.                     doSecondOrderNewtonRaphson = true;
  407.                 }
  408.             }
  409.             if (doSecondOrderNewtonRaphson) {
  410.                 delta_epw = f / (fdot + 0.5 * esinE * delta_epw);
  411.             }
  412.             epw += delta_epw;
  413.         }

  414.         // Short period preliminary quantities
  415.         temp = 1.0 - elsq;
  416.         final double pl = a * temp;
  417.         final double r = a * (1.0 - ecosE);
  418.         double temp2 = a / r;
  419.         final double betal = FastMath.sqrt(temp);
  420.         temp = esinE / (1.0 + betal);
  421.         final double cosu = temp2 * (cosEPW - axn + ayn * temp);
  422.         final double sinu = temp2 * (sinEPW - ayn - axn * temp);
  423.         final double u = FastMath.atan2(sinu, cosu);
  424.         final double sin2u = 2.0 * sinu * cosu;
  425.         final double cos2u = 2.0 * cosu * cosu - 1.0;
  426.         final double temp1 = TLEConstants.CK2 / pl;
  427.         temp2 = temp1 / pl;

  428.         // Update for short periodics
  429.         final double rk = r * (1.0 - 1.5 * temp2 * betal * x3thm1) + 0.5 * temp1 * x1mth2 * cos2u;
  430.         final double uk = u - 0.25 * temp2 * x7thm1 * sin2u;
  431.         final double xnodek = xnode + 1.5 * temp2 * cosi0 * sin2u;
  432.         final double xinck = i + 1.5 * temp2 * cosi0 * sini0 * cos2u;

  433.         // Orientation vectors
  434.         final SinCos scuk   = FastMath.sinCos(uk);
  435.         final SinCos scik   = FastMath.sinCos(xinck);
  436.         final SinCos scnok  = FastMath.sinCos(xnodek);
  437.         final double sinuk  = scuk.sin();
  438.         final double cosuk  = scuk.cos();
  439.         final double sinik  = scik.sin();
  440.         final double cosik  = scik.cos();
  441.         final double sinnok = scnok.sin();
  442.         final double cosnok = scnok.cos();
  443.         final double xmx = -sinnok * cosik;
  444.         final double xmy = cosnok * cosik;
  445.         final double ux  = xmx * sinuk + cosnok * cosuk;
  446.         final double uy  = xmy * sinuk + sinnok * cosuk;
  447.         final double uz  = sinik * sinuk;

  448.         // Position and velocity
  449.         final double cr = 1000 * rk * TLEConstants.EARTH_RADIUS;
  450.         final Vector3D pos = new Vector3D(cr * ux, cr * uy, cr * uz);

  451.         final double rdot   = TLEConstants.XKE * FastMath.sqrt(a) * esinE / r;
  452.         final double rfdot  = TLEConstants.XKE * FastMath.sqrt(pl) / r;
  453.         final double xn     = TLEConstants.XKE / (a * FastMath.sqrt(a));
  454.         final double rdotk  = rdot - xn * temp1 * x1mth2 * sin2u;
  455.         final double rfdotk = rfdot + xn * temp1 * (x1mth2 * cos2u + 1.5 * x3thm1);
  456.         final double vx     = xmx * cosuk - cosnok * sinuk;
  457.         final double vy     = xmy * cosuk - sinnok * sinuk;
  458.         final double vz     = sinik * cosuk;

  459.         final double cv = 1000.0 * TLEConstants.EARTH_RADIUS / 60.0;
  460.         final Vector3D vel = new Vector3D(cv * (rdotk * ux + rfdotk * vx),
  461.                                           cv * (rdotk * uy + rfdotk * vy),
  462.                                           cv * (rdotk * uz + rfdotk * vz));

  463.         return new PVCoordinates(pos, vel);

  464.     }

  465.     /** Initialization proper to each propagator (SGP or SDP).
  466.      */
  467.     protected abstract void sxpInitialize();

  468.     /** Propagation proper to each propagator (SGP or SDP).
  469.      * @param t the offset from initial epoch (min)
  470.      */
  471.     protected abstract void sxpPropagate(double t);

  472.     /** {@inheritDoc}
  473.      * <p>
  474.      * For TLE propagator, calling this method is only recommended
  475.      * for covariance propagation when the new <code>state</code>
  476.      * differs from the previous one by only adding the additional
  477.      * state containing the derivatives.
  478.      * </p>
  479.      */
  480.     public void resetInitialState(final SpacecraftState state) {
  481.         super.resetInitialState(state);
  482.         resetTle(state);
  483.         masses = new HashMap<>();
  484.         masses.put(tle, state.getMass());
  485.         tles = new TimeSpanMap<>(tle);
  486.     }

  487.     /** {@inheritDoc} */
  488.     protected void resetIntermediateState(final SpacecraftState state, final boolean forward) {
  489.         resetTle(state);
  490.         if (forward) {
  491.             tles.addValidAfter(tle, state.getDate(), false);
  492.         } else {
  493.             tles.addValidBefore(tle, state.getDate(), false);
  494.         }
  495.         stateChanged(state);
  496.         masses.put(tle, state.getMass());
  497.     }

  498.     /** Reset internal TLE from a SpacecraftState.
  499.      * @param state spacecraft state on which to base new TLE
  500.      */
  501.     private void resetTle(final SpacecraftState state) {
  502.         final TleGenerationAlgorithm algorithm = getDefaultTleGenerationAlgorithm(utc, teme);
  503.         final TLE newTle = algorithm.generate(state, tle);
  504.         initializeTle(newTle);
  505.     }

  506.     /** Initialize internal TLE.
  507.      * @param newTle tle to replace current one
  508.      */
  509.     private void initializeTle(final TLE newTle) {
  510.         tle = newTle;
  511.         initializeCommons();
  512.         sxpInitialize();
  513.     }

  514.     /** {@inheritDoc} */
  515.     protected double getMass(final AbsoluteDate date) {
  516.         return masses.get(tles.get(date));
  517.     }

  518.     /** {@inheritDoc} */
  519.     protected Orbit propagateOrbit(final AbsoluteDate date) {
  520.         final TLE closestTle = tles.get(date);
  521.         if (!tle.equals(closestTle)) {
  522.             initializeTle(closestTle);
  523.         }
  524.         return new CartesianOrbit(getPVCoordinates(date), teme, date, TLEConstants.MU);
  525.     }

  526.     /** Get the underlying TLE.
  527.      * If there has been calls to #resetInitialState or #resetIntermediateState,
  528.      * it will not be the same as given to the constructor.
  529.      * @return underlying TLE
  530.      */
  531.     public TLE getTLE() {
  532.         return tle;
  533.     }

  534.     /** {@inheritDoc} */
  535.     public Frame getFrame() {
  536.         return teme;
  537.     }

  538.     /** {@inheritDoc} */
  539.     @Override
  540.     protected AbstractMatricesHarvester createHarvester(final String stmName, final RealMatrix initialStm,
  541.                                                         final DoubleArrayDictionary initialJacobianColumns) {
  542.         // Create the harvester
  543.         final TLEHarvester harvester = new TLEHarvester(this, stmName, initialStm, initialJacobianColumns);
  544.         // Update the list of additional state provider
  545.         addAdditionalStateProvider(harvester);
  546.         // Return the configured harvester
  547.         return harvester;
  548.     }

  549.     /**
  550.      * Get the names of the parameters in the matrix returned by {@link MatricesHarvester#getParametersJacobian}.
  551.      * @return names of the parameters (i.e. columns) of the Jacobian matrix
  552.      */
  553.     protected List<String> getJacobiansColumnsNames() {
  554.         final List<String> columnsNames = new ArrayList<>();
  555.         for (final ParameterDriver driver : tle.getParametersDrivers()) {

  556.             if (driver.isSelected() && !columnsNames.contains(driver.getNamesSpanMap().getFirstSpan().getData())) {
  557.                 // As driver with same name should have same NamesSpanMap we only check the if condition on the
  558.                 // first span map and then if the condition is OK all the span names are added to the jacobian column names
  559.                 for (Span<String> span = driver.getNamesSpanMap().getFirstSpan(); span != null; span = span.next()) {
  560.                     columnsNames.add(span.getData());
  561.                 }
  562.             }
  563.         }
  564.         Collections.sort(columnsNames);
  565.         return columnsNames;
  566.     }

  567.     /**
  568.      * Get the default TLE generation algorithm.
  569.      * @param utc UTC time scale
  570.      * @param teme TEME frame
  571.      * @return a TLE generation algorithm
  572.      * @since 12.0
  573.      */
  574.     public static TleGenerationAlgorithm getDefaultTleGenerationAlgorithm(final TimeScale utc, final Frame teme) {
  575.         return new FixedPointTleGenerationAlgorithm(FixedPointTleGenerationAlgorithm.EPSILON_DEFAULT,
  576.                                                     FixedPointTleGenerationAlgorithm.MAX_ITERATIONS_DEFAULT,
  577.                                                     FixedPointTleGenerationAlgorithm.SCALE_DEFAULT, utc, teme);
  578.     }

  579. }