TLEPropagator.java

  1. /* Copyright 2002-2022 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.List;

  21. import org.hipparchus.geometry.euclidean.threed.Vector3D;
  22. import org.hipparchus.linear.RealMatrix;
  23. import org.hipparchus.util.FastMath;
  24. import org.hipparchus.util.MathUtils;
  25. import org.hipparchus.util.SinCos;
  26. import org.orekit.annotation.DefaultDataContext;
  27. import org.orekit.attitudes.Attitude;
  28. import org.orekit.attitudes.AttitudeProvider;
  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.CartesianOrbit;
  36. import org.orekit.orbits.Orbit;
  37. import org.orekit.propagation.AbstractMatricesHarvester;
  38. import org.orekit.propagation.MatricesHarvester;
  39. import org.orekit.propagation.SpacecraftState;
  40. import org.orekit.propagation.analytical.AbstractAnalyticalPropagator;
  41. import org.orekit.time.AbsoluteDate;
  42. import org.orekit.time.TimeScale;
  43. import org.orekit.utils.DoubleArrayDictionary;
  44. import org.orekit.utils.PVCoordinates;
  45. import org.orekit.utils.ParameterDriver;


  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(TLE)} 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.  * @see TLE
  70.  */
  71. public abstract class TLEPropagator extends AbstractAnalyticalPropagator {

  72.     // CHECKSTYLE: stop VisibilityModifier check

  73.     /** Initial state. */
  74.     protected TLE tle;

  75.     /** UTC time scale. */
  76.     protected final TimeScale utc;

  77.     /** final RAAN. */
  78.     protected double xnode;

  79.     /** final semi major axis. */
  80.     protected double a;

  81.     /** final eccentricity. */
  82.     protected double e;

  83.     /** final inclination. */
  84.     protected double i;

  85.     /** final perigee argument. */
  86.     protected double omega;

  87.     /** L from SPTRCK #3. */
  88.     protected double xl;

  89.     /** original recovered semi major axis. */
  90.     protected double a0dp;

  91.     /** original recovered mean motion. */
  92.     protected double xn0dp;

  93.     /** cosinus original inclination. */
  94.     protected double cosi0;

  95.     /** cos io squared. */
  96.     protected double theta2;

  97.     /** sinus original inclination. */
  98.     protected double sini0;

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

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

  103.     /** common parameter for raan (OMEGA) computation. */
  104.     protected double xnodot;

  105.     /** original eccentricity squared. */
  106.     protected double e0sq;
  107.     /** 1 - e2. */
  108.     protected double beta02;

  109.     /** sqrt (1 - e2). */
  110.     protected double beta0;

  111.     /** perigee, expressed in KM and ALTITUDE. */
  112.     protected double perige;

  113.     /** eta squared. */
  114.     protected double etasq;

  115.     /** original eccentricity * eta. */
  116.     protected double eeta;

  117.     /** s* new value for the contant s. */
  118.     protected double s4;

  119.     /** tsi from SPTRCK #3. */
  120.     protected double tsi;

  121.     /** eta from SPTRCK #3. */
  122.     protected double eta;

  123.     /** coef for SGP C3 computation. */
  124.     protected double coef;

  125.     /** coef for SGP C5 computation. */
  126.     protected double coef1;

  127.     /** C1 from SPTRCK #3. */
  128.     protected double c1;

  129.     /** C2 from SPTRCK #3. */
  130.     protected double c2;

  131.     /** C4 from SPTRCK #3. */
  132.     protected double c4;

  133.     /** common parameter for raan (OMEGA) computation. */
  134.     protected double xnodcf;

  135.     /** 3/2 * C1. */
  136.     protected double t2cof;

  137.     // CHECKSTYLE: resume VisibilityModifier check

  138.     /** TLE frame. */
  139.     private final Frame teme;

  140.     /** Spacecraft mass (kg). */
  141.     private final double mass;

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

  157.     /** Protected constructor for derived classes.
  158.      * @param initialTLE the unique TLE to propagate
  159.      * @param attitudeProvider provider for attitude computation
  160.      * @param mass spacecraft mass (kg)
  161.      * @param teme the TEME frame to use for propagation.
  162.      * @since 10.1
  163.      */
  164.     protected TLEPropagator(final TLE initialTLE,
  165.                             final AttitudeProvider attitudeProvider,
  166.                             final double mass,
  167.                             final Frame teme) {
  168.         super(attitudeProvider);
  169.         setStartDate(initialTLE.getDate());
  170.         this.tle       = initialTLE;
  171.         this.teme      = teme;
  172.         this.mass      = mass;
  173.         this.utc       = initialTLE.getUtc();

  174.         initializeCommons();
  175.         sxpInitialize();
  176.         // set the initial state
  177.         final Orbit orbit = propagateOrbit(initialTLE.getDate());
  178.         final Attitude attitude = attitudeProvider.getAttitude(orbit, orbit.getDate(), orbit.getFrame());
  179.         super.resetInitialState(new SpacecraftState(orbit, attitude, mass));
  180.     }

  181.     /** Selects the extrapolator to use with the selected TLE.
  182.      *
  183.      * <p>This method uses the {@link DataContext#getDefault() default data context}.
  184.      *
  185.      * @param tle the TLE to propagate.
  186.      * @return the correct propagator.
  187.      * @see #selectExtrapolator(TLE, Frames)
  188.      */
  189.     @DefaultDataContext
  190.     public static TLEPropagator selectExtrapolator(final TLE tle) {
  191.         return selectExtrapolator(tle, DataContext.getDefault().getFrames());
  192.     }

  193.     /** Selects the extrapolator to use with the selected TLE.
  194.      * @param tle the TLE to propagate.
  195.      * @param frames set of Frames to use in the propagator.
  196.      * @return the correct propagator.
  197.      * @since 10.1
  198.      */
  199.     public static TLEPropagator selectExtrapolator(final TLE tle, final Frames frames) {
  200.         return selectExtrapolator(
  201.                 tle,
  202.                 InertialProvider.of(frames.getTEME()),
  203.                 DEFAULT_MASS,
  204.                 frames.getTEME());
  205.     }

  206.     /** Selects the extrapolator to use with the selected TLE.
  207.      *
  208.      * <p>This method uses the {@link DataContext#getDefault() default data context}.
  209.      *
  210.      * @param tle the TLE to propagate.
  211.      * @param attitudeProvider provider for attitude computation
  212.      * @param mass spacecraft mass (kg)
  213.      * @return the correct propagator.
  214.      * @see #selectExtrapolator(TLE, AttitudeProvider, double, Frame)
  215.      */
  216.     @DefaultDataContext
  217.     public static TLEPropagator selectExtrapolator(final TLE tle, final AttitudeProvider attitudeProvider,
  218.                                                    final double mass) {
  219.         return selectExtrapolator(tle, attitudeProvider, mass,
  220.                 DataContext.getDefault().getFrames().getTEME());
  221.     }

  222.     /** Selects the extrapolator to use with the selected TLE.
  223.      * @param tle the TLE to propagate.
  224.      * @param attitudeProvider provider for attitude computation
  225.      * @param mass spacecraft mass (kg)
  226.      * @param teme the TEME frame to use for propagation.
  227.      * @return the correct propagator.
  228.      * @since 10.1
  229.      */
  230.     public static TLEPropagator selectExtrapolator(final TLE tle,
  231.                                                    final AttitudeProvider attitudeProvider,
  232.                                                    final double mass,
  233.                                                    final Frame teme) {

  234.         final double a1 = FastMath.pow( TLEConstants.XKE / (tle.getMeanMotion() * 60.0), TLEConstants.TWO_THIRD);
  235.         final double cosi0 = FastMath.cos(tle.getI());
  236.         final double temp = TLEConstants.CK2 * 1.5 * (3 * cosi0 * cosi0 - 1.0) *
  237.                             FastMath.pow(1.0 - tle.getE() * tle.getE(), -1.5);
  238.         final double delta1 = temp / (a1 * a1);
  239.         final double a0 = a1 * (1.0 - delta1 * (TLEConstants.ONE_THIRD + delta1 * (delta1 * 134.0 / 81.0 + 1.0)));
  240.         final double delta0 = temp / (a0 * a0);

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

  243.         // Period >= 225 minutes is deep space
  244.         if (MathUtils.TWO_PI / (xn0dp * TLEConstants.MINUTES_PER_DAY) >= (1.0 / 6.4)) {
  245.             return new DeepSDP4(tle, attitudeProvider, mass, teme);
  246.         } else {
  247.             return new SGP4(tle, attitudeProvider, mass, teme);
  248.         }
  249.     }

  250.     /** Get the Earth gravity coefficient used for TLE propagation.
  251.      * @return the Earth gravity coefficient.
  252.      */
  253.     public static double getMU() {
  254.         return TLEConstants.MU;
  255.     }

  256.     /** Get the extrapolated position and velocity from an initial TLE.
  257.      * @param date the final date
  258.      * @return the final PVCoordinates
  259.      */
  260.     public PVCoordinates getPVCoordinates(final AbsoluteDate date) {

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

  262.         // Compute PV with previous calculated parameters
  263.         return computePVCoordinates();
  264.     }

  265.     /** Computation of the first commons parameters.
  266.      */
  267.     private void initializeCommons() {

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

  270.         final double a1 = FastMath.pow(TLEConstants.XKE / (tle.getMeanMotion() * 60.0), TLEConstants.TWO_THIRD);
  271.         cosi0 = scI0.cos();
  272.         theta2 = cosi0 * cosi0;
  273.         final double x3thm1 = 3.0 * theta2 - 1.0;
  274.         e0sq = tle.getE() * tle.getE();
  275.         beta02 = 1.0 - e0sq;
  276.         beta0 = FastMath.sqrt(beta02);
  277.         final double tval = TLEConstants.CK2 * 1.5 * x3thm1 / (beta0 * beta02);
  278.         final double delta1 = tval / (a1 * a1);
  279.         final double a0 = a1 * (1.0 - delta1 * (TLEConstants.ONE_THIRD + delta1 * (1.0 + 134.0 / 81.0 * delta1)));
  280.         final double delta0 = tval / (a0 * a0);

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

  284.         // Values of s and qms2t :
  285.         s4 = TLEConstants.S;  // unmodified value for s
  286.         double q0ms24 = TLEConstants.QOMS2T; // unmodified value for q0ms2T

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

  288.         //  For perigee below 156 km, the values of s and qoms2t are changed :
  289.         if (perige < 156.0) {
  290.             if (perige <= 98.0) {
  291.                 s4 = 20.0;
  292.             } else {
  293.                 s4 = perige - 78.0;
  294.             }
  295.             final double temp_val = (120.0 - s4) * TLEConstants.NORMALIZED_EQUATORIAL_RADIUS / TLEConstants.EARTH_RADIUS;
  296.             final double temp_val_squared = temp_val * temp_val;
  297.             q0ms24 = temp_val_squared * temp_val_squared;
  298.             s4 = s4 / TLEConstants.EARTH_RADIUS + TLEConstants.NORMALIZED_EQUATORIAL_RADIUS; // new value for q0ms2T and s
  299.         }

  300.         final double pinv = 1.0 / (a0dp * beta02);
  301.         final double pinvsq = pinv * pinv;
  302.         tsi = 1.0 / (a0dp - s4);
  303.         eta = a0dp * tle.getE() * tsi;
  304.         etasq = eta * eta;
  305.         eeta = tle.getE() * eta;

  306.         final double psisq = FastMath.abs(1.0 - etasq); // abs because pow 3.5 needs positive value
  307.         final double tsi_squared = tsi * tsi;
  308.         coef = q0ms24 * tsi_squared * tsi_squared;
  309.         coef1 = coef / FastMath.pow(psisq, 3.5);

  310.         // C2 and C1 coefficients computation :
  311.         c2 = coef1 * xn0dp * (a0dp * (1.0 + 1.5 * etasq + eeta * (4.0 + etasq)) +
  312.              0.75 * TLEConstants.CK2 * tsi / psisq * x3thm1 * (8.0 + 3.0 * etasq * (8.0 + etasq)));
  313.         c1 = tle.getBStar() * c2;
  314.         sini0 = scI0.sin();

  315.         final double x1mth2 = 1.0 - theta2;

  316.         // C4 coefficient computation :
  317.         c4 = 2.0 * xn0dp * coef1 * a0dp * beta02 * (eta * (2.0 + 0.5 * etasq) +
  318.              tle.getE() * (0.5 + 2.0 * etasq) -
  319.              2 * TLEConstants.CK2 * tsi / (a0dp * psisq) *
  320.              (-3.0 * x3thm1 * (1.0 - 2.0 * eeta + etasq * (1.5 - 0.5 * eeta)) +
  321.               0.75 * x1mth2 * (2.0 * etasq - eeta * (1.0 + etasq)) * FastMath.cos(2.0 * tle.getPerigeeArgument())));

  322.         final double theta4 = theta2 * theta2;
  323.         final double temp1 = 3 * TLEConstants.CK2 * pinvsq * xn0dp;
  324.         final double temp2 = temp1 * TLEConstants.CK2 * pinvsq;
  325.         final double temp3 = 1.25 * TLEConstants.CK4 * pinvsq * pinvsq * xn0dp;

  326.         // atmospheric and gravitation coefs :(Mdf and OMEGAdf)
  327.         xmdot = xn0dp +
  328.                 0.5 * temp1 * beta0 * x3thm1 +
  329.                 0.0625 * temp2 * beta0 * (13.0 - 78.0 * theta2 + 137.0 * theta4);

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

  331.         omgdot = -0.5 * temp1 * x1m5th +
  332.                  0.0625 * temp2 * (7.0 - 114.0 * theta2 + 395.0 * theta4) +
  333.                  temp3 * (3.0 - 36.0 * theta2 + 49.0 * theta4);

  334.         final double xhdot1 = -temp1 * cosi0;

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

  338.     }

  339.     /** Retrieves the position and velocity.
  340.      * @return the computed PVCoordinates.
  341.      */
  342.     private PVCoordinates computePVCoordinates() {

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

  345.         // Long period periodics
  346.         final double axn = e * scOmega.cos();
  347.         double temp = 1.0 / (a * (1.0 - e * e));
  348.         final double xlcof = 0.125 * TLEConstants.A3OVK2 * sini0 * (3.0 + 5.0 * cosi0) / (1.0 + cosi0);
  349.         final double aycof = 0.25 * TLEConstants.A3OVK2 * sini0;
  350.         final double xll = temp * xlcof * axn;
  351.         final double aynl = temp * aycof;
  352.         final double xlt = xl + xll;
  353.         final double ayn = e * scOmega.sin() + aynl;
  354.         final double elsq = axn * axn + ayn * ayn;
  355.         final double capu = MathUtils.normalizeAngle(xlt - xnode, FastMath.PI);
  356.         double epw = capu;
  357.         double ecosE = 0;
  358.         double esinE = 0;
  359.         double sinEPW = 0;
  360.         double cosEPW = 0;

  361.         // Dundee changes:  items dependent on cosio get recomputed:
  362.         final double cosi0Sq = cosi0 * cosi0;
  363.         final double x3thm1 = 3.0 * cosi0Sq - 1.0;
  364.         final double x1mth2 = 1.0 - cosi0Sq;
  365.         final double x7thm1 = 7.0 * cosi0Sq - 1.0;

  366.         if (e > (1 - 1e-6)) {
  367.             throw new OrekitException(OrekitMessages.TOO_LARGE_ECCENTRICITY_FOR_PROPAGATION_MODEL, e);
  368.         }

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

  372.             boolean doSecondOrderNewtonRaphson = true;

  373.             final SinCos scEPW = FastMath.sinCos(epw);
  374.             sinEPW = scEPW.sin();
  375.             cosEPW = scEPW.cos();
  376.             ecosE = axn * cosEPW + ayn * sinEPW;
  377.             esinE = axn * sinEPW - ayn * cosEPW;
  378.             final double f = capu - epw + esinE;
  379.             if (FastMath.abs(f) < newtonRaphsonEpsilon) {
  380.                 break;
  381.             }
  382.             final double fdot = 1.0 - ecosE;
  383.             double delta_epw = f / fdot;
  384.             if (j == 0) {
  385.                 final double maxNewtonRaphson = 1.25 * FastMath.abs(e);
  386.                 doSecondOrderNewtonRaphson = false;
  387.                 if (delta_epw > maxNewtonRaphson) {
  388.                     delta_epw = maxNewtonRaphson;
  389.                 } else if (delta_epw < -maxNewtonRaphson) {
  390.                     delta_epw = -maxNewtonRaphson;
  391.                 } else {
  392.                     doSecondOrderNewtonRaphson = true;
  393.                 }
  394.             }
  395.             if (doSecondOrderNewtonRaphson) {
  396.                 delta_epw = f / (fdot + 0.5 * esinE * delta_epw);
  397.             }
  398.             epw += delta_epw;
  399.         }

  400.         // Short period preliminary quantities
  401.         temp = 1.0 - elsq;
  402.         final double pl = a * temp;
  403.         final double r = a * (1.0 - ecosE);
  404.         double temp2 = a / r;
  405.         final double betal = FastMath.sqrt(temp);
  406.         temp = esinE / (1.0 + betal);
  407.         final double cosu = temp2 * (cosEPW - axn + ayn * temp);
  408.         final double sinu = temp2 * (sinEPW - ayn - axn * temp);
  409.         final double u = FastMath.atan2(sinu, cosu);
  410.         final double sin2u = 2.0 * sinu * cosu;
  411.         final double cos2u = 2.0 * cosu * cosu - 1.0;
  412.         final double temp1 = TLEConstants.CK2 / pl;
  413.         temp2 = temp1 / pl;

  414.         // Update for short periodics
  415.         final double rk = r * (1.0 - 1.5 * temp2 * betal * x3thm1) + 0.5 * temp1 * x1mth2 * cos2u;
  416.         final double uk = u - 0.25 * temp2 * x7thm1 * sin2u;
  417.         final double xnodek = xnode + 1.5 * temp2 * cosi0 * sin2u;
  418.         final double xinck = i + 1.5 * temp2 * cosi0 * sini0 * cos2u;

  419.         // Orientation vectors
  420.         final SinCos scuk   = FastMath.sinCos(uk);
  421.         final SinCos scik   = FastMath.sinCos(xinck);
  422.         final SinCos scnok  = FastMath.sinCos(xnodek);
  423.         final double sinuk  = scuk.sin();
  424.         final double cosuk  = scuk.cos();
  425.         final double sinik  = scik.sin();
  426.         final double cosik  = scik.cos();
  427.         final double sinnok = scnok.sin();
  428.         final double cosnok = scnok.cos();
  429.         final double xmx = -sinnok * cosik;
  430.         final double xmy = cosnok * cosik;
  431.         final double ux  = xmx * sinuk + cosnok * cosuk;
  432.         final double uy  = xmy * sinuk + sinnok * cosuk;
  433.         final double uz  = sinik * sinuk;

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

  437.         final double rdot   = TLEConstants.XKE * FastMath.sqrt(a) * esinE / r;
  438.         final double rfdot  = TLEConstants.XKE * FastMath.sqrt(pl) / r;
  439.         final double xn     = TLEConstants.XKE / (a * FastMath.sqrt(a));
  440.         final double rdotk  = rdot - xn * temp1 * x1mth2 * sin2u;
  441.         final double rfdotk = rfdot + xn * temp1 * (x1mth2 * cos2u + 1.5 * x3thm1);
  442.         final double vx     = xmx * cosuk - cosnok * sinuk;
  443.         final double vy     = xmy * cosuk - sinnok * sinuk;
  444.         final double vz     = sinik * cosuk;

  445.         final double cv = 1000.0 * TLEConstants.EARTH_RADIUS / 60.0;
  446.         final Vector3D vel = new Vector3D(cv * (rdotk * ux + rfdotk * vx),
  447.                                           cv * (rdotk * uy + rfdotk * vy),
  448.                                           cv * (rdotk * uz + rfdotk * vz));

  449.         return new PVCoordinates(pos, vel);

  450.     }

  451.     /** Initialization proper to each propagator (SGP or SDP).
  452.      */
  453.     protected abstract void sxpInitialize();

  454.     /** Propagation proper to each propagator (SGP or SDP).
  455.      * @param t the offset from initial epoch (min)
  456.      */
  457.     protected abstract void sxpPropagate(double t);

  458.     /** {@inheritDoc}
  459.      * <p>
  460.      * For TLE propagator, calling this method is only recommended
  461.      * for covariance propagation when the new <code>state</code>
  462.      * differs from the previous one by only adding the additional
  463.      * state containing the derivatives.
  464.      * </p>
  465.      */
  466.     public void resetInitialState(final SpacecraftState state) {
  467.         super.resetInitialState(state);
  468.         super.setStartDate(state.getDate());
  469.         final TLE newTLE = TLE.stateToTLE(state, tle, utc, teme);
  470.         this.tle = newTLE;
  471.         initializeCommons();
  472.         sxpInitialize();
  473.     }

  474.     /** {@inheritDoc} */
  475.     protected void resetIntermediateState(final SpacecraftState state, final boolean forward) {
  476.         throw new OrekitException(OrekitMessages.NON_RESETABLE_STATE);
  477.     }

  478.     /** {@inheritDoc} */
  479.     protected double getMass(final AbsoluteDate date) {
  480.         return mass;
  481.     }

  482.     /** {@inheritDoc} */
  483.     protected Orbit propagateOrbit(final AbsoluteDate date) {
  484.         return new CartesianOrbit(getPVCoordinates(date), teme, date, TLEConstants.MU);
  485.     }

  486.     /** Get the underlying TLE.
  487.      * @return underlying TLE
  488.      */
  489.     public TLE getTLE() {
  490.         return tle;
  491.     }

  492.     /** {@inheritDoc} */
  493.     public Frame getFrame() {
  494.         return teme;
  495.     }

  496.     /** {@inheritDoc} */
  497.     @Override
  498.     protected AbstractMatricesHarvester createHarvester(final String stmName, final RealMatrix initialStm,
  499.                                                         final DoubleArrayDictionary initialJacobianColumns) {
  500.         // Create the harvester
  501.         final TLEHarvester harvester = new TLEHarvester(this, stmName, initialStm, initialJacobianColumns);
  502.         // Update the list of additional state provider
  503.         addAdditionalStateProvider(harvester);
  504.         // Return the configured harvester
  505.         return harvester;
  506.     }

  507.     /**
  508.      * Get the names of the parameters in the matrix returned by {@link MatricesHarvester#getParametersJacobian}.
  509.      * @return names of the parameters (i.e. columns) of the Jacobian matrix
  510.      */
  511.     protected List<String> getJacobiansColumnsNames() {
  512.         final List<String> columnsNames = new ArrayList<>();
  513.         for (final ParameterDriver driver : tle.getParametersDrivers()) {
  514.             if (driver.isSelected() && !columnsNames.contains(driver.getName())) {
  515.                 columnsNames.add(driver.getName());
  516.             }
  517.         }
  518.         Collections.sort(columnsNames);
  519.         return columnsNames;
  520.     }

  521. }