1   /* Copyright 2002-2021 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.semianalytical.dsst.forces;
18  
19  import java.lang.reflect.Array;
20  import java.util.ArrayList;
21  import java.util.Arrays;
22  import java.util.Collections;
23  import java.util.HashMap;
24  import java.util.List;
25  import java.util.Map;
26  import java.util.Set;
27  import java.util.SortedMap;
28  import java.util.TreeMap;
29  
30  import org.hipparchus.Field;
31  import org.hipparchus.CalculusFieldElement;
32  import org.hipparchus.analysis.differentiation.FieldGradient;
33  import org.hipparchus.analysis.differentiation.Gradient;
34  import org.hipparchus.exception.LocalizedCoreFormats;
35  import org.hipparchus.geometry.euclidean.threed.FieldVector3D;
36  import org.hipparchus.geometry.euclidean.threed.Vector3D;
37  import org.hipparchus.util.FastMath;
38  import org.hipparchus.util.FieldSinCos;
39  import org.hipparchus.util.MathArrays;
40  import org.hipparchus.util.MathUtils;
41  import org.hipparchus.util.SinCos;
42  import org.orekit.attitudes.AttitudeProvider;
43  import org.orekit.errors.OrekitException;
44  import org.orekit.errors.OrekitInternalError;
45  import org.orekit.forces.gravity.potential.UnnormalizedSphericalHarmonicsProvider;
46  import org.orekit.forces.gravity.potential.UnnormalizedSphericalHarmonicsProvider.UnnormalizedSphericalHarmonics;
47  import org.orekit.frames.FieldTransform;
48  import org.orekit.frames.Frame;
49  import org.orekit.frames.Transform;
50  import org.orekit.orbits.FieldOrbit;
51  import org.orekit.orbits.Orbit;
52  import org.orekit.propagation.FieldSpacecraftState;
53  import org.orekit.propagation.PropagationType;
54  import org.orekit.propagation.SpacecraftState;
55  import org.orekit.propagation.events.EventDetector;
56  import org.orekit.propagation.events.FieldEventDetector;
57  import org.orekit.propagation.semianalytical.dsst.utilities.AuxiliaryElements;
58  import org.orekit.propagation.semianalytical.dsst.utilities.CoefficientsFactory;
59  import org.orekit.propagation.semianalytical.dsst.utilities.FieldAuxiliaryElements;
60  import org.orekit.propagation.semianalytical.dsst.utilities.FieldGHmsjPolynomials;
61  import org.orekit.propagation.semianalytical.dsst.utilities.FieldGammaMnsFunction;
62  import org.orekit.propagation.semianalytical.dsst.utilities.FieldShortPeriodicsInterpolatedCoefficient;
63  import org.orekit.propagation.semianalytical.dsst.utilities.GHmsjPolynomials;
64  import org.orekit.propagation.semianalytical.dsst.utilities.GammaMnsFunction;
65  import org.orekit.propagation.semianalytical.dsst.utilities.JacobiPolynomials;
66  import org.orekit.propagation.semianalytical.dsst.utilities.ShortPeriodicsInterpolatedCoefficient;
67  import org.orekit.propagation.semianalytical.dsst.utilities.hansen.FieldHansenTesseralLinear;
68  import org.orekit.propagation.semianalytical.dsst.utilities.hansen.HansenTesseralLinear;
69  import org.orekit.time.AbsoluteDate;
70  import org.orekit.time.FieldAbsoluteDate;
71  import org.orekit.utils.FieldTimeSpanMap;
72  import org.orekit.utils.ParameterDriver;
73  import org.orekit.utils.TimeSpanMap;
74  
75  /** Tesseral contribution to the central body gravitational perturbation.
76   *  <p>
77   *  Only resonant tesserals are considered.
78   *  </p>
79   *
80   *  @author Romain Di Costanzo
81   *  @author Pascal Parraud
82   *  @author Bryan Cazabonne (field translation)
83   */
84  public class DSSTTesseral implements DSSTForceModel {
85  
86      /**  Name of the prefix for short period coefficients keys. */
87      public static final String SHORT_PERIOD_PREFIX = "DSST-central-body-tesseral-";
88  
89      /** Identifier for cMm coefficients. */
90      public static final String CM_COEFFICIENTS = "cM";
91  
92      /** Identifier for sMm coefficients. */
93      public static final String SM_COEFFICIENTS = "sM";
94  
95      /** Retrograde factor I.
96       *  <p>
97       *  DSST model needs equinoctial orbit as internal representation.
98       *  Classical equinoctial elements have discontinuities when inclination
99       *  is close to zero. In this representation, I = +1. <br>
100      *  To avoid this discontinuity, another representation exists and equinoctial
101      *  elements can be expressed in a different way, called "retrograde" orbit.
102      *  This implies I = -1. <br>
103      *  As Orekit doesn't implement the retrograde orbit, I is always set to +1.
104      *  But for the sake of consistency with the theory, the retrograde factor
105      *  has been kept in the formulas.
106      *  </p>
107      */
108     private static final int I = 1;
109 
110     /** Central attraction scaling factor.
111      * <p>
112      * We use a power of 2 to avoid numeric noise introduction
113      * in the multiplications/divisions sequences.
114      * </p>
115      */
116     private static final double MU_SCALE = FastMath.scalb(1.0, 32);
117 
118     /** Minimum period for analytically averaged high-order resonant
119      *  central body spherical harmonics in seconds.
120      */
121     private static final double MIN_PERIOD_IN_SECONDS = 864000.;
122 
123     /** Minimum period for analytically averaged high-order resonant
124      *  central body spherical harmonics in satellite revolutions.
125      */
126     private static final double MIN_PERIOD_IN_SAT_REV = 10.;
127 
128     /** Number of points for interpolation. */
129     private static final int INTERPOLATION_POINTS = 3;
130 
131     /** Provider for spherical harmonics. */
132     private final UnnormalizedSphericalHarmonicsProvider provider;
133 
134     /** Central body rotating frame. */
135     private final Frame bodyFrame;
136 
137     /** Central body rotation rate (rad/s). */
138     private final double centralBodyRotationRate;
139 
140     /** Central body rotation period (seconds). */
141     private final double bodyPeriod;
142 
143     /** Maximal degree to consider for harmonics potential. */
144     private final int maxDegree;
145 
146     /** Maximal degree to consider for short periodics tesseral harmonics potential (without m-daily). */
147     private final int maxDegreeTesseralSP;
148 
149     /** Maximal degree to consider for short periodics m-daily tesseral harmonics potential . */
150     private final int maxDegreeMdailyTesseralSP;
151 
152     /** Maximal order to consider for harmonics potential. */
153     private final int maxOrder;
154 
155     /** Maximal order to consider for short periodics tesseral harmonics potential (without m-daily). */
156     private final int maxOrderTesseralSP;
157 
158     /** Maximal order to consider for short periodics m-daily tesseral harmonics potential . */
159     private final int maxOrderMdailyTesseralSP;
160 
161     /** Maximum power of the eccentricity to use in summation over s for
162      * short periodic tesseral harmonics (without m-daily). */
163     private final int maxEccPowTesseralSP;
164 
165     /** Maximum power of the eccentricity to use in summation over s for
166      * m-daily tesseral harmonics. */
167     private final int maxEccPowMdailyTesseralSP;
168 
169     /** Maximum value for j. */
170     private final int maxFrequencyShortPeriodics;
171 
172     /** Maximum power of the eccentricity to use in summation over s. */
173     private int maxEccPow;
174 
175     /** Maximum power of the eccentricity to use in Hansen coefficient Kernel expansion. */
176     private int maxHansen;
177 
178     /** Maximum value between maxOrderMdailyTesseralSP and maxOrderTesseralSP. */
179     private int mMax;
180 
181     /** List of non resonant orders with j != 0. */
182     private final SortedMap<Integer, List<Integer> > nonResOrders;
183 
184     /** List of resonant orders. */
185     private final List<Integer> resOrders;
186 
187     /** Short period terms. */
188     private TesseralShortPeriodicCoefficients shortPeriodTerms;
189 
190     /** Short period terms. */
191     private Map<Field<?>, FieldTesseralShortPeriodicCoefficients<?>> fieldShortPeriodTerms;
192 
193     /** Driver for gravitational parameter. */
194     private final ParameterDriver gmParameterDriver;
195 
196     /** Hansen objects. */
197     private HansenObjects hansen;
198 
199     /** Hansen objects for field elements. */
200     private Map<Field<?>, FieldHansenObjects<?>> fieldHansen;
201 
202     /** Flag for force model initialization with field elements. */
203     private boolean pendingInitialization;
204 
205     /** Simple constructor with default reference values.
206      * <p>
207      * When this constructor is used, maximum allowed values are used
208      * for the short periodic coefficients:
209      * </p>
210      * <ul>
211      *    <li> {@link #maxDegreeTesseralSP} is set to {@code provider.getMaxDegree()} </li>
212      *    <li> {@link #maxOrderTesseralSP} is set to {@code provider.getMaxOrder()}. </li>
213      *    <li> {@link #maxEccPowTesseralSP} is set to 4 </li>
214      *    <li> {@link #maxFrequencyShortPeriodics} is set to {@code min(provider.getMaxDegree() + 4, 12)}.
215      *         This parameter should not exceed 12 as higher values will exceed computer capacity </li>
216      *    <li> {@link #maxDegreeMdailyTesseralSP} is set to {@code provider.getMaxDegree()} </li>
217      *    <li> {@link #maxOrderMdailyTesseralSP} is set to {@code provider.getMaxOrder()} </li>
218      *    <li> {@link #maxEccPowMdailyTesseralSP} is set to min(provider.getMaxDegree() - 2, 4).
219      *         This parameter should not exceed 4 as higher values will exceed computer capacity </li>
220      * </ul>
221      * @param centralBodyFrame rotating body frame
222      * @param centralBodyRotationRate central body rotation rate (rad/s)
223      * @param provider provider for spherical harmonics
224      * @since 10.1
225      */
226     public DSSTTesseral(final Frame centralBodyFrame,
227                         final double centralBodyRotationRate,
228                         final UnnormalizedSphericalHarmonicsProvider provider) {
229         this(centralBodyFrame, centralBodyRotationRate, provider, provider.getMaxDegree(),
230              provider.getMaxOrder(), 4,  FastMath.min(12, provider.getMaxDegree() + 4),
231              provider.getMaxDegree(), provider.getMaxOrder(), FastMath.min(4, provider.getMaxDegree() - 2));
232     }
233 
234     /** Simple constructor.
235      * @param centralBodyFrame rotating body frame
236      * @param centralBodyRotationRate central body rotation rate (rad/s)
237      * @param provider provider for spherical harmonics
238      * @param maxDegreeTesseralSP maximal degree to consider for short periodics tesseral harmonics potential
239      *  (must be between 2 and {@code provider.getMaxDegree()})
240      * @param maxOrderTesseralSP maximal order to consider for short periodics tesseral harmonics potential
241      *  (must be between 0 and {@code provider.getMaxOrder()})
242      * @param maxEccPowTesseralSP maximum power of the eccentricity to use in summation over s for
243      * short periodic tesseral harmonics (without m-daily), should typically not exceed 4 as higher
244      * values will exceed computer capacity
245      * @param maxFrequencyShortPeriodics maximum frequency in mean longitude for short periodic computations
246      * (typically {@code maxDegreeTesseralSP} + {@code maxEccPowTesseralSP and no more than 12})
247      * @param maxDegreeMdailyTesseralSP maximal degree to consider for short periodics m-daily tesseral harmonics potential
248      *  (must be between 2 and {@code provider.getMaxDegree()})
249      * @param maxOrderMdailyTesseralSP maximal order to consider for short periodics m-daily tesseral harmonics potential
250      *  (must be between 0 and {@code provider.getMaxOrder()})
251      * @param maxEccPowMdailyTesseralSP maximum power of the eccentricity to use in summation over s for
252      * m-daily tesseral harmonics, (must be between 0 and {@code maxDegreeMdailyTesseralSP - 2},
253      * but should typically not exceed 4 as higher values will exceed computer capacity)
254      * @since 7.2
255      */
256     public DSSTTesseral(final Frame centralBodyFrame,
257                         final double centralBodyRotationRate,
258                         final UnnormalizedSphericalHarmonicsProvider provider,
259                         final int maxDegreeTesseralSP, final int maxOrderTesseralSP,
260                         final int maxEccPowTesseralSP, final int maxFrequencyShortPeriodics,
261                         final int maxDegreeMdailyTesseralSP, final int maxOrderMdailyTesseralSP,
262                         final int maxEccPowMdailyTesseralSP) {
263 
264         gmParameterDriver = new ParameterDriver(DSSTNewtonianAttraction.CENTRAL_ATTRACTION_COEFFICIENT,
265                                                 provider.getMu(), MU_SCALE,
266                                                 0.0, Double.POSITIVE_INFINITY);
267 
268         // Central body rotating frame
269         this.bodyFrame = centralBodyFrame;
270 
271         //Save the rotation rate
272         this.centralBodyRotationRate = centralBodyRotationRate;
273 
274         // Central body rotation period in seconds
275         this.bodyPeriod = MathUtils.TWO_PI / centralBodyRotationRate;
276 
277         // Provider for spherical harmonics
278         this.provider      = provider;
279         this.maxDegree     = provider.getMaxDegree();
280         this.maxOrder      = provider.getMaxOrder();
281 
282         //set the maximum degree order for short periodics
283         checkIndexRange(maxDegreeTesseralSP, 2, maxDegree);
284         this.maxDegreeTesseralSP       = maxDegreeTesseralSP;
285 
286         checkIndexRange(maxDegreeMdailyTesseralSP, 2, maxDegree);
287         this.maxDegreeMdailyTesseralSP = maxDegreeMdailyTesseralSP;
288 
289         checkIndexRange(maxOrderTesseralSP, 0, maxOrder);
290         this.maxOrderTesseralSP        = maxOrderTesseralSP;
291 
292         checkIndexRange(maxOrderMdailyTesseralSP, 0, maxOrder);
293         this.maxOrderMdailyTesseralSP  = maxOrderMdailyTesseralSP;
294 
295         // set the maximum value for eccentricity power
296         this.maxEccPowTesseralSP       = maxEccPowTesseralSP;
297 
298         checkIndexRange(maxEccPowMdailyTesseralSP, 0, maxDegreeMdailyTesseralSP - 2);
299         this.maxEccPowMdailyTesseralSP = maxEccPowMdailyTesseralSP;
300 
301         // set the maximum value for frequency
302         this.maxFrequencyShortPeriodics = maxFrequencyShortPeriodics;
303 
304         // Initialize default values
305         this.resOrders    = new ArrayList<Integer>();
306         this.nonResOrders = new TreeMap<Integer, List <Integer>>();
307 
308         pendingInitialization = true;
309 
310         // Initialize default values
311         this.fieldShortPeriodTerms = new HashMap<>();
312         this.fieldHansen           = new HashMap<>();
313         this.maxEccPow             = 0;
314         this.maxHansen             = 0;
315 
316     }
317 
318     /** Check an index range.
319      * @param index index value
320      * @param min minimum value for index
321      * @param max maximum value for index
322      */
323     private void checkIndexRange(final int index, final int min, final int max) {
324         if (index < min || index > max) {
325             throw new OrekitException(LocalizedCoreFormats.OUT_OF_RANGE_SIMPLE, index, min, max);
326         }
327     }
328 
329     /** {@inheritDoc} */
330     @Override
331     public List<ShortPeriodTerms> initializeShortPeriodTerms(final AuxiliaryElements auxiliaryElements,
332                                              final PropagationType type,
333                                              final double[] parameters) {
334 
335         // Initializes specific parameters.
336         final DSSTTesseralContext context = initializeStep(auxiliaryElements, parameters);
337 
338         // Set the highest power of the eccentricity in the analytical power
339         // series expansion for the averaged high order resonant central body
340         // spherical harmonic perturbation
341         maxEccPow = getMaxEccPow(auxiliaryElements.getEcc());
342 
343         // Set the maximum power of the eccentricity to use in Hansen coefficient Kernel expansion.
344         maxHansen = maxEccPow / 2;
345 
346         // The following terms are only used for hansen objects initialization
347         final double ratio = context.getRatio();
348 
349         // Compute the non resonant tesseral harmonic terms if not set by the user
350         getResonantAndNonResonantTerms(type, context.getOrbitPeriod(), ratio);
351 
352         hansen = new HansenObjects(ratio, type);
353 
354         mMax = FastMath.max(maxOrderTesseralSP, maxOrderMdailyTesseralSP);
355 
356         shortPeriodTerms = new TesseralShortPeriodicCoefficients(bodyFrame, maxOrderMdailyTesseralSP,
357                                                                  maxDegreeTesseralSP < 0, nonResOrders,
358                                                                  mMax, maxFrequencyShortPeriodics, INTERPOLATION_POINTS,
359                                                                  new TimeSpanMap<Slot>(new Slot(mMax, maxFrequencyShortPeriodics,
360                                                                                                 INTERPOLATION_POINTS)));
361 
362         final List<ShortPeriodTerms> list = new ArrayList<ShortPeriodTerms>();
363         list.add(shortPeriodTerms);
364         return list;
365 
366     }
367 
368     /** {@inheritDoc} */
369     @Override
370     public <T extends CalculusFieldElement<T>> List<FieldShortPeriodTerms<T>> initializeShortPeriodTerms(final FieldAuxiliaryElements<T> auxiliaryElements,
371                                                                                      final PropagationType type,
372                                                                                      final T[] parameters) {
373 
374         // Field used by default
375         final Field<T> field = auxiliaryElements.getDate().getField();
376 
377         if (pendingInitialization == true) {
378 
379             // Initializes specific parameters.
380             final FieldDSSTTesseralContext<T> context = initializeStep(auxiliaryElements, parameters);
381 
382             // Set the highest power of the eccentricity in the analytical power
383             // series expansion for the averaged high order resonant central body
384             // spherical harmonic perturbation
385             maxEccPow = getMaxEccPow(auxiliaryElements.getEcc().getReal());
386 
387             // Set the maximum power of the eccentricity to use in Hansen coefficient Kernel expansion.
388             maxHansen = maxEccPow / 2;
389 
390             // The following terms are only used for hansen objects initialization
391             final T ratio = context.getRatio();
392 
393             // Compute the non resonant tesseral harmonic terms if not set by the user
394             // Field information is not important here
395             getResonantAndNonResonantTerms(type, context.getOrbitPeriod().getReal(), ratio.getReal());
396 
397             mMax = FastMath.max(maxOrderTesseralSP, maxOrderMdailyTesseralSP);
398 
399             fieldHansen.put(field, new FieldHansenObjects<>(ratio, type));
400 
401             pendingInitialization = false;
402         }
403 
404         final FieldTesseralShortPeriodicCoefficients<T> ftspc =
405                         new FieldTesseralShortPeriodicCoefficients<>(bodyFrame, maxOrderMdailyTesseralSP,
406                                                                      maxDegreeTesseralSP < 0, nonResOrders,
407                                                                      mMax, maxFrequencyShortPeriodics, INTERPOLATION_POINTS,
408                                                                      new FieldTimeSpanMap<>(new FieldSlot<>(mMax,
409                                                                                                             maxFrequencyShortPeriodics,
410                                                                                                             INTERPOLATION_POINTS),
411                                                                                             field));
412 
413         fieldShortPeriodTerms.put(field, ftspc);
414         return Collections.singletonList(ftspc);
415 
416     }
417 
418     /**
419      * Get the maximum power of the eccentricity to use in summation over s.
420      * @param e eccentricity
421      * @return the maximum power of the eccentricity
422      */
423     private int getMaxEccPow(final double e) {
424         // maxEccPow depends on satellite eccentricity
425         if (e <= 0.005) {
426             return 3;
427         } else if (e <= 0.02) {
428             return 4;
429         } else if (e <= 0.1) {
430             return 7;
431         } else if (e <= 0.2) {
432             return 10;
433         } else if (e <= 0.3) {
434             return 12;
435         } else if (e <= 0.4) {
436             return 15;
437         } else {
438             return 20;
439         }
440     }
441 
442     /** Performs initialization at each integration step for the current force model.
443      *  <p>
444      *  This method aims at being called before mean elements rates computation.
445      *  </p>
446      *  @param auxiliaryElements auxiliary elements related to the current orbit
447      *  @param parameters values of the force model parameters
448      *  @return new force model context
449      */
450     private DSSTTesseralContext initializeStep(final AuxiliaryElements auxiliaryElements, final double[] parameters) {
451         return new DSSTTesseralContext(auxiliaryElements, bodyFrame, provider, maxFrequencyShortPeriodics, bodyPeriod, parameters);
452     }
453 
454     /** Performs initialization at each integration step for the current force model.
455      *  <p>
456      *  This method aims at being called before mean elements rates computation.
457      *  </p>
458      *  @param <T> type of the elements
459      *  @param auxiliaryElements auxiliary elements related to the current orbit
460      *  @param parameters values of the force model parameters
461      *  @return new force model context
462      */
463     private <T extends CalculusFieldElement<T>> FieldDSSTTesseralContext<T> initializeStep(final FieldAuxiliaryElements<T> auxiliaryElements,
464                                                                                        final T[] parameters) {
465         return new FieldDSSTTesseralContext<>(auxiliaryElements, bodyFrame, provider, maxFrequencyShortPeriodics, bodyPeriod, parameters);
466     }
467 
468     /** {@inheritDoc} */
469     @Override
470     public double[] getMeanElementRate(final SpacecraftState spacecraftState,
471                                        final AuxiliaryElements auxiliaryElements, final double[] parameters) {
472 
473         // Container for attributes
474         final DSSTTesseralContext context = initializeStep(auxiliaryElements, parameters);
475 
476         // Access to potential U derivatives
477         final UAnddU udu = new UAnddU(spacecraftState.getDate(), context, hansen);
478 
479         // Compute the cross derivative operator :
480         final double UAlphaGamma   = auxiliaryElements.getAlpha() * udu.getdUdGa() - auxiliaryElements.getGamma() * udu.getdUdAl();
481         final double UAlphaBeta    = auxiliaryElements.getAlpha() * udu.getdUdBe() - auxiliaryElements.getBeta()  * udu.getdUdAl();
482         final double UBetaGamma    = auxiliaryElements.getBeta() * udu.getdUdGa() - auxiliaryElements.getGamma() * udu.getdUdBe();
483         final double Uhk           = auxiliaryElements.getH() * udu.getdUdk()  - auxiliaryElements.getK() * udu.getdUdh();
484         final double pUagmIqUbgoAB = (auxiliaryElements.getP() * UAlphaGamma - I * auxiliaryElements.getQ() * UBetaGamma) * context.getOoAB();
485         final double UhkmUabmdUdl  = Uhk - UAlphaBeta - udu.getdUdl();
486 
487         final double da =  context.getAx2oA() * udu.getdUdl();
488         final double dh =  context.getBoA() * udu.getdUdk() + auxiliaryElements.getK() * pUagmIqUbgoAB - auxiliaryElements.getH() * context.getBoABpo() * udu.getdUdl();
489         final double dk =  -(context.getBoA() * udu.getdUdh() + auxiliaryElements.getH() * pUagmIqUbgoAB + auxiliaryElements.getK() * context.getBoABpo() * udu.getdUdl());
490         final double dp =  context.getCo2AB() * (auxiliaryElements.getP() * UhkmUabmdUdl - UBetaGamma);
491         final double dq =  context.getCo2AB() * (auxiliaryElements.getQ() * UhkmUabmdUdl - I * UAlphaGamma);
492         final double dM = -context.getAx2oA() * udu.getdUda() + context.getBoABpo() * (auxiliaryElements.getH() * udu.getdUdh() + auxiliaryElements.getK() * udu.getdUdk()) + pUagmIqUbgoAB;
493 
494         return new double[] {da, dk, dh, dq, dp, dM};
495     }
496 
497     /** {@inheritDoc} */
498     @Override
499     public <T extends CalculusFieldElement<T>> T[] getMeanElementRate(final FieldSpacecraftState<T> spacecraftState,
500                                                                   final FieldAuxiliaryElements<T> auxiliaryElements,
501                                                                   final T[] parameters) {
502 
503         // Field used by default
504         final Field<T> field = auxiliaryElements.getDate().getField();
505 
506         // Container for attributes
507         final FieldDSSTTesseralContext<T> context = initializeStep(auxiliaryElements, parameters);
508 
509         @SuppressWarnings("unchecked")
510         final FieldHansenObjects<T> fho = (FieldHansenObjects<T>) fieldHansen.get(field);
511         // Access to potential U derivatives
512         final FieldUAnddU<T> udu = new FieldUAnddU<>(spacecraftState.getDate(), context, fho);
513 
514         // Compute the cross derivative operator :
515         final T UAlphaGamma   = udu.getdUdGa().multiply(auxiliaryElements.getAlpha()).subtract(udu.getdUdAl().multiply(auxiliaryElements.getGamma()));
516         final T UAlphaBeta    = udu.getdUdBe().multiply(auxiliaryElements.getAlpha()).subtract(udu.getdUdAl().multiply(auxiliaryElements.getBeta()));
517         final T UBetaGamma    = udu.getdUdGa().multiply(auxiliaryElements.getBeta()).subtract(udu.getdUdBe().multiply(auxiliaryElements.getGamma()));
518         final T Uhk           = udu.getdUdk().multiply(auxiliaryElements.getH()).subtract(udu.getdUdh().multiply(auxiliaryElements.getK()));
519         final T pUagmIqUbgoAB = (UAlphaGamma.multiply(auxiliaryElements.getP()).subtract(UBetaGamma.multiply(auxiliaryElements.getQ()).multiply(I))).multiply(context.getOoAB());
520         final T UhkmUabmdUdl  = Uhk.subtract(UAlphaBeta).subtract(udu.getdUdl());
521 
522         final T da = udu.getdUdl().multiply(context.getAx2oA());
523         final T dh = udu.getdUdk().multiply(context.getBoA()).add(pUagmIqUbgoAB.multiply(auxiliaryElements.getK())).subtract(udu.getdUdl().multiply(auxiliaryElements.getH()).multiply(context.getBoABpo()));
524         final T dk = (udu.getdUdh().multiply(context.getBoA()).add(pUagmIqUbgoAB.multiply(auxiliaryElements.getH())).add(udu.getdUdl().multiply(context.getBoABpo()).multiply(auxiliaryElements.getK()))).negate();
525         final T dp = context.getCo2AB().multiply(auxiliaryElements.getP().multiply(UhkmUabmdUdl).subtract(UBetaGamma));
526         final T dq = context.getCo2AB().multiply(auxiliaryElements.getQ().multiply(UhkmUabmdUdl).subtract(UAlphaGamma.multiply(I)));
527         final T dM = pUagmIqUbgoAB.add(udu.getdUda().multiply(context.getAx2oA()).negate()).add((udu.getdUdh().multiply(auxiliaryElements.getH()).add(udu.getdUdk().multiply(auxiliaryElements.getK()))).multiply(context.getBoABpo()));
528 
529         final T[] elements = MathArrays.buildArray(field, 6);
530         elements[0] = da;
531         elements[1] = dk;
532         elements[2] = dh;
533         elements[3] = dq;
534         elements[4] = dp;
535         elements[5] = dM;
536 
537         return elements;
538 
539     }
540 
541     /** {@inheritDoc} */
542     @Override
543     public void updateShortPeriodTerms(final double[] parameters, final SpacecraftState... meanStates) {
544 
545         final Slot slot = shortPeriodTerms.createSlot(meanStates);
546 
547         for (final SpacecraftState meanState : meanStates) {
548 
549             final AuxiliaryElements auxiliaryElements = new AuxiliaryElements(meanState.getOrbit(), I);
550 
551             final DSSTTesseralContext context = initializeStep(auxiliaryElements, parameters);
552 
553             // Initialise the Hansen coefficients
554             for (int s = -maxDegree; s <= maxDegree; s++) {
555                 // coefficients with j == 0 are always needed
556                 hansen.computeHansenObjectsInitValues(context, s + maxDegree, 0);
557                 if (maxDegreeTesseralSP >= 0) {
558                     // initialize other objects only if required
559                     for (int j = 1; j <= maxFrequencyShortPeriodics; j++) {
560                         hansen.computeHansenObjectsInitValues(context, s + maxDegree, j);
561                     }
562                 }
563             }
564 
565             final FourierCjSjCoefficients cjsjFourier = new FourierCjSjCoefficients(maxFrequencyShortPeriodics, mMax);
566 
567             // Compute coefficients
568             // Compute only if there is at least one non-resonant tesseral
569             if (!nonResOrders.isEmpty() || maxDegreeTesseralSP < 0) {
570                 // Generate the fourrier coefficients
571                 cjsjFourier.generateCoefficients(meanState.getDate(), context, hansen);
572 
573                 // the coefficient 3n / 2a
574                 final double tnota = 1.5 * context.getMeanMotion() / auxiliaryElements.getSma();
575 
576                 // build the mDaily coefficients
577                 for (int m = 1; m <= maxOrderMdailyTesseralSP; m++) {
578                     // build the coefficients
579                     buildCoefficients(cjsjFourier, meanState.getDate(), slot, m, 0, tnota, context);
580                 }
581 
582                 if (maxDegreeTesseralSP >= 0) {
583                     // generate the other coefficients, if required
584                     for (final Map.Entry<Integer, List<Integer>> entry : nonResOrders.entrySet()) {
585 
586                         for (int j : entry.getValue()) {
587                             // build the coefficients
588                             buildCoefficients(cjsjFourier, meanState.getDate(), slot, entry.getKey(), j, tnota, context);
589                         }
590                     }
591                 }
592             }
593 
594         }
595 
596     }
597 
598     /** {@inheritDoc} */
599     @Override
600     @SuppressWarnings("unchecked")
601     public <T extends CalculusFieldElement<T>> void updateShortPeriodTerms(final T[] parameters,
602                                                                        final FieldSpacecraftState<T>... meanStates) {
603 
604         // Field used by default
605         final Field<T> field = meanStates[0].getDate().getField();
606 
607         final FieldTesseralShortPeriodicCoefficients<T> ftspc = (FieldTesseralShortPeriodicCoefficients<T>) fieldShortPeriodTerms.get(field);
608         final FieldSlot<T> slot = ftspc.createSlot(meanStates);
609 
610         for (final FieldSpacecraftState<T> meanState : meanStates) {
611 
612             final FieldAuxiliaryElements<T> auxiliaryElements = new FieldAuxiliaryElements<>(meanState.getOrbit(), I);
613 
614             final FieldDSSTTesseralContext<T> context = initializeStep(auxiliaryElements, parameters);
615 
616             final FieldHansenObjects<T> fho = (FieldHansenObjects<T>) fieldHansen.get(field);
617             // Initialise the Hansen coefficients
618             for (int s = -maxDegree; s <= maxDegree; s++) {
619                 // coefficients with j == 0 are always needed
620                 fho.computeHansenObjectsInitValues(context, s + maxDegree, 0);
621                 if (maxDegreeTesseralSP >= 0) {
622                     // initialize other objects only if required
623                     for (int j = 1; j <= maxFrequencyShortPeriodics; j++) {
624                         fho.computeHansenObjectsInitValues(context, s + maxDegree, j);
625                     }
626                 }
627             }
628 
629             final FieldFourierCjSjCoefficients<T> cjsjFourier = new FieldFourierCjSjCoefficients<>(maxFrequencyShortPeriodics, mMax, field);
630 
631             // Compute coefficients
632             // Compute only if there is at least one non-resonant tesseral
633             if (!nonResOrders.isEmpty() || maxDegreeTesseralSP < 0) {
634                 // Generate the fourrier coefficients
635                 cjsjFourier.generateCoefficients(meanState.getDate(), context, fho, field);
636 
637                 // the coefficient 3n / 2a
638                 final T tnota = context.getMeanMotion().multiply(1.5).divide(auxiliaryElements.getSma());
639 
640                 // build the mDaily coefficients
641                 for (int m = 1; m <= maxOrderMdailyTesseralSP; m++) {
642                     // build the coefficients
643                     buildCoefficients(cjsjFourier, meanState.getDate(), slot, m, 0, tnota, context, field);
644                 }
645 
646                 if (maxDegreeTesseralSP >= 0) {
647                     // generate the other coefficients, if required
648                     for (final Map.Entry<Integer, List<Integer>> entry : nonResOrders.entrySet()) {
649 
650                         for (int j : entry.getValue()) {
651                             // build the coefficients
652                             buildCoefficients(cjsjFourier, meanState.getDate(), slot, entry.getKey(), j, tnota, context, field);
653                         }
654                     }
655                 }
656             }
657 
658         }
659 
660     }
661 
662     /** {@inheritDoc} */
663     public List<ParameterDriver> getParametersDrivers() {
664         return Collections.singletonList(gmParameterDriver);
665     }
666 
667     /** Build a set of coefficients.
668      * @param cjsjFourier the fourier coefficients C<sub>i</sub><sup>j</sup> and the S<sub>i</sub><sup>j</sup>
669      * @param date the current date
670      * @param slot slot to which the coefficients belong
671      * @param m m index
672      * @param j j index
673      * @param tnota 3n/2a
674      * @param context container for attributes
675      */
676     private void buildCoefficients(final FourierCjSjCoefficients cjsjFourier,
677                                    final AbsoluteDate date, final Slot slot,
678                                    final int m, final int j, final double tnota, final DSSTTesseralContext context) {
679 
680         // Create local arrays
681         final double[] currentCijm = new double[] {0., 0., 0., 0., 0., 0.};
682         final double[] currentSijm = new double[] {0., 0., 0., 0., 0., 0.};
683 
684         // compute the term 1 / (jn - mθ<sup>.</sup>)
685         final double oojnmt = 1. / (j * context.getMeanMotion() - m * centralBodyRotationRate);
686 
687         // initialise the coeficients
688         for (int i = 0; i < 6; i++) {
689             currentCijm[i] = -cjsjFourier.getSijm(i, j, m);
690             currentSijm[i] = cjsjFourier.getCijm(i, j, m);
691         }
692         // Add the separate part for δ<sub>6i</sub>
693         currentCijm[5] += tnota * oojnmt * cjsjFourier.getCijm(0, j, m);
694         currentSijm[5] += tnota * oojnmt * cjsjFourier.getSijm(0, j, m);
695 
696         //Multiply by 1 / (jn - mθ<sup>.</sup>)
697         for (int i = 0; i < 6; i++) {
698             currentCijm[i] *= oojnmt;
699             currentSijm[i] *= oojnmt;
700         }
701 
702         // Add the coefficients to the interpolation grid
703         slot.cijm[m][j + maxFrequencyShortPeriodics].addGridPoint(date, currentCijm);
704         slot.sijm[m][j + maxFrequencyShortPeriodics].addGridPoint(date, currentSijm);
705 
706     }
707 
708      /** Build a set of coefficients.
709      * @param <T> the type of the field elements
710      * @param cjsjFourier the fourier coefficients C<sub>i</sub><sup>j</sup> and the S<sub>i</sub><sup>j</sup>
711      * @param date the current date
712      * @param slot slot to which the coefficients belong
713      * @param m m index
714      * @param j j index
715      * @param tnota 3n/2a
716      * @param context container for attributes
717      * @param field field used by default
718      */
719     private <T extends CalculusFieldElement<T>> void buildCoefficients(final FieldFourierCjSjCoefficients<T> cjsjFourier,
720                                                                    final FieldAbsoluteDate<T> date,
721                                                                    final FieldSlot<T> slot,
722                                                                    final int m, final int j, final T tnota,
723                                                                    final FieldDSSTTesseralContext<T> context,
724                                                                    final Field<T> field) {
725 
726         // Zero
727         final T zero = field.getZero();
728 
729         // Create local arrays
730         final T[] currentCijm = MathArrays.buildArray(field, 6);
731         final T[] currentSijm = MathArrays.buildArray(field, 6);
732 
733         Arrays.fill(currentCijm, zero);
734         Arrays.fill(currentSijm, zero);
735 
736         // compute the term 1 / (jn - mθ<sup>.</sup>)
737         final T oojnmt = (context.getMeanMotion().multiply(j).subtract(m * centralBodyRotationRate)).reciprocal();
738 
739         // initialise the coeficients
740         for (int i = 0; i < 6; i++) {
741             currentCijm[i] = cjsjFourier.getSijm(i, j, m).negate();
742             currentSijm[i] = cjsjFourier.getCijm(i, j, m);
743         }
744         // Add the separate part for δ<sub>6i</sub>
745         currentCijm[5] = currentCijm[5].add(tnota.multiply(oojnmt).multiply(cjsjFourier.getCijm(0, j, m)));
746         currentSijm[5] = currentSijm[5].add(tnota.multiply(oojnmt).multiply(cjsjFourier.getSijm(0, j, m)));
747 
748         //Multiply by 1 / (jn - mθ<sup>.</sup>)
749         for (int i = 0; i < 6; i++) {
750             currentCijm[i] = currentCijm[i].multiply(oojnmt);
751             currentSijm[i] = currentSijm[i].multiply(oojnmt);
752         }
753 
754         // Add the coefficients to the interpolation grid
755         slot.cijm[m][j + maxFrequencyShortPeriodics].addGridPoint(date, currentCijm);
756         slot.sijm[m][j + maxFrequencyShortPeriodics].addGridPoint(date, currentSijm);
757 
758     }
759 
760     /** {@inheritDoc} */
761     @Override
762     public EventDetector[] getEventsDetectors() {
763         return null;
764     }
765 
766     /** {@inheritDoc} */
767     @Override
768     public <T extends CalculusFieldElement<T>> FieldEventDetector<T>[] getFieldEventsDetectors(final Field<T> field) {
769         return null;
770     }
771 
772      /**
773       * Get the resonant and non-resonant tesseral terms in the central body spherical harmonic field.
774       *
775       * @param type type of the elements used during the propagation
776       * @param orbitPeriod Keplerian period
777       * @param ratio ratio of satellite period to central body rotation period
778       */
779     private void getResonantAndNonResonantTerms(final PropagationType type, final double orbitPeriod,
780                                                 final double ratio) {
781 
782         // Compute natural resonant terms
783         final double tolerance = 1. / FastMath.max(MIN_PERIOD_IN_SAT_REV,
784                                                    MIN_PERIOD_IN_SECONDS / orbitPeriod);
785 
786         // Search the resonant orders in the tesseral harmonic field
787         resOrders.clear();
788         nonResOrders.clear();
789         for (int m = 1; m <= maxOrder; m++) {
790             final double resonance = ratio * m;
791             int jRes = 0;
792             final int jComputedRes = (int) FastMath.round(resonance);
793             if (jComputedRes > 0 && jComputedRes <= maxFrequencyShortPeriodics && FastMath.abs(resonance - jComputedRes) <= tolerance) {
794                 // Store each resonant index and order
795                 resOrders.add(m);
796                 jRes = jComputedRes;
797             }
798 
799             if (type == PropagationType.OSCULATING && maxDegreeTesseralSP >= 0 && m <= maxOrderTesseralSP) {
800                 //compute non resonant orders in the tesseral harmonic field
801                 final List<Integer> listJofM = new ArrayList<Integer>();
802                 //for the moment we take only the pairs (j,m) with |j| <= maxDegree + maxEccPow (from |s-j| <= maxEccPow and |s| <= maxDegree)
803                 for (int j = -maxFrequencyShortPeriodics; j <= maxFrequencyShortPeriodics; j++) {
804                     if (j != 0 && j != jRes) {
805                         listJofM.add(j);
806                     }
807                 }
808 
809                 nonResOrders.put(m, listJofM);
810             }
811         }
812     }
813 
814     /** Compute the n-SUM for potential derivatives components.
815      *  @param date current date
816      *  @param j resonant index <i>j</i>
817      *  @param m resonant order <i>m</i>
818      *  @param s d'Alembert characteristic <i>s</i>
819      *  @param maxN maximum possible value for <i>n</i> index
820      *  @param roaPow powers of R/a up to degree <i>n</i>
821      *  @param ghMSJ G<sup>j</sup><sub>m,s</sub> and H<sup>j</sup><sub>m,s</sub> polynomials
822      *  @param gammaMNS &Gamma;<sup>m</sup><sub>n,s</sub>(γ) function
823      *  @param context container for attributes
824      *  @param hansenObjects initialization of hansen objects
825      *  @return Components of U<sub>n</sub> derivatives for fixed j, m, s
826      */
827     private double[][] computeNSum(final AbsoluteDate date,
828                                    final int j, final int m, final int s, final int maxN, final double[] roaPow,
829                                    final GHmsjPolynomials ghMSJ, final GammaMnsFunction gammaMNS, final DSSTTesseralContext context,
830                                    final HansenObjects hansenObjects) {
831 
832         // Auxiliary elements related to the current orbit
833         final AuxiliaryElements auxiliaryElements = context.getAuxiliaryElements();
834 
835         //spherical harmonics
836         final UnnormalizedSphericalHarmonics harmonics = provider.onDate(date);
837 
838         // Potential derivatives components
839         double dUdaCos  = 0.;
840         double dUdaSin  = 0.;
841         double dUdhCos  = 0.;
842         double dUdhSin  = 0.;
843         double dUdkCos  = 0.;
844         double dUdkSin  = 0.;
845         double dUdlCos  = 0.;
846         double dUdlSin  = 0.;
847         double dUdAlCos = 0.;
848         double dUdAlSin = 0.;
849         double dUdBeCos = 0.;
850         double dUdBeSin = 0.;
851         double dUdGaCos = 0.;
852         double dUdGaSin = 0.;
853 
854         // I^m
855         @SuppressWarnings("unused")
856         final int Im = I > 0 ? 1 : (m % 2 == 0 ? 1 : -1);
857 
858         // jacobi v, w, indices from 2.7.1-(15)
859         final int v = FastMath.abs(m - s);
860         final int w = FastMath.abs(m + s);
861 
862         // Initialise lower degree nmin = (Max(2, m, |s|)) for summation over n
863         final int nmin = FastMath.max(FastMath.max(2, m), FastMath.abs(s));
864 
865         //Get the corresponding Hansen object
866         final int sIndex = maxDegree + (j < 0 ? -s : s);
867         final int jIndex = FastMath.abs(j);
868         final HansenTesseralLinear hans = hansenObjects.getHansenObjects()[sIndex][jIndex];
869 
870         // n-SUM from nmin to N
871         for (int n = nmin; n <= maxN; n++) {
872             // If (n - s) is odd, the contribution is null because of Vmns
873             if ((n - s) % 2 == 0) {
874 
875                 // Vmns coefficient
876                 final double vMNS   = CoefficientsFactory.getVmns(m, n, s);
877 
878                 // Inclination function Gamma and derivative
879                 final double gaMNS  = gammaMNS.getValue(m, n, s);
880                 final double dGaMNS = gammaMNS.getDerivative(m, n, s);
881 
882                 // Hansen kernel value and derivative
883                 final double kJNS   = hans.getValue(-n - 1, context.getChi());
884                 final double dkJNS  = hans.getDerivative(-n - 1, context.getChi());
885 
886                 // Gjms, Hjms polynomials and derivatives
887                 final double gMSJ   = ghMSJ.getGmsj(m, s, j);
888                 final double hMSJ   = ghMSJ.getHmsj(m, s, j);
889                 final double dGdh   = ghMSJ.getdGmsdh(m, s, j);
890                 final double dGdk   = ghMSJ.getdGmsdk(m, s, j);
891                 final double dGdA   = ghMSJ.getdGmsdAlpha(m, s, j);
892                 final double dGdB   = ghMSJ.getdGmsdBeta(m, s, j);
893                 final double dHdh   = ghMSJ.getdHmsdh(m, s, j);
894                 final double dHdk   = ghMSJ.getdHmsdk(m, s, j);
895                 final double dHdA   = ghMSJ.getdHmsdAlpha(m, s, j);
896                 final double dHdB   = ghMSJ.getdHmsdBeta(m, s, j);
897 
898                 // Jacobi l-index from 2.7.1-(15)
899                 final int l = FastMath.min(n - m, n - FastMath.abs(s));
900                 // Jacobi polynomial and derivative
901                 final Gradient jacobi =
902                         JacobiPolynomials.getValue(l, v, w, Gradient.variable(1, 0, auxiliaryElements.getGamma()));
903 
904                 // Geopotential coefficients
905                 final double cnm = harmonics.getUnnormalizedCnm(n, m);
906                 final double snm = harmonics.getUnnormalizedSnm(n, m);
907 
908                 // Common factors from expansion of equations 3.3-4
909                 final double cf_0      = roaPow[n] * Im * vMNS;
910                 final double cf_1      = cf_0 * gaMNS * jacobi.getValue();
911                 final double cf_2      = cf_1 * kJNS;
912                 final double gcPhs     = gMSJ * cnm + hMSJ * snm;
913                 final double gsMhc     = gMSJ * snm - hMSJ * cnm;
914                 final double dKgcPhsx2 = 2. * dkJNS * gcPhs;
915                 final double dKgsMhcx2 = 2. * dkJNS * gsMhc;
916                 final double dUdaCoef  = (n + 1) * cf_2;
917                 final double dUdlCoef  = j * cf_2;
918                 final double dUdGaCoef = cf_0 * kJNS * (jacobi.getValue() * dGaMNS + gaMNS * jacobi.getGradient()[0]);
919 
920                 // dU / da components
921                 dUdaCos  += dUdaCoef * gcPhs;
922                 dUdaSin  += dUdaCoef * gsMhc;
923 
924                 // dU / dh components
925                 dUdhCos  += cf_1 * (kJNS * (cnm * dGdh + snm * dHdh) + auxiliaryElements.getH() * dKgcPhsx2);
926                 dUdhSin  += cf_1 * (kJNS * (snm * dGdh - cnm * dHdh) + auxiliaryElements.getH() * dKgsMhcx2);
927 
928                 // dU / dk components
929                 dUdkCos  += cf_1 * (kJNS * (cnm * dGdk + snm * dHdk) + auxiliaryElements.getK() * dKgcPhsx2);
930                 dUdkSin  += cf_1 * (kJNS * (snm * dGdk - cnm * dHdk) + auxiliaryElements.getK() * dKgsMhcx2);
931 
932                 // dU / dLambda components
933                 dUdlCos  +=  dUdlCoef * gsMhc;
934                 dUdlSin  += -dUdlCoef * gcPhs;
935 
936                 // dU / alpha components
937                 dUdAlCos += cf_2 * (dGdA * cnm + dHdA * snm);
938                 dUdAlSin += cf_2 * (dGdA * snm - dHdA * cnm);
939 
940                 // dU / dBeta components
941                 dUdBeCos += cf_2 * (dGdB * cnm + dHdB * snm);
942                 dUdBeSin += cf_2 * (dGdB * snm - dHdB * cnm);
943 
944                 // dU / dGamma components
945                 dUdGaCos += dUdGaCoef * gcPhs;
946                 dUdGaSin += dUdGaCoef * gsMhc;
947             }
948         }
949 
950         return new double[][] { { dUdaCos,  dUdaSin  },
951                                 { dUdhCos,  dUdhSin  },
952                                 { dUdkCos,  dUdkSin  },
953                                 { dUdlCos,  dUdlSin  },
954                                 { dUdAlCos, dUdAlSin },
955                                 { dUdBeCos, dUdBeSin },
956                                 { dUdGaCos, dUdGaSin } };
957     }
958 
959     /** Compute the n-SUM for potential derivatives components.
960      *  @param <T> the type of the field elements
961      *  @param date current date
962      *  @param j resonant index <i>j</i>
963      *  @param m resonant order <i>m</i>
964      *  @param s d'Alembert characteristic <i>s</i>
965      *  @param maxN maximum possible value for <i>n</i> index
966      *  @param roaPow powers of R/a up to degree <i>n</i>
967      *  @param ghMSJ G<sup>j</sup><sub>m,s</sub> and H<sup>j</sup><sub>m,s</sub> polynomials
968      *  @param gammaMNS &Gamma;<sup>m</sup><sub>n,s</sub>(γ) function
969      *  @param context container for attributes
970      *  @param hansenObjects initialization of hansen objects
971      *  @return Components of U<sub>n</sub> derivatives for fixed j, m, s
972      */
973     private <T extends CalculusFieldElement<T>> T[][] computeNSum(final FieldAbsoluteDate<T> date,
974                                                               final int j, final int m, final int s, final int maxN,
975                                                               final T[] roaPow,
976                                                               final FieldGHmsjPolynomials<T> ghMSJ,
977                                                               final FieldGammaMnsFunction<T> gammaMNS,
978                                                               final FieldDSSTTesseralContext<T> context,
979                                                               final FieldHansenObjects<T> hansenObjects) {
980 
981         // Auxiliary elements related to the current orbit
982         final FieldAuxiliaryElements<T> auxiliaryElements = context.getFieldAuxiliaryElements();
983         // Zero for initialization
984         final Field<T> field = date.getField();
985         final T zero = field.getZero();
986 
987         //spherical harmonics
988         final UnnormalizedSphericalHarmonics harmonics = provider.onDate(date.toAbsoluteDate());
989 
990         // Potential derivatives components
991         T dUdaCos  = zero;
992         T dUdaSin  = zero;
993         T dUdhCos  = zero;
994         T dUdhSin  = zero;
995         T dUdkCos  = zero;
996         T dUdkSin  = zero;
997         T dUdlCos  = zero;
998         T dUdlSin  = zero;
999         T dUdAlCos = zero;
1000         T dUdAlSin = zero;
1001         T dUdBeCos = zero;
1002         T dUdBeSin = zero;
1003         T dUdGaCos = zero;
1004         T dUdGaSin = zero;
1005 
1006         // I^m
1007         @SuppressWarnings("unused")
1008         final int Im = I > 0 ? 1 : (m % 2 == 0 ? 1 : -1);
1009 
1010         // jacobi v, w, indices from 2.7.1-(15)
1011         final int v = FastMath.abs(m - s);
1012         final int w = FastMath.abs(m + s);
1013 
1014         // Initialise lower degree nmin = (Max(2, m, |s|)) for summation over n
1015         final int nmin = FastMath.max(FastMath.max(2, m), FastMath.abs(s));
1016 
1017         //Get the corresponding Hansen object
1018         final int sIndex = maxDegree + (j < 0 ? -s : s);
1019         final int jIndex = FastMath.abs(j);
1020         final FieldHansenTesseralLinear<T> hans = hansenObjects.getHansenObjects()[sIndex][jIndex];
1021 
1022         // n-SUM from nmin to N
1023         for (int n = nmin; n <= maxN; n++) {
1024             // If (n - s) is odd, the contribution is null because of Vmns
1025             if ((n - s) % 2 == 0) {
1026 
1027                 // Vmns coefficient
1028                 final T vMNS   = zero.add(CoefficientsFactory.getVmns(m, n, s));
1029 
1030                 // Inclination function Gamma and derivative
1031                 final T gaMNS  = gammaMNS.getValue(m, n, s);
1032                 final T dGaMNS = gammaMNS.getDerivative(m, n, s);
1033 
1034                 // Hansen kernel value and derivative
1035                 final T kJNS   = hans.getValue(-n - 1, context.getChi());
1036                 final T dkJNS  = hans.getDerivative(-n - 1, context.getChi());
1037 
1038                 // Gjms, Hjms polynomials and derivatives
1039                 final T gMSJ   = ghMSJ.getGmsj(m, s, j);
1040                 final T hMSJ   = ghMSJ.getHmsj(m, s, j);
1041                 final T dGdh   = ghMSJ.getdGmsdh(m, s, j);
1042                 final T dGdk   = ghMSJ.getdGmsdk(m, s, j);
1043                 final T dGdA   = ghMSJ.getdGmsdAlpha(m, s, j);
1044                 final T dGdB   = ghMSJ.getdGmsdBeta(m, s, j);
1045                 final T dHdh   = ghMSJ.getdHmsdh(m, s, j);
1046                 final T dHdk   = ghMSJ.getdHmsdk(m, s, j);
1047                 final T dHdA   = ghMSJ.getdHmsdAlpha(m, s, j);
1048                 final T dHdB   = ghMSJ.getdHmsdBeta(m, s, j);
1049 
1050                 // Jacobi l-index from 2.7.1-(15)
1051                 final int l = FastMath.min(n - m, n - FastMath.abs(s));
1052                 // Jacobi polynomial and derivative
1053                 final FieldGradient<T> jacobi =
1054                         JacobiPolynomials.getValue(l, v, w, FieldGradient.variable(1, 0, auxiliaryElements.getGamma()));
1055 
1056                 // Geopotential coefficients
1057                 final T cnm = zero.add(harmonics.getUnnormalizedCnm(n, m));
1058                 final T snm = zero.add(harmonics.getUnnormalizedSnm(n, m));
1059 
1060                 // Common factors from expansion of equations 3.3-4
1061                 final T cf_0      = roaPow[n].multiply(Im).multiply(vMNS);
1062                 final T cf_1      = cf_0.multiply(gaMNS).multiply(jacobi.getValue());
1063                 final T cf_2      = cf_1.multiply(kJNS);
1064                 final T gcPhs     = gMSJ.multiply(cnm).add(hMSJ.multiply(snm));
1065                 final T gsMhc     = gMSJ.multiply(snm).subtract(hMSJ.multiply(cnm));
1066                 final T dKgcPhsx2 = dkJNS.multiply(gcPhs).multiply(2.);
1067                 final T dKgsMhcx2 = dkJNS.multiply(gsMhc).multiply(2.);
1068                 final T dUdaCoef  = cf_2.multiply(n + 1);
1069                 final T dUdlCoef  = cf_2.multiply(j);
1070                 final T dUdGaCoef = cf_0.multiply(kJNS).multiply(dGaMNS.multiply(jacobi.getValue()).add(gaMNS.multiply(jacobi.getGradient()[0])));
1071 
1072                 // dU / da components
1073                 dUdaCos  = dUdaCos.add(dUdaCoef.multiply(gcPhs));
1074                 dUdaSin  = dUdaSin.add(dUdaCoef.multiply(gsMhc));
1075 
1076                 // dU / dh components
1077                 dUdhCos  = dUdhCos.add(cf_1.multiply(kJNS.multiply(cnm.multiply(dGdh).add(snm.multiply(dHdh))).add(dKgcPhsx2.multiply(auxiliaryElements.getH()))));
1078                 dUdhSin  = dUdhSin.add(cf_1.multiply(kJNS.multiply(snm.multiply(dGdh).subtract(cnm.multiply(dHdh))).add(dKgsMhcx2.multiply(auxiliaryElements.getH()))));
1079 
1080                 // dU / dk components
1081                 dUdkCos  = dUdkCos.add(cf_1.multiply(kJNS.multiply(cnm.multiply(dGdk).add(snm.multiply(dHdk))).add(dKgcPhsx2.multiply(auxiliaryElements.getK()))));
1082                 dUdkSin  = dUdkSin.add(cf_1.multiply(kJNS.multiply(snm.multiply(dGdk).subtract(cnm.multiply(dHdk))).add(dKgsMhcx2.multiply(auxiliaryElements.getK()))));
1083 
1084                 // dU / dLambda components
1085                 dUdlCos  = dUdlCos.add(dUdlCoef.multiply(gsMhc));
1086                 dUdlSin  = dUdlSin.add(dUdlCoef.multiply(gcPhs).negate());
1087 
1088                 // dU / alpha components
1089                 dUdAlCos = dUdAlCos.add(cf_2.multiply(dGdA.multiply(cnm).add(dHdA.multiply(snm))));
1090                 dUdAlSin = dUdAlSin.add(cf_2.multiply(dGdA.multiply(snm).subtract(dHdA.multiply(cnm))));
1091 
1092                 // dU / dBeta components
1093                 dUdBeCos = dUdBeCos.add(cf_2.multiply(dGdB.multiply(cnm).add(dHdB.multiply(snm))));
1094                 dUdBeSin = dUdBeSin.add(cf_2.multiply(dGdB.multiply(snm).subtract(dHdB.multiply(cnm))));
1095 
1096                 // dU / dGamma components
1097                 dUdGaCos = dUdGaCos.add(dUdGaCoef.multiply(gcPhs));
1098                 dUdGaSin = dUdGaSin.add(dUdGaCoef.multiply(gsMhc));
1099             }
1100         }
1101 
1102         final T[][] derivatives = MathArrays.buildArray(field, 7, 2);
1103         derivatives[0][0] = dUdaCos;
1104         derivatives[0][1] = dUdaSin;
1105         derivatives[1][0] = dUdhCos;
1106         derivatives[1][1] = dUdhSin;
1107         derivatives[2][0] = dUdkCos;
1108         derivatives[2][1] = dUdkSin;
1109         derivatives[3][0] = dUdlCos;
1110         derivatives[3][1] = dUdlSin;
1111         derivatives[4][0] = dUdAlCos;
1112         derivatives[4][1] = dUdAlSin;
1113         derivatives[5][0] = dUdBeCos;
1114         derivatives[5][1] = dUdBeSin;
1115         derivatives[6][0] = dUdGaCos;
1116         derivatives[6][1] = dUdGaSin;
1117 
1118         return derivatives;
1119 
1120     }
1121 
1122     /** {@inheritDoc} */
1123     @Override
1124     public void registerAttitudeProvider(final AttitudeProvider attitudeProvider) {
1125         //nothing is done since this contribution is not sensitive to attitude
1126     }
1127 
1128     /** Compute the C<sup>j</sup> and the S<sup>j</sup> coefficients.
1129      *  <p>
1130      *  Those coefficients are given in Danielson paper by substituting the
1131      *  disturbing function (2.7.1-16) with m != 0 into (2.2-10)
1132      *  </p>
1133      */
1134     private class FourierCjSjCoefficients {
1135 
1136         /** Absolute limit for j ( -jMax <= j <= jMax ).  */
1137         private final int jMax;
1138 
1139         /** The C<sub>i</sub><sup>jm</sup> coefficients.
1140          * <p>
1141          * The index order is [m][j][i] <br/>
1142          * The i index corresponds to the C<sub>i</sub><sup>jm</sup> coefficients used to
1143          * compute the following: <br/>
1144          * - da/dt <br/>
1145          * - dk/dt <br/>
1146          * - dh/dt / dk <br/>
1147          * - dq/dt <br/>
1148          * - dp/dt / dα <br/>
1149          * - dλ/dt / dβ <br/>
1150          * </p>
1151          */
1152         private final double[][][] cCoef;
1153 
1154         /** The S<sub>i</sub><sup>jm</sup> coefficients.
1155          * <p>
1156          * The index order is [m][j][i] <br/>
1157          * The i index corresponds to the C<sub>i</sub><sup>jm</sup> coefficients used to
1158          * compute the following: <br/>
1159          * - da/dt <br/>
1160          * - dk/dt <br/>
1161          * - dh/dt / dk <br/>
1162          * - dq/dt <br/>
1163          * - dp/dt / dα <br/>
1164          * - dλ/dt / dβ <br/>
1165          * </p>
1166          */
1167         private final double[][][] sCoef;
1168 
1169         /** G<sub>ms</sub><sup>j</sup> and H<sub>ms</sub><sup>j</sup> polynomials. */
1170         private GHmsjPolynomials ghMSJ;
1171 
1172         /** &Gamma;<sub>ns</sub><sup>m</sup> function. */
1173         private GammaMnsFunction gammaMNS;
1174 
1175         /** R / a up to power degree. */
1176         private final double[] roaPow;
1177 
1178         /** Create a set of C<sub>i</sub><sup>jm</sup> and S<sub>i</sub><sup>jm</sup> coefficients.
1179          *  @param jMax absolute limit for j ( -jMax <= j <= jMax )
1180          *  @param mMax maximum value for m
1181          */
1182         FourierCjSjCoefficients(final int jMax, final int mMax) {
1183             // initialise fields
1184             final int rows    = mMax + 1;
1185             final int columns = 2 * jMax + 1;
1186             this.jMax         = jMax;
1187             this.cCoef        = new double[rows][columns][6];
1188             this.sCoef        = new double[rows][columns][6];
1189             this.roaPow       = new double[maxDegree + 1];
1190             roaPow[0] = 1.;
1191         }
1192 
1193         /**
1194          * Generate the coefficients.
1195          * @param date the current date
1196          * @param context container for attributes
1197          * @param hansenObjects initialization of hansen objects
1198          */
1199         public void generateCoefficients(final AbsoluteDate date, final DSSTTesseralContext context,
1200                                          final HansenObjects hansenObjects) {
1201 
1202             final AuxiliaryElements auxiliaryElements = context.getAuxiliaryElements();
1203 
1204             // Compute only if there is at least one non-resonant tesseral
1205             if (!nonResOrders.isEmpty() || maxDegreeTesseralSP < 0) {
1206                 // Gmsj and Hmsj polynomials
1207                 ghMSJ = new GHmsjPolynomials(auxiliaryElements.getK(), auxiliaryElements.getH(), auxiliaryElements.getAlpha(), auxiliaryElements.getBeta(), I);
1208 
1209                 // GAMMAmns function
1210                 gammaMNS = new GammaMnsFunction(maxDegree, auxiliaryElements.getGamma(), I);
1211 
1212                 final int maxRoaPower = FastMath.max(maxDegreeTesseralSP, maxDegreeMdailyTesseralSP);
1213 
1214                 // R / a up to power degree
1215                 for (int i = 1; i <= maxRoaPower; i++) {
1216                     roaPow[i] = context.getRoa() * roaPow[i - 1];
1217                 }
1218 
1219                 //generate the m-daily coefficients
1220                 for (int m = 1; m <= maxOrderMdailyTesseralSP; m++) {
1221                     buildFourierCoefficients(date, m, 0, maxDegreeMdailyTesseralSP, context, hansenObjects);
1222                 }
1223 
1224                 // generate the other coefficients only if required
1225                 if (maxDegreeTesseralSP >= 0) {
1226                     for (int m: nonResOrders.keySet()) {
1227                         final List<Integer> listJ = nonResOrders.get(m);
1228 
1229                         for (int j: listJ) {
1230                             buildFourierCoefficients(date, m, j, maxDegreeTesseralSP, context, hansenObjects);
1231                         }
1232                     }
1233                 }
1234             }
1235         }
1236 
1237         /** Build a set of fourier coefficients for a given m and j.
1238          *
1239          * @param date the date of the coefficients
1240          * @param m m index
1241          * @param j j index
1242          * @param maxN  maximum value for n index
1243          * @param context container for attributes
1244          * @param hansenObjects initialization of hansen objects
1245          */
1246         private void buildFourierCoefficients(final AbsoluteDate date,
1247                final int m, final int j, final int maxN, final DSSTTesseralContext context,
1248                final HansenObjects hansenObjects) {
1249 
1250             final AuxiliaryElements auxiliaryElements = context.getAuxiliaryElements();
1251 
1252             // Potential derivatives components for a given non-resonant pair {j,m}
1253             double dRdaCos  = 0.;
1254             double dRdaSin  = 0.;
1255             double dRdhCos  = 0.;
1256             double dRdhSin  = 0.;
1257             double dRdkCos  = 0.;
1258             double dRdkSin  = 0.;
1259             double dRdlCos  = 0.;
1260             double dRdlSin  = 0.;
1261             double dRdAlCos = 0.;
1262             double dRdAlSin = 0.;
1263             double dRdBeCos = 0.;
1264             double dRdBeSin = 0.;
1265             double dRdGaCos = 0.;
1266             double dRdGaSin = 0.;
1267 
1268             // s-SUM from -sMin to sMax
1269             final int sMin = j == 0 ? maxEccPowMdailyTesseralSP : maxEccPowTesseralSP;
1270             final int sMax = j == 0 ? maxEccPowMdailyTesseralSP : maxEccPowTesseralSP;
1271             for (int s = 0; s <= sMax; s++) {
1272 
1273                 // n-SUM for s positive
1274                 final double[][] nSumSpos = computeNSum(date, j, m, s, maxN,
1275                                                         roaPow, ghMSJ, gammaMNS, context, hansenObjects);
1276                 dRdaCos  += nSumSpos[0][0];
1277                 dRdaSin  += nSumSpos[0][1];
1278                 dRdhCos  += nSumSpos[1][0];
1279                 dRdhSin  += nSumSpos[1][1];
1280                 dRdkCos  += nSumSpos[2][0];
1281                 dRdkSin  += nSumSpos[2][1];
1282                 dRdlCos  += nSumSpos[3][0];
1283                 dRdlSin  += nSumSpos[3][1];
1284                 dRdAlCos += nSumSpos[4][0];
1285                 dRdAlSin += nSumSpos[4][1];
1286                 dRdBeCos += nSumSpos[5][0];
1287                 dRdBeSin += nSumSpos[5][1];
1288                 dRdGaCos += nSumSpos[6][0];
1289                 dRdGaSin += nSumSpos[6][1];
1290 
1291                 // n-SUM for s negative
1292                 if (s > 0 && s <= sMin) {
1293                     final double[][] nSumSneg = computeNSum(date, j, m, -s, maxN,
1294                                                             roaPow, ghMSJ, gammaMNS, context, hansenObjects);
1295                     dRdaCos  += nSumSneg[0][0];
1296                     dRdaSin  += nSumSneg[0][1];
1297                     dRdhCos  += nSumSneg[1][0];
1298                     dRdhSin  += nSumSneg[1][1];
1299                     dRdkCos  += nSumSneg[2][0];
1300                     dRdkSin  += nSumSneg[2][1];
1301                     dRdlCos  += nSumSneg[3][0];
1302                     dRdlSin  += nSumSneg[3][1];
1303                     dRdAlCos += nSumSneg[4][0];
1304                     dRdAlSin += nSumSneg[4][1];
1305                     dRdBeCos += nSumSneg[5][0];
1306                     dRdBeSin += nSumSneg[5][1];
1307                     dRdGaCos += nSumSneg[6][0];
1308                     dRdGaSin += nSumSneg[6][1];
1309                 }
1310             }
1311             dRdaCos  *= -context.getMoa() / auxiliaryElements.getSma();
1312             dRdaSin  *= -context.getMoa() / auxiliaryElements.getSma();
1313             dRdhCos  *=  context.getMoa();
1314             dRdhSin  *=  context.getMoa();
1315             dRdkCos  *=  context.getMoa();
1316             dRdkSin  *=  context.getMoa();
1317             dRdlCos  *=  context.getMoa();
1318             dRdlSin  *=  context.getMoa();
1319             dRdAlCos *=  context.getMoa();
1320             dRdAlSin *=  context.getMoa();
1321             dRdBeCos *=  context.getMoa();
1322             dRdBeSin *=  context.getMoa();
1323             dRdGaCos *=  context.getMoa();
1324             dRdGaSin *=  context.getMoa();
1325 
1326             // Compute the cross derivative operator :
1327             final double RAlphaGammaCos   = auxiliaryElements.getAlpha() * dRdGaCos - auxiliaryElements.getGamma() * dRdAlCos;
1328             final double RAlphaGammaSin   = auxiliaryElements.getAlpha() * dRdGaSin - auxiliaryElements.getGamma() * dRdAlSin;
1329             final double RAlphaBetaCos    = auxiliaryElements.getAlpha() * dRdBeCos - auxiliaryElements.getBeta()  * dRdAlCos;
1330             final double RAlphaBetaSin    = auxiliaryElements.getAlpha() * dRdBeSin - auxiliaryElements.getBeta()  * dRdAlSin;
1331             final double RBetaGammaCos    =  auxiliaryElements.getBeta() * dRdGaCos - auxiliaryElements.getGamma() * dRdBeCos;
1332             final double RBetaGammaSin    =  auxiliaryElements.getBeta() * dRdGaSin - auxiliaryElements.getGamma() * dRdBeSin;
1333             final double RhkCos           =     auxiliaryElements.getH() * dRdkCos  -     auxiliaryElements.getK() * dRdhCos;
1334             final double RhkSin           =     auxiliaryElements.getH() * dRdkSin  -     auxiliaryElements.getK() * dRdhSin;
1335             final double pRagmIqRbgoABCos = (auxiliaryElements.getP() * RAlphaGammaCos - I * auxiliaryElements.getQ() * RBetaGammaCos) * context.getOoAB();
1336             final double pRagmIqRbgoABSin = (auxiliaryElements.getP() * RAlphaGammaSin - I * auxiliaryElements.getQ() * RBetaGammaSin) * context.getOoAB();
1337             final double RhkmRabmdRdlCos  =  RhkCos - RAlphaBetaCos - dRdlCos;
1338             final double RhkmRabmdRdlSin  =  RhkSin - RAlphaBetaSin - dRdlSin;
1339 
1340             // da/dt
1341             cCoef[m][j + jMax][0] = context.getAx2oA() * dRdlCos;
1342             sCoef[m][j + jMax][0] = context.getAx2oA() * dRdlSin;
1343 
1344             // dk/dt
1345             cCoef[m][j + jMax][1] = -(context.getBoA() * dRdhCos + auxiliaryElements.getH() * pRagmIqRbgoABCos + auxiliaryElements.getK() * context.getBoABpo() * dRdlCos);
1346             sCoef[m][j + jMax][1] = -(context.getBoA() * dRdhSin + auxiliaryElements.getH() * pRagmIqRbgoABSin + auxiliaryElements.getK() * context.getBoABpo() * dRdlSin);
1347 
1348             // dh/dt
1349             cCoef[m][j + jMax][2] = context.getBoA() * dRdkCos + auxiliaryElements.getK() * pRagmIqRbgoABCos - auxiliaryElements.getH() * context.getBoABpo() * dRdlCos;
1350             sCoef[m][j + jMax][2] = context.getBoA() * dRdkSin + auxiliaryElements.getK() * pRagmIqRbgoABSin - auxiliaryElements.getH() * context.getBoABpo() * dRdlSin;
1351 
1352             // dq/dt
1353             cCoef[m][j + jMax][3] = context.getCo2AB() * (auxiliaryElements.getQ() * RhkmRabmdRdlCos - I * RAlphaGammaCos);
1354             sCoef[m][j + jMax][3] = context.getCo2AB() * (auxiliaryElements.getQ() * RhkmRabmdRdlSin - I * RAlphaGammaSin);
1355 
1356             // dp/dt
1357             cCoef[m][j + jMax][4] = context.getCo2AB() * (auxiliaryElements.getP() * RhkmRabmdRdlCos - RBetaGammaCos);
1358             sCoef[m][j + jMax][4] = context.getCo2AB() * (auxiliaryElements.getP() * RhkmRabmdRdlSin - RBetaGammaSin);
1359 
1360             // dλ/dt
1361             cCoef[m][j + jMax][5] = -context.getAx2oA() * dRdaCos + context.getBoABpo() * (auxiliaryElements.getH() * dRdhCos + auxiliaryElements.getK() * dRdkCos) + pRagmIqRbgoABCos;
1362             sCoef[m][j + jMax][5] = -context.getAx2oA() * dRdaSin + context.getBoABpo() * (auxiliaryElements.getH() * dRdhSin + auxiliaryElements.getK() * dRdkSin) + pRagmIqRbgoABSin;
1363         }
1364 
1365         /** Get the coefficient C<sub>i</sub><sup>jm</sup>.
1366          * @param i i index - corresponds to the required variation
1367          * @param j j index
1368          * @param m m index
1369          * @return the coefficient C<sub>i</sub><sup>jm</sup>
1370          */
1371         public double getCijm(final int i, final int j, final int m) {
1372             return cCoef[m][j + jMax][i];
1373         }
1374 
1375         /** Get the coefficient S<sub>i</sub><sup>jm</sup>.
1376          * @param i i index - corresponds to the required variation
1377          * @param j j index
1378          * @param m m index
1379          * @return the coefficient S<sub>i</sub><sup>jm</sup>
1380          */
1381         public double getSijm(final int i, final int j, final int m) {
1382             return sCoef[m][j + jMax][i];
1383         }
1384     }
1385 
1386     /** Compute the C<sup>j</sup> and the S<sup>j</sup> coefficients.
1387      *  <p>
1388      *  Those coefficients are given in Danielson paper by substituting the
1389      *  disturbing function (2.7.1-16) with m != 0 into (2.2-10)
1390      *  </p>
1391      */
1392     private class FieldFourierCjSjCoefficients <T extends CalculusFieldElement<T>> {
1393 
1394         /** Absolute limit for j ( -jMax <= j <= jMax ).  */
1395         private final int jMax;
1396 
1397         /** The C<sub>i</sub><sup>jm</sup> coefficients.
1398          * <p>
1399          * The index order is [m][j][i] <br/>
1400          * The i index corresponds to the C<sub>i</sub><sup>jm</sup> coefficients used to
1401          * compute the following: <br/>
1402          * - da/dt <br/>
1403          * - dk/dt <br/>
1404          * - dh/dt / dk <br/>
1405          * - dq/dt <br/>
1406          * - dp/dt / dα <br/>
1407          * - dλ/dt / dβ <br/>
1408          * </p>
1409          */
1410         private final T[][][] cCoef;
1411 
1412         /** The S<sub>i</sub><sup>jm</sup> coefficients.
1413          * <p>
1414          * The index order is [m][j][i] <br/>
1415          * The i index corresponds to the C<sub>i</sub><sup>jm</sup> coefficients used to
1416          * compute the following: <br/>
1417          * - da/dt <br/>
1418          * - dk/dt <br/>
1419          * - dh/dt / dk <br/>
1420          * - dq/dt <br/>
1421          * - dp/dt / dα <br/>
1422          * - dλ/dt / dβ <br/>
1423          * </p>
1424          */
1425         private final T[][][] sCoef;
1426 
1427         /** G<sub>ms</sub><sup>j</sup> and H<sub>ms</sub><sup>j</sup> polynomials. */
1428         private FieldGHmsjPolynomials<T> ghMSJ;
1429 
1430         /** &Gamma;<sub>ns</sub><sup>m</sup> function. */
1431         private FieldGammaMnsFunction<T> gammaMNS;
1432 
1433         /** R / a up to power degree. */
1434         private final T[] roaPow;
1435 
1436         /** Create a set of C<sub>i</sub><sup>jm</sup> and S<sub>i</sub><sup>jm</sup> coefficients.
1437          *  @param jMax absolute limit for j ( -jMax <= j <= jMax )
1438          *  @param mMax maximum value for m
1439          *  @param field field used by default
1440          */
1441         FieldFourierCjSjCoefficients(final int jMax, final int mMax, final Field<T> field) {
1442             // initialise fields
1443             final T zero = field.getZero();
1444             final int rows    = mMax + 1;
1445             final int columns = 2 * jMax + 1;
1446             this.jMax         = jMax;
1447             this.cCoef        = MathArrays.buildArray(field, rows, columns, 6);
1448             this.sCoef        = MathArrays.buildArray(field, rows, columns, 6);
1449             this.roaPow       = MathArrays.buildArray(field, maxDegree + 1);
1450             roaPow[0] = zero.add(1.);
1451         }
1452 
1453         /**
1454          * Generate the coefficients.
1455          * @param date the current date
1456          * @param context container for attributes
1457          * @param hansenObjects initialization of hansen objects
1458          * @param field field used by default
1459          */
1460         public void generateCoefficients(final FieldAbsoluteDate<T> date,
1461                                          final FieldDSSTTesseralContext<T> context,
1462                                          final FieldHansenObjects<T> hansenObjects,
1463                                          final Field<T> field) {
1464 
1465             final FieldAuxiliaryElements<T> auxiliaryElements = context.getFieldAuxiliaryElements();
1466             // Compute only if there is at least one non-resonant tesseral
1467             if (!nonResOrders.isEmpty() || maxDegreeTesseralSP < 0) {
1468                 // Gmsj and Hmsj polynomials
1469                 ghMSJ = new FieldGHmsjPolynomials<>(auxiliaryElements.getK(), auxiliaryElements.getH(), auxiliaryElements.getAlpha(), auxiliaryElements.getBeta(), I, field);
1470 
1471                 // GAMMAmns function
1472                 gammaMNS = new FieldGammaMnsFunction<>(maxDegree, auxiliaryElements.getGamma(), I, field);
1473 
1474                 final int maxRoaPower = FastMath.max(maxDegreeTesseralSP, maxDegreeMdailyTesseralSP);
1475 
1476                 // R / a up to power degree
1477                 for (int i = 1; i <= maxRoaPower; i++) {
1478                     roaPow[i] = context.getRoa().multiply(roaPow[i - 1]);
1479                 }
1480 
1481                 //generate the m-daily coefficients
1482                 for (int m = 1; m <= maxOrderMdailyTesseralSP; m++) {
1483                     buildFourierCoefficients(date, m, 0, maxDegreeMdailyTesseralSP, context, hansenObjects, field);
1484                 }
1485 
1486                 // generate the other coefficients only if required
1487                 if (maxDegreeTesseralSP >= 0) {
1488                     for (int m: nonResOrders.keySet()) {
1489                         final List<Integer> listJ = nonResOrders.get(m);
1490 
1491                         for (int j: listJ) {
1492                             buildFourierCoefficients(date, m, j, maxDegreeTesseralSP, context, hansenObjects, field);
1493                         }
1494                     }
1495                 }
1496             }
1497         }
1498 
1499         /** Build a set of fourier coefficients for a given m and j.
1500          *
1501          * @param date the date of the coefficients
1502          * @param m m index
1503          * @param j j index
1504          * @param maxN  maximum value for n index
1505          * @param context container for attributes
1506          * @param hansenObjects initialization of hansen objects
1507          * @param field field used by default
1508          */
1509         private void buildFourierCoefficients(final FieldAbsoluteDate<T> date,
1510                                               final int m, final int j, final int maxN,
1511                                               final FieldDSSTTesseralContext<T> context,
1512                                               final FieldHansenObjects<T> hansenObjects,
1513                                               final Field<T> field) {
1514 
1515             // Zero
1516             final T zero = field.getZero();
1517             // Common parameters
1518             final FieldAuxiliaryElements<T> auxiliaryElements = context.getFieldAuxiliaryElements();
1519 
1520             // Potential derivatives components for a given non-resonant pair {j,m}
1521             T dRdaCos  = zero;
1522             T dRdaSin  = zero;
1523             T dRdhCos  = zero;
1524             T dRdhSin  = zero;
1525             T dRdkCos  = zero;
1526             T dRdkSin  = zero;
1527             T dRdlCos  = zero;
1528             T dRdlSin  = zero;
1529             T dRdAlCos = zero;
1530             T dRdAlSin = zero;
1531             T dRdBeCos = zero;
1532             T dRdBeSin = zero;
1533             T dRdGaCos = zero;
1534             T dRdGaSin = zero;
1535 
1536             // s-SUM from -sMin to sMax
1537             final int sMin = j == 0 ? maxEccPowMdailyTesseralSP : maxEccPowTesseralSP;
1538             final int sMax = j == 0 ? maxEccPowMdailyTesseralSP : maxEccPowTesseralSP;
1539             for (int s = 0; s <= sMax; s++) {
1540 
1541                 // n-SUM for s positive
1542                 final T[][] nSumSpos = computeNSum(date, j, m, s, maxN,
1543                                                         roaPow, ghMSJ, gammaMNS, context, hansenObjects);
1544                 dRdaCos  =  dRdaCos.add(nSumSpos[0][0]);
1545                 dRdaSin  =  dRdaSin.add(nSumSpos[0][1]);
1546                 dRdhCos  =  dRdhCos.add(nSumSpos[1][0]);
1547                 dRdhSin  =  dRdhSin.add(nSumSpos[1][1]);
1548                 dRdkCos  =  dRdkCos.add(nSumSpos[2][0]);
1549                 dRdkSin  =  dRdkSin.add(nSumSpos[2][1]);
1550                 dRdlCos  =  dRdlCos.add(nSumSpos[3][0]);
1551                 dRdlSin  =  dRdlSin.add(nSumSpos[3][1]);
1552                 dRdAlCos = dRdAlCos.add(nSumSpos[4][0]);
1553                 dRdAlSin = dRdAlSin.add(nSumSpos[4][1]);
1554                 dRdBeCos = dRdBeCos.add(nSumSpos[5][0]);
1555                 dRdBeSin = dRdBeSin.add(nSumSpos[5][1]);
1556                 dRdGaCos = dRdGaCos.add(nSumSpos[6][0]);
1557                 dRdGaSin = dRdGaSin.add(nSumSpos[6][1]);
1558 
1559                 // n-SUM for s negative
1560                 if (s > 0 && s <= sMin) {
1561                     final T[][] nSumSneg = computeNSum(date, j, m, -s, maxN,
1562                                                             roaPow, ghMSJ, gammaMNS, context, hansenObjects);
1563                     dRdaCos  =  dRdaCos.add(nSumSneg[0][0]);
1564                     dRdaSin  =  dRdaSin.add(nSumSneg[0][1]);
1565                     dRdhCos  =  dRdhCos.add(nSumSneg[1][0]);
1566                     dRdhSin  =  dRdhSin.add(nSumSneg[1][1]);
1567                     dRdkCos  =  dRdkCos.add(nSumSneg[2][0]);
1568                     dRdkSin  =  dRdkSin.add(nSumSneg[2][1]);
1569                     dRdlCos  =  dRdlCos.add(nSumSneg[3][0]);
1570                     dRdlSin  =  dRdlSin.add(nSumSneg[3][1]);
1571                     dRdAlCos = dRdAlCos.add(nSumSneg[4][0]);
1572                     dRdAlSin = dRdAlSin.add(nSumSneg[4][1]);
1573                     dRdBeCos = dRdBeCos.add(nSumSneg[5][0]);
1574                     dRdBeSin = dRdBeSin.add(nSumSneg[5][1]);
1575                     dRdGaCos = dRdGaCos.add(nSumSneg[6][0]);
1576                     dRdGaSin = dRdGaSin.add(nSumSneg[6][1]);
1577                 }
1578             }
1579             dRdaCos  =  dRdaCos.multiply(context.getMoa().negate().divide(auxiliaryElements.getSma()));
1580             dRdaSin  =  dRdaSin.multiply(context.getMoa().negate().divide(auxiliaryElements.getSma()));
1581             dRdhCos  =  dRdhCos.multiply(context.getMoa());
1582             dRdhSin  =  dRdhSin.multiply(context.getMoa());
1583             dRdkCos  =  dRdkCos.multiply(context.getMoa());
1584             dRdkSin  =  dRdkSin.multiply(context.getMoa());
1585             dRdlCos  =  dRdlCos.multiply(context.getMoa());
1586             dRdlSin  =  dRdlSin.multiply(context.getMoa());
1587             dRdAlCos = dRdAlCos.multiply(context.getMoa());
1588             dRdAlSin = dRdAlSin.multiply(context.getMoa());
1589             dRdBeCos = dRdBeCos.multiply(context.getMoa());
1590             dRdBeSin = dRdBeSin.multiply(context.getMoa());
1591             dRdGaCos = dRdGaCos.multiply(context.getMoa());
1592             dRdGaSin = dRdGaSin.multiply(context.getMoa());
1593 
1594             // Compute the cross derivative operator :
1595             final T RAlphaGammaCos   = auxiliaryElements.getAlpha().multiply(dRdGaCos).subtract(auxiliaryElements.getGamma().multiply(dRdAlCos));
1596             final T RAlphaGammaSin   = auxiliaryElements.getAlpha().multiply(dRdGaSin).subtract(auxiliaryElements.getGamma().multiply(dRdAlSin));
1597             final T RAlphaBetaCos    = auxiliaryElements.getAlpha().multiply(dRdBeCos).subtract(auxiliaryElements.getBeta().multiply(dRdAlCos));
1598             final T RAlphaBetaSin    = auxiliaryElements.getAlpha().multiply(dRdBeSin).subtract(auxiliaryElements.getBeta().multiply(dRdAlSin));
1599             final T RBetaGammaCos    =  auxiliaryElements.getBeta().multiply(dRdGaCos).subtract(auxiliaryElements.getGamma().multiply(dRdBeCos));
1600             final T RBetaGammaSin    =  auxiliaryElements.getBeta().multiply(dRdGaSin).subtract(auxiliaryElements.getGamma().multiply(dRdBeSin));
1601             final T RhkCos           =     auxiliaryElements.getH().multiply(dRdkCos).subtract(auxiliaryElements.getK().multiply(dRdhCos));
1602             final T RhkSin           =     auxiliaryElements.getH().multiply(dRdkSin).subtract(auxiliaryElements.getK().multiply(dRdhSin));
1603             final T pRagmIqRbgoABCos = (auxiliaryElements.getP().multiply(RAlphaGammaCos).subtract(auxiliaryElements.getQ().multiply(RBetaGammaCos).multiply(I))).multiply(context.getOoAB());
1604             final T pRagmIqRbgoABSin = (auxiliaryElements.getP().multiply(RAlphaGammaSin).subtract(auxiliaryElements.getQ().multiply(RBetaGammaSin).multiply(I))).multiply(context.getOoAB());
1605             final T RhkmRabmdRdlCos  =  RhkCos.subtract(RAlphaBetaCos).subtract(dRdlCos);
1606             final T RhkmRabmdRdlSin  =  RhkSin.subtract(RAlphaBetaSin).subtract(dRdlSin);
1607 
1608             // da/dt
1609             cCoef[m][j + jMax][0] = context.getAx2oA().multiply(dRdlCos);
1610             sCoef[m][j + jMax][0] = context.getAx2oA().multiply(dRdlSin);
1611 
1612             // dk/dt
1613             cCoef[m][j + jMax][1] = (context.getBoA().multiply(dRdhCos).add(auxiliaryElements.getH().multiply(pRagmIqRbgoABCos)).add(auxiliaryElements.getK().multiply(context.getBoABpo()).multiply(dRdlCos))).negate();
1614             sCoef[m][j + jMax][1] = (context.getBoA().multiply(dRdhSin).add(auxiliaryElements.getH().multiply(pRagmIqRbgoABSin)).add(auxiliaryElements.getK().multiply(context.getBoABpo()).multiply(dRdlSin))).negate();
1615 
1616             // dh/dt
1617             cCoef[m][j + jMax][2] = context.getBoA().multiply(dRdkCos).add(auxiliaryElements.getK().multiply(pRagmIqRbgoABCos)).subtract(auxiliaryElements.getH().multiply(context.getBoABpo()).multiply(dRdlCos));
1618             sCoef[m][j + jMax][2] = context.getBoA().multiply(dRdkSin).add(auxiliaryElements.getK().multiply(pRagmIqRbgoABSin)).subtract(auxiliaryElements.getH().multiply(context.getBoABpo()).multiply(dRdlSin));
1619 
1620             // dq/dt
1621             cCoef[m][j + jMax][3] = context.getCo2AB().multiply(auxiliaryElements.getQ().multiply(RhkmRabmdRdlCos).subtract(RAlphaGammaCos.multiply(I)));
1622             sCoef[m][j + jMax][3] = context.getCo2AB().multiply(auxiliaryElements.getQ().multiply(RhkmRabmdRdlSin).subtract(RAlphaGammaSin.multiply(I)));
1623 
1624             // dp/dt
1625             cCoef[m][j + jMax][4] = context.getCo2AB().multiply(auxiliaryElements.getP().multiply(RhkmRabmdRdlCos).subtract(RBetaGammaCos));
1626             sCoef[m][j + jMax][4] = context.getCo2AB().multiply(auxiliaryElements.getP().multiply(RhkmRabmdRdlSin).subtract(RBetaGammaSin));
1627 
1628             // dλ/dt
1629             cCoef[m][j + jMax][5] = context.getAx2oA().negate().multiply(dRdaCos).add(context.getBoABpo().multiply(auxiliaryElements.getH().multiply(dRdhCos).add(auxiliaryElements.getK().multiply(dRdkCos)))).add(pRagmIqRbgoABCos);
1630             sCoef[m][j + jMax][5] = context.getAx2oA().negate().multiply(dRdaSin).add(context.getBoABpo().multiply(auxiliaryElements.getH().multiply(dRdhSin).add(auxiliaryElements.getK().multiply(dRdkSin)))).add(pRagmIqRbgoABSin);
1631         }
1632 
1633         /** Get the coefficient C<sub>i</sub><sup>jm</sup>.
1634          * @param i i index - corresponds to the required variation
1635          * @param j j index
1636          * @param m m index
1637          * @return the coefficient C<sub>i</sub><sup>jm</sup>
1638          */
1639         public T getCijm(final int i, final int j, final int m) {
1640             return cCoef[m][j + jMax][i];
1641         }
1642 
1643         /** Get the coefficient S<sub>i</sub><sup>jm</sup>.
1644          * @param i i index - corresponds to the required variation
1645          * @param j j index
1646          * @param m m index
1647          * @return the coefficient S<sub>i</sub><sup>jm</sup>
1648          */
1649         public T getSijm(final int i, final int j, final int m) {
1650             return sCoef[m][j + jMax][i];
1651         }
1652     }
1653 
1654     /** The C<sup>i</sup><sub>m</sub><sub>t</sub> and S<sup>i</sup><sub>m</sub><sub>t</sub> coefficients used to compute
1655      * the short-periodic zonal contribution.
1656      *   <p>
1657      *  Those coefficients are given by expression 2.5.4-5 from the Danielson paper.
1658      *   </p>
1659      *
1660      * @author Sorin Scortan
1661      *
1662      */
1663     private static class TesseralShortPeriodicCoefficients implements ShortPeriodTerms {
1664 
1665         /** Retrograde factor I.
1666          *  <p>
1667          *  DSST model needs equinoctial orbit as internal representation.
1668          *  Classical equinoctial elements have discontinuities when inclination
1669          *  is close to zero. In this representation, I = +1. <br>
1670          *  To avoid this discontinuity, another representation exists and equinoctial
1671          *  elements can be expressed in a different way, called "retrograde" orbit.
1672          *  This implies I = -1. <br>
1673          *  As Orekit doesn't implement the retrograde orbit, I is always set to +1.
1674          *  But for the sake of consistency with the theory, the retrograde factor
1675          *  has been kept in the formulas.
1676          *  </p>
1677          */
1678         private static final int I = 1;
1679 
1680         /** Central body rotating frame. */
1681         private final Frame bodyFrame;
1682 
1683         /** Maximal order to consider for short periodics m-daily tesseral harmonics potential. */
1684         private final int maxOrderMdailyTesseralSP;
1685 
1686         /** Flag to take into account only M-dailies harmonic tesserals for short periodic perturbations.  */
1687         private final boolean mDailiesOnly;
1688 
1689         /** List of non resonant orders with j != 0. */
1690         private final SortedMap<Integer, List<Integer> > nonResOrders;
1691 
1692         /** Maximum value for m index. */
1693         private final int mMax;
1694 
1695         /** Maximum value for j index. */
1696         private final int jMax;
1697 
1698         /** Number of points used in the interpolation process. */
1699         private final int interpolationPoints;
1700 
1701         /** All coefficients slots. */
1702         private final transient TimeSpanMap<Slot> slots;
1703 
1704         /** Constructor.
1705          * @param bodyFrame central body rotating frame
1706          * @param maxOrderMdailyTesseralSP maximal order to consider for short periodics m-daily tesseral harmonics potential
1707          * @param mDailiesOnly flag to take into account only M-dailies harmonic tesserals for short periodic perturbations
1708          * @param nonResOrders lst of non resonant orders with j != 0
1709          * @param mMax maximum value for m index
1710          * @param jMax maximum value for j index
1711          * @param interpolationPoints number of points used in the interpolation process
1712          * @param slots all coefficients slots
1713          */
1714         TesseralShortPeriodicCoefficients(final Frame bodyFrame, final int maxOrderMdailyTesseralSP,
1715                                           final boolean mDailiesOnly, final SortedMap<Integer, List<Integer> > nonResOrders,
1716                                           final int mMax, final int jMax, final int interpolationPoints,
1717                                           final TimeSpanMap<Slot> slots) {
1718             this.bodyFrame                = bodyFrame;
1719             this.maxOrderMdailyTesseralSP = maxOrderMdailyTesseralSP;
1720             this.mDailiesOnly             = mDailiesOnly;
1721             this.nonResOrders             = nonResOrders;
1722             this.mMax                     = mMax;
1723             this.jMax                     = jMax;
1724             this.interpolationPoints      = interpolationPoints;
1725             this.slots                    = slots;
1726         }
1727 
1728         /** Get the slot valid for some date.
1729          * @param meanStates mean states defining the slot
1730          * @return slot valid at the specified date
1731          */
1732         public Slot createSlot(final SpacecraftState... meanStates) {
1733             final Slot         slot  = new Slot(mMax, jMax, interpolationPoints);
1734             final AbsoluteDate first = meanStates[0].getDate();
1735             final AbsoluteDate last  = meanStates[meanStates.length - 1].getDate();
1736             final int compare = first.compareTo(last);
1737             if (compare < 0) {
1738                 slots.addValidAfter(slot, first);
1739             } else if (compare > 0) {
1740                 slots.addValidBefore(slot, first);
1741             } else {
1742                 // single date, valid for all time
1743                 slots.addValidAfter(slot, AbsoluteDate.PAST_INFINITY);
1744             }
1745             return slot;
1746         }
1747 
1748         /** {@inheritDoc} */
1749         @Override
1750         public double[] value(final Orbit meanOrbit) {
1751 
1752             // select the coefficients slot
1753             final Slot slot = slots.get(meanOrbit.getDate());
1754 
1755             // Initialise the short periodic variations
1756             final double[] shortPeriodicVariation = new double[6];
1757 
1758             // Compute only if there is at least one non-resonant tesseral or
1759             // only the m-daily tesseral should be taken into account
1760             if (!nonResOrders.isEmpty() || mDailiesOnly) {
1761 
1762                 //Build an auxiliary object
1763                 final AuxiliaryElements auxiliaryElements = new AuxiliaryElements(meanOrbit, I);
1764 
1765                 // Central body rotation angle from equation 2.7.1-(3)(4).
1766                 final Transform t = bodyFrame.getTransformTo(auxiliaryElements.getFrame(), auxiliaryElements.getDate());
1767                 final Vector3D xB = t.transformVector(Vector3D.PLUS_I);
1768                 final Vector3D yB = t.transformVector(Vector3D.PLUS_J);
1769                 final Vector3D  f = auxiliaryElements.getVectorF();
1770                 final Vector3D  g = auxiliaryElements.getVectorG();
1771                 final double currentTheta = FastMath.atan2(-f.dotProduct(yB) + I * g.dotProduct(xB),
1772                                                             f.dotProduct(xB) + I * g.dotProduct(yB));
1773 
1774                 //Add the m-daily contribution
1775                 for (int m = 1; m <= maxOrderMdailyTesseralSP; m++) {
1776                     // Phase angle
1777                     final double jlMmt  = -m * currentTheta;
1778                     final SinCos scPhi  = FastMath.sinCos(jlMmt);
1779                     final double sinPhi = scPhi.sin();
1780                     final double cosPhi = scPhi.cos();
1781 
1782                     // compute contribution for each element
1783                     final double[] c = slot.getCijm(0, m, meanOrbit.getDate());
1784                     final double[] s = slot.getSijm(0, m, meanOrbit.getDate());
1785                     for (int i = 0; i < 6; i++) {
1786                         shortPeriodicVariation[i] += c[i] * cosPhi + s[i] * sinPhi;
1787                     }
1788                 }
1789 
1790                 // loop through all non-resonant (j,m) pairs
1791                 for (final Map.Entry<Integer, List<Integer>> entry : nonResOrders.entrySet()) {
1792                     final int           m     = entry.getKey();
1793                     final List<Integer> listJ = entry.getValue();
1794 
1795                     for (int j : listJ) {
1796                         // Phase angle
1797                         final double jlMmt  = j * meanOrbit.getLM() - m * currentTheta;
1798                         final SinCos scPhi  = FastMath.sinCos(jlMmt);
1799                         final double sinPhi = scPhi.sin();
1800                         final double cosPhi = scPhi.cos();
1801 
1802                         // compute contribution for each element
1803                         final double[] c = slot.getCijm(j, m, meanOrbit.getDate());
1804                         final double[] s = slot.getSijm(j, m, meanOrbit.getDate());
1805                         for (int i = 0; i < 6; i++) {
1806                             shortPeriodicVariation[i] += c[i] * cosPhi + s[i] * sinPhi;
1807                         }
1808 
1809                     }
1810                 }
1811             }
1812 
1813             return shortPeriodicVariation;
1814 
1815         }
1816 
1817         /** {@inheritDoc} */
1818         @Override
1819         public String getCoefficientsKeyPrefix() {
1820             return DSSTTesseral.SHORT_PERIOD_PREFIX;
1821         }
1822 
1823         /** {@inheritDoc}
1824          * <p>
1825          * For tesseral terms contributions,there are maxOrderMdailyTesseralSP
1826          * m-daily cMm coefficients, maxOrderMdailyTesseralSP m-daily sMm
1827          * coefficients, nbNonResonant cjm coefficients and nbNonResonant
1828          * sjm coefficients, where maxOrderMdailyTesseralSP and nbNonResonant both
1829          * depend on the orbit. The j index is the integer multiplier for the true
1830          * longitude argument and the m index is the integer multiplier for m-dailies.
1831          * </p>
1832          */
1833         @Override
1834         public Map<String, double[]> getCoefficients(final AbsoluteDate date, final Set<String> selected) {
1835 
1836             // select the coefficients slot
1837             final Slot slot = slots.get(date);
1838 
1839             if (!nonResOrders.isEmpty() || mDailiesOnly) {
1840                 final Map<String, double[]> coefficients =
1841                                 new HashMap<String, double[]>(12 * maxOrderMdailyTesseralSP +
1842                                                               12 * nonResOrders.size());
1843 
1844                 for (int m = 1; m <= maxOrderMdailyTesseralSP; m++) {
1845                     storeIfSelected(coefficients, selected, slot.getCijm(0, m, date), DSSTTesseral.CM_COEFFICIENTS, m);
1846                     storeIfSelected(coefficients, selected, slot.getSijm(0, m, date), DSSTTesseral.SM_COEFFICIENTS, m);
1847                 }
1848 
1849                 for (final Map.Entry<Integer, List<Integer>> entry : nonResOrders.entrySet()) {
1850                     final int           m     = entry.getKey();
1851                     final List<Integer> listJ = entry.getValue();
1852 
1853                     for (int j : listJ) {
1854                         for (int i = 0; i < 6; ++i) {
1855                             storeIfSelected(coefficients, selected, slot.getCijm(j, m, date), "c", j, m);
1856                             storeIfSelected(coefficients, selected, slot.getSijm(j, m, date), "s", j, m);
1857                         }
1858                     }
1859                 }
1860 
1861                 return coefficients;
1862 
1863             } else {
1864                 return Collections.emptyMap();
1865             }
1866 
1867         }
1868 
1869         /** Put a coefficient in a map if selected.
1870          * @param map map to populate
1871          * @param selected set of coefficients that should be put in the map
1872          * (empty set means all coefficients are selected)
1873          * @param value coefficient value
1874          * @param id coefficient identifier
1875          * @param indices list of coefficient indices
1876          */
1877         private void storeIfSelected(final Map<String, double[]> map, final Set<String> selected,
1878                                      final double[] value, final String id, final int... indices) {
1879             final StringBuilder keyBuilder = new StringBuilder(getCoefficientsKeyPrefix());
1880             keyBuilder.append(id);
1881             for (int index : indices) {
1882                 keyBuilder.append('[').append(index).append(']');
1883             }
1884             final String key = keyBuilder.toString();
1885             if (selected.isEmpty() || selected.contains(key)) {
1886                 map.put(key, value);
1887             }
1888         }
1889 
1890     }
1891 
1892     /** The C<sup>i</sup><sub>m</sub><sub>t</sub> and S<sup>i</sup><sub>m</sub><sub>t</sub> coefficients used to compute
1893      * the short-periodic zonal contribution.
1894      *   <p>
1895      *  Those coefficients are given by expression 2.5.4-5 from the Danielson paper.
1896      *   </p>
1897      *
1898      * @author Sorin Scortan
1899      *
1900      */
1901     private static class FieldTesseralShortPeriodicCoefficients <T extends CalculusFieldElement<T>> implements FieldShortPeriodTerms<T> {
1902 
1903         /** Retrograde factor I.
1904          *  <p>
1905          *  DSST model needs equinoctial orbit as internal representation.
1906          *  Classical equinoctial elements have discontinuities when inclination
1907          *  is close to zero. In this representation, I = +1. <br>
1908          *  To avoid this discontinuity, another representation exists and equinoctial
1909          *  elements can be expressed in a different way, called "retrograde" orbit.
1910          *  This implies I = -1. <br>
1911          *  As Orekit doesn't implement the retrograde orbit, I is always set to +1.
1912          *  But for the sake of consistency with the theory, the retrograde factor
1913          *  has been kept in the formulas.
1914          *  </p>
1915          */
1916         private static final int I = 1;
1917 
1918         /** Central body rotating frame. */
1919         private final Frame bodyFrame;
1920 
1921         /** Maximal order to consider for short periodics m-daily tesseral harmonics potential. */
1922         private final int maxOrderMdailyTesseralSP;
1923 
1924         /** Flag to take into account only M-dailies harmonic tesserals for short periodic perturbations.  */
1925         private final boolean mDailiesOnly;
1926 
1927         /** List of non resonant orders with j != 0. */
1928         private final SortedMap<Integer, List<Integer> > nonResOrders;
1929 
1930         /** Maximum value for m index. */
1931         private final int mMax;
1932 
1933         /** Maximum value for j index. */
1934         private final int jMax;
1935 
1936         /** Number of points used in the interpolation process. */
1937         private final int interpolationPoints;
1938 
1939         /** All coefficients slots. */
1940         private final transient FieldTimeSpanMap<FieldSlot<T>, T> slots;
1941 
1942         /** Constructor.
1943          * @param bodyFrame central body rotating frame
1944          * @param maxOrderMdailyTesseralSP maximal order to consider for short periodics m-daily tesseral harmonics potential
1945          * @param mDailiesOnly flag to take into account only M-dailies harmonic tesserals for short periodic perturbations
1946          * @param nonResOrders lst of non resonant orders with j != 0
1947          * @param mMax maximum value for m index
1948          * @param jMax maximum value for j index
1949          * @param interpolationPoints number of points used in the interpolation process
1950          * @param slots all coefficients slots
1951          */
1952         FieldTesseralShortPeriodicCoefficients(final Frame bodyFrame, final int maxOrderMdailyTesseralSP,
1953                                                final boolean mDailiesOnly, final SortedMap<Integer, List<Integer> > nonResOrders,
1954                                                final int mMax, final int jMax, final int interpolationPoints,
1955                                                final FieldTimeSpanMap<FieldSlot<T>, T> slots) {
1956             this.bodyFrame                = bodyFrame;
1957             this.maxOrderMdailyTesseralSP = maxOrderMdailyTesseralSP;
1958             this.mDailiesOnly             = mDailiesOnly;
1959             this.nonResOrders             = nonResOrders;
1960             this.mMax                     = mMax;
1961             this.jMax                     = jMax;
1962             this.interpolationPoints      = interpolationPoints;
1963             this.slots                    = slots;
1964         }
1965 
1966         /** Get the slot valid for some date.
1967          * @param meanStates mean states defining the slot
1968          * @return slot valid at the specified date
1969          */
1970         @SuppressWarnings("unchecked")
1971         public FieldSlot<T> createSlot(final FieldSpacecraftState<T>... meanStates) {
1972             final FieldSlot<T>         slot  = new FieldSlot<>(mMax, jMax, interpolationPoints);
1973             final FieldAbsoluteDate<T> first = meanStates[0].getDate();
1974             final FieldAbsoluteDate<T> last  = meanStates[meanStates.length - 1].getDate();
1975             if (first.compareTo(last) <= 0) {
1976                 slots.addValidAfter(slot, first);
1977             } else {
1978                 slots.addValidBefore(slot, first);
1979             }
1980             return slot;
1981         }
1982 
1983         /** {@inheritDoc} */
1984         @Override
1985         public T[] value(final FieldOrbit<T> meanOrbit) {
1986 
1987             // select the coefficients slot
1988             final FieldSlot<T> slot = slots.get(meanOrbit.getDate());
1989 
1990             // Initialise the short periodic variations
1991             final T[] shortPeriodicVariation = MathArrays.buildArray(meanOrbit.getDate().getField(), 6);
1992 
1993             // Compute only if there is at least one non-resonant tesseral or
1994             // only the m-daily tesseral should be taken into account
1995             if (!nonResOrders.isEmpty() || mDailiesOnly) {
1996 
1997                 //Build an auxiliary object
1998                 final FieldAuxiliaryElements<T> auxiliaryElements = new FieldAuxiliaryElements<>(meanOrbit, I);
1999 
2000                 // Central body rotation angle from equation 2.7.1-(3)(4).
2001                 final FieldTransform<T> t = bodyFrame.getTransformTo(auxiliaryElements.getFrame(), auxiliaryElements.getDate());
2002                 final FieldVector3D<T> xB = t.transformVector(Vector3D.PLUS_I);
2003                 final FieldVector3D<T> yB = t.transformVector(Vector3D.PLUS_J);
2004                 final FieldVector3D<T>  f = auxiliaryElements.getVectorF();
2005                 final FieldVector3D<T>  g = auxiliaryElements.getVectorG();
2006                 final T currentTheta = FastMath.atan2(f.dotProduct(yB).negate().add(g.dotProduct(xB).multiply(I)),
2007                                                       f.dotProduct(xB).add(g.dotProduct(yB).multiply(I)));
2008 
2009                 //Add the m-daily contribution
2010                 for (int m = 1; m <= maxOrderMdailyTesseralSP; m++) {
2011                     // Phase angle
2012                     final T jlMmt              = currentTheta.multiply(-m);
2013                     final FieldSinCos<T> scPhi = FastMath.sinCos(jlMmt);
2014                     final T sinPhi             = scPhi.sin();
2015                     final T cosPhi             = scPhi.cos();
2016 
2017                     // compute contribution for each element
2018                     final T[] c = slot.getCijm(0, m, meanOrbit.getDate());
2019                     final T[] s = slot.getSijm(0, m, meanOrbit.getDate());
2020                     for (int i = 0; i < 6; i++) {
2021                         shortPeriodicVariation[i] = shortPeriodicVariation[i].add(c[i].multiply(cosPhi).add(s[i].multiply(sinPhi)));
2022                     }
2023                 }
2024 
2025                 // loop through all non-resonant (j,m) pairs
2026                 for (final Map.Entry<Integer, List<Integer>> entry : nonResOrders.entrySet()) {
2027                     final int           m     = entry.getKey();
2028                     final List<Integer> listJ = entry.getValue();
2029 
2030                     for (int j : listJ) {
2031                         // Phase angle
2032                         final T jlMmt              = meanOrbit.getLM().multiply(j).subtract(currentTheta.multiply(m));
2033                         final FieldSinCos<T> scPhi = FastMath.sinCos(jlMmt);
2034                         final T sinPhi             = scPhi.sin();
2035                         final T cosPhi             = scPhi.cos();
2036 
2037                         // compute contribution for each element
2038                         final T[] c = slot.getCijm(j, m, meanOrbit.getDate());
2039                         final T[] s = slot.getSijm(j, m, meanOrbit.getDate());
2040                         for (int i = 0; i < 6; i++) {
2041                             shortPeriodicVariation[i] = shortPeriodicVariation[i].add(c[i].multiply(cosPhi).add(s[i].multiply(sinPhi)));
2042                         }
2043 
2044                     }
2045                 }
2046             }
2047 
2048             return shortPeriodicVariation;
2049 
2050         }
2051 
2052         /** {@inheritDoc} */
2053         @Override
2054         public String getCoefficientsKeyPrefix() {
2055             return DSSTTesseral.SHORT_PERIOD_PREFIX;
2056         }
2057 
2058         /** {@inheritDoc}
2059          * <p>
2060          * For tesseral terms contributions,there are maxOrderMdailyTesseralSP
2061          * m-daily cMm coefficients, maxOrderMdailyTesseralSP m-daily sMm
2062          * coefficients, nbNonResonant cjm coefficients and nbNonResonant
2063          * sjm coefficients, where maxOrderMdailyTesseralSP and nbNonResonant both
2064          * depend on the orbit. The j index is the integer multiplier for the true
2065          * longitude argument and the m index is the integer multiplier for m-dailies.
2066          * </p>
2067          */
2068         @Override
2069         public Map<String, T[]> getCoefficients(final FieldAbsoluteDate<T> date, final Set<String> selected) {
2070 
2071             // select the coefficients slot
2072             final FieldSlot<T> slot = slots.get(date);
2073 
2074             if (!nonResOrders.isEmpty() || mDailiesOnly) {
2075                 final Map<String, T[]> coefficients =
2076                                 new HashMap<String, T[]>(12 * maxOrderMdailyTesseralSP +
2077                                                          12 * nonResOrders.size());
2078 
2079                 for (int m = 1; m <= maxOrderMdailyTesseralSP; m++) {
2080                     storeIfSelected(coefficients, selected, slot.getCijm(0, m, date), DSSTTesseral.CM_COEFFICIENTS, m);
2081                     storeIfSelected(coefficients, selected, slot.getSijm(0, m, date), DSSTTesseral.SM_COEFFICIENTS, m);
2082                 }
2083 
2084                 for (final Map.Entry<Integer, List<Integer>> entry : nonResOrders.entrySet()) {
2085                     final int           m     = entry.getKey();
2086                     final List<Integer> listJ = entry.getValue();
2087 
2088                     for (int j : listJ) {
2089                         for (int i = 0; i < 6; ++i) {
2090                             storeIfSelected(coefficients, selected, slot.getCijm(j, m, date), "c", j, m);
2091                             storeIfSelected(coefficients, selected, slot.getSijm(j, m, date), "s", j, m);
2092                         }
2093                     }
2094                 }
2095 
2096                 return coefficients;
2097 
2098             } else {
2099                 return Collections.emptyMap();
2100             }
2101 
2102         }
2103 
2104         /** Put a coefficient in a map if selected.
2105          * @param map map to populate
2106          * @param selected set of coefficients that should be put in the map
2107          * (empty set means all coefficients are selected)
2108          * @param value coefficient value
2109          * @param id coefficient identifier
2110          * @param indices list of coefficient indices
2111          */
2112         private void storeIfSelected(final Map<String, T[]> map, final Set<String> selected,
2113                                      final T[] value, final String id, final int... indices) {
2114             final StringBuilder keyBuilder = new StringBuilder(getCoefficientsKeyPrefix());
2115             keyBuilder.append(id);
2116             for (int index : indices) {
2117                 keyBuilder.append('[').append(index).append(']');
2118             }
2119             final String key = keyBuilder.toString();
2120             if (selected.isEmpty() || selected.contains(key)) {
2121                 map.put(key, value);
2122             }
2123         }
2124     }
2125 
2126     /** Coefficients valid for one time slot. */
2127     private static class Slot {
2128 
2129         /** The coefficients C<sub>i</sub><sup>j</sup><sup>m</sup>.
2130          * <p>
2131          * The index order is cijm[m][j][i] <br/>
2132          * i corresponds to the equinoctial element, as follows: <br/>
2133          * - i=0 for a <br/>
2134          * - i=1 for k <br/>
2135          * - i=2 for h <br/>
2136          * - i=3 for q <br/>
2137          * - i=4 for p <br/>
2138          * - i=5 for λ <br/>
2139          * </p>
2140          */
2141         private final ShortPeriodicsInterpolatedCoefficient[][] cijm;
2142 
2143         /** The coefficients S<sub>i</sub><sup>j</sup><sup>m</sup>.
2144          * <p>
2145          * The index order is sijm[m][j][i] <br/>
2146          * i corresponds to the equinoctial element, as follows: <br/>
2147          * - i=0 for a <br/>
2148          * - i=1 for k <br/>
2149          * - i=2 for h <br/>
2150          * - i=3 for q <br/>
2151          * - i=4 for p <br/>
2152          * - i=5 for λ <br/>
2153          * </p>
2154          */
2155         private final ShortPeriodicsInterpolatedCoefficient[][] sijm;
2156 
2157         /** Simple constructor.
2158          *  @param mMax maximum value for m index
2159          *  @param jMax maximum value for j index
2160          *  @param interpolationPoints number of points used in the interpolation process
2161          */
2162         Slot(final int mMax, final int jMax, final int interpolationPoints) {
2163 
2164             final int rows    = mMax + 1;
2165             final int columns = 2 * jMax + 1;
2166             cijm = new ShortPeriodicsInterpolatedCoefficient[rows][columns];
2167             sijm = new ShortPeriodicsInterpolatedCoefficient[rows][columns];
2168             for (int m = 1; m <= mMax; m++) {
2169                 for (int j = -jMax; j <= jMax; j++) {
2170                     cijm[m][j + jMax] = new ShortPeriodicsInterpolatedCoefficient(interpolationPoints);
2171                     sijm[m][j + jMax] = new ShortPeriodicsInterpolatedCoefficient(interpolationPoints);
2172                 }
2173             }
2174 
2175         }
2176 
2177         /** Get C<sub>i</sub><sup>j</sup><sup>m</sup>.
2178          *
2179          * @param j j index
2180          * @param m m index
2181          * @param date the date
2182          * @return C<sub>i</sub><sup>j</sup><sup>m</sup>
2183          */
2184         double[] getCijm(final int j, final int m, final AbsoluteDate date) {
2185             final int jMax = (cijm[m].length - 1) / 2;
2186             return cijm[m][j + jMax].value(date);
2187         }
2188 
2189         /** Get S<sub>i</sub><sup>j</sup><sup>m</sup>.
2190          *
2191          * @param j j index
2192          * @param m m index
2193          * @param date the date
2194          * @return S<sub>i</sub><sup>j</sup><sup>m</sup>
2195          */
2196         double[] getSijm(final int j, final int m, final AbsoluteDate date) {
2197             final int jMax = (cijm[m].length - 1) / 2;
2198             return sijm[m][j + jMax].value(date);
2199         }
2200 
2201     }
2202 
2203     /** Coefficients valid for one time slot. */
2204     private static class FieldSlot <T extends CalculusFieldElement<T>> {
2205 
2206         /** The coefficients C<sub>i</sub><sup>j</sup><sup>m</sup>.
2207          * <p>
2208          * The index order is cijm[m][j][i] <br/>
2209          * i corresponds to the equinoctial element, as follows: <br/>
2210          * - i=0 for a <br/>
2211          * - i=1 for k <br/>
2212          * - i=2 for h <br/>
2213          * - i=3 for q <br/>
2214          * - i=4 for p <br/>
2215          * - i=5 for λ <br/>
2216          * </p>
2217          */
2218         private final FieldShortPeriodicsInterpolatedCoefficient<T>[][] cijm;
2219 
2220         /** The coefficients S<sub>i</sub><sup>j</sup><sup>m</sup>.
2221          * <p>
2222          * The index order is sijm[m][j][i] <br/>
2223          * i corresponds to the equinoctial element, as follows: <br/>
2224          * - i=0 for a <br/>
2225          * - i=1 for k <br/>
2226          * - i=2 for h <br/>
2227          * - i=3 for q <br/>
2228          * - i=4 for p <br/>
2229          * - i=5 for λ <br/>
2230          * </p>
2231          */
2232         private final FieldShortPeriodicsInterpolatedCoefficient<T>[][] sijm;
2233 
2234         /** Simple constructor.
2235          *  @param mMax maximum value for m index
2236          *  @param jMax maximum value for j index
2237          *  @param interpolationPoints number of points used in the interpolation process
2238          */
2239         @SuppressWarnings("unchecked")
2240         FieldSlot(final int mMax, final int jMax, final int interpolationPoints) {
2241 
2242             final int rows    = mMax + 1;
2243             final int columns = 2 * jMax + 1;
2244             cijm = (FieldShortPeriodicsInterpolatedCoefficient<T>[][]) Array.newInstance(FieldShortPeriodicsInterpolatedCoefficient.class, rows, columns);
2245             sijm = (FieldShortPeriodicsInterpolatedCoefficient<T>[][]) Array.newInstance(FieldShortPeriodicsInterpolatedCoefficient.class, rows, columns);
2246             for (int m = 1; m <= mMax; m++) {
2247                 for (int j = -jMax; j <= jMax; j++) {
2248                     cijm[m][j + jMax] = new FieldShortPeriodicsInterpolatedCoefficient<>(interpolationPoints);
2249                     sijm[m][j + jMax] = new FieldShortPeriodicsInterpolatedCoefficient<>(interpolationPoints);
2250                 }
2251             }
2252 
2253         }
2254 
2255         /** Get C<sub>i</sub><sup>j</sup><sup>m</sup>.
2256          *
2257          * @param j j index
2258          * @param m m index
2259          * @param date the date
2260          * @return C<sub>i</sub><sup>j</sup><sup>m</sup>
2261          */
2262         T[] getCijm(final int j, final int m, final FieldAbsoluteDate<T> date) {
2263             final int jMax = (cijm[m].length - 1) / 2;
2264             return cijm[m][j + jMax].value(date);
2265         }
2266 
2267         /** Get S<sub>i</sub><sup>j</sup><sup>m</sup>.
2268          *
2269          * @param j j index
2270          * @param m m index
2271          * @param date the date
2272          * @return S<sub>i</sub><sup>j</sup><sup>m</sup>
2273          */
2274         T[] getSijm(final int j, final int m, final FieldAbsoluteDate<T> date) {
2275             final int jMax = (cijm[m].length - 1) / 2;
2276             return sijm[m][j + jMax].value(date);
2277         }
2278 
2279     }
2280 
2281     /** Compute potential and potential derivatives with respect to orbital parameters.
2282      *  <p>The following elements are computed from expression 3.3 - (4).
2283      *  <pre>
2284      *  dU / da
2285      *  dU / dh
2286      *  dU / dk
2287      *  dU / dλ
2288      *  dU / dα
2289      *  dU / dβ
2290      *  dU / dγ
2291      *  </pre>
2292      *  </p>
2293      */
2294     private class UAnddU {
2295 
2296         /** dU / da. */
2297         private  double dUda;
2298 
2299         /** dU / dk. */
2300         private double dUdk;
2301 
2302         /** dU / dh. */
2303         private double dUdh;
2304 
2305         /** dU / dl. */
2306         private double dUdl;
2307 
2308         /** dU / dAlpha. */
2309         private double dUdAl;
2310 
2311         /** dU / dBeta. */
2312         private double dUdBe;
2313 
2314         /** dU / dGamma. */
2315         private double dUdGa;
2316 
2317         /** Simple constuctor.
2318          * @param date current date
2319          * @param context container for attributes
2320          * @param hansen hansen objects
2321          */
2322         UAnddU(final AbsoluteDate date, final DSSTTesseralContext context, final HansenObjects hansen) {
2323 
2324             // Auxiliary elements related to the current orbit
2325             final AuxiliaryElements auxiliaryElements = context.getAuxiliaryElements();
2326 
2327             // Potential derivatives
2328             dUda  = 0.;
2329             dUdh  = 0.;
2330             dUdk  = 0.;
2331             dUdl  = 0.;
2332             dUdAl = 0.;
2333             dUdBe = 0.;
2334             dUdGa = 0.;
2335 
2336             // Compute only if there is at least one resonant tesseral
2337             if (!resOrders.isEmpty()) {
2338                 // Gmsj and Hmsj polynomials
2339                 final GHmsjPolynomials ghMSJ = new GHmsjPolynomials(auxiliaryElements.getK(), auxiliaryElements.getH(), auxiliaryElements.getAlpha(), auxiliaryElements.getBeta(), I);
2340 
2341                 // GAMMAmns function
2342                 final GammaMnsFunction gammaMNS = new GammaMnsFunction(maxDegree, auxiliaryElements.getGamma(), I);
2343 
2344                 // R / a up to power degree
2345                 final double[] roaPow = new double[maxDegree + 1];
2346                 roaPow[0] = 1.;
2347                 for (int i = 1; i <= maxDegree; i++) {
2348                     roaPow[i] = context.getRoa() * roaPow[i - 1];
2349                 }
2350 
2351                 // SUM over resonant terms {j,m}
2352                 for (int m : resOrders) {
2353 
2354                     // Resonant index for the current resonant order
2355                     final int j = FastMath.max(1, (int) FastMath.round(context.getRatio() * m));
2356 
2357                     // Phase angle
2358                     final double jlMmt  = j * auxiliaryElements.getLM() - m * context.getTheta();
2359                     final SinCos scPhi  = FastMath.sinCos(jlMmt);
2360                     final double sinPhi = scPhi.sin();
2361                     final double cosPhi = scPhi.cos();
2362 
2363                     // Potential derivatives components for a given resonant pair {j,m}
2364                     double dUdaCos  = 0.;
2365                     double dUdaSin  = 0.;
2366                     double dUdhCos  = 0.;
2367                     double dUdhSin  = 0.;
2368                     double dUdkCos  = 0.;
2369                     double dUdkSin  = 0.;
2370                     double dUdlCos  = 0.;
2371                     double dUdlSin  = 0.;
2372                     double dUdAlCos = 0.;
2373                     double dUdAlSin = 0.;
2374                     double dUdBeCos = 0.;
2375                     double dUdBeSin = 0.;
2376                     double dUdGaCos = 0.;
2377                     double dUdGaSin = 0.;
2378 
2379                     // s-SUM from -sMin to sMax
2380                     final int sMin = FastMath.min(maxEccPow - j, maxDegree);
2381                     final int sMax = FastMath.min(maxEccPow + j, maxDegree);
2382                     for (int s = 0; s <= sMax; s++) {
2383 
2384                         //Compute the initial values for Hansen coefficients using newComb operators
2385                         hansen.computeHansenObjectsInitValues(context, s + maxDegree, j);
2386 
2387                         // n-SUM for s positive
2388                         final double[][] nSumSpos = computeNSum(date, j, m, s, maxDegree,
2389                                                                 roaPow, ghMSJ, gammaMNS, context, hansen);
2390                         dUdaCos  += nSumSpos[0][0];
2391                         dUdaSin  += nSumSpos[0][1];
2392                         dUdhCos  += nSumSpos[1][0];
2393                         dUdhSin  += nSumSpos[1][1];
2394                         dUdkCos  += nSumSpos[2][0];
2395                         dUdkSin  += nSumSpos[2][1];
2396                         dUdlCos  += nSumSpos[3][0];
2397                         dUdlSin  += nSumSpos[3][1];
2398                         dUdAlCos += nSumSpos[4][0];
2399                         dUdAlSin += nSumSpos[4][1];
2400                         dUdBeCos += nSumSpos[5][0];
2401                         dUdBeSin += nSumSpos[5][1];
2402                         dUdGaCos += nSumSpos[6][0];
2403                         dUdGaSin += nSumSpos[6][1];
2404 
2405                         // n-SUM for s negative
2406                         if (s > 0 && s <= sMin) {
2407                             //Compute the initial values for Hansen coefficients using newComb operators
2408                             hansen.computeHansenObjectsInitValues(context, maxDegree - s, j);
2409 
2410                             final double[][] nSumSneg = computeNSum(date, j, m, -s, maxDegree,
2411                                                                     roaPow, ghMSJ, gammaMNS, context, hansen);
2412                             dUdaCos  += nSumSneg[0][0];
2413                             dUdaSin  += nSumSneg[0][1];
2414                             dUdhCos  += nSumSneg[1][0];
2415                             dUdhSin  += nSumSneg[1][1];
2416                             dUdkCos  += nSumSneg[2][0];
2417                             dUdkSin  += nSumSneg[2][1];
2418                             dUdlCos  += nSumSneg[3][0];
2419                             dUdlSin  += nSumSneg[3][1];
2420                             dUdAlCos += nSumSneg[4][0];
2421                             dUdAlSin += nSumSneg[4][1];
2422                             dUdBeCos += nSumSneg[5][0];
2423                             dUdBeSin += nSumSneg[5][1];
2424                             dUdGaCos += nSumSneg[6][0];
2425                             dUdGaSin += nSumSneg[6][1];
2426                         }
2427                     }
2428 
2429                     // Assembly of potential derivatives componants
2430                     dUda  += cosPhi * dUdaCos  + sinPhi * dUdaSin;
2431                     dUdh  += cosPhi * dUdhCos  + sinPhi * dUdhSin;
2432                     dUdk  += cosPhi * dUdkCos  + sinPhi * dUdkSin;
2433                     dUdl  += cosPhi * dUdlCos  + sinPhi * dUdlSin;
2434                     dUdAl += cosPhi * dUdAlCos + sinPhi * dUdAlSin;
2435                     dUdBe += cosPhi * dUdBeCos + sinPhi * dUdBeSin;
2436                     dUdGa += cosPhi * dUdGaCos + sinPhi * dUdGaSin;
2437                 }
2438 
2439                 this.dUda  = dUda * (-context.getMoa() / auxiliaryElements.getSma());
2440                 this.dUdh  = dUdh * context.getMoa();
2441                 this.dUdk  = dUdk * context.getMoa();
2442                 this.dUdl  = dUdl * context.getMoa();
2443                 this.dUdAl = dUdAl * context.getMoa();
2444                 this.dUdBe = dUdBe * context.getMoa();
2445                 this.dUdGa = dUdGa * context.getMoa();
2446             }
2447 
2448         }
2449 
2450         /** Return value of dU / da.
2451          * @return dUda
2452          */
2453         public double getdUda() {
2454             return dUda;
2455         }
2456 
2457         /** Return value of dU / dk.
2458          * @return dUdk
2459          */
2460         public double getdUdk() {
2461             return dUdk;
2462         }
2463 
2464         /** Return value of dU / dh.
2465          * @return dUdh
2466          */
2467         public double getdUdh() {
2468             return dUdh;
2469         }
2470 
2471         /** Return value of dU / dl.
2472          * @return dUdl
2473          */
2474         public double getdUdl() {
2475             return dUdl;
2476         }
2477 
2478         /** Return value of dU / dAlpha.
2479          * @return dUdAl
2480          */
2481         public double getdUdAl() {
2482             return dUdAl;
2483         }
2484 
2485         /** Return value of dU / dBeta.
2486          * @return dUdBe
2487          */
2488         public double getdUdBe() {
2489             return dUdBe;
2490         }
2491 
2492         /** Return value of dU / dGamma.
2493          * @return dUdGa
2494          */
2495         public double getdUdGa() {
2496             return dUdGa;
2497         }
2498 
2499     }
2500 
2501     /**  Computes the potential U derivatives.
2502      *  <p>The following elements are computed from expression 3.3 - (4).
2503      *  <pre>
2504      *  dU / da
2505      *  dU / dh
2506      *  dU / dk
2507      *  dU / dλ
2508      *  dU / dα
2509      *  dU / dβ
2510      *  dU / dγ
2511      *  </pre>
2512      *  </p>
2513      */
2514     private class FieldUAnddU <T extends CalculusFieldElement<T>> {
2515 
2516         /** dU / da. */
2517         private T dUda;
2518 
2519         /** dU / dk. */
2520         private T dUdk;
2521 
2522         /** dU / dh. */
2523         private T dUdh;
2524 
2525         /** dU / dl. */
2526         private T dUdl;
2527 
2528         /** dU / dAlpha. */
2529         private T dUdAl;
2530 
2531         /** dU / dBeta. */
2532         private T dUdBe;
2533 
2534         /** dU / dGamma. */
2535         private T dUdGa;
2536 
2537         /** Simple constuctor.
2538          * @param date current date
2539          * @param context container for attributes
2540          * @param hansen hansen objects
2541          */
2542         FieldUAnddU(final FieldAbsoluteDate<T> date, final FieldDSSTTesseralContext<T> context,
2543                     final FieldHansenObjects<T> hansen) {
2544 
2545             // Auxiliary elements related to the current orbit
2546             final FieldAuxiliaryElements<T> auxiliaryElements = context.getFieldAuxiliaryElements();
2547 
2548             // Zero for initialization
2549             final Field<T> field = date.getField();
2550             final T zero = field.getZero();
2551 
2552             // Potential derivatives
2553             dUda  = zero;
2554             dUdh  = zero;
2555             dUdk  = zero;
2556             dUdl  = zero;
2557             dUdAl = zero;
2558             dUdBe = zero;
2559             dUdGa = zero;
2560 
2561             // Compute only if there is at least one resonant tesseral
2562             if (!resOrders.isEmpty()) {
2563                 // Gmsj and Hmsj polynomials
2564                 final FieldGHmsjPolynomials<T> ghMSJ = new FieldGHmsjPolynomials<>(auxiliaryElements.getK(), auxiliaryElements.getH(), auxiliaryElements.getAlpha(), auxiliaryElements.getBeta(), I, field);
2565 
2566                 // GAMMAmns function
2567                 final FieldGammaMnsFunction<T> gammaMNS = new FieldGammaMnsFunction<>(maxDegree, auxiliaryElements.getGamma(), I, field);
2568 
2569                 // R / a up to power degree
2570                 final T[] roaPow = MathArrays.buildArray(field, maxDegree + 1);
2571                 roaPow[0] = zero.add(1.);
2572                 for (int i = 1; i <= maxDegree; i++) {
2573                     roaPow[i] = roaPow[i - 1].multiply(context.getRoa());
2574                 }
2575 
2576                 // SUM over resonant terms {j,m}
2577                 for (int m : resOrders) {
2578 
2579                     // Resonant index for the current resonant order
2580                     final int j = FastMath.max(1, (int) FastMath.round(context.getRatio().multiply(m)));
2581 
2582                     // Phase angle
2583                     final T jlMmt              = auxiliaryElements.getLM().multiply(j).subtract(context.getTheta().multiply(m));
2584                     final FieldSinCos<T> scPhi = FastMath.sinCos(jlMmt);
2585                     final T sinPhi             = scPhi.sin();
2586                     final T cosPhi             = scPhi.cos();
2587 
2588                     // Potential derivatives components for a given resonant pair {j,m}
2589                     T dUdaCos  = zero;
2590                     T dUdaSin  = zero;
2591                     T dUdhCos  = zero;
2592                     T dUdhSin  = zero;
2593                     T dUdkCos  = zero;
2594                     T dUdkSin  = zero;
2595                     T dUdlCos  = zero;
2596                     T dUdlSin  = zero;
2597                     T dUdAlCos = zero;
2598                     T dUdAlSin = zero;
2599                     T dUdBeCos = zero;
2600                     T dUdBeSin = zero;
2601                     T dUdGaCos = zero;
2602                     T dUdGaSin = zero;
2603 
2604                     // s-SUM from -sMin to sMax
2605                     final int sMin = FastMath.min(maxEccPow - j, maxDegree);
2606                     final int sMax = FastMath.min(maxEccPow + j, maxDegree);
2607                     for (int s = 0; s <= sMax; s++) {
2608 
2609                         //Compute the initial values for Hansen coefficients using newComb operators
2610                         hansen.computeHansenObjectsInitValues(context, s + maxDegree, j);
2611 
2612                         // n-SUM for s positive
2613                         final T[][] nSumSpos = computeNSum(date, j, m, s, maxDegree,
2614                                                                 roaPow, ghMSJ, gammaMNS, context, hansen);
2615                         dUdaCos  = dUdaCos.add(nSumSpos[0][0]);
2616                         dUdaSin  = dUdaSin.add(nSumSpos[0][1]);
2617                         dUdhCos  = dUdhCos.add(nSumSpos[1][0]);
2618                         dUdhSin  = dUdhSin.add(nSumSpos[1][1]);
2619                         dUdkCos  = dUdkCos.add(nSumSpos[2][0]);
2620                         dUdkSin  = dUdkSin.add(nSumSpos[2][1]);
2621                         dUdlCos  = dUdlCos.add(nSumSpos[3][0]);
2622                         dUdlSin  = dUdlSin.add(nSumSpos[3][1]);
2623                         dUdAlCos = dUdAlCos.add(nSumSpos[4][0]);
2624                         dUdAlSin = dUdAlSin.add(nSumSpos[4][1]);
2625                         dUdBeCos = dUdBeCos.add(nSumSpos[5][0]);
2626                         dUdBeSin = dUdBeSin.add(nSumSpos[5][1]);
2627                         dUdGaCos = dUdGaCos.add(nSumSpos[6][0]);
2628                         dUdGaSin = dUdGaSin.add(nSumSpos[6][1]);
2629 
2630                         // n-SUM for s negative
2631                         if (s > 0 && s <= sMin) {
2632                             //Compute the initial values for Hansen coefficients using newComb operators
2633                             hansen.computeHansenObjectsInitValues(context, maxDegree - s, j);
2634 
2635                             final T[][] nSumSneg = computeNSum(date, j, m, -s, maxDegree,
2636                                                                     roaPow, ghMSJ, gammaMNS, context, hansen);
2637                             dUdaCos  = dUdaCos.add(nSumSneg[0][0]);
2638                             dUdaSin  = dUdaSin.add(nSumSneg[0][1]);
2639                             dUdhCos  = dUdhCos.add(nSumSneg[1][0]);
2640                             dUdhSin  = dUdhSin.add(nSumSneg[1][1]);
2641                             dUdkCos  = dUdkCos.add(nSumSneg[2][0]);
2642                             dUdkSin  = dUdkSin.add(nSumSneg[2][1]);
2643                             dUdlCos  = dUdlCos.add(nSumSneg[3][0]);
2644                             dUdlSin  = dUdlSin.add(nSumSneg[3][1]);
2645                             dUdAlCos = dUdAlCos.add(nSumSneg[4][0]);
2646                             dUdAlSin = dUdAlSin.add(nSumSneg[4][1]);
2647                             dUdBeCos = dUdBeCos.add(nSumSneg[5][0]);
2648                             dUdBeSin = dUdBeSin.add(nSumSneg[5][1]);
2649                             dUdGaCos = dUdGaCos.add(nSumSneg[6][0]);
2650                             dUdGaSin = dUdGaSin.add(nSumSneg[6][1]);
2651                         }
2652                     }
2653 
2654                     // Assembly of potential derivatives componants
2655                     dUda  = dUda.add(dUdaCos.multiply(cosPhi).add(dUdaSin.multiply(sinPhi)));
2656                     dUdh  = dUdh.add(dUdhCos.multiply(cosPhi).add(dUdhSin.multiply(sinPhi)));
2657                     dUdk  = dUdk.add(dUdkCos.multiply(cosPhi).add(dUdkSin.multiply(sinPhi)));
2658                     dUdl  = dUdl.add(dUdlCos.multiply(cosPhi).add(dUdlSin.multiply(sinPhi)));
2659                     dUdAl = dUdAl.add(dUdAlCos.multiply(cosPhi).add(dUdAlSin.multiply(sinPhi)));
2660                     dUdBe = dUdBe.add(dUdBeCos.multiply(cosPhi).add(dUdBeSin.multiply(sinPhi)));
2661                     dUdGa = dUdGa.add(dUdGaCos.multiply(cosPhi).add(dUdGaSin.multiply(sinPhi)));
2662                 }
2663 
2664                 dUda  =  dUda.multiply(context.getMoa().divide(auxiliaryElements.getSma())).negate();
2665                 dUdh  =  dUdh.multiply(context.getMoa());
2666                 dUdk  =  dUdk.multiply(context.getMoa());
2667                 dUdl  =  dUdl.multiply(context.getMoa());
2668                 dUdAl =  dUdAl.multiply(context.getMoa());
2669                 dUdBe =  dUdBe.multiply(context.getMoa());
2670                 dUdGa =  dUdGa.multiply(context.getMoa());
2671 
2672             }
2673 
2674         }
2675 
2676         /** Return value of dU / da.
2677          * @return dUda
2678          */
2679         public T getdUda() {
2680             return dUda;
2681         }
2682 
2683         /** Return value of dU / dk.
2684          * @return dUdk
2685          */
2686         public T getdUdk() {
2687             return dUdk;
2688         }
2689 
2690         /** Return value of dU / dh.
2691          * @return dUdh
2692          */
2693         public T getdUdh() {
2694             return dUdh;
2695         }
2696 
2697         /** Return value of dU / dl.
2698          * @return dUdl
2699          */
2700         public T getdUdl() {
2701             return dUdl;
2702         }
2703 
2704         /** Return value of dU / dAlpha.
2705          * @return dUdAl
2706          */
2707         public T getdUdAl() {
2708             return dUdAl;
2709         }
2710 
2711         /** Return value of dU / dBeta.
2712          * @return dUdBe
2713          */
2714         public T getdUdBe() {
2715             return dUdBe;
2716         }
2717 
2718         /** Return value of dU / dGamma.
2719          * @return dUdGa
2720          */
2721         public T getdUdGa() {
2722             return dUdGa;
2723         }
2724 
2725     }
2726 
2727     /** Computes init values of the Hansen Objects. */
2728     private class HansenObjects {
2729 
2730         /** A two dimensional array that contains the objects needed to build the Hansen coefficients. <br/>
2731          * The indexes are s + maxDegree and j */
2732         private HansenTesseralLinear[][] hansenObjects;
2733 
2734         /** Simple constructor.
2735          * @param ratio Ratio of satellite period to central body rotation period
2736          * @param type type of the elements used during the propagation
2737          */
2738         HansenObjects(final double ratio,
2739                       final PropagationType type) {
2740 
2741             //Allocate the two dimensional array
2742             final int rows     = 2 * maxDegree + 1;
2743             final int columns  = maxFrequencyShortPeriodics + 1;
2744             this.hansenObjects = new HansenTesseralLinear[rows][columns];
2745 
2746             switch (type) {
2747                 case MEAN:
2748                     // loop through the resonant orders
2749                     for (int m : resOrders) {
2750                         //Compute the corresponding j term
2751                         final int j = FastMath.max(1, (int) FastMath.round(ratio * m));
2752 
2753                         //Compute the sMin and sMax values
2754                         final int sMin = FastMath.min(maxEccPow - j, maxDegree);
2755                         final int sMax = FastMath.min(maxEccPow + j, maxDegree);
2756 
2757                         //loop through the s values
2758                         for (int s = 0; s <= sMax; s++) {
2759                             //Compute the n0 value
2760                             final int n0 = FastMath.max(FastMath.max(2, m), s);
2761 
2762                             //Create the object for the pair j, s
2763                             this.hansenObjects[s + maxDegree][j] = new HansenTesseralLinear(maxDegree, s, j, n0, maxHansen);
2764 
2765                             if (s > 0 && s <= sMin) {
2766                                 //Also create the object for the pair j, -s
2767                                 this.hansenObjects[maxDegree - s][j] =  new HansenTesseralLinear(maxDegree, -s, j, n0, maxHansen);
2768                             }
2769                         }
2770                     }
2771                     break;
2772 
2773                 case OSCULATING:
2774                     // create all objects
2775                     for (int j = 0; j <= maxFrequencyShortPeriodics; j++) {
2776                         for (int s = -maxDegree; s <= maxDegree; s++) {
2777                             //Compute the n0 value
2778                             final int n0 = FastMath.max(2, FastMath.abs(s));
2779                             this.hansenObjects[s + maxDegree][j] = new HansenTesseralLinear(maxDegree, s, j, n0, maxHansen);
2780                         }
2781                     }
2782                     break;
2783 
2784                 default:
2785                     throw new OrekitInternalError(null);
2786             }
2787 
2788         }
2789 
2790         /** Compute init values for hansen objects.
2791          * @param context container for attributes
2792          * @param rows number of rows of the hansen matrix
2793          * @param columns columns number of columns of the hansen matrix
2794          */
2795         public void computeHansenObjectsInitValues(final DSSTTesseralContext context, final int rows, final int columns) {
2796             hansenObjects[rows][columns].computeInitValues(context.getE2(), context.getChi(), context.getChi2());
2797         }
2798 
2799         /** Get hansen object.
2800          * @return hansenObjects
2801          */
2802         public HansenTesseralLinear[][] getHansenObjects() {
2803             return hansenObjects;
2804         }
2805 
2806     }
2807 
2808     /** Computes init values of the Hansen Objects. */
2809     private class FieldHansenObjects<T extends CalculusFieldElement<T>> {
2810 
2811         /** A two dimensional array that contains the objects needed to build the Hansen coefficients. <br/>
2812          * The indexes are s + maxDegree and j */
2813         private FieldHansenTesseralLinear<T>[][] hansenObjects;
2814 
2815         /** Simple constructor.
2816          * @param ratio Ratio of satellite period to central body rotation period
2817          * @param type type of the elements used during the propagation
2818          */
2819         @SuppressWarnings("unchecked")
2820         FieldHansenObjects(final T ratio,
2821                            final PropagationType type) {
2822 
2823             // Set the maximum power of the eccentricity to use in Hansen coefficient Kernel expansion.
2824             maxHansen = maxEccPow / 2;
2825 
2826             //Allocate the two dimensional array
2827             final int rows     = 2 * maxDegree + 1;
2828             final int columns  = maxFrequencyShortPeriodics + 1;
2829             this.hansenObjects = (FieldHansenTesseralLinear<T>[][]) Array.newInstance(FieldHansenTesseralLinear.class, rows, columns);
2830 
2831             switch (type) {
2832                 case MEAN:
2833                  // loop through the resonant orders
2834                     for (int m : resOrders) {
2835                         //Compute the corresponding j term
2836                         final int j = FastMath.max(1, (int) FastMath.round(ratio.multiply(m)));
2837 
2838                         //Compute the sMin and sMax values
2839                         final int sMin = FastMath.min(maxEccPow - j, maxDegree);
2840                         final int sMax = FastMath.min(maxEccPow + j, maxDegree);
2841 
2842                         //loop through the s values
2843                         for (int s = 0; s <= sMax; s++) {
2844                             //Compute the n0 value
2845                             final int n0 = FastMath.max(FastMath.max(2, m), s);
2846 
2847                             //Create the object for the pair j, s
2848                             this.hansenObjects[s + maxDegree][j] = new FieldHansenTesseralLinear<>(maxDegree, s, j, n0, maxHansen, ratio.getField());
2849 
2850                             if (s > 0 && s <= sMin) {
2851                                 //Also create the object for the pair j, -s
2852                                 this.hansenObjects[maxDegree - s][j] =  new FieldHansenTesseralLinear<>(maxDegree, -s, j, n0, maxHansen, ratio.getField());
2853                             }
2854                         }
2855                     }
2856                     break;
2857 
2858                 case OSCULATING:
2859                     // create all objects
2860                     for (int j = 0; j <= maxFrequencyShortPeriodics; j++) {
2861                         for (int s = -maxDegree; s <= maxDegree; s++) {
2862                             //Compute the n0 value
2863                             final int n0 = FastMath.max(2, FastMath.abs(s));
2864                             this.hansenObjects[s + maxDegree][j] = new FieldHansenTesseralLinear<>(maxDegree, s, j, n0, maxHansen, ratio.getField());
2865                         }
2866                     }
2867                     break;
2868 
2869                 default:
2870                     throw new OrekitInternalError(null);
2871             }
2872 
2873         }
2874 
2875         /** Compute init values for hansen objects.
2876          * @param context container for attributes
2877          * @param rows number of rows of the hansen matrix
2878          * @param columns columns number of columns of the hansen matrix
2879          */
2880         public void computeHansenObjectsInitValues(final FieldDSSTTesseralContext<T> context,
2881                                                    final int rows, final int columns) {
2882             hansenObjects[rows][columns].computeInitValues(context.getE2(), context.getChi(), context.getChi2());
2883         }
2884 
2885         /** Get hansen object.
2886          * @return hansenObjects
2887          */
2888         public FieldHansenTesseralLinear<T>[][] getHansenObjects() {
2889             return hansenObjects;
2890         }
2891 
2892     }
2893 
2894 }