TLEPropagator.java

  1. /* Copyright 2002-2020 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 org.hipparchus.geometry.euclidean.threed.Vector3D;
  19. import org.hipparchus.util.FastMath;
  20. import org.hipparchus.util.MathUtils;
  21. import org.orekit.annotation.DefaultDataContext;
  22. import org.orekit.attitudes.Attitude;
  23. import org.orekit.attitudes.AttitudeProvider;
  24. import org.orekit.data.DataContext;
  25. import org.orekit.errors.OrekitException;
  26. import org.orekit.errors.OrekitMessages;
  27. import org.orekit.frames.Frame;
  28. import org.orekit.frames.Frames;
  29. import org.orekit.orbits.CartesianOrbit;
  30. import org.orekit.orbits.Orbit;
  31. import org.orekit.propagation.Propagator;
  32. import org.orekit.propagation.SpacecraftState;
  33. import org.orekit.propagation.analytical.AbstractAnalyticalPropagator;
  34. import org.orekit.time.AbsoluteDate;
  35. import org.orekit.time.TimeScale;
  36. import org.orekit.utils.PVCoordinates;


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

  63.     // CHECKSTYLE: stop VisibilityModifier check

  64.     /** Initial state. */
  65.     protected final TLE tle;

  66.     /** UTC time scale. */
  67.     protected final TimeScale utc;

  68.     /** final RAAN. */
  69.     protected double xnode;

  70.     /** final semi major axis. */
  71.     protected double a;

  72.     /** final eccentricity. */
  73.     protected double e;

  74.     /** final inclination. */
  75.     protected double i;

  76.     /** final perigee argument. */
  77.     protected double omega;

  78.     /** L from SPTRCK #3. */
  79.     protected double xl;

  80.     /** original recovered semi major axis. */
  81.     protected double a0dp;

  82.     /** original recovered mean motion. */
  83.     protected double xn0dp;

  84.     /** cosinus original inclination. */
  85.     protected double cosi0;

  86.     /** cos io squared. */
  87.     protected double theta2;

  88.     /** sinus original inclination. */
  89.     protected double sini0;

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

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

  94.     /** common parameter for raan (OMEGA) computation. */
  95.     protected double xnodot;

  96.     /** original eccentricity squared. */
  97.     protected double e0sq;
  98.     /** 1 - e2. */
  99.     protected double beta02;

  100.     /** sqrt (1 - e2). */
  101.     protected double beta0;

  102.     /** perigee, expressed in KM and ALTITUDE. */
  103.     protected double perige;

  104.     /** eta squared. */
  105.     protected double etasq;

  106.     /** original eccentricity * eta. */
  107.     protected double eeta;

  108.     /** s* new value for the contant s. */
  109.     protected double s4;

  110.     /** tsi from SPTRCK #3. */
  111.     protected double tsi;

  112.     /** eta from SPTRCK #3. */
  113.     protected double eta;

  114.     /** coef for SGP C3 computation. */
  115.     protected double coef;

  116.     /** coef for SGP C5 computation. */
  117.     protected double coef1;

  118.     /** C1 from SPTRCK #3. */
  119.     protected double c1;

  120.     /** C2 from SPTRCK #3. */
  121.     protected double c2;

  122.     /** C4 from SPTRCK #3. */
  123.     protected double c4;

  124.     /** common parameter for raan (OMEGA) computation. */
  125.     protected double xnodcf;

  126.     /** 3/2 * C1. */
  127.     protected double t2cof;

  128.     // CHECKSTYLE: resume VisibilityModifier check

  129.     /** TLE frame. */
  130.     private final Frame teme;

  131.     /** Spacecraft mass (kg). */
  132.     private final double mass;

  133.     /** Protected constructor for derived classes.
  134.      *
  135.      * <p>This constructor uses the {@link DataContext#getDefault() default data context}.
  136.      *
  137.      * @param initialTLE the unique TLE to propagate
  138.      * @param attitudeProvider provider for attitude computation
  139.      * @param mass spacecraft mass (kg)
  140.      * @see #TLEPropagator(TLE, AttitudeProvider, double, Frame)
  141.      */
  142.     @DefaultDataContext
  143.     protected TLEPropagator(final TLE initialTLE, final AttitudeProvider attitudeProvider,
  144.                             final double mass) {
  145.         this(initialTLE, attitudeProvider, mass,
  146.                 DataContext.getDefault().getFrames().getTEME());
  147.     }

  148.     /** Protected constructor for derived classes.
  149.      * @param initialTLE the unique TLE to propagate
  150.      * @param attitudeProvider provider for attitude computation
  151.      * @param mass spacecraft mass (kg)
  152.      * @param teme the TEME frame to use for propagation.
  153.      * @since 10.1
  154.      */
  155.     protected TLEPropagator(final TLE initialTLE,
  156.                             final AttitudeProvider attitudeProvider,
  157.                             final double mass,
  158.                             final Frame teme) {
  159.         super(attitudeProvider);
  160.         setStartDate(initialTLE.getDate());
  161.         this.tle  = initialTLE;
  162.         this.teme = teme;
  163.         this.mass = mass;
  164.         this.utc = initialTLE.getUtc();
  165.         initializeCommons();
  166.         sxpInitialize();
  167.         // set the initial state
  168.         final Orbit orbit = propagateOrbit(initialTLE.getDate());
  169.         final Attitude attitude = attitudeProvider.getAttitude(orbit, orbit.getDate(), orbit.getFrame());
  170.         super.resetInitialState(new SpacecraftState(orbit, attitude, mass));
  171.     }

  172.     /** Selects the extrapolator to use with the selected TLE.
  173.      *
  174.      * <p>This method uses the {@link DataContext#getDefault() default data context}.
  175.      *
  176.      * @param tle the TLE to propagate.
  177.      * @return the correct propagator.
  178.      * @see #selectExtrapolator(TLE, Frames)
  179.      */
  180.     @DefaultDataContext
  181.     public static TLEPropagator selectExtrapolator(final TLE tle) {
  182.         return selectExtrapolator(tle, DataContext.getDefault().getFrames());
  183.     }

  184.     /** Selects the extrapolator to use with the selected TLE.
  185.      * @param tle the TLE to propagate.
  186.      * @param frames set of Frames to use in the propagator.
  187.      * @return the correct propagator.
  188.      * @since 10.1
  189.      */
  190.     public static TLEPropagator selectExtrapolator(final TLE tle, final Frames frames) {
  191.         return selectExtrapolator(
  192.                 tle,
  193.                 Propagator.getDefaultLaw(frames),
  194.                 DEFAULT_MASS,
  195.                 frames.getTEME());
  196.     }

  197.     /** Selects the extrapolator to use with the selected TLE.
  198.      *
  199.      * <p>This method uses the {@link DataContext#getDefault() default data context}.
  200.      *
  201.      * @param tle the TLE to propagate.
  202.      * @param attitudeProvider provider for attitude computation
  203.      * @param mass spacecraft mass (kg)
  204.      * @return the correct propagator.
  205.      * @see #selectExtrapolator(TLE, AttitudeProvider, double, Frame)
  206.      */
  207.     @DefaultDataContext
  208.     public static TLEPropagator selectExtrapolator(final TLE tle, final AttitudeProvider attitudeProvider,
  209.                                                    final double mass) {
  210.         return selectExtrapolator(tle, attitudeProvider, mass,
  211.                 DataContext.getDefault().getFrames().getTEME());
  212.     }

  213.     /** Selects the extrapolator to use with the selected TLE.
  214.      * @param tle the TLE to propagate.
  215.      * @param attitudeProvider provider for attitude computation
  216.      * @param mass spacecraft mass (kg)
  217.      * @param teme the TEME frame to use for propagation.
  218.      * @return the correct propagator.
  219.      * @since 10.1
  220.      */
  221.     public static TLEPropagator selectExtrapolator(final TLE tle,
  222.                                                    final AttitudeProvider attitudeProvider,
  223.                                                    final double mass,
  224.                                                    final Frame teme) {

  225.         final double a1 = FastMath.pow( TLEConstants.XKE / (tle.getMeanMotion() * 60.0), TLEConstants.TWO_THIRD);
  226.         final double cosi0 = FastMath.cos(tle.getI());
  227.         final double temp = TLEConstants.CK2 * 1.5 * (3 * cosi0 * cosi0 - 1.0) *
  228.                             FastMath.pow(1.0 - tle.getE() * tle.getE(), -1.5);
  229.         final double delta1 = temp / (a1 * a1);
  230.         final double a0 = a1 * (1.0 - delta1 * (TLEConstants.ONE_THIRD + delta1 * (delta1 * 134.0 / 81.0 + 1.0)));
  231.         final double delta0 = temp / (a0 * a0);

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

  234.         // Period >= 225 minutes is deep space
  235.         if (MathUtils.TWO_PI / (xn0dp * TLEConstants.MINUTES_PER_DAY) >= (1.0 / 6.4)) {
  236.             return new DeepSDP4(tle, attitudeProvider, mass, teme);
  237.         } else {
  238.             return new SGP4(tle, attitudeProvider, mass, teme);
  239.         }
  240.     }

  241.     /** Get the Earth gravity coefficient used for TLE propagation.
  242.      * @return the Earth gravity coefficient.
  243.      */
  244.     public static double getMU() {
  245.         return TLEConstants.MU;
  246.     }

  247.     /** Get the extrapolated position and velocity from an initial TLE.
  248.      * @param date the final date
  249.      * @return the final PVCoordinates
  250.      */
  251.     public PVCoordinates getPVCoordinates(final AbsoluteDate date) {

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

  253.         // Compute PV with previous calculated parameters
  254.         return computePVCoordinates();
  255.     }

  256.     /** Computation of the first commons parameters.
  257.      */
  258.     private void initializeCommons() {

  259.         final double a1 = FastMath.pow(TLEConstants.XKE / (tle.getMeanMotion() * 60.0), TLEConstants.TWO_THIRD);
  260.         cosi0 = FastMath.cos(tle.getI());
  261.         theta2 = cosi0 * cosi0;
  262.         final double x3thm1 = 3.0 * theta2 - 1.0;
  263.         e0sq = tle.getE() * tle.getE();
  264.         beta02 = 1.0 - e0sq;
  265.         beta0 = FastMath.sqrt(beta02);
  266.         final double tval = TLEConstants.CK2 * 1.5 * x3thm1 / (beta0 * beta02);
  267.         final double delta1 = tval / (a1 * a1);
  268.         final double a0 = a1 * (1.0 - delta1 * (TLEConstants.ONE_THIRD + delta1 * (1.0 + 134.0 / 81.0 * delta1)));
  269.         final double delta0 = tval / (a0 * a0);

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

  273.         // Values of s and qms2t :
  274.         s4 = TLEConstants.S;  // unmodified value for s
  275.         double q0ms24 = TLEConstants.QOMS2T; // unmodified value for q0ms2T

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

  277.         //  For perigee below 156 km, the values of s and qoms2t are changed :
  278.         if (perige < 156.0) {
  279.             if (perige <= 98.0) {
  280.                 s4 = 20.0;
  281.             } else {
  282.                 s4 = perige - 78.0;
  283.             }
  284.             final double temp_val = (120.0 - s4) * TLEConstants.NORMALIZED_EQUATORIAL_RADIUS / TLEConstants.EARTH_RADIUS;
  285.             final double temp_val_squared = temp_val * temp_val;
  286.             q0ms24 = temp_val_squared * temp_val_squared;
  287.             s4 = s4 / TLEConstants.EARTH_RADIUS + TLEConstants.NORMALIZED_EQUATORIAL_RADIUS; // new value for q0ms2T and s
  288.         }

  289.         final double pinv = 1.0 / (a0dp * beta02);
  290.         final double pinvsq = pinv * pinv;
  291.         tsi = 1.0 / (a0dp - s4);
  292.         eta = a0dp * tle.getE() * tsi;
  293.         etasq = eta * eta;
  294.         eeta = tle.getE() * eta;

  295.         final double psisq = FastMath.abs(1.0 - etasq); // abs because pow 3.5 needs positive value
  296.         final double tsi_squared = tsi * tsi;
  297.         coef = q0ms24 * tsi_squared * tsi_squared;
  298.         coef1 = coef / FastMath.pow(psisq, 3.5);

  299.         // C2 and C1 coefficients computation :
  300.         c2 = coef1 * xn0dp * (a0dp * (1.0 + 1.5 * etasq + eeta * (4.0 + etasq)) +
  301.              0.75 * TLEConstants.CK2 * tsi / psisq * x3thm1 * (8.0 + 3.0 * etasq * (8.0 + etasq)));
  302.         c1 = tle.getBStar() * c2;
  303.         sini0 = FastMath.sin(tle.getI());

  304.         final double x1mth2 = 1.0 - theta2;

  305.         // C4 coefficient computation :
  306.         c4 = 2.0 * xn0dp * coef1 * a0dp * beta02 * (eta * (2.0 + 0.5 * etasq) +
  307.              tle.getE() * (0.5 + 2.0 * etasq) -
  308.              2 * TLEConstants.CK2 * tsi / (a0dp * psisq) *
  309.              (-3.0 * x3thm1 * (1.0 - 2.0 * eeta + etasq * (1.5 - 0.5 * eeta)) +
  310.               0.75 * x1mth2 * (2.0 * etasq - eeta * (1.0 + etasq)) * FastMath.cos(2.0 * tle.getPerigeeArgument())));

  311.         final double theta4 = theta2 * theta2;
  312.         final double temp1 = 3 * TLEConstants.CK2 * pinvsq * xn0dp;
  313.         final double temp2 = temp1 * TLEConstants.CK2 * pinvsq;
  314.         final double temp3 = 1.25 * TLEConstants.CK4 * pinvsq * pinvsq * xn0dp;

  315.         // atmospheric and gravitation coefs :(Mdf and OMEGAdf)
  316.         xmdot = xn0dp +
  317.                 0.5 * temp1 * beta0 * x3thm1 +
  318.                 0.0625 * temp2 * beta0 * (13.0 - 78.0 * theta2 + 137.0 * theta4);

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

  320.         omgdot = -0.5 * temp1 * x1m5th +
  321.                  0.0625 * temp2 * (7.0 - 114.0 * theta2 + 395.0 * theta4) +
  322.                  temp3 * (3.0 - 36.0 * theta2 + 49.0 * theta4);

  323.         final double xhdot1 = -temp1 * cosi0;

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

  327.     }

  328.     /** Retrieves the position and velocity.
  329.      * @return the computed PVCoordinates.
  330.      */
  331.     private PVCoordinates computePVCoordinates() {

  332.         // Long period periodics
  333.         final double axn = e * FastMath.cos(omega);
  334.         double temp = 1.0 / (a * (1.0 - e * e));
  335.         final double xlcof = 0.125 * TLEConstants.A3OVK2 * sini0 * (3.0 + 5.0 * cosi0) / (1.0 + cosi0);
  336.         final double aycof = 0.25 * TLEConstants.A3OVK2 * sini0;
  337.         final double xll = temp * xlcof * axn;
  338.         final double aynl = temp * aycof;
  339.         final double xlt = xl + xll;
  340.         final double ayn = e * FastMath.sin(omega) + aynl;
  341.         final double elsq = axn * axn + ayn * ayn;
  342.         final double capu = MathUtils.normalizeAngle(xlt - xnode, FastMath.PI);
  343.         double epw = capu;
  344.         double ecosE = 0;
  345.         double esinE = 0;
  346.         double sinEPW = 0;
  347.         double cosEPW = 0;

  348.         // Dundee changes:  items dependent on cosio get recomputed:
  349.         final double cosi0Sq = cosi0 * cosi0;
  350.         final double x3thm1 = 3.0 * cosi0Sq - 1.0;
  351.         final double x1mth2 = 1.0 - cosi0Sq;
  352.         final double x7thm1 = 7.0 * cosi0Sq - 1.0;

  353.         if (e > (1 - 1e-6)) {
  354.             throw new OrekitException(OrekitMessages.TOO_LARGE_ECCENTRICITY_FOR_PROPAGATION_MODEL, e);
  355.         }

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

  359.             boolean doSecondOrderNewtonRaphson = true;

  360.             sinEPW = FastMath.sin( epw);
  361.             cosEPW = FastMath.cos( epw);
  362.             ecosE = axn * cosEPW + ayn * sinEPW;
  363.             esinE = axn * sinEPW - ayn * cosEPW;
  364.             final double f = capu - epw + esinE;
  365.             if (FastMath.abs(f) < newtonRaphsonEpsilon) {
  366.                 break;
  367.             }
  368.             final double fdot = 1.0 - ecosE;
  369.             double delta_epw = f / fdot;
  370.             if (j == 0) {
  371.                 final double maxNewtonRaphson = 1.25 * FastMath.abs(e);
  372.                 doSecondOrderNewtonRaphson = false;
  373.                 if (delta_epw > maxNewtonRaphson) {
  374.                     delta_epw = maxNewtonRaphson;
  375.                 } else if (delta_epw < -maxNewtonRaphson) {
  376.                     delta_epw = -maxNewtonRaphson;
  377.                 } else {
  378.                     doSecondOrderNewtonRaphson = true;
  379.                 }
  380.             }
  381.             if (doSecondOrderNewtonRaphson) {
  382.                 delta_epw = f / (fdot + 0.5 * esinE * delta_epw);
  383.             }
  384.             epw += delta_epw;
  385.         }

  386.         // Short period preliminary quantities
  387.         temp = 1.0 - elsq;
  388.         final double pl = a * temp;
  389.         final double r = a * (1.0 - ecosE);
  390.         double temp2 = a / r;
  391.         final double betal = FastMath.sqrt(temp);
  392.         temp = esinE / (1.0 + betal);
  393.         final double cosu = temp2 * (cosEPW - axn + ayn * temp);
  394.         final double sinu = temp2 * (sinEPW - ayn - axn * temp);
  395.         final double u = FastMath.atan2(sinu, cosu);
  396.         final double sin2u = 2.0 * sinu * cosu;
  397.         final double cos2u = 2.0 * cosu * cosu - 1.0;
  398.         final double temp1 = TLEConstants.CK2 / pl;
  399.         temp2 = temp1 / pl;

  400.         // Update for short periodics
  401.         final double rk = r * (1.0 - 1.5 * temp2 * betal * x3thm1) + 0.5 * temp1 * x1mth2 * cos2u;
  402.         final double uk = u - 0.25 * temp2 * x7thm1 * sin2u;
  403.         final double xnodek = xnode + 1.5 * temp2 * cosi0 * sin2u;
  404.         final double xinck = i + 1.5 * temp2 * cosi0 * sini0 * cos2u;

  405.         // Orientation vectors
  406.         final double sinuk = FastMath.sin(uk);
  407.         final double cosuk = FastMath.cos(uk);
  408.         final double sinik = FastMath.sin(xinck);
  409.         final double cosik = FastMath.cos(xinck);
  410.         final double sinnok = FastMath.sin(xnodek);
  411.         final double cosnok = FastMath.cos(xnodek);
  412.         final double xmx = -sinnok * cosik;
  413.         final double xmy = cosnok * cosik;
  414.         final double ux = xmx * sinuk + cosnok * cosuk;
  415.         final double uy = xmy * sinuk + sinnok * cosuk;
  416.         final double uz = sinik * sinuk;

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

  420.         final double rdot   = TLEConstants.XKE * FastMath.sqrt(a) * esinE / r;
  421.         final double rfdot  = TLEConstants.XKE * FastMath.sqrt(pl) / r;
  422.         final double xn     = TLEConstants.XKE / (a * FastMath.sqrt(a));
  423.         final double rdotk  = rdot - xn * temp1 * x1mth2 * sin2u;
  424.         final double rfdotk = rfdot + xn * temp1 * (x1mth2 * cos2u + 1.5 * x3thm1);
  425.         final double vx     = xmx * cosuk - cosnok * sinuk;
  426.         final double vy     = xmy * cosuk - sinnok * sinuk;
  427.         final double vz     = sinik * cosuk;

  428.         final double cv = 1000.0 * TLEConstants.EARTH_RADIUS / 60.0;
  429.         final Vector3D vel = new Vector3D(cv * (rdotk * ux + rfdotk * vx),
  430.                                           cv * (rdotk * uy + rfdotk * vy),
  431.                                           cv * (rdotk * uz + rfdotk * vz));

  432.         return new PVCoordinates(pos, vel);

  433.     }

  434.     /** Initialization proper to each propagator (SGP or SDP).
  435.      */
  436.     protected abstract void sxpInitialize();

  437.     /** Propagation proper to each propagator (SGP or SDP).
  438.      * @param t the offset from initial epoch (min)
  439.      */
  440.     protected abstract void sxpPropagate(double t);

  441.     /** {@inheritDoc} */
  442.     public void resetInitialState(final SpacecraftState state) {
  443.         throw new OrekitException(OrekitMessages.NON_RESETABLE_STATE);
  444.     }

  445.     /** {@inheritDoc} */
  446.     protected void resetIntermediateState(final SpacecraftState state, final boolean forward) {
  447.         throw new OrekitException(OrekitMessages.NON_RESETABLE_STATE);
  448.     }

  449.     /** {@inheritDoc} */
  450.     protected double getMass(final AbsoluteDate date) {
  451.         return mass;
  452.     }

  453.     /** {@inheritDoc} */
  454.     protected Orbit propagateOrbit(final AbsoluteDate date) {
  455.         return new CartesianOrbit(getPVCoordinates(date), teme, date, TLEConstants.MU);
  456.     }

  457.     /** Get the underlying TLE.
  458.      * @return underlying TLE
  459.      */
  460.     public TLE getTLE() {
  461.         return tle;
  462.     }

  463.     /** {@inheritDoc} */
  464.     public Frame getFrame() {
  465.         return teme;
  466.     }

  467. }