HolmesFeatherstoneAttractionModel.java

  1. /* Copyright 2002-2018 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.      * @exception OrekitException if position cannot be converted to central body frame
  166.      */
  167.     public double value(final AbsoluteDate date, final Vector3D position,
  168.                         final double mu)
  169.         throws OrekitException {
  170.         return mu / position.getNorm() + nonCentralPart(date, position, mu);
  171.     }

  172.     /** Compute the non-central part of the gravity field.
  173.      * @param date current date
  174.      * @param position position at which gravity field is desired in body frame
  175.      * @param mu central attraction coefficient to use
  176.      * @return value of the non-central part of the gravity field
  177.      * @exception OrekitException if position cannot be converted to central body frame
  178.      */
  179.     public double nonCentralPart(final AbsoluteDate date, final Vector3D position, final double mu)
  180.         throws OrekitException {

  181.         final int degree = provider.getMaxDegree();
  182.         final int order  = provider.getMaxOrder();
  183.         final NormalizedSphericalHarmonics harmonics = provider.onDate(date);

  184.         // allocate the columns for recursion
  185.         double[] pnm0Plus2 = new double[degree + 1];
  186.         double[] pnm0Plus1 = new double[degree + 1];
  187.         double[] pnm0      = new double[degree + 1];

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

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

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

  205.         // outer summation over order
  206.         int    index = 0;
  207.         double value = 0;
  208.         for (int m = degree; m >= 0; --m) {

  209.             // compute tesseral terms without derivatives
  210.             index = computeTesseral(m, degree, index, t, u, tOu,
  211.                                     pnm0Plus2, pnm0Plus1, null, pnm0, null, null);

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

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

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

  223.             }

  224.             // rotate the recursion arrays
  225.             final double[] tmp = pnm0Plus2;
  226.             pnm0Plus2 = pnm0Plus1;
  227.             pnm0Plus1 = pnm0;
  228.             pnm0      = tmp;

  229.         }

  230.         // scale back
  231.         value = FastMath.scalb(value, SCALING);

  232.         // apply the global mu/r factor
  233.         return mu * value / r;

  234.     }

  235.     /** Compute the gradient of the non-central part of the gravity field.
  236.      * @param date current date
  237.      * @param position position at which gravity field is desired in body frame
  238.      * @param mu central attraction coefficient to use
  239.      * @return gradient of the non-central part of the gravity field
  240.      * @exception OrekitException if position cannot be converted to central body frame
  241.      */
  242.     public double[] gradient(final AbsoluteDate date, final Vector3D position, final double mu)
  243.         throws OrekitException {

  244.         final int degree = provider.getMaxDegree();
  245.         final int order  = provider.getMaxOrder();
  246.         final NormalizedSphericalHarmonics harmonics = provider.onDate(date);

  247.         // allocate the columns for recursion
  248.         double[] pnm0Plus2  = new double[degree + 1];
  249.         double[] pnm0Plus1  = new double[degree + 1];
  250.         double[] pnm0       = new double[degree + 1];
  251.         final double[] pnm1 = new double[degree + 1];

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

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

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

  270.         // outer summation over order
  271.         int    index = 0;
  272.         double value = 0;
  273.         final double[] gradient = new double[3];
  274.         for (int m = degree; m >= 0; --m) {

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

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

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

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

  312.             }

  313.             // rotate the recursion arrays
  314.             final double[] tmp = pnm0Plus2;
  315.             pnm0Plus2 = pnm0Plus1;
  316.             pnm0Plus1 = pnm0;
  317.             pnm0      = tmp;

  318.         }

  319.         // scale back
  320.         value       = FastMath.scalb(value,       SCALING);
  321.         gradient[0] = FastMath.scalb(gradient[0], SCALING);
  322.         gradient[1] = FastMath.scalb(gradient[1], SCALING);
  323.         gradient[2] = FastMath.scalb(gradient[2], SCALING);

  324.         // apply the global mu/r factor
  325.         final double muOr = mu / r;
  326.         value            *= muOr;
  327.         gradient[0]       = muOr * gradient[0] - value / r;
  328.         gradient[1]      *= muOr;
  329.         gradient[2]      *= muOr;

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

  332.     }

  333.     /** Compute the gradient of the non-central part of the gravity field.
  334.      * @param date current date
  335.      * @param position position at which gravity field is desired in body frame
  336.      * @param mu central attraction coefficient to use
  337.      * @param <T> type of field used
  338.      * @return gradient of the non-central part of the gravity field
  339.      * @exception OrekitException if position cannot be converted to central body frame
  340.      */
  341.     public <T extends RealFieldElement<T>> T[] gradient(final FieldAbsoluteDate<T> date, final FieldVector3D<T> position,
  342.                                                         final T mu)
  343.         throws OrekitException {

  344.         final int degree = provider.getMaxDegree();
  345.         final int order  = provider.getMaxOrder();
  346.         final NormalizedSphericalHarmonics harmonics = provider.onDate(date.toAbsoluteDate());
  347.         final T zero = date.getField().getZero();
  348.         // allocate the columns for recursion
  349.         T[] pnm0Plus2  = MathArrays.buildArray(date.getField(), degree + 1);
  350.         T[] pnm0Plus1  = MathArrays.buildArray(date.getField(), degree + 1);
  351.         T[] pnm0       = MathArrays.buildArray(date.getField(), degree + 1);
  352.         final T[] pnm1 = MathArrays.buildArray(date.getField(), degree + 1);

  353.         // compute polar coordinates
  354.         final T x    = position.getX();
  355.         final T y    = position.getY();
  356.         final T z    = position.getZ();
  357.         final T x2   = x.multiply(x);
  358.         final T y2   = y.multiply(y);
  359.         final T z2   = z.multiply(z);
  360.         final T r2   = x2.add(y2).add(z2);
  361.         final T r    = r2.sqrt();
  362.         final T rho2 = x2.add(y2);
  363.         final T rho  = rho2.sqrt();
  364.         final T t    = z.divide(r);   // cos(theta), where theta is the polar angle
  365.         final T u    = rho.divide(r); // sin(theta), where theta is the polar angle
  366.         final T tOu  = z.divide(rho);

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

  369.         // compute longitude cosines/sines
  370.         final T[][] cosSinLambda = createCosSinArrays(rho.reciprocal().multiply(position.getX()), rho.reciprocal().multiply(position.getY()));
  371.         // outer summation over order
  372.         int    index = 0;
  373.         T value = zero;
  374.         final T[] gradient = MathArrays.buildArray(zero.getField(), 3);
  375.         for (int m = degree; m >= 0; --m) {

  376.             // compute tesseral terms with derivatives
  377.             index = computeTesseral(m, degree, index, t, u, tOu,
  378.                                     pnm0Plus2, pnm0Plus1, null, pnm0, pnm1, null);
  379.             if (m <= order) {
  380.                 // compute contribution of current order to field (equation 5 of the paper)

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

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

  419.         }
  420.         // scale back
  421.         value       = value.scalb(SCALING);
  422.         gradient[0] = gradient[0].scalb(SCALING);
  423.         gradient[1] = gradient[1].scalb(SCALING);
  424.         gradient[2] = gradient[2].scalb(SCALING);

  425.         // apply the global mu/r factor
  426.         final T muOr = r.reciprocal().multiply(mu);
  427.         value            = value.multiply(muOr);
  428.         gradient[0]       = muOr.multiply(gradient[0]).subtract(value.divide(r));
  429.         gradient[1]      = gradient[1].multiply(muOr);
  430.         gradient[2]      = gradient[2].multiply(muOr);

  431.         // convert gradient from spherical to Cartesian
  432.         // Cartesian coordinates
  433.         // remaining spherical coordinates
  434.         final T rPos     = position.getNorm();
  435.         // intermediate variables
  436.         final T xPos    = position.getX();
  437.         final T yPos    = position.getY();
  438.         final T zPos    = position.getZ();
  439.         final T rho2Pos = x.multiply(x).add(y.multiply(y));
  440.         final T rhoPos  = rho2.sqrt();
  441.         final T r2Pos   = rho2.add(z.multiply(z));

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

  443.         // row representing the gradient of r
  444.         jacobianPos[0][0] = xPos.divide(rPos);
  445.         jacobianPos[0][1] = yPos.divide(rPos);
  446.         jacobianPos[0][2] = zPos.divide(rPos);

  447.         // row representing the gradient of theta
  448.         jacobianPos[1][0] =  yPos.negate().divide(rho2Pos);
  449.         jacobianPos[1][1] =  xPos.divide(rho2Pos);
  450.         // jacobian[1][2] is already set to 0 at allocation time

  451.         // row representing the gradient of phi
  452.         jacobianPos[2][0] = xPos.multiply(zPos).divide(rhoPos.multiply(r2Pos));
  453.         jacobianPos[2][1] = yPos.multiply(zPos).divide(rhoPos.multiply(r2Pos));
  454.         jacobianPos[2][2] = rhoPos.negate().divide(r2Pos);
  455.         final T[] cartGradPos = MathArrays.buildArray(zero.getField(), 3);
  456.         cartGradPos[0] = gradient[0].multiply(jacobianPos[0][0]).add(gradient[1].multiply(jacobianPos[1][0])).add(gradient[2].multiply(jacobianPos[2][0]));
  457.         cartGradPos[1] = gradient[0].multiply(jacobianPos[0][1]).add(gradient[1].multiply(jacobianPos[1][1])).add(gradient[2].multiply(jacobianPos[2][1]));
  458.         cartGradPos[2] = gradient[0].multiply(jacobianPos[0][2])                                      .add(gradient[2].multiply(jacobianPos[2][2]));
  459.         return cartGradPos;

  460.     }

  461.     /** Compute both the gradient and the hessian of the non-central part of the gravity field.
  462.      * @param date current date
  463.      * @param position position at which gravity field is desired in body frame
  464.      * @param mu central attraction coefficient to use
  465.      * @return gradient and hessian of the non-central part of the gravity field
  466.      * @exception OrekitException if position cannot be converted to central body frame
  467.      */
  468.     private GradientHessian gradientHessian(final AbsoluteDate date, final Vector3D position, final double mu)
  469.         throws OrekitException {

  470.         final int degree = provider.getMaxDegree();
  471.         final int order  = provider.getMaxOrder();
  472.         final NormalizedSphericalHarmonics harmonics = provider.onDate(date);

  473.         // allocate the columns for recursion
  474.         double[] pnm0Plus2  = new double[degree + 1];
  475.         double[] pnm0Plus1  = new double[degree + 1];
  476.         double[] pnm0       = new double[degree + 1];
  477.         double[] pnm1Plus1  = new double[degree + 1];
  478.         double[] pnm1       = new double[degree + 1];
  479.         final double[] pnm2 = new double[degree + 1];

  480.         // compute polar coordinates
  481.         final double x    = position.getX();
  482.         final double y    = position.getY();
  483.         final double z    = position.getZ();
  484.         final double x2   = x * x;
  485.         final double y2   = y * y;
  486.         final double z2   = z * z;
  487.         final double r2   = x2 + y2 + z2;
  488.         final double r    = FastMath.sqrt (r2);
  489.         final double rho2 = x2 + y2;
  490.         final double rho  = FastMath.sqrt(rho2);
  491.         final double t    = z / r;   // cos(theta), where theta is the polar angle
  492.         final double u    = rho / r; // sin(theta), where theta is the polar angle
  493.         final double tOu  = z / rho;

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

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

  498.         // outer summation over order
  499.         int    index = 0;
  500.         double value = 0;
  501.         final double[]   gradient = new double[3];
  502.         final double[][] hessian  = new double[3][3];
  503.         for (int m = degree; m >= 0; --m) {

  504.             // compute tesseral terms
  505.             index = computeTesseral(m, degree, index, t, u, tOu,
  506.                                     pnm0Plus2, pnm0Plus1, pnm1Plus1, pnm0, pnm1, pnm2);

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

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

  546.                 // contribution to outer summation over order
  547.                 final double sML = cosSinLambda[1][m];
  548.                 final double cML = cosSinLambda[0][m];
  549.                 value            = value         * u + sML * sumDegreeS + cML * sumDegreeC;
  550.                 gradient[0]      = gradient[0]   * u + sML * dSumDegreeSdR + cML * dSumDegreeCdR;
  551.                 gradient[1]      = gradient[1]   * u + m * (cML * sumDegreeS - sML * sumDegreeC);
  552.                 gradient[2]      = gradient[2]   * u + sML * dSumDegreeSdTheta + cML * dSumDegreeCdTheta;
  553.                 hessian[0][0]    = hessian[0][0] * u + sML * d2SumDegreeSdRdR + cML * d2SumDegreeCdRdR;
  554.                 hessian[1][0]    = hessian[1][0] * u + m * (cML * dSumDegreeSdR - sML * dSumDegreeCdR);
  555.                 hessian[2][0]    = hessian[2][0] * u + sML * d2SumDegreeSdRdTheta + cML * d2SumDegreeCdRdTheta;
  556.                 hessian[1][1]    = hessian[1][1] * u - m * m * (sML * sumDegreeS + cML * sumDegreeC);
  557.                 hessian[2][1]    = hessian[2][1] * u + m * (cML * dSumDegreeSdTheta - sML * dSumDegreeCdTheta);
  558.                 hessian[2][2]    = hessian[2][2] * u + sML * d2SumDegreeSdThetadTheta + cML * d2SumDegreeCdThetadTheta;

  559.             }

  560.             // rotate the recursion arrays
  561.             final double[] tmp0 = pnm0Plus2;
  562.             pnm0Plus2 = pnm0Plus1;
  563.             pnm0Plus1 = pnm0;
  564.             pnm0      = tmp0;
  565.             final double[] tmp1 = pnm1Plus1;
  566.             pnm1Plus1 = pnm1;
  567.             pnm1      = tmp1;

  568.         }

  569.         // scale back
  570.         value = FastMath.scalb(value, SCALING);
  571.         for (int i = 0; i < 3; ++i) {
  572.             gradient[i] = FastMath.scalb(gradient[i], SCALING);
  573.             for (int j = 0; j <= i; ++j) {
  574.                 hessian[i][j] = FastMath.scalb(hessian[i][j], SCALING);
  575.             }
  576.         }


  577.         // apply the global mu/r factor
  578.         final double muOr = mu / r;
  579.         value         *= muOr;
  580.         gradient[0]    = muOr * gradient[0] - value / r;
  581.         gradient[1]   *= muOr;
  582.         gradient[2]   *= muOr;
  583.         hessian[0][0]  = muOr * hessian[0][0] - 2 * gradient[0] / r;
  584.         hessian[1][0]  = muOr * hessian[1][0] -     gradient[1] / r;
  585.         hessian[2][0]  = muOr * hessian[2][0] -     gradient[2] / r;
  586.         hessian[1][1] *= muOr;
  587.         hessian[2][1] *= muOr;
  588.         hessian[2][2] *= muOr;

  589.         // convert gradient and Hessian from spherical to Cartesian
  590.         final SphericalCoordinates sc = new SphericalCoordinates(position);
  591.         return new GradientHessian(sc.toCartesianGradient(gradient),
  592.                                    sc.toCartesianHessian(hessian, gradient));


  593.     }

  594.     /** Container for gradient and Hessian. */
  595.     private static class GradientHessian {

  596.         /** Gradient. */
  597.         private final double[] gradient;

  598.         /** Hessian. */
  599.         private final double[][] hessian;

  600.         /** Simple constructor.
  601.          * <p>
  602.          * A reference to the arrays is stored, they are <strong>not</strong> cloned.
  603.          * </p>
  604.          * @param gradient gradient
  605.          * @param hessian hessian
  606.          */
  607.         GradientHessian(final double[] gradient, final double[][] hessian) {
  608.             this.gradient = gradient;
  609.             this.hessian  = hessian;
  610.         }

  611.         /** Get a reference to the gradient.
  612.          * @return gradient (a reference to the internal array is returned)
  613.          */
  614.         public double[] getGradient() {
  615.             return gradient;
  616.         }

  617.         /** Get a reference to the Hessian.
  618.          * @return Hessian (a reference to the internal array is returned)
  619.          */
  620.         public double[][] getHessian() {
  621.             return hessian;
  622.         }

  623.     }

  624.     /** Compute a/r powers array.
  625.      * @param aOr a/r
  626.      * @return array containing (a/r)<sup>n</sup>
  627.      */
  628.     private double[] createDistancePowersArray(final double aOr) {

  629.         // initialize array
  630.         final double[] aOrN = new double[provider.getMaxDegree() + 1];
  631.         aOrN[0] = 1;
  632.         aOrN[1] = aOr;

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

  639.         return aOrN;

  640.     }
  641.     /** Compute a/r powers array.
  642.      * @param aOr a/r
  643.      * @param <T> type of field used
  644.      * @return array containing (a/r)<sup>n</sup>
  645.      */
  646.     private <T extends RealFieldElement<T>> T[] createDistancePowersArray(final T aOr) {

  647.         // initialize array
  648.         final T[] aOrN = MathArrays.buildArray(aOr.getField(), provider.getMaxDegree() + 1);
  649.         aOrN[0] = aOr.getField().getOne();
  650.         aOrN[1] = aOr;

  651.         // fill up array
  652.         for (int n = 2; n < aOrN.length; ++n) {
  653.             final int p = n / 2;
  654.             final int q = n - p;
  655.             aOrN[n] = aOrN[p].multiply(aOrN[q]);
  656.         }

  657.         return aOrN;

  658.     }

  659.     /** Compute longitude cosines and sines.
  660.      * @param cosLambda cos(λ)
  661.      * @param sinLambda sin(λ)
  662.      * @return array containing cos(m &times; λ) in row 0
  663.      * and sin(m &times; λ) in row 1
  664.      */
  665.     private double[][] createCosSinArrays(final double cosLambda, final double sinLambda) {

  666.         // initialize arrays
  667.         final double[][] cosSin = new double[2][provider.getMaxOrder() + 1];
  668.         cosSin[0][0] = 1;
  669.         cosSin[1][0] = 0;
  670.         if (provider.getMaxOrder() > 0) {
  671.             cosSin[0][1] = cosLambda;
  672.             cosSin[1][1] = sinLambda;

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

  675.                 // m * lambda is split as p * lambda + q * lambda, trying to avoid
  676.                 // p or q being much larger than the other. This reduces the number of
  677.                 // intermediate results reused to compute each value, and hence should limit
  678.                 // as much as possible roundoff error accumulation
  679.                 // (this does not change the number of floating point operations)
  680.                 final int p = m / 2;
  681.                 final int q = m - p;

  682.                 cosSin[0][m] = cosSin[0][p] * cosSin[0][q] - cosSin[1][p] * cosSin[1][q];
  683.                 cosSin[1][m] = cosSin[1][p] * cosSin[0][q] + cosSin[0][p] * cosSin[1][q];
  684.             }
  685.         }

  686.         return cosSin;

  687.     }

  688.     /** Compute longitude cosines and sines.
  689.      * @param cosLambda cos(λ)
  690.      * @param sinLambda sin(λ)
  691.      * @param <T> type of field used
  692.      * @return array containing cos(m &times; λ) in row 0
  693.      * and sin(m &times; λ) in row 1
  694.      */
  695.     private <T extends RealFieldElement<T>> T[][] createCosSinArrays(final T cosLambda, final T sinLambda) {

  696.         final T one = cosLambda.getField().getOne();
  697.         final T zero = cosLambda.getField().getZero();
  698.         // initialize arrays
  699.         final T[][] cosSin = MathArrays.buildArray(one.getField(), 2, provider.getMaxOrder() + 1);
  700.         cosSin[0][0] = one;
  701.         cosSin[1][0] = zero;
  702.         if (provider.getMaxOrder() > 0) {
  703.             cosSin[0][1] = cosLambda;
  704.             cosSin[1][1] = sinLambda;

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

  707.                 // m * lambda is split as p * lambda + q * lambda, trying to avoid
  708.                 // p or q being much larger than the other. This reduces the number of
  709.                 // intermediate results reused to compute each value, and hence should limit
  710.                 // as much as possible roundoff error accumulation
  711.                 // (this does not change the number of floating point operations)
  712.                 final int p = m / 2;
  713.                 final int q = m - p;

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

  716.             }
  717.         }

  718.         return cosSin;

  719.     }

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

  745.         final double u2 = u * u;

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

  752.         // compute tesseral values
  753.         int localIndex = index;
  754.         while (n <= degree) {

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

  757.             ++localIndex;
  758.             ++n;

  759.         }

  760.         if (pnm1 != null) {

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

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

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

  772.                 ++localIndex;
  773.                 ++n;

  774.             }

  775.             if (pnm2 != null) {

  776.                 // initialize recursion from sectorial terms
  777.                 n = FastMath.max(2, m);
  778.                 if (n == m) {
  779.                     pnm2[n] = m * (tOu * pnm1[n] - pnm0[n] / u2);
  780.                     ++n;
  781.                 }

  782.                 // compute tesseral values and derivatives with respect to polar angle
  783.                 localIndex = index;
  784.                 while (n <= degree) {

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

  787.                     ++localIndex;
  788.                     ++n;

  789.                 }

  790.             }

  791.         }

  792.         return localIndex;

  793.     }

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

  820.         final T u2 = u.multiply(u);
  821.         final T zero = u.getField().getZero();
  822.         // initialize recursion from sectorial terms
  823.         int n = FastMath.max(2, m);
  824.         if (n == m) {
  825.             pnm0[n] = zero.add(sectorial[n]);
  826.             ++n;
  827.         }

  828.         // compute tesseral values
  829.         int localIndex = index;
  830.         while (n <= degree) {

  831.             // value (equation 27 of the paper)
  832.             pnm0[n] = t.multiply(gnmOj[localIndex]).multiply(pnm0Plus1[n]).subtract(u2.multiply(pnm0Plus2[n]).multiply(hnmOj[localIndex]));
  833.             ++localIndex;
  834.             ++n;

  835.         }
  836.         if (pnm1 != null) {

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

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

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

  848.                 ++localIndex;
  849.                 ++n;

  850.             }

  851.             if (pnm2 != null) {

  852.                 // initialize recursion from sectorial terms
  853.                 n = FastMath.max(2, m);
  854.                 if (n == m) {
  855.                     pnm2[n] =   tOu.multiply(pnm1[n]).subtract(pnm0[n].divide(u2)).multiply(m);
  856.                     ++n;
  857.                 }

  858.                 // compute tesseral values and derivatives with respect to polar angle
  859.                 localIndex = index;
  860.                 while (n <= degree) {

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

  865.                 }

  866.             }

  867.         }
  868.         return localIndex;

  869.     }

  870.     /** {@inheritDoc} */
  871.     @Override
  872.     public Vector3D acceleration(final SpacecraftState s, final double[] parameters)
  873.         throws OrekitException {

  874.         final double mu = parameters[0];

  875.         // get the position in body frame
  876.         final AbsoluteDate date       = s.getDate();
  877.         final Transform fromBodyFrame = bodyFrame.getTransformTo(s.getFrame(), date);
  878.         final Transform toBodyFrame   = fromBodyFrame.getInverse();
  879.         final Vector3D position       = toBodyFrame.transformPosition(s.getPVCoordinates().getPosition());

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

  882.     }

  883.     /** {@inheritDoc} */
  884.     public <T extends RealFieldElement<T>> FieldVector3D<T> acceleration(final FieldSpacecraftState<T> s,
  885.                                                                          final T[] parameters)
  886.         throws OrekitException {

  887.         final T mu = parameters[0];

  888.         // check for faster computation dedicated to derivatives with respect to state
  889.         if (isStateDerivative(s)) {
  890.             @SuppressWarnings("unchecked")
  891.             final FieldVector3D<DerivativeStructure> p = (FieldVector3D<DerivativeStructure>) s.getPVCoordinates().getPosition();
  892.             @SuppressWarnings("unchecked")
  893.             final FieldVector3D<T> a = (FieldVector3D<T>) accelerationWrtState(s.getDate().toAbsoluteDate(),
  894.                                                                                s.getFrame(), p,
  895.                                                                                (DerivativeStructure) mu);
  896.             return a;
  897.         }

  898.         // get the position in body frame
  899.         final FieldAbsoluteDate<T> date          = s.getDate();
  900.         final Transform            fromBodyFrame = bodyFrame.getTransformTo(s.getFrame(), date.toAbsoluteDate());
  901.         final Transform            toBodyFrame   = fromBodyFrame.getInverse();
  902.         final FieldVector3D<T>     position      = toBodyFrame.transformPosition(s.getPVCoordinates().getPosition());

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

  905.     }

  906.     /** {@inheritDoc} */
  907.     public Stream<EventDetector> getEventsDetectors() {
  908.         return Stream.empty();
  909.     }

  910.     @Override
  911.     /** {@inheritDoc} */
  912.     public <T extends RealFieldElement<T>> Stream<FieldEventDetector<T>> getFieldEventsDetectors(final Field<T> field) {
  913.         return Stream.empty();
  914.     }

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

  938.     /** Check if a derivative represents a specified variable.
  939.      * @param ds derivative to check
  940.      * @param index index of the variable
  941.      * @return true if the derivative represents a specified variable
  942.      * @since 9.0
  943.      */
  944.     private boolean isVariable(final DerivativeStructure ds, final int index) {
  945.         final double[] derivatives = ds.getAllDerivatives();
  946.         boolean check = true;
  947.         for (int i = 1; i < derivatives.length; ++i) {
  948.             check &= derivatives[i] == ((index + 1 == i) ? 1.0 : 0.0);
  949.         }
  950.         return check;
  951.     }

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

  984.         // get the position in body frame
  985.         final Transform fromBodyFrame = bodyFrame.getTransformTo(frame, date);
  986.         final Transform toBodyFrame   = fromBodyFrame.getInverse();
  987.         final Vector3D positionBody   = toBodyFrame.transformPosition(position.toVector3D());

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

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

  992.         // Hessian of the non-central part of the gravity field
  993.         final RealMatrix hBody     = new Array2DRowRealMatrix(gh.getHessian(), false);
  994.         final RealMatrix rot       = new Array2DRowRealMatrix(toBodyFrame.getRotation().getMatrix());
  995.         final RealMatrix hInertial = rot.transpose().multiply(hBody).multiply(rot);

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

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

  1002.             // next three elements are one row of the Jacobian of acceleration (i.e. Hessian of field)
  1003.             derivatives[1] = hInertial.getEntry(i, 0);
  1004.             derivatives[2] = hInertial.getEntry(i, 1);
  1005.             derivatives[3] = hInertial.getEntry(i, 2);

  1006.             // next element is derivative with respect to parameter mu
  1007.             if (derivatives.length > 4 && isVariable(mu, 3)) {
  1008.                 derivatives[4] = gInertial[i] / mu.getReal();
  1009.             }

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

  1011.         }

  1012.         return new FieldVector3D<>(accDer);

  1013.     }

  1014.     /** {@inheritDoc} */
  1015.     public ParameterDriver[] getParametersDrivers() {
  1016.         return new ParameterDriver[] {
  1017.             gmParameterDriver
  1018.         };
  1019.     }

  1020. }