SGP4.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 org.hipparchus.util.FastMath;
  19. import org.hipparchus.util.SinCos;
  20. import org.orekit.annotation.DefaultDataContext;
  21. import org.orekit.attitudes.AttitudeProvider;
  22. import org.orekit.data.DataContext;
  23. import org.orekit.frames.Frame;

  24. /** This class contains methods to compute propagated coordinates with the SGP4 model.
  25.  * <p>
  26.  * The user should not bother in this class since it is handled internaly by the
  27.  * {@link TLEPropagator}.
  28.  * </p>
  29.  * <p>This implementation is largely inspired from the paper and source code <a
  30.  * href="https://www.celestrak.com/publications/AIAA/2006-6753/">Revisiting Spacetrack
  31.  * Report #3</a> and is fully compliant with its results and tests cases.</p>
  32.  * @author Felix R. Hoots, Ronald L. Roehrich, December 1980 (original fortran)
  33.  * @author David A. Vallado, Paul Crawford, Richard Hujsak, T.S. Kelso (C++ translation and improvements)
  34.  * @author Fabien Maussion (java translation)
  35.  */
  36. public class SGP4 extends TLEPropagator {

  37.     /** If perige is less than 220 km, some calculus are avoided. */
  38.     private boolean lessThan220;

  39.     /** (1 + eta * cos(M0))³. */
  40.     private double delM0;

  41.     // CHECKSTYLE: stop JavadocVariable check
  42.     private double d2;
  43.     private double d3;
  44.     private double d4;
  45.     private double t3cof;
  46.     private double t4cof;
  47.     private double t5cof;
  48.     private double sinM0;
  49.     private double omgcof;
  50.     private double xmcof;
  51.     private double c5;
  52.     // CHECKSTYLE: resume JavadocVariable check

  53.     /** Constructor for a unique initial TLE.
  54.      *
  55.      * <p>This constructor uses the {@link DataContext#getDefault() default data context}.
  56.      *
  57.      * @param initialTLE the TLE to propagate.
  58.      * @param attitudeProvider provider for attitude computation
  59.      * @param mass spacecraft mass (kg)
  60.      * @see #SGP4(TLE, AttitudeProvider, double, Frame)
  61.      */
  62.     @DefaultDataContext
  63.     public SGP4(final TLE initialTLE, final AttitudeProvider attitudeProvider,
  64.                 final double mass) {
  65.         this(initialTLE, attitudeProvider, mass,
  66.                 DataContext.getDefault().getFrames().getTEME());
  67.     }

  68.     /** Constructor for a unique initial TLE.
  69.      * @param initialTLE the TLE to propagate.
  70.      * @param attitudeProvider provider for attitude computation
  71.      * @param mass spacecraft mass (kg)
  72.      * @param teme the TEME frame to use for propagation.
  73.      * @since 10.1
  74.      */
  75.     public SGP4(final TLE initialTLE,
  76.                 final AttitudeProvider attitudeProvider,
  77.                 final double mass,
  78.                 final Frame teme) {
  79.         super(initialTLE, attitudeProvider, mass, teme);
  80.     }

  81.     /** Initialization proper to each propagator (SGP or SDP).
  82.      */
  83.     protected void sxpInitialize() {

  84.         // For perigee less than 220 kilometers, the equations are truncated to
  85.         // linear variation in sqrt a and quadratic variation in mean anomaly.
  86.         // Also, the c3 term, the delta omega term, and the delta m term are dropped.
  87.         lessThan220 = perige < 220;
  88.         if (!lessThan220) {
  89.             final SinCos scM0 = FastMath.sinCos(tle.getMeanAnomaly());
  90.             final double c1sq = c1 * c1;
  91.             delM0 = 1.0 + eta * scM0.cos();
  92.             delM0 *= delM0 * delM0;
  93.             d2 = 4 * a0dp * tsi * c1sq;
  94.             final double temp = d2 * tsi * c1 / 3.0;
  95.             d3 = (17 * a0dp + s4) * temp;
  96.             d4 = 0.5 * temp * a0dp * tsi * (221 * a0dp + 31 * s4) * c1;
  97.             t3cof = d2 + 2 * c1sq;
  98.             t4cof = 0.25 * (3 * d3 + c1 * (12 * d2 + 10 * c1sq));
  99.             t5cof = 0.2 * (3 * d4 + 12 * c1 * d3 + 6 * d2 * d2 + 15 * c1sq * (2 * d2 + c1sq));
  100.             sinM0 = scM0.sin();
  101.             if (tle.getE() < 1e-4) {
  102.                 omgcof = 0.;
  103.                 xmcof = 0.;
  104.             } else  {
  105.                 final double c3 = coef * tsi * TLEConstants.A3OVK2 * xn0dp *
  106.                                   TLEConstants.NORMALIZED_EQUATORIAL_RADIUS * sini0 / tle.getE();
  107.                 xmcof = -TLEConstants.TWO_THIRD * coef * tle.getBStar() *
  108.                         TLEConstants.NORMALIZED_EQUATORIAL_RADIUS / eeta;
  109.                 omgcof = tle.getBStar() * c3 * FastMath.cos(tle.getPerigeeArgument());
  110.             }
  111.         }

  112.         c5 = 2 * coef1 * a0dp * beta02 * (1 + 2.75 * (etasq + eeta) + eeta * etasq);
  113.         // initialized
  114.     }

  115.     /** Propagation proper to each propagator (SGP or SDP).
  116.      * @param tSince the offset from initial epoch (min)
  117.      */
  118.     protected void sxpPropagate(final double tSince) {

  119.         // Update for secular gravity and atmospheric drag.
  120.         final double xmdf = tle.getMeanAnomaly() + xmdot * tSince;
  121.         final double omgadf = tle.getPerigeeArgument() + omgdot * tSince;
  122.         final double xn0ddf = tle.getRaan() + xnodot * tSince;
  123.         omega = omgadf;
  124.         double xmp = xmdf;
  125.         final double tsq = tSince * tSince;
  126.         xnode = xn0ddf + xnodcf * tsq;
  127.         double tempa = 1 - c1 * tSince;
  128.         double tempe = tle.getBStar(tle.getDate().shiftedBy(tSince)) * c4 * tSince;
  129.         double templ = t2cof * tsq;

  130.         if (!lessThan220) {
  131.             final double delomg = omgcof * tSince;
  132.             double delm = 1. + eta * FastMath.cos(xmdf);
  133.             delm = xmcof * (delm * delm * delm - delM0);
  134.             final double temp = delomg + delm;
  135.             xmp = xmdf + temp;
  136.             omega = omgadf - temp;
  137.             final double tcube = tsq * tSince;
  138.             final double tfour = tSince * tcube;
  139.             tempa = tempa - d2 * tsq - d3 * tcube - d4 * tfour;
  140.             tempe = tempe + tle.getBStar(tle.getDate().shiftedBy(tSince)) * c5 * (FastMath.sin(xmp) - sinM0);
  141.             templ = templ + t3cof * tcube + tfour * (t4cof + tSince * t5cof);
  142.         }

  143.         a = a0dp * tempa * tempa;
  144.         e = tle.getE() - tempe;

  145.         // A highly arbitrary lower limit on e,  of 1e-6:
  146.         if (e < 1e-6) {
  147.             e = 1e-6;
  148.         }

  149.         xl = xmp + omega + xnode + xn0dp * templ;

  150.         i = tle.getI();

  151.     }

  152. }