AbstractGaussianContribution.java

  1. /* Copyright 2002-2021 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.propagation.semianalytical.dsst.forces;

  18. import java.lang.reflect.Array;
  19. import java.util.ArrayList;
  20. import java.util.Collections;
  21. import java.util.HashMap;
  22. import java.util.List;
  23. import java.util.Map;
  24. import java.util.Set;

  25. import org.hipparchus.Field;
  26. import org.hipparchus.CalculusFieldElement;
  27. import org.hipparchus.analysis.CalculusFieldUnivariateVectorFunction;
  28. import org.hipparchus.analysis.UnivariateVectorFunction;
  29. import org.hipparchus.geometry.euclidean.threed.FieldVector3D;
  30. import org.hipparchus.geometry.euclidean.threed.Vector3D;
  31. import org.hipparchus.util.FastMath;
  32. import org.hipparchus.util.FieldSinCos;
  33. import org.hipparchus.util.MathArrays;
  34. import org.hipparchus.util.SinCos;
  35. import org.orekit.attitudes.Attitude;
  36. import org.orekit.attitudes.AttitudeProvider;
  37. import org.orekit.attitudes.FieldAttitude;
  38. import org.orekit.forces.ForceModel;
  39. import org.orekit.orbits.EquinoctialOrbit;
  40. import org.orekit.orbits.FieldEquinoctialOrbit;
  41. import org.orekit.orbits.FieldOrbit;
  42. import org.orekit.orbits.Orbit;
  43. import org.orekit.orbits.OrbitType;
  44. import org.orekit.orbits.PositionAngle;
  45. import org.orekit.propagation.FieldSpacecraftState;
  46. import org.orekit.propagation.PropagationType;
  47. import org.orekit.propagation.SpacecraftState;
  48. import org.orekit.propagation.semianalytical.dsst.utilities.AuxiliaryElements;
  49. import org.orekit.propagation.semianalytical.dsst.utilities.CjSjCoefficient;
  50. import org.orekit.propagation.semianalytical.dsst.utilities.FieldAuxiliaryElements;
  51. import org.orekit.propagation.semianalytical.dsst.utilities.FieldCjSjCoefficient;
  52. import org.orekit.propagation.semianalytical.dsst.utilities.FieldShortPeriodicsInterpolatedCoefficient;
  53. import org.orekit.propagation.semianalytical.dsst.utilities.ShortPeriodicsInterpolatedCoefficient;
  54. import org.orekit.time.AbsoluteDate;
  55. import org.orekit.time.FieldAbsoluteDate;
  56. import org.orekit.utils.FieldTimeSpanMap;
  57. import org.orekit.utils.ParameterDriver;
  58. import org.orekit.utils.TimeSpanMap;

  59. /**
  60.  * Common handling of {@link DSSTForceModel} methods for Gaussian contributions
  61.  * to DSST propagation.
  62.  * <p>
  63.  * This abstract class allows to provide easily a subset of
  64.  * {@link DSSTForceModel} methods for specific Gaussian contributions.
  65.  * </p>
  66.  * <p>
  67.  * This class implements the notion of numerical averaging of the DSST theory.
  68.  * Numerical averaging is mainly used for non-conservative disturbing forces
  69.  * such as atmospheric drag and solar radiation pressure.
  70.  * </p>
  71.  * <p>
  72.  * Gaussian contributions can be expressed as: da<sub>i</sub>/dt =
  73.  * δa<sub>i</sub>/δv . q<br>
  74.  * where:
  75.  * <ul>
  76.  * <li>a<sub>i</sub> are the six equinoctial elements</li>
  77.  * <li>v is the velocity vector</li>
  78.  * <li>q is the perturbing acceleration due to the considered force</li>
  79.  * </ul>
  80.  *
  81.  * <p>
  82.  * The averaging process and other considerations lead to integrate this
  83.  * contribution over the true longitude L possibly taking into account some
  84.  * limits.
  85.  *
  86.  * <p>
  87.  * To create a numerically averaged contribution, one needs only to provide a
  88.  * {@link ForceModel} and to implement in the derived class the methods:
  89.  * {@link #getLLimits(SpacecraftState, AuxiliaryElements)} and
  90.  * {@link #getParametersDriversWithoutMu()}.
  91.  * </p>
  92.  * @author Pascal Parraud
  93.  * @author Bryan Cazabonne (field translation)
  94.  */
  95. public abstract class AbstractGaussianContribution implements DSSTForceModel {

  96.     /**
  97.      * Retrograde factor I.
  98.      * <p>
  99.      * DSST model needs equinoctial orbit as internal representation. Classical
  100.      * equinoctial elements have discontinuities when inclination is close to zero.
  101.      * In this representation, I = +1. <br>
  102.      * To avoid this discontinuity, another representation exists and equinoctial
  103.      * elements can be expressed in a different way, called "retrograde" orbit. This
  104.      * implies I = -1. <br>
  105.      * As Orekit doesn't implement the retrograde orbit, I is always set to +1. But
  106.      * for the sake of consistency with the theory, the retrograde factor has been
  107.      * kept in the formulas.
  108.      * </p>
  109.      */
  110.     private static final int I = 1;

  111.     /**
  112.      * Central attraction scaling factor.
  113.      * <p>
  114.      * We use a power of 2 to avoid numeric noise introduction in the
  115.      * multiplications/divisions sequences.
  116.      * </p>
  117.      */
  118.     private static final double MU_SCALE = FastMath.scalb(1.0, 32);

  119.     /** Available orders for Gauss quadrature. */
  120.     private static final int[] GAUSS_ORDER = { 12, 16, 20, 24, 32, 40, 48 };

  121.     /** Max rank in Gauss quadrature orders array. */
  122.     private static final int MAX_ORDER_RANK = GAUSS_ORDER.length - 1;

  123.     /** Number of points for interpolation. */
  124.     private static final int INTERPOLATION_POINTS = 3;

  125.     /** Maximum value for j index. */
  126.     private static final int JMAX = 12;

  127.     /** Contribution to be numerically averaged. */
  128.     private final ForceModel contribution;

  129.     /** Gauss integrator. */
  130.     private final double threshold;

  131.     /** Gauss integrator. */
  132.     private GaussQuadrature integrator;

  133.     /** Flag for Gauss order computation. */
  134.     private boolean isDirty;

  135.     /** Attitude provider. */
  136.     private AttitudeProvider attitudeProvider;

  137.     /** Prefix for coefficients keys. */
  138.     private final String coefficientsKeyPrefix;

  139.     /** Short period terms. */
  140.     private GaussianShortPeriodicCoefficients gaussianSPCoefs;

  141.     /** Short period terms. */
  142.     private Map<Field<?>, FieldGaussianShortPeriodicCoefficients<?>> gaussianFieldSPCoefs;

  143.     /** Driver for gravitational parameter. */
  144.     private final ParameterDriver gmParameterDriver;

  145.     /**
  146.      * Build a new instance.
  147.      * @param coefficientsKeyPrefix prefix for coefficients keys
  148.      * @param threshold             tolerance for the choice of the Gauss quadrature
  149.      *                              order
  150.      * @param contribution          the {@link ForceModel} to be numerically
  151.      *                              averaged
  152.      * @param mu                    central attraction coefficient
  153.      */
  154.     protected AbstractGaussianContribution(final String coefficientsKeyPrefix, final double threshold,
  155.             final ForceModel contribution, final double mu) {

  156.         gmParameterDriver = new ParameterDriver(DSSTNewtonianAttraction.CENTRAL_ATTRACTION_COEFFICIENT, mu, MU_SCALE,
  157.                 0.0, Double.POSITIVE_INFINITY);

  158.         this.coefficientsKeyPrefix = coefficientsKeyPrefix;
  159.         this.contribution = contribution;
  160.         this.threshold = threshold;
  161.         this.integrator = new GaussQuadrature(GAUSS_ORDER[MAX_ORDER_RANK]);
  162.         this.isDirty = true;

  163.         gaussianFieldSPCoefs = new HashMap<>();
  164.     }

  165.     /** {@inheritDoc} */
  166.     @Override
  167.     public void init(final SpacecraftState initialState, final AbsoluteDate target) {
  168.         // Initialize the numerical force model
  169.         contribution.init(initialState, target);
  170.     }

  171.     /** {@inheritDoc} */
  172.     @Override
  173.     public List<ParameterDriver> getParametersDrivers() {

  174.         // Parameter drivers
  175.         final List<ParameterDriver> drivers = new ArrayList<>();

  176.         // Loop on drivers (without central attraction coefficient driver)
  177.         for (final ParameterDriver driver : getParametersDriversWithoutMu()) {
  178.             drivers.add(driver);
  179.         }

  180.         // We put central attraction coefficient driver at the end of the array
  181.         drivers.add(gmParameterDriver);
  182.         return drivers;

  183.     }

  184.     /**
  185.      * Get the drivers for force model parameters except the one for the central
  186.      * attraction coefficient.
  187.      * <p>
  188.      * The driver for central attraction coefficient is automatically added at the
  189.      * last element of the {@link ParameterDriver} array into
  190.      * {@link #getParametersDrivers()} method.
  191.      * </p>
  192.      * @return drivers for force model parameters
  193.      */
  194.     protected abstract List<ParameterDriver> getParametersDriversWithoutMu();

  195.     /** {@inheritDoc} */
  196.     @Override
  197.     public List<ShortPeriodTerms> initializeShortPeriodTerms(final AuxiliaryElements auxiliaryElements, final PropagationType type,
  198.             final double[] parameters) {

  199.         final List<ShortPeriodTerms> list = new ArrayList<ShortPeriodTerms>();
  200.         gaussianSPCoefs = new GaussianShortPeriodicCoefficients(coefficientsKeyPrefix, JMAX, INTERPOLATION_POINTS,
  201.                 new TimeSpanMap<Slot>(new Slot(JMAX, INTERPOLATION_POINTS)));
  202.         list.add(gaussianSPCoefs);
  203.         return list;

  204.     }

  205.     /** {@inheritDoc} */
  206.     @Override
  207.     public <T extends CalculusFieldElement<T>> List<FieldShortPeriodTerms<T>> initializeShortPeriodTerms(
  208.             final FieldAuxiliaryElements<T> auxiliaryElements, final PropagationType type, final T[] parameters) {

  209.         final Field<T> field = auxiliaryElements.getDate().getField();

  210.         final FieldGaussianShortPeriodicCoefficients<T> fgspc = new FieldGaussianShortPeriodicCoefficients<>(
  211.                 coefficientsKeyPrefix, JMAX, INTERPOLATION_POINTS,
  212.                 new FieldTimeSpanMap<>(new FieldSlot<>(JMAX, INTERPOLATION_POINTS), field));
  213.         gaussianFieldSPCoefs.put(field, fgspc);
  214.         return Collections.singletonList(fgspc);
  215.     }

  216.     /**
  217.      * Performs initialization at each integration step for the current force model.
  218.      * <p>
  219.      * This method aims at being called before mean elements rates computation.
  220.      * </p>
  221.      * @param auxiliaryElements auxiliary elements related to the current orbit
  222.      * @param parameters        parameters values of the force model parameters
  223.      * @return new force model context
  224.      */
  225.     private AbstractGaussianContributionContext initializeStep(final AuxiliaryElements auxiliaryElements,
  226.             final double[] parameters) {
  227.         return new AbstractGaussianContributionContext(auxiliaryElements, parameters);
  228.     }

  229.     /**
  230.      * Performs initialization at each integration step for the current force model.
  231.      * <p>
  232.      * This method aims at being called before mean elements rates computation.
  233.      * </p>
  234.      * @param <T>               type of the elements
  235.      * @param auxiliaryElements auxiliary elements related to the current orbit
  236.      * @param parameters        parameters values of the force model parameters
  237.      * @return new force model context
  238.      */
  239.     private <T extends CalculusFieldElement<T>> FieldAbstractGaussianContributionContext<T> initializeStep(
  240.             final FieldAuxiliaryElements<T> auxiliaryElements, final T[] parameters) {
  241.         return new FieldAbstractGaussianContributionContext<>(auxiliaryElements, parameters);
  242.     }

  243.     /** {@inheritDoc} */
  244.     @Override
  245.     public double[] getMeanElementRate(final SpacecraftState state, final AuxiliaryElements auxiliaryElements,
  246.             final double[] parameters) {

  247.         // Container for attributes
  248.         final AbstractGaussianContributionContext context = initializeStep(auxiliaryElements, parameters);
  249.         double[] meanElementRate = new double[6];
  250.         // Computes the limits for the integral
  251.         final double[] ll = getLLimits(state, auxiliaryElements);
  252.         // Computes integrated mean element rates if Llow < Lhigh
  253.         if (ll[0] < ll[1]) {
  254.             meanElementRate = getMeanElementRate(state, integrator, ll[0], ll[1], context, parameters);
  255.             if (isDirty) {
  256.                 boolean next = true;
  257.                 for (int i = 0; i < MAX_ORDER_RANK && next; i++) {
  258.                     final double[] meanRates = getMeanElementRate(state, new GaussQuadrature(GAUSS_ORDER[i]), ll[0],
  259.                             ll[1], context, parameters);
  260.                     if (getRatesDiff(meanElementRate, meanRates, context) < threshold) {
  261.                         integrator = new GaussQuadrature(GAUSS_ORDER[i]);
  262.                         next = false;
  263.                     }
  264.                 }
  265.                 isDirty = false;
  266.             }
  267.         }
  268.         return meanElementRate;
  269.     }

  270.     /** {@inheritDoc} */
  271.     @Override
  272.     public <T extends CalculusFieldElement<T>> T[] getMeanElementRate(final FieldSpacecraftState<T> state,
  273.             final FieldAuxiliaryElements<T> auxiliaryElements, final T[] parameters) {

  274.         // Container for attributes
  275.         final FieldAbstractGaussianContributionContext<T> context = initializeStep(auxiliaryElements, parameters);
  276.         final Field<T> field = state.getDate().getField();

  277.         T[] meanElementRate = MathArrays.buildArray(field, 6);
  278.         // Computes the limits for the integral
  279.         final T[] ll = getLLimits(state, auxiliaryElements);
  280.         // Computes integrated mean element rates if Llow < Lhigh
  281.         if (ll[0].getReal() < ll[1].getReal()) {
  282.             meanElementRate = getMeanElementRate(state, integrator, ll[0], ll[1], context, parameters);
  283.             if (isDirty) {
  284.                 boolean next = true;
  285.                 for (int i = 0; i < MAX_ORDER_RANK && next; i++) {
  286.                     final T[] meanRates = getMeanElementRate(state, new GaussQuadrature(GAUSS_ORDER[i]), ll[0], ll[1],
  287.                             context, parameters);
  288.                     if (getRatesDiff(meanElementRate, meanRates, context).getReal() < threshold) {
  289.                         integrator = new GaussQuadrature(GAUSS_ORDER[i]);
  290.                         next = false;
  291.                     }
  292.                 }
  293.                 isDirty = false;
  294.             }
  295.         }

  296.         return meanElementRate;
  297.     }

  298.     /**
  299.      * Compute the limits in L, the true longitude, for integration.
  300.      *
  301.      * @param state             current state information: date, kinematics,
  302.      *                          attitude
  303.      * @param auxiliaryElements auxiliary elements related to the current orbit
  304.      * @return the integration limits in L
  305.      */
  306.     protected abstract double[] getLLimits(SpacecraftState state, AuxiliaryElements auxiliaryElements);

  307.     /**
  308.      * Compute the limits in L, the true longitude, for integration.
  309.      *
  310.      * @param <T>               type of the elements
  311.      * @param state             current state information: date, kinematics,
  312.      *                          attitude
  313.      * @param auxiliaryElements auxiliary elements related to the current orbit
  314.      * @return the integration limits in L
  315.      */
  316.     protected abstract <T extends CalculusFieldElement<T>> T[] getLLimits(FieldSpacecraftState<T> state,
  317.             FieldAuxiliaryElements<T> auxiliaryElements);

  318.     /**
  319.      * Computes the mean equinoctial elements rates da<sub>i</sub> / dt.
  320.      *
  321.      * @param state      current state
  322.      * @param gauss      Gauss quadrature
  323.      * @param low        lower bound of the integral interval
  324.      * @param high       upper bound of the integral interval
  325.      * @param context    container for attributes
  326.      * @param parameters values of the force model parameters
  327.      * @return the mean element rates
  328.      */
  329.     protected double[] getMeanElementRate(final SpacecraftState state, final GaussQuadrature gauss, final double low,
  330.             final double high, final AbstractGaussianContributionContext context, final double[] parameters) {

  331.         // Auxiliary elements related to the current orbit
  332.         final AuxiliaryElements auxiliaryElements = context.getAuxiliaryElements();

  333.         final double[] meanElementRate = gauss.integrate(new IntegrableFunction(state, true, 0, parameters), low, high);

  334.         // Constant multiplier for integral
  335.         final double coef = 1. / (2. * FastMath.PI * auxiliaryElements.getB());
  336.         // Corrects mean element rates
  337.         for (int i = 0; i < 6; i++) {
  338.             meanElementRate[i] *= coef;
  339.         }
  340.         return meanElementRate;
  341.     }

  342.     /**
  343.      * Computes the mean equinoctial elements rates da<sub>i</sub> / dt.
  344.      *
  345.      * @param <T>        type of the elements
  346.      * @param state      current state
  347.      * @param gauss      Gauss quadrature
  348.      * @param low        lower bound of the integral interval
  349.      * @param high       upper bound of the integral interval
  350.      * @param context    container for attributes
  351.      * @param parameters values of the force model parameters
  352.      * @return the mean element rates
  353.      */
  354.     protected <T extends CalculusFieldElement<T>> T[] getMeanElementRate(final FieldSpacecraftState<T> state,
  355.             final GaussQuadrature gauss, final T low, final T high,
  356.             final FieldAbstractGaussianContributionContext<T> context, final T[] parameters) {

  357.         // Field
  358.         final Field<T> field = context.getA().getField();

  359.         // Auxiliary elements related to the current orbit
  360.         final FieldAuxiliaryElements<T> auxiliaryElements = context.getFieldAuxiliaryElements();

  361.         final T[] meanElementRate = gauss.integrate(new FieldIntegrableFunction<>(state, true, 0, parameters, field),
  362.                 low, high, field);
  363.         // Constant multiplier for integral
  364.         final T coef = auxiliaryElements.getB().multiply(low.getPi()).multiply(2.).reciprocal();
  365.         // Corrects mean element rates
  366.         for (int i = 0; i < 6; i++) {
  367.             meanElementRate[i] = meanElementRate[i].multiply(coef);
  368.         }
  369.         return meanElementRate;
  370.     }

  371.     /**
  372.      * Estimates the weighted magnitude of the difference between 2 sets of
  373.      * equinoctial elements rates.
  374.      *
  375.      * @param meanRef reference rates
  376.      * @param meanCur current rates
  377.      * @param context container for attributes
  378.      * @return estimated magnitude of weighted differences
  379.      */
  380.     private double getRatesDiff(final double[] meanRef, final double[] meanCur,
  381.             final AbstractGaussianContributionContext context) {

  382.         // Auxiliary elements related to the current orbit
  383.         final AuxiliaryElements auxiliaryElements = context.getAuxiliaryElements();

  384.         double maxDiff = FastMath.abs(meanRef[0] - meanCur[0]) / auxiliaryElements.getSma();
  385.         // Corrects mean element rates
  386.         for (int i = 1; i < meanRef.length; i++) {
  387.             maxDiff = FastMath.max(maxDiff, FastMath.abs(meanRef[i] - meanCur[i]));
  388.         }
  389.         return maxDiff;
  390.     }

  391.     /**
  392.      * Estimates the weighted magnitude of the difference between 2 sets of
  393.      * equinoctial elements rates.
  394.      *
  395.      * @param <T>     type of the elements
  396.      * @param meanRef reference rates
  397.      * @param meanCur current rates
  398.      * @param context container for attributes
  399.      * @return estimated magnitude of weighted differences
  400.      */
  401.     private <T extends CalculusFieldElement<T>> T getRatesDiff(final T[] meanRef, final T[] meanCur,
  402.             final FieldAbstractGaussianContributionContext<T> context) {

  403.         // Auxiliary elements related to the current orbit
  404.         final FieldAuxiliaryElements<T> auxiliaryElements = context.getFieldAuxiliaryElements();

  405.         T maxDiff = FastMath.abs(meanRef[0].subtract(meanCur[0])).divide(auxiliaryElements.getSma());
  406.         ;
  407.         // Corrects mean element rates
  408.         for (int i = 1; i < meanRef.length; i++) {
  409.             maxDiff = FastMath.max(maxDiff, FastMath.abs(meanRef[i].subtract(meanCur[i])));
  410.         }
  411.         return maxDiff;
  412.     }

  413.     /** {@inheritDoc} */
  414.     @Override
  415.     public void registerAttitudeProvider(final AttitudeProvider provider) {
  416.         this.attitudeProvider = provider;
  417.     }

  418.     /** {@inheritDoc} */
  419.     @Override
  420.     public void updateShortPeriodTerms(final double[] parameters, final SpacecraftState... meanStates) {

  421.         final Slot slot = gaussianSPCoefs.createSlot(meanStates);
  422.         for (final SpacecraftState meanState : meanStates) {

  423.             // Auxiliary elements related to the current orbit
  424.             final AuxiliaryElements auxiliaryElements = new AuxiliaryElements(meanState.getOrbit(), I);

  425.             // Container of attributes
  426.             final AbstractGaussianContributionContext context = initializeStep(auxiliaryElements, parameters);

  427.             // Compute rhoj and sigmaj
  428.             final double[][] currentRhoSigmaj = computeRhoSigmaCoefficients(meanState.getDate(), auxiliaryElements);

  429.             // Generate the Cij and Sij coefficients
  430.             final FourierCjSjCoefficients fourierCjSj = new FourierCjSjCoefficients(meanState, JMAX, auxiliaryElements,
  431.                     parameters);

  432.             // Generate the Uij and Vij coefficients
  433.             final UijVijCoefficients uijvij = new UijVijCoefficients(currentRhoSigmaj, fourierCjSj, JMAX);

  434.             gaussianSPCoefs.computeCoefficients(meanState, slot, fourierCjSj, uijvij, context.getMeanMotion(),
  435.                     auxiliaryElements.getSma());

  436.         }

  437.     }

  438.     /** {@inheritDoc} */
  439.     @Override
  440.     @SuppressWarnings("unchecked")
  441.     public <T extends CalculusFieldElement<T>> void updateShortPeriodTerms(final T[] parameters,
  442.             final FieldSpacecraftState<T>... meanStates) {

  443.         // Field used by default
  444.         final Field<T> field = meanStates[0].getDate().getField();

  445.         final FieldGaussianShortPeriodicCoefficients<T> fgspc = (FieldGaussianShortPeriodicCoefficients<T>) gaussianFieldSPCoefs
  446.                 .get(field);
  447.         final FieldSlot<T> slot = fgspc.createSlot(meanStates);
  448.         for (final FieldSpacecraftState<T> meanState : meanStates) {

  449.             // Auxiliary elements related to the current orbit
  450.             final FieldAuxiliaryElements<T> auxiliaryElements = new FieldAuxiliaryElements<>(meanState.getOrbit(), I);

  451.             // Container of attributes
  452.             final FieldAbstractGaussianContributionContext<T> context = initializeStep(auxiliaryElements, parameters);

  453.             // Compute rhoj and sigmaj
  454.             final T[][] currentRhoSigmaj = computeRhoSigmaCoefficients(meanState.getDate(), context, field);

  455.             // Generate the Cij and Sij coefficients
  456.             final FieldFourierCjSjCoefficients<T> fourierCjSj = new FieldFourierCjSjCoefficients<>(meanState, JMAX,
  457.                     auxiliaryElements, parameters, field);

  458.             // Generate the Uij and Vij coefficients
  459.             final FieldUijVijCoefficients<T> uijvij = new FieldUijVijCoefficients<>(currentRhoSigmaj, fourierCjSj, JMAX,
  460.                     field);

  461.             fgspc.computeCoefficients(meanState, slot, fourierCjSj, uijvij, context.getMeanMotion(),
  462.                     auxiliaryElements.getSma(), field);

  463.         }

  464.     }

  465.     /**
  466.      * Compute the auxiliary quantities ρ<sub>j</sub> and σ<sub>j</sub>.
  467.      * <p>
  468.      * The expressions used are equations 2.5.3-(4) from the Danielson paper. <br/>
  469.      * ρ<sub>j</sub> = (1+jB)(-b)<sup>j</sup>C<sub>j</sub>(k, h) <br/>
  470.      * σ<sub>j</sub> = (1+jB)(-b)<sup>j</sup>S<sub>j</sub>(k, h) <br/>
  471.      * </p>
  472.      * @param date              current date
  473.      * @param auxiliaryElements auxiliary elements related to the current orbit
  474.      * @return computed coefficients
  475.      */
  476.     private double[][] computeRhoSigmaCoefficients(final AbsoluteDate date, final AuxiliaryElements auxiliaryElements) {
  477.         final double[][] currentRhoSigmaj = new double[2][3 * JMAX + 1];
  478.         final CjSjCoefficient cjsjKH = new CjSjCoefficient(auxiliaryElements.getK(), auxiliaryElements.getH());
  479.         final double b = 1. / (1 + auxiliaryElements.getB());

  480.         // (-b)<sup>j</sup>
  481.         double mbtj = 1;

  482.         for (int j = 1; j <= 3 * JMAX; j++) {

  483.             // Compute current rho and sigma;
  484.             mbtj *= -b;
  485.             final double coef = (1 + j * auxiliaryElements.getB()) * mbtj;
  486.             currentRhoSigmaj[0][j] = coef * cjsjKH.getCj(j);
  487.             currentRhoSigmaj[1][j] = coef * cjsjKH.getSj(j);
  488.         }
  489.         return currentRhoSigmaj;
  490.     }

  491.     /**
  492.      * Compute the auxiliary quantities ρ<sub>j</sub> and σ<sub>j</sub>.
  493.      * <p>
  494.      * The expressions used are equations 2.5.3-(4) from the Danielson paper. <br/>
  495.      * ρ<sub>j</sub> = (1+jB)(-b)<sup>j</sup>C<sub>j</sub>(k, h) <br/>
  496.      * σ<sub>j</sub> = (1+jB)(-b)<sup>j</sup>S<sub>j</sub>(k, h) <br/>
  497.      * </p>
  498.      * @param <T>     type of the elements
  499.      * @param date    current date
  500.      * @param context container for attributes
  501.      * @param field   field used by default
  502.      * @return computed coefficients
  503.      */
  504.     private <T extends CalculusFieldElement<T>> T[][] computeRhoSigmaCoefficients(final FieldAbsoluteDate<T> date,
  505.             final FieldAbstractGaussianContributionContext<T> context, final Field<T> field) {
  506.         // zero
  507.         final T zero = field.getZero();

  508.         final FieldAuxiliaryElements<T> auxiliaryElements = context.getFieldAuxiliaryElements();
  509.         final T[][] currentRhoSigmaj = MathArrays.buildArray(field, 2, 3 * JMAX + 1);
  510.         final FieldCjSjCoefficient<T> cjsjKH = new FieldCjSjCoefficient<>(auxiliaryElements.getK(),
  511.                 auxiliaryElements.getH(), field);
  512.         final T b = auxiliaryElements.getB().add(1.).reciprocal();

  513.         // (-b)<sup>j</sup>
  514.         T mbtj = zero.add(1.);

  515.         for (int j = 1; j <= 3 * JMAX; j++) {

  516.             // Compute current rho and sigma;
  517.             mbtj = mbtj.multiply(b.negate());
  518.             final T coef = mbtj.multiply(auxiliaryElements.getB().multiply(j).add(1.));
  519.             currentRhoSigmaj[0][j] = coef.multiply(cjsjKH.getCj(j));
  520.             currentRhoSigmaj[1][j] = coef.multiply(cjsjKH.getSj(j));
  521.         }
  522.         return currentRhoSigmaj;
  523.     }

  524.     /**
  525.      * Internal class for numerical quadrature.
  526.      * <p>
  527.      * This class is a rewrite of {@link IntegrableFunction} for field elements
  528.      * </p>
  529.      */
  530.     protected class FieldIntegrableFunction<T extends CalculusFieldElement<T>>
  531.             implements CalculusFieldUnivariateVectorFunction<T> {

  532.         /** Current state. */
  533.         private final FieldSpacecraftState<T> state;

  534.         /**
  535.          * Signal that this class is used to compute the values required by the mean
  536.          * element variations or by the short periodic element variations.
  537.          */
  538.         private final boolean meanMode;

  539.         /**
  540.          * The j index.
  541.          * <p>
  542.          * Used only for short periodic variation. Ignored for mean elements variation.
  543.          * </p>
  544.          */
  545.         private final int j;

  546.         /** Container for attributes. */
  547.         private final FieldAbstractGaussianContributionContext<T> context;

  548.         /** Auxiliary Elements. */
  549.         private final FieldAuxiliaryElements<T> auxiliaryElements;

  550.         /** Drivers for solar radiation and atmospheric drag forces. */
  551.         private final T[] parameters;

  552.         /**
  553.          * Build a new instance with a new field.
  554.          * @param state      current state information: date, kinematics, attitude
  555.          * @param meanMode   if true return the value associated to the mean elements
  556.          *                   variation, if false return the values associated to the
  557.          *                   short periodic elements variation
  558.          * @param j          the j index. used only for short periodic variation.
  559.          *                   Ignored for mean elements variation.
  560.          * @param parameters values of the force model parameters
  561.          * @param field      field utilized by default
  562.          */
  563.         public FieldIntegrableFunction(final FieldSpacecraftState<T> state, final boolean meanMode, final int j,
  564.                 final T[] parameters, final Field<T> field) {

  565.             this.meanMode = meanMode;
  566.             this.j = j;
  567.             this.parameters = parameters.clone();
  568.             this.auxiliaryElements = new FieldAuxiliaryElements<>(state.getOrbit(), I);
  569.             this.context = new FieldAbstractGaussianContributionContext<>(auxiliaryElements, this.parameters);
  570.             // remove derivatives from state
  571.             final T[] stateVector = MathArrays.buildArray(field, 6);
  572.             OrbitType.EQUINOCTIAL.mapOrbitToArray(state.getOrbit(), PositionAngle.TRUE, stateVector, null);
  573.             final FieldOrbit<T> fixedOrbit = OrbitType.EQUINOCTIAL.mapArrayToOrbit(stateVector, null,
  574.                     PositionAngle.TRUE, state.getDate(), context.getMu(), state.getFrame());
  575.             this.state = new FieldSpacecraftState<>(fixedOrbit, state.getAttitude(), state.getMass());
  576.         }

  577.         /** {@inheritDoc} */
  578.         @Override
  579.         public T[] value(final T x) {

  580.             // Parameters for array building
  581.             final Field<T> field = auxiliaryElements.getDate().getField();
  582.             final int dimension = 6;

  583.             // Compute the time difference from the true longitude difference
  584.             final T shiftedLm = trueToMean(x);
  585.             final T dLm = shiftedLm.subtract(auxiliaryElements.getLM());
  586.             final T dt = dLm.divide(context.getMeanMotion());

  587.             final FieldSinCos<T> scL = FastMath.sinCos(x);
  588.             final T cosL = scL.cos();
  589.             final T sinL = scL.sin();
  590.             final T roa  = auxiliaryElements.getB().multiply(auxiliaryElements.getB()).divide(auxiliaryElements.getH().multiply(sinL).add(auxiliaryElements.getK().multiply(cosL)).add(1.));
  591.             final T roa2 = roa.multiply(roa);
  592.             final T r = auxiliaryElements.getSma().multiply(roa);
  593.             final T X = r.multiply(cosL);
  594.             final T Y = r.multiply(sinL);
  595.             final T naob = context.getMeanMotion().multiply(auxiliaryElements.getSma())
  596.                     .divide(auxiliaryElements.getB());
  597.             final T Xdot = naob.multiply(auxiliaryElements.getH().add(sinL)).negate();
  598.             final T Ydot = naob.multiply(auxiliaryElements.getK().add(cosL));
  599.             final FieldVector3D<T> vel = new FieldVector3D<>(Xdot, auxiliaryElements.getVectorF(), Ydot,
  600.                     auxiliaryElements.getVectorG());

  601.             // Compute acceleration
  602.             FieldVector3D<T> acc = FieldVector3D.getZero(field);

  603.             // shift the orbit to dt
  604.             final FieldOrbit<T> shiftedOrbit = state.getOrbit().shiftedBy(dt);

  605.             // Recompose an orbit with time held fixed to be compliant with DSST theory
  606.             final FieldOrbit<T> recomposedOrbit = new FieldEquinoctialOrbit<>(shiftedOrbit.getA(),
  607.                     shiftedOrbit.getEquinoctialEx(), shiftedOrbit.getEquinoctialEy(), shiftedOrbit.getHx(),
  608.                     shiftedOrbit.getHy(), shiftedOrbit.getLv(), PositionAngle.TRUE, shiftedOrbit.getFrame(),
  609.                     state.getDate(), context.getMu());

  610.             // Get the corresponding attitude
  611.             final FieldAttitude<T> recomposedAttitude = attitudeProvider.getAttitude(recomposedOrbit,
  612.                     recomposedOrbit.getDate(), recomposedOrbit.getFrame());

  613.             // create shifted SpacecraftState with attitude at specified time
  614.             final FieldSpacecraftState<T> shiftedState = new FieldSpacecraftState<>(recomposedOrbit, recomposedAttitude,
  615.                     state.getMass());

  616.             acc = contribution.acceleration(shiftedState, parameters);

  617.             // Compute the derivatives of the elements by the speed
  618.             final T[] deriv = MathArrays.buildArray(field, dimension);
  619.             // da/dv
  620.             deriv[0] = getAoV(vel).dotProduct(acc);
  621.             // dex/dv
  622.             deriv[1] = getKoV(X, Y, Xdot, Ydot).dotProduct(acc);
  623.             // dey/dv
  624.             deriv[2] = getHoV(X, Y, Xdot, Ydot).dotProduct(acc);
  625.             // dhx/dv
  626.             deriv[3] = getQoV(X).dotProduct(acc);
  627.             // dhy/dv
  628.             deriv[4] = getPoV(Y).dotProduct(acc);
  629.             // dλ/dv
  630.             deriv[5] = getLoV(X, Y, Xdot, Ydot).dotProduct(acc);

  631.             // Compute mean elements rates
  632.             T[] val = null;
  633.             if (meanMode) {
  634.                 val = MathArrays.buildArray(field, dimension);
  635.                 for (int i = 0; i < 6; i++) {
  636.                     // da<sub>i</sub>/dt
  637.                     val[i] = deriv[i].multiply(roa2);
  638.                 }
  639.             } else {
  640.                 val = MathArrays.buildArray(field, dimension * 2);
  641.                 //Compute cos(j*L) and sin(j*L);
  642.                 final FieldSinCos<T> scjL = FastMath.sinCos(x.multiply(j));
  643.                 final T cosjL = j == 1 ? cosL : scjL.cos();
  644.                 final T sinjL = j == 1 ? sinL : scjL.sin();

  645.                 for (int i = 0; i < 6; i++) {
  646.                     // da<sub>i</sub>/dv * cos(jL)
  647.                     val[i] = deriv[i].multiply(cosjL);
  648.                     // da<sub>i</sub>/dv * sin(jL)
  649.                     val[i + 6] = deriv[i].multiply(sinjL);
  650.                 }
  651.             }

  652.             return val;
  653.         }

  654.         /**
  655.          * Converts true longitude to mean longitude.
  656.          * @param x True longitude
  657.          * @return Eccentric longitude
  658.          */
  659.         private T trueToMean(final T x) {
  660.             return eccentricToMean(trueToEccentric(x));
  661.         }

  662.         /**
  663.          * Converts true longitude to eccentric longitude.
  664.          * @param lv True longitude
  665.          * @return Eccentric longitude
  666.          */
  667.         private T trueToEccentric (final T lv) {
  668.             final FieldSinCos<T> sclV = FastMath.sinCos(lv);
  669.             final T cosLv   = sclV.cos();
  670.             final T sinLv   = sclV.sin();
  671.             final T num     = auxiliaryElements.getH().multiply(cosLv).subtract(auxiliaryElements.getK().multiply(sinLv));
  672.             final T den     = auxiliaryElements.getB().add(auxiliaryElements.getK().multiply(cosLv)).add(auxiliaryElements.getH().multiply(sinLv)).add(1.);
  673.             return FastMath.atan(num.divide(den)).multiply(2.).add(lv);
  674.         }

  675.         /**
  676.          * Converts eccentric longitude to mean longitude.
  677.          * @param le Eccentric longitude
  678.          * @return Mean longitude
  679.          */
  680.         private T eccentricToMean (final T le) {
  681.             final FieldSinCos<T> scle = FastMath.sinCos(le);
  682.             return le.subtract(auxiliaryElements.getK().multiply(scle.sin())).add(auxiliaryElements.getH().multiply(scle.cos()));
  683.         }

  684.         /**
  685.          * Compute δa/δv.
  686.          * @param vel satellite velocity
  687.          * @return δa/δv
  688.          */
  689.         private FieldVector3D<T> getAoV(final FieldVector3D<T> vel) {
  690.             return new FieldVector3D<>(context.getTon2a(), vel);
  691.         }

  692.         /**
  693.          * Compute δh/δv.
  694.          * @param X    satellite position component along f, equinoctial reference frame
  695.          *             1st vector
  696.          * @param Y    satellite position component along g, equinoctial reference frame
  697.          *             2nd vector
  698.          * @param Xdot satellite velocity component along f, equinoctial reference frame
  699.          *             1st vector
  700.          * @param Ydot satellite velocity component along g, equinoctial reference frame
  701.          *             2nd vector
  702.          * @return δh/δv
  703.          */
  704.         private FieldVector3D<T> getHoV(final T X, final T Y, final T Xdot, final T Ydot) {
  705.             final T kf = (Xdot.multiply(Y).multiply(2.).subtract(X.multiply(Ydot))).multiply(context.getOoMU());
  706.             final T kg = X.multiply(Xdot).multiply(context.getOoMU());
  707.             final T kw = auxiliaryElements.getK().multiply(
  708.                     auxiliaryElements.getQ().multiply(Y).multiply(I).subtract(auxiliaryElements.getP().multiply(X)))
  709.                     .multiply(context.getOOAB());
  710.             return new FieldVector3D<>(kf, auxiliaryElements.getVectorF(), kg.negate(), auxiliaryElements.getVectorG(),
  711.                     kw, auxiliaryElements.getVectorW());
  712.         }

  713.         /**
  714.          * Compute δk/δv.
  715.          * @param X    satellite position component along f, equinoctial reference frame
  716.          *             1st vector
  717.          * @param Y    satellite position component along g, equinoctial reference frame
  718.          *             2nd vector
  719.          * @param Xdot satellite velocity component along f, equinoctial reference frame
  720.          *             1st vector
  721.          * @param Ydot satellite velocity component along g, equinoctial reference frame
  722.          *             2nd vector
  723.          * @return δk/δv
  724.          */
  725.         private FieldVector3D<T> getKoV(final T X, final T Y, final T Xdot, final T Ydot) {
  726.             final T kf = Y.multiply(Ydot).multiply(context.getOoMU());
  727.             final T kg = (X.multiply(Ydot).multiply(2.).subtract(Xdot.multiply(Y))).multiply(context.getOoMU());
  728.             final T kw = auxiliaryElements.getH().multiply(
  729.                     auxiliaryElements.getQ().multiply(Y).multiply(I).subtract(auxiliaryElements.getP().multiply(X)))
  730.                     .multiply(context.getOOAB());
  731.             return new FieldVector3D<>(kf.negate(), auxiliaryElements.getVectorF(), kg, auxiliaryElements.getVectorG(),
  732.                     kw.negate(), auxiliaryElements.getVectorW());
  733.         }

  734.         /**
  735.          * Compute δp/δv.
  736.          * @param Y satellite position component along g, equinoctial reference frame
  737.          *          2nd vector
  738.          * @return δp/δv
  739.          */
  740.         private FieldVector3D<T> getPoV(final T Y) {
  741.             return new FieldVector3D<>(context.getCo2AB().multiply(Y), auxiliaryElements.getVectorW());
  742.         }

  743.         /**
  744.          * Compute δq/δv.
  745.          * @param X satellite position component along f, equinoctial reference frame
  746.          *          1st vector
  747.          * @return δq/δv
  748.          */
  749.         private FieldVector3D<T> getQoV(final T X) {
  750.             return new FieldVector3D<>(context.getCo2AB().multiply(X).multiply(I), auxiliaryElements.getVectorW());
  751.         }

  752.         /**
  753.          * Compute δλ/δv.
  754.          * @param X    satellite position component along f, equinoctial reference frame
  755.          *             1st vector
  756.          * @param Y    satellite position component along g, equinoctial reference frame
  757.          *             2nd vector
  758.          * @param Xdot satellite velocity component along f, equinoctial reference frame
  759.          *             1st vector
  760.          * @param Ydot satellite velocity component along g, equinoctial reference frame
  761.          *             2nd vector
  762.          * @return δλ/δv
  763.          */
  764.         private FieldVector3D<T> getLoV(final T X, final T Y, final T Xdot, final T Ydot) {
  765.             final FieldVector3D<T> pos = new FieldVector3D<>(X, auxiliaryElements.getVectorF(), Y,
  766.                     auxiliaryElements.getVectorG());
  767.             final FieldVector3D<T> v2 = new FieldVector3D<>(auxiliaryElements.getK(), getHoV(X, Y, Xdot, Ydot),
  768.                     auxiliaryElements.getH().negate(), getKoV(X, Y, Xdot, Ydot));
  769.             return new FieldVector3D<>(context.getOOA().multiply(-2.), pos, context.getOoBpo(), v2,
  770.                     context.getOOA().multiply(auxiliaryElements.getQ().multiply(Y).multiply(I)
  771.                             .subtract(auxiliaryElements.getP().multiply(X))),
  772.                     auxiliaryElements.getVectorW());
  773.         }

  774.     }

  775.     /** Internal class for numerical quadrature. */
  776.     protected class IntegrableFunction implements UnivariateVectorFunction {

  777.         /** Current state. */
  778.         private final SpacecraftState state;

  779.         /**
  780.          * Signal that this class is used to compute the values required by the mean
  781.          * element variations or by the short periodic element variations.
  782.          */
  783.         private final boolean meanMode;

  784.         /**
  785.          * The j index.
  786.          * <p>
  787.          * Used only for short periodic variation. Ignored for mean elements variation.
  788.          * </p>
  789.          */
  790.         private final int j;

  791.         /** Container for attributes. */
  792.         private final AbstractGaussianContributionContext context;

  793.         /** Auxiliary Elements. */
  794.         private final AuxiliaryElements auxiliaryElements;

  795.         /** Drivers for solar radiation and atmospheric drag forces. */
  796.         private final double[] parameters;

  797.         /**
  798.          * Build a new instance.
  799.          * @param state      current state information: date, kinematics, attitude
  800.          * @param meanMode   if true return the value associated to the mean elements
  801.          *                   variation, if false return the values associated to the
  802.          *                   short periodic elements variation
  803.          * @param j          the j index. used only for short periodic variation.
  804.          *                   Ignored for mean elements variation.
  805.          * @param parameters values of the force model parameters
  806.          */
  807.         IntegrableFunction(final SpacecraftState state, final boolean meanMode, final int j,
  808.                 final double[] parameters) {

  809.             this.meanMode = meanMode;
  810.             this.j = j;
  811.             this.parameters = parameters.clone();
  812.             this.auxiliaryElements = new AuxiliaryElements(state.getOrbit(), I);
  813.             this.context = new AbstractGaussianContributionContext(auxiliaryElements, this.parameters);
  814.             // remove derivatives from state
  815.             final double[] stateVector = new double[6];
  816.             OrbitType.EQUINOCTIAL.mapOrbitToArray(state.getOrbit(), PositionAngle.TRUE, stateVector, null);
  817.             final Orbit fixedOrbit = OrbitType.EQUINOCTIAL.mapArrayToOrbit(stateVector, null, PositionAngle.TRUE,
  818.                     state.getDate(), context.getMu(), state.getFrame());
  819.             this.state = new SpacecraftState(fixedOrbit, state.getAttitude(), state.getMass());
  820.         }

  821.         /** {@inheritDoc} */
  822.         @Override
  823.         public double[] value(final double x) {

  824.             // Compute the time difference from the true longitude difference
  825.             final double shiftedLm = trueToMean(x);
  826.             final double dLm = shiftedLm - auxiliaryElements.getLM();
  827.             final double dt = dLm / context.getMeanMotion();

  828.             final SinCos scL  = FastMath.sinCos(x);
  829.             final double cosL = scL.cos();
  830.             final double sinL = scL.sin();
  831.             final double roa  = auxiliaryElements.getB() * auxiliaryElements.getB() / (1. + auxiliaryElements.getH() * sinL + auxiliaryElements.getK() * cosL);
  832.             final double roa2 = roa * roa;
  833.             final double r = auxiliaryElements.getSma() * roa;
  834.             final double X = r * cosL;
  835.             final double Y = r * sinL;
  836.             final double naob = context.getMeanMotion() * auxiliaryElements.getSma() / auxiliaryElements.getB();
  837.             final double Xdot = -naob * (auxiliaryElements.getH() + sinL);
  838.             final double Ydot = naob * (auxiliaryElements.getK() + cosL);
  839.             final Vector3D vel = new Vector3D(Xdot, auxiliaryElements.getVectorF(), Ydot,
  840.                     auxiliaryElements.getVectorG());

  841.             // Compute acceleration
  842.             Vector3D acc = Vector3D.ZERO;

  843.             // shift the orbit to dt
  844.             final Orbit shiftedOrbit = state.getOrbit().shiftedBy(dt);

  845.             // Recompose an orbit with time held fixed to be compliant with DSST theory
  846.             final Orbit recomposedOrbit = new EquinoctialOrbit(shiftedOrbit.getA(), shiftedOrbit.getEquinoctialEx(),
  847.                     shiftedOrbit.getEquinoctialEy(), shiftedOrbit.getHx(), shiftedOrbit.getHy(), shiftedOrbit.getLv(),
  848.                     PositionAngle.TRUE, shiftedOrbit.getFrame(), state.getDate(), context.getMu());

  849.             // Get the corresponding attitude
  850.             final Attitude recomposedAttitude = attitudeProvider.getAttitude(recomposedOrbit, recomposedOrbit.getDate(),
  851.                     recomposedOrbit.getFrame());

  852.             // create shifted SpacecraftState with attitude at specified time
  853.             final SpacecraftState shiftedState = new SpacecraftState(recomposedOrbit, recomposedAttitude,
  854.                     state.getMass());

  855.             acc = contribution.acceleration(shiftedState, parameters);

  856.             // Compute the derivatives of the elements by the speed
  857.             final double[] deriv = new double[6];
  858.             // da/dv
  859.             deriv[0] = getAoV(vel).dotProduct(acc);
  860.             // dex/dv
  861.             deriv[1] = getKoV(X, Y, Xdot, Ydot).dotProduct(acc);
  862.             // dey/dv
  863.             deriv[2] = getHoV(X, Y, Xdot, Ydot).dotProduct(acc);
  864.             // dhx/dv
  865.             deriv[3] = getQoV(X).dotProduct(acc);
  866.             // dhy/dv
  867.             deriv[4] = getPoV(Y).dotProduct(acc);
  868.             // dλ/dv
  869.             deriv[5] = getLoV(X, Y, Xdot, Ydot).dotProduct(acc);

  870.             // Compute mean elements rates
  871.             double[] val = null;
  872.             if (meanMode) {
  873.                 val = new double[6];
  874.                 for (int i = 0; i < 6; i++) {
  875.                     // da<sub>i</sub>/dt
  876.                     val[i] = roa2 * deriv[i];
  877.                 }
  878.             } else {
  879.                 val = new double[12];
  880.                 //Compute cos(j*L) and sin(j*L);
  881.                 final SinCos scjL  = FastMath.sinCos(j * x);
  882.                 final double cosjL = j == 1 ? cosL : scjL.cos();
  883.                 final double sinjL = j == 1 ? sinL : scjL.sin();

  884.                 for (int i = 0; i < 6; i++) {
  885.                     // da<sub>i</sub>/dv * cos(jL)
  886.                     val[i] = cosjL * deriv[i];
  887.                     // da<sub>i</sub>/dv * sin(jL)
  888.                     val[i + 6] = sinjL * deriv[i];
  889.                 }
  890.             }
  891.             return val;
  892.         }

  893.         /**
  894.          * Converts true longitude to eccentric longitude.
  895.          * @param lv True longitude
  896.          * @return Eccentric longitude
  897.          */
  898.         private double trueToEccentric (final double lv) {
  899.             final SinCos scLv    = FastMath.sinCos(lv);
  900.             final double num     = auxiliaryElements.getH() * scLv.cos() - auxiliaryElements.getK() * scLv.sin();
  901.             final double den     = auxiliaryElements.getB() + 1. + auxiliaryElements.getK() * scLv.cos() + auxiliaryElements.getH() * scLv.sin();
  902.             return lv + 2. * FastMath.atan(num / den);
  903.         }

  904.         /**
  905.          * Converts eccentric longitude to mean longitude.
  906.          * @param le Eccentric longitude
  907.          * @return Mean longitude
  908.          */
  909.         private double eccentricToMean (final double le) {
  910.             final SinCos scLe = FastMath.sinCos(le);
  911.             return le - auxiliaryElements.getK() * scLe.sin() + auxiliaryElements.getH() * scLe.cos();
  912.         }

  913.         /**
  914.          * Converts true longitude to mean longitude.
  915.          * @param lv True longitude
  916.          * @return Eccentric longitude
  917.          */
  918.         private double trueToMean(final double lv) {
  919.             return eccentricToMean(trueToEccentric(lv));
  920.         }

  921.         /**
  922.          * Compute δa/δv.
  923.          * @param vel satellite velocity
  924.          * @return δa/δv
  925.          */
  926.         private Vector3D getAoV(final Vector3D vel) {
  927.             return new Vector3D(context.getTon2a(), vel);
  928.         }

  929.         /**
  930.          * Compute δh/δv.
  931.          * @param X    satellite position component along f, equinoctial reference frame
  932.          *             1st vector
  933.          * @param Y    satellite position component along g, equinoctial reference frame
  934.          *             2nd vector
  935.          * @param Xdot satellite velocity component along f, equinoctial reference frame
  936.          *             1st vector
  937.          * @param Ydot satellite velocity component along g, equinoctial reference frame
  938.          *             2nd vector
  939.          * @return δh/δv
  940.          */
  941.         private Vector3D getHoV(final double X, final double Y, final double Xdot, final double Ydot) {
  942.             final double kf = (2. * Xdot * Y - X * Ydot) * context.getOoMU();
  943.             final double kg = X * Xdot * context.getOoMU();
  944.             final double kw = auxiliaryElements.getK() *
  945.                     (I * auxiliaryElements.getQ() * Y - auxiliaryElements.getP() * X) * context.getOOAB();
  946.             return new Vector3D(kf, auxiliaryElements.getVectorF(), -kg, auxiliaryElements.getVectorG(), kw,
  947.                     auxiliaryElements.getVectorW());
  948.         }

  949.         /**
  950.          * Compute δk/δv.
  951.          * @param X    satellite position component along f, equinoctial reference frame
  952.          *             1st vector
  953.          * @param Y    satellite position component along g, equinoctial reference frame
  954.          *             2nd vector
  955.          * @param Xdot satellite velocity component along f, equinoctial reference frame
  956.          *             1st vector
  957.          * @param Ydot satellite velocity component along g, equinoctial reference frame
  958.          *             2nd vector
  959.          * @return δk/δv
  960.          */
  961.         private Vector3D getKoV(final double X, final double Y, final double Xdot, final double Ydot) {
  962.             final double kf = Y * Ydot * context.getOoMU();
  963.             final double kg = (2. * X * Ydot - Xdot * Y) * context.getOoMU();
  964.             final double kw = auxiliaryElements.getH() *
  965.                     (I * auxiliaryElements.getQ() * Y - auxiliaryElements.getP() * X) * context.getOOAB();
  966.             return new Vector3D(-kf, auxiliaryElements.getVectorF(), kg, auxiliaryElements.getVectorG(), -kw,
  967.                     auxiliaryElements.getVectorW());
  968.         }

  969.         /**
  970.          * Compute δp/δv.
  971.          * @param Y satellite position component along g, equinoctial reference frame
  972.          *          2nd vector
  973.          * @return δp/δv
  974.          */
  975.         private Vector3D getPoV(final double Y) {
  976.             return new Vector3D(context.getCo2AB() * Y, auxiliaryElements.getVectorW());
  977.         }

  978.         /**
  979.          * Compute δq/δv.
  980.          * @param X satellite position component along f, equinoctial reference frame
  981.          *          1st vector
  982.          * @return δq/δv
  983.          */
  984.         private Vector3D getQoV(final double X) {
  985.             return new Vector3D(I * context.getCo2AB() * X, auxiliaryElements.getVectorW());
  986.         }

  987.         /**
  988.          * Compute δλ/δv.
  989.          * @param X    satellite position component along f, equinoctial reference frame
  990.          *             1st vector
  991.          * @param Y    satellite position component along g, equinoctial reference frame
  992.          *             2nd vector
  993.          * @param Xdot satellite velocity component along f, equinoctial reference frame
  994.          *             1st vector
  995.          * @param Ydot satellite velocity component along g, equinoctial reference frame
  996.          *             2nd vector
  997.          * @return δλ/δv
  998.          */
  999.         private Vector3D getLoV(final double X, final double Y, final double Xdot, final double Ydot) {
  1000.             final Vector3D pos = new Vector3D(X, auxiliaryElements.getVectorF(), Y, auxiliaryElements.getVectorG());
  1001.             final Vector3D v2 = new Vector3D(auxiliaryElements.getK(), getHoV(X, Y, Xdot, Ydot),
  1002.                     -auxiliaryElements.getH(), getKoV(X, Y, Xdot, Ydot));
  1003.             return new Vector3D(-2. * context.getOOA(), pos, context.getOoBpo(), v2,
  1004.                     (I * auxiliaryElements.getQ() * Y - auxiliaryElements.getP() * X) * context.getOOA(),
  1005.                     auxiliaryElements.getVectorW());
  1006.         }

  1007.     }

  1008.     /**
  1009.      * Class used to {@link #integrate(UnivariateVectorFunction, double, double)
  1010.      * integrate} a {@link org.hipparchus.analysis.UnivariateVectorFunction
  1011.      * function} of the orbital elements using the Gaussian quadrature rule to get
  1012.      * the acceleration.
  1013.      */
  1014.     protected static class GaussQuadrature {

  1015.         // Points and weights for the available quadrature orders

  1016.         /** Points for quadrature of order 12. */
  1017.         private static final double[] P_12 = { -0.98156063424671910000, -0.90411725637047490000,
  1018.             -0.76990267419430470000, -0.58731795428661740000, -0.36783149899818024000, -0.12523340851146890000,
  1019.             0.12523340851146890000, 0.36783149899818024000, 0.58731795428661740000, 0.76990267419430470000,
  1020.             0.90411725637047490000, 0.98156063424671910000 };

  1021.         /** Weights for quadrature of order 12. */
  1022.         private static final double[] W_12 = { 0.04717533638651220000, 0.10693932599531830000, 0.16007832854334633000,
  1023.             0.20316742672306584000, 0.23349253653835478000, 0.24914704581340286000, 0.24914704581340286000,
  1024.             0.23349253653835478000, 0.20316742672306584000, 0.16007832854334633000, 0.10693932599531830000,
  1025.             0.04717533638651220000 };

  1026.         /** Points for quadrature of order 16. */
  1027.         private static final double[] P_16 = { -0.98940093499164990000, -0.94457502307323260000,
  1028.             -0.86563120238783160000, -0.75540440835500310000, -0.61787624440264380000, -0.45801677765722737000,
  1029.             -0.28160355077925890000, -0.09501250983763745000, 0.09501250983763745000, 0.28160355077925890000,
  1030.             0.45801677765722737000, 0.61787624440264380000, 0.75540440835500310000, 0.86563120238783160000,
  1031.             0.94457502307323260000, 0.98940093499164990000 };

  1032.         /** Weights for quadrature of order 16. */
  1033.         private static final double[] W_16 = { 0.02715245941175405800, 0.06225352393864777000, 0.09515851168249283000,
  1034.             0.12462897125553388000, 0.14959598881657685000, 0.16915651939500256000, 0.18260341504492360000,
  1035.             0.18945061045506847000, 0.18945061045506847000, 0.18260341504492360000, 0.16915651939500256000,
  1036.             0.14959598881657685000, 0.12462897125553388000, 0.09515851168249283000, 0.06225352393864777000,
  1037.             0.02715245941175405800 };

  1038.         /** Points for quadrature of order 20. */
  1039.         private static final double[] P_20 = { -0.99312859918509490000, -0.96397192727791390000,
  1040.             -0.91223442825132600000, -0.83911697182221890000, -0.74633190646015080000, -0.63605368072651510000,
  1041.             -0.51086700195082700000, -0.37370608871541955000, -0.22778585114164507000, -0.07652652113349734000,
  1042.             0.07652652113349734000, 0.22778585114164507000, 0.37370608871541955000, 0.51086700195082700000,
  1043.             0.63605368072651510000, 0.74633190646015080000, 0.83911697182221890000, 0.91223442825132600000,
  1044.             0.96397192727791390000, 0.99312859918509490000 };

  1045.         /** Weights for quadrature of order 20. */
  1046.         private static final double[] W_20 = { 0.01761400713915226400, 0.04060142980038684000, 0.06267204833410904000,
  1047.             0.08327674157670477000, 0.10193011981724048000, 0.11819453196151844000, 0.13168863844917678000,
  1048.             0.14209610931838212000, 0.14917298647260380000, 0.15275338713072600000, 0.15275338713072600000,
  1049.             0.14917298647260380000, 0.14209610931838212000, 0.13168863844917678000, 0.11819453196151844000,
  1050.             0.10193011981724048000, 0.08327674157670477000, 0.06267204833410904000, 0.04060142980038684000,
  1051.             0.01761400713915226400 };

  1052.         /** Points for quadrature of order 24. */
  1053.         private static final double[] P_24 = { -0.99518721999702130000, -0.97472855597130950000,
  1054.             -0.93827455200273270000, -0.88641552700440100000, -0.82000198597390300000, -0.74012419157855440000,
  1055.             -0.64809365193697550000, -0.54542147138883950000, -0.43379350762604520000, -0.31504267969616340000,
  1056.             -0.19111886747361634000, -0.06405689286260563000, 0.06405689286260563000, 0.19111886747361634000,
  1057.             0.31504267969616340000, 0.43379350762604520000, 0.54542147138883950000, 0.64809365193697550000,
  1058.             0.74012419157855440000, 0.82000198597390300000, 0.88641552700440100000, 0.93827455200273270000,
  1059.             0.97472855597130950000, 0.99518721999702130000 };

  1060.         /** Weights for quadrature of order 24. */
  1061.         private static final double[] W_24 = { 0.01234122979998733500, 0.02853138862893380600, 0.04427743881741981000,
  1062.             0.05929858491543691500, 0.07334648141108027000, 0.08619016153195320000, 0.09761865210411391000,
  1063.             0.10744427011596558000, 0.11550566805372553000, 0.12167047292780335000, 0.12583745634682825000,
  1064.             0.12793819534675221000, 0.12793819534675221000, 0.12583745634682825000, 0.12167047292780335000,
  1065.             0.11550566805372553000, 0.10744427011596558000, 0.09761865210411391000, 0.08619016153195320000,
  1066.             0.07334648141108027000, 0.05929858491543691500, 0.04427743881741981000, 0.02853138862893380600,
  1067.             0.01234122979998733500 };

  1068.         /** Points for quadrature of order 32. */
  1069.         private static final double[] P_32 = { -0.99726386184948160000, -0.98561151154526840000,
  1070.             -0.96476225558750640000, -0.93490607593773970000, -0.89632115576605220000, -0.84936761373256990000,
  1071.             -0.79448379596794250000, -0.73218211874028970000, -0.66304426693021520000, -0.58771575724076230000,
  1072.             -0.50689990893222950000, -0.42135127613063540000, -0.33186860228212767000, -0.23928736225213710000,
  1073.             -0.14447196158279646000, -0.04830766568773831000, 0.04830766568773831000, 0.14447196158279646000,
  1074.             0.23928736225213710000, 0.33186860228212767000, 0.42135127613063540000, 0.50689990893222950000,
  1075.             0.58771575724076230000, 0.66304426693021520000, 0.73218211874028970000, 0.79448379596794250000,
  1076.             0.84936761373256990000, 0.89632115576605220000, 0.93490607593773970000, 0.96476225558750640000,
  1077.             0.98561151154526840000, 0.99726386184948160000 };

  1078.         /** Weights for quadrature of order 32. */
  1079.         private static final double[] W_32 = { 0.00701861000947013600, 0.01627439473090571200, 0.02539206530926214200,
  1080.             0.03427386291302141000, 0.04283589802222658600, 0.05099805926237621600, 0.05868409347853559000,
  1081.             0.06582222277636193000, 0.07234579410884862000, 0.07819389578707042000, 0.08331192422694673000,
  1082.             0.08765209300440380000, 0.09117387869576390000, 0.09384439908080441000, 0.09563872007927487000,
  1083.             0.09654008851472784000, 0.09654008851472784000, 0.09563872007927487000, 0.09384439908080441000,
  1084.             0.09117387869576390000, 0.08765209300440380000, 0.08331192422694673000, 0.07819389578707042000,
  1085.             0.07234579410884862000, 0.06582222277636193000, 0.05868409347853559000, 0.05099805926237621600,
  1086.             0.04283589802222658600, 0.03427386291302141000, 0.02539206530926214200, 0.01627439473090571200,
  1087.             0.00701861000947013600 };

  1088.         /** Points for quadrature of order 40. */
  1089.         private static final double[] P_40 = { -0.99823770971055930000, -0.99072623869945710000,
  1090.             -0.97725994998377420000, -0.95791681921379170000, -0.93281280827867660000, -0.90209880696887420000,
  1091.             -0.86595950321225960000, -0.82461223083331170000, -0.77830565142651940000, -0.72731825518992710000,
  1092.             -0.67195668461417960000, -0.61255388966798030000, -0.54946712509512820000, -0.48307580168617870000,
  1093.             -0.41377920437160500000, -0.34199409082575850000, -0.26815218500725370000, -0.19269758070137110000,
  1094.             -0.11608407067525522000, -0.03877241750605081600, 0.03877241750605081600, 0.11608407067525522000,
  1095.             0.19269758070137110000, 0.26815218500725370000, 0.34199409082575850000, 0.41377920437160500000,
  1096.             0.48307580168617870000, 0.54946712509512820000, 0.61255388966798030000, 0.67195668461417960000,
  1097.             0.72731825518992710000, 0.77830565142651940000, 0.82461223083331170000, 0.86595950321225960000,
  1098.             0.90209880696887420000, 0.93281280827867660000, 0.95791681921379170000, 0.97725994998377420000,
  1099.             0.99072623869945710000, 0.99823770971055930000 };

  1100.         /** Weights for quadrature of order 40. */
  1101.         private static final double[] W_40 = { 0.00452127709853309800, 0.01049828453115270400, 0.01642105838190797300,
  1102.             0.02224584919416689000, 0.02793700698002338000, 0.03346019528254786500, 0.03878216797447199000,
  1103.             0.04387090818567333000, 0.04869580763507221000, 0.05322784698393679000, 0.05743976909939157000,
  1104.             0.06130624249292891000, 0.06480401345660108000, 0.06791204581523394000, 0.07061164739128681000,
  1105.             0.07288658239580408000, 0.07472316905796833000, 0.07611036190062619000, 0.07703981816424793000,
  1106.             0.07750594797842482000, 0.07750594797842482000, 0.07703981816424793000, 0.07611036190062619000,
  1107.             0.07472316905796833000, 0.07288658239580408000, 0.07061164739128681000, 0.06791204581523394000,
  1108.             0.06480401345660108000, 0.06130624249292891000, 0.05743976909939157000, 0.05322784698393679000,
  1109.             0.04869580763507221000, 0.04387090818567333000, 0.03878216797447199000, 0.03346019528254786500,
  1110.             0.02793700698002338000, 0.02224584919416689000, 0.01642105838190797300, 0.01049828453115270400,
  1111.             0.00452127709853309800 };

  1112.         /** Points for quadrature of order 48. */
  1113.         private static final double[] P_48 = { -0.99877100725242610000, -0.99353017226635080000,
  1114.             -0.98412458372282700000, -0.97059159254624720000, -0.95298770316043080000, -0.93138669070655440000,
  1115.             -0.90587913671556960000, -0.87657202027424800000, -0.84358826162439350000, -0.80706620402944250000,
  1116.             -0.76715903251574020000, -0.72403413092381470000, -0.67787237963266400000, -0.62886739677651370000,
  1117.             -0.57722472608397270000, -0.52316097472223300000, -0.46690290475095840000, -0.40868648199071680000,
  1118.             -0.34875588629216070000, -0.28736248735545555000, -0.22476379039468908000, -0.16122235606889174000,
  1119.             -0.09700469920946270000, -0.03238017096286937000, 0.03238017096286937000, 0.09700469920946270000,
  1120.             0.16122235606889174000, 0.22476379039468908000, 0.28736248735545555000, 0.34875588629216070000,
  1121.             0.40868648199071680000, 0.46690290475095840000, 0.52316097472223300000, 0.57722472608397270000,
  1122.             0.62886739677651370000, 0.67787237963266400000, 0.72403413092381470000, 0.76715903251574020000,
  1123.             0.80706620402944250000, 0.84358826162439350000, 0.87657202027424800000, 0.90587913671556960000,
  1124.             0.93138669070655440000, 0.95298770316043080000, 0.97059159254624720000, 0.98412458372282700000,
  1125.             0.99353017226635080000, 0.99877100725242610000 };

  1126.         /** Weights for quadrature of order 48. */
  1127.         private static final double[] W_48 = { 0.00315334605230596250, 0.00732755390127620800, 0.01147723457923446900,
  1128.             0.01557931572294386600, 0.01961616045735556700, 0.02357076083932435600, 0.02742650970835688000,
  1129.             0.03116722783279807000, 0.03477722256477045000, 0.03824135106583080600, 0.04154508294346483000,
  1130.             0.04467456085669424000, 0.04761665849249054000, 0.05035903555385448000, 0.05289018948519365000,
  1131.             0.05519950369998416500, 0.05727729210040315000, 0.05911483969839566000, 0.06070443916589384000,
  1132.             0.06203942315989268000, 0.06311419228625403000, 0.06392423858464817000, 0.06446616443595010000,
  1133.             0.06473769681268386000, 0.06473769681268386000, 0.06446616443595010000, 0.06392423858464817000,
  1134.             0.06311419228625403000, 0.06203942315989268000, 0.06070443916589384000, 0.05911483969839566000,
  1135.             0.05727729210040315000, 0.05519950369998416500, 0.05289018948519365000, 0.05035903555385448000,
  1136.             0.04761665849249054000, 0.04467456085669424000, 0.04154508294346483000, 0.03824135106583080600,
  1137.             0.03477722256477045000, 0.03116722783279807000, 0.02742650970835688000, 0.02357076083932435600,
  1138.             0.01961616045735556700, 0.01557931572294386600, 0.01147723457923446900, 0.00732755390127620800,
  1139.             0.00315334605230596250 };

  1140.         /** Node points. */
  1141.         private final double[] nodePoints;

  1142.         /** Node weights. */
  1143.         private final double[] nodeWeights;

  1144.         /** Number of points. */
  1145.         private final int numberOfPoints;

  1146.         /**
  1147.          * Creates a Gauss integrator of the given order.
  1148.          *
  1149.          * @param numberOfPoints Order of the integration rule.
  1150.          */
  1151.         GaussQuadrature(final int numberOfPoints) {

  1152.             this.numberOfPoints = numberOfPoints;

  1153.             switch (numberOfPoints) {
  1154.                 case 12:
  1155.                     this.nodePoints = P_12.clone();
  1156.                     this.nodeWeights = W_12.clone();
  1157.                     break;
  1158.                 case 16:
  1159.                     this.nodePoints = P_16.clone();
  1160.                     this.nodeWeights = W_16.clone();
  1161.                     break;
  1162.                 case 20:
  1163.                     this.nodePoints = P_20.clone();
  1164.                     this.nodeWeights = W_20.clone();
  1165.                     break;
  1166.                 case 24:
  1167.                     this.nodePoints = P_24.clone();
  1168.                     this.nodeWeights = W_24.clone();
  1169.                     break;
  1170.                 case 32:
  1171.                     this.nodePoints = P_32.clone();
  1172.                     this.nodeWeights = W_32.clone();
  1173.                     break;
  1174.                 case 40:
  1175.                     this.nodePoints = P_40.clone();
  1176.                     this.nodeWeights = W_40.clone();
  1177.                     break;
  1178.                 case 48:
  1179.                 default:
  1180.                     this.nodePoints = P_48.clone();
  1181.                     this.nodeWeights = W_48.clone();
  1182.                     break;
  1183.             }

  1184.         }

  1185.         /**
  1186.          * Integrates a given function on the given interval.
  1187.          *
  1188.          * @param f          Function to integrate.
  1189.          * @param lowerBound Lower bound of the integration interval.
  1190.          * @param upperBound Upper bound of the integration interval.
  1191.          * @return the integral of the weighted function.
  1192.          */
  1193.         public double[] integrate(final UnivariateVectorFunction f, final double lowerBound, final double upperBound) {

  1194.             final double[] adaptedPoints = nodePoints.clone();
  1195.             final double[] adaptedWeights = nodeWeights.clone();
  1196.             transform(adaptedPoints, adaptedWeights, lowerBound, upperBound);
  1197.             return basicIntegrate(f, adaptedPoints, adaptedWeights);
  1198.         }

  1199.         /**
  1200.          * Integrates a given function on the given interval.
  1201.          *
  1202.          * @param <T>        the type of the field elements
  1203.          * @param f          Function to integrate.
  1204.          * @param lowerBound Lower bound of the integration interval.
  1205.          * @param upperBound Upper bound of the integration interval.
  1206.          * @param field      field utilized by default
  1207.          * @return the integral of the weighted function.
  1208.          */
  1209.         public <T extends CalculusFieldElement<T>> T[] integrate(final CalculusFieldUnivariateVectorFunction<T> f,
  1210.                 final T lowerBound, final T upperBound, final Field<T> field) {

  1211.             final T zero = field.getZero();

  1212.             final T[] adaptedPoints = MathArrays.buildArray(field, numberOfPoints);
  1213.             final T[] adaptedWeights = MathArrays.buildArray(field, numberOfPoints);

  1214.             for (int i = 0; i < numberOfPoints; i++) {
  1215.                 adaptedPoints[i] = zero.add(nodePoints[i]);
  1216.                 adaptedWeights[i] = zero.add(nodeWeights[i]);
  1217.             }

  1218.             transform(adaptedPoints, adaptedWeights, lowerBound, upperBound);
  1219.             return basicIntegrate(f, adaptedPoints, adaptedWeights, field);
  1220.         }

  1221.         /**
  1222.          * Performs a change of variable so that the integration can be performed on an
  1223.          * arbitrary interval {@code [a, b]}.
  1224.          * <p>
  1225.          * It is assumed that the natural interval is {@code [-1, 1]}.
  1226.          * </p>
  1227.          *
  1228.          * @param points  Points to adapt to the new interval.
  1229.          * @param weights Weights to adapt to the new interval.
  1230.          * @param a       Lower bound of the integration interval.
  1231.          * @param b       Lower bound of the integration interval.
  1232.          */
  1233.         private void transform(final double[] points, final double[] weights, final double a, final double b) {
  1234.             // Scaling
  1235.             final double scale = (b - a) / 2;
  1236.             final double shift = a + scale;
  1237.             for (int i = 0; i < points.length; i++) {
  1238.                 points[i] = points[i] * scale + shift;
  1239.                 weights[i] *= scale;
  1240.             }
  1241.         }

  1242.         /**
  1243.          * Performs a change of variable so that the integration can be performed on an
  1244.          * arbitrary interval {@code [a, b]}.
  1245.          * <p>
  1246.          * It is assumed that the natural interval is {@code [-1, 1]}.
  1247.          * </p>
  1248.          * @param <T>     the type of the field elements
  1249.          * @param points  Points to adapt to the new interval.
  1250.          * @param weights Weights to adapt to the new interval.
  1251.          * @param a       Lower bound of the integration interval.
  1252.          * @param b       Lower bound of the integration interval
  1253.          */
  1254.         private <T extends CalculusFieldElement<T>> void transform(final T[] points, final T[] weights, final T a,
  1255.                 final T b) {
  1256.             // Scaling
  1257.             final T scale = (b.subtract(a)).divide(2.);
  1258.             final T shift = a.add(scale);
  1259.             for (int i = 0; i < points.length; i++) {
  1260.                 points[i] = scale.multiply(points[i]).add(shift);
  1261.                 weights[i] = scale.multiply(weights[i]);
  1262.             }
  1263.         }

  1264.         /**
  1265.          * Returns an estimate of the integral of {@code f(x) * w(x)}, where {@code w}
  1266.          * is a weight function that depends on the actual flavor of the Gauss
  1267.          * integration scheme.
  1268.          *
  1269.          * @param f       Function to integrate.
  1270.          * @param points  Nodes.
  1271.          * @param weights Nodes weights.
  1272.          * @return the integral of the weighted function.
  1273.          */
  1274.         private double[] basicIntegrate(final UnivariateVectorFunction f, final double[] points,
  1275.                 final double[] weights) {
  1276.             double x = points[0];
  1277.             double w = weights[0];
  1278.             double[] v = f.value(x);
  1279.             final double[] y = new double[v.length];
  1280.             for (int j = 0; j < v.length; j++) {
  1281.                 y[j] = w * v[j];
  1282.             }
  1283.             final double[] t = y.clone();
  1284.             final double[] c = new double[v.length];
  1285.             final double[] s = t.clone();
  1286.             for (int i = 1; i < points.length; i++) {
  1287.                 x = points[i];
  1288.                 w = weights[i];
  1289.                 v = f.value(x);
  1290.                 for (int j = 0; j < v.length; j++) {
  1291.                     y[j] = w * v[j] - c[j];
  1292.                     t[j] = s[j] + y[j];
  1293.                     c[j] = (t[j] - s[j]) - y[j];
  1294.                     s[j] = t[j];
  1295.                 }
  1296.             }
  1297.             return s;
  1298.         }

  1299.         /**
  1300.          * Returns an estimate of the integral of {@code f(x) * w(x)}, where {@code w}
  1301.          * is a weight function that depends on the actual flavor of the Gauss
  1302.          * integration scheme.
  1303.          *
  1304.          * @param <T>     the type of the field elements.
  1305.          * @param f       Function to integrate.
  1306.          * @param points  Nodes.
  1307.          * @param weights Nodes weight
  1308.          * @param field   field utilized by default
  1309.          * @return the integral of the weighted function.
  1310.          */
  1311.         private <T extends CalculusFieldElement<T>> T[] basicIntegrate(final CalculusFieldUnivariateVectorFunction<T> f,
  1312.                 final T[] points, final T[] weights, final Field<T> field) {

  1313.             T x = points[0];
  1314.             T w = weights[0];
  1315.             T[] v = f.value(x);

  1316.             final T[] y = MathArrays.buildArray(field, v.length);
  1317.             for (int j = 0; j < v.length; j++) {
  1318.                 y[j] = v[j].multiply(w);
  1319.             }
  1320.             final T[] t = y.clone();
  1321.             final T[] c = MathArrays.buildArray(field, v.length);
  1322.             ;
  1323.             final T[] s = t.clone();
  1324.             for (int i = 1; i < points.length; i++) {
  1325.                 x = points[i];
  1326.                 w = weights[i];
  1327.                 v = f.value(x);
  1328.                 for (int j = 0; j < v.length; j++) {
  1329.                     y[j] = v[j].multiply(w).subtract(c[j]);
  1330.                     t[j] = y[j].add(s[j]);
  1331.                     c[j] = (t[j].subtract(s[j])).subtract(y[j]);
  1332.                     s[j] = t[j];
  1333.                 }
  1334.             }
  1335.             return s;
  1336.         }

  1337.     }

  1338.     /**
  1339.      * Compute the C<sub>i</sub><sup>j</sup> and the S<sub>i</sub><sup>j</sup>
  1340.      * coefficients.
  1341.      * <p>
  1342.      * Those coefficients are given in Danielson paper by expression 4.4-(6)
  1343.      * </p>
  1344.      * @author Petre Bazavan
  1345.      * @author Lucian Barbulescu
  1346.      */
  1347.     protected class FourierCjSjCoefficients {

  1348.         /** Maximum possible value for j. */
  1349.         private final int jMax;

  1350.         /**
  1351.          * The C<sub>i</sub><sup>j</sup> coefficients.
  1352.          * <p>
  1353.          * the index i corresponds to the following elements: <br/>
  1354.          * - 0 for a <br>
  1355.          * - 1 for k <br>
  1356.          * - 2 for h <br>
  1357.          * - 3 for q <br>
  1358.          * - 4 for p <br>
  1359.          * - 5 for λ <br>
  1360.          * </p>
  1361.          */
  1362.         private final double[][] cCoef;

  1363.         /**
  1364.          * The C<sub>i</sub><sup>j</sup> coefficients.
  1365.          * <p>
  1366.          * the index i corresponds to the following elements: <br/>
  1367.          * - 0 for a <br>
  1368.          * - 1 for k <br>
  1369.          * - 2 for h <br>
  1370.          * - 3 for q <br>
  1371.          * - 4 for p <br>
  1372.          * - 5 for λ <br>
  1373.          * </p>
  1374.          */
  1375.         private final double[][] sCoef;

  1376.         /**
  1377.          * Standard constructor.
  1378.          * @param state             the current state
  1379.          * @param jMax              maximum value for j
  1380.          * @param auxiliaryElements auxiliary elements related to the current orbit
  1381.          * @param parameters        values of the force model parameters
  1382.          */
  1383.         FourierCjSjCoefficients(final SpacecraftState state, final int jMax, final AuxiliaryElements auxiliaryElements,
  1384.                 final double[] parameters) {

  1385.             // Initialise the fields
  1386.             this.jMax = jMax;

  1387.             // Allocate the arrays
  1388.             final int rows = jMax + 1;
  1389.             cCoef = new double[rows][6];
  1390.             sCoef = new double[rows][6];

  1391.             // Compute the coefficients
  1392.             computeCoefficients(state, auxiliaryElements, parameters);
  1393.         }

  1394.         /**
  1395.          * Compute the Fourrier coefficients.
  1396.          * <p>
  1397.          * Only the C<sub>i</sub><sup>j</sup> and S<sub>i</sub><sup>j</sup> coefficients
  1398.          * need to be computed as D<sub>i</sub><sup>m</sup> is always 0.
  1399.          * </p>
  1400.          * @param state             the current state
  1401.          * @param auxiliaryElements auxiliary elements related to the current orbit
  1402.          * @param parameters        values of the force model parameters
  1403.          */
  1404.         private void computeCoefficients(final SpacecraftState state, final AuxiliaryElements auxiliaryElements,
  1405.                 final double[] parameters) {

  1406.             // Computes the limits for the integral
  1407.             final double[] ll = getLLimits(state, auxiliaryElements);
  1408.             // Computes integrated mean element rates if Llow < Lhigh
  1409.             if (ll[0] < ll[1]) {
  1410.                 // Compute 1 / PI
  1411.                 final double ooPI = 1 / FastMath.PI;

  1412.                 // loop through all values of j
  1413.                 for (int j = 0; j <= jMax; j++) {
  1414.                     final double[] curentCoefficients = integrator
  1415.                             .integrate(new IntegrableFunction(state, false, j, parameters), ll[0], ll[1]);

  1416.                     // divide by PI and set the values for the coefficients
  1417.                     for (int i = 0; i < 6; i++) {
  1418.                         cCoef[j][i] = ooPI * curentCoefficients[i];
  1419.                         sCoef[j][i] = ooPI * curentCoefficients[i + 6];
  1420.                     }
  1421.                 }
  1422.             }
  1423.         }

  1424.         /**
  1425.          * Get the coefficient C<sub>i</sub><sup>j</sup>.
  1426.          * @param i i index - corresponds to the required variation
  1427.          * @param j j index
  1428.          * @return the coefficient C<sub>i</sub><sup>j</sup>
  1429.          */
  1430.         public double getCij(final int i, final int j) {
  1431.             return cCoef[j][i];
  1432.         }

  1433.         /**
  1434.          * Get the coefficient S<sub>i</sub><sup>j</sup>.
  1435.          * @param i i index - corresponds to the required variation
  1436.          * @param j j index
  1437.          * @return the coefficient S<sub>i</sub><sup>j</sup>
  1438.          */
  1439.         public double getSij(final int i, final int j) {
  1440.             return sCoef[j][i];
  1441.         }
  1442.     }

  1443.     /**
  1444.      * Compute the C<sub>i</sub><sup>j</sup> and the S<sub>i</sub><sup>j</sup>
  1445.      * coefficients with field elements.
  1446.      * <p>
  1447.      * Those coefficients are given in Danielson paper by expression 4.4-(6)
  1448.      * </p>
  1449.      * @author Petre Bazavan
  1450.      * @author Lucian Barbulescu
  1451.      */
  1452.     protected class FieldFourierCjSjCoefficients<T extends CalculusFieldElement<T>> {

  1453.         /** Maximum possible value for j. */
  1454.         private final int jMax;

  1455.         /**
  1456.          * The C<sub>i</sub><sup>j</sup> coefficients.
  1457.          * <p>
  1458.          * the index i corresponds to the following elements: <br/>
  1459.          * - 0 for a <br>
  1460.          * - 1 for k <br>
  1461.          * - 2 for h <br>
  1462.          * - 3 for q <br>
  1463.          * - 4 for p <br>
  1464.          * - 5 for λ <br>
  1465.          * </p>
  1466.          */
  1467.         private final T[][] cCoef;

  1468.         /**
  1469.          * The C<sub>i</sub><sup>j</sup> coefficients.
  1470.          * <p>
  1471.          * the index i corresponds to the following elements: <br/>
  1472.          * - 0 for a <br>
  1473.          * - 1 for k <br>
  1474.          * - 2 for h <br>
  1475.          * - 3 for q <br>
  1476.          * - 4 for p <br>
  1477.          * - 5 for λ <br>
  1478.          * </p>
  1479.          */
  1480.         private final T[][] sCoef;

  1481.         /**
  1482.          * Standard constructor.
  1483.          * @param state             the current state
  1484.          * @param jMax              maximum value for j
  1485.          * @param auxiliaryElements auxiliary elements related to the current orbit
  1486.          * @param parameters        values of the force model parameters
  1487.          * @param field             field used by default
  1488.          */
  1489.         FieldFourierCjSjCoefficients(final FieldSpacecraftState<T> state, final int jMax,
  1490.                 final FieldAuxiliaryElements<T> auxiliaryElements, final T[] parameters, final Field<T> field) {
  1491.             // Initialise the fields
  1492.             this.jMax = jMax;

  1493.             // Allocate the arrays
  1494.             final int rows = jMax + 1;
  1495.             cCoef = MathArrays.buildArray(field, rows, 6);
  1496.             sCoef = MathArrays.buildArray(field, rows, 6);

  1497.             // Compute the coefficients
  1498.             computeCoefficients(state, auxiliaryElements, parameters, field);
  1499.         }

  1500.         /**
  1501.          * Compute the Fourrier coefficients.
  1502.          * <p>
  1503.          * Only the C<sub>i</sub><sup>j</sup> and S<sub>i</sub><sup>j</sup> coefficients
  1504.          * need to be computed as D<sub>i</sub><sup>m</sup> is always 0.
  1505.          * </p>
  1506.          * @param state             the current state
  1507.          * @param auxiliaryElements auxiliary elements related to the current orbit
  1508.          * @param parameters        values of the force model parameters
  1509.          * @param field             field used by default
  1510.          */
  1511.         private void computeCoefficients(final FieldSpacecraftState<T> state,
  1512.                 final FieldAuxiliaryElements<T> auxiliaryElements, final T[] parameters, final Field<T> field) {
  1513.             // Zero
  1514.             final T zero = field.getZero();
  1515.             // Computes the limits for the integral
  1516.             final T[] ll = getLLimits(state, auxiliaryElements);
  1517.             // Computes integrated mean element rates if Llow < Lhigh
  1518.             if (ll[0].getReal() < ll[1].getReal()) {
  1519.                 // Compute 1 / PI
  1520.                 final T ooPI = zero.getPi().reciprocal();

  1521.                 // loop through all values of j
  1522.                 for (int j = 0; j <= jMax; j++) {
  1523.                     final T[] curentCoefficients = integrator.integrate(
  1524.                             new FieldIntegrableFunction<>(state, false, j, parameters, field), ll[0], ll[1], field);

  1525.                     // divide by PI and set the values for the coefficients
  1526.                     for (int i = 0; i < 6; i++) {
  1527.                         cCoef[j][i] = curentCoefficients[i].multiply(ooPI);
  1528.                         sCoef[j][i] = curentCoefficients[i + 6].multiply(ooPI);
  1529.                     }
  1530.                 }
  1531.             }
  1532.         }

  1533.         /**
  1534.          * Get the coefficient C<sub>i</sub><sup>j</sup>.
  1535.          * @param i i index - corresponds to the required variation
  1536.          * @param j j index
  1537.          * @return the coefficient C<sub>i</sub><sup>j</sup>
  1538.          */
  1539.         public T getCij(final int i, final int j) {
  1540.             return cCoef[j][i];
  1541.         }

  1542.         /**
  1543.          * Get the coefficient S<sub>i</sub><sup>j</sup>.
  1544.          * @param i i index - corresponds to the required variation
  1545.          * @param j j index
  1546.          * @return the coefficient S<sub>i</sub><sup>j</sup>
  1547.          */
  1548.         public T getSij(final int i, final int j) {
  1549.             return sCoef[j][i];
  1550.         }
  1551.     }

  1552.     /**
  1553.      * This class handles the short periodic coefficients described in Danielson
  1554.      * 2.5.3-26.
  1555.      *
  1556.      * <p>
  1557.      * The value of M is 0. Also, since the values of the Fourier coefficient
  1558.      * D<sub>i</sub><sup>m</sup> is 0 then the values of the coefficients
  1559.      * D<sub>i</sub><sup>m</sup> for m &gt; 2 are also 0.
  1560.      * </p>
  1561.      * @author Petre Bazavan
  1562.      * @author Lucian Barbulescu
  1563.      *
  1564.      */
  1565.     protected static class GaussianShortPeriodicCoefficients implements ShortPeriodTerms {

  1566.         /** Maximum value for j index. */
  1567.         private final int jMax;

  1568.         /** Number of points used in the interpolation process. */
  1569.         private final int interpolationPoints;

  1570.         /** Prefix for coefficients keys. */
  1571.         private final String coefficientsKeyPrefix;

  1572.         /** All coefficients slots. */
  1573.         private final transient TimeSpanMap<Slot> slots;

  1574.         /**
  1575.          * Constructor.
  1576.          * @param coefficientsKeyPrefix prefix for coefficients keys
  1577.          * @param jMax                  maximum value for j index
  1578.          * @param interpolationPoints   number of points used in the interpolation
  1579.          *                              process
  1580.          * @param slots                 all coefficients slots
  1581.          */
  1582.         GaussianShortPeriodicCoefficients(final String coefficientsKeyPrefix, final int jMax,
  1583.                 final int interpolationPoints, final TimeSpanMap<Slot> slots) {
  1584.             // Initialize fields
  1585.             this.jMax = jMax;
  1586.             this.interpolationPoints = interpolationPoints;
  1587.             this.coefficientsKeyPrefix = coefficientsKeyPrefix;
  1588.             this.slots = slots;
  1589.         }

  1590.         /**
  1591.          * Get the slot valid for some date.
  1592.          * @param meanStates mean states defining the slot
  1593.          * @return slot valid at the specified date
  1594.          */
  1595.         public Slot createSlot(final SpacecraftState... meanStates) {
  1596.             final Slot slot = new Slot(jMax, interpolationPoints);
  1597.             final AbsoluteDate first = meanStates[0].getDate();
  1598.             final AbsoluteDate last = meanStates[meanStates.length - 1].getDate();
  1599.             final int compare = first.compareTo(last);
  1600.             if (compare < 0) {
  1601.                 slots.addValidAfter(slot, first);
  1602.             } else if (compare > 0) {
  1603.                 slots.addValidBefore(slot, first);
  1604.             } else {
  1605.                 // single date, valid for all time
  1606.                 slots.addValidAfter(slot, AbsoluteDate.PAST_INFINITY);
  1607.             }
  1608.             return slot;
  1609.         }

  1610.         /**
  1611.          * Compute the short periodic coefficients.
  1612.          *
  1613.          * @param state       current state information: date, kinematics, attitude
  1614.          * @param slot        coefficients slot
  1615.          * @param fourierCjSj Fourier coefficients
  1616.          * @param uijvij      U and V coefficients
  1617.          * @param n           Keplerian mean motion
  1618.          * @param a           semi major axis
  1619.          */
  1620.         private void computeCoefficients(final SpacecraftState state, final Slot slot,
  1621.                 final FourierCjSjCoefficients fourierCjSj, final UijVijCoefficients uijvij, final double n,
  1622.                 final double a) {

  1623.             // get the current date
  1624.             final AbsoluteDate date = state.getDate();

  1625.             // compute the k₂⁰ coefficient
  1626.             final double k20 = computeK20(jMax, uijvij.currentRhoSigmaj);

  1627.             // 1. / n
  1628.             final double oon = 1. / n;
  1629.             // 3. / (2 * a * n)
  1630.             final double to2an = 1.5 * oon / a;
  1631.             // 3. / (4 * a * n)
  1632.             final double to4an = to2an / 2;

  1633.             // Compute the coefficients for each element
  1634.             final int size = jMax + 1;
  1635.             final double[] di1 = new double[6];
  1636.             final double[] di2 = new double[6];
  1637.             final double[][] currentCij = new double[size][6];
  1638.             final double[][] currentSij = new double[size][6];
  1639.             for (int i = 0; i < 6; i++) {

  1640.                 // compute D<sub>i</sub>¹ and D<sub>i</sub>² (all others are 0)
  1641.                 di1[i] = -oon * fourierCjSj.getCij(i, 0);
  1642.                 if (i == 5) {
  1643.                     di1[i] += to2an * uijvij.getU1(0, 0);
  1644.                 }
  1645.                 di2[i] = 0.;
  1646.                 if (i == 5) {
  1647.                     di2[i] += -to4an * fourierCjSj.getCij(0, 0);
  1648.                 }

  1649.                 // the C<sub>i</sub>⁰ is computed based on all others
  1650.                 currentCij[0][i] = -di2[i] * k20;

  1651.                 for (int j = 1; j <= jMax; j++) {
  1652.                     // compute the current C<sub>i</sub><sup>j</sup> and S<sub>i</sub><sup>j</sup>
  1653.                     currentCij[j][i] = oon * uijvij.getU1(j, i);
  1654.                     if (i == 5) {
  1655.                         currentCij[j][i] += -to2an * uijvij.getU2(j);
  1656.                     }
  1657.                     currentSij[j][i] = oon * uijvij.getV1(j, i);
  1658.                     if (i == 5) {
  1659.                         currentSij[j][i] += -to2an * uijvij.getV2(j);
  1660.                     }

  1661.                     // add the computed coefficients to C<sub>i</sub>⁰
  1662.                     currentCij[0][i] += -(currentCij[j][i] * uijvij.currentRhoSigmaj[0][j] +
  1663.                             currentSij[j][i] * uijvij.currentRhoSigmaj[1][j]);
  1664.                 }

  1665.             }

  1666.             // add the values to the interpolators
  1667.             slot.cij[0].addGridPoint(date, currentCij[0]);
  1668.             slot.dij[1].addGridPoint(date, di1);
  1669.             slot.dij[2].addGridPoint(date, di2);
  1670.             for (int j = 1; j <= jMax; j++) {
  1671.                 slot.cij[j].addGridPoint(date, currentCij[j]);
  1672.                 slot.sij[j].addGridPoint(date, currentSij[j]);
  1673.             }

  1674.         }

  1675.         /**
  1676.          * Compute the coefficient k₂⁰ by using the equation 2.5.3-(9a) from Danielson.
  1677.          * <p>
  1678.          * After inserting 2.5.3-(8) into 2.5.3-(9a) the result becomes:<br>
  1679.          * k₂⁰ = &Sigma;<sub>k=1</sub><sup>kMax</sup>[(2 / k²) * (σ<sub>k</sub>² +
  1680.          * ρ<sub>k</sub>²)]
  1681.          * </p>
  1682.          * @param kMax             max value fot k index
  1683.          * @param currentRhoSigmaj the current computed values for the ρ<sub>j</sub> and
  1684.          *                         σ<sub>j</sub> coefficients
  1685.          * @return the coefficient k₂⁰
  1686.          */
  1687.         private double computeK20(final int kMax, final double[][] currentRhoSigmaj) {
  1688.             double k20 = 0.;

  1689.             for (int kIndex = 1; kIndex <= kMax; kIndex++) {
  1690.                 // After inserting 2.5.3-(8) into 2.5.3-(9a) the result becomes:
  1691.                 // k₂⁰ = &Sigma;<sub>k=1</sub><sup>kMax</sup>[(2 / k²) * (σ<sub>k</sub>² +
  1692.                 // ρ<sub>k</sub>²)]
  1693.                 double currentTerm = currentRhoSigmaj[1][kIndex] * currentRhoSigmaj[1][kIndex] +
  1694.                         currentRhoSigmaj[0][kIndex] * currentRhoSigmaj[0][kIndex];

  1695.                 // multiply by 2 / k²
  1696.                 currentTerm *= 2. / (kIndex * kIndex);

  1697.                 // add the term to the result
  1698.                 k20 += currentTerm;
  1699.             }

  1700.             return k20;
  1701.         }

  1702.         /** {@inheritDoc} */
  1703.         @Override
  1704.         public double[] value(final Orbit meanOrbit) {

  1705.             // select the coefficients slot
  1706.             final Slot slot = slots.get(meanOrbit.getDate());

  1707.             // Get the True longitude L
  1708.             final double L = meanOrbit.getLv();

  1709.             // Compute the center (l - λ)
  1710.             final double center = L - meanOrbit.getLM();
  1711.             // Compute (l - λ)²
  1712.             final double center2 = center * center;

  1713.             // Initialize short periodic variations
  1714.             final double[] shortPeriodicVariation = slot.cij[0].value(meanOrbit.getDate());
  1715.             final double[] d1 = slot.dij[1].value(meanOrbit.getDate());
  1716.             final double[] d2 = slot.dij[2].value(meanOrbit.getDate());
  1717.             for (int i = 0; i < 6; i++) {
  1718.                 shortPeriodicVariation[i] += center * d1[i] + center2 * d2[i];
  1719.             }

  1720.             for (int j = 1; j <= JMAX; j++) {
  1721.                 final double[] c = slot.cij[j].value(meanOrbit.getDate());
  1722.                 final double[] s = slot.sij[j].value(meanOrbit.getDate());
  1723.                 final SinCos sc  = FastMath.sinCos(j * L);
  1724.                 final double cos = sc.cos();
  1725.                 final double sin = sc.sin();
  1726.                 for (int i = 0; i < 6; i++) {
  1727.                     // add corresponding term to the short periodic variation
  1728.                     shortPeriodicVariation[i] += c[i] * cos;
  1729.                     shortPeriodicVariation[i] += s[i] * sin;
  1730.                 }
  1731.             }

  1732.             return shortPeriodicVariation;

  1733.         }

  1734.         /** {@inheritDoc} */
  1735.         public String getCoefficientsKeyPrefix() {
  1736.             return coefficientsKeyPrefix;
  1737.         }

  1738.         /**
  1739.          * {@inheritDoc}
  1740.          * <p>
  1741.          * For Gaussian forces, there are JMAX cj coefficients, JMAX sj coefficients and
  1742.          * 3 dj coefficients. As JMAX = 12, this sums up to 27 coefficients. The j index
  1743.          * is the integer multiplier for the true longitude argument in the cj and sj
  1744.          * coefficients and to the degree in the polynomial dj coefficients.
  1745.          * </p>
  1746.          */
  1747.         @Override
  1748.         public Map<String, double[]> getCoefficients(final AbsoluteDate date, final Set<String> selected) {

  1749.             // select the coefficients slot
  1750.             final Slot slot = slots.get(date);

  1751.             final Map<String, double[]> coefficients = new HashMap<String, double[]>(2 * JMAX + 3);
  1752.             storeIfSelected(coefficients, selected, slot.cij[0].value(date), "d", 0);
  1753.             storeIfSelected(coefficients, selected, slot.dij[1].value(date), "d", 1);
  1754.             storeIfSelected(coefficients, selected, slot.dij[2].value(date), "d", 2);
  1755.             for (int j = 1; j <= JMAX; j++) {
  1756.                 storeIfSelected(coefficients, selected, slot.cij[j].value(date), "c", j);
  1757.                 storeIfSelected(coefficients, selected, slot.sij[j].value(date), "s", j);
  1758.             }

  1759.             return coefficients;

  1760.         }

  1761.         /**
  1762.          * Put a coefficient in a map if selected.
  1763.          * @param map      map to populate
  1764.          * @param selected set of coefficients that should be put in the map (empty set
  1765.          *                 means all coefficients are selected)
  1766.          * @param value    coefficient value
  1767.          * @param id       coefficient identifier
  1768.          * @param indices  list of coefficient indices
  1769.          */
  1770.         private void storeIfSelected(final Map<String, double[]> map, final Set<String> selected, final double[] value,
  1771.                 final String id, final int... indices) {
  1772.             final StringBuilder keyBuilder = new StringBuilder(getCoefficientsKeyPrefix());
  1773.             keyBuilder.append(id);
  1774.             for (int index : indices) {
  1775.                 keyBuilder.append('[').append(index).append(']');
  1776.             }
  1777.             final String key = keyBuilder.toString();
  1778.             if (selected.isEmpty() || selected.contains(key)) {
  1779.                 map.put(key, value);
  1780.             }
  1781.         }

  1782.     }

  1783.     /**
  1784.      * This class handles the short periodic coefficients described in Danielson
  1785.      * 2.5.3-26.
  1786.      *
  1787.      * <p>
  1788.      * The value of M is 0. Also, since the values of the Fourier coefficient
  1789.      * D<sub>i</sub><sup>m</sup> is 0 then the values of the coefficients
  1790.      * D<sub>i</sub><sup>m</sup> for m &gt; 2 are also 0.
  1791.      * </p>
  1792.      * @author Petre Bazavan
  1793.      * @author Lucian Barbulescu
  1794.      *
  1795.      */
  1796.     protected static class FieldGaussianShortPeriodicCoefficients<T extends CalculusFieldElement<T>>
  1797.             implements FieldShortPeriodTerms<T> {

  1798.         /** Maximum value for j index. */
  1799.         private final int jMax;

  1800.         /** Number of points used in the interpolation process. */
  1801.         private final int interpolationPoints;

  1802.         /** Prefix for coefficients keys. */
  1803.         private final String coefficientsKeyPrefix;

  1804.         /** All coefficients slots. */
  1805.         private final transient FieldTimeSpanMap<FieldSlot<T>, T> slots;

  1806.         /**
  1807.          * Constructor.
  1808.          * @param coefficientsKeyPrefix prefix for coefficients keys
  1809.          * @param jMax                  maximum value for j index
  1810.          * @param interpolationPoints   number of points used in the interpolation
  1811.          *                              process
  1812.          * @param slots                 all coefficients slots
  1813.          */
  1814.         FieldGaussianShortPeriodicCoefficients(final String coefficientsKeyPrefix, final int jMax,
  1815.                 final int interpolationPoints, final FieldTimeSpanMap<FieldSlot<T>, T> slots) {
  1816.             // Initialize fields
  1817.             this.jMax = jMax;
  1818.             this.interpolationPoints = interpolationPoints;
  1819.             this.coefficientsKeyPrefix = coefficientsKeyPrefix;
  1820.             this.slots = slots;
  1821.         }

  1822.         /**
  1823.          * Get the slot valid for some date.
  1824.          * @param meanStates mean states defining the slot
  1825.          * @return slot valid at the specified date
  1826.          */
  1827.         @SuppressWarnings("unchecked")
  1828.         public FieldSlot<T> createSlot(final FieldSpacecraftState<T>... meanStates) {
  1829.             final FieldSlot<T> slot = new FieldSlot<>(jMax, interpolationPoints);
  1830.             final FieldAbsoluteDate<T> first = meanStates[0].getDate();
  1831.             final FieldAbsoluteDate<T> last = meanStates[meanStates.length - 1].getDate();
  1832.             if (first.compareTo(last) <= 0) {
  1833.                 slots.addValidAfter(slot, first);
  1834.             } else {
  1835.                 slots.addValidBefore(slot, first);
  1836.             }
  1837.             return slot;
  1838.         }

  1839.         /**
  1840.          * Compute the short periodic coefficients.
  1841.          *
  1842.          * @param state       current state information: date, kinematics, attitude
  1843.          * @param slot        coefficients slot
  1844.          * @param fourierCjSj Fourier coefficients
  1845.          * @param uijvij      U and V coefficients
  1846.          * @param n           Keplerian mean motion
  1847.          * @param a           semi major axis
  1848.          * @param field       field used by default
  1849.          */
  1850.         private void computeCoefficients(final FieldSpacecraftState<T> state, final FieldSlot<T> slot,
  1851.                 final FieldFourierCjSjCoefficients<T> fourierCjSj, final FieldUijVijCoefficients<T> uijvij, final T n,
  1852.                 final T a, final Field<T> field) {

  1853.             // Zero
  1854.             final T zero = field.getZero();

  1855.             // get the current date
  1856.             final FieldAbsoluteDate<T> date = state.getDate();

  1857.             // compute the k₂⁰ coefficient
  1858.             final T k20 = computeK20(jMax, uijvij.currentRhoSigmaj, field);

  1859.             // 1. / n
  1860.             final T oon = n.reciprocal();
  1861.             // 3. / (2 * a * n)
  1862.             final T to2an = oon.multiply(1.5).divide(a);
  1863.             // 3. / (4 * a * n)
  1864.             final T to4an = to2an.divide(2.);

  1865.             // Compute the coefficients for each element
  1866.             final int size = jMax + 1;
  1867.             final T[] di1 = MathArrays.buildArray(field, 6);
  1868.             final T[] di2 = MathArrays.buildArray(field, 6);
  1869.             final T[][] currentCij = MathArrays.buildArray(field, size, 6);
  1870.             final T[][] currentSij = MathArrays.buildArray(field, size, 6);
  1871.             for (int i = 0; i < 6; i++) {

  1872.                 // compute D<sub>i</sub>¹ and D<sub>i</sub>² (all others are 0)
  1873.                 di1[i] = oon.negate().multiply(fourierCjSj.getCij(i, 0));
  1874.                 if (i == 5) {
  1875.                     di1[i] = di1[i].add(to2an.multiply(uijvij.getU1(0, 0)));
  1876.                 }
  1877.                 di2[i] = zero;
  1878.                 if (i == 5) {
  1879.                     di2[i] = di2[i].add(to4an.negate().multiply(fourierCjSj.getCij(0, 0)));
  1880.                 }

  1881.                 // the C<sub>i</sub>⁰ is computed based on all others
  1882.                 currentCij[0][i] = di2[i].negate().multiply(k20);

  1883.                 for (int j = 1; j <= jMax; j++) {
  1884.                     // compute the current C<sub>i</sub><sup>j</sup> and S<sub>i</sub><sup>j</sup>
  1885.                     currentCij[j][i] = oon.multiply(uijvij.getU1(j, i));
  1886.                     if (i == 5) {
  1887.                         currentCij[j][i] = currentCij[j][i].add(to2an.negate().multiply(uijvij.getU2(j)));
  1888.                     }
  1889.                     currentSij[j][i] = oon.multiply(uijvij.getV1(j, i));
  1890.                     if (i == 5) {
  1891.                         currentSij[j][i] = currentSij[j][i].add(to2an.negate().multiply(uijvij.getV2(j)));
  1892.                     }

  1893.                     // add the computed coefficients to C<sub>i</sub>⁰
  1894.                     currentCij[0][i] = currentCij[0][i].add(currentCij[j][i].multiply(uijvij.currentRhoSigmaj[0][j])
  1895.                             .add(currentSij[j][i].multiply(uijvij.currentRhoSigmaj[1][j])).negate());
  1896.                 }

  1897.             }

  1898.             // add the values to the interpolators
  1899.             slot.cij[0].addGridPoint(date, currentCij[0]);
  1900.             slot.dij[1].addGridPoint(date, di1);
  1901.             slot.dij[2].addGridPoint(date, di2);
  1902.             for (int j = 1; j <= jMax; j++) {
  1903.                 slot.cij[j].addGridPoint(date, currentCij[j]);
  1904.                 slot.sij[j].addGridPoint(date, currentSij[j]);
  1905.             }

  1906.         }

  1907.         /**
  1908.          * Compute the coefficient k₂⁰ by using the equation 2.5.3-(9a) from Danielson.
  1909.          * <p>
  1910.          * After inserting 2.5.3-(8) into 2.5.3-(9a) the result becomes:<br>
  1911.          * k₂⁰ = &Sigma;<sub>k=1</sub><sup>kMax</sup>[(2 / k²) * (σ<sub>k</sub>² +
  1912.          * ρ<sub>k</sub>²)]
  1913.          * </p>
  1914.          * @param kMax             max value fot k index
  1915.          * @param currentRhoSigmaj the current computed values for the ρ<sub>j</sub> and
  1916.          *                         σ<sub>j</sub> coefficients
  1917.          * @param field            field used by default
  1918.          * @return the coefficient k₂⁰
  1919.          */
  1920.         private T computeK20(final int kMax, final T[][] currentRhoSigmaj, final Field<T> field) {
  1921.             final T zero = field.getZero();
  1922.             T k20 = zero;

  1923.             for (int kIndex = 1; kIndex <= kMax; kIndex++) {
  1924.                 // After inserting 2.5.3-(8) into 2.5.3-(9a) the result becomes:
  1925.                 // k₂⁰ = &Sigma;<sub>k=1</sub><sup>kMax</sup>[(2 / k²) * (σ<sub>k</sub>² +
  1926.                 // ρ<sub>k</sub>²)]
  1927.                 T currentTerm = currentRhoSigmaj[1][kIndex].multiply(currentRhoSigmaj[1][kIndex])
  1928.                         .add(currentRhoSigmaj[0][kIndex].multiply(currentRhoSigmaj[0][kIndex]));

  1929.                 // multiply by 2 / k²
  1930.                 currentTerm = currentTerm.multiply(2. / (kIndex * kIndex));

  1931.                 // add the term to the result
  1932.                 k20 = k20.add(currentTerm);
  1933.             }

  1934.             return k20;
  1935.         }

  1936.         /** {@inheritDoc} */
  1937.         @Override
  1938.         public T[] value(final FieldOrbit<T> meanOrbit) {

  1939.             // select the coefficients slot
  1940.             final FieldSlot<T> slot = slots.get(meanOrbit.getDate());

  1941.             // Get the True longitude L
  1942.             final T L = meanOrbit.getLv();

  1943.             // Compute the center (l - λ)
  1944.             final T center = L.subtract(meanOrbit.getLM());
  1945.             // Compute (l - λ)²
  1946.             final T center2 = center.multiply(center);

  1947.             // Initialize short periodic variations
  1948.             final T[] shortPeriodicVariation = slot.cij[0].value(meanOrbit.getDate());
  1949.             final T[] d1 = slot.dij[1].value(meanOrbit.getDate());
  1950.             final T[] d2 = slot.dij[2].value(meanOrbit.getDate());
  1951.             for (int i = 0; i < 6; i++) {
  1952.                 shortPeriodicVariation[i] = shortPeriodicVariation[i]
  1953.                         .add(center.multiply(d1[i]).add(center2.multiply(d2[i])));
  1954.             }

  1955.             for (int j = 1; j <= JMAX; j++) {
  1956.                 final T[] c = slot.cij[j].value(meanOrbit.getDate());
  1957.                 final T[] s = slot.sij[j].value(meanOrbit.getDate());
  1958.                 final FieldSinCos<T> sc = FastMath.sinCos(L.multiply(j));
  1959.                 final T cos = sc.cos();
  1960.                 final T sin = sc.sin();
  1961.                 for (int i = 0; i < 6; i++) {
  1962.                     // add corresponding term to the short periodic variation
  1963.                     shortPeriodicVariation[i] = shortPeriodicVariation[i].add(c[i].multiply(cos));
  1964.                     shortPeriodicVariation[i] = shortPeriodicVariation[i].add(s[i].multiply(sin));
  1965.                 }
  1966.             }

  1967.             return shortPeriodicVariation;

  1968.         }

  1969.         /** {@inheritDoc} */
  1970.         public String getCoefficientsKeyPrefix() {
  1971.             return coefficientsKeyPrefix;
  1972.         }

  1973.         /**
  1974.          * {@inheritDoc}
  1975.          * <p>
  1976.          * For Gaussian forces, there are JMAX cj coefficients, JMAX sj coefficients and
  1977.          * 3 dj coefficients. As JMAX = 12, this sums up to 27 coefficients. The j index
  1978.          * is the integer multiplier for the true longitude argument in the cj and sj
  1979.          * coefficients and to the degree in the polynomial dj coefficients.
  1980.          * </p>
  1981.          */
  1982.         @Override
  1983.         public Map<String, T[]> getCoefficients(final FieldAbsoluteDate<T> date, final Set<String> selected) {

  1984.             // select the coefficients slot
  1985.             final FieldSlot<T> slot = slots.get(date);

  1986.             final Map<String, T[]> coefficients = new HashMap<String, T[]>(2 * JMAX + 3);
  1987.             storeIfSelected(coefficients, selected, slot.cij[0].value(date), "d", 0);
  1988.             storeIfSelected(coefficients, selected, slot.dij[1].value(date), "d", 1);
  1989.             storeIfSelected(coefficients, selected, slot.dij[2].value(date), "d", 2);
  1990.             for (int j = 1; j <= JMAX; j++) {
  1991.                 storeIfSelected(coefficients, selected, slot.cij[j].value(date), "c", j);
  1992.                 storeIfSelected(coefficients, selected, slot.sij[j].value(date), "s", j);
  1993.             }

  1994.             return coefficients;

  1995.         }

  1996.         /**
  1997.          * Put a coefficient in a map if selected.
  1998.          * @param map      map to populate
  1999.          * @param selected set of coefficients that should be put in the map (empty set
  2000.          *                 means all coefficients are selected)
  2001.          * @param value    coefficient value
  2002.          * @param id       coefficient identifier
  2003.          * @param indices  list of coefficient indices
  2004.          */
  2005.         private void storeIfSelected(final Map<String, T[]> map, final Set<String> selected, final T[] value,
  2006.                 final String id, final int... indices) {
  2007.             final StringBuilder keyBuilder = new StringBuilder(getCoefficientsKeyPrefix());
  2008.             keyBuilder.append(id);
  2009.             for (int index : indices) {
  2010.                 keyBuilder.append('[').append(index).append(']');
  2011.             }
  2012.             final String key = keyBuilder.toString();
  2013.             if (selected.isEmpty() || selected.contains(key)) {
  2014.                 map.put(key, value);
  2015.             }
  2016.         }

  2017.     }

  2018.     /**
  2019.      * The U<sub>i</sub><sup>j</sup> and V<sub>i</sub><sup>j</sup> coefficients
  2020.      * described by equations 2.5.3-(21) and 2.5.3-(22) from Danielson.
  2021.      * <p>
  2022.      * The index i takes only the values 1 and 2<br>
  2023.      * For U only the index 0 for j is used.
  2024.      * </p>
  2025.      *
  2026.      * @author Petre Bazavan
  2027.      * @author Lucian Barbulescu
  2028.      */
  2029.     protected static class UijVijCoefficients {

  2030.         /**
  2031.          * The U₁<sup>j</sup> coefficients.
  2032.          * <p>
  2033.          * The first index identifies the Fourier coefficients used<br>
  2034.          * Those coefficients are computed for all Fourier C<sub>i</sub><sup>j</sup> and
  2035.          * S<sub>i</sub><sup>j</sup><br>
  2036.          * The only exception is when j = 0 when only the coefficient for fourier index
  2037.          * = 1 (i == 0) is needed.<br>
  2038.          * Also, for fourier index = 1 (i == 0), the coefficients up to 2 * jMax are
  2039.          * computed, because are required to compute the coefficients U₂<sup>j</sup>
  2040.          * </p>
  2041.          */
  2042.         private final double[][] u1ij;

  2043.         /**
  2044.          * The V₁<sup>j</sup> coefficients.
  2045.          * <p>
  2046.          * The first index identifies the Fourier coefficients used<br>
  2047.          * Those coefficients are computed for all Fourier C<sub>i</sub><sup>j</sup> and
  2048.          * S<sub>i</sub><sup>j</sup><br>
  2049.          * for fourier index = 1 (i == 0), the coefficients up to 2 * jMax are computed,
  2050.          * because are required to compute the coefficients V₂<sup>j</sup>
  2051.          * </p>
  2052.          */
  2053.         private final double[][] v1ij;

  2054.         /**
  2055.          * The U₂<sup>j</sup> coefficients.
  2056.          * <p>
  2057.          * Only the coefficients that use the Fourier index = 1 (i == 0) are computed as
  2058.          * they are the only ones required.
  2059.          * </p>
  2060.          */
  2061.         private final double[] u2ij;

  2062.         /**
  2063.          * The V₂<sup>j</sup> coefficients.
  2064.          * <p>
  2065.          * Only the coefficients that use the Fourier index = 1 (i == 0) are computed as
  2066.          * they are the only ones required.
  2067.          * </p>
  2068.          */
  2069.         private final double[] v2ij;

  2070.         /**
  2071.          * The current computed values for the ρ<sub>j</sub> and σ<sub>j</sub>
  2072.          * coefficients.
  2073.          */
  2074.         private final double[][] currentRhoSigmaj;

  2075.         /**
  2076.          * The C<sub>i</sub><sup>j</sup> and the S<sub>i</sub><sup>j</sup> Fourier
  2077.          * coefficients.
  2078.          */
  2079.         private final FourierCjSjCoefficients fourierCjSj;

  2080.         /** The maximum value for j index. */
  2081.         private final int jMax;

  2082.         /**
  2083.          * Constructor.
  2084.          * @param currentRhoSigmaj the current computed values for the ρ<sub>j</sub> and
  2085.          *                         σ<sub>j</sub> coefficients
  2086.          * @param fourierCjSj      the fourier coefficients C<sub>i</sub><sup>j</sup>
  2087.          *                         and the S<sub>i</sub><sup>j</sup>
  2088.          * @param jMax             maximum value for j index
  2089.          */
  2090.         UijVijCoefficients(final double[][] currentRhoSigmaj, final FourierCjSjCoefficients fourierCjSj,
  2091.                 final int jMax) {
  2092.             this.currentRhoSigmaj = currentRhoSigmaj;
  2093.             this.fourierCjSj = fourierCjSj;
  2094.             this.jMax = jMax;

  2095.             // initialize the internal arrays.
  2096.             this.u1ij = new double[6][2 * jMax + 1];
  2097.             this.v1ij = new double[6][2 * jMax + 1];
  2098.             this.u2ij = new double[jMax + 1];
  2099.             this.v2ij = new double[jMax + 1];

  2100.             // compute the coefficients
  2101.             computeU1V1Coefficients();
  2102.             computeU2V2Coefficients();
  2103.         }

  2104.         /** Build the U₁<sup>j</sup> and V₁<sup>j</sup> coefficients. */
  2105.         private void computeU1V1Coefficients() {
  2106.             // generate the U₁<sup>j</sup> and V₁<sup>j</sup> coefficients
  2107.             // for j >= 1
  2108.             // also the U₁⁰ for Fourier index = 1 (i == 0) coefficient will be computed
  2109.             u1ij[0][0] = 0;
  2110.             for (int j = 1; j <= jMax; j++) {
  2111.                 // compute 1 / j
  2112.                 final double ooj = 1. / j;

  2113.                 for (int i = 0; i < 6; i++) {
  2114.                     // j is aready between 1 and J
  2115.                     u1ij[i][j] = fourierCjSj.getSij(i, j);
  2116.                     v1ij[i][j] = fourierCjSj.getCij(i, j);

  2117.                     // 1 - δ<sub>1j</sub> is 1 for all j > 1
  2118.                     if (j > 1) {
  2119.                         // k starts with 1 because j-J is less than or equal to 0
  2120.                         for (int kIndex = 1; kIndex <= j - 1; kIndex++) {
  2121.                             // C<sub>i</sub><sup>j-k</sup> * σ<sub>k</sub> +
  2122.                             // S<sub>i</sub><sup>j-k</sup> * ρ<sub>k</sub>
  2123.                             u1ij[i][j] += fourierCjSj.getCij(i, j - kIndex) * currentRhoSigmaj[1][kIndex] +
  2124.                                     fourierCjSj.getSij(i, j - kIndex) * currentRhoSigmaj[0][kIndex];

  2125.                             // C<sub>i</sub><sup>j-k</sup> * ρ<sub>k</sub> -
  2126.                             // S<sub>i</sub><sup>j-k</sup> * σ<sub>k</sub>
  2127.                             v1ij[i][j] += fourierCjSj.getCij(i, j - kIndex) * currentRhoSigmaj[0][kIndex] -
  2128.                                     fourierCjSj.getSij(i, j - kIndex) * currentRhoSigmaj[1][kIndex];
  2129.                         }
  2130.                     }

  2131.                     // since j must be between 1 and J-1 and is already between 1 and J
  2132.                     // the following sum is skiped only for j = jMax
  2133.                     if (j != jMax) {
  2134.                         for (int kIndex = 1; kIndex <= jMax - j; kIndex++) {
  2135.                             // -C<sub>i</sub><sup>j+k</sup> * σ<sub>k</sub> +
  2136.                             // S<sub>i</sub><sup>j+k</sup> * ρ<sub>k</sub>
  2137.                             u1ij[i][j] += -fourierCjSj.getCij(i, j + kIndex) * currentRhoSigmaj[1][kIndex] +
  2138.                                     fourierCjSj.getSij(i, j + kIndex) * currentRhoSigmaj[0][kIndex];

  2139.                             // C<sub>i</sub><sup>j+k</sup> * ρ<sub>k</sub> +
  2140.                             // S<sub>i</sub><sup>j+k</sup> * σ<sub>k</sub>
  2141.                             v1ij[i][j] += fourierCjSj.getCij(i, j + kIndex) * currentRhoSigmaj[0][kIndex] +
  2142.                                     fourierCjSj.getSij(i, j + kIndex) * currentRhoSigmaj[1][kIndex];
  2143.                         }
  2144.                     }

  2145.                     for (int kIndex = 1; kIndex <= jMax; kIndex++) {
  2146.                         // C<sub>i</sub><sup>k</sup> * σ<sub>j+k</sub> -
  2147.                         // S<sub>i</sub><sup>k</sup> * ρ<sub>j+k</sub>
  2148.                         u1ij[i][j] += -fourierCjSj.getCij(i, kIndex) * currentRhoSigmaj[1][j + kIndex] -
  2149.                                 fourierCjSj.getSij(i, kIndex) * currentRhoSigmaj[0][j + kIndex];

  2150.                         // C<sub>i</sub><sup>k</sup> * ρ<sub>j+k</sub> +
  2151.                         // S<sub>i</sub><sup>k</sup> * σ<sub>j+k</sub>
  2152.                         v1ij[i][j] += fourierCjSj.getCij(i, kIndex) * currentRhoSigmaj[0][j + kIndex] +
  2153.                                 fourierCjSj.getSij(i, kIndex) * currentRhoSigmaj[1][j + kIndex];
  2154.                     }

  2155.                     // divide by 1 / j
  2156.                     u1ij[i][j] *= -ooj;
  2157.                     v1ij[i][j] *= ooj;

  2158.                     // if index = 1 (i == 0) add the computed terms to U₁⁰
  2159.                     if (i == 0) {
  2160.                         // - (U₁<sup>j</sup> * ρ<sub>j</sub> + V₁<sup>j</sup> * σ<sub>j</sub>
  2161.                         u1ij[0][0] += -u1ij[0][j] * currentRhoSigmaj[0][j] - v1ij[0][j] * currentRhoSigmaj[1][j];
  2162.                     }
  2163.                 }
  2164.             }

  2165.             // Terms with j > jMax are required only when computing the coefficients
  2166.             // U₂<sup>j</sup> and V₂<sup>j</sup>
  2167.             // and those coefficients are only required for Fourier index = 1 (i == 0).
  2168.             for (int j = jMax + 1; j <= 2 * jMax; j++) {
  2169.                 // compute 1 / j
  2170.                 final double ooj = 1. / j;
  2171.                 // the value of i is 0
  2172.                 u1ij[0][j] = 0.;
  2173.                 v1ij[0][j] = 0.;

  2174.                 // k starts from j-J as it is always greater than or equal to 1
  2175.                 for (int kIndex = j - jMax; kIndex <= j - 1; kIndex++) {
  2176.                     // C<sub>i</sub><sup>j-k</sup> * σ<sub>k</sub> +
  2177.                     // S<sub>i</sub><sup>j-k</sup> * ρ<sub>k</sub>
  2178.                     u1ij[0][j] += fourierCjSj.getCij(0, j - kIndex) * currentRhoSigmaj[1][kIndex] +
  2179.                             fourierCjSj.getSij(0, j - kIndex) * currentRhoSigmaj[0][kIndex];

  2180.                     // C<sub>i</sub><sup>j-k</sup> * ρ<sub>k</sub> -
  2181.                     // S<sub>i</sub><sup>j-k</sup> * σ<sub>k</sub>
  2182.                     v1ij[0][j] += fourierCjSj.getCij(0, j - kIndex) * currentRhoSigmaj[0][kIndex] -
  2183.                             fourierCjSj.getSij(0, j - kIndex) * currentRhoSigmaj[1][kIndex];
  2184.                 }
  2185.                 for (int kIndex = 1; kIndex <= jMax; kIndex++) {
  2186.                     // C<sub>i</sub><sup>k</sup> * σ<sub>j+k</sub> -
  2187.                     // S<sub>i</sub><sup>k</sup> * ρ<sub>j+k</sub>
  2188.                     u1ij[0][j] += -fourierCjSj.getCij(0, kIndex) * currentRhoSigmaj[1][j + kIndex] -
  2189.                             fourierCjSj.getSij(0, kIndex) * currentRhoSigmaj[0][j + kIndex];

  2190.                     // C<sub>i</sub><sup>k</sup> * ρ<sub>j+k</sub> +
  2191.                     // S<sub>i</sub><sup>k</sup> * σ<sub>j+k</sub>
  2192.                     v1ij[0][j] += fourierCjSj.getCij(0, kIndex) * currentRhoSigmaj[0][j + kIndex] +
  2193.                             fourierCjSj.getSij(0, kIndex) * currentRhoSigmaj[1][j + kIndex];
  2194.                 }

  2195.                 // divide by 1 / j
  2196.                 u1ij[0][j] *= -ooj;
  2197.                 v1ij[0][j] *= ooj;
  2198.             }
  2199.         }

  2200.         /**
  2201.          * Build the U₁<sup>j</sup> and V₁<sup>j</sup> coefficients.
  2202.          * <p>
  2203.          * Only the coefficients for Fourier index = 1 (i == 0) are required.
  2204.          * </p>
  2205.          */
  2206.         private void computeU2V2Coefficients() {
  2207.             for (int j = 1; j <= jMax; j++) {
  2208.                 // compute 1 / j
  2209.                 final double ooj = 1. / j;

  2210.                 // only the values for i == 0 are computed
  2211.                 u2ij[j] = v1ij[0][j];
  2212.                 v2ij[j] = u1ij[0][j];

  2213.                 // 1 - δ<sub>1j</sub> is 1 for all j > 1
  2214.                 if (j > 1) {
  2215.                     for (int l = 1; l <= j - 1; l++) {
  2216.                         // U₁<sup>j-l</sup> * σ<sub>l</sub> +
  2217.                         // V₁<sup>j-l</sup> * ρ<sub>l</sub>
  2218.                         u2ij[j] += u1ij[0][j - l] * currentRhoSigmaj[1][l] + v1ij[0][j - l] * currentRhoSigmaj[0][l];

  2219.                         // U₁<sup>j-l</sup> * ρ<sub>l</sub> -
  2220.                         // V₁<sup>j-l</sup> * σ<sub>l</sub>
  2221.                         v2ij[j] += u1ij[0][j - l] * currentRhoSigmaj[0][l] - v1ij[0][j - l] * currentRhoSigmaj[1][l];
  2222.                     }
  2223.                 }

  2224.                 for (int l = 1; l <= jMax; l++) {
  2225.                     // -U₁<sup>j+l</sup> * σ<sub>l</sub> +
  2226.                     // U₁<sup>l</sup> * σ<sub>j+l</sub> +
  2227.                     // V₁<sup>j+l</sup> * ρ<sub>l</sub> -
  2228.                     // V₁<sup>l</sup> * ρ<sub>j+l</sub>
  2229.                     u2ij[j] += -u1ij[0][j + l] * currentRhoSigmaj[1][l] + u1ij[0][l] * currentRhoSigmaj[1][j + l] +
  2230.                             v1ij[0][j + l] * currentRhoSigmaj[0][l] - v1ij[0][l] * currentRhoSigmaj[0][j + l];

  2231.                     // U₁<sup>j+l</sup> * ρ<sub>l</sub> +
  2232.                     // U₁<sup>l</sup> * ρ<sub>j+l</sub> +
  2233.                     // V₁<sup>j+l</sup> * σ<sub>l</sub> +
  2234.                     // V₁<sup>l</sup> * σ<sub>j+l</sub>
  2235.                     u2ij[j] += u1ij[0][j + l] * currentRhoSigmaj[0][l] + u1ij[0][l] * currentRhoSigmaj[0][j + l] +
  2236.                             v1ij[0][j + l] * currentRhoSigmaj[1][l] + v1ij[0][l] * currentRhoSigmaj[1][j + l];
  2237.                 }

  2238.                 // divide by 1 / j
  2239.                 u2ij[j] *= -ooj;
  2240.                 v2ij[j] *= ooj;
  2241.             }
  2242.         }

  2243.         /**
  2244.          * Get the coefficient U₁<sup>j</sup> for Fourier index i.
  2245.          *
  2246.          * @param j j index
  2247.          * @param i Fourier index (starts at 0)
  2248.          * @return the coefficient U₁<sup>j</sup> for the given Fourier index i
  2249.          */
  2250.         public double getU1(final int j, final int i) {
  2251.             return u1ij[i][j];
  2252.         }

  2253.         /**
  2254.          * Get the coefficient V₁<sup>j</sup> for Fourier index i.
  2255.          *
  2256.          * @param j j index
  2257.          * @param i Fourier index (starts at 0)
  2258.          * @return the coefficient V₁<sup>j</sup> for the given Fourier index i
  2259.          */
  2260.         public double getV1(final int j, final int i) {
  2261.             return v1ij[i][j];
  2262.         }

  2263.         /**
  2264.          * Get the coefficient U₂<sup>j</sup> for Fourier index = 1 (i == 0).
  2265.          *
  2266.          * @param j j index
  2267.          * @return the coefficient U₂<sup>j</sup> for Fourier index = 1 (i == 0)
  2268.          */
  2269.         public double getU2(final int j) {
  2270.             return u2ij[j];
  2271.         }

  2272.         /**
  2273.          * Get the coefficient V₂<sup>j</sup> for Fourier index = 1 (i == 0).
  2274.          *
  2275.          * @param j j index
  2276.          * @return the coefficient V₂<sup>j</sup> for Fourier index = 1 (i == 0)
  2277.          */
  2278.         public double getV2(final int j) {
  2279.             return v2ij[j];
  2280.         }
  2281.     }

  2282.     /**
  2283.      * The U<sub>i</sub><sup>j</sup> and V<sub>i</sub><sup>j</sup> coefficients
  2284.      * described by equations 2.5.3-(21) and 2.5.3-(22) from Danielson.
  2285.      * <p>
  2286.      * The index i takes only the values 1 and 2<br>
  2287.      * For U only the index 0 for j is used.
  2288.      * </p>
  2289.      *
  2290.      * @author Petre Bazavan
  2291.      * @author Lucian Barbulescu
  2292.      */
  2293.     protected static class FieldUijVijCoefficients<T extends CalculusFieldElement<T>> {

  2294.         /**
  2295.          * The U₁<sup>j</sup> coefficients.
  2296.          * <p>
  2297.          * The first index identifies the Fourier coefficients used<br>
  2298.          * Those coefficients are computed for all Fourier C<sub>i</sub><sup>j</sup> and
  2299.          * S<sub>i</sub><sup>j</sup><br>
  2300.          * The only exception is when j = 0 when only the coefficient for fourier index
  2301.          * = 1 (i == 0) is needed.<br>
  2302.          * Also, for fourier index = 1 (i == 0), the coefficients up to 2 * jMax are
  2303.          * computed, because are required to compute the coefficients U₂<sup>j</sup>
  2304.          * </p>
  2305.          */
  2306.         private final T[][] u1ij;

  2307.         /**
  2308.          * The V₁<sup>j</sup> coefficients.
  2309.          * <p>
  2310.          * The first index identifies the Fourier coefficients used<br>
  2311.          * Those coefficients are computed for all Fourier C<sub>i</sub><sup>j</sup> and
  2312.          * S<sub>i</sub><sup>j</sup><br>
  2313.          * for fourier index = 1 (i == 0), the coefficients up to 2 * jMax are computed,
  2314.          * because are required to compute the coefficients V₂<sup>j</sup>
  2315.          * </p>
  2316.          */
  2317.         private final T[][] v1ij;

  2318.         /**
  2319.          * The U₂<sup>j</sup> coefficients.
  2320.          * <p>
  2321.          * Only the coefficients that use the Fourier index = 1 (i == 0) are computed as
  2322.          * they are the only ones required.
  2323.          * </p>
  2324.          */
  2325.         private final T[] u2ij;

  2326.         /**
  2327.          * The V₂<sup>j</sup> coefficients.
  2328.          * <p>
  2329.          * Only the coefficients that use the Fourier index = 1 (i == 0) are computed as
  2330.          * they are the only ones required.
  2331.          * </p>
  2332.          */
  2333.         private final T[] v2ij;

  2334.         /**
  2335.          * The current computed values for the ρ<sub>j</sub> and σ<sub>j</sub>
  2336.          * coefficients.
  2337.          */
  2338.         private final T[][] currentRhoSigmaj;

  2339.         /**
  2340.          * The C<sub>i</sub><sup>j</sup> and the S<sub>i</sub><sup>j</sup> Fourier
  2341.          * coefficients.
  2342.          */
  2343.         private final FieldFourierCjSjCoefficients<T> fourierCjSj;

  2344.         /** The maximum value for j index. */
  2345.         private final int jMax;

  2346.         /**
  2347.          * Constructor.
  2348.          * @param currentRhoSigmaj the current computed values for the ρ<sub>j</sub> and
  2349.          *                         σ<sub>j</sub> coefficients
  2350.          * @param fourierCjSj      the fourier coefficients C<sub>i</sub><sup>j</sup>
  2351.          *                         and the S<sub>i</sub><sup>j</sup>
  2352.          * @param jMax             maximum value for j index
  2353.          * @param field            field used by default
  2354.          */
  2355.         FieldUijVijCoefficients(final T[][] currentRhoSigmaj, final FieldFourierCjSjCoefficients<T> fourierCjSj,
  2356.                 final int jMax, final Field<T> field) {
  2357.             this.currentRhoSigmaj = currentRhoSigmaj;
  2358.             this.fourierCjSj = fourierCjSj;
  2359.             this.jMax = jMax;

  2360.             // initialize the internal arrays.
  2361.             this.u1ij = MathArrays.buildArray(field, 6, 2 * jMax + 1);
  2362.             this.v1ij = MathArrays.buildArray(field, 6, 2 * jMax + 1);
  2363.             this.u2ij = MathArrays.buildArray(field, jMax + 1);
  2364.             this.v2ij = MathArrays.buildArray(field, jMax + 1);

  2365.             // compute the coefficients
  2366.             computeU1V1Coefficients(field);
  2367.             computeU2V2Coefficients(field);
  2368.         }

  2369.         /**
  2370.          * Build the U₁<sup>j</sup> and V₁<sup>j</sup> coefficients.
  2371.          * @param field field used by default
  2372.          */
  2373.         private void computeU1V1Coefficients(final Field<T> field) {
  2374.             // Zero
  2375.             final T zero = field.getZero();

  2376.             // generate the U₁<sup>j</sup> and V₁<sup>j</sup> coefficients
  2377.             // for j >= 1
  2378.             // also the U₁⁰ for Fourier index = 1 (i == 0) coefficient will be computed
  2379.             u1ij[0][0] = zero;
  2380.             for (int j = 1; j <= jMax; j++) {
  2381.                 // compute 1 / j
  2382.                 final double ooj = 1. / j;

  2383.                 for (int i = 0; i < 6; i++) {
  2384.                     // j is aready between 1 and J
  2385.                     u1ij[i][j] = fourierCjSj.getSij(i, j);
  2386.                     v1ij[i][j] = fourierCjSj.getCij(i, j);

  2387.                     // 1 - δ<sub>1j</sub> is 1 for all j > 1
  2388.                     if (j > 1) {
  2389.                         // k starts with 1 because j-J is less than or equal to 0
  2390.                         for (int kIndex = 1; kIndex <= j - 1; kIndex++) {
  2391.                             // C<sub>i</sub><sup>j-k</sup> * σ<sub>k</sub> +
  2392.                             // S<sub>i</sub><sup>j-k</sup> * ρ<sub>k</sub>
  2393.                             u1ij[i][j] = u1ij[i][j]
  2394.                                     .add(fourierCjSj.getCij(i, j - kIndex).multiply(currentRhoSigmaj[1][kIndex]).add(
  2395.                                             fourierCjSj.getSij(i, j - kIndex).multiply(currentRhoSigmaj[0][kIndex])));

  2396.                             // C<sub>i</sub><sup>j-k</sup> * ρ<sub>k</sub> -
  2397.                             // S<sub>i</sub><sup>j-k</sup> * σ<sub>k</sub>
  2398.                             v1ij[i][j] = v1ij[i][j].add(
  2399.                                     fourierCjSj.getCij(i, j - kIndex).multiply(currentRhoSigmaj[0][kIndex]).subtract(
  2400.                                             fourierCjSj.getSij(i, j - kIndex).multiply(currentRhoSigmaj[1][kIndex])));
  2401.                         }
  2402.                     }

  2403.                     // since j must be between 1 and J-1 and is already between 1 and J
  2404.                     // the following sum is skiped only for j = jMax
  2405.                     if (j != jMax) {
  2406.                         for (int kIndex = 1; kIndex <= jMax - j; kIndex++) {
  2407.                             // -C<sub>i</sub><sup>j+k</sup> * σ<sub>k</sub> +
  2408.                             // S<sub>i</sub><sup>j+k</sup> * ρ<sub>k</sub>
  2409.                             u1ij[i][j] = u1ij[i][j].add(fourierCjSj.getCij(i, j + kIndex).negate()
  2410.                                     .multiply(currentRhoSigmaj[1][kIndex])
  2411.                                     .add(fourierCjSj.getSij(i, j + kIndex).multiply(currentRhoSigmaj[0][kIndex])));

  2412.                             // C<sub>i</sub><sup>j+k</sup> * ρ<sub>k</sub> +
  2413.                             // S<sub>i</sub><sup>j+k</sup> * σ<sub>k</sub>
  2414.                             v1ij[i][j] = v1ij[i][j]
  2415.                                     .add(fourierCjSj.getCij(i, j + kIndex).multiply(currentRhoSigmaj[0][kIndex]).add(
  2416.                                             fourierCjSj.getSij(i, j + kIndex).multiply(currentRhoSigmaj[1][kIndex])));
  2417.                         }
  2418.                     }

  2419.                     for (int kIndex = 1; kIndex <= jMax; kIndex++) {
  2420.                         // C<sub>i</sub><sup>k</sup> * σ<sub>j+k</sub> -
  2421.                         // S<sub>i</sub><sup>k</sup> * ρ<sub>j+k</sub>
  2422.                         u1ij[i][j] = u1ij[i][j].add(fourierCjSj.getCij(i, kIndex).negate()
  2423.                                 .multiply(currentRhoSigmaj[1][j + kIndex])
  2424.                                 .subtract(fourierCjSj.getSij(i, kIndex).multiply(currentRhoSigmaj[0][j + kIndex])));

  2425.                         // C<sub>i</sub><sup>k</sup> * ρ<sub>j+k</sub> +
  2426.                         // S<sub>i</sub><sup>k</sup> * σ<sub>j+k</sub>
  2427.                         v1ij[i][j] = v1ij[i][j]
  2428.                                 .add(fourierCjSj.getCij(i, kIndex).multiply(currentRhoSigmaj[0][j + kIndex])
  2429.                                         .add(fourierCjSj.getSij(i, kIndex).multiply(currentRhoSigmaj[1][j + kIndex])));
  2430.                     }

  2431.                     // divide by 1 / j
  2432.                     u1ij[i][j] = u1ij[i][j].multiply(-ooj);
  2433.                     v1ij[i][j] = v1ij[i][j].multiply(ooj);

  2434.                     // if index = 1 (i == 0) add the computed terms to U₁⁰
  2435.                     if (i == 0) {
  2436.                         // - (U₁<sup>j</sup> * ρ<sub>j</sub> + V₁<sup>j</sup> * σ<sub>j</sub>
  2437.                         u1ij[0][0] = u1ij[0][0].add(u1ij[0][j].negate().multiply(currentRhoSigmaj[0][j])
  2438.                                 .subtract(v1ij[0][j].multiply(currentRhoSigmaj[1][j])));
  2439.                     }
  2440.                 }
  2441.             }

  2442.             // Terms with j > jMax are required only when computing the coefficients
  2443.             // U₂<sup>j</sup> and V₂<sup>j</sup>
  2444.             // and those coefficients are only required for Fourier index = 1 (i == 0).
  2445.             for (int j = jMax + 1; j <= 2 * jMax; j++) {
  2446.                 // compute 1 / j
  2447.                 final double ooj = 1. / j;
  2448.                 // the value of i is 0
  2449.                 u1ij[0][j] = zero;
  2450.                 v1ij[0][j] = zero;

  2451.                 // k starts from j-J as it is always greater than or equal to 1
  2452.                 for (int kIndex = j - jMax; kIndex <= j - 1; kIndex++) {
  2453.                     // C<sub>i</sub><sup>j-k</sup> * σ<sub>k</sub> +
  2454.                     // S<sub>i</sub><sup>j-k</sup> * ρ<sub>k</sub>
  2455.                     u1ij[0][j] = u1ij[0][j].add(fourierCjSj.getCij(0, j - kIndex).multiply(currentRhoSigmaj[1][kIndex])
  2456.                             .add(fourierCjSj.getSij(0, j - kIndex).multiply(currentRhoSigmaj[0][kIndex])));

  2457.                     // C<sub>i</sub><sup>j-k</sup> * ρ<sub>k</sub> -
  2458.                     // S<sub>i</sub><sup>j-k</sup> * σ<sub>k</sub>
  2459.                     v1ij[0][j] = v1ij[0][j].add(fourierCjSj.getCij(0, j - kIndex).multiply(currentRhoSigmaj[0][kIndex])
  2460.                             .subtract(fourierCjSj.getSij(0, j - kIndex).multiply(currentRhoSigmaj[1][kIndex])));
  2461.                 }
  2462.                 for (int kIndex = 1; kIndex <= jMax; kIndex++) {
  2463.                     // C<sub>i</sub><sup>k</sup> * σ<sub>j+k</sub> -
  2464.                     // S<sub>i</sub><sup>k</sup> * ρ<sub>j+k</sub>
  2465.                     u1ij[0][j] = u1ij[0][j]
  2466.                             .add(fourierCjSj.getCij(0, kIndex).negate().multiply(currentRhoSigmaj[1][j + kIndex])
  2467.                                     .subtract(fourierCjSj.getSij(0, kIndex).multiply(currentRhoSigmaj[0][j + kIndex])));

  2468.                     // C<sub>i</sub><sup>k</sup> * ρ<sub>j+k</sub> +
  2469.                     // S<sub>i</sub><sup>k</sup> * σ<sub>j+k</sub>
  2470.                     v1ij[0][j] = v1ij[0][j].add(fourierCjSj.getCij(0, kIndex).multiply(currentRhoSigmaj[0][j + kIndex])
  2471.                             .add(fourierCjSj.getSij(0, kIndex).multiply(currentRhoSigmaj[1][j + kIndex])));
  2472.                 }

  2473.                 // divide by 1 / j
  2474.                 u1ij[0][j] = u1ij[0][j].multiply(-ooj);
  2475.                 v1ij[0][j] = v1ij[0][j].multiply(ooj);
  2476.             }
  2477.         }

  2478.         /**
  2479.          * Build the U₁<sup>j</sup> and V₁<sup>j</sup> coefficients.
  2480.          * <p>
  2481.          * Only the coefficients for Fourier index = 1 (i == 0) are required.
  2482.          * </p>
  2483.          * @param field field used by default
  2484.          */
  2485.         private void computeU2V2Coefficients(final Field<T> field) {
  2486.             for (int j = 1; j <= jMax; j++) {
  2487.                 // compute 1 / j
  2488.                 final double ooj = 1. / j;

  2489.                 // only the values for i == 0 are computed
  2490.                 u2ij[j] = v1ij[0][j];
  2491.                 v2ij[j] = u1ij[0][j];

  2492.                 // 1 - δ<sub>1j</sub> is 1 for all j > 1
  2493.                 if (j > 1) {
  2494.                     for (int l = 1; l <= j - 1; l++) {
  2495.                         // U₁<sup>j-l</sup> * σ<sub>l</sub> +
  2496.                         // V₁<sup>j-l</sup> * ρ<sub>l</sub>
  2497.                         u2ij[j] = u2ij[j].add(u1ij[0][j - l].multiply(currentRhoSigmaj[1][l])
  2498.                                 .add(v1ij[0][j - l].multiply(currentRhoSigmaj[0][l])));

  2499.                         // U₁<sup>j-l</sup> * ρ<sub>l</sub> -
  2500.                         // V₁<sup>j-l</sup> * σ<sub>l</sub>
  2501.                         v2ij[j] = v2ij[j].add(u1ij[0][j - l].multiply(currentRhoSigmaj[0][l])
  2502.                                 .subtract(v1ij[0][j - l].multiply(currentRhoSigmaj[1][l])));
  2503.                     }
  2504.                 }

  2505.                 for (int l = 1; l <= jMax; l++) {
  2506.                     // -U₁<sup>j+l</sup> * σ<sub>l</sub> +
  2507.                     // U₁<sup>l</sup> * σ<sub>j+l</sub> +
  2508.                     // V₁<sup>j+l</sup> * ρ<sub>l</sub> -
  2509.                     // V₁<sup>l</sup> * ρ<sub>j+l</sub>
  2510.                     u2ij[j] = u2ij[j].add(u1ij[0][j + l].negate().multiply(currentRhoSigmaj[1][l])
  2511.                             .add(u1ij[0][l].multiply(currentRhoSigmaj[1][j + l]))
  2512.                             .add(v1ij[0][j + l].multiply(currentRhoSigmaj[0][l]))
  2513.                             .subtract(v1ij[0][l].multiply(currentRhoSigmaj[0][j + l])));

  2514.                     // U₁<sup>j+l</sup> * ρ<sub>l</sub> +
  2515.                     // U₁<sup>l</sup> * ρ<sub>j+l</sub> +
  2516.                     // V₁<sup>j+l</sup> * σ<sub>l</sub> +
  2517.                     // V₁<sup>l</sup> * σ<sub>j+l</sub>
  2518.                     u2ij[j] = u2ij[j].add(u1ij[0][j + l].multiply(currentRhoSigmaj[0][l])
  2519.                             .add(u1ij[0][l].multiply(currentRhoSigmaj[0][j + l]))
  2520.                             .add(v1ij[0][j + l].multiply(currentRhoSigmaj[1][l]))
  2521.                             .add(v1ij[0][l].multiply(currentRhoSigmaj[1][j + l])));
  2522.                 }

  2523.                 // divide by 1 / j
  2524.                 u2ij[j] = u2ij[j].multiply(-ooj);
  2525.                 v2ij[j] = v2ij[j].multiply(ooj);
  2526.             }
  2527.         }

  2528.         /**
  2529.          * Get the coefficient U₁<sup>j</sup> for Fourier index i.
  2530.          *
  2531.          * @param j j index
  2532.          * @param i Fourier index (starts at 0)
  2533.          * @return the coefficient U₁<sup>j</sup> for the given Fourier index i
  2534.          */
  2535.         public T getU1(final int j, final int i) {
  2536.             return u1ij[i][j];
  2537.         }

  2538.         /**
  2539.          * Get the coefficient V₁<sup>j</sup> for Fourier index i.
  2540.          *
  2541.          * @param j j index
  2542.          * @param i Fourier index (starts at 0)
  2543.          * @return the coefficient V₁<sup>j</sup> for the given Fourier index i
  2544.          */
  2545.         public T getV1(final int j, final int i) {
  2546.             return v1ij[i][j];
  2547.         }

  2548.         /**
  2549.          * Get the coefficient U₂<sup>j</sup> for Fourier index = 1 (i == 0).
  2550.          *
  2551.          * @param j j index
  2552.          * @return the coefficient U₂<sup>j</sup> for Fourier index = 1 (i == 0)
  2553.          */
  2554.         public T getU2(final int j) {
  2555.             return u2ij[j];
  2556.         }

  2557.         /**
  2558.          * Get the coefficient V₂<sup>j</sup> for Fourier index = 1 (i == 0).
  2559.          *
  2560.          * @param j j index
  2561.          * @return the coefficient V₂<sup>j</sup> for Fourier index = 1 (i == 0)
  2562.          */
  2563.         public T getV2(final int j) {
  2564.             return v2ij[j];
  2565.         }
  2566.     }

  2567.     /** Coefficients valid for one time slot. */
  2568.     protected static class Slot {

  2569.         /**
  2570.          * The coefficients D<sub>i</sub><sup>j</sup>.
  2571.          * <p>
  2572.          * Only for j = 1 and j = 2 the coefficients are not 0. <br>
  2573.          * i corresponds to the equinoctial element, as follows: - i=0 for a <br/>
  2574.          * - i=1 for k <br/>
  2575.          * - i=2 for h <br/>
  2576.          * - i=3 for q <br/>
  2577.          * - i=4 for p <br/>
  2578.          * - i=5 for λ <br/>
  2579.          * </p>
  2580.          */
  2581.         private final ShortPeriodicsInterpolatedCoefficient[] dij;

  2582.         /**
  2583.          * The coefficients C<sub>i</sub><sup>j</sup>.
  2584.          * <p>
  2585.          * The index order is cij[j][i] <br/>
  2586.          * i corresponds to the equinoctial element, as follows: <br/>
  2587.          * - i=0 for a <br/>
  2588.          * - i=1 for k <br/>
  2589.          * - i=2 for h <br/>
  2590.          * - i=3 for q <br/>
  2591.          * - i=4 for p <br/>
  2592.          * - i=5 for λ <br/>
  2593.          * </p>
  2594.          */
  2595.         private final ShortPeriodicsInterpolatedCoefficient[] cij;

  2596.         /**
  2597.          * The coefficients S<sub>i</sub><sup>j</sup>.
  2598.          * <p>
  2599.          * The index order is sij[j][i] <br/>
  2600.          * i corresponds to the equinoctial element, as follows: <br/>
  2601.          * - i=0 for a <br/>
  2602.          * - i=1 for k <br/>
  2603.          * - i=2 for h <br/>
  2604.          * - i=3 for q <br/>
  2605.          * - i=4 for p <br/>
  2606.          * - i=5 for λ <br/>
  2607.          * </p>
  2608.          */
  2609.         private final ShortPeriodicsInterpolatedCoefficient[] sij;

  2610.         /**
  2611.          * Simple constructor.
  2612.          * @param jMax                maximum value for j index
  2613.          * @param interpolationPoints number of points used in the interpolation process
  2614.          */
  2615.         Slot(final int jMax, final int interpolationPoints) {

  2616.             dij = new ShortPeriodicsInterpolatedCoefficient[3];
  2617.             cij = new ShortPeriodicsInterpolatedCoefficient[jMax + 1];
  2618.             sij = new ShortPeriodicsInterpolatedCoefficient[jMax + 1];

  2619.             // Initialize the C<sub>i</sub><sup>j</sup>, S<sub>i</sub><sup>j</sup> and
  2620.             // D<sub>i</sub><sup>j</sup> coefficients
  2621.             for (int j = 0; j <= jMax; j++) {
  2622.                 cij[j] = new ShortPeriodicsInterpolatedCoefficient(interpolationPoints);
  2623.                 if (j > 0) {
  2624.                     sij[j] = new ShortPeriodicsInterpolatedCoefficient(interpolationPoints);
  2625.                 }
  2626.                 // Initialize only the non-zero D<sub>i</sub><sup>j</sup> coefficients
  2627.                 if (j == 1 || j == 2) {
  2628.                     dij[j] = new ShortPeriodicsInterpolatedCoefficient(interpolationPoints);
  2629.                 }
  2630.             }

  2631.         }

  2632.     }

  2633.     /** Coefficients valid for one time slot. */
  2634.     protected static class FieldSlot<T extends CalculusFieldElement<T>> {

  2635.         /**
  2636.          * The coefficients D<sub>i</sub><sup>j</sup>.
  2637.          * <p>
  2638.          * Only for j = 1 and j = 2 the coefficients are not 0. <br>
  2639.          * i corresponds to the equinoctial element, as follows: - i=0 for a <br/>
  2640.          * - i=1 for k <br/>
  2641.          * - i=2 for h <br/>
  2642.          * - i=3 for q <br/>
  2643.          * - i=4 for p <br/>
  2644.          * - i=5 for λ <br/>
  2645.          * </p>
  2646.          */
  2647.         private final FieldShortPeriodicsInterpolatedCoefficient<T>[] dij;

  2648.         /**
  2649.          * The coefficients C<sub>i</sub><sup>j</sup>.
  2650.          * <p>
  2651.          * The index order is cij[j][i] <br/>
  2652.          * i corresponds to the equinoctial element, as follows: <br/>
  2653.          * - i=0 for a <br/>
  2654.          * - i=1 for k <br/>
  2655.          * - i=2 for h <br/>
  2656.          * - i=3 for q <br/>
  2657.          * - i=4 for p <br/>
  2658.          * - i=5 for λ <br/>
  2659.          * </p>
  2660.          */
  2661.         private final FieldShortPeriodicsInterpolatedCoefficient<T>[] cij;

  2662.         /**
  2663.          * The coefficients S<sub>i</sub><sup>j</sup>.
  2664.          * <p>
  2665.          * The index order is sij[j][i] <br/>
  2666.          * i corresponds to the equinoctial element, as follows: <br/>
  2667.          * - i=0 for a <br/>
  2668.          * - i=1 for k <br/>
  2669.          * - i=2 for h <br/>
  2670.          * - i=3 for q <br/>
  2671.          * - i=4 for p <br/>
  2672.          * - i=5 for λ <br/>
  2673.          * </p>
  2674.          */
  2675.         private final FieldShortPeriodicsInterpolatedCoefficient<T>[] sij;

  2676.         /**
  2677.          * Simple constructor.
  2678.          * @param jMax                maximum value for j index
  2679.          * @param interpolationPoints number of points used in the interpolation process
  2680.          */
  2681.         @SuppressWarnings("unchecked")
  2682.         FieldSlot(final int jMax, final int interpolationPoints) {

  2683.             dij = (FieldShortPeriodicsInterpolatedCoefficient<T>[]) Array
  2684.                     .newInstance(FieldShortPeriodicsInterpolatedCoefficient.class, 3);
  2685.             cij = (FieldShortPeriodicsInterpolatedCoefficient<T>[]) Array
  2686.                     .newInstance(FieldShortPeriodicsInterpolatedCoefficient.class, jMax + 1);
  2687.             sij = (FieldShortPeriodicsInterpolatedCoefficient<T>[]) Array
  2688.                     .newInstance(FieldShortPeriodicsInterpolatedCoefficient.class, jMax + 1);

  2689.             // Initialize the C<sub>i</sub><sup>j</sup>, S<sub>i</sub><sup>j</sup> and
  2690.             // D<sub>i</sub><sup>j</sup> coefficients
  2691.             for (int j = 0; j <= jMax; j++) {
  2692.                 cij[j] = new FieldShortPeriodicsInterpolatedCoefficient<>(interpolationPoints);
  2693.                 if (j > 0) {
  2694.                     sij[j] = new FieldShortPeriodicsInterpolatedCoefficient<>(interpolationPoints);
  2695.                 }
  2696.                 // Initialize only the non-zero D<sub>i</sub><sup>j</sup> coefficients
  2697.                 if (j == 1 || j == 2) {
  2698.                     dij[j] = new FieldShortPeriodicsInterpolatedCoefficient<>(interpolationPoints);
  2699.                 }
  2700.             }

  2701.         }

  2702.     }

  2703. }