HolmesFeatherstoneAttractionModel.java

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


  18. import java.util.Collections;
  19. import java.util.List;
  20. import java.util.stream.Stream;

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

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

  73. public class HolmesFeatherstoneAttractionModel extends AbstractForceModel implements TideSystemProvider {

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

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

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

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

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

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

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

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

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

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

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

  109.         this.provider  = provider;
  110.         this.bodyFrame = centralBodyFrame;

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

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

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

  140.     }

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

  215.             }

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

  221.         }

  222.         // scale back
  223.         value = FastMath.scalb(value, SCALING);

  224.         // apply the global mu/r factor
  225.         return mu * value / r;

  226.     }

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

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

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

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

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

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

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

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

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

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

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

  302.             }

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

  308.         }

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

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

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

  322.     }

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

  448.     }

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

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

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

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

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

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

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

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

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

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

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

  545.             }

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

  554.         }

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


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

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


  579.     }

  580.     /** Container for gradient and Hessian. */
  581.     private static class GradientHessian {

  582.         /** Gradient. */
  583.         private final double[] gradient;

  584.         /** Hessian. */
  585.         private final double[][] hessian;

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

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

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

  609.     }

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

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

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

  625.         return aOrN;

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

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

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

  643.         return aOrN;

  644.     }

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

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

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

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

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

  672.         return cosSin;

  673.     }

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

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

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

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

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

  702.             }
  703.         }

  704.         return cosSin;

  705.     }

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

  731.         final double u2 = u * u;

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

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

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

  743.             ++localIndex;
  744.             ++n;

  745.         }

  746.         if (pnm1 != null) {

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

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

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

  758.                 ++localIndex;
  759.                 ++n;

  760.             }

  761.             if (pnm2 != null) {

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

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

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

  773.                     ++localIndex;
  774.                     ++n;

  775.                 }

  776.             }

  777.         }

  778.         return localIndex;

  779.     }

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

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

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

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

  821.         }
  822.         if (pnm1 != null) {

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

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

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

  834.                 ++localIndex;
  835.                 ++n;

  836.             }

  837.             if (pnm2 != null) {

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

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

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

  851.                 }

  852.             }

  853.         }
  854.         return localIndex;

  855.     }

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

  859.         final double mu = parameters[0];

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

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

  867.     }

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

  871.         final T mu = parameters[0];

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

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

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

  897.     }

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

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

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

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

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

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

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

  1010.         // free parameters
  1011.         final int freeParameters = mu.getFreeParameters();

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

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

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

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

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

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

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

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

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

  1039.         }

  1040.         return new FieldVector3D<>(accDer);

  1041.     }

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

  1072.         // free parameters
  1073.         final int freeParameters = mu.getFreeParameters();

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

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

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

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

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

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

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

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

  1099.         }

  1100.         return new FieldVector3D<>(accDer);

  1101.     }

  1102.     /** {@inheritDoc} */
  1103.     public List<ParameterDriver> getParametersDrivers() {
  1104.         return Collections.singletonList(gmParameterDriver);
  1105.     }

  1106. }