HolmesFeatherstoneAttractionModel.java

  1. /* Copyright 2002-2020 CS GROUP
  2.  * Licensed to CS GROUP (CS) under one or more
  3.  * contributor license agreements.  See the NOTICE file distributed with
  4.  * this work for additional information regarding copyright ownership.
  5.  * CS licenses this file to You under the Apache License, Version 2.0
  6.  * (the "License"); you may not use this file except in compliance with
  7.  * the License.  You may obtain a copy of the License at
  8.  *
  9.  *   http://www.apache.org/licenses/LICENSE-2.0
  10.  *
  11.  * Unless required by applicable law or agreed to in writing, software
  12.  * distributed under the License is distributed on an "AS IS" BASIS,
  13.  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
  14.  * See the License for the specific language governing permissions and
  15.  * limitations under the License.
  16.  */
  17. package org.orekit.forces.gravity;


  18. import java.util.stream.Stream;

  19. import org.hipparchus.Field;
  20. import org.hipparchus.RealFieldElement;
  21. import org.hipparchus.analysis.differentiation.DerivativeStructure;
  22. import org.hipparchus.analysis.differentiation.Gradient;
  23. import org.hipparchus.geometry.euclidean.threed.FieldVector3D;
  24. import org.hipparchus.geometry.euclidean.threed.SphericalCoordinates;
  25. import org.hipparchus.geometry.euclidean.threed.Vector3D;
  26. import org.hipparchus.linear.Array2DRowRealMatrix;
  27. import org.hipparchus.linear.RealMatrix;
  28. import org.hipparchus.util.FastMath;
  29. import org.hipparchus.util.MathArrays;
  30. import org.orekit.forces.AbstractForceModel;
  31. import org.orekit.forces.gravity.potential.NormalizedSphericalHarmonicsProvider;
  32. import org.orekit.forces.gravity.potential.NormalizedSphericalHarmonicsProvider.NormalizedSphericalHarmonics;
  33. import org.orekit.forces.gravity.potential.TideSystem;
  34. import org.orekit.forces.gravity.potential.TideSystemProvider;
  35. import org.orekit.frames.Frame;
  36. import org.orekit.frames.Transform;
  37. import org.orekit.propagation.FieldSpacecraftState;
  38. import org.orekit.propagation.SpacecraftState;
  39. import org.orekit.propagation.events.EventDetector;
  40. import org.orekit.propagation.events.FieldEventDetector;
  41. import org.orekit.time.AbsoluteDate;
  42. import org.orekit.time.FieldAbsoluteDate;
  43. import org.orekit.utils.FieldPVCoordinates;
  44. import org.orekit.utils.ParameterDriver;

  45. /** This class represents the gravitational field of a celestial body.
  46.  * <p>
  47.  * The algorithm implemented in this class has been designed by S. A. Holmes
  48.  * and W. E. Featherstone from Department of Spatial Sciences, Curtin University
  49.  * of Technology, Perth, Australia. It is described in their 2002 paper: <a
  50.  * 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">
  51.  * A unified approach to he Clenshaw summation and the recursive computation of
  52.  * very high degree and order normalised associated Legendre functions</a>
  53.  * (Journal of Geodesy (2002) 76: 279–299).
  54.  * </p>
  55.  * <p>
  56.  * This model directly uses normalized coefficients and stable recursion algorithms
  57.  * so it is more suited to high degree gravity fields than the classical Cunningham
  58.  * Droziner models which use un-normalized coefficients.
  59.  * </p>
  60.  * <p>
  61.  * Among the different algorithms presented in Holmes and Featherstone paper, this
  62.  * class implements the <em>modified forward row method</em>. All recursion coefficients
  63.  * are precomputed and stored for greater performance. This caching was suggested in the
  64.  * paper but not used due to the large memory requirements. Since 2002, even low end
  65.  * computers and mobile devices do have sufficient memory so this caching has become
  66.  * feasible nowadays.
  67.  * <p>
  68.  * @author Luc Maisonobe
  69.  * @since 6.0
  70.  */

  71. public class HolmesFeatherstoneAttractionModel extends AbstractForceModel implements TideSystemProvider {

  72.     /** Exponent scaling to avoid floating point overflow.
  73.      * <p>The paper uses 10^280, we prefer a power of two to preserve accuracy thanks to
  74.      * {@link FastMath#scalb(double, int)}, so we use 2^930 which has the same order of magnitude.
  75.      */
  76.     private static final int SCALING = 930;

  77.     /** Central attraction scaling factor.
  78.      * <p>
  79.      * We use a power of 2 to avoid numeric noise introduction
  80.      * in the multiplications/divisions sequences.
  81.      * </p>
  82.      */
  83.     private static final double MU_SCALE = FastMath.scalb(1.0, 32);

  84.     /** Driver for gravitational parameter. */
  85.     private final ParameterDriver gmParameterDriver;

  86.     /** Provider for the spherical harmonics. */
  87.     private final NormalizedSphericalHarmonicsProvider provider;

  88.     /** Rotating body. */
  89.     private final Frame bodyFrame;

  90.     /** Recursion coefficients g<sub>n,m</sub>/√j. */
  91.     private final double[] gnmOj;

  92.     /** Recursion coefficients h<sub>n,m</sub>/√j. */
  93.     private final double[] hnmOj;

  94.     /** Recursion coefficients e<sub>n,m</sub>. */
  95.     private final double[] enm;

  96.     /** Scaled sectorial Pbar<sub>m,m</sub>/u<sup>m</sup> &times; 2<sup>-SCALING</sup>. */
  97.     private final double[] sectorial;

  98.     /** Creates a new instance.
  99.      * @param centralBodyFrame rotating body frame
  100.      * @param provider provider for spherical harmonics
  101.      * @since 6.0
  102.      */
  103.     public HolmesFeatherstoneAttractionModel(final Frame centralBodyFrame,
  104.                                              final NormalizedSphericalHarmonicsProvider provider) {

  105.         gmParameterDriver = new ParameterDriver(NewtonianAttraction.CENTRAL_ATTRACTION_COEFFICIENT,
  106.                                                 provider.getMu(), MU_SCALE, 0.0, Double.POSITIVE_INFINITY);

  107.         this.provider  = provider;
  108.         this.bodyFrame = centralBodyFrame;

  109.         // the pre-computed arrays hold coefficients from triangular arrays in a single
  110.         // storing neither diagonal elements (n = m) nor the non-diagonal element n=1, m=0
  111.         final int degree = provider.getMaxDegree();
  112.         final int size = FastMath.max(0, degree * (degree + 1) / 2 - 1);
  113.         gnmOj = new double[size];
  114.         hnmOj = new double[size];
  115.         enm   = new double[size];

  116.         // pre-compute the recursion coefficients corresponding to equations 19 and 22
  117.         // from Holmes and Featherstone paper
  118.         // for cache efficiency, elements are stored in the same order they will be used
  119.         // later on, i.e. from rightmost column to leftmost column
  120.         int index = 0;
  121.         for (int m = degree; m >= 0; --m) {
  122.             final int j = (m == 0) ? 2 : 1;
  123.             for (int n = FastMath.max(2, m + 1); n <= degree; ++n) {
  124.                 final double f = (n - m) * (n + m + 1);
  125.                 gnmOj[index] = 2 * (m + 1) / FastMath.sqrt(j * f);
  126.                 hnmOj[index] = FastMath.sqrt((n + m + 2) * (n - m - 1) / (j * f));
  127.                 enm[index]   = FastMath.sqrt(f / j);
  128.                 ++index;
  129.             }
  130.         }

  131.         // scaled sectorial terms corresponding to equation 28 in Holmes and Featherstone paper
  132.         sectorial    = new double[degree + 1];
  133.         sectorial[0] = FastMath.scalb(1.0, -SCALING);
  134.         sectorial[1] = FastMath.sqrt(3) * sectorial[0];
  135.         for (int m = 2; m < sectorial.length; ++m) {
  136.             sectorial[m] = FastMath.sqrt((2 * m + 1) / (2.0 * m)) * sectorial[m - 1];
  137.         }

  138.     }

  139.     /** {@inheritDoc} */
  140.     @Override
  141.     public boolean dependsOnPositionOnly() {
  142.         return true;
  143.     }

  144.     /** {@inheritDoc} */
  145.     public TideSystem getTideSystem() {
  146.         return provider.getTideSystem();
  147.     }

  148.     /** Get the central attraction coefficient μ.
  149.      * @return mu central attraction coefficient (m³/s²)
  150.      */
  151.     public double getMu() {
  152.         return gmParameterDriver.getValue();
  153.     }

  154.     /** Compute the value of the gravity field.
  155.      * @param date current date
  156.      * @param position position at which gravity field is desired in body frame
  157.      * @param mu central attraction coefficient to use
  158.      * @return value of the gravity field (central and non-central parts summed together)
  159.      */
  160.     public double value(final AbsoluteDate date, final Vector3D position,
  161.                         final double mu) {
  162.         return mu / position.getNorm() + nonCentralPart(date, position, mu);
  163.     }

  164.     /** Compute the non-central part of the gravity field.
  165.      * @param date current date
  166.      * @param position position at which gravity field is desired in body frame
  167.      * @param mu central attraction coefficient to use
  168.      * @return value of the non-central part of the gravity field
  169.      */
  170.     public double nonCentralPart(final AbsoluteDate date, final Vector3D position, final double mu) {

  171.         final int degree = provider.getMaxDegree();
  172.         final int order  = provider.getMaxOrder();
  173.         final NormalizedSphericalHarmonics harmonics = provider.onDate(date);

  174.         // allocate the columns for recursion
  175.         double[] pnm0Plus2 = new double[degree + 1];
  176.         double[] pnm0Plus1 = new double[degree + 1];
  177.         double[] pnm0      = new double[degree + 1];

  178.         // compute polar coordinates
  179.         final double x   = position.getX();
  180.         final double y   = position.getY();
  181.         final double z   = position.getZ();
  182.         final double x2  = x * x;
  183.         final double y2  = y * y;
  184.         final double z2  = z * z;
  185.         final double r2  = x2 + y2 + z2;
  186.         final double r   = FastMath.sqrt (r2);
  187.         final double rho = FastMath.sqrt(x2 + y2);
  188.         final double t   = z / r;   // cos(theta), where theta is the polar angle
  189.         final double u   = rho / r; // sin(theta), where theta is the polar angle
  190.         final double tOu = z / rho;

  191.         // compute distance powers
  192.         final double[] aOrN = createDistancePowersArray(provider.getAe() / r);

  193.         // compute longitude cosines/sines
  194.         final double[][] cosSinLambda = createCosSinArrays(position.getX() / rho, position.getY() / rho);

  195.         // outer summation over order
  196.         int    index = 0;
  197.         double value = 0;
  198.         for (int m = degree; m >= 0; --m) {

  199.             // compute tesseral terms without derivatives
  200.             index = computeTesseral(m, degree, index, t, u, tOu,
  201.                                     pnm0Plus2, pnm0Plus1, null, pnm0, null, null);

  202.             if (m <= order) {
  203.                 // compute contribution of current order to field (equation 5 of the paper)

  204.                 // inner summation over degree, for fixed order
  205.                 double sumDegreeS        = 0;
  206.                 double sumDegreeC        = 0;
  207.                 for (int n = FastMath.max(2, m); n <= degree; ++n) {
  208.                     sumDegreeS += pnm0[n] * aOrN[n] * harmonics.getNormalizedSnm(n, m);
  209.                     sumDegreeC += pnm0[n] * aOrN[n] * harmonics.getNormalizedCnm(n, m);
  210.                 }

  211.                 // contribution to outer summation over order
  212.                 value = value * u + cosSinLambda[1][m] * sumDegreeS + cosSinLambda[0][m] * sumDegreeC;

  213.             }

  214.             // rotate the recursion arrays
  215.             final double[] tmp = pnm0Plus2;
  216.             pnm0Plus2 = pnm0Plus1;
  217.             pnm0Plus1 = pnm0;
  218.             pnm0      = tmp;

  219.         }

  220.         // scale back
  221.         value = FastMath.scalb(value, SCALING);

  222.         // apply the global mu/r factor
  223.         return mu * value / r;

  224.     }

  225.     /** Compute the gradient of the non-central part of the gravity field.
  226.      * @param date current date
  227.      * @param position position at which gravity field is desired in body frame
  228.      * @param mu central attraction coefficient to use
  229.      * @return gradient of the non-central part of the gravity field
  230.      */
  231.     public double[] gradient(final AbsoluteDate date, final Vector3D position, final double mu) {

  232.         final int degree = provider.getMaxDegree();
  233.         final int order  = provider.getMaxOrder();
  234.         final NormalizedSphericalHarmonics harmonics = provider.onDate(date);

  235.         // allocate the columns for recursion
  236.         double[] pnm0Plus2  = new double[degree + 1];
  237.         double[] pnm0Plus1  = new double[degree + 1];
  238.         double[] pnm0       = new double[degree + 1];
  239.         final double[] pnm1 = new double[degree + 1];

  240.         // compute polar coordinates
  241.         final double x    = position.getX();
  242.         final double y    = position.getY();
  243.         final double z    = position.getZ();
  244.         final double x2   = x * x;
  245.         final double y2   = y * y;
  246.         final double z2   = z * z;
  247.         final double r2   = x2 + y2 + z2;
  248.         final double r    = FastMath.sqrt (r2);
  249.         final double rho2 = x2 + y2;
  250.         final double rho  = FastMath.sqrt(rho2);
  251.         final double t    = z / r;   // cos(theta), where theta is the polar angle
  252.         final double u    = rho / r; // sin(theta), where theta is the polar angle
  253.         final double tOu  = z / rho;

  254.         // compute distance powers
  255.         final double[] aOrN = createDistancePowersArray(provider.getAe() / r);

  256.         // compute longitude cosines/sines
  257.         final double[][] cosSinLambda = createCosSinArrays(position.getX() / rho, position.getY() / rho);

  258.         // outer summation over order
  259.         int    index = 0;
  260.         double value = 0;
  261.         final double[] gradient = new double[3];
  262.         for (int m = degree; m >= 0; --m) {

  263.             // compute tesseral terms with derivatives
  264.             index = computeTesseral(m, degree, index, t, u, tOu,
  265.                                     pnm0Plus2, pnm0Plus1, null, pnm0, pnm1, null);

  266.             if (m <= order) {
  267.                 // compute contribution of current order to field (equation 5 of the paper)

  268.                 // inner summation over degree, for fixed order
  269.                 double sumDegreeS        = 0;
  270.                 double sumDegreeC        = 0;
  271.                 double dSumDegreeSdR     = 0;
  272.                 double dSumDegreeCdR     = 0;
  273.                 double dSumDegreeSdTheta = 0;
  274.                 double dSumDegreeCdTheta = 0;
  275.                 for (int n = FastMath.max(2, m); n <= degree; ++n) {
  276.                     final double qSnm  = aOrN[n] * harmonics.getNormalizedSnm(n, m);
  277.                     final double qCnm  = aOrN[n] * harmonics.getNormalizedCnm(n, m);
  278.                     final double nOr   = n / r;
  279.                     final double s0    = pnm0[n] * qSnm;
  280.                     final double c0    = pnm0[n] * qCnm;
  281.                     final double s1    = pnm1[n] * qSnm;
  282.                     final double c1    = pnm1[n] * qCnm;
  283.                     sumDegreeS        += s0;
  284.                     sumDegreeC        += c0;
  285.                     dSumDegreeSdR     -= nOr * s0;
  286.                     dSumDegreeCdR     -= nOr * c0;
  287.                     dSumDegreeSdTheta += s1;
  288.                     dSumDegreeCdTheta += c1;
  289.                 }

  290.                 // contribution to outer summation over order
  291.                 // beware that we need to order gradient using the mathematical conventions
  292.                 // compliant with the SphericalCoordinates class, so our lambda is its theta
  293.                 // (and hence at index 1) and our theta is its phi (and hence at index 2)
  294.                 final double sML = cosSinLambda[1][m];
  295.                 final double cML = cosSinLambda[0][m];
  296.                 value            = value       * u + sML * sumDegreeS        + cML * sumDegreeC;
  297.                 gradient[0]      = gradient[0] * u + sML * dSumDegreeSdR     + cML * dSumDegreeCdR;
  298.                 gradient[1]      = gradient[1] * u + m * (cML * sumDegreeS - sML * sumDegreeC);
  299.                 gradient[2]      = gradient[2] * u + sML * dSumDegreeSdTheta + cML * dSumDegreeCdTheta;

  300.             }

  301.             // rotate the recursion arrays
  302.             final double[] tmp = pnm0Plus2;
  303.             pnm0Plus2 = pnm0Plus1;
  304.             pnm0Plus1 = pnm0;
  305.             pnm0      = tmp;

  306.         }

  307.         // scale back
  308.         value       = FastMath.scalb(value,       SCALING);
  309.         gradient[0] = FastMath.scalb(gradient[0], SCALING);
  310.         gradient[1] = FastMath.scalb(gradient[1], SCALING);
  311.         gradient[2] = FastMath.scalb(gradient[2], SCALING);

  312.         // apply the global mu/r factor
  313.         final double muOr = mu / r;
  314.         value            *= muOr;
  315.         gradient[0]       = muOr * gradient[0] - value / r;
  316.         gradient[1]      *= muOr;
  317.         gradient[2]      *= muOr;

  318.         // convert gradient from spherical to Cartesian
  319.         return new SphericalCoordinates(position).toCartesianGradient(gradient);

  320.     }

  321.     /** Compute the gradient of the non-central part of the gravity field.
  322.      * @param date current date
  323.      * @param position position at which gravity field is desired in body frame
  324.      * @param mu central attraction coefficient to use
  325.      * @param <T> type of field used
  326.      * @return gradient of the non-central part of the gravity field
  327.      */
  328.     public <T extends RealFieldElement<T>> T[] gradient(final FieldAbsoluteDate<T> date, final FieldVector3D<T> position,
  329.                                                         final T mu) {

  330.         final int degree = provider.getMaxDegree();
  331.         final int order  = provider.getMaxOrder();
  332.         final NormalizedSphericalHarmonics harmonics = provider.onDate(date.toAbsoluteDate());
  333.         final T zero = date.getField().getZero();
  334.         // allocate the columns for recursion
  335.         T[] pnm0Plus2  = MathArrays.buildArray(date.getField(), degree + 1);
  336.         T[] pnm0Plus1  = MathArrays.buildArray(date.getField(), degree + 1);
  337.         T[] pnm0       = MathArrays.buildArray(date.getField(), degree + 1);
  338.         final T[] pnm1 = MathArrays.buildArray(date.getField(), degree + 1);

  339.         // compute polar coordinates
  340.         final T x    = position.getX();
  341.         final T y    = position.getY();
  342.         final T z    = position.getZ();
  343.         final T x2   = x.multiply(x);
  344.         final T y2   = y.multiply(y);
  345.         final T z2   = z.multiply(z);
  346.         final T r2   = x2.add(y2).add(z2);
  347.         final T r    = r2.sqrt();
  348.         final T rho2 = x2.add(y2);
  349.         final T rho  = rho2.sqrt();
  350.         final T t    = z.divide(r);   // cos(theta), where theta is the polar angle
  351.         final T u    = rho.divide(r); // sin(theta), where theta is the polar angle
  352.         final T tOu  = z.divide(rho);

  353.         // compute distance powers
  354.         final T[] aOrN = createDistancePowersArray(r.reciprocal().multiply(provider.getAe()));

  355.         // compute longitude cosines/sines
  356.         final T[][] cosSinLambda = createCosSinArrays(rho.reciprocal().multiply(position.getX()), rho.reciprocal().multiply(position.getY()));
  357.         // outer summation over order
  358.         int    index = 0;
  359.         T value = zero;
  360.         final T[] gradient = MathArrays.buildArray(zero.getField(), 3);
  361.         for (int m = degree; m >= 0; --m) {

  362.             // compute tesseral terms with derivatives
  363.             index = computeTesseral(m, degree, index, t, u, tOu,
  364.                                     pnm0Plus2, pnm0Plus1, null, pnm0, pnm1, null);
  365.             if (m <= order) {
  366.                 // compute contribution of current order to field (equation 5 of the paper)

  367.                 // inner summation over degree, for fixed order
  368.                 T sumDegreeS        = zero;
  369.                 T sumDegreeC        = zero;
  370.                 T dSumDegreeSdR     = zero;
  371.                 T dSumDegreeCdR     = zero;
  372.                 T dSumDegreeSdTheta = zero;
  373.                 T dSumDegreeCdTheta = zero;
  374.                 for (int n = FastMath.max(2, m); n <= degree; ++n) {
  375.                     final T qSnm  = aOrN[n].multiply(harmonics.getNormalizedSnm(n, m));
  376.                     final T qCnm  = aOrN[n].multiply(harmonics.getNormalizedCnm(n, m));
  377.                     final T nOr   = r.reciprocal().multiply(n);
  378.                     final T s0    = pnm0[n].multiply(qSnm);
  379.                     final T c0    = pnm0[n].multiply(qCnm);
  380.                     final T s1    = pnm1[n].multiply(qSnm);
  381.                     final T c1    = pnm1[n].multiply(qCnm);
  382.                     sumDegreeS        = sumDegreeS       .add(s0);
  383.                     sumDegreeC        = sumDegreeC       .add(c0);
  384.                     dSumDegreeSdR     = dSumDegreeSdR    .subtract(nOr.multiply(s0));
  385.                     dSumDegreeCdR     = dSumDegreeCdR    .subtract(nOr.multiply(c0));
  386.                     dSumDegreeSdTheta = dSumDegreeSdTheta.add(s1);
  387.                     dSumDegreeCdTheta = dSumDegreeCdTheta.add(c1);
  388.                 }

  389.                 // contribution to outer summation over order
  390.                 // beware that we need to order gradient using the mathematical conventions
  391.                 // compliant with the SphericalCoordinates class, so our lambda is its theta
  392.                 // (and hence at index 1) and our theta is its phi (and hence at index 2)
  393.                 final T sML = cosSinLambda[1][m];
  394.                 final T cML = cosSinLambda[0][m];
  395.                 value            = value      .multiply(u).add(sML.multiply(sumDegreeS   )).add(cML.multiply(sumDegreeC));
  396.                 gradient[0]      = gradient[0].multiply(u).add(sML.multiply(dSumDegreeSdR)).add(cML.multiply(dSumDegreeCdR));
  397.                 gradient[1]      = gradient[1].multiply(u).add(cML.multiply(sumDegreeS).subtract(sML.multiply(sumDegreeC)).multiply(m));
  398.                 gradient[2]      = gradient[2].multiply(u).add(sML.multiply(dSumDegreeSdTheta)).add(cML.multiply(dSumDegreeCdTheta));
  399.             }
  400.             // rotate the recursion arrays
  401.             final T[] tmp = pnm0Plus2;
  402.             pnm0Plus2 = pnm0Plus1;
  403.             pnm0Plus1 = pnm0;
  404.             pnm0      = tmp;

  405.         }
  406.         // scale back
  407.         value       = value.scalb(SCALING);
  408.         gradient[0] = gradient[0].scalb(SCALING);
  409.         gradient[1] = gradient[1].scalb(SCALING);
  410.         gradient[2] = gradient[2].scalb(SCALING);

  411.         // apply the global mu/r factor
  412.         final T muOr = r.reciprocal().multiply(mu);
  413.         value            = value.multiply(muOr);
  414.         gradient[0]       = muOr.multiply(gradient[0]).subtract(value.divide(r));
  415.         gradient[1]      = gradient[1].multiply(muOr);
  416.         gradient[2]      = gradient[2].multiply(muOr);

  417.         // convert gradient from spherical to Cartesian
  418.         // Cartesian coordinates
  419.         // remaining spherical coordinates
  420.         final T rPos     = position.getNorm();
  421.         // intermediate variables
  422.         final T xPos    = position.getX();
  423.         final T yPos    = position.getY();
  424.         final T zPos    = position.getZ();
  425.         final T rho2Pos = x.multiply(x).add(y.multiply(y));
  426.         final T rhoPos  = rho2.sqrt();
  427.         final T r2Pos   = rho2.add(z.multiply(z));

  428.         final T[][] jacobianPos = MathArrays.buildArray(zero.getField(), 3, 3);

  429.         // row representing the gradient of r
  430.         jacobianPos[0][0] = xPos.divide(rPos);
  431.         jacobianPos[0][1] = yPos.divide(rPos);
  432.         jacobianPos[0][2] = zPos.divide(rPos);

  433.         // row representing the gradient of theta
  434.         jacobianPos[1][0] =  yPos.negate().divide(rho2Pos);
  435.         jacobianPos[1][1] =  xPos.divide(rho2Pos);
  436.         // jacobian[1][2] is already set to 0 at allocation time

  437.         // row representing the gradient of phi
  438.         jacobianPos[2][0] = xPos.multiply(zPos).divide(rhoPos.multiply(r2Pos));
  439.         jacobianPos[2][1] = yPos.multiply(zPos).divide(rhoPos.multiply(r2Pos));
  440.         jacobianPos[2][2] = rhoPos.negate().divide(r2Pos);
  441.         final T[] cartGradPos = MathArrays.buildArray(zero.getField(), 3);
  442.         cartGradPos[0] = gradient[0].multiply(jacobianPos[0][0]).add(gradient[1].multiply(jacobianPos[1][0])).add(gradient[2].multiply(jacobianPos[2][0]));
  443.         cartGradPos[1] = gradient[0].multiply(jacobianPos[0][1]).add(gradient[1].multiply(jacobianPos[1][1])).add(gradient[2].multiply(jacobianPos[2][1]));
  444.         cartGradPos[2] = gradient[0].multiply(jacobianPos[0][2])                                      .add(gradient[2].multiply(jacobianPos[2][2]));
  445.         return cartGradPos;

  446.     }

  447.     /** Compute both the gradient and the hessian of the non-central part of the gravity field.
  448.      * @param date current date
  449.      * @param position position at which gravity field is desired in body frame
  450.      * @param mu central attraction coefficient to use
  451.      * @return gradient and hessian of the non-central part of the gravity field
  452.      */
  453.     private GradientHessian gradientHessian(final AbsoluteDate date, final Vector3D position, final double mu) {

  454.         final int degree = provider.getMaxDegree();
  455.         final int order  = provider.getMaxOrder();
  456.         final NormalizedSphericalHarmonics harmonics = provider.onDate(date);

  457.         // allocate the columns for recursion
  458.         double[] pnm0Plus2  = new double[degree + 1];
  459.         double[] pnm0Plus1  = new double[degree + 1];
  460.         double[] pnm0       = new double[degree + 1];
  461.         double[] pnm1Plus1  = new double[degree + 1];
  462.         double[] pnm1       = new double[degree + 1];
  463.         final double[] pnm2 = new double[degree + 1];

  464.         // compute polar coordinates
  465.         final double x    = position.getX();
  466.         final double y    = position.getY();
  467.         final double z    = position.getZ();
  468.         final double x2   = x * x;
  469.         final double y2   = y * y;
  470.         final double z2   = z * z;
  471.         final double r2   = x2 + y2 + z2;
  472.         final double r    = FastMath.sqrt (r2);
  473.         final double rho2 = x2 + y2;
  474.         final double rho  = FastMath.sqrt(rho2);
  475.         final double t    = z / r;   // cos(theta), where theta is the polar angle
  476.         final double u    = rho / r; // sin(theta), where theta is the polar angle
  477.         final double tOu  = z / rho;

  478.         // compute distance powers
  479.         final double[] aOrN = createDistancePowersArray(provider.getAe() / r);

  480.         // compute longitude cosines/sines
  481.         final double[][] cosSinLambda = createCosSinArrays(position.getX() / rho, position.getY() / rho);

  482.         // outer summation over order
  483.         int    index = 0;
  484.         double value = 0;
  485.         final double[]   gradient = new double[3];
  486.         final double[][] hessian  = new double[3][3];
  487.         for (int m = degree; m >= 0; --m) {

  488.             // compute tesseral terms
  489.             index = computeTesseral(m, degree, index, t, u, tOu,
  490.                                     pnm0Plus2, pnm0Plus1, pnm1Plus1, pnm0, pnm1, pnm2);

  491.             if (m <= order) {
  492.                 // compute contribution of current order to field (equation 5 of the paper)

  493.                 // inner summation over degree, for fixed order
  494.                 double sumDegreeS               = 0;
  495.                 double sumDegreeC               = 0;
  496.                 double dSumDegreeSdR            = 0;
  497.                 double dSumDegreeCdR            = 0;
  498.                 double dSumDegreeSdTheta        = 0;
  499.                 double dSumDegreeCdTheta        = 0;
  500.                 double d2SumDegreeSdRdR         = 0;
  501.                 double d2SumDegreeSdRdTheta     = 0;
  502.                 double d2SumDegreeSdThetadTheta = 0;
  503.                 double d2SumDegreeCdRdR         = 0;
  504.                 double d2SumDegreeCdRdTheta     = 0;
  505.                 double d2SumDegreeCdThetadTheta = 0;
  506.                 for (int n = FastMath.max(2, m); n <= degree; ++n) {
  507.                     final double qSnm         = aOrN[n] * harmonics.getNormalizedSnm(n, m);
  508.                     final double qCnm         = aOrN[n] * harmonics.getNormalizedCnm(n, m);
  509.                     final double nOr          = n / r;
  510.                     final double nnP1Or2      = nOr * (n + 1) / r;
  511.                     final double s0           = pnm0[n] * qSnm;
  512.                     final double c0           = pnm0[n] * qCnm;
  513.                     final double s1           = pnm1[n] * qSnm;
  514.                     final double c1           = pnm1[n] * qCnm;
  515.                     final double s2           = pnm2[n] * qSnm;
  516.                     final double c2           = pnm2[n] * qCnm;
  517.                     sumDegreeS               += s0;
  518.                     sumDegreeC               += c0;
  519.                     dSumDegreeSdR            -= nOr * s0;
  520.                     dSumDegreeCdR            -= nOr * c0;
  521.                     dSumDegreeSdTheta        += s1;
  522.                     dSumDegreeCdTheta        += c1;
  523.                     d2SumDegreeSdRdR         += nnP1Or2 * s0;
  524.                     d2SumDegreeSdRdTheta     -= nOr * s1;
  525.                     d2SumDegreeSdThetadTheta += s2;
  526.                     d2SumDegreeCdRdR         += nnP1Or2 * c0;
  527.                     d2SumDegreeCdRdTheta     -= nOr * c1;
  528.                     d2SumDegreeCdThetadTheta += c2;
  529.                 }

  530.                 // contribution to outer summation over order
  531.                 final double sML = cosSinLambda[1][m];
  532.                 final double cML = cosSinLambda[0][m];
  533.                 value            = value         * u + sML * sumDegreeS + cML * sumDegreeC;
  534.                 gradient[0]      = gradient[0]   * u + sML * dSumDegreeSdR + cML * dSumDegreeCdR;
  535.                 gradient[1]      = gradient[1]   * u + m * (cML * sumDegreeS - sML * sumDegreeC);
  536.                 gradient[2]      = gradient[2]   * u + sML * dSumDegreeSdTheta + cML * dSumDegreeCdTheta;
  537.                 hessian[0][0]    = hessian[0][0] * u + sML * d2SumDegreeSdRdR + cML * d2SumDegreeCdRdR;
  538.                 hessian[1][0]    = hessian[1][0] * u + m * (cML * dSumDegreeSdR - sML * dSumDegreeCdR);
  539.                 hessian[2][0]    = hessian[2][0] * u + sML * d2SumDegreeSdRdTheta + cML * d2SumDegreeCdRdTheta;
  540.                 hessian[1][1]    = hessian[1][1] * u - m * m * (sML * sumDegreeS + cML * sumDegreeC);
  541.                 hessian[2][1]    = hessian[2][1] * u + m * (cML * dSumDegreeSdTheta - sML * dSumDegreeCdTheta);
  542.                 hessian[2][2]    = hessian[2][2] * u + sML * d2SumDegreeSdThetadTheta + cML * d2SumDegreeCdThetadTheta;

  543.             }

  544.             // rotate the recursion arrays
  545.             final double[] tmp0 = pnm0Plus2;
  546.             pnm0Plus2 = pnm0Plus1;
  547.             pnm0Plus1 = pnm0;
  548.             pnm0      = tmp0;
  549.             final double[] tmp1 = pnm1Plus1;
  550.             pnm1Plus1 = pnm1;
  551.             pnm1      = tmp1;

  552.         }

  553.         // scale back
  554.         value = FastMath.scalb(value, SCALING);
  555.         for (int i = 0; i < 3; ++i) {
  556.             gradient[i] = FastMath.scalb(gradient[i], SCALING);
  557.             for (int j = 0; j <= i; ++j) {
  558.                 hessian[i][j] = FastMath.scalb(hessian[i][j], SCALING);
  559.             }
  560.         }


  561.         // apply the global mu/r factor
  562.         final double muOr = mu / r;
  563.         value         *= muOr;
  564.         gradient[0]    = muOr * gradient[0] - value / r;
  565.         gradient[1]   *= muOr;
  566.         gradient[2]   *= muOr;
  567.         hessian[0][0]  = muOr * hessian[0][0] - 2 * gradient[0] / r;
  568.         hessian[1][0]  = muOr * hessian[1][0] -     gradient[1] / r;
  569.         hessian[2][0]  = muOr * hessian[2][0] -     gradient[2] / r;
  570.         hessian[1][1] *= muOr;
  571.         hessian[2][1] *= muOr;
  572.         hessian[2][2] *= muOr;

  573.         // convert gradient and Hessian from spherical to Cartesian
  574.         final SphericalCoordinates sc = new SphericalCoordinates(position);
  575.         return new GradientHessian(sc.toCartesianGradient(gradient),
  576.                                    sc.toCartesianHessian(hessian, gradient));


  577.     }

  578.     /** Container for gradient and Hessian. */
  579.     private static class GradientHessian {

  580.         /** Gradient. */
  581.         private final double[] gradient;

  582.         /** Hessian. */
  583.         private final double[][] hessian;

  584.         /** Simple constructor.
  585.          * <p>
  586.          * A reference to the arrays is stored, they are <strong>not</strong> cloned.
  587.          * </p>
  588.          * @param gradient gradient
  589.          * @param hessian hessian
  590.          */
  591.         GradientHessian(final double[] gradient, final double[][] hessian) {
  592.             this.gradient = gradient;
  593.             this.hessian  = hessian;
  594.         }

  595.         /** Get a reference to the gradient.
  596.          * @return gradient (a reference to the internal array is returned)
  597.          */
  598.         public double[] getGradient() {
  599.             return gradient;
  600.         }

  601.         /** Get a reference to the Hessian.
  602.          * @return Hessian (a reference to the internal array is returned)
  603.          */
  604.         public double[][] getHessian() {
  605.             return hessian;
  606.         }

  607.     }

  608.     /** Compute a/r powers array.
  609.      * @param aOr a/r
  610.      * @return array containing (a/r)<sup>n</sup>
  611.      */
  612.     private double[] createDistancePowersArray(final double aOr) {

  613.         // initialize array
  614.         final double[] aOrN = new double[provider.getMaxDegree() + 1];
  615.         aOrN[0] = 1;
  616.         aOrN[1] = aOr;

  617.         // fill up array
  618.         for (int n = 2; n < aOrN.length; ++n) {
  619.             final int p = n / 2;
  620.             final int q = n - p;
  621.             aOrN[n] = aOrN[p] * aOrN[q];
  622.         }

  623.         return aOrN;

  624.     }
  625.     /** Compute a/r powers array.
  626.      * @param aOr a/r
  627.      * @param <T> type of field used
  628.      * @return array containing (a/r)<sup>n</sup>
  629.      */
  630.     private <T extends RealFieldElement<T>> T[] createDistancePowersArray(final T aOr) {

  631.         // initialize array
  632.         final T[] aOrN = MathArrays.buildArray(aOr.getField(), provider.getMaxDegree() + 1);
  633.         aOrN[0] = aOr.getField().getOne();
  634.         aOrN[1] = aOr;

  635.         // fill up array
  636.         for (int n = 2; n < aOrN.length; ++n) {
  637.             final int p = n / 2;
  638.             final int q = n - p;
  639.             aOrN[n] = aOrN[p].multiply(aOrN[q]);
  640.         }

  641.         return aOrN;

  642.     }

  643.     /** Compute longitude cosines and sines.
  644.      * @param cosLambda cos(λ)
  645.      * @param sinLambda sin(λ)
  646.      * @return array containing cos(m &times; λ) in row 0
  647.      * and sin(m &times; λ) in row 1
  648.      */
  649.     private double[][] createCosSinArrays(final double cosLambda, final double sinLambda) {

  650.         // initialize arrays
  651.         final double[][] cosSin = new double[2][provider.getMaxOrder() + 1];
  652.         cosSin[0][0] = 1;
  653.         cosSin[1][0] = 0;
  654.         if (provider.getMaxOrder() > 0) {
  655.             cosSin[0][1] = cosLambda;
  656.             cosSin[1][1] = sinLambda;

  657.             // fill up array
  658.             for (int m = 2; m < cosSin[0].length; ++m) {

  659.                 // m * lambda is split as p * lambda + q * lambda, trying to avoid
  660.                 // p or q being much larger than the other. This reduces the number of
  661.                 // intermediate results reused to compute each value, and hence should limit
  662.                 // as much as possible roundoff error accumulation
  663.                 // (this does not change the number of floating point operations)
  664.                 final int p = m / 2;
  665.                 final int q = m - p;

  666.                 cosSin[0][m] = cosSin[0][p] * cosSin[0][q] - cosSin[1][p] * cosSin[1][q];
  667.                 cosSin[1][m] = cosSin[1][p] * cosSin[0][q] + cosSin[0][p] * cosSin[1][q];
  668.             }
  669.         }

  670.         return cosSin;

  671.     }

  672.     /** Compute longitude cosines and sines.
  673.      * @param cosLambda cos(λ)
  674.      * @param sinLambda sin(λ)
  675.      * @param <T> type of field used
  676.      * @return array containing cos(m &times; λ) in row 0
  677.      * and sin(m &times; λ) in row 1
  678.      */
  679.     private <T extends RealFieldElement<T>> T[][] createCosSinArrays(final T cosLambda, final T sinLambda) {

  680.         final T one = cosLambda.getField().getOne();
  681.         final T zero = cosLambda.getField().getZero();
  682.         // initialize arrays
  683.         final T[][] cosSin = MathArrays.buildArray(one.getField(), 2, provider.getMaxOrder() + 1);
  684.         cosSin[0][0] = one;
  685.         cosSin[1][0] = zero;
  686.         if (provider.getMaxOrder() > 0) {
  687.             cosSin[0][1] = cosLambda;
  688.             cosSin[1][1] = sinLambda;

  689.             // fill up array
  690.             for (int m = 2; m < cosSin[0].length; ++m) {

  691.                 // m * lambda is split as p * lambda + q * lambda, trying to avoid
  692.                 // p or q being much larger than the other. This reduces the number of
  693.                 // intermediate results reused to compute each value, and hence should limit
  694.                 // as much as possible roundoff error accumulation
  695.                 // (this does not change the number of floating point operations)
  696.                 final int p = m / 2;
  697.                 final int q = m - p;

  698.                 cosSin[0][m] = cosSin[0][p].multiply(cosSin[0][q]).subtract(cosSin[1][p].multiply(cosSin[1][q]));
  699.                 cosSin[1][m] = cosSin[1][p].multiply(cosSin[0][q]).add(cosSin[0][p].multiply(cosSin[1][q]));

  700.             }
  701.         }

  702.         return cosSin;

  703.     }

  704.     /** Compute one order of tesseral terms.
  705.      * <p>
  706.      * This corresponds to equations 27 and 30 of the paper.
  707.      * </p>
  708.      * @param m current order
  709.      * @param degree max degree
  710.      * @param index index in the flattened array
  711.      * @param t cos(θ), where θ is the polar angle
  712.      * @param u sin(θ), where θ is the polar angle
  713.      * @param tOu t/u
  714.      * @param pnm0Plus2 array containing scaled P<sub>n,m+2</sub>/u<sup>m+2</sup>
  715.      * @param pnm0Plus1 array containing scaled P<sub>n,m+1</sub>/u<sup>m+1</sup>
  716.      * @param pnm1Plus1 array containing scaled dP<sub>n,m+1</sub>/u<sup>m+1</sup>
  717.      * (may be null if second derivatives are not needed)
  718.      * @param pnm0 array to fill with scaled P<sub>n,m</sub>/u<sup>m</sup>
  719.      * @param pnm1 array to fill with scaled dP<sub>n,m</sub>/u<sup>m</sup>
  720.      * (may be null if first derivatives are not needed)
  721.      * @param pnm2 array to fill with scaled d²P<sub>n,m</sub>/u<sup>m</sup>
  722.      * (may be null if second derivatives are not needed)
  723.      * @return new value for index
  724.      */
  725.     private int computeTesseral(final int m, final int degree, final int index,
  726.                                 final double t, final double u, final double tOu,
  727.                                 final double[] pnm0Plus2, final double[] pnm0Plus1, final double[] pnm1Plus1,
  728.                                 final double[] pnm0, final double[] pnm1, final double[] pnm2) {

  729.         final double u2 = u * u;

  730.         // initialize recursion from sectorial terms
  731.         int n = FastMath.max(2, m);
  732.         if (n == m) {
  733.             pnm0[n] = sectorial[n];
  734.             ++n;
  735.         }

  736.         // compute tesseral values
  737.         int localIndex = index;
  738.         while (n <= degree) {

  739.             // value (equation 27 of the paper)
  740.             pnm0[n] = gnmOj[localIndex] * t * pnm0Plus1[n] - hnmOj[localIndex] * u2 * pnm0Plus2[n];

  741.             ++localIndex;
  742.             ++n;

  743.         }

  744.         if (pnm1 != null) {

  745.             // initialize recursion from sectorial terms
  746.             n = FastMath.max(2, m);
  747.             if (n == m) {
  748.                 pnm1[n] = m * tOu * pnm0[n];
  749.                 ++n;
  750.             }

  751.             // compute tesseral values and derivatives with respect to polar angle
  752.             localIndex = index;
  753.             while (n <= degree) {

  754.                 // first derivative (equation 30 of the paper)
  755.                 pnm1[n] = m * tOu * pnm0[n] - enm[localIndex] * u * pnm0Plus1[n];

  756.                 ++localIndex;
  757.                 ++n;

  758.             }

  759.             if (pnm2 != null) {

  760.                 // initialize recursion from sectorial terms
  761.                 n = FastMath.max(2, m);
  762.                 if (n == m) {
  763.                     pnm2[n] = m * (tOu * pnm1[n] - pnm0[n] / u2);
  764.                     ++n;
  765.                 }

  766.                 // compute tesseral values and derivatives with respect to polar angle
  767.                 localIndex = index;
  768.                 while (n <= degree) {

  769.                     // second derivative (differential of equation 30 with respect to theta)
  770.                     pnm2[n] = m * (tOu * pnm1[n] - pnm0[n] / u2) - enm[localIndex] * u * pnm1Plus1[n];

  771.                     ++localIndex;
  772.                     ++n;

  773.                 }

  774.             }

  775.         }

  776.         return localIndex;

  777.     }

  778.     /** Compute one order of tesseral terms.
  779.      * <p>
  780.      * This corresponds to equations 27 and 30 of the paper.
  781.      * </p>
  782.      * @param m current order
  783.      * @param degree max degree
  784.      * @param index index in the flattened array
  785.      * @param t cos(θ), where θ is the polar angle
  786.      * @param u sin(θ), where θ is the polar angle
  787.      * @param tOu t/u
  788.      * @param pnm0Plus2 array containing scaled P<sub>n,m+2</sub>/u<sup>m+2</sup>
  789.      * @param pnm0Plus1 array containing scaled P<sub>n,m+1</sub>/u<sup>m+1</sup>
  790.      * @param pnm1Plus1 array containing scaled dP<sub>n,m+1</sub>/u<sup>m+1</sup>
  791.      * (may be null if second derivatives are not needed)
  792.      * @param pnm0 array to fill with scaled P<sub>n,m</sub>/u<sup>m</sup>
  793.      * @param pnm1 array to fill with scaled dP<sub>n,m</sub>/u<sup>m</sup>
  794.      * (may be null if first derivatives are not needed)
  795.      * @param pnm2 array to fill with scaled d²P<sub>n,m</sub>/u<sup>m</sup>
  796.      * (may be null if second derivatives are not needed)
  797.      * @param <T> instance of field element
  798.      * @return new value for index
  799.      */
  800.     private <T extends RealFieldElement<T>> int computeTesseral(final int m, final int degree, final int index,
  801.                                                                 final T t, final T u, final T tOu,
  802.                                                                 final T[] pnm0Plus2, final T[] pnm0Plus1, final T[] pnm1Plus1,
  803.                                                                 final T[] pnm0, final T[] pnm1, final T[] pnm2) {

  804.         final T u2 = u.multiply(u);
  805.         final T zero = u.getField().getZero();
  806.         // initialize recursion from sectorial terms
  807.         int n = FastMath.max(2, m);
  808.         if (n == m) {
  809.             pnm0[n] = zero.add(sectorial[n]);
  810.             ++n;
  811.         }

  812.         // compute tesseral values
  813.         int localIndex = index;
  814.         while (n <= degree) {

  815.             // value (equation 27 of the paper)
  816.             pnm0[n] = t.multiply(gnmOj[localIndex]).multiply(pnm0Plus1[n]).subtract(u2.multiply(pnm0Plus2[n]).multiply(hnmOj[localIndex]));
  817.             ++localIndex;
  818.             ++n;

  819.         }
  820.         if (pnm1 != null) {

  821.             // initialize recursion from sectorial terms
  822.             n = FastMath.max(2, m);
  823.             if (n == m) {
  824.                 pnm1[n] = tOu.multiply(m).multiply(pnm0[n]);
  825.                 ++n;
  826.             }

  827.             // compute tesseral values and derivatives with respect to polar angle
  828.             localIndex = index;
  829.             while (n <= degree) {

  830.                 // first derivative (equation 30 of the paper)
  831.                 pnm1[n] = tOu.multiply(m).multiply(pnm0[n]).subtract(u.multiply(enm[localIndex]).multiply(pnm0Plus1[n]));

  832.                 ++localIndex;
  833.                 ++n;

  834.             }

  835.             if (pnm2 != null) {

  836.                 // initialize recursion from sectorial terms
  837.                 n = FastMath.max(2, m);
  838.                 if (n == m) {
  839.                     pnm2[n] =   tOu.multiply(pnm1[n]).subtract(pnm0[n].divide(u2)).multiply(m);
  840.                     ++n;
  841.                 }

  842.                 // compute tesseral values and derivatives with respect to polar angle
  843.                 localIndex = index;
  844.                 while (n <= degree) {

  845.                     // second derivative (differential of equation 30 with respect to theta)
  846.                     pnm2[n] = tOu.multiply(pnm1[n]).subtract(pnm0[n].divide(u2)).multiply(m).subtract(u.multiply(pnm1Plus1[n]).multiply(enm[localIndex]));
  847.                     ++localIndex;
  848.                     ++n;

  849.                 }

  850.             }

  851.         }
  852.         return localIndex;

  853.     }

  854.     /** {@inheritDoc} */
  855.     @Override
  856.     public Vector3D acceleration(final SpacecraftState s, final double[] parameters) {

  857.         final double mu = parameters[0];

  858.         // get the position in body frame
  859.         final AbsoluteDate date       = s.getDate();
  860.         final Transform fromBodyFrame = bodyFrame.getTransformTo(s.getFrame(), date);
  861.         final Transform toBodyFrame   = fromBodyFrame.getInverse();
  862.         final Vector3D position       = toBodyFrame.transformPosition(s.getPVCoordinates().getPosition());

  863.         // gradient of the non-central part of the gravity field
  864.         return fromBodyFrame.transformVector(new Vector3D(gradient(date, position, mu)));

  865.     }

  866.     /** {@inheritDoc} */
  867.     public <T extends RealFieldElement<T>> FieldVector3D<T> acceleration(final FieldSpacecraftState<T> s,
  868.                                                                          final T[] parameters) {

  869.         final T mu = parameters[0];

  870.         // check for faster computation dedicated to derivatives with respect to state
  871.         if (isGradientStateDerivative(s)) {
  872.             @SuppressWarnings("unchecked")
  873.             final FieldVector3D<Gradient> p = (FieldVector3D<Gradient>) s.getPVCoordinates().getPosition();
  874.             @SuppressWarnings("unchecked")
  875.             final FieldVector3D<T> a = (FieldVector3D<T>) accelerationWrtState(s.getDate().toAbsoluteDate(),
  876.                                                                                s.getFrame(), p,
  877.                                                                                (Gradient) mu);
  878.             return a;
  879.         } else if (isDSStateDerivative(s)) {
  880.             @SuppressWarnings("unchecked")
  881.             final FieldVector3D<DerivativeStructure> p = (FieldVector3D<DerivativeStructure>) s.getPVCoordinates().getPosition();
  882.             @SuppressWarnings("unchecked")
  883.             final FieldVector3D<T> a = (FieldVector3D<T>) accelerationWrtState(s.getDate().toAbsoluteDate(),
  884.                                                                                s.getFrame(), p,
  885.                                                                                (DerivativeStructure) mu);
  886.             return a;
  887.         }

  888.         // get the position in body frame
  889.         final FieldAbsoluteDate<T> date          = s.getDate();
  890.         final Transform            fromBodyFrame = bodyFrame.getTransformTo(s.getFrame(), date.toAbsoluteDate());
  891.         final Transform            toBodyFrame   = fromBodyFrame.getInverse();
  892.         final FieldVector3D<T>     position      = toBodyFrame.transformPosition(s.getPVCoordinates().getPosition());

  893.         // gradient of the non-central part of the gravity field
  894.         return fromBodyFrame.transformVector(new FieldVector3D<>(gradient(date, position, mu)));

  895.     }

  896.     /** {@inheritDoc} */
  897.     public Stream<EventDetector> getEventsDetectors() {
  898.         return Stream.empty();
  899.     }

  900.     @Override
  901.     /** {@inheritDoc} */
  902.     public <T extends RealFieldElement<T>> Stream<FieldEventDetector<T>> getFieldEventsDetectors(final Field<T> field) {
  903.         return Stream.empty();
  904.     }

  905.     /** Check if a field state corresponds to derivatives with respect to state.
  906.      * @param state state to check
  907.      * @param <T> type of the filed elements
  908.      * @return true if state corresponds to derivatives with respect to state
  909.      * @since 9.0
  910.      */
  911.     private <T extends RealFieldElement<T>> boolean isDSStateDerivative(final FieldSpacecraftState<T> state) {
  912.         try {
  913.             final DerivativeStructure dsMass = (DerivativeStructure) state.getMass();
  914.             final int o = dsMass.getOrder();
  915.             final int p = dsMass.getFreeParameters();
  916.             if (o != 1 || p < 3) {
  917.                 return false;
  918.             }
  919.             @SuppressWarnings("unchecked")
  920.             final FieldPVCoordinates<DerivativeStructure> pv = (FieldPVCoordinates<DerivativeStructure>) state.getPVCoordinates();
  921.             return isVariable(pv.getPosition().getX(), 0) &&
  922.                    isVariable(pv.getPosition().getY(), 1) &&
  923.                    isVariable(pv.getPosition().getZ(), 2);
  924.         } catch (ClassCastException cce) {
  925.             return false;
  926.         }
  927.     }

  928.     /** Check if a field state corresponds to derivatives with respect to state.
  929.      * @param state state to check
  930.      * @param <T> type of the filed elements
  931.      * @return true if state corresponds to derivatives with respect to state
  932.      * @since 10.2
  933.      */
  934.     private <T extends RealFieldElement<T>> boolean isGradientStateDerivative(final FieldSpacecraftState<T> state) {
  935.         try {
  936.             final Gradient gMass = (Gradient) state.getMass();
  937.             final int p = gMass.getFreeParameters();
  938.             if (p < 3) {
  939.                 return false;
  940.             }
  941.             @SuppressWarnings("unchecked")
  942.             final FieldPVCoordinates<Gradient> pv = (FieldPVCoordinates<Gradient>) state.getPVCoordinates();
  943.             return isVariable(pv.getPosition().getX(), 0) &&
  944.                    isVariable(pv.getPosition().getY(), 1) &&
  945.                    isVariable(pv.getPosition().getZ(), 2);
  946.         } catch (ClassCastException cce) {
  947.             return false;
  948.         }
  949.     }

  950.     /** Check if a derivative represents a specified variable.
  951.      * @param ds derivative to check
  952.      * @param index index of the variable
  953.      * @return true if the derivative represents a specified variable
  954.      * @since 9.0
  955.      */
  956.     private boolean isVariable(final DerivativeStructure ds, final int index) {
  957.         final double[] derivatives = ds.getAllDerivatives();
  958.         boolean check = true;
  959.         for (int i = 1; i < derivatives.length; ++i) {
  960.             check &= derivatives[i] == ((index + 1 == i) ? 1.0 : 0.0);
  961.         }
  962.         return check;
  963.     }

  964.     /** Check if a derivative represents a specified variable.
  965.      * @param g derivative to check
  966.      * @param index index of the variable
  967.      * @return true if the derivative represents a specified variable
  968.      * @since 10.2
  969.      */
  970.     private boolean isVariable(final Gradient g, final int index) {
  971.         final double[] derivatives = g.getGradient();
  972.         boolean check = true;
  973.         for (int i = 0; i < derivatives.length; ++i) {
  974.             check &= derivatives[i] == ((index == i) ? 1.0 : 0.0);
  975.         }
  976.         return check;
  977.     }

  978.     /** Compute acceleration derivatives with respect to state parameters.
  979.      * <p>
  980.      * From a theoretical point of view, this method computes the same values
  981.      * as {@link #acceleration(FieldSpacecraftState, RealFieldElement[])} in the
  982.      * specific case of {@link DerivativeStructure} with respect to state, so
  983.      * it is less general. However, it is *much* faster in this important case.
  984.      * <p>
  985.      * <p>
  986.      * The derivatives should be computed with respect to position. The input
  987.      * parameters already take into account the free parameters (6 or 7 depending
  988.      * on derivation with respect to mass being considered or not) and order
  989.      * (always 1). Free parameters at indices 0, 1 and 2 correspond to derivatives
  990.      * with respect to position. Free parameters at indices 3, 4 and 5 correspond
  991.      * to derivatives with respect to velocity (these derivatives will remain zero
  992.      * as acceleration due to gravity does not depend on velocity). Free parameter
  993.      * at index 6 (if present) corresponds to to derivatives with respect to mass
  994.      * (this derivative will remain zero as acceleration due to gravity does not
  995.      * depend on mass).
  996.      * </p>
  997.      * @param date current date
  998.      * @param frame inertial reference frame for state (both orbit and attitude)
  999.      * @param position position of spacecraft in inertial frame
  1000.      * @param mu central attraction coefficient to use
  1001.      * @return acceleration with all derivatives specified by the input parameters
  1002.      * own derivatives
  1003.      * @since 6.0
  1004.      */
  1005.     private FieldVector3D<DerivativeStructure> accelerationWrtState(final AbsoluteDate date, final Frame frame,
  1006.                                                                     final FieldVector3D<DerivativeStructure> position,
  1007.                                                                     final DerivativeStructure mu) {

  1008.         // free parameters
  1009.         final int freeParameters = mu.getFreeParameters();

  1010.         // get the position in body frame
  1011.         final Transform fromBodyFrame = bodyFrame.getTransformTo(frame, date);
  1012.         final Transform toBodyFrame   = fromBodyFrame.getInverse();
  1013.         final Vector3D positionBody   = toBodyFrame.transformPosition(position.toVector3D());

  1014.         // compute gradient and Hessian
  1015.         final GradientHessian gh   = gradientHessian(date, positionBody, mu.getReal());

  1016.         // gradient of the non-central part of the gravity field
  1017.         final double[] gInertial = fromBodyFrame.transformVector(new Vector3D(gh.getGradient())).toArray();

  1018.         // Hessian of the non-central part of the gravity field
  1019.         final RealMatrix hBody     = new Array2DRowRealMatrix(gh.getHessian(), false);
  1020.         final RealMatrix rot       = new Array2DRowRealMatrix(toBodyFrame.getRotation().getMatrix());
  1021.         final RealMatrix hInertial = rot.transpose().multiply(hBody).multiply(rot);

  1022.         // distribute all partial derivatives in a compact acceleration vector
  1023.         final double[] derivatives = new double[freeParameters + 1];
  1024.         final DerivativeStructure[] accDer = new DerivativeStructure[3];
  1025.         for (int i = 0; i < 3; ++i) {

  1026.             // first element is value of acceleration (i.e. gradient of field)
  1027.             derivatives[0] = gInertial[i];

  1028.             // Jacobian of acceleration (i.e. Hessian of field)
  1029.             derivatives[1] = hInertial.getEntry(i, 0);
  1030.             derivatives[2] = hInertial.getEntry(i, 1);
  1031.             derivatives[3] = hInertial.getEntry(i, 2);

  1032.             // next element is derivative with respect to parameter mu
  1033.             if (derivatives.length > 4 && isVariable(mu, 3)) {
  1034.                 derivatives[4] = gInertial[i] / mu.getReal();
  1035.             }

  1036.             accDer[i] = position.getX().getFactory().build(derivatives);

  1037.         }

  1038.         return new FieldVector3D<>(accDer);

  1039.     }

  1040.     /** Compute acceleration derivatives with respect to state parameters.
  1041.      * <p>
  1042.      * From a theoretical point of view, this method computes the same values
  1043.      * as {@link #acceleration(FieldSpacecraftState, RealFieldElement[])} in the
  1044.      * specific case of {@link DerivativeStructure} with respect to state, so
  1045.      * it is less general. However, it is *much* faster in this important case.
  1046.      * <p>
  1047.      * <p>
  1048.      * The derivatives should be computed with respect to position. The input
  1049.      * parameters already take into account the free parameters (6 or 7 depending
  1050.      * on derivation with respect to mass being considered or not) and order
  1051.      * (always 1). Free parameters at indices 0, 1 and 2 correspond to derivatives
  1052.      * with respect to position. Free parameters at indices 3, 4 and 5 correspond
  1053.      * to derivatives with respect to velocity (these derivatives will remain zero
  1054.      * as acceleration due to gravity does not depend on velocity). Free parameter
  1055.      * at index 6 (if present) corresponds to to derivatives with respect to mass
  1056.      * (this derivative will remain zero as acceleration due to gravity does not
  1057.      * depend on mass).
  1058.      * </p>
  1059.      * @param date current date
  1060.      * @param frame inertial reference frame for state (both orbit and attitude)
  1061.      * @param position position of spacecraft in inertial frame
  1062.      * @param mu central attraction coefficient to use
  1063.      * @return acceleration with all derivatives specified by the input parameters
  1064.      * own derivatives
  1065.      * @since 10.2
  1066.      */
  1067.     private FieldVector3D<Gradient> accelerationWrtState(final AbsoluteDate date, final Frame frame,
  1068.                                                          final FieldVector3D<Gradient> position,
  1069.                                                          final Gradient mu) {

  1070.         // free parameters
  1071.         final int freeParameters = mu.getFreeParameters();

  1072.         // get the position in body frame
  1073.         final Transform fromBodyFrame = bodyFrame.getTransformTo(frame, date);
  1074.         final Transform toBodyFrame   = fromBodyFrame.getInverse();
  1075.         final Vector3D positionBody   = toBodyFrame.transformPosition(position.toVector3D());

  1076.         // compute gradient and Hessian
  1077.         final GradientHessian gh   = gradientHessian(date, positionBody, mu.getReal());

  1078.         // gradient of the non-central part of the gravity field
  1079.         final double[] gInertial = fromBodyFrame.transformVector(new Vector3D(gh.getGradient())).toArray();

  1080.         // Hessian of the non-central part of the gravity field
  1081.         final RealMatrix hBody     = new Array2DRowRealMatrix(gh.getHessian(), false);
  1082.         final RealMatrix rot       = new Array2DRowRealMatrix(toBodyFrame.getRotation().getMatrix());
  1083.         final RealMatrix hInertial = rot.transpose().multiply(hBody).multiply(rot);

  1084.         // distribute all partial derivatives in a compact acceleration vector
  1085.         final double[] derivatives = new double[freeParameters];
  1086.         final Gradient[] accDer = new Gradient[3];
  1087.         for (int i = 0; i < 3; ++i) {

  1088.             // Jacobian of acceleration (i.e. Hessian of field)
  1089.             derivatives[0] = hInertial.getEntry(i, 0);
  1090.             derivatives[1] = hInertial.getEntry(i, 1);
  1091.             derivatives[2] = hInertial.getEntry(i, 2);

  1092.             // next element is derivative with respect to parameter mu
  1093.             if (derivatives.length > 3 && isVariable(mu, 3)) {
  1094.                 derivatives[3] = gInertial[i] / mu.getReal();
  1095.             }

  1096.             accDer[i] = new Gradient(gInertial[i], derivatives);

  1097.         }

  1098.         return new FieldVector3D<>(accDer);

  1099.     }

  1100.     /** {@inheritDoc} */
  1101.     public ParameterDriver[] getParametersDrivers() {
  1102.         return new ParameterDriver[] {
  1103.             gmParameterDriver
  1104.         };
  1105.     }

  1106. }