DeepSDP4.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.util.FastMath;
  19. import org.hipparchus.util.MathUtils;
  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. import org.orekit.time.DateTimeComponents;
  25. import org.orekit.utils.Constants;


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

  39.     // CHECKSTYLE: stop JavadocVariable check

  40.     // Internal constants
  41.     private static final double ZNS      = 1.19459E-5;
  42.     private static final double ZES      = 0.01675;
  43.     private static final double ZNL      = 1.5835218E-4;
  44.     private static final double ZEL      = 0.05490;
  45.     private static final double THDT     = 4.3752691E-3;
  46.     private static final double C1SS     =  2.9864797E-6;
  47.     private static final double C1L      = 4.7968065E-7;

  48.     private static final double ROOT22   = 1.7891679E-6;
  49.     private static final double ROOT32   = 3.7393792E-7;
  50.     private static final double ROOT44   = 7.3636953E-9;
  51.     private static final double ROOT52   = 1.1428639E-7;
  52.     private static final double ROOT54   = 2.1765803E-9;

  53.     private static final double Q22      =  1.7891679E-6;
  54.     private static final double Q31      =  2.1460748E-6;
  55.     private static final double Q33      =  2.2123015E-7;

  56.     private static final double C_FASX2  =  0.99139134268488593;
  57.     private static final double S_FASX2  =  0.13093206501640101;
  58.     private static final double C_2FASX4 =  0.87051638752972937;
  59.     private static final double S_2FASX4 = -0.49213943048915526;
  60.     private static final double C_3FASX6 =  0.43258117585763334;
  61.     private static final double S_3FASX6 =  0.90159499016666422;

  62.     private static final double C_G22    =  0.87051638752972937;
  63.     private static final double S_G22    = -0.49213943048915526;
  64.     private static final double C_G32    =  0.57972190187001149;
  65.     private static final double S_G32    =  0.81481440616389245;
  66.     private static final double C_G44    = -0.22866241528815548;
  67.     private static final double S_G44    =  0.97350577801807991;
  68.     private static final double C_G52    =  0.49684831179884198;
  69.     private static final double S_G52    =  0.86783740128127729;
  70.     private static final double C_G54    = -0.29695209575316894;
  71.     private static final double S_G54    = -0.95489237761529999;

  72.     /** Integration step (seconds). */
  73.     private static final double SECULAR_INTEGRATION_STEP  = 720.0;

  74.     /** Intermediate values. */
  75.     private double thgr;
  76.     private double xnq;
  77.     private double omegaq;
  78.     private double zcosil;
  79.     private double zsinil;
  80.     private double zsinhl;
  81.     private double zcoshl;
  82.     private double zmol;
  83.     private double zcosgl;
  84.     private double zsingl;
  85.     private double zmos;
  86.     private double savtsn;

  87.     private double ee2;
  88.     private double e3;
  89.     private double xi2;
  90.     private double xi3;
  91.     private double xl2;
  92.     private double xl3;
  93.     private double xl4;
  94.     private double xgh2;
  95.     private double xgh3;
  96.     private double xgh4;
  97.     private double xh2;
  98.     private double xh3;

  99.     private double d2201;
  100.     private double d2211;
  101.     private double d3210;
  102.     private double d3222;
  103.     private double d4410;
  104.     private double d4422;
  105.     private double d5220;
  106.     private double d5232;
  107.     private double d5421;
  108.     private double d5433;
  109.     private double xlamo;

  110.     private double sse;
  111.     private double ssi;
  112.     private double ssl;
  113.     private double ssh;
  114.     private double ssg;
  115.     private double se2;
  116.     private double si2;
  117.     private double sl2;
  118.     private double sgh2;
  119.     private double sh2;
  120.     private double se3;
  121.     private double si3;
  122.     private double sl3;
  123.     private double sgh3;
  124.     private double sh3;
  125.     private double sl4;
  126.     private double sgh4;

  127.     private double del1;
  128.     private double del2;
  129.     private double del3;
  130.     private double xfact;
  131.     private double xli;
  132.     private double xni;
  133.     private double atime;

  134.     private double pe;
  135.     private double pinc;
  136.     private double pl;
  137.     private double pgh;
  138.     private double ph;

  139.     private double[] derivs;

  140.     // CHECKSTYLE: resume JavadocVariable check

  141.     /** Flag for resonant orbits. */
  142.     private boolean resonant;

  143.     /** Flag for synchronous orbits. */
  144.     private boolean synchronous;

  145.     /** Flag for compliance with Dundee modifications. */
  146.     private boolean isDundeeCompliant = true;

  147.     /** Constructor for a unique initial TLE.
  148.      *
  149.      * <p>This constructor uses the {@link DataContext#getDefault() default data context}.
  150.      *
  151.      * @param initialTLE the TLE to propagate.
  152.      * @param attitudeProvider provider for attitude computation
  153.      * @param mass spacecraft mass (kg)
  154.      * @see #DeepSDP4(TLE, AttitudeProvider, double, Frame)
  155.      */
  156.     @DefaultDataContext
  157.     public DeepSDP4(final TLE initialTLE, final AttitudeProvider attitudeProvider,
  158.                     final double mass) {
  159.         this(initialTLE, attitudeProvider, mass,
  160.                 DataContext.getDefault().getFrames().getTEME());
  161.     }

  162.     /** Constructor for a unique initial TLE.
  163.      * @param initialTLE the TLE to propagate.
  164.      * @param attitudeProvider provider for attitude computation
  165.      * @param mass spacecraft mass (kg)
  166.      * @param teme the TEME frame to use for propagation.
  167.      * @since 10.1
  168.      */
  169.     public DeepSDP4(final TLE initialTLE,
  170.                     final AttitudeProvider attitudeProvider,
  171.                     final double mass,
  172.                     final Frame teme) {
  173.         super(initialTLE, attitudeProvider, mass, teme);
  174.     }

  175.     /** Computes luni - solar terms from initial coordinates and epoch.
  176.      */
  177.     protected void luniSolarTermsComputation() {

  178.         final double sing = FastMath.sin(tle.getPerigeeArgument());
  179.         final double cosg = FastMath.cos(tle.getPerigeeArgument());

  180.         final double sinq = FastMath.sin(tle.getRaan());
  181.         final double cosq = FastMath.cos(tle.getRaan());
  182.         final double aqnv = 1.0 / a0dp;

  183.         // Compute julian days since 1900
  184.         final double daysSince1900 = (tle.getDate()
  185.                 .getComponents(utc)
  186.                 .offsetFrom(DateTimeComponents.JULIAN_EPOCH)) /
  187.                 Constants.JULIAN_DAY - 2415020;

  188.         double cc = C1SS;
  189.         double ze = ZES;
  190.         double zn = ZNS;
  191.         double zsinh = sinq;
  192.         double zcosh = cosq;

  193.         thgr = thetaG(tle.getDate());
  194.         xnq = xn0dp;
  195.         omegaq = tle.getPerigeeArgument();

  196.         final double xnodce = 4.5236020 - 9.2422029e-4 * daysSince1900;
  197.         final double stem = FastMath.sin(xnodce);
  198.         final double ctem = FastMath.cos(xnodce);
  199.         final double c_minus_gam = 0.228027132 * daysSince1900 - 1.1151842;
  200.         final double gam = 5.8351514 + 0.0019443680 * daysSince1900;

  201.         zcosil = 0.91375164 - 0.03568096 * ctem;
  202.         zsinil = FastMath.sqrt(1.0 - zcosil * zcosil);
  203.         zsinhl = 0.089683511 * stem / zsinil;
  204.         zcoshl = FastMath.sqrt(1.0 - zsinhl * zsinhl);
  205.         zmol = MathUtils.normalizeAngle(c_minus_gam, FastMath.PI);

  206.         double zx = 0.39785416 * stem / zsinil;
  207.         final double zy = zcoshl * ctem + 0.91744867 * zsinhl * stem;
  208.         zx = FastMath.atan2( zx, zy) + gam - xnodce;
  209.         zcosgl = FastMath.cos( zx);
  210.         zsingl = FastMath.sin( zx);
  211.         zmos = MathUtils.normalizeAngle(6.2565837 + 0.017201977 * daysSince1900, FastMath.PI);

  212.         // Do solar terms
  213.         savtsn = 1e20;

  214.         double zcosi =  0.91744867;
  215.         double zsini =  0.39785416;
  216.         double zsing = -0.98088458;
  217.         double zcosg =  0.1945905;

  218.         double se = 0;
  219.         double sgh = 0;
  220.         double sh = 0;
  221.         double si = 0;
  222.         double sl = 0;

  223.         // There was previously some convoluted logic here, but it boils
  224.         // down to this:  we compute the solar terms,  then the lunar terms.
  225.         // On a second pass,  we recompute the solar terms, taking advantage
  226.         // of the improved data that resulted from computing lunar terms.
  227.         for (int iteration = 0; iteration < 2; ++iteration) {
  228.             final double a1 = zcosg * zcosh + zsing * zcosi * zsinh;
  229.             final double a3 = -zsing * zcosh + zcosg * zcosi * zsinh;
  230.             final double a7 = -zcosg * zsinh + zsing * zcosi * zcosh;
  231.             final double a8 = zsing * zsini;
  232.             final double a9 = zsing * zsinh + zcosg * zcosi * zcosh;
  233.             final double a10 = zcosg * zsini;
  234.             final double a2 = cosi0 * a7 + sini0 * a8;
  235.             final double a4 = cosi0 * a9 + sini0 * a10;
  236.             final double a5 = -sini0 * a7 + cosi0 * a8;
  237.             final double a6 = -sini0 * a9 + cosi0 * a10;
  238.             final double x1 = a1 * cosg + a2 * sing;
  239.             final double x2 = a3 * cosg + a4 * sing;
  240.             final double x3 = -a1 * sing + a2 * cosg;
  241.             final double x4 = -a3 * sing + a4 * cosg;
  242.             final double x5 = a5 * sing;
  243.             final double x6 = a6 * sing;
  244.             final double x7 = a5 * cosg;
  245.             final double x8 = a6 * cosg;
  246.             final double z31 = 12 * x1 * x1 - 3 * x3 * x3;
  247.             final double z32 = 24 * x1 * x2 - 6 * x3 * x4;
  248.             final double z33 = 12 * x2 * x2 - 3 * x4 * x4;
  249.             final double z11 = -6 * a1 * a5 + e0sq * (-24 * x1 * x7 - 6 * x3 * x5);
  250.             final double z12 = -6 * (a1 * a6 + a3 * a5) +
  251.                                e0sq * (-24 * (x2 * x7 + x1 * x8) - 6 * (x3 * x6 + x4 * x5));
  252.             final double z13 = -6 * a3 * a6 + e0sq * (-24 * x2 * x8 - 6 * x4 * x6);
  253.             final double z21 = 6 * a2 * a5 + e0sq * (24 * x1 * x5 - 6 * x3 * x7);
  254.             final double z22 = 6 * (a4 * a5 + a2 * a6) +
  255.                                e0sq * (24 * (x2 * x5 + x1 * x6) - 6 * (x4 * x7 + x3 * x8));
  256.             final double z23 = 6 * a4 * a6 + e0sq * (24 * x2 * x6 - 6 * x4 * x8);
  257.             final double s3 = cc / xnq;
  258.             final double s2 = -0.5 * s3 / beta0;
  259.             final double s4 = s3 * beta0;
  260.             final double s1 = -15 * tle.getE() * s4;
  261.             final double s5 = x1 * x3 + x2 * x4;
  262.             final double s6 = x2 * x3 + x1 * x4;
  263.             final double s7 = x2 * x4 - x1 * x3;
  264.             double z1 = 3 * (a1 * a1 + a2 * a2) + z31 * e0sq;
  265.             double z2 = 6 * (a1 * a3 + a2 * a4) + z32 * e0sq;
  266.             double z3 = 3 * (a3 * a3 + a4 * a4) + z33 * e0sq;

  267.             z1 = z1 + z1 + beta02 * z31;
  268.             z2 = z2 + z2 + beta02 * z32;
  269.             z3 = z3 + z3 + beta02 * z33;
  270.             se = s1 * zn * s5;
  271.             si = s2 * zn * (z11 + z13);
  272.             sl = -zn * s3 * (z1 + z3 - 14 - 6 * e0sq);
  273.             sgh = s4 * zn * (z31 + z33 - 6);
  274.             if (tle.getI() < (FastMath.PI / 60.0)) {
  275.                 // inclination smaller than 3 degrees
  276.                 sh = 0;
  277.             } else {
  278.                 sh = -zn * s2 * (z21 + z23);
  279.             }
  280.             ee2  =  2 * s1 * s6;
  281.             e3   =  2 * s1 * s7;
  282.             xi2  =  2 * s2 * z12;
  283.             xi3  =  2 * s2 * (z13 - z11);
  284.             xl2  = -2 * s3 * z2;
  285.             xl3  = -2 * s3 * (z3 - z1);
  286.             xl4  = -2 * s3 * (-21 - 9 * e0sq) * ze;
  287.             xgh2 =  2 * s4 * z32;
  288.             xgh3 =  2 * s4 * (z33 - z31);
  289.             xgh4 = -18 * s4 * ze;
  290.             xh2  = -2 * s2 * z22;
  291.             xh3  = -2 * s2 * (z23 - z21);

  292.             if (iteration == 0) { // we compute lunar terms only on the first pass:
  293.                 sse = se;
  294.                 ssi = si;
  295.                 ssl = sl;
  296.                 ssh = (tle.getI() < (FastMath.PI / 60.0)) ? 0 : sh / sini0;
  297.                 ssg = sgh - cosi0 * ssh;
  298.                 se2 = ee2;
  299.                 si2 = xi2;
  300.                 sl2 = xl2;
  301.                 sgh2 = xgh2;
  302.                 sh2 = xh2;
  303.                 se3 = e3;
  304.                 si3 = xi3;
  305.                 sl3 = xl3;
  306.                 sgh3 = xgh3;
  307.                 sh3 = xh3;
  308.                 sl4 = xl4;
  309.                 sgh4 = xgh4;
  310.                 zcosg = zcosgl;
  311.                 zsing = zsingl;
  312.                 zcosi = zcosil;
  313.                 zsini = zsinil;
  314.                 zcosh = zcoshl * cosq + zsinhl * sinq;
  315.                 zsinh = sinq * zcoshl - cosq * zsinhl;
  316.                 zn = ZNL;
  317.                 cc = C1L;
  318.                 ze = ZEL;
  319.             }
  320.         } // end of solar - lunar - solar terms computation

  321.         sse += se;
  322.         ssi += si;
  323.         ssl += sl;
  324.         ssg += sgh - ((tle.getI() < (FastMath.PI / 60.0)) ? 0 : (cosi0 / sini0 * sh));
  325.         ssh += (tle.getI() < (FastMath.PI / 60.0)) ? 0 : sh / sini0;



  326.         //        Start the resonant-synchronous tests and initialization

  327.         double bfact = 0;

  328.         // if mean motion is 1.893053 to 2.117652 revs/day, and eccentricity >= 0.5,
  329.         // start of the 12-hour orbit, e > 0.5 section
  330.         if ((xnq >= 0.00826) && (xnq <= 0.00924) && (tle.getE() >= 0.5)) {

  331.             final double g201 = -0.306 - (tle.getE() - 0.64) * 0.440;
  332.             final double eoc = tle.getE() * e0sq;
  333.             final double sini2 = sini0 * sini0;
  334.             final double f220 = 0.75 * (1 + 2 * cosi0 + theta2);
  335.             final double f221 = 1.5 * sini2;
  336.             final double f321 =  1.875 * sini0 * (1 - 2 * cosi0 - 3 * theta2);
  337.             final double f322 = -1.875 * sini0 * (1 + 2 * cosi0 - 3 * theta2);
  338.             final double f441 = 35 * sini2 * f220;
  339.             final double f442 = 39.3750 * sini2 * sini2;
  340.             final double f522 = 9.84375 * sini0 * (sini2 * (1 - 2 * cosi0 - 5 * theta2) +
  341.                                                    0.33333333 * (-2 + 4 * cosi0 + 6 * theta2));
  342.             final double f523 = sini0 * (4.92187512 * sini2 * (-2 - 4 * cosi0 + 10 * theta2) +
  343.                                          6.56250012 * (1 + 2 * cosi0 - 3 * theta2));
  344.             final double f542 = 29.53125 * sini0 * (2 - 8 * cosi0 + theta2 * (-12 + 8 * cosi0 + 10 * theta2));
  345.             final double f543 = 29.53125 * sini0 * (-2 - 8 * cosi0 + theta2 * (12 + 8 * cosi0 - 10 * theta2));
  346.             final double g211;
  347.             final double g310;
  348.             final double g322;
  349.             final double g410;
  350.             final double g422;
  351.             final double g520;

  352.             resonant = true;       // it is resonant...
  353.             synchronous = false;     // but it's not synchronous

  354.             // Geopotential resonance initialization for 12 hour orbits :
  355.             if (tle.getE() <= 0.65) {
  356.                 g211 =    3.616  -   13.247  * tle.getE() +   16.290  * e0sq;
  357.                 g310 =  -19.302  +  117.390  * tle.getE() -  228.419  * e0sq +  156.591  * eoc;
  358.                 g322 =  -18.9068 +  109.7927 * tle.getE() -  214.6334 * e0sq +  146.5816 * eoc;
  359.                 g410 =  -41.122  +  242.694  * tle.getE() -  471.094  * e0sq +  313.953  * eoc;
  360.                 g422 = -146.407  +  841.880  * tle.getE() - 1629.014  * e0sq + 1083.435  * eoc;
  361.                 g520 = -532.114  + 3017.977  * tle.getE() - 5740.032  * e0sq + 3708.276  * eoc;
  362.             } else  {
  363.                 g211 =   -72.099 +   331.819 * tle.getE() -   508.738 * e0sq +   266.724 * eoc;
  364.                 g310 =  -346.844 +  1582.851 * tle.getE() -  2415.925 * e0sq +  1246.113 * eoc;
  365.                 g322 =  -342.585 +  1554.908 * tle.getE() -  2366.899 * e0sq +  1215.972 * eoc;
  366.                 g410 = -1052.797 +  4758.686 * tle.getE() -  7193.992 * e0sq +  3651.957 * eoc;
  367.                 g422 = -3581.69  + 16178.11  * tle.getE() - 24462.77  * e0sq + 12422.52  * eoc;
  368.                 if (tle.getE() <= 0.715) {
  369.                     g520 = 1464.74 - 4664.75 * tle.getE() + 3763.64 * e0sq;
  370.                 } else {
  371.                     g520 = -5149.66 + 29936.92 * tle.getE() - 54087.36 * e0sq + 31324.56 * eoc;
  372.                 }
  373.             }

  374.             final double g533;
  375.             final double g521;
  376.             final double g532;
  377.             if (tle.getE() < 0.7) {
  378.                 g533 = -919.2277  + 4988.61   * tle.getE() - 9064.77   * e0sq + 5542.21  * eoc;
  379.                 g521 = -822.71072 + 4568.6173 * tle.getE() - 8491.4146 * e0sq + 5337.524 * eoc;
  380.                 g532 = -853.666   + 4690.25   * tle.getE() - 8624.77   * e0sq + 5341.4   * eoc;
  381.             } else {
  382.                 g533 = -37995.78  + 161616.52 * tle.getE() - 229838.2  * e0sq + 109377.94 * eoc;
  383.                 g521 = -51752.104 + 218913.95 * tle.getE() - 309468.16 * e0sq + 146349.42 * eoc;
  384.                 g532 = -40023.88  + 170470.89 * tle.getE() - 242699.48 * e0sq + 115605.82 * eoc;
  385.             }

  386.             double temp1 = 3 * xnq * xnq * aqnv * aqnv;
  387.             double temp = temp1 * ROOT22;
  388.             d2201 = temp * f220 * g201;
  389.             d2211 = temp * f221 * g211;
  390.             temp1 *= aqnv;
  391.             temp = temp1 * ROOT32;
  392.             d3210 = temp * f321 * g310;
  393.             d3222 = temp * f322 * g322;
  394.             temp1 *= aqnv;
  395.             temp = 2 * temp1 * ROOT44;
  396.             d4410 = temp * f441 * g410;
  397.             d4422 = temp * f442 * g422;
  398.             temp1 *= aqnv;
  399.             temp = temp1 * ROOT52;
  400.             d5220 = temp * f522 * g520;
  401.             d5232 = temp * f523 * g532;
  402.             temp = 2 * temp1 * ROOT54;
  403.             d5421 = temp * f542 * g521;
  404.             d5433 = temp * f543 * g533;
  405.             xlamo = tle.getMeanAnomaly() + tle.getRaan() + tle.getRaan() - thgr - thgr;
  406.             bfact = xmdot + xnodot + xnodot - THDT - THDT;
  407.             bfact += ssl + ssh + ssh;
  408.         } else if ((xnq < 0.0052359877) && (xnq > 0.0034906585)) {
  409.             // if mean motion is .8 to 1.2 revs/day : (geosynch)

  410.             final double cosio_plus_1 = 1.0 + cosi0;
  411.             final double g200 = 1 + e0sq * (-2.5 + 0.8125  * e0sq);
  412.             final double g300 = 1 + e0sq * (-6   + 6.60937 * e0sq);
  413.             final double f311 = 0.9375 * sini0 * sini0 * (1 + 3 * cosi0) - 0.75 * cosio_plus_1;
  414.             final double g310 = 1 + 2 * e0sq;
  415.             final double f220 = 0.75 * cosio_plus_1 * cosio_plus_1;
  416.             final double f330 = 2.5 * f220 * cosio_plus_1;

  417.             resonant = true;
  418.             synchronous = true;

  419.             // Synchronous resonance terms initialization
  420.             del1 = 3 * xnq * xnq * aqnv * aqnv;
  421.             del2 = 2 * del1 * f220 * g200 * Q22;
  422.             del3 = 3 * del1 * f330 * g300 * Q33 * aqnv;
  423.             del1 = del1 * f311 * g310 * Q31 * aqnv;
  424.             xlamo = tle.getMeanAnomaly() + tle.getRaan() + tle.getPerigeeArgument() - thgr;
  425.             bfact = xmdot + omgdot + xnodot - THDT;
  426.             bfact = bfact + ssl + ssg + ssh;
  427.         } else {
  428.             // it's neither a high-e 12-hours orbit nor a geosynchronous:
  429.             resonant = false;
  430.             synchronous = false;
  431.         }

  432.         if (resonant) {
  433.             xfact = bfact - xnq;

  434.             // Initialize integrator
  435.             xli   = xlamo;
  436.             xni   = xnq;
  437.             atime = 0;
  438.         }
  439.         derivs = new double[2];
  440.     }

  441.     /** Computes secular terms from current coordinates and epoch.
  442.      * @param t offset from initial epoch (minutes)
  443.      */
  444.     protected void deepSecularEffects(final double t)  {

  445.         xll    += ssl * t;
  446.         omgadf += ssg * t;
  447.         xnode  += ssh * t;
  448.         em      = tle.getE() + sse * t;
  449.         xinc    = tle.getI() + ssi * t;

  450.         if (resonant) {
  451.             // If we're closer to t = 0 than to the currently-stored data
  452.             // from the previous call to this function,  then we're
  453.             // better off "restarting",  going back to the initial data.
  454.             // The Dundee code rigs things up to _always_ take 720-minute
  455.             // steps from epoch to end time,  except for the final step.
  456.             // Easiest way to arrange similar behavior in this code is
  457.             // just to always do a restart,  if we're in Dundee-compliant
  458.             // mode.
  459.             if (FastMath.abs(t) < FastMath.abs(t - atime) || isDundeeCompliant)  {
  460.                 // Epoch restart
  461.                 atime = 0;
  462.                 xni = xnq;
  463.                 xli = xlamo;
  464.             }
  465.             boolean lastIntegrationStep = false;
  466.             // if |step|>|step max| then do one step at step max
  467.             while (!lastIntegrationStep) {
  468.                 double delt = t - atime;
  469.                 if (delt > SECULAR_INTEGRATION_STEP) {
  470.                     delt = SECULAR_INTEGRATION_STEP;
  471.                 } else if (delt < -SECULAR_INTEGRATION_STEP) {
  472.                     delt = -SECULAR_INTEGRATION_STEP;
  473.                 } else {
  474.                     lastIntegrationStep = true;
  475.                 }

  476.                 computeSecularDerivs();

  477.                 final double xldot = xni + xfact;

  478.                 double xlpow = 1.;
  479.                 xli += delt * xldot;
  480.                 xni += delt * derivs[0];
  481.                 double delt_factor = delt;
  482.                 xlpow *= xldot;
  483.                 derivs[1] *= xlpow;
  484.                 delt_factor *= delt / 2;
  485.                 xli += delt_factor * derivs[0];
  486.                 xni += delt_factor * derivs[1];
  487.                 atime += delt;
  488.             }
  489.             xn = xni;
  490.             final double temp = -xnode + thgr + t * THDT;
  491.             xll = xli + temp + (synchronous ? -omgadf : temp);
  492.         }
  493.     }

  494.     /** Computes periodic terms from current coordinates and epoch.
  495.      * @param t offset from initial epoch (min)
  496.      */
  497.     protected void deepPeriodicEffects(final double t)  {

  498.         // If the time didn't change by more than 30 minutes,
  499.         // there's no good reason to recompute the perturbations;
  500.         // they don't change enough over so short a time span.
  501.         // However,  the Dundee code _always_ recomputes,  so if
  502.         // we're attempting to replicate its results,  we've gotta
  503.         // recompute everything,  too.
  504.         if ((FastMath.abs(savtsn - t) >= 30.0) || isDundeeCompliant)  {

  505.             savtsn = t;

  506.             // Update solar perturbations for time T
  507.             double zm = zmos + ZNS * t;
  508.             double zf = zm + 2 * ZES * FastMath.sin(zm);
  509.             double sinzf = FastMath.sin(zf);
  510.             double f2 = 0.5 * sinzf * sinzf - 0.25;
  511.             double f3 = -0.5 * sinzf * FastMath.cos(zf);
  512.             final double ses = se2 * f2 + se3 * f3;
  513.             final double sis = si2 * f2 + si3 * f3;
  514.             final double sls = sl2 * f2 + sl3 * f3 + sl4 * sinzf;
  515.             final double sghs = sgh2 * f2 + sgh3 * f3 + sgh4 * sinzf;
  516.             final double shs = sh2 * f2 + sh3 * f3;

  517.             // Update lunar perturbations for time T
  518.             zm = zmol + ZNL * t;
  519.             zf = zm + 2 * ZEL * FastMath.sin(zm);
  520.             sinzf = FastMath.sin(zf);
  521.             f2 =  0.5 * sinzf * sinzf - 0.25;
  522.             f3 = -0.5 * sinzf * FastMath.cos(zf);
  523.             final double sel = ee2 * f2 + e3 * f3;
  524.             final double sil = xi2 * f2 + xi3 * f3;
  525.             final double sll = xl2 * f2 + xl3 * f3 + xl4 * sinzf;
  526.             final double sghl = xgh2 * f2 + xgh3 * f3 + xgh4 * sinzf;
  527.             final double sh1 = xh2 * f2 + xh3 * f3;

  528.             // Sum the solar and lunar contributions
  529.             pe   = ses  + sel;
  530.             pinc = sis  + sil;
  531.             pl   = sls  + sll;
  532.             pgh  = sghs + sghl;
  533.             ph   = shs  + sh1;
  534.         }

  535.         xinc += pinc;

  536.         final double sinis = FastMath.sin( xinc);
  537.         final double cosis = FastMath.cos( xinc);

  538.         /* Add solar/lunar perturbation correction to eccentricity: */
  539.         em     += pe;
  540.         xll    += pl;
  541.         omgadf += pgh;
  542.         xinc    = MathUtils.normalizeAngle(xinc, 0);

  543.         if (FastMath.abs(xinc) >= 0.2) {
  544.             // Apply periodics directly
  545.             final double temp_val = ph / sinis;
  546.             omgadf -= cosis * temp_val;
  547.             xnode += temp_val;
  548.         } else {
  549.             // Apply periodics with Lyddane modification
  550.             final double sinok = FastMath.sin(xnode);
  551.             final double cosok = FastMath.cos(xnode);
  552.             final double alfdp =  ph * cosok + (pinc * cosis + sinis) * sinok;
  553.             final double betdp = -ph * sinok + (pinc * cosis + sinis) * cosok;
  554.             final double delta_xnode = MathUtils.normalizeAngle(FastMath.atan2(alfdp, betdp) - xnode, 0);
  555.             final double dls = -xnode * sinis * pinc;
  556.             omgadf += dls - cosis * delta_xnode;
  557.             xnode  += delta_xnode;
  558.         }
  559.     }

  560.     /** Computes internal secular derivs. */
  561.     private void computeSecularDerivs() {

  562.         final double sin_li = FastMath.sin(xli);
  563.         final double cos_li = FastMath.cos(xli);
  564.         final double sin_2li = 2. * sin_li * cos_li;
  565.         final double cos_2li = 2. * cos_li * cos_li - 1.;

  566.         // Dot terms calculated :
  567.         if (synchronous)  {
  568.             final double sin_3li = sin_2li * cos_li + cos_2li * sin_li;
  569.             final double cos_3li = cos_2li * cos_li - sin_2li * sin_li;
  570.             final double term1a = del1 * (sin_li  * C_FASX2  - cos_li  * S_FASX2);
  571.             final double term2a = del2 * (sin_2li * C_2FASX4 - cos_2li * S_2FASX4);
  572.             final double term3a = del3 * (sin_3li * C_3FASX6 - cos_3li * S_3FASX6);
  573.             final double term1b = del1 * (cos_li  * C_FASX2  + sin_li  * S_FASX2);
  574.             final double term2b = 2.0 * del2 * (cos_2li * C_2FASX4 + sin_2li * S_2FASX4);
  575.             final double term3b = 3.0 * del3 * (cos_3li * C_3FASX6 + sin_3li * S_3FASX6);
  576.             derivs[0] = term1a + term2a + term3a;
  577.             derivs[1] = term1b + term2b + term3b;
  578.         } else {
  579.             // orbit is a 12-hour resonant one
  580.             final double xomi = omegaq + omgdot * atime;
  581.             final double sin_omi = FastMath.sin(xomi);
  582.             final double cos_omi = FastMath.cos(xomi);
  583.             final double sin_li_m_omi = sin_li * cos_omi - sin_omi * cos_li;
  584.             final double sin_li_p_omi = sin_li * cos_omi + sin_omi * cos_li;
  585.             final double cos_li_m_omi = cos_li * cos_omi + sin_omi * sin_li;
  586.             final double cos_li_p_omi = cos_li * cos_omi - sin_omi * sin_li;
  587.             final double sin_2omi = 2. * sin_omi * cos_omi;
  588.             final double cos_2omi = 2. * cos_omi * cos_omi - 1.;
  589.             final double sin_2li_m_omi = sin_2li * cos_omi - sin_omi * cos_2li;
  590.             final double sin_2li_p_omi = sin_2li * cos_omi + sin_omi * cos_2li;
  591.             final double cos_2li_m_omi = cos_2li * cos_omi + sin_omi * sin_2li;
  592.             final double cos_2li_p_omi = cos_2li * cos_omi - sin_omi * sin_2li;
  593.             final double sin_2li_p_2omi = sin_2li * cos_2omi + sin_2omi * cos_2li;
  594.             final double cos_2li_p_2omi = cos_2li * cos_2omi - sin_2omi * sin_2li;
  595.             final double sin_2omi_p_li = sin_li * cos_2omi + sin_2omi * cos_li;
  596.             final double cos_2omi_p_li = cos_li * cos_2omi - sin_2omi * sin_li;
  597.             final double term1a = d2201 * (sin_2omi_p_li * C_G22 - cos_2omi_p_li * S_G22) +
  598.                                   d2211 * (sin_li * C_G22 - cos_li * S_G22) +
  599.                                   d3210 * (sin_li_p_omi * C_G32 - cos_li_p_omi * S_G32) +
  600.                                   d3222 * (sin_li_m_omi * C_G32 - cos_li_m_omi * S_G32) +
  601.                                   d5220 * (sin_li_p_omi * C_G52 - cos_li_p_omi * S_G52) +
  602.                                   d5232 * (sin_li_m_omi * C_G52 - cos_li_m_omi * S_G52);
  603.             final double term2a = d4410 * (sin_2li_p_2omi * C_G44 - cos_2li_p_2omi * S_G44) +
  604.                                   d4422 * (sin_2li * C_G44 - cos_2li * S_G44) +
  605.                                   d5421 * (sin_2li_p_omi * C_G54 - cos_2li_p_omi * S_G54) +
  606.                                   d5433 * (sin_2li_m_omi * C_G54 - cos_2li_m_omi * S_G54);
  607.             final double term1b = d2201 * (cos_2omi_p_li * C_G22 + sin_2omi_p_li * S_G22) +
  608.                                   d2211 * (cos_li * C_G22 + sin_li * S_G22) +
  609.                                   d3210 * (cos_li_p_omi * C_G32 + sin_li_p_omi * S_G32) +
  610.                                   d3222 * (cos_li_m_omi * C_G32 + sin_li_m_omi * S_G32) +
  611.                                   d5220 * (cos_li_p_omi * C_G52 + sin_li_p_omi * S_G52) +
  612.                                   d5232 * (cos_li_m_omi * C_G52 + sin_li_m_omi * S_G52);
  613.             final double term2b = 2.0 * (d4410 * (cos_2li_p_2omi * C_G44 + sin_2li_p_2omi * S_G44) +
  614.                                          d4422 * (cos_2li * C_G44 + sin_2li * S_G44) +
  615.                                          d5421 * (cos_2li_p_omi * C_G54 + sin_2li_p_omi * S_G54) +
  616.                                          d5433 * (cos_2li_m_omi * C_G54 + sin_2li_m_omi * S_G54));

  617.             derivs[0] = term1a + term2a;
  618.             derivs[1] = term1b + term2b;

  619.         }
  620.     }

  621. }