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  
19  import org.hipparchus.util.FastMath;
20  import org.hipparchus.util.MathUtils;
21  import org.orekit.attitudes.AttitudeProvider;
22  import org.orekit.time.AbsoluteDate;
23  import org.orekit.time.TimeScalesFactory;
24  import org.orekit.utils.Constants;
25  
26  
27  /** This class contains the methods that compute deep space perturbation terms.
28   * <p>
29   * The user should not bother in this class since it is handled internaly by the
30   * {@link TLEPropagator}.
31   * </p>
32   * <p>This implementation is largely inspired from the paper and source code <a
33   * href="http://www.celestrak.com/publications/AIAA/2006-6753/">Revisiting Spacetrack
34   * Report #3</a> and is fully compliant with its results and tests cases.</p>
35   * @author Felix R. Hoots, Ronald L. Roehrich, December 1980 (original fortran)
36   * @author David A. Vallado, Paul Crawford, Richard Hujsak, T.S. Kelso (C++ translation and improvements)
37   * @author Fabien Maussion (java translation)
38   */
39  public class DeepSDP4 extends SDP4 {
40  
41      // CHECKSTYLE: stop JavadocVariable check
42  
43      // Internal constants
44      private static final double ZNS      = 1.19459E-5;
45      private static final double ZES      = 0.01675;
46      private static final double ZNL      = 1.5835218E-4;
47      private static final double ZEL      = 0.05490;
48      private static final double THDT     = 4.3752691E-3;
49      private static final double C1SS     =  2.9864797E-6;
50      private static final double C1L      = 4.7968065E-7;
51  
52      private static final double ROOT22   = 1.7891679E-6;
53      private static final double ROOT32   = 3.7393792E-7;
54      private static final double ROOT44   = 7.3636953E-9;
55      private static final double ROOT52   = 1.1428639E-7;
56      private static final double ROOT54   = 2.1765803E-9;
57  
58      private static final double Q22      =  1.7891679E-6;
59      private static final double Q31      =  2.1460748E-6;
60      private static final double Q33      =  2.2123015E-7;
61  
62      private static final double C_FASX2  =  0.99139134268488593;
63      private static final double S_FASX2  =  0.13093206501640101;
64      private static final double C_2FASX4 =  0.87051638752972937;
65      private static final double S_2FASX4 = -0.49213943048915526;
66      private static final double C_3FASX6 =  0.43258117585763334;
67      private static final double S_3FASX6 =  0.90159499016666422;
68  
69      private static final double C_G22    =  0.87051638752972937;
70      private static final double S_G22    = -0.49213943048915526;
71      private static final double C_G32    =  0.57972190187001149;
72      private static final double S_G32    =  0.81481440616389245;
73      private static final double C_G44    = -0.22866241528815548;
74      private static final double S_G44    =  0.97350577801807991;
75      private static final double C_G52    =  0.49684831179884198;
76      private static final double S_G52    =  0.86783740128127729;
77      private static final double C_G54    = -0.29695209575316894;
78      private static final double S_G54    = -0.95489237761529999;
79  
80      /** Integration step (seconds). */
81      private static final double SECULAR_INTEGRATION_STEP  = 720.0;
82  
83      /** Intermediate values. */
84      private double thgr;
85      private double xnq;
86      private double omegaq;
87      private double zcosil;
88      private double zsinil;
89      private double zsinhl;
90      private double zcoshl;
91      private double zmol;
92      private double zcosgl;
93      private double zsingl;
94      private double zmos;
95      private double savtsn;
96  
97      private double ee2;
98      private double e3;
99      private double xi2;
100     private double xi3;
101     private double xl2;
102     private double xl3;
103     private double xl4;
104     private double xgh2;
105     private double xgh3;
106     private double xgh4;
107     private double xh2;
108     private double xh3;
109 
110     private double d2201;
111     private double d2211;
112     private double d3210;
113     private double d3222;
114     private double d4410;
115     private double d4422;
116     private double d5220;
117     private double d5232;
118     private double d5421;
119     private double d5433;
120     private double xlamo;
121 
122     private double sse;
123     private double ssi;
124     private double ssl;
125     private double ssh;
126     private double ssg;
127     private double se2;
128     private double si2;
129     private double sl2;
130     private double sgh2;
131     private double sh2;
132     private double se3;
133     private double si3;
134     private double sl3;
135     private double sgh3;
136     private double sh3;
137     private double sl4;
138     private double sgh4;
139 
140     private double del1;
141     private double del2;
142     private double del3;
143     private double xfact;
144     private double xli;
145     private double xni;
146     private double atime;
147 
148     private double pe;
149     private double pinc;
150     private double pl;
151     private double pgh;
152     private double ph;
153 
154     private double[] derivs;
155 
156     // CHECKSTYLE: resume JavadocVariable check
157 
158     /** Flag for resonant orbits. */
159     private boolean resonant;
160 
161     /** Flag for synchronous orbits. */
162     private boolean synchronous;
163 
164     /** Flag for compliance with Dundee modifications. */
165     private boolean isDundeeCompliant = true;
166 
167     /** Constructor for a unique initial TLE.
168      * @param initialTLE the TLE to propagate.
169      * @param attitudeProvider provider for attitude computation
170      * @param mass spacecraft mass (kg)
171      */
172     public DeepSDP4(final TLE initialTLE, final AttitudeProvider attitudeProvider,
173                        final double mass) {
174         super(initialTLE, attitudeProvider, mass);
175     }
176 
177     /** Computes luni - solar terms from initial coordinates and epoch.
178      */
179     protected void luniSolarTermsComputation() {
180 
181         final double sing = FastMath.sin(tle.getPerigeeArgument());
182         final double cosg = FastMath.cos(tle.getPerigeeArgument());
183 
184         final double sinq = FastMath.sin(tle.getRaan());
185         final double cosq = FastMath.cos(tle.getRaan());
186         final double aqnv = 1.0 / a0dp;
187 
188         // Compute julian days since 1900
189         final double daysSince1900 =
190             (tle.getDate().durationFrom(AbsoluteDate.JULIAN_EPOCH) +
191              tle.getDate().timeScalesOffset(TimeScalesFactory.getUTC(), TimeScalesFactory.getTT())) / Constants.JULIAN_DAY - 2415020;
192 
193 
194         double cc = C1SS;
195         double ze = ZES;
196         double zn = ZNS;
197         double zsinh = sinq;
198         double zcosh = cosq;
199 
200         thgr = thetaG(tle.getDate());
201         xnq = xn0dp;
202         omegaq = tle.getPerigeeArgument();
203 
204         final double xnodce = 4.5236020 - 9.2422029e-4 * daysSince1900;
205         final double stem = FastMath.sin(xnodce);
206         final double ctem = FastMath.cos(xnodce);
207         final double c_minus_gam = 0.228027132 * daysSince1900 - 1.1151842;
208         final double gam = 5.8351514 + 0.0019443680 * daysSince1900;
209 
210         zcosil = 0.91375164 - 0.03568096 * ctem;
211         zsinil = FastMath.sqrt(1.0 - zcosil * zcosil);
212         zsinhl = 0.089683511 * stem / zsinil;
213         zcoshl = FastMath.sqrt(1.0 - zsinhl * zsinhl);
214         zmol = MathUtils.normalizeAngle(c_minus_gam, FastMath.PI);
215 
216         double zx = 0.39785416 * stem / zsinil;
217         final double zy = zcoshl * ctem + 0.91744867 * zsinhl * stem;
218         zx = FastMath.atan2( zx, zy) + gam - xnodce;
219         zcosgl = FastMath.cos( zx);
220         zsingl = FastMath.sin( zx);
221         zmos = MathUtils.normalizeAngle(6.2565837 + 0.017201977 * daysSince1900, FastMath.PI);
222 
223         // Do solar terms
224         savtsn = 1e20;
225 
226         double zcosi =  0.91744867;
227         double zsini =  0.39785416;
228         double zsing = -0.98088458;
229         double zcosg =  0.1945905;
230 
231         double se = 0;
232         double sgh = 0;
233         double sh = 0;
234         double si = 0;
235         double sl = 0;
236 
237         // There was previously some convoluted logic here, but it boils
238         // down to this:  we compute the solar terms,  then the lunar terms.
239         // On a second pass,  we recompute the solar terms, taking advantage
240         // of the improved data that resulted from computing lunar terms.
241         for (int iteration = 0; iteration < 2; ++iteration) {
242             final double a1 = zcosg * zcosh + zsing * zcosi * zsinh;
243             final double a3 = -zsing * zcosh + zcosg * zcosi * zsinh;
244             final double a7 = -zcosg * zsinh + zsing * zcosi * zcosh;
245             final double a8 = zsing * zsini;
246             final double a9 = zsing * zsinh + zcosg * zcosi * zcosh;
247             final double a10 = zcosg * zsini;
248             final double a2 = cosi0 * a7 + sini0 * a8;
249             final double a4 = cosi0 * a9 + sini0 * a10;
250             final double a5 = -sini0 * a7 + cosi0 * a8;
251             final double a6 = -sini0 * a9 + cosi0 * a10;
252             final double x1 = a1 * cosg + a2 * sing;
253             final double x2 = a3 * cosg + a4 * sing;
254             final double x3 = -a1 * sing + a2 * cosg;
255             final double x4 = -a3 * sing + a4 * cosg;
256             final double x5 = a5 * sing;
257             final double x6 = a6 * sing;
258             final double x7 = a5 * cosg;
259             final double x8 = a6 * cosg;
260             final double z31 = 12 * x1 * x1 - 3 * x3 * x3;
261             final double z32 = 24 * x1 * x2 - 6 * x3 * x4;
262             final double z33 = 12 * x2 * x2 - 3 * x4 * x4;
263             final double z11 = -6 * a1 * a5 + e0sq * (-24 * x1 * x7 - 6 * x3 * x5);
264             final double z12 = -6 * (a1 * a6 + a3 * a5) +
265                                e0sq * (-24 * (x2 * x7 + x1 * x8) - 6 * (x3 * x6 + x4 * x5));
266             final double z13 = -6 * a3 * a6 + e0sq * (-24 * x2 * x8 - 6 * x4 * x6);
267             final double z21 = 6 * a2 * a5 + e0sq * (24 * x1 * x5 - 6 * x3 * x7);
268             final double z22 = 6 * (a4 * a5 + a2 * a6) +
269                                e0sq * (24 * (x2 * x5 + x1 * x6) - 6 * (x4 * x7 + x3 * x8));
270             final double z23 = 6 * a4 * a6 + e0sq * (24 * x2 * x6 - 6 * x4 * x8);
271             final double s3 = cc / xnq;
272             final double s2 = -0.5 * s3 / beta0;
273             final double s4 = s3 * beta0;
274             final double s1 = -15 * tle.getE() * s4;
275             final double s5 = x1 * x3 + x2 * x4;
276             final double s6 = x2 * x3 + x1 * x4;
277             final double s7 = x2 * x4 - x1 * x3;
278             double z1 = 3 * (a1 * a1 + a2 * a2) + z31 * e0sq;
279             double z2 = 6 * (a1 * a3 + a2 * a4) + z32 * e0sq;
280             double z3 = 3 * (a3 * a3 + a4 * a4) + z33 * e0sq;
281 
282             z1 = z1 + z1 + beta02 * z31;
283             z2 = z2 + z2 + beta02 * z32;
284             z3 = z3 + z3 + beta02 * z33;
285             se = s1 * zn * s5;
286             si = s2 * zn * (z11 + z13);
287             sl = -zn * s3 * (z1 + z3 - 14 - 6 * e0sq);
288             sgh = s4 * zn * (z31 + z33 - 6);
289             if (tle.getI() < (FastMath.PI / 60.0)) {
290                 // inclination smaller than 3 degrees
291                 sh = 0;
292             } else {
293                 sh = -zn * s2 * (z21 + z23);
294             }
295             ee2  =  2 * s1 * s6;
296             e3   =  2 * s1 * s7;
297             xi2  =  2 * s2 * z12;
298             xi3  =  2 * s2 * (z13 - z11);
299             xl2  = -2 * s3 * z2;
300             xl3  = -2 * s3 * (z3 - z1);
301             xl4  = -2 * s3 * (-21 - 9 * e0sq) * ze;
302             xgh2 =  2 * s4 * z32;
303             xgh3 =  2 * s4 * (z33 - z31);
304             xgh4 = -18 * s4 * ze;
305             xh2  = -2 * s2 * z22;
306             xh3  = -2 * s2 * (z23 - z21);
307 
308             if (iteration == 0) { // we compute lunar terms only on the first pass:
309                 sse = se;
310                 ssi = si;
311                 ssl = sl;
312                 ssh = (tle.getI() < (FastMath.PI / 60.0)) ? 0 : sh / sini0;
313                 ssg = sgh - cosi0 * ssh;
314                 se2 = ee2;
315                 si2 = xi2;
316                 sl2 = xl2;
317                 sgh2 = xgh2;
318                 sh2 = xh2;
319                 se3 = e3;
320                 si3 = xi3;
321                 sl3 = xl3;
322                 sgh3 = xgh3;
323                 sh3 = xh3;
324                 sl4 = xl4;
325                 sgh4 = xgh4;
326                 zcosg = zcosgl;
327                 zsing = zsingl;
328                 zcosi = zcosil;
329                 zsini = zsinil;
330                 zcosh = zcoshl * cosq + zsinhl * sinq;
331                 zsinh = sinq * zcoshl - cosq * zsinhl;
332                 zn = ZNL;
333                 cc = C1L;
334                 ze = ZEL;
335             }
336         } // end of solar - lunar - solar terms computation
337 
338         sse += se;
339         ssi += si;
340         ssl += sl;
341         ssg += sgh - ((tle.getI() < (FastMath.PI / 60.0)) ? 0 : (cosi0 / sini0 * sh));
342         ssh += (tle.getI() < (FastMath.PI / 60.0)) ? 0 : sh / sini0;
343 
344 
345 
346         //        Start the resonant-synchronous tests and initialization
347 
348         double bfact = 0;
349 
350         // if mean motion is 1.893053 to 2.117652 revs/day, and eccentricity >= 0.5,
351         // start of the 12-hour orbit, e > 0.5 section
352         if ((xnq >= 0.00826) && (xnq <= 0.00924) && (tle.getE() >= 0.5)) {
353 
354             final double g201 = -0.306 - (tle.getE() - 0.64) * 0.440;
355             final double eoc = tle.getE() * e0sq;
356             final double sini2 = sini0 * sini0;
357             final double f220 = 0.75 * (1 + 2 * cosi0 + theta2);
358             final double f221 = 1.5 * sini2;
359             final double f321 =  1.875 * sini0 * (1 - 2 * cosi0 - 3 * theta2);
360             final double f322 = -1.875 * sini0 * (1 + 2 * cosi0 - 3 * theta2);
361             final double f441 = 35 * sini2 * f220;
362             final double f442 = 39.3750 * sini2 * sini2;
363             final double f522 = 9.84375 * sini0 * (sini2 * (1 - 2 * cosi0 - 5 * theta2) +
364                                                    0.33333333 * (-2 + 4 * cosi0 + 6 * theta2));
365             final double f523 = sini0 * (4.92187512 * sini2 * (-2 - 4 * cosi0 + 10 * theta2) +
366                                          6.56250012 * (1 + 2 * cosi0 - 3 * theta2));
367             final double f542 = 29.53125 * sini0 * (2 - 8 * cosi0 + theta2 * (-12 + 8 * cosi0 + 10 * theta2));
368             final double f543 = 29.53125 * sini0 * (-2 - 8 * cosi0 + theta2 * (12 + 8 * cosi0 - 10 * theta2));
369             final double g211;
370             final double g310;
371             final double g322;
372             final double g410;
373             final double g422;
374             final double g520;
375 
376             resonant = true;       // it is resonant...
377             synchronous = false;     // but it's not synchronous
378 
379             // Geopotential resonance initialization for 12 hour orbits :
380             if (tle.getE() <= 0.65) {
381                 g211 =    3.616  -   13.247  * tle.getE() +   16.290  * e0sq;
382                 g310 =  -19.302  +  117.390  * tle.getE() -  228.419  * e0sq +  156.591  * eoc;
383                 g322 =  -18.9068 +  109.7927 * tle.getE() -  214.6334 * e0sq +  146.5816 * eoc;
384                 g410 =  -41.122  +  242.694  * tle.getE() -  471.094  * e0sq +  313.953  * eoc;
385                 g422 = -146.407  +  841.880  * tle.getE() - 1629.014  * e0sq + 1083.435  * eoc;
386                 g520 = -532.114  + 3017.977  * tle.getE() - 5740.032  * e0sq + 3708.276  * eoc;
387             } else  {
388                 g211 =   -72.099 +   331.819 * tle.getE() -   508.738 * e0sq +   266.724 * eoc;
389                 g310 =  -346.844 +  1582.851 * tle.getE() -  2415.925 * e0sq +  1246.113 * eoc;
390                 g322 =  -342.585 +  1554.908 * tle.getE() -  2366.899 * e0sq +  1215.972 * eoc;
391                 g410 = -1052.797 +  4758.686 * tle.getE() -  7193.992 * e0sq +  3651.957 * eoc;
392                 g422 = -3581.69  + 16178.11  * tle.getE() - 24462.77  * e0sq + 12422.52  * eoc;
393                 if (tle.getE() <= 0.715) {
394                     g520 = 1464.74 - 4664.75 * tle.getE() + 3763.64 * e0sq;
395                 } else {
396                     g520 = -5149.66 + 29936.92 * tle.getE() - 54087.36 * e0sq + 31324.56 * eoc;
397                 }
398             }
399 
400             final double g533;
401             final double g521;
402             final double g532;
403             if (tle.getE() < 0.7) {
404                 g533 = -919.2277  + 4988.61   * tle.getE() - 9064.77   * e0sq + 5542.21  * eoc;
405                 g521 = -822.71072 + 4568.6173 * tle.getE() - 8491.4146 * e0sq + 5337.524 * eoc;
406                 g532 = -853.666   + 4690.25   * tle.getE() - 8624.77   * e0sq + 5341.4   * eoc;
407             } else {
408                 g533 = -37995.78  + 161616.52 * tle.getE() - 229838.2  * e0sq + 109377.94 * eoc;
409                 g521 = -51752.104 + 218913.95 * tle.getE() - 309468.16 * e0sq + 146349.42 * eoc;
410                 g532 = -40023.88  + 170470.89 * tle.getE() - 242699.48 * e0sq + 115605.82 * eoc;
411             }
412 
413             double temp1 = 3 * xnq * xnq * aqnv * aqnv;
414             double temp = temp1 * ROOT22;
415             d2201 = temp * f220 * g201;
416             d2211 = temp * f221 * g211;
417             temp1 *= aqnv;
418             temp = temp1 * ROOT32;
419             d3210 = temp * f321 * g310;
420             d3222 = temp * f322 * g322;
421             temp1 *= aqnv;
422             temp = 2 * temp1 * ROOT44;
423             d4410 = temp * f441 * g410;
424             d4422 = temp * f442 * g422;
425             temp1 *= aqnv;
426             temp = temp1 * ROOT52;
427             d5220 = temp * f522 * g520;
428             d5232 = temp * f523 * g532;
429             temp = 2 * temp1 * ROOT54;
430             d5421 = temp * f542 * g521;
431             d5433 = temp * f543 * g533;
432             xlamo = tle.getMeanAnomaly() + tle.getRaan() + tle.getRaan() - thgr - thgr;
433             bfact = xmdot + xnodot + xnodot - THDT - THDT;
434             bfact += ssl + ssh + ssh;
435         } else if ((xnq < 0.0052359877) && (xnq > 0.0034906585)) {
436             // if mean motion is .8 to 1.2 revs/day : (geosynch)
437 
438             final double cosio_plus_1 = 1.0 + cosi0;
439             final double g200 = 1 + e0sq * (-2.5 + 0.8125  * e0sq);
440             final double g300 = 1 + e0sq * (-6   + 6.60937 * e0sq);
441             final double f311 = 0.9375 * sini0 * sini0 * (1 + 3 * cosi0) - 0.75 * cosio_plus_1;
442             final double g310 = 1 + 2 * e0sq;
443             final double f220 = 0.75 * cosio_plus_1 * cosio_plus_1;
444             final double f330 = 2.5 * f220 * cosio_plus_1;
445 
446             resonant = true;
447             synchronous = true;
448 
449             // Synchronous resonance terms initialization
450             del1 = 3 * xnq * xnq * aqnv * aqnv;
451             del2 = 2 * del1 * f220 * g200 * Q22;
452             del3 = 3 * del1 * f330 * g300 * Q33 * aqnv;
453             del1 = del1 * f311 * g310 * Q31 * aqnv;
454             xlamo = tle.getMeanAnomaly() + tle.getRaan() + tle.getPerigeeArgument() - thgr;
455             bfact = xmdot + omgdot + xnodot - THDT;
456             bfact = bfact + ssl + ssg + ssh;
457         } else {
458             // it's neither a high-e 12-hours orbit nor a geosynchronous:
459             resonant = false;
460             synchronous = false;
461         }
462 
463         if (resonant) {
464             xfact = bfact - xnq;
465 
466             // Initialize integrator
467             xli   = xlamo;
468             xni   = xnq;
469             atime = 0;
470         }
471         derivs = new double[2];
472     }
473 
474     /** Computes secular terms from current coordinates and epoch.
475      * @param t offset from initial epoch (minutes)
476      */
477     protected void deepSecularEffects(final double t)  {
478 
479         xll    += ssl * t;
480         omgadf += ssg * t;
481         xnode  += ssh * t;
482         em      = tle.getE() + sse * t;
483         xinc    = tle.getI() + ssi * t;
484 
485         if (resonant) {
486             // If we're closer to t = 0 than to the currently-stored data
487             // from the previous call to this function,  then we're
488             // better off "restarting",  going back to the initial data.
489             // The Dundee code rigs things up to _always_ take 720-minute
490             // steps from epoch to end time,  except for the final step.
491             // Easiest way to arrange similar behavior in this code is
492             // just to always do a restart,  if we're in Dundee-compliant
493             // mode.
494             if (FastMath.abs(t) < FastMath.abs(t - atime) || isDundeeCompliant)  {
495                 // Epoch restart
496                 atime = 0;
497                 xni = xnq;
498                 xli = xlamo;
499             }
500             boolean lastIntegrationStep = false;
501             // if |step|>|step max| then do one step at step max
502             while (!lastIntegrationStep) {
503                 double delt = t - atime;
504                 if (delt > SECULAR_INTEGRATION_STEP) {
505                     delt = SECULAR_INTEGRATION_STEP;
506                 } else if (delt < -SECULAR_INTEGRATION_STEP) {
507                     delt = -SECULAR_INTEGRATION_STEP;
508                 } else {
509                     lastIntegrationStep = true;
510                 }
511 
512                 computeSecularDerivs();
513 
514                 final double xldot = xni + xfact;
515 
516                 double xlpow = 1.;
517                 xli += delt * xldot;
518                 xni += delt * derivs[0];
519                 double delt_factor = delt;
520                 xlpow *= xldot;
521                 derivs[1] *= xlpow;
522                 delt_factor *= delt / 2;
523                 xli += delt_factor * derivs[0];
524                 xni += delt_factor * derivs[1];
525                 atime += delt;
526             }
527             xn = xni;
528             final double temp = -xnode + thgr + t * THDT;
529             xll = xli + temp + (synchronous ? -omgadf : temp);
530         }
531     }
532 
533     /** Computes periodic terms from current coordinates and epoch.
534      * @param t offset from initial epoch (min)
535      */
536     protected void deepPeriodicEffects(final double t)  {
537 
538         // If the time didn't change by more than 30 minutes,
539         // there's no good reason to recompute the perturbations;
540         // they don't change enough over so short a time span.
541         // However,  the Dundee code _always_ recomputes,  so if
542         // we're attempting to replicate its results,  we've gotta
543         // recompute everything,  too.
544         if ((FastMath.abs(savtsn - t) >= 30.0) || isDundeeCompliant)  {
545 
546             savtsn = t;
547 
548             // Update solar perturbations for time T
549             double zm = zmos + ZNS * t;
550             double zf = zm + 2 * ZES * FastMath.sin(zm);
551             double sinzf = FastMath.sin(zf);
552             double f2 = 0.5 * sinzf * sinzf - 0.25;
553             double f3 = -0.5 * sinzf * FastMath.cos(zf);
554             final double ses = se2 * f2 + se3 * f3;
555             final double sis = si2 * f2 + si3 * f3;
556             final double sls = sl2 * f2 + sl3 * f3 + sl4 * sinzf;
557             final double sghs = sgh2 * f2 + sgh3 * f3 + sgh4 * sinzf;
558             final double shs = sh2 * f2 + sh3 * f3;
559 
560             // Update lunar perturbations for time T
561             zm = zmol + ZNL * t;
562             zf = zm + 2 * ZEL * FastMath.sin(zm);
563             sinzf = FastMath.sin(zf);
564             f2 =  0.5 * sinzf * sinzf - 0.25;
565             f3 = -0.5 * sinzf * FastMath.cos(zf);
566             final double sel = ee2 * f2 + e3 * f3;
567             final double sil = xi2 * f2 + xi3 * f3;
568             final double sll = xl2 * f2 + xl3 * f3 + xl4 * sinzf;
569             final double sghl = xgh2 * f2 + xgh3 * f3 + xgh4 * sinzf;
570             final double sh1 = xh2 * f2 + xh3 * f3;
571 
572             // Sum the solar and lunar contributions
573             pe   = ses  + sel;
574             pinc = sis  + sil;
575             pl   = sls  + sll;
576             pgh  = sghs + sghl;
577             ph   = shs  + sh1;
578         }
579 
580         xinc += pinc;
581 
582         final double sinis = FastMath.sin( xinc);
583         final double cosis = FastMath.cos( xinc);
584 
585         /* Add solar/lunar perturbation correction to eccentricity: */
586         em     += pe;
587         xll    += pl;
588         omgadf += pgh;
589         xinc    = MathUtils.normalizeAngle(xinc, 0);
590 
591         if (FastMath.abs(xinc) >= 0.2) {
592             // Apply periodics directly
593             final double temp_val = ph / sinis;
594             omgadf -= cosis * temp_val;
595             xnode += temp_val;
596         } else {
597             // Apply periodics with Lyddane modification
598             final double sinok = FastMath.sin(xnode);
599             final double cosok = FastMath.cos(xnode);
600             final double alfdp =  ph * cosok + (pinc * cosis + sinis) * sinok;
601             final double betdp = -ph * sinok + (pinc * cosis + sinis) * cosok;
602             final double delta_xnode = MathUtils.normalizeAngle(FastMath.atan2(alfdp, betdp) - xnode, 0);
603             final double dls = -xnode * sinis * pinc;
604             omgadf += dls - cosis * delta_xnode;
605             xnode  += delta_xnode;
606         }
607     }
608 
609     /** Computes internal secular derivs. */
610     private void computeSecularDerivs() {
611 
612         final double sin_li = FastMath.sin(xli);
613         final double cos_li = FastMath.cos(xli);
614         final double sin_2li = 2. * sin_li * cos_li;
615         final double cos_2li = 2. * cos_li * cos_li - 1.;
616 
617         // Dot terms calculated :
618         if (synchronous)  {
619             final double sin_3li = sin_2li * cos_li + cos_2li * sin_li;
620             final double cos_3li = cos_2li * cos_li - sin_2li * sin_li;
621             final double term1a = del1 * (sin_li  * C_FASX2  - cos_li  * S_FASX2);
622             final double term2a = del2 * (sin_2li * C_2FASX4 - cos_2li * S_2FASX4);
623             final double term3a = del3 * (sin_3li * C_3FASX6 - cos_3li * S_3FASX6);
624             final double term1b = del1 * (cos_li  * C_FASX2  + sin_li  * S_FASX2);
625             final double term2b = 2.0 * del2 * (cos_2li * C_2FASX4 + sin_2li * S_2FASX4);
626             final double term3b = 3.0 * del3 * (cos_3li * C_3FASX6 + sin_3li * S_3FASX6);
627             derivs[0] = term1a + term2a + term3a;
628             derivs[1] = term1b + term2b + term3b;
629         } else {
630             // orbit is a 12-hour resonant one
631             final double xomi = omegaq + omgdot * atime;
632             final double sin_omi = FastMath.sin(xomi);
633             final double cos_omi = FastMath.cos(xomi);
634             final double sin_li_m_omi = sin_li * cos_omi - sin_omi * cos_li;
635             final double sin_li_p_omi = sin_li * cos_omi + sin_omi * cos_li;
636             final double cos_li_m_omi = cos_li * cos_omi + sin_omi * sin_li;
637             final double cos_li_p_omi = cos_li * cos_omi - sin_omi * sin_li;
638             final double sin_2omi = 2. * sin_omi * cos_omi;
639             final double cos_2omi = 2. * cos_omi * cos_omi - 1.;
640             final double sin_2li_m_omi = sin_2li * cos_omi - sin_omi * cos_2li;
641             final double sin_2li_p_omi = sin_2li * cos_omi + sin_omi * cos_2li;
642             final double cos_2li_m_omi = cos_2li * cos_omi + sin_omi * sin_2li;
643             final double cos_2li_p_omi = cos_2li * cos_omi - sin_omi * sin_2li;
644             final double sin_2li_p_2omi = sin_2li * cos_2omi + sin_2omi * cos_2li;
645             final double cos_2li_p_2omi = cos_2li * cos_2omi - sin_2omi * sin_2li;
646             final double sin_2omi_p_li = sin_li * cos_2omi + sin_2omi * cos_li;
647             final double cos_2omi_p_li = cos_li * cos_2omi - sin_2omi * sin_li;
648             final double term1a = d2201 * (sin_2omi_p_li * C_G22 - cos_2omi_p_li * S_G22) +
649                                   d2211 * (sin_li * C_G22 - cos_li * S_G22) +
650                                   d3210 * (sin_li_p_omi * C_G32 - cos_li_p_omi * S_G32) +
651                                   d3222 * (sin_li_m_omi * C_G32 - cos_li_m_omi * S_G32) +
652                                   d5220 * (sin_li_p_omi * C_G52 - cos_li_p_omi * S_G52) +
653                                   d5232 * (sin_li_m_omi * C_G52 - cos_li_m_omi * S_G52);
654             final double term2a = d4410 * (sin_2li_p_2omi * C_G44 - cos_2li_p_2omi * S_G44) +
655                                   d4422 * (sin_2li * C_G44 - cos_2li * S_G44) +
656                                   d5421 * (sin_2li_p_omi * C_G54 - cos_2li_p_omi * S_G54) +
657                                   d5433 * (sin_2li_m_omi * C_G54 - cos_2li_m_omi * S_G54);
658             final double term1b = d2201 * (cos_2omi_p_li * C_G22 + sin_2omi_p_li * S_G22) +
659                                   d2211 * (cos_li * C_G22 + sin_li * S_G22) +
660                                   d3210 * (cos_li_p_omi * C_G32 + sin_li_p_omi * S_G32) +
661                                   d3222 * (cos_li_m_omi * C_G32 + sin_li_m_omi * S_G32) +
662                                   d5220 * (cos_li_p_omi * C_G52 + sin_li_p_omi * S_G52) +
663                                   d5232 * (cos_li_m_omi * C_G52 + sin_li_m_omi * S_G52);
664             final double term2b = 2.0 * (d4410 * (cos_2li_p_2omi * C_G44 + sin_2li_p_2omi * S_G44) +
665                                          d4422 * (cos_2li * C_G44 + sin_2li * S_G44) +
666                                          d5421 * (cos_2li_p_omi * C_G54 + sin_2li_p_omi * S_G54) +
667                                          d5433 * (cos_2li_m_omi * C_G54 + sin_2li_m_omi * S_G54));
668 
669             derivs[0] = term1a + term2a;
670             derivs[1] = term1b + term2b;
671 
672         }
673     }
674 
675 }