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