1   /* Copyright 2002-2025 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.forces.gravity;
18  
19  
20  import java.util.Collections;
21  import java.util.List;
22  
23  import org.hipparchus.CalculusFieldElement;
24  import org.hipparchus.analysis.differentiation.DerivativeStructure;
25  import org.hipparchus.analysis.differentiation.Gradient;
26  import org.hipparchus.geometry.euclidean.threed.FieldVector3D;
27  import org.hipparchus.geometry.euclidean.threed.SphericalCoordinates;
28  import org.hipparchus.geometry.euclidean.threed.Vector3D;
29  import org.hipparchus.linear.Array2DRowRealMatrix;
30  import org.hipparchus.linear.RealMatrix;
31  import org.hipparchus.util.FastMath;
32  import org.hipparchus.util.MathArrays;
33  import org.orekit.forces.ForceModel;
34  import org.orekit.forces.gravity.potential.NormalizedSphericalHarmonicsProvider;
35  import org.orekit.forces.gravity.potential.NormalizedSphericalHarmonicsProvider.NormalizedSphericalHarmonics;
36  import org.orekit.forces.gravity.potential.TideSystem;
37  import org.orekit.forces.gravity.potential.TideSystemProvider;
38  import org.orekit.frames.FieldStaticTransform;
39  import org.orekit.frames.Frame;
40  import org.orekit.frames.StaticTransform;
41  import org.orekit.propagation.FieldSpacecraftState;
42  import org.orekit.propagation.SpacecraftState;
43  import org.orekit.time.AbsoluteDate;
44  import org.orekit.time.FieldAbsoluteDate;
45  import org.orekit.utils.FieldPVCoordinates;
46  import org.orekit.utils.ParameterDriver;
47  
48  /** This class represents the gravitational field of a celestial body.
49   * <p>
50   * The algorithm implemented in this class has been designed by S. A. Holmes
51   * and W. E. Featherstone from Department of Spatial Sciences, Curtin University
52   * of Technology, Perth, Australia. It is described in their 2002 paper: <a
53   * href="https://www.researchgate.net/publication/226460594_A_unified_approach_to_the_Clenshaw_summation_and_the_recursive_computation_of_very_high_degree_and_order_normalised_associated_Legendre_functions">
54   * A unified approach to he Clenshaw summation and the recursive computation of
55   * very high degree and order normalised associated Legendre functions</a>
56   * (Journal of Geodesy (2002) 76: 279–299).
57   * </p>
58   * <p>
59   * This model directly uses normalized coefficients and stable recursion algorithms
60   * so it is more suited to high degree gravity fields than the classical Cunningham
61   * Droziner models which use un-normalized coefficients.
62   * </p>
63   * <p>
64   * Among the different algorithms presented in Holmes and Featherstone paper, this
65   * class implements the <em>modified forward row method</em>. All recursion coefficients
66   * are precomputed and stored for greater performance. This caching was suggested in the
67   * paper but not used due to the large memory requirements. Since 2002, even low end
68   * computers and mobile devices do have sufficient memory so this caching has become
69   * feasible nowadays.
70   * </p>
71   * @author Luc Maisonobe
72   * @since 6.0
73   */
74  
75  public class HolmesFeatherstoneAttractionModel implements ForceModel, TideSystemProvider {
76  
77      /** Exponent scaling to avoid floating point overflow.
78       * <p>The paper uses 10^280, we prefer a power of two to preserve accuracy thanks to
79       * {@link FastMath#scalb(double, int)}, so we use 2^930 which has the same order of magnitude.
80       */
81      private static final int SCALING = 930;
82  
83      /** Central attraction scaling factor.
84       * <p>
85       * We use a power of 2 to avoid numeric noise introduction
86       * in the multiplications/divisions sequences.
87       * </p>
88       */
89      private static final double MU_SCALE = FastMath.scalb(1.0, 32);
90  
91      /** Driver for gravitational parameter. */
92      private final ParameterDriver gmParameterDriver;
93  
94      /** Provider for the spherical harmonics. */
95      private final NormalizedSphericalHarmonicsProvider provider;
96  
97      /** Rotating body. */
98      private final Frame bodyFrame;
99  
100     /** Recursion coefficients g<sub>n,m</sub>/√j. */
101     private final double[] gnmOj;
102 
103     /** Recursion coefficients h<sub>n,m</sub>/√j. */
104     private final double[] hnmOj;
105 
106     /** Recursion coefficients e<sub>n,m</sub>. */
107     private final double[] enm;
108 
109     /** Scaled sectorial Pbar<sub>m,m</sub>/u<sup>m</sup> &times; 2<sup>-SCALING</sup>. */
110     private final double[] sectorial;
111 
112     /** Creates a new instance.
113      * @param centralBodyFrame rotating body frame
114      * @param provider provider for spherical harmonics
115      * @since 6.0
116      */
117     public HolmesFeatherstoneAttractionModel(final Frame centralBodyFrame,
118                                              final NormalizedSphericalHarmonicsProvider provider) {
119 
120         gmParameterDriver = new ParameterDriver(NewtonianAttraction.CENTRAL_ATTRACTION_COEFFICIENT,
121                                                 provider.getMu(), MU_SCALE, 0.0, Double.POSITIVE_INFINITY);
122 
123         this.provider  = provider;
124         this.bodyFrame = centralBodyFrame;
125 
126         // the pre-computed arrays hold coefficients from triangular arrays in a single
127         // storing neither diagonal elements (n = m) nor the non-diagonal element n=1, m=0
128         final int degree = provider.getMaxDegree();
129         final int size = FastMath.max(0, degree * (degree + 1) / 2 - 1);
130         gnmOj = new double[size];
131         hnmOj = new double[size];
132         enm   = new double[size];
133 
134         // pre-compute the recursion coefficients corresponding to equations 19 and 22
135         // from Holmes and Featherstone paper
136         // for cache efficiency, elements are stored in the same order they will be used
137         // later on, i.e. from rightmost column to leftmost column
138         int index = 0;
139         for (int m = degree; m >= 0; --m) {
140             final int j = (m == 0) ? 2 : 1;
141             for (int n = FastMath.max(2, m + 1); n <= degree; ++n) {
142                 final double f = (n - m) * (n + m + 1);
143                 gnmOj[index] = 2 * (m + 1) / FastMath.sqrt(j * f);
144                 hnmOj[index] = FastMath.sqrt((n + m + 2) * (n - m - 1) / (j * f));
145                 enm[index]   = FastMath.sqrt(f / j);
146                 ++index;
147             }
148         }
149 
150         // scaled sectorial terms corresponding to equation 28 in Holmes and Featherstone paper
151         sectorial    = new double[degree + 1];
152         sectorial[0] = FastMath.scalb(1.0, -SCALING);
153         if (degree > 0) {
154             sectorial[1] = FastMath.sqrt(3) * sectorial[0];
155         }
156         for (int m = 2; m < sectorial.length; ++m) {
157             sectorial[m] = FastMath.sqrt((2 * m + 1) / (2.0 * m)) * sectorial[m - 1];
158         }
159 
160     }
161 
162     /** {@inheritDoc} */
163     @Override
164     public boolean dependsOnPositionOnly() {
165         return true;
166     }
167 
168     /** {@inheritDoc} */
169     public TideSystem getTideSystem() {
170         return provider.getTideSystem();
171     }
172 
173     /** Get the central attraction coefficient μ.
174      * @return mu central attraction coefficient (m³/s²),
175      * will throw an exception if gm PDriver has several
176      * values driven (in this case the method
177      * {@link #getMu(AbsoluteDate)} must be used.
178      */
179     public double getMu() {
180         return gmParameterDriver.getValue();
181     }
182 
183     /** Get the central attraction coefficient μ.
184      * @param date date at which mu wants to be known
185      * @return mu central attraction coefficient (m³/s²)
186      */
187     public double getMu(final AbsoluteDate date) {
188         return gmParameterDriver.getValue(date);
189     }
190 
191     /** Compute the value of the gravity field.
192      * @param date current date
193      * @param position position at which gravity field is desired in body frame
194      * @param mu central attraction coefficient to use
195      * @return value of the gravity field (central and non-central parts summed together)
196      */
197     public double value(final AbsoluteDate date, final Vector3D position,
198                         final double mu) {
199         return mu / position.getNorm() + nonCentralPart(date, position, mu);
200     }
201 
202     /** Compute the non-central part of the gravity field.
203      * @param date current date
204      * @param position position at which gravity field is desired in body frame
205      * @param mu central attraction coefficient to use
206      * @return value of the non-central part of the gravity field
207      */
208     public double nonCentralPart(final AbsoluteDate date, final Vector3D position, final double mu) {
209 
210         final int degree = provider.getMaxDegree();
211         final int order  = provider.getMaxOrder();
212         final NormalizedSphericalHarmonics harmonics = provider.onDate(date);
213 
214         // allocate the columns for recursion
215         double[] pnm0Plus2 = new double[degree + 1];
216         double[] pnm0Plus1 = new double[degree + 1];
217         double[] pnm0      = new double[degree + 1];
218 
219         // compute polar coordinates
220         final double x    = position.getX();
221         final double y    = position.getY();
222         final double z    = position.getZ();
223         final double x2   = x * x;
224         final double y2   = y * y;
225         final double z2   = z * z;
226         final double rho2 = x2 + y2;
227         final double r2   = rho2 + z2;
228         final double r    = FastMath.sqrt(r2);
229         final double rho  = FastMath.sqrt(rho2);
230         final double t    = z / r;   // cos(theta), where theta is the polar angle
231         final double u    = rho / r; // sin(theta), where theta is the polar angle
232         final double tOu  = z / rho;
233 
234         // compute distance powers
235         final double[] aOrN = createDistancePowersArray(provider.getAe() / r);
236 
237         // compute longitude cosines/sines
238         final double[][] cosSinLambda = createCosSinArrays(x / rho, y / rho);
239 
240         // outer summation over order
241         int    index = 0;
242         double value = 0;
243         for (int m = degree; m >= 0; --m) {
244 
245             // compute tesseral terms without derivatives
246             index = computeTesseral(m, degree, index, t, u, tOu,
247                                     pnm0Plus2, pnm0Plus1, null, pnm0, null, null);
248 
249             if (m <= order) {
250                 // compute contribution of current order to field (equation 5 of the paper)
251 
252                 // inner summation over degree, for fixed order
253                 double sumDegreeS        = 0;
254                 double sumDegreeC        = 0;
255                 for (int n = FastMath.max(2, m); n <= degree; ++n) {
256                     sumDegreeS += pnm0[n] * aOrN[n] * harmonics.getNormalizedSnm(n, m);
257                     sumDegreeC += pnm0[n] * aOrN[n] * harmonics.getNormalizedCnm(n, m);
258                 }
259 
260                 // contribution to outer summation over order
261                 value = value * u + cosSinLambda[1][m] * sumDegreeS + cosSinLambda[0][m] * sumDegreeC;
262 
263             }
264 
265             // rotate the recursion arrays
266             final double[] tmp = pnm0Plus2;
267             pnm0Plus2 = pnm0Plus1;
268             pnm0Plus1 = pnm0;
269             pnm0      = tmp;
270 
271         }
272 
273         // scale back
274         value = FastMath.scalb(value, SCALING);
275 
276         // apply the global mu/r factor
277         return mu * value / r;
278 
279     }
280 
281     /** Compute the gradient of the non-central part of the gravity field.
282      * <p>
283      * If U represents the non-central part of the gravity field,
284      * this method returns the negative gradient of U (-grad(U)).
285      * </p>
286      * @param date current date
287      * @param position position at which gravity field is desired in body frame
288      * @param mu central attraction coefficient to use
289      * @return gradient of the non-central part of the gravity field
290      */
291     public double[] gradient(final AbsoluteDate date, final Vector3D position, final double mu) {
292 
293         final int degree = provider.getMaxDegree();
294         final int order  = provider.getMaxOrder();
295         final NormalizedSphericalHarmonics harmonics = provider.onDate(date);
296 
297         // allocate the columns for recursion
298         double[] pnm0Plus2  = new double[degree + 1];
299         double[] pnm0Plus1  = new double[degree + 1];
300         double[] pnm0       = new double[degree + 1];
301         final double[] pnm1 = new double[degree + 1];
302 
303         // compute polar coordinates
304         final double x    = position.getX();
305         final double y    = position.getY();
306         final double z    = position.getZ();
307         final double x2   = x * x;
308         final double y2   = y * y;
309         final double z2   = z * z;
310         final double r2   = x2 + y2 + z2;
311         final double r    = FastMath.sqrt (r2);
312         final double rho2 = x2 + y2;
313         final double rho  = FastMath.sqrt(rho2);
314         final double t    = z / r;   // cos(theta), where theta is the polar angle
315         final double u    = rho / r; // sin(theta), where theta is the polar angle
316         final double tOu  = z / rho;
317 
318         // compute distance powers
319         final double[] aOrN = createDistancePowersArray(provider.getAe() / r);
320 
321         // compute longitude cosines/sines
322         final double[][] cosSinLambda = createCosSinArrays(x / rho, y / rho);
323 
324         // outer summation over order
325         int    index = 0;
326         double value = 0;
327         final double[] gradient = new double[3];
328         for (int m = degree; m >= 0; --m) {
329 
330             // compute tesseral terms with derivatives
331             index = computeTesseral(m, degree, index, t, u, tOu,
332                                     pnm0Plus2, pnm0Plus1, null, pnm0, pnm1, null);
333 
334             if (m <= order) {
335                 // compute contribution of current order to field (equation 5 of the paper)
336 
337                 // inner summation over degree, for fixed order
338                 double sumDegreeS        = 0;
339                 double sumDegreeC        = 0;
340                 double dSumDegreeSdR     = 0;
341                 double dSumDegreeCdR     = 0;
342                 double dSumDegreeSdTheta = 0;
343                 double dSumDegreeCdTheta = 0;
344                 for (int n = FastMath.max(2, m); n <= degree; ++n) {
345                     final double qSnm  = aOrN[n] * harmonics.getNormalizedSnm(n, m);
346                     final double qCnm  = aOrN[n] * harmonics.getNormalizedCnm(n, m);
347                     final double nOr   = n / r;
348                     final double s0    = pnm0[n] * qSnm;
349                     final double c0    = pnm0[n] * qCnm;
350                     final double s1    = pnm1[n] * qSnm;
351                     final double c1    = pnm1[n] * qCnm;
352                     sumDegreeS        += s0;
353                     sumDegreeC        += c0;
354                     dSumDegreeSdR     -= nOr * s0;
355                     dSumDegreeCdR     -= nOr * c0;
356                     dSumDegreeSdTheta += s1;
357                     dSumDegreeCdTheta += c1;
358                 }
359 
360                 // contribution to outer summation over order
361                 // beware that we need to order gradient using the mathematical conventions
362                 // compliant with the SphericalCoordinates class, so our lambda is its theta
363                 // (and hence at index 1) and our theta is its phi (and hence at index 2)
364                 final double sML = cosSinLambda[1][m];
365                 final double cML = cosSinLambda[0][m];
366                 value            = value       * u + sML * sumDegreeS        + cML * sumDegreeC;
367                 gradient[0]      = gradient[0] * u + sML * dSumDegreeSdR     + cML * dSumDegreeCdR;
368                 gradient[1]      = gradient[1] * u + m * (cML * sumDegreeS - sML * sumDegreeC);
369                 gradient[2]      = gradient[2] * u + sML * dSumDegreeSdTheta + cML * dSumDegreeCdTheta;
370 
371             }
372 
373             // rotate the recursion arrays
374             final double[] tmp = pnm0Plus2;
375             pnm0Plus2 = pnm0Plus1;
376             pnm0Plus1 = pnm0;
377             pnm0      = tmp;
378 
379         }
380 
381         // scale back
382         value       = FastMath.scalb(value,       SCALING);
383         gradient[0] = FastMath.scalb(gradient[0], SCALING);
384         gradient[1] = FastMath.scalb(gradient[1], SCALING);
385         gradient[2] = FastMath.scalb(gradient[2], SCALING);
386 
387         // apply the global mu/r factor
388         final double muOr = mu / r;
389         value            *= muOr;
390         gradient[0]       = muOr * gradient[0] - value / r;
391         gradient[1]      *= muOr;
392         gradient[2]      *= muOr;
393 
394         // convert gradient from spherical to Cartesian
395         return new SphericalCoordinates(position).toCartesianGradient(gradient);
396 
397     }
398 
399     /** Compute the gradient of the non-central part of the gravity field.
400      * <p>
401      * If U represents the non-central part of the gravity field,
402      * this method returns the negative gradient of U (-grad(U)).
403      * </p>
404      * @param date current date
405      * @param position position at which gravity field is desired in body frame
406      * @param mu central attraction coefficient to use
407      * @param <T> type of field used
408      * @return gradient of the non-central part of the gravity field
409      */
410     public <T extends CalculusFieldElement<T>> T[] gradient(final FieldAbsoluteDate<T> date, final FieldVector3D<T> position,
411                                                         final T mu) {
412 
413         final int degree = provider.getMaxDegree();
414         final int order  = provider.getMaxOrder();
415         final NormalizedSphericalHarmonics harmonics = provider.onDate(date.toAbsoluteDate());
416         final T zero = date.getField().getZero();
417         // allocate the columns for recursion
418         T[] pnm0Plus2  = MathArrays.buildArray(date.getField(), degree + 1);
419         T[] pnm0Plus1  = MathArrays.buildArray(date.getField(), degree + 1);
420         T[] pnm0       = MathArrays.buildArray(date.getField(), degree + 1);
421         final T[] pnm1 = MathArrays.buildArray(date.getField(), degree + 1);
422 
423         // compute polar coordinates
424         final T x    = position.getX();
425         final T y    = position.getY();
426         final T z    = position.getZ();
427         final T x2   = x.square();
428         final T y2   = y.square();
429         final T rho2 = x2.add(y2);
430         final T rho  = rho2.sqrt();
431         final T z2   = z.square();
432         final T r2   = rho2.add(z2);
433         final T r    = r2.sqrt();
434         final T t    = z.divide(r);   // cos(theta), where theta is the polar angle
435         final T u    = rho.divide(r); // sin(theta), where theta is the polar angle
436         final T tOu  = z.divide(rho);
437 
438         // compute distance powers
439         final T[] aOrN = createDistancePowersArray(r.reciprocal().multiply(provider.getAe()));
440 
441         // compute longitude cosines/sines
442         final T[][] cosSinLambda = createCosSinArrays(x.divide(rho), y.divide(rho));
443         // outer summation over order
444         int    index = 0;
445         T value = zero;
446         final T[] gradient = MathArrays.buildArray(zero.getField(), 3);
447         for (int m = degree; m >= 0; --m) {
448 
449             // compute tesseral terms with derivatives
450             index = computeTesseral(m, degree, index, t, u, tOu,
451                                     pnm0Plus2, pnm0Plus1, null, pnm0, pnm1, null);
452             if (m <= order) {
453                 // compute contribution of current order to field (equation 5 of the paper)
454 
455                 // inner summation over degree, for fixed order
456                 T sumDegreeS        = zero;
457                 T sumDegreeC        = zero;
458                 T dSumDegreeSdR     = zero;
459                 T dSumDegreeCdR     = zero;
460                 T dSumDegreeSdTheta = zero;
461                 T dSumDegreeCdTheta = zero;
462                 for (int n = FastMath.max(2, m); n <= degree; ++n) {
463                     final T qSnm  = aOrN[n].multiply(harmonics.getNormalizedSnm(n, m));
464                     final T qCnm  = aOrN[n].multiply(harmonics.getNormalizedCnm(n, m));
465                     final T nOr   = r.reciprocal().multiply(n);
466                     final T s0    = pnm0[n].multiply(qSnm);
467                     final T c0    = pnm0[n].multiply(qCnm);
468                     final T s1    = pnm1[n].multiply(qSnm);
469                     final T c1    = pnm1[n].multiply(qCnm);
470                     sumDegreeS        = sumDegreeS       .add(s0);
471                     sumDegreeC        = sumDegreeC       .add(c0);
472                     dSumDegreeSdR     = dSumDegreeSdR    .subtract(nOr.multiply(s0));
473                     dSumDegreeCdR     = dSumDegreeCdR    .subtract(nOr.multiply(c0));
474                     dSumDegreeSdTheta = dSumDegreeSdTheta.add(s1);
475                     dSumDegreeCdTheta = dSumDegreeCdTheta.add(c1);
476                 }
477 
478                 // contribution to outer summation over order
479                 // beware that we need to order gradient using the mathematical conventions
480                 // compliant with the SphericalCoordinates class, so our lambda is its theta
481                 // (and hence at index 1) and our theta is its phi (and hence at index 2)
482                 final T sML = cosSinLambda[1][m];
483                 final T cML = cosSinLambda[0][m];
484                 value            = value      .multiply(u).add(sML.multiply(sumDegreeS   )).add(cML.multiply(sumDegreeC));
485                 gradient[0]      = gradient[0].multiply(u).add(sML.multiply(dSumDegreeSdR)).add(cML.multiply(dSumDegreeCdR));
486                 gradient[1]      = gradient[1].multiply(u).add(cML.multiply(sumDegreeS).subtract(sML.multiply(sumDegreeC)).multiply(m));
487                 gradient[2]      = gradient[2].multiply(u).add(sML.multiply(dSumDegreeSdTheta)).add(cML.multiply(dSumDegreeCdTheta));
488             }
489             // rotate the recursion arrays
490             final T[] tmp = pnm0Plus2;
491             pnm0Plus2 = pnm0Plus1;
492             pnm0Plus1 = pnm0;
493             pnm0      = tmp;
494 
495         }
496         // scale back
497         value       = value.scalb(SCALING);
498         gradient[0] = gradient[0].scalb(SCALING);
499         gradient[1] = gradient[1].scalb(SCALING);
500         gradient[2] = gradient[2].scalb(SCALING);
501 
502         // apply the global mu/r factor
503         final T muOr = r.reciprocal().multiply(mu);
504         value            = value.multiply(muOr);
505         gradient[0]      = muOr.multiply(gradient[0]).subtract(value.divide(r));
506         gradient[1]      = gradient[1].multiply(muOr);
507         gradient[2]      = gradient[2].multiply(muOr);
508 
509         // convert gradient from spherical to Cartesian
510         // Cartesian coordinates
511         // remaining spherical coordinates
512 
513         // intermediate variables
514         final T xPos    = position.getX();
515         final T yPos    = position.getY();
516         final T zPos    = position.getZ();
517         final T rho2Pos = x.square().add(y.square());
518         final T rhoPos  = rho2.sqrt();
519         final T r2Pos   = rho2.add(z.square());
520         final T rPos    = r2Pos.sqrt();
521 
522         final T[][] jacobianPos = MathArrays.buildArray(zero.getField(), 3, 3);
523 
524         // row representing the gradient of r
525         jacobianPos[0][0] = xPos.divide(rPos);
526         jacobianPos[0][1] = yPos.divide(rPos);
527         jacobianPos[0][2] = zPos.divide(rPos);
528 
529         // row representing the gradient of theta
530         jacobianPos[1][0] =  yPos.negate().divide(rho2Pos);
531         jacobianPos[1][1] =  xPos.divide(rho2Pos);
532         // jacobian[1][2] is already set to 0 at allocation time
533 
534         // row representing the gradient of phi
535         final T rhoPosTimesR2Pos = rhoPos.multiply(r2Pos);
536         jacobianPos[2][0] = xPos.multiply(zPos).divide(rhoPosTimesR2Pos);
537         jacobianPos[2][1] = yPos.multiply(zPos).divide(rhoPosTimesR2Pos);
538         jacobianPos[2][2] = rhoPos.negate().divide(r2Pos);
539         final T[] cartGradPos = MathArrays.buildArray(zero.getField(), 3);
540         cartGradPos[0] = gradient[0].multiply(jacobianPos[0][0]).add(gradient[1].multiply(jacobianPos[1][0])).add(gradient[2].multiply(jacobianPos[2][0]));
541         cartGradPos[1] = gradient[0].multiply(jacobianPos[0][1]).add(gradient[1].multiply(jacobianPos[1][1])).add(gradient[2].multiply(jacobianPos[2][1]));
542         cartGradPos[2] = gradient[0].multiply(jacobianPos[0][2])                                      .add(gradient[2].multiply(jacobianPos[2][2]));
543         return cartGradPos;
544 
545     }
546 
547     /** Compute both the gradient and the hessian of the non-central part of the gravity field.
548      * <p>
549      * If U represents the non-central part of the gravity field,
550      * this method returns the negative gradient and hessian of U (-grad(U) and -hessian(U)).
551      * </p>
552      * @param date current date
553      * @param position position at which gravity field is desired in body frame
554      * @param mu central attraction coefficient to use
555      * @return gradient and hessian of the non-central part of the gravity field
556      */
557     private GradientHessian gradientHessian(final AbsoluteDate date, final Vector3D position, final double mu) {
558 
559         final int degree = provider.getMaxDegree();
560         final int order  = provider.getMaxOrder();
561         final NormalizedSphericalHarmonics harmonics = provider.onDate(date);
562 
563         // allocate the columns for recursion
564         double[] pnm0Plus2  = new double[degree + 1];
565         double[] pnm0Plus1  = new double[degree + 1];
566         double[] pnm0       = new double[degree + 1];
567         double[] pnm1Plus1  = new double[degree + 1];
568         double[] pnm1       = new double[degree + 1];
569         final double[] pnm2 = new double[degree + 1];
570 
571         // compute polar coordinates
572         final double x    = position.getX();
573         final double y    = position.getY();
574         final double z    = position.getZ();
575         final double x2   = x * x;
576         final double y2   = y * y;
577         final double z2   = z * z;
578         final double rho2 = x2 + y2;
579         final double rho  = FastMath.sqrt(rho2);
580         final double r2   = rho2 + z2;
581         final double r    = FastMath.sqrt(r2);
582         final double t    = z / r;   // cos(theta), where theta is the polar angle
583         final double u    = rho / r; // sin(theta), where theta is the polar angle
584         final double tOu  = z / rho;
585 
586         // compute distance powers
587         final double[] aOrN = createDistancePowersArray(provider.getAe() / r);
588 
589         // compute longitude cosines/sines
590         final double[][] cosSinLambda = createCosSinArrays(x / rho, y / rho);
591 
592         // outer summation over order
593         int    index = 0;
594         double value = 0;
595         final double[]   gradient = new double[3];
596         final double[][] hessian  = new double[3][3];
597         for (int m = degree; m >= 0; --m) {
598 
599             // compute tesseral terms
600             index = computeTesseral(m, degree, index, t, u, tOu,
601                                     pnm0Plus2, pnm0Plus1, pnm1Plus1, pnm0, pnm1, pnm2);
602 
603             if (m <= order) {
604                 // compute contribution of current order to field (equation 5 of the paper)
605 
606                 // inner summation over degree, for fixed order
607                 double sumDegreeS               = 0;
608                 double sumDegreeC               = 0;
609                 double dSumDegreeSdR            = 0;
610                 double dSumDegreeCdR            = 0;
611                 double dSumDegreeSdTheta        = 0;
612                 double dSumDegreeCdTheta        = 0;
613                 double d2SumDegreeSdRdR         = 0;
614                 double d2SumDegreeSdRdTheta     = 0;
615                 double d2SumDegreeSdThetadTheta = 0;
616                 double d2SumDegreeCdRdR         = 0;
617                 double d2SumDegreeCdRdTheta     = 0;
618                 double d2SumDegreeCdThetadTheta = 0;
619                 for (int n = FastMath.max(2, m); n <= degree; ++n) {
620                     final double qSnm         = aOrN[n] * harmonics.getNormalizedSnm(n, m);
621                     final double qCnm         = aOrN[n] * harmonics.getNormalizedCnm(n, m);
622                     final double nOr          = n / r;
623                     final double nnP1Or2      = nOr * (n + 1) / r;
624                     final double s0           = pnm0[n] * qSnm;
625                     final double c0           = pnm0[n] * qCnm;
626                     final double s1           = pnm1[n] * qSnm;
627                     final double c1           = pnm1[n] * qCnm;
628                     final double s2           = pnm2[n] * qSnm;
629                     final double c2           = pnm2[n] * qCnm;
630                     sumDegreeS               += s0;
631                     sumDegreeC               += c0;
632                     dSumDegreeSdR            -= nOr * s0;
633                     dSumDegreeCdR            -= nOr * c0;
634                     dSumDegreeSdTheta        += s1;
635                     dSumDegreeCdTheta        += c1;
636                     d2SumDegreeSdRdR         += nnP1Or2 * s0;
637                     d2SumDegreeSdRdTheta     -= nOr * s1;
638                     d2SumDegreeSdThetadTheta += s2;
639                     d2SumDegreeCdRdR         += nnP1Or2 * c0;
640                     d2SumDegreeCdRdTheta     -= nOr * c1;
641                     d2SumDegreeCdThetadTheta += c2;
642                 }
643 
644                 // contribution to outer summation over order
645                 final double sML = cosSinLambda[1][m];
646                 final double cML = cosSinLambda[0][m];
647                 value            = value         * u + sML * sumDegreeS + cML * sumDegreeC;
648                 gradient[0]      = gradient[0]   * u + sML * dSumDegreeSdR + cML * dSumDegreeCdR;
649                 gradient[1]      = gradient[1]   * u + m * (cML * sumDegreeS - sML * sumDegreeC);
650                 gradient[2]      = gradient[2]   * u + sML * dSumDegreeSdTheta + cML * dSumDegreeCdTheta;
651                 hessian[0][0]    = hessian[0][0] * u + sML * d2SumDegreeSdRdR + cML * d2SumDegreeCdRdR;
652                 hessian[1][0]    = hessian[1][0] * u + m * (cML * dSumDegreeSdR - sML * dSumDegreeCdR);
653                 hessian[2][0]    = hessian[2][0] * u + sML * d2SumDegreeSdRdTheta + cML * d2SumDegreeCdRdTheta;
654                 hessian[1][1]    = hessian[1][1] * u - m * m * (sML * sumDegreeS + cML * sumDegreeC);
655                 hessian[2][1]    = hessian[2][1] * u + m * (cML * dSumDegreeSdTheta - sML * dSumDegreeCdTheta);
656                 hessian[2][2]    = hessian[2][2] * u + sML * d2SumDegreeSdThetadTheta + cML * d2SumDegreeCdThetadTheta;
657 
658             }
659 
660             // rotate the recursion arrays
661             final double[] tmp0 = pnm0Plus2;
662             pnm0Plus2 = pnm0Plus1;
663             pnm0Plus1 = pnm0;
664             pnm0      = tmp0;
665             final double[] tmp1 = pnm1Plus1;
666             pnm1Plus1 = pnm1;
667             pnm1      = tmp1;
668 
669         }
670 
671         // scale back
672         value = FastMath.scalb(value, SCALING);
673         for (int i = 0; i < 3; ++i) {
674             gradient[i] = FastMath.scalb(gradient[i], SCALING);
675             for (int j = 0; j <= i; ++j) {
676                 hessian[i][j] = FastMath.scalb(hessian[i][j], SCALING);
677             }
678         }
679 
680 
681         // apply the global mu/r factor
682         final double muOr = mu / r;
683         value         *= muOr;
684         gradient[0]    = muOr * gradient[0] - value / r;
685         gradient[1]   *= muOr;
686         gradient[2]   *= muOr;
687         hessian[0][0]  = muOr * hessian[0][0] - 2 * gradient[0] / r;
688         hessian[1][0]  = muOr * hessian[1][0] -     gradient[1] / r;
689         hessian[2][0]  = muOr * hessian[2][0] -     gradient[2] / r;
690         hessian[1][1] *= muOr;
691         hessian[2][1] *= muOr;
692         hessian[2][2] *= muOr;
693 
694         // convert gradient and Hessian from spherical to Cartesian
695         final SphericalCoordinates sc = new SphericalCoordinates(position);
696         return new GradientHessian(sc.toCartesianGradient(gradient),
697                                    sc.toCartesianHessian(hessian, gradient));
698 
699 
700     }
701 
702     /** Container for gradient and Hessian. */
703     private static class GradientHessian {
704 
705         /** Gradient. */
706         private final double[] gradient;
707 
708         /** Hessian. */
709         private final double[][] hessian;
710 
711         /** Simple constructor.
712          * <p>
713          * A reference to the arrays is stored, they are <strong>not</strong> cloned.
714          * </p>
715          * @param gradient gradient
716          * @param hessian hessian
717          */
718         GradientHessian(final double[] gradient, final double[][] hessian) {
719             this.gradient = gradient;
720             this.hessian  = hessian;
721         }
722 
723         /** Get a reference to the gradient.
724          * @return gradient (a reference to the internal array is returned)
725          */
726         public double[] getGradient() {
727             return gradient;
728         }
729 
730         /** Get a reference to the Hessian.
731          * @return Hessian (a reference to the internal array is returned)
732          */
733         public double[][] getHessian() {
734             return hessian;
735         }
736 
737     }
738 
739     /** Compute a/r powers array.
740      * @param aOr a/r
741      * @return array containing (a/r)<sup>n</sup>
742      */
743     private double[] createDistancePowersArray(final double aOr) {
744 
745         // initialize array
746         final double[] aOrN = new double[provider.getMaxDegree() + 1];
747         aOrN[0] = 1;
748         if (provider.getMaxDegree() > 0) {
749             aOrN[1] = aOr;
750         }
751 
752         // fill up array
753         for (int n = 2; n < aOrN.length; ++n) {
754             final int p = n / 2;
755             final int q = n - p;
756             aOrN[n] = aOrN[p] * aOrN[q];
757         }
758 
759         return aOrN;
760 
761     }
762     /** Compute a/r powers array.
763      * @param aOr a/r
764      * @param <T> type of field used
765      * @return array containing (a/r)<sup>n</sup>
766      */
767     private <T extends CalculusFieldElement<T>> T[] createDistancePowersArray(final T aOr) {
768 
769         // initialize array
770         final T[] aOrN = MathArrays.buildArray(aOr.getField(), provider.getMaxDegree() + 1);
771         aOrN[0] = aOr.getField().getOne();
772         if (provider.getMaxDegree() > 0) {
773             aOrN[1] = aOr;
774         }
775 
776         // fill up array
777         for (int n = 2; n < aOrN.length; ++n) {
778             final int p = n / 2;
779             final int q = n - p;
780             aOrN[n] = aOrN[p].multiply(aOrN[q]);
781         }
782 
783         return aOrN;
784 
785     }
786 
787     /** Compute longitude cosines and sines.
788      * @param cosLambda cos(λ)
789      * @param sinLambda sin(λ)
790      * @return array containing cos(m &times; λ) in row 0
791      * and sin(m &times; λ) in row 1
792      */
793     private double[][] createCosSinArrays(final double cosLambda, final double sinLambda) {
794 
795         // initialize arrays
796         final double[][] cosSin = new double[2][provider.getMaxOrder() + 1];
797         cosSin[0][0] = 1;
798         cosSin[1][0] = 0;
799         if (provider.getMaxOrder() > 0) {
800             cosSin[0][1] = cosLambda;
801             cosSin[1][1] = sinLambda;
802 
803             // fill up array
804             for (int m = 2; m < cosSin[0].length; ++m) {
805 
806                 // m * lambda is split as p * lambda + q * lambda, trying to avoid
807                 // p or q being much larger than the other. This reduces the number of
808                 // intermediate results reused to compute each value, and hence should limit
809                 // as much as possible roundoff error accumulation
810                 // (this does not change the number of floating point operations)
811                 final int p = m / 2;
812                 final int q = m - p;
813 
814                 cosSin[0][m] = cosSin[0][p] * cosSin[0][q] - cosSin[1][p] * cosSin[1][q];
815                 cosSin[1][m] = cosSin[1][p] * cosSin[0][q] + cosSin[0][p] * cosSin[1][q];
816             }
817         }
818 
819         return cosSin;
820 
821     }
822 
823     /** Compute longitude cosines and sines.
824      * @param cosLambda cos(λ)
825      * @param sinLambda sin(λ)
826      * @param <T> type of field used
827      * @return array containing cos(m &times; λ) in row 0
828      * and sin(m &times; λ) in row 1
829      */
830     private <T extends CalculusFieldElement<T>> T[][] createCosSinArrays(final T cosLambda, final T sinLambda) {
831 
832         final T one = cosLambda.getField().getOne();
833         final T zero = cosLambda.getField().getZero();
834         // initialize arrays
835         final T[][] cosSin = MathArrays.buildArray(one.getField(), 2, provider.getMaxOrder() + 1);
836         cosSin[0][0] = one;
837         cosSin[1][0] = zero;
838         if (provider.getMaxOrder() > 0) {
839             cosSin[0][1] = cosLambda;
840             cosSin[1][1] = sinLambda;
841 
842             // fill up array
843             for (int m = 2; m < cosSin[0].length; ++m) {
844 
845                 // m * lambda is split as p * lambda + q * lambda, trying to avoid
846                 // p or q being much larger than the other. This reduces the number of
847                 // intermediate results reused to compute each value, and hence should limit
848                 // as much as possible roundoff error accumulation
849                 // (this does not change the number of floating point operations)
850                 final int p = m / 2;
851                 final int q = m - p;
852 
853                 cosSin[0][m] = cosSin[0][p].multiply(cosSin[0][q]).subtract(cosSin[1][p].multiply(cosSin[1][q]));
854                 cosSin[1][m] = cosSin[1][p].multiply(cosSin[0][q]).add(cosSin[0][p].multiply(cosSin[1][q]));
855 
856             }
857         }
858 
859         return cosSin;
860 
861     }
862 
863     /** Compute one order of tesseral terms.
864      * <p>
865      * This corresponds to equations 27 and 30 of the paper.
866      * </p>
867      * @param m current order
868      * @param degree max degree
869      * @param index index in the flattened array
870      * @param t cos(θ), where θ is the polar angle
871      * @param u sin(θ), where θ is the polar angle
872      * @param tOu t/u
873      * @param pnm0Plus2 array containing scaled P<sub>n,m+2</sub>/u<sup>m+2</sup>
874      * @param pnm0Plus1 array containing scaled P<sub>n,m+1</sub>/u<sup>m+1</sup>
875      * @param pnm1Plus1 array containing scaled dP<sub>n,m+1</sub>/u<sup>m+1</sup>
876      * (may be null if second derivatives are not needed)
877      * @param pnm0 array to fill with scaled P<sub>n,m</sub>/u<sup>m</sup>
878      * @param pnm1 array to fill with scaled dP<sub>n,m</sub>/u<sup>m</sup>
879      * (may be null if first derivatives are not needed)
880      * @param pnm2 array to fill with scaled d²P<sub>n,m</sub>/u<sup>m</sup>
881      * (may be null if second derivatives are not needed)
882      * @return new value for index
883      */
884     private int computeTesseral(final int m, final int degree, final int index,
885                                 final double t, final double u, final double tOu,
886                                 final double[] pnm0Plus2, final double[] pnm0Plus1, final double[] pnm1Plus1,
887                                 final double[] pnm0, final double[] pnm1, final double[] pnm2) {
888 
889         final double u2 = u * u;
890 
891         // initialize recursion from sectorial terms
892         int n = FastMath.max(2, m);
893         if (n == m) {
894             pnm0[n] = sectorial[n];
895             ++n;
896         }
897 
898         // compute tesseral values
899         int localIndex = index;
900         while (n <= degree) {
901 
902             // value (equation 27 of the paper)
903             pnm0[n] = gnmOj[localIndex] * t * pnm0Plus1[n] - hnmOj[localIndex] * u2 * pnm0Plus2[n];
904 
905             ++localIndex;
906             ++n;
907 
908         }
909 
910         if (pnm1 != null) {
911 
912             // initialize recursion from sectorial terms
913             n = FastMath.max(2, m);
914             if (n == m) {
915                 pnm1[n] = m * tOu * pnm0[n];
916                 ++n;
917             }
918 
919             // compute tesseral values and derivatives with respect to polar angle
920             localIndex = index;
921             while (n <= degree) {
922 
923                 // first derivative (equation 30 of the paper)
924                 pnm1[n] = m * tOu * pnm0[n] - enm[localIndex] * u * pnm0Plus1[n];
925 
926                 ++localIndex;
927                 ++n;
928 
929             }
930 
931             if (pnm2 != null) {
932 
933                 // initialize recursion from sectorial terms
934                 n = FastMath.max(2, m);
935                 if (n == m) {
936                     pnm2[n] = m * (tOu * pnm1[n] - pnm0[n] / u2);
937                     ++n;
938                 }
939 
940                 // compute tesseral values and derivatives with respect to polar angle
941                 localIndex = index;
942                 while (n <= degree) {
943 
944                     // second derivative (differential of equation 30 with respect to theta)
945                     pnm2[n] = m * (tOu * pnm1[n] - pnm0[n] / u2) - enm[localIndex] * u * pnm1Plus1[n];
946 
947                     ++localIndex;
948                     ++n;
949 
950                 }
951 
952             }
953 
954         }
955 
956         return localIndex;
957 
958     }
959 
960     /** Compute one order of tesseral terms.
961      * <p>
962      * This corresponds to equations 27 and 30 of the paper.
963      * </p>
964      * @param m current order
965      * @param degree max degree
966      * @param index index in the flattened array
967      * @param t cos(θ), where θ is the polar angle
968      * @param u sin(θ), where θ is the polar angle
969      * @param tOu t/u
970      * @param pnm0Plus2 array containing scaled P<sub>n,m+2</sub>/u<sup>m+2</sup>
971      * @param pnm0Plus1 array containing scaled P<sub>n,m+1</sub>/u<sup>m+1</sup>
972      * @param pnm1Plus1 array containing scaled dP<sub>n,m+1</sub>/u<sup>m+1</sup>
973      * (may be null if second derivatives are not needed)
974      * @param pnm0 array to fill with scaled P<sub>n,m</sub>/u<sup>m</sup>
975      * @param pnm1 array to fill with scaled dP<sub>n,m</sub>/u<sup>m</sup>
976      * (may be null if first derivatives are not needed)
977      * @param pnm2 array to fill with scaled d²P<sub>n,m</sub>/u<sup>m</sup>
978      * (may be null if second derivatives are not needed)
979      * @param <T> instance of field element
980      * @return new value for index
981      */
982     private <T extends CalculusFieldElement<T>> int computeTesseral(final int m, final int degree, final int index,
983                                                                 final T t, final T u, final T tOu,
984                                                                 final T[] pnm0Plus2, final T[] pnm0Plus1, final T[] pnm1Plus1,
985                                                                 final T[] pnm0, final T[] pnm1, final T[] pnm2) {
986 
987         final T u2 = u.square();
988         final T zero = u.getField().getZero();
989         // initialize recursion from sectorial terms
990         int n = FastMath.max(2, m);
991         if (n == m) {
992             pnm0[n] = zero.newInstance(sectorial[n]);
993             ++n;
994         }
995 
996         // compute tesseral values
997         int localIndex = index;
998         while (n <= degree) {
999 
1000             // value (equation 27 of the paper)
1001             pnm0[n] = t.multiply(gnmOj[localIndex]).multiply(pnm0Plus1[n]).subtract(u2.multiply(pnm0Plus2[n]).multiply(hnmOj[localIndex]));
1002             ++localIndex;
1003             ++n;
1004 
1005         }
1006         if (pnm1 != null) {
1007 
1008             // initialize recursion from sectorial terms
1009             n = FastMath.max(2, m);
1010             if (n == m) {
1011                 pnm1[n] = tOu.multiply(m).multiply(pnm0[n]);
1012                 ++n;
1013             }
1014 
1015             // compute tesseral values and derivatives with respect to polar angle
1016             localIndex = index;
1017             while (n <= degree) {
1018 
1019                 // first derivative (equation 30 of the paper)
1020                 pnm1[n] = tOu.multiply(m).multiply(pnm0[n]).subtract(u.multiply(enm[localIndex]).multiply(pnm0Plus1[n]));
1021 
1022                 ++localIndex;
1023                 ++n;
1024 
1025             }
1026 
1027             if (pnm2 != null) {
1028 
1029                 // initialize recursion from sectorial terms
1030                 n = FastMath.max(2, m);
1031                 if (n == m) {
1032                     pnm2[n] =   tOu.multiply(pnm1[n]).subtract(pnm0[n].divide(u2)).multiply(m);
1033                     ++n;
1034                 }
1035 
1036                 // compute tesseral values and derivatives with respect to polar angle
1037                 localIndex = index;
1038                 while (n <= degree) {
1039 
1040                     // second derivative (differential of equation 30 with respect to theta)
1041                     pnm2[n] = tOu.multiply(pnm1[n]).subtract(pnm0[n].divide(u2)).multiply(m).subtract(u.multiply(pnm1Plus1[n]).multiply(enm[localIndex]));
1042                     ++localIndex;
1043                     ++n;
1044 
1045                 }
1046 
1047             }
1048 
1049         }
1050         return localIndex;
1051 
1052     }
1053 
1054     /** {@inheritDoc} */
1055     @Override
1056     public Vector3D acceleration(final SpacecraftState s, final double[] parameters) {
1057 
1058         final double mu = parameters[0];
1059 
1060         // get the position in body frame
1061         final AbsoluteDate date       = s.getDate();
1062         final StaticTransform fromBodyFrame = bodyFrame.getStaticTransformTo(s.getFrame(), date);
1063         final StaticTransform toBodyFrame   = fromBodyFrame.getInverse();
1064         final Vector3D position       = toBodyFrame.transformPosition(s.getPosition());
1065 
1066         // gradient of the non-central part of the gravity field
1067         return fromBodyFrame.transformVector(new Vector3D(gradient(date, position, mu)));
1068 
1069     }
1070 
1071     /** {@inheritDoc} */
1072     public <T extends CalculusFieldElement<T>> FieldVector3D<T> acceleration(final FieldSpacecraftState<T> s,
1073                                                                          final T[] parameters) {
1074 
1075         final T mu = parameters[0];
1076 
1077         // check for faster computation dedicated to derivatives with respect to state and gravitational parameter
1078         if (s.getDate().hasZeroField()) {
1079             if (isGradientStateDerivative(s) && isGradientConstantOrMuDerivative(mu)) {
1080                 @SuppressWarnings("unchecked")
1081                 final FieldVector3D<Gradient> p = (FieldVector3D<Gradient>) s.getPosition();
1082                 @SuppressWarnings("unchecked")
1083                 final FieldVector3D<T> a = (FieldVector3D<T>) accelerationWrtState(s.getDate().toAbsoluteDate(),
1084                         s.getFrame(), p,
1085                         (Gradient) mu);
1086                 return a;
1087             } else if (isDSStateDerivative(s) && isDSConstantOrMuDerivative(mu)) {
1088                 @SuppressWarnings("unchecked")
1089                 final FieldVector3D<DerivativeStructure> p = (FieldVector3D<DerivativeStructure>) s.getPosition();
1090                 @SuppressWarnings("unchecked")
1091                 final FieldVector3D<T> a = (FieldVector3D<T>) accelerationWrtState(s.getDate().toAbsoluteDate(),
1092                         s.getFrame(), p,
1093                         (DerivativeStructure) mu);
1094                 return a;
1095             }
1096         }
1097 
1098         // get the position in body frame
1099         final FieldAbsoluteDate<T> date             = s.getDate();
1100         final FieldStaticTransform<T> fromBodyFrame = bodyFrame.getStaticTransformTo(s.getFrame(), date);
1101         final FieldStaticTransform<T> toBodyFrame   = fromBodyFrame.getInverse();
1102         final FieldVector3D<T> position             = toBodyFrame.transformPosition(s.getPosition());
1103 
1104         // gradient of the non-central part of the gravity field
1105         return fromBodyFrame.transformVector(new FieldVector3D<>(gradient(date, position, mu)));
1106 
1107     }
1108 
1109     /** Check if a field state corresponds to derivatives with respect to state.
1110      * @param state state to check
1111      * @param <T> type of the filed elements
1112      * @return true if state corresponds to derivatives with respect to state
1113      * @since 9.0
1114      */
1115     private <T extends CalculusFieldElement<T>> boolean isDSStateDerivative(final FieldSpacecraftState<T> state) {
1116         try {
1117             final DerivativeStructure dsMass = (DerivativeStructure) state.getMass();
1118             final int o = dsMass.getOrder();
1119             final int p = dsMass.getFreeParameters();
1120             if (o != 1 || p < 3) {
1121                 return false;
1122             }
1123             @SuppressWarnings("unchecked")
1124             final FieldPVCoordinates<DerivativeStructure> pv = (FieldPVCoordinates<DerivativeStructure>) state.getPVCoordinates();
1125             return isVariable(pv.getPosition().getX(), 0) &&
1126                    isVariable(pv.getPosition().getY(), 1) &&
1127                    isVariable(pv.getPosition().getZ(), 2);
1128         } catch (ClassCastException cce) {
1129             return false;
1130         }
1131     }
1132 
1133     /** Check if a field state corresponds to derivatives with respect to state.
1134      * @param state state to check
1135      * @param <T> type of the filed elements
1136      * @return true if state corresponds to derivatives with respect to state
1137      * @since 10.2
1138      */
1139     private <T extends CalculusFieldElement<T>> boolean isGradientStateDerivative(final FieldSpacecraftState<T> state) {
1140         try {
1141             final Gradient gMass = (Gradient) state.getMass();
1142             final int p = gMass.getFreeParameters();
1143             if (p < 3) {
1144                 return false;
1145             }
1146             @SuppressWarnings("unchecked")
1147             final FieldPVCoordinates<Gradient> pv = (FieldPVCoordinates<Gradient>) state.getPVCoordinates();
1148             return isVariable(pv.getPosition().getX(), 0) &&
1149                    isVariable(pv.getPosition().getY(), 1) &&
1150                    isVariable(pv.getPosition().getZ(), 2);
1151         } catch (ClassCastException cce) {
1152             return false;
1153         }
1154     }
1155 
1156     /**
1157      * Check if a field gravitational parameter is constant or represents an independent variable.
1158      * @param mu field gravitational parameter
1159      * @return true if the field gravitational parameter is constant or represents an independent variable
1160      * @param <T> type of field
1161      * @since 13.1.3
1162      */
1163     private <T extends CalculusFieldElement<T>> boolean isDSConstantOrMuDerivative(final T mu) {
1164         try {
1165             final double[] derivatives = ((DerivativeStructure) mu).getAllDerivatives();
1166             boolean check = true;
1167             for (int i = 1; i < derivatives.length; i++) {
1168                 if (i == 4) {
1169                     check &= derivatives[i] == 0.0 || derivatives[i] == 1.0;
1170                 } else {
1171                     check &= derivatives[i] == 0.0;
1172                 }
1173             }
1174             return check;
1175         } catch (ClassCastException cce) {
1176             return false;
1177         }
1178     }
1179 
1180     /** Check if a derivative represents a specified variable.
1181      * @param ds derivative to check
1182      * @param index index of the variable
1183      * @return true if the derivative represents a specified variable
1184      * @since 9.0
1185      */
1186     private boolean isVariable(final DerivativeStructure ds, final int index) {
1187         final double[] derivatives = ds.getAllDerivatives();
1188         boolean check = true;
1189         for (int i = 1; i < derivatives.length; ++i) {
1190             check &= derivatives[i] == ((index + 1 == i) ? 1.0 : 0.0);
1191         }
1192         return check;
1193     }
1194 
1195     /**
1196      * Check if a field gravitational parameter is constant or represents an independent variable.
1197      * @param mu field gravitational parameter
1198      * @return true if the field gravitational parameter is constant or represents an independent variable
1199      * @param <T> type of field
1200      * @since 13.1.3
1201      */
1202     private <T extends CalculusFieldElement<T>> boolean isGradientConstantOrMuDerivative(final T mu) {
1203         try {
1204             final double[] derivatives = ((Gradient) mu).getGradient();
1205             boolean check = true;
1206             for (int i = 0; i < derivatives.length; i++) {
1207                 if (i == 3) {
1208                     check &= derivatives[i] == 0.0 || derivatives[i] == 1.0;
1209                 } else {
1210                     check &= derivatives[i] == 0.0;
1211                 }
1212             }
1213             return check;
1214         } catch (ClassCastException cce) {
1215             return false;
1216         }
1217     }
1218 
1219     /** Check if a derivative represents a specified variable.
1220      * @param g derivative to check
1221      * @param index index of the variable
1222      * @return true if the derivative represents a specified variable
1223      * @since 10.2
1224      */
1225     private boolean isVariable(final Gradient g, final int index) {
1226         final double[] derivatives = g.getGradient();
1227         boolean check = true;
1228         for (int i = 0; i < derivatives.length; ++i) {
1229             check &= derivatives[i] == ((index == i) ? 1.0 : 0.0);
1230         }
1231         return check;
1232     }
1233 
1234     /** Compute acceleration derivatives with respect to state parameters.
1235      * <p>
1236      * From a theoretical point of view, this method computes the same values
1237      * as {@link #acceleration(FieldSpacecraftState, CalculusFieldElement[])} in the
1238      * specific case of {@link DerivativeStructure} with respect to state, so
1239      * it is less general. However, it is *much* faster in this important case.
1240      * <p>
1241      * <p>
1242      * The derivatives should be computed with respect to position. The input
1243      * parameters already take into account the free parameters (6 or 7 depending
1244      * on derivation with respect to mass being considered or not) and order
1245      * (always 1). Free parameters at indices 0, 1 and 2 correspond to derivatives
1246      * with respect to position. Free parameters at indices 3, 4 and 5 correspond
1247      * to derivatives with respect to velocity (these derivatives will remain zero
1248      * as acceleration due to gravity does not depend on velocity). Free parameter
1249      * at index 6 (if present) corresponds to to derivatives with respect to mass
1250      * (this derivative will remain zero as acceleration due to gravity does not
1251      * depend on mass).
1252      * </p>
1253      * @param date current date
1254      * @param frame inertial reference frame for state (both orbit and attitude)
1255      * @param position position of spacecraft in inertial frame
1256      * @param mu central attraction coefficient to use
1257      * @return acceleration with all derivatives specified by the input parameters
1258      * own derivatives
1259      * @since 6.0
1260      */
1261     private FieldVector3D<DerivativeStructure> accelerationWrtState(final AbsoluteDate date, final Frame frame,
1262                                                                     final FieldVector3D<DerivativeStructure> position,
1263                                                                     final DerivativeStructure mu) {
1264 
1265         // free parameters
1266         final int freeParameters = mu.getFreeParameters();
1267 
1268         // get the position in body frame
1269         final StaticTransform fromBodyFrame = bodyFrame.getStaticTransformTo(frame, date);
1270         final StaticTransform toBodyFrame   = fromBodyFrame.getInverse();
1271         final Vector3D positionBody   = toBodyFrame.transformPosition(position.toVector3D());
1272 
1273         // compute gradient and Hessian
1274         final GradientHessian gh   = gradientHessian(date, positionBody, mu.getReal());
1275 
1276         // gradient of the non-central part of the gravity field
1277         final double[] gInertial = fromBodyFrame.transformVector(new Vector3D(gh.getGradient())).toArray();
1278 
1279         // Hessian of the non-central part of the gravity field
1280         final RealMatrix hBody     = new Array2DRowRealMatrix(gh.getHessian(), false);
1281         final RealMatrix rot       = new Array2DRowRealMatrix(toBodyFrame.getRotation().getMatrix());
1282         final RealMatrix hInertial = rot.transposeMultiply(hBody).multiply(rot);
1283 
1284         // distribute all partial derivatives in a compact acceleration vector
1285         final double[] derivatives = new double[freeParameters + 1];
1286         final DerivativeStructure[] accDer = new DerivativeStructure[3];
1287         for (int i = 0; i < 3; ++i) {
1288 
1289             // first element is value of acceleration (i.e. gradient of field)
1290             derivatives[0] = gInertial[i];
1291 
1292             // Jacobian of acceleration (i.e. Hessian of field)
1293             derivatives[1] = hInertial.getEntry(i, 0);
1294             derivatives[2] = hInertial.getEntry(i, 1);
1295             derivatives[3] = hInertial.getEntry(i, 2);
1296 
1297             // next element is derivative with respect to parameter mu
1298             if (derivatives.length > 4 && isVariable(mu, 3)) {
1299                 derivatives[4] = gInertial[i] / mu.getReal();
1300             }
1301 
1302             accDer[i] = position.getX().getFactory().build(derivatives);
1303 
1304         }
1305 
1306         return new FieldVector3D<>(accDer);
1307 
1308     }
1309 
1310     /** Compute acceleration derivatives with respect to state parameters.
1311      * <p>
1312      * From a theoretical point of view, this method computes the same values
1313      * as {@link #acceleration(FieldSpacecraftState, CalculusFieldElement[])} in the
1314      * specific case of {@link DerivativeStructure} with respect to state, so
1315      * it is less general. However, it is *much* faster in this important case.
1316      * <p>
1317      * <p>
1318      * The derivatives should be computed with respect to position. The input
1319      * parameters already take into account the free parameters (6 or 7 depending
1320      * on derivation with respect to mass being considered or not) and order
1321      * (always 1). Free parameters at indices 0, 1 and 2 correspond to derivatives
1322      * with respect to position. Free parameters at indices 3, 4 and 5 correspond
1323      * to derivatives with respect to velocity (these derivatives will remain zero
1324      * as acceleration due to gravity does not depend on velocity). Free parameter
1325      * at index 6 (if present) corresponds to to derivatives with respect to mass
1326      * (this derivative will remain zero as acceleration due to gravity does not
1327      * depend on mass).
1328      * </p>
1329      * @param date current date
1330      * @param frame inertial reference frame for state (both orbit and attitude)
1331      * @param position position of spacecraft in inertial frame
1332      * @param mu central attraction coefficient to use
1333      * @return acceleration with all derivatives specified by the input parameters
1334      * own derivatives
1335      * @since 10.2
1336      */
1337     private FieldVector3D<Gradient> accelerationWrtState(final AbsoluteDate date, final Frame frame,
1338                                                          final FieldVector3D<Gradient> position,
1339                                                          final Gradient mu) {
1340 
1341         // free parameters
1342         final int freeParameters = mu.getFreeParameters();
1343 
1344         // get the position in body frame
1345         final StaticTransform fromBodyFrame = bodyFrame.getStaticTransformTo(frame, date);
1346         final StaticTransform toBodyFrame   = fromBodyFrame.getInverse();
1347         final Vector3D positionBody   = toBodyFrame.transformPosition(position.toVector3D());
1348 
1349         // compute gradient and Hessian
1350         final GradientHessian gh   = gradientHessian(date, positionBody, mu.getReal());
1351 
1352         // gradient of the non-central part of the gravity field
1353         final double[] gInertial = fromBodyFrame.transformVector(new Vector3D(gh.getGradient())).toArray();
1354 
1355         // Hessian of the non-central part of the gravity field
1356         final RealMatrix hBody     = new Array2DRowRealMatrix(gh.getHessian(), false);
1357         final RealMatrix rot       = new Array2DRowRealMatrix(toBodyFrame.getRotation().getMatrix());
1358         final RealMatrix hInertial = rot.transposeMultiply(hBody).multiply(rot);
1359 
1360         // distribute all partial derivatives in a compact acceleration vector
1361         final double[] derivatives = new double[freeParameters];
1362         final Gradient[] accDer = new Gradient[3];
1363         for (int i = 0; i < 3; ++i) {
1364 
1365             // Jacobian of acceleration (i.e. Hessian of field)
1366             derivatives[0] = hInertial.getEntry(i, 0);
1367             derivatives[1] = hInertial.getEntry(i, 1);
1368             derivatives[2] = hInertial.getEntry(i, 2);
1369 
1370             // next element is derivative with respect to parameter mu
1371             if (derivatives.length > 3 && isVariable(mu, 3)) {
1372                 derivatives[3] = gInertial[i] / mu.getReal();
1373             }
1374 
1375             accDer[i] = new Gradient(gInertial[i], derivatives);
1376 
1377         }
1378 
1379         return new FieldVector3D<>(accDer);
1380 
1381     }
1382 
1383     /** {@inheritDoc} */
1384     public List<ParameterDriver> getParametersDrivers() {
1385         return Collections.singletonList(gmParameterDriver);
1386     }
1387 
1388 }