SGP4.java

  1. /* Copyright 2002-2013 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.apache.commons.math3.util.FastMath;
  19. import org.orekit.attitudes.AttitudeProvider;
  20. import org.orekit.errors.OrekitException;

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

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

  36.     /** (1 + eta * cos(M0))<sup>3</sup>. */
  37.     private double delM0;

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

  50.     /** Constructor for a unique initial TLE.
  51.      * @param initialTLE the TLE to propagate.
  52.      * @param attitudeProvider provider for attitude computation
  53.      * @param mass spacecraft mass (kg)
  54.      * @exception OrekitException if some specific error occurs
  55.      */
  56.     protected SGP4(final TLE initialTLE, final AttitudeProvider attitudeProvider,
  57.                        final double mass) throws OrekitException {
  58.         super(initialTLE, attitudeProvider, mass);
  59.     }

  60.     /** Initialization proper to each propagator (SGP or SDP).
  61.      */
  62.     protected void sxpInitialize() {

  63.         // For perigee less than 220 kilometers, the equations are truncated to
  64.         // linear variation in sqrt a and quadratic variation in mean anomaly.
  65.         // Also, the c3 term, the delta omega term, and the delta m term are dropped.
  66.         lessThan220 = perige < 220;
  67.         if (!lessThan220) {
  68.             final double c1sq = c1 * c1;
  69.             delM0 = 1.0 + eta * FastMath.cos(tle.getMeanAnomaly());
  70.             delM0 *= delM0 * delM0;
  71.             d2 = 4 * a0dp * tsi * c1sq;
  72.             final double temp = d2 * tsi * c1 / 3.0;
  73.             d3 = (17 * a0dp + s4) * temp;
  74.             d4 = 0.5 * temp * a0dp * tsi * (221 * a0dp + 31 * s4) * c1;
  75.             t3cof = d2 + 2 * c1sq;
  76.             t4cof = 0.25 * (3 * d3 + c1 * (12 * d2 + 10 * c1sq));
  77.             t5cof = 0.2 * (3 * d4 + 12 * c1 * d3 + 6 * d2 * d2 + 15 * c1sq * (2 * d2 + c1sq));
  78.             sinM0 = FastMath.sin(tle.getMeanAnomaly());
  79.             if (tle.getE() < 1e-4) {
  80.                 omgcof = 0.;
  81.                 xmcof = 0.;
  82.             } else  {
  83.                 final double c3 = coef * tsi * TLEConstants.A3OVK2 * xn0dp *
  84.                                   TLEConstants.NORMALIZED_EQUATORIAL_RADIUS * sini0 / tle.getE();
  85.                 xmcof = -TLEConstants.TWO_THIRD * coef * tle.getBStar() *
  86.                         TLEConstants.NORMALIZED_EQUATORIAL_RADIUS / eeta;
  87.                 omgcof = tle.getBStar() * c3 * FastMath.cos(tle.getPerigeeArgument());
  88.             }
  89.         }

  90.         c5 = 2 * coef1 * a0dp * beta02 * (1 + 2.75 * (etasq + eeta) + eeta * etasq);
  91.         // initialized
  92.     }

  93.     /** Propagation proper to each propagator (SGP or SDP).
  94.      * @param tSince the offset from initial epoch (min)
  95.      */
  96.     protected void sxpPropagate(final double tSince) {

  97.         // Update for secular gravity and atmospheric drag.
  98.         final double xmdf = tle.getMeanAnomaly() + xmdot * tSince;
  99.         final double omgadf = tle.getPerigeeArgument() + omgdot * tSince;
  100.         final double xn0ddf = tle.getRaan() + xnodot * tSince;
  101.         omega = omgadf;
  102.         double xmp = xmdf;
  103.         final double tsq = tSince * tSince;
  104.         xnode = xn0ddf + xnodcf * tsq;
  105.         double tempa = 1 - c1 * tSince;
  106.         double tempe = tle.getBStar() * c4 * tSince;
  107.         double templ = t2cof * tsq;

  108.         if (!lessThan220) {
  109.             final double delomg = omgcof * tSince;
  110.             double delm = 1. + eta * FastMath.cos(xmdf);
  111.             delm = xmcof * (delm * delm * delm - delM0);
  112.             final double temp = delomg + delm;
  113.             xmp = xmdf + temp;
  114.             omega = omgadf - temp;
  115.             final double tcube = tsq * tSince;
  116.             final double tfour = tSince * tcube;
  117.             tempa = tempa - d2 * tsq - d3 * tcube - d4 * tfour;
  118.             tempe = tempe + tle.getBStar() * c5 * (FastMath.sin(xmp) - sinM0);
  119.             templ = templ + t3cof * tcube + tfour * (t4cof + tSince * t5cof);
  120.         }

  121.         a = a0dp * tempa * tempa;
  122.         e = tle.getE() - tempe;

  123.         // A highly arbitrary lower limit on e,  of 1e-6:
  124.         if (e < 1e-6) {
  125.             e = 1e-6;
  126.         }

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

  128.         i = tle.getI();

  129.     }

  130. }