HolmesFeatherstoneAttractionModel.java

  1. /* Copyright 2002-2019 CS Systèmes d'Information
  2.  * Licensed to CS Systèmes d'Information (CS) under one or more
  3.  * contributor license agreements.  See the NOTICE file distributed with
  4.  * this work for additional information regarding copyright ownership.
  5.  * CS licenses this file to You under the Apache License, Version 2.0
  6.  * (the "License"); you may not use this file except in compliance with
  7.  * the License.  You may obtain a copy of the License at
  8.  *
  9.  *   http://www.apache.org/licenses/LICENSE-2.0
  10.  *
  11.  * Unless required by applicable law or agreed to in writing, software
  12.  * distributed under the License is distributed on an "AS IS" BASIS,
  13.  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
  14.  * See the License for the specific language governing permissions and
  15.  * limitations under the License.
  16.  */
  17. package org.orekit.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.geometry.euclidean.threed.FieldVector3D;
  23. import org.hipparchus.geometry.euclidean.threed.SphericalCoordinates;
  24. import org.hipparchus.geometry.euclidean.threed.Vector3D;
  25. import org.hipparchus.linear.Array2DRowRealMatrix;
  26. import org.hipparchus.linear.RealMatrix;
  27. import org.hipparchus.util.FastMath;
  28. import org.hipparchus.util.MathArrays;
  29. import org.orekit.errors.OrekitException;
  30. import org.orekit.errors.OrekitInternalError;
  31. import org.orekit.forces.AbstractForceModel;
  32. import org.orekit.forces.gravity.potential.NormalizedSphericalHarmonicsProvider;
  33. import org.orekit.forces.gravity.potential.NormalizedSphericalHarmonicsProvider.NormalizedSphericalHarmonics;
  34. import org.orekit.forces.gravity.potential.TideSystem;
  35. import org.orekit.forces.gravity.potential.TideSystemProvider;
  36. import org.orekit.frames.Frame;
  37. import org.orekit.frames.Transform;
  38. import org.orekit.propagation.FieldSpacecraftState;
  39. import org.orekit.propagation.SpacecraftState;
  40. import org.orekit.propagation.events.EventDetector;
  41. import org.orekit.propagation.events.FieldEventDetector;
  42. import org.orekit.time.AbsoluteDate;
  43. import org.orekit.time.FieldAbsoluteDate;
  44. import org.orekit.utils.FieldPVCoordinates;
  45. import org.orekit.utils.ParameterDriver;

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

  72. public class HolmesFeatherstoneAttractionModel extends AbstractForceModel implements TideSystemProvider {

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

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

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

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

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

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

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

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

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

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

  106.         try {
  107.             gmParameterDriver = new ParameterDriver(NewtonianAttraction.CENTRAL_ATTRACTION_COEFFICIENT,
  108.                                                     provider.getMu(), MU_SCALE, 0.0, Double.POSITIVE_INFINITY);
  109.         } catch (OrekitException oe) {
  110.             // this should never occur as valueChanged above never throws an exception
  111.             throw new OrekitInternalError(oe);
  112.         }

  113.         this.provider  = provider;
  114.         this.bodyFrame = centralBodyFrame;

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

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

  137.         // scaled sectorial terms corresponding to equation 28 in Holmes and Featherstone paper
  138.         sectorial    = new double[degree + 1];
  139.         sectorial[0] = FastMath.scalb(1.0, -SCALING);
  140.         sectorial[1] = FastMath.sqrt(3) * sectorial[0];
  141.         for (int m = 2; m < sectorial.length; ++m) {
  142.             sectorial[m] = FastMath.sqrt((2 * m + 1) / (2.0 * m)) * sectorial[m - 1];
  143.         }

  144.     }

  145.     /** {@inheritDoc} */
  146.     @Override
  147.     public boolean dependsOnPositionOnly() {
  148.         return true;
  149.     }

  150.     /** {@inheritDoc} */
  151.     public TideSystem getTideSystem() {
  152.         return provider.getTideSystem();
  153.     }

  154.     /** Get the central attraction coefficient μ.
  155.      * @return mu central attraction coefficient (m³/s²)
  156.      */
  157.     public double getMu() {
  158.         return gmParameterDriver.getValue();
  159.     }

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

  170.     /** Compute the non-central part of the gravity field.
  171.      * @param date current date
  172.      * @param position position at which gravity field is desired in body frame
  173.      * @param mu central attraction coefficient to use
  174.      * @return value of the non-central part of the gravity field
  175.      */
  176.     public double nonCentralPart(final AbsoluteDate date, final Vector3D position, final double mu) {

  177.         final int degree = provider.getMaxDegree();
  178.         final int order  = provider.getMaxOrder();
  179.         final NormalizedSphericalHarmonics harmonics = provider.onDate(date);

  180.         // allocate the columns for recursion
  181.         double[] pnm0Plus2 = new double[degree + 1];
  182.         double[] pnm0Plus1 = new double[degree + 1];
  183.         double[] pnm0      = new double[degree + 1];

  184.         // compute polar coordinates
  185.         final double x   = position.getX();
  186.         final double y   = position.getY();
  187.         final double z   = position.getZ();
  188.         final double x2  = x * x;
  189.         final double y2  = y * y;
  190.         final double z2  = z * z;
  191.         final double r2  = x2 + y2 + z2;
  192.         final double r   = FastMath.sqrt (r2);
  193.         final double rho = FastMath.sqrt(x2 + y2);
  194.         final double t   = z / r;   // cos(theta), where theta is the polar angle
  195.         final double u   = rho / r; // sin(theta), where theta is the polar angle
  196.         final double tOu = z / rho;

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

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

  201.         // outer summation over order
  202.         int    index = 0;
  203.         double value = 0;
  204.         for (int m = degree; m >= 0; --m) {

  205.             // compute tesseral terms without derivatives
  206.             index = computeTesseral(m, degree, index, t, u, tOu,
  207.                                     pnm0Plus2, pnm0Plus1, null, pnm0, null, null);

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

  210.                 // inner summation over degree, for fixed order
  211.                 double sumDegreeS        = 0;
  212.                 double sumDegreeC        = 0;
  213.                 for (int n = FastMath.max(2, m); n <= degree; ++n) {
  214.                     sumDegreeS += pnm0[n] * aOrN[n] * harmonics.getNormalizedSnm(n, m);
  215.                     sumDegreeC += pnm0[n] * aOrN[n] * harmonics.getNormalizedCnm(n, m);
  216.                 }

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

  219.             }

  220.             // rotate the recursion arrays
  221.             final double[] tmp = pnm0Plus2;
  222.             pnm0Plus2 = pnm0Plus1;
  223.             pnm0Plus1 = pnm0;
  224.             pnm0      = tmp;

  225.         }

  226.         // scale back
  227.         value = FastMath.scalb(value, SCALING);

  228.         // apply the global mu/r factor
  229.         return mu * value / r;

  230.     }

  231.     /** Compute the gradient of the non-central part of the gravity field.
  232.      * @param date current date
  233.      * @param position position at which gravity field is desired in body frame
  234.      * @param mu central attraction coefficient to use
  235.      * @return gradient of the non-central part of the gravity field
  236.      */
  237.     public double[] gradient(final AbsoluteDate date, final Vector3D position, final double mu) {

  238.         final int degree = provider.getMaxDegree();
  239.         final int order  = provider.getMaxOrder();
  240.         final NormalizedSphericalHarmonics harmonics = provider.onDate(date);

  241.         // allocate the columns for recursion
  242.         double[] pnm0Plus2  = new double[degree + 1];
  243.         double[] pnm0Plus1  = new double[degree + 1];
  244.         double[] pnm0       = new double[degree + 1];
  245.         final double[] pnm1 = new double[degree + 1];

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

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

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

  264.         // outer summation over order
  265.         int    index = 0;
  266.         double value = 0;
  267.         final double[] gradient = new double[3];
  268.         for (int m = degree; m >= 0; --m) {

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

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

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

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

  306.             }

  307.             // rotate the recursion arrays
  308.             final double[] tmp = pnm0Plus2;
  309.             pnm0Plus2 = pnm0Plus1;
  310.             pnm0Plus1 = pnm0;
  311.             pnm0      = tmp;

  312.         }

  313.         // scale back
  314.         value       = FastMath.scalb(value,       SCALING);
  315.         gradient[0] = FastMath.scalb(gradient[0], SCALING);
  316.         gradient[1] = FastMath.scalb(gradient[1], SCALING);
  317.         gradient[2] = FastMath.scalb(gradient[2], SCALING);

  318.         // apply the global mu/r factor
  319.         final double muOr = mu / r;
  320.         value            *= muOr;
  321.         gradient[0]       = muOr * gradient[0] - value / r;
  322.         gradient[1]      *= muOr;
  323.         gradient[2]      *= muOr;

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

  326.     }

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

  336.         final int degree = provider.getMaxDegree();
  337.         final int order  = provider.getMaxOrder();
  338.         final NormalizedSphericalHarmonics harmonics = provider.onDate(date.toAbsoluteDate());
  339.         final T zero = date.getField().getZero();
  340.         // allocate the columns for recursion
  341.         T[] pnm0Plus2  = MathArrays.buildArray(date.getField(), degree + 1);
  342.         T[] pnm0Plus1  = MathArrays.buildArray(date.getField(), degree + 1);
  343.         T[] pnm0       = MathArrays.buildArray(date.getField(), degree + 1);
  344.         final T[] pnm1 = MathArrays.buildArray(date.getField(), degree + 1);

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

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

  361.         // compute longitude cosines/sines
  362.         final T[][] cosSinLambda = createCosSinArrays(rho.reciprocal().multiply(position.getX()), rho.reciprocal().multiply(position.getY()));
  363.         // outer summation over order
  364.         int    index = 0;
  365.         T value = zero;
  366.         final T[] gradient = MathArrays.buildArray(zero.getField(), 3);
  367.         for (int m = degree; m >= 0; --m) {

  368.             // compute tesseral terms with derivatives
  369.             index = computeTesseral(m, degree, index, t, u, tOu,
  370.                                     pnm0Plus2, pnm0Plus1, null, pnm0, pnm1, null);
  371.             if (m <= order) {
  372.                 // compute contribution of current order to field (equation 5 of the paper)

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

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

  411.         }
  412.         // scale back
  413.         value       = value.scalb(SCALING);
  414.         gradient[0] = gradient[0].scalb(SCALING);
  415.         gradient[1] = gradient[1].scalb(SCALING);
  416.         gradient[2] = gradient[2].scalb(SCALING);

  417.         // apply the global mu/r factor
  418.         final T muOr = r.reciprocal().multiply(mu);
  419.         value            = value.multiply(muOr);
  420.         gradient[0]       = muOr.multiply(gradient[0]).subtract(value.divide(r));
  421.         gradient[1]      = gradient[1].multiply(muOr);
  422.         gradient[2]      = gradient[2].multiply(muOr);

  423.         // convert gradient from spherical to Cartesian
  424.         // Cartesian coordinates
  425.         // remaining spherical coordinates
  426.         final T rPos     = position.getNorm();
  427.         // intermediate variables
  428.         final T xPos    = position.getX();
  429.         final T yPos    = position.getY();
  430.         final T zPos    = position.getZ();
  431.         final T rho2Pos = x.multiply(x).add(y.multiply(y));
  432.         final T rhoPos  = rho2.sqrt();
  433.         final T r2Pos   = rho2.add(z.multiply(z));

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

  435.         // row representing the gradient of r
  436.         jacobianPos[0][0] = xPos.divide(rPos);
  437.         jacobianPos[0][1] = yPos.divide(rPos);
  438.         jacobianPos[0][2] = zPos.divide(rPos);

  439.         // row representing the gradient of theta
  440.         jacobianPos[1][0] =  yPos.negate().divide(rho2Pos);
  441.         jacobianPos[1][1] =  xPos.divide(rho2Pos);
  442.         // jacobian[1][2] is already set to 0 at allocation time

  443.         // row representing the gradient of phi
  444.         jacobianPos[2][0] = xPos.multiply(zPos).divide(rhoPos.multiply(r2Pos));
  445.         jacobianPos[2][1] = yPos.multiply(zPos).divide(rhoPos.multiply(r2Pos));
  446.         jacobianPos[2][2] = rhoPos.negate().divide(r2Pos);
  447.         final T[] cartGradPos = MathArrays.buildArray(zero.getField(), 3);
  448.         cartGradPos[0] = gradient[0].multiply(jacobianPos[0][0]).add(gradient[1].multiply(jacobianPos[1][0])).add(gradient[2].multiply(jacobianPos[2][0]));
  449.         cartGradPos[1] = gradient[0].multiply(jacobianPos[0][1]).add(gradient[1].multiply(jacobianPos[1][1])).add(gradient[2].multiply(jacobianPos[2][1]));
  450.         cartGradPos[2] = gradient[0].multiply(jacobianPos[0][2])                                      .add(gradient[2].multiply(jacobianPos[2][2]));
  451.         return cartGradPos;

  452.     }

  453.     /** Compute both the gradient and the hessian of the non-central part of the gravity field.
  454.      * @param date current date
  455.      * @param position position at which gravity field is desired in body frame
  456.      * @param mu central attraction coefficient to use
  457.      * @return gradient and hessian of the non-central part of the gravity field
  458.      */
  459.     private GradientHessian gradientHessian(final AbsoluteDate date, final Vector3D position, final double mu) {

  460.         final int degree = provider.getMaxDegree();
  461.         final int order  = provider.getMaxOrder();
  462.         final NormalizedSphericalHarmonics harmonics = provider.onDate(date);

  463.         // allocate the columns for recursion
  464.         double[] pnm0Plus2  = new double[degree + 1];
  465.         double[] pnm0Plus1  = new double[degree + 1];
  466.         double[] pnm0       = new double[degree + 1];
  467.         double[] pnm1Plus1  = new double[degree + 1];
  468.         double[] pnm1       = new double[degree + 1];
  469.         final double[] pnm2 = new double[degree + 1];

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

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

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

  488.         // outer summation over order
  489.         int    index = 0;
  490.         double value = 0;
  491.         final double[]   gradient = new double[3];
  492.         final double[][] hessian  = new double[3][3];
  493.         for (int m = degree; m >= 0; --m) {

  494.             // compute tesseral terms
  495.             index = computeTesseral(m, degree, index, t, u, tOu,
  496.                                     pnm0Plus2, pnm0Plus1, pnm1Plus1, pnm0, pnm1, pnm2);

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

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

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

  549.             }

  550.             // rotate the recursion arrays
  551.             final double[] tmp0 = pnm0Plus2;
  552.             pnm0Plus2 = pnm0Plus1;
  553.             pnm0Plus1 = pnm0;
  554.             pnm0      = tmp0;
  555.             final double[] tmp1 = pnm1Plus1;
  556.             pnm1Plus1 = pnm1;
  557.             pnm1      = tmp1;

  558.         }

  559.         // scale back
  560.         value = FastMath.scalb(value, SCALING);
  561.         for (int i = 0; i < 3; ++i) {
  562.             gradient[i] = FastMath.scalb(gradient[i], SCALING);
  563.             for (int j = 0; j <= i; ++j) {
  564.                 hessian[i][j] = FastMath.scalb(hessian[i][j], SCALING);
  565.             }
  566.         }


  567.         // apply the global mu/r factor
  568.         final double muOr = mu / r;
  569.         value         *= muOr;
  570.         gradient[0]    = muOr * gradient[0] - value / r;
  571.         gradient[1]   *= muOr;
  572.         gradient[2]   *= muOr;
  573.         hessian[0][0]  = muOr * hessian[0][0] - 2 * gradient[0] / r;
  574.         hessian[1][0]  = muOr * hessian[1][0] -     gradient[1] / r;
  575.         hessian[2][0]  = muOr * hessian[2][0] -     gradient[2] / r;
  576.         hessian[1][1] *= muOr;
  577.         hessian[2][1] *= muOr;
  578.         hessian[2][2] *= muOr;

  579.         // convert gradient and Hessian from spherical to Cartesian
  580.         final SphericalCoordinates sc = new SphericalCoordinates(position);
  581.         return new GradientHessian(sc.toCartesianGradient(gradient),
  582.                                    sc.toCartesianHessian(hessian, gradient));


  583.     }

  584.     /** Container for gradient and Hessian. */
  585.     private static class GradientHessian {

  586.         /** Gradient. */
  587.         private final double[] gradient;

  588.         /** Hessian. */
  589.         private final double[][] hessian;

  590.         /** Simple constructor.
  591.          * <p>
  592.          * A reference to the arrays is stored, they are <strong>not</strong> cloned.
  593.          * </p>
  594.          * @param gradient gradient
  595.          * @param hessian hessian
  596.          */
  597.         GradientHessian(final double[] gradient, final double[][] hessian) {
  598.             this.gradient = gradient;
  599.             this.hessian  = hessian;
  600.         }

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

  607.         /** Get a reference to the Hessian.
  608.          * @return Hessian (a reference to the internal array is returned)
  609.          */
  610.         public double[][] getHessian() {
  611.             return hessian;
  612.         }

  613.     }

  614.     /** Compute a/r powers array.
  615.      * @param aOr a/r
  616.      * @return array containing (a/r)<sup>n</sup>
  617.      */
  618.     private double[] createDistancePowersArray(final double aOr) {

  619.         // initialize array
  620.         final double[] aOrN = new double[provider.getMaxDegree() + 1];
  621.         aOrN[0] = 1;
  622.         aOrN[1] = aOr;

  623.         // fill up array
  624.         for (int n = 2; n < aOrN.length; ++n) {
  625.             final int p = n / 2;
  626.             final int q = n - p;
  627.             aOrN[n] = aOrN[p] * aOrN[q];
  628.         }

  629.         return aOrN;

  630.     }
  631.     /** Compute a/r powers array.
  632.      * @param aOr a/r
  633.      * @param <T> type of field used
  634.      * @return array containing (a/r)<sup>n</sup>
  635.      */
  636.     private <T extends RealFieldElement<T>> T[] createDistancePowersArray(final T aOr) {

  637.         // initialize array
  638.         final T[] aOrN = MathArrays.buildArray(aOr.getField(), provider.getMaxDegree() + 1);
  639.         aOrN[0] = aOr.getField().getOne();
  640.         aOrN[1] = aOr;

  641.         // fill up array
  642.         for (int n = 2; n < aOrN.length; ++n) {
  643.             final int p = n / 2;
  644.             final int q = n - p;
  645.             aOrN[n] = aOrN[p].multiply(aOrN[q]);
  646.         }

  647.         return aOrN;

  648.     }

  649.     /** Compute longitude cosines and sines.
  650.      * @param cosLambda cos(λ)
  651.      * @param sinLambda sin(λ)
  652.      * @return array containing cos(m &times; λ) in row 0
  653.      * and sin(m &times; λ) in row 1
  654.      */
  655.     private double[][] createCosSinArrays(final double cosLambda, final double sinLambda) {

  656.         // initialize arrays
  657.         final double[][] cosSin = new double[2][provider.getMaxOrder() + 1];
  658.         cosSin[0][0] = 1;
  659.         cosSin[1][0] = 0;
  660.         if (provider.getMaxOrder() > 0) {
  661.             cosSin[0][1] = cosLambda;
  662.             cosSin[1][1] = sinLambda;

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

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

  672.                 cosSin[0][m] = cosSin[0][p] * cosSin[0][q] - cosSin[1][p] * cosSin[1][q];
  673.                 cosSin[1][m] = cosSin[1][p] * cosSin[0][q] + cosSin[0][p] * cosSin[1][q];
  674.             }
  675.         }

  676.         return cosSin;

  677.     }

  678.     /** Compute longitude cosines and sines.
  679.      * @param cosLambda cos(λ)
  680.      * @param sinLambda sin(λ)
  681.      * @param <T> type of field used
  682.      * @return array containing cos(m &times; λ) in row 0
  683.      * and sin(m &times; λ) in row 1
  684.      */
  685.     private <T extends RealFieldElement<T>> T[][] createCosSinArrays(final T cosLambda, final T sinLambda) {

  686.         final T one = cosLambda.getField().getOne();
  687.         final T zero = cosLambda.getField().getZero();
  688.         // initialize arrays
  689.         final T[][] cosSin = MathArrays.buildArray(one.getField(), 2, provider.getMaxOrder() + 1);
  690.         cosSin[0][0] = one;
  691.         cosSin[1][0] = zero;
  692.         if (provider.getMaxOrder() > 0) {
  693.             cosSin[0][1] = cosLambda;
  694.             cosSin[1][1] = sinLambda;

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

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

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

  706.             }
  707.         }

  708.         return cosSin;

  709.     }

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

  735.         final double u2 = u * u;

  736.         // initialize recursion from sectorial terms
  737.         int n = FastMath.max(2, m);
  738.         if (n == m) {
  739.             pnm0[n] = sectorial[n];
  740.             ++n;
  741.         }

  742.         // compute tesseral values
  743.         int localIndex = index;
  744.         while (n <= degree) {

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

  747.             ++localIndex;
  748.             ++n;

  749.         }

  750.         if (pnm1 != null) {

  751.             // initialize recursion from sectorial terms
  752.             n = FastMath.max(2, m);
  753.             if (n == m) {
  754.                 pnm1[n] = m * tOu * pnm0[n];
  755.                 ++n;
  756.             }

  757.             // compute tesseral values and derivatives with respect to polar angle
  758.             localIndex = index;
  759.             while (n <= degree) {

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

  762.                 ++localIndex;
  763.                 ++n;

  764.             }

  765.             if (pnm2 != null) {

  766.                 // initialize recursion from sectorial terms
  767.                 n = FastMath.max(2, m);
  768.                 if (n == m) {
  769.                     pnm2[n] = m * (tOu * pnm1[n] - pnm0[n] / u2);
  770.                     ++n;
  771.                 }

  772.                 // compute tesseral values and derivatives with respect to polar angle
  773.                 localIndex = index;
  774.                 while (n <= degree) {

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

  777.                     ++localIndex;
  778.                     ++n;

  779.                 }

  780.             }

  781.         }

  782.         return localIndex;

  783.     }

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

  810.         final T u2 = u.multiply(u);
  811.         final T zero = u.getField().getZero();
  812.         // initialize recursion from sectorial terms
  813.         int n = FastMath.max(2, m);
  814.         if (n == m) {
  815.             pnm0[n] = zero.add(sectorial[n]);
  816.             ++n;
  817.         }

  818.         // compute tesseral values
  819.         int localIndex = index;
  820.         while (n <= degree) {

  821.             // value (equation 27 of the paper)
  822.             pnm0[n] = t.multiply(gnmOj[localIndex]).multiply(pnm0Plus1[n]).subtract(u2.multiply(pnm0Plus2[n]).multiply(hnmOj[localIndex]));
  823.             ++localIndex;
  824.             ++n;

  825.         }
  826.         if (pnm1 != null) {

  827.             // initialize recursion from sectorial terms
  828.             n = FastMath.max(2, m);
  829.             if (n == m) {
  830.                 pnm1[n] = tOu.multiply(m).multiply(pnm0[n]);
  831.                 ++n;
  832.             }

  833.             // compute tesseral values and derivatives with respect to polar angle
  834.             localIndex = index;
  835.             while (n <= degree) {

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

  838.                 ++localIndex;
  839.                 ++n;

  840.             }

  841.             if (pnm2 != null) {

  842.                 // initialize recursion from sectorial terms
  843.                 n = FastMath.max(2, m);
  844.                 if (n == m) {
  845.                     pnm2[n] =   tOu.multiply(pnm1[n]).subtract(pnm0[n].divide(u2)).multiply(m);
  846.                     ++n;
  847.                 }

  848.                 // compute tesseral values and derivatives with respect to polar angle
  849.                 localIndex = index;
  850.                 while (n <= degree) {

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

  855.                 }

  856.             }

  857.         }
  858.         return localIndex;

  859.     }

  860.     /** {@inheritDoc} */
  861.     @Override
  862.     public Vector3D acceleration(final SpacecraftState s, final double[] parameters) {

  863.         final double mu = parameters[0];

  864.         // get the position in body frame
  865.         final AbsoluteDate date       = s.getDate();
  866.         final Transform fromBodyFrame = bodyFrame.getTransformTo(s.getFrame(), date);
  867.         final Transform toBodyFrame   = fromBodyFrame.getInverse();
  868.         final Vector3D position       = toBodyFrame.transformPosition(s.getPVCoordinates().getPosition());

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

  871.     }

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

  875.         final T mu = parameters[0];

  876.         // check for faster computation dedicated to derivatives with respect to state
  877.         if (isStateDerivative(s)) {
  878.             @SuppressWarnings("unchecked")
  879.             final FieldVector3D<DerivativeStructure> p = (FieldVector3D<DerivativeStructure>) s.getPVCoordinates().getPosition();
  880.             @SuppressWarnings("unchecked")
  881.             final FieldVector3D<T> a = (FieldVector3D<T>) accelerationWrtState(s.getDate().toAbsoluteDate(),
  882.                                                                                s.getFrame(), p,
  883.                                                                                (DerivativeStructure) mu);
  884.             return a;
  885.         }

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

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

  893.     }

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

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

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

  926.     /** Check if a derivative represents a specified variable.
  927.      * @param ds derivative to check
  928.      * @param index index of the variable
  929.      * @return true if the derivative represents a specified variable
  930.      * @since 9.0
  931.      */
  932.     private boolean isVariable(final DerivativeStructure ds, final int index) {
  933.         final double[] derivatives = ds.getAllDerivatives();
  934.         boolean check = true;
  935.         for (int i = 1; i < derivatives.length; ++i) {
  936.             check &= derivatives[i] == ((index + 1 == i) ? 1.0 : 0.0);
  937.         }
  938.         return check;
  939.     }

  940.     /** Compute acceleration derivatives with respect to state parameters.
  941.      * <p>
  942.      * From a theoretical point of view, this method computes the same values
  943.      * as {@link #acceleration(FieldSpacecraftState, RealFieldElement[])} in the
  944.      * specific case of {@link DerivativeStructure} with respect to state, so
  945.      * it is less general. However, it is *much* faster in this important case.
  946.      * <p>
  947.      * <p>
  948.      * The derivatives should be computed with respect to position. The input
  949.      * parameters already take into account the free parameters (6 or 7 depending
  950.      * on derivation with respect to mass being considered or not) and order
  951.      * (always 1). Free parameters at indices 0, 1 and 2 correspond to derivatives
  952.      * with respect to position. Free parameters at indices 3, 4 and 5 correspond
  953.      * to derivatives with respect to velocity (these derivatives will remain zero
  954.      * as acceleration due to gravity does not depend on velocity). Free parameter
  955.      * at index 6 (if present) corresponds to to derivatives with respect to mass
  956.      * (this derivative will remain zero as acceleration due to gravity does not
  957.      * depend on mass).
  958.      * </p>
  959.      * @param date current date
  960.     * @param frame inertial reference frame for state (both orbit and attitude)
  961.      * @param position position of spacecraft in inertial frame
  962.      * @param mu central attraction coefficient to use
  963.      * @return acceleration with all derivatives specified by the input parameters
  964.      * own derivatives
  965.           * @since 6.0
  966.      */
  967.     private FieldVector3D<DerivativeStructure> accelerationWrtState(final AbsoluteDate date, final Frame frame,
  968.                                                                     final FieldVector3D<DerivativeStructure> position,
  969.                                                                     final DerivativeStructure mu) {

  970.         // get the position in body frame
  971.         final Transform fromBodyFrame = bodyFrame.getTransformTo(frame, date);
  972.         final Transform toBodyFrame   = fromBodyFrame.getInverse();
  973.         final Vector3D positionBody   = toBodyFrame.transformPosition(position.toVector3D());

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

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

  978.         // Hessian of the non-central part of the gravity field
  979.         final RealMatrix hBody     = new Array2DRowRealMatrix(gh.getHessian(), false);
  980.         final RealMatrix rot       = new Array2DRowRealMatrix(toBodyFrame.getRotation().getMatrix());
  981.         final RealMatrix hInertial = rot.transpose().multiply(hBody).multiply(rot);

  982.         // distribute all partial derivatives in a compact acceleration vector
  983.         final double[] derivatives = new double[1 + position.getX().getFreeParameters()];
  984.         final DerivativeStructure[] accDer = new DerivativeStructure[3];
  985.         for (int i = 0; i < 3; ++i) {

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

  988.             // next three elements are one row of the Jacobian of acceleration (i.e. Hessian of field)
  989.             derivatives[1] = hInertial.getEntry(i, 0);
  990.             derivatives[2] = hInertial.getEntry(i, 1);
  991.             derivatives[3] = hInertial.getEntry(i, 2);

  992.             // next element is derivative with respect to parameter mu
  993.             if (derivatives.length > 4 && isVariable(mu, 3)) {
  994.                 derivatives[4] = gInertial[i] / mu.getReal();
  995.             }

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

  997.         }

  998.         return new FieldVector3D<>(accDer);

  999.     }

  1000.     /** {@inheritDoc} */
  1001.     public ParameterDriver[] getParametersDrivers() {
  1002.         return new ParameterDriver[] {
  1003.             gmParameterDriver
  1004.         };
  1005.     }

  1006. }