TLEPropagator.java

  1. /* Copyright 2002-2019 CS Systèmes d'Information
  2.  * Licensed to CS Systèmes d'Information (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.attitudes.Attitude;
  22. import org.orekit.attitudes.AttitudeProvider;
  23. import org.orekit.errors.OrekitException;
  24. import org.orekit.errors.OrekitMessages;
  25. import org.orekit.frames.Frame;
  26. import org.orekit.frames.FramesFactory;
  27. import org.orekit.orbits.CartesianOrbit;
  28. import org.orekit.orbits.Orbit;
  29. import org.orekit.propagation.SpacecraftState;
  30. import org.orekit.propagation.analytical.AbstractAnalyticalPropagator;
  31. import org.orekit.time.AbsoluteDate;
  32. import org.orekit.utils.PVCoordinates;


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

  59.     // CHECKSTYLE: stop VisibilityModifier check

  60.     /** Initial state. */
  61.     protected final TLE tle;

  62.     /** final RAAN. */
  63.     protected double xnode;

  64.     /** final semi major axis. */
  65.     protected double a;

  66.     /** final eccentricity. */
  67.     protected double e;

  68.     /** final inclination. */
  69.     protected double i;

  70.     /** final perigee argument. */
  71.     protected double omega;

  72.     /** L from SPTRCK #3. */
  73.     protected double xl;

  74.     /** original recovered semi major axis. */
  75.     protected double a0dp;

  76.     /** original recovered mean motion. */
  77.     protected double xn0dp;

  78.     /** cosinus original inclination. */
  79.     protected double cosi0;

  80.     /** cos io squared. */
  81.     protected double theta2;

  82.     /** sinus original inclination. */
  83.     protected double sini0;

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

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

  88.     /** common parameter for raan (OMEGA) computation. */
  89.     protected double xnodot;

  90.     /** original eccentricity squared. */
  91.     protected double e0sq;
  92.     /** 1 - e2. */
  93.     protected double beta02;

  94.     /** sqrt (1 - e2). */
  95.     protected double beta0;

  96.     /** perigee, expressed in KM and ALTITUDE. */
  97.     protected double perige;

  98.     /** eta squared. */
  99.     protected double etasq;

  100.     /** original eccentricity * eta. */
  101.     protected double eeta;

  102.     /** s* new value for the contant s. */
  103.     protected double s4;

  104.     /** tsi from SPTRCK #3. */
  105.     protected double tsi;

  106.     /** eta from SPTRCK #3. */
  107.     protected double eta;

  108.     /** coef for SGP C3 computation. */
  109.     protected double coef;

  110.     /** coef for SGP C5 computation. */
  111.     protected double coef1;

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

  114.     /** C2 from SPTRCK #3. */
  115.     protected double c2;

  116.     /** C4 from SPTRCK #3. */
  117.     protected double c4;

  118.     /** common parameter for raan (OMEGA) computation. */
  119.     protected double xnodcf;

  120.     /** 3/2 * C1. */
  121.     protected double t2cof;

  122.     // CHECKSTYLE: resume VisibilityModifier check

  123.     /** TLE frame. */
  124.     private final Frame teme;

  125.     /** Spacecraft mass (kg). */
  126.     private final double mass;

  127.     /** Protected constructor for derived classes.
  128.      * @param initialTLE the unique TLE to propagate
  129.      * @param attitudeProvider provider for attitude computation
  130.      * @param mass spacecraft mass (kg)
  131.      */
  132.     protected TLEPropagator(final TLE initialTLE, final AttitudeProvider attitudeProvider,
  133.                             final double mass) {
  134.         super(attitudeProvider);
  135.         setStartDate(initialTLE.getDate());
  136.         this.tle  = initialTLE;
  137.         this.teme = FramesFactory.getTEME();
  138.         this.mass = mass;
  139.         initializeCommons();
  140.         sxpInitialize();
  141.         // set the initial state
  142.         final Orbit orbit = propagateOrbit(initialTLE.getDate());
  143.         final Attitude attitude = attitudeProvider.getAttitude(orbit, orbit.getDate(), orbit.getFrame());
  144.         super.resetInitialState(new SpacecraftState(orbit, attitude, mass));
  145.     }

  146.     /** Selects the extrapolator to use with the selected TLE.
  147.      * @param tle the TLE to propagate.
  148.      * @return the correct propagator.
  149.      */
  150.     public static TLEPropagator selectExtrapolator(final TLE tle) {
  151.         return selectExtrapolator(tle, DEFAULT_LAW, DEFAULT_MASS);
  152.     }

  153.     /** Selects the extrapolator to use with the selected TLE.
  154.      * @param tle the TLE to propagate.
  155.      * @param attitudeProvider provider for attitude computation
  156.      * @param mass spacecraft mass (kg)
  157.      * @return the correct propagator.
  158.      */
  159.     public static TLEPropagator selectExtrapolator(final TLE tle, final AttitudeProvider attitudeProvider,
  160.                                                    final double mass) {

  161.         final double a1 = FastMath.pow( TLEConstants.XKE / (tle.getMeanMotion() * 60.0), TLEConstants.TWO_THIRD);
  162.         final double cosi0 = FastMath.cos(tle.getI());
  163.         final double temp = TLEConstants.CK2 * 1.5 * (3 * cosi0 * cosi0 - 1.0) *
  164.                             FastMath.pow(1.0 - tle.getE() * tle.getE(), -1.5);
  165.         final double delta1 = temp / (a1 * a1);
  166.         final double a0 = a1 * (1.0 - delta1 * (TLEConstants.ONE_THIRD + delta1 * (delta1 * 134.0 / 81.0 + 1.0)));
  167.         final double delta0 = temp / (a0 * a0);

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

  170.         // Period >= 225 minutes is deep space
  171.         if (MathUtils.TWO_PI / (xn0dp * TLEConstants.MINUTES_PER_DAY) >= (1.0 / 6.4)) {
  172.             return new DeepSDP4(tle, attitudeProvider, mass);
  173.         } else {
  174.             return new SGP4(tle, attitudeProvider, mass);
  175.         }
  176.     }

  177.     /** Get the Earth gravity coefficient used for TLE propagation.
  178.      * @return the Earth gravity coefficient.
  179.      */
  180.     public static double getMU() {
  181.         return TLEConstants.MU;
  182.     }

  183.     /** Get the extrapolated position and velocity from an initial TLE.
  184.      * @param date the final date
  185.      * @return the final PVCoordinates
  186.      */
  187.     public PVCoordinates getPVCoordinates(final AbsoluteDate date) {

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

  189.         // Compute PV with previous calculated parameters
  190.         return computePVCoordinates();
  191.     }

  192.     /** Computation of the first commons parameters.
  193.      */
  194.     private void initializeCommons() {

  195.         final double a1 = FastMath.pow(TLEConstants.XKE / (tle.getMeanMotion() * 60.0), TLEConstants.TWO_THIRD);
  196.         cosi0 = FastMath.cos(tle.getI());
  197.         theta2 = cosi0 * cosi0;
  198.         final double x3thm1 = 3.0 * theta2 - 1.0;
  199.         e0sq = tle.getE() * tle.getE();
  200.         beta02 = 1.0 - e0sq;
  201.         beta0 = FastMath.sqrt(beta02);
  202.         final double tval = TLEConstants.CK2 * 1.5 * x3thm1 / (beta0 * beta02);
  203.         final double delta1 = tval / (a1 * a1);
  204.         final double a0 = a1 * (1.0 - delta1 * (TLEConstants.ONE_THIRD + delta1 * (1.0 + 134.0 / 81.0 * delta1)));
  205.         final double delta0 = tval / (a0 * a0);

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

  209.         // Values of s and qms2t :
  210.         s4 = TLEConstants.S;  // unmodified value for s
  211.         double q0ms24 = TLEConstants.QOMS2T; // unmodified value for q0ms2T

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

  213.         //  For perigee below 156 km, the values of s and qoms2t are changed :
  214.         if (perige < 156.0) {
  215.             if (perige <= 98.0) {
  216.                 s4 = 20.0;
  217.             } else {
  218.                 s4 = perige - 78.0;
  219.             }
  220.             final double temp_val = (120.0 - s4) * TLEConstants.NORMALIZED_EQUATORIAL_RADIUS / TLEConstants.EARTH_RADIUS;
  221.             final double temp_val_squared = temp_val * temp_val;
  222.             q0ms24 = temp_val_squared * temp_val_squared;
  223.             s4 = s4 / TLEConstants.EARTH_RADIUS + TLEConstants.NORMALIZED_EQUATORIAL_RADIUS; // new value for q0ms2T and s
  224.         }

  225.         final double pinv = 1.0 / (a0dp * beta02);
  226.         final double pinvsq = pinv * pinv;
  227.         tsi = 1.0 / (a0dp - s4);
  228.         eta = a0dp * tle.getE() * tsi;
  229.         etasq = eta * eta;
  230.         eeta = tle.getE() * eta;

  231.         final double psisq = FastMath.abs(1.0 - etasq); // abs because pow 3.5 needs positive value
  232.         final double tsi_squared = tsi * tsi;
  233.         coef = q0ms24 * tsi_squared * tsi_squared;
  234.         coef1 = coef / FastMath.pow(psisq, 3.5);

  235.         // C2 and C1 coefficients computation :
  236.         c2 = coef1 * xn0dp * (a0dp * (1.0 + 1.5 * etasq + eeta * (4.0 + etasq)) +
  237.              0.75 * TLEConstants.CK2 * tsi / psisq * x3thm1 * (8.0 + 3.0 * etasq * (8.0 + etasq)));
  238.         c1 = tle.getBStar() * c2;
  239.         sini0 = FastMath.sin(tle.getI());

  240.         final double x1mth2 = 1.0 - theta2;

  241.         // C4 coefficient computation :
  242.         c4 = 2.0 * xn0dp * coef1 * a0dp * beta02 * (eta * (2.0 + 0.5 * etasq) +
  243.              tle.getE() * (0.5 + 2.0 * etasq) -
  244.              2 * TLEConstants.CK2 * tsi / (a0dp * psisq) *
  245.              (-3.0 * x3thm1 * (1.0 - 2.0 * eeta + etasq * (1.5 - 0.5 * eeta)) +
  246.               0.75 * x1mth2 * (2.0 * etasq - eeta * (1.0 + etasq)) * FastMath.cos(2.0 * tle.getPerigeeArgument())));

  247.         final double theta4 = theta2 * theta2;
  248.         final double temp1 = 3 * TLEConstants.CK2 * pinvsq * xn0dp;
  249.         final double temp2 = temp1 * TLEConstants.CK2 * pinvsq;
  250.         final double temp3 = 1.25 * TLEConstants.CK4 * pinvsq * pinvsq * xn0dp;

  251.         // atmospheric and gravitation coefs :(Mdf and OMEGAdf)
  252.         xmdot = xn0dp +
  253.                 0.5 * temp1 * beta0 * x3thm1 +
  254.                 0.0625 * temp2 * beta0 * (13.0 - 78.0 * theta2 + 137.0 * theta4);

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

  256.         omgdot = -0.5 * temp1 * x1m5th +
  257.                  0.0625 * temp2 * (7.0 - 114.0 * theta2 + 395.0 * theta4) +
  258.                  temp3 * (3.0 - 36.0 * theta2 + 49.0 * theta4);

  259.         final double xhdot1 = -temp1 * cosi0;

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

  263.     }

  264.     /** Retrieves the position and velocity.
  265.      * @return the computed PVCoordinates.
  266.      */
  267.     private PVCoordinates computePVCoordinates() {

  268.         // Long period periodics
  269.         final double axn = e * FastMath.cos(omega);
  270.         double temp = 1.0 / (a * (1.0 - e * e));
  271.         final double xlcof = 0.125 * TLEConstants.A3OVK2 * sini0 * (3.0 + 5.0 * cosi0) / (1.0 + cosi0);
  272.         final double aycof = 0.25 * TLEConstants.A3OVK2 * sini0;
  273.         final double xll = temp * xlcof * axn;
  274.         final double aynl = temp * aycof;
  275.         final double xlt = xl + xll;
  276.         final double ayn = e * FastMath.sin(omega) + aynl;
  277.         final double elsq = axn * axn + ayn * ayn;
  278.         final double capu = MathUtils.normalizeAngle(xlt - xnode, FastMath.PI);
  279.         double epw = capu;
  280.         double ecosE = 0;
  281.         double esinE = 0;
  282.         double sinEPW = 0;
  283.         double cosEPW = 0;

  284.         // Dundee changes:  items dependent on cosio get recomputed:
  285.         final double cosi0Sq = cosi0 * cosi0;
  286.         final double x3thm1 = 3.0 * cosi0Sq - 1.0;
  287.         final double x1mth2 = 1.0 - cosi0Sq;
  288.         final double x7thm1 = 7.0 * cosi0Sq - 1.0;

  289.         if (e > (1 - 1e-6)) {
  290.             throw new OrekitException(OrekitMessages.TOO_LARGE_ECCENTRICITY_FOR_PROPAGATION_MODEL, e);
  291.         }

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

  295.             boolean doSecondOrderNewtonRaphson = true;

  296.             sinEPW = FastMath.sin( epw);
  297.             cosEPW = FastMath.cos( epw);
  298.             ecosE = axn * cosEPW + ayn * sinEPW;
  299.             esinE = axn * sinEPW - ayn * cosEPW;
  300.             final double f = capu - epw + esinE;
  301.             if (FastMath.abs(f) < newtonRaphsonEpsilon) {
  302.                 break;
  303.             }
  304.             final double fdot = 1.0 - ecosE;
  305.             double delta_epw = f / fdot;
  306.             if (j == 0) {
  307.                 final double maxNewtonRaphson = 1.25 * FastMath.abs(e);
  308.                 doSecondOrderNewtonRaphson = false;
  309.                 if (delta_epw > maxNewtonRaphson) {
  310.                     delta_epw = maxNewtonRaphson;
  311.                 } else if (delta_epw < -maxNewtonRaphson) {
  312.                     delta_epw = -maxNewtonRaphson;
  313.                 } else {
  314.                     doSecondOrderNewtonRaphson = true;
  315.                 }
  316.             }
  317.             if (doSecondOrderNewtonRaphson) {
  318.                 delta_epw = f / (fdot + 0.5 * esinE * delta_epw);
  319.             }
  320.             epw += delta_epw;
  321.         }

  322.         // Short period preliminary quantities
  323.         temp = 1.0 - elsq;
  324.         final double pl = a * temp;
  325.         final double r = a * (1.0 - ecosE);
  326.         double temp2 = a / r;
  327.         final double betal = FastMath.sqrt(temp);
  328.         temp = esinE / (1.0 + betal);
  329.         final double cosu = temp2 * (cosEPW - axn + ayn * temp);
  330.         final double sinu = temp2 * (sinEPW - ayn - axn * temp);
  331.         final double u = FastMath.atan2(sinu, cosu);
  332.         final double sin2u = 2.0 * sinu * cosu;
  333.         final double cos2u = 2.0 * cosu * cosu - 1.0;
  334.         final double temp1 = TLEConstants.CK2 / pl;
  335.         temp2 = temp1 / pl;

  336.         // Update for short periodics
  337.         final double rk = r * (1.0 - 1.5 * temp2 * betal * x3thm1) + 0.5 * temp1 * x1mth2 * cos2u;
  338.         final double uk = u - 0.25 * temp2 * x7thm1 * sin2u;
  339.         final double xnodek = xnode + 1.5 * temp2 * cosi0 * sin2u;
  340.         final double xinck = i + 1.5 * temp2 * cosi0 * sini0 * cos2u;

  341.         // Orientation vectors
  342.         final double sinuk = FastMath.sin(uk);
  343.         final double cosuk = FastMath.cos(uk);
  344.         final double sinik = FastMath.sin(xinck);
  345.         final double cosik = FastMath.cos(xinck);
  346.         final double sinnok = FastMath.sin(xnodek);
  347.         final double cosnok = FastMath.cos(xnodek);
  348.         final double xmx = -sinnok * cosik;
  349.         final double xmy = cosnok * cosik;
  350.         final double ux = xmx * sinuk + cosnok * cosuk;
  351.         final double uy = xmy * sinuk + sinnok * cosuk;
  352.         final double uz = sinik * sinuk;

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

  356.         final double rdot   = TLEConstants.XKE * FastMath.sqrt(a) * esinE / r;
  357.         final double rfdot  = TLEConstants.XKE * FastMath.sqrt(pl) / r;
  358.         final double xn     = TLEConstants.XKE / (a * FastMath.sqrt(a));
  359.         final double rdotk  = rdot - xn * temp1 * x1mth2 * sin2u;
  360.         final double rfdotk = rfdot + xn * temp1 * (x1mth2 * cos2u + 1.5 * x3thm1);
  361.         final double vx     = xmx * cosuk - cosnok * sinuk;
  362.         final double vy     = xmy * cosuk - sinnok * sinuk;
  363.         final double vz     = sinik * cosuk;

  364.         final double cv = 1000.0 * TLEConstants.EARTH_RADIUS / 60.0;
  365.         final Vector3D vel = new Vector3D(cv * (rdotk * ux + rfdotk * vx),
  366.                                           cv * (rdotk * uy + rfdotk * vy),
  367.                                           cv * (rdotk * uz + rfdotk * vz));

  368.         return new PVCoordinates(pos, vel);

  369.     }

  370.     /** Initialization proper to each propagator (SGP or SDP).
  371.      */
  372.     protected abstract void sxpInitialize();

  373.     /** Propagation proper to each propagator (SGP or SDP).
  374.      * @param t the offset from initial epoch (min)
  375.      */
  376.     protected abstract void sxpPropagate(double t);

  377.     /** {@inheritDoc} */
  378.     public void resetInitialState(final SpacecraftState state) {
  379.         throw new OrekitException(OrekitMessages.NON_RESETABLE_STATE);
  380.     }

  381.     /** {@inheritDoc} */
  382.     protected void resetIntermediateState(final SpacecraftState state, final boolean forward) {
  383.         throw new OrekitException(OrekitMessages.NON_RESETABLE_STATE);
  384.     }

  385.     /** {@inheritDoc} */
  386.     protected double getMass(final AbsoluteDate date) {
  387.         return mass;
  388.     }

  389.     /** {@inheritDoc} */
  390.     protected Orbit propagateOrbit(final AbsoluteDate date) {
  391.         return new CartesianOrbit(getPVCoordinates(date), teme, date, TLEConstants.MU);
  392.     }

  393.     /** Get the underlying TLE.
  394.      * @return underlying TLE
  395.      */
  396.     public TLE getTLE() {
  397.         return tle;
  398.     }

  399.     /** {@inheritDoc} */
  400.     public Frame getFrame() {
  401.         return teme;
  402.     }

  403. }