AbstractGaussianContribution.java

  1. /* Copyright 2002-2023 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.PositionAngleType;
  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 <T extends CalculusFieldElement<T>> void init(final FieldSpacecraftState<T> initialState, final FieldAbsoluteDate<T> target) {
  174.         // Initialize the numerical force model
  175.         contribution.init(initialState, target);
  176.     }

  177.     /** {@inheritDoc} */
  178.     @Override
  179.     public List<ParameterDriver> getParametersDrivers() {

  180.         // Parameter drivers
  181.         final List<ParameterDriver> drivers = new ArrayList<>();

  182.         // Loop on drivers (without central attraction coefficient driver)
  183.         for (final ParameterDriver driver : getParametersDriversWithoutMu()) {
  184.             drivers.add(driver);
  185.         }

  186.         // We put central attraction coefficient driver at the end of the array
  187.         drivers.add(gmParameterDriver);
  188.         return drivers;

  189.     }

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

  201.     /** {@inheritDoc} */
  202.     @Override
  203.     public List<ShortPeriodTerms> initializeShortPeriodTerms(final AuxiliaryElements auxiliaryElements, final PropagationType type,
  204.             final double[] parameters) {

  205.         final List<ShortPeriodTerms> list = new ArrayList<ShortPeriodTerms>();
  206.         gaussianSPCoefs = new GaussianShortPeriodicCoefficients(coefficientsKeyPrefix, JMAX, INTERPOLATION_POINTS,
  207.                 new TimeSpanMap<Slot>(new Slot(JMAX, INTERPOLATION_POINTS)));
  208.         list.add(gaussianSPCoefs);
  209.         return list;

  210.     }

  211.     /** {@inheritDoc} */
  212.     @Override
  213.     public <T extends CalculusFieldElement<T>> List<FieldShortPeriodTerms<T>> initializeShortPeriodTerms(
  214.             final FieldAuxiliaryElements<T> auxiliaryElements, final PropagationType type, final T[] parameters) {

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

  216.         final FieldGaussianShortPeriodicCoefficients<T> fgspc = new FieldGaussianShortPeriodicCoefficients<>(
  217.                 coefficientsKeyPrefix, JMAX, INTERPOLATION_POINTS,
  218.                 new FieldTimeSpanMap<>(new FieldSlot<>(JMAX, INTERPOLATION_POINTS), field));
  219.         gaussianFieldSPCoefs.put(field, fgspc);
  220.         return Collections.singletonList(fgspc);
  221.     }

  222.     /**
  223.      * Performs initialization at each integration step for the current force model.
  224.      * <p>
  225.      * This method aims at being called before mean elements rates computation.
  226.      * </p>
  227.      * @param auxiliaryElements auxiliary elements related to the current orbit
  228.      * @param parameters        parameters values of the force model parameters
  229.      *                          only 1 value for each parameterDriver
  230.      * @return new force model context
  231.      */
  232.     private AbstractGaussianContributionContext initializeStep(final AuxiliaryElements auxiliaryElements,
  233.             final double[] parameters) {
  234.         return new AbstractGaussianContributionContext(auxiliaryElements, parameters);
  235.     }

  236.     /**
  237.      * Performs initialization at each integration step for the current force model.
  238.      * <p>
  239.      * This method aims at being called before mean elements rates computation.
  240.      * </p>
  241.      * @param <T>               type of the elements
  242.      * @param auxiliaryElements auxiliary elements related to the current orbit
  243.      * @param parameters        parameters values of the force model parameters
  244.      *                          (only 1 values for each parameters corresponding
  245.      *                          to state date) by getting the parameters for a specific date.
  246.      * @return new force model context
  247.      */
  248.     private <T extends CalculusFieldElement<T>> FieldAbstractGaussianContributionContext<T> initializeStep(
  249.             final FieldAuxiliaryElements<T> auxiliaryElements, final T[] parameters) {
  250.         return new FieldAbstractGaussianContributionContext<>(auxiliaryElements, parameters);
  251.     }

  252.     /** {@inheritDoc} */
  253.     @Override
  254.     public double[] getMeanElementRate(final SpacecraftState state, final AuxiliaryElements auxiliaryElements,
  255.             final double[] parameters) {

  256.         // Container for attributes

  257.         final AbstractGaussianContributionContext context = initializeStep(auxiliaryElements, parameters);
  258.         double[] meanElementRate = new double[6];
  259.         // Computes the limits for the integral
  260.         final double[] ll = getLLimits(state, auxiliaryElements);
  261.         // Computes integrated mean element rates if Llow < Lhigh
  262.         if (ll[0] < ll[1]) {
  263.             meanElementRate = getMeanElementRate(state, integrator, ll[0], ll[1], context, parameters);
  264.             if (isDirty) {
  265.                 boolean next = true;
  266.                 for (int i = 0; i < MAX_ORDER_RANK && next; i++) {
  267.                     final double[] meanRates = getMeanElementRate(state, new GaussQuadrature(GAUSS_ORDER[i]), ll[0],
  268.                             ll[1], context, parameters);
  269.                     if (getRatesDiff(meanElementRate, meanRates, context) < threshold) {
  270.                         integrator = new GaussQuadrature(GAUSS_ORDER[i]);
  271.                         next = false;
  272.                     }
  273.                 }
  274.                 isDirty = false;
  275.             }
  276.         }
  277.         return meanElementRate;
  278.     }

  279.     /** {@inheritDoc} */
  280.     @Override
  281.     public <T extends CalculusFieldElement<T>> T[] getMeanElementRate(final FieldSpacecraftState<T> state,
  282.             final FieldAuxiliaryElements<T> auxiliaryElements, final T[] parameters) {

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

  286.         T[] meanElementRate = MathArrays.buildArray(field, 6);
  287.         // Computes the limits for the integral
  288.         final T[] ll = getLLimits(state, auxiliaryElements);
  289.         // Computes integrated mean element rates if Llow < Lhigh
  290.         if (ll[0].getReal() < ll[1].getReal()) {
  291.             meanElementRate = getMeanElementRate(state, integrator, ll[0], ll[1], context, parameters);
  292.             if (isDirty) {
  293.                 boolean next = true;
  294.                 for (int i = 0; i < MAX_ORDER_RANK && next; i++) {
  295.                     final T[] meanRates = getMeanElementRate(state, new GaussQuadrature(GAUSS_ORDER[i]), ll[0], ll[1],
  296.                             context, parameters);
  297.                     if (getRatesDiff(meanElementRate, meanRates, context).getReal() < threshold) {
  298.                         integrator = new GaussQuadrature(GAUSS_ORDER[i]);
  299.                         next = false;
  300.                     }
  301.                 }
  302.                 isDirty = false;
  303.             }
  304.         }

  305.         return meanElementRate;
  306.     }

  307.     /**
  308.      * Compute the limits in L, the true longitude, for integration.
  309.      *
  310.      * @param state             current state information: date, kinematics,
  311.      *                          attitude
  312.      * @param auxiliaryElements auxiliary elements related to the current orbit
  313.      * @return the integration limits in L
  314.      */
  315.     protected abstract double[] getLLimits(SpacecraftState state, AuxiliaryElements auxiliaryElements);

  316.     /**
  317.      * Compute the limits in L, the true longitude, for integration.
  318.      *
  319.      * @param <T>               type of the elements
  320.      * @param state             current state information: date, kinematics,
  321.      *                          attitude
  322.      * @param auxiliaryElements auxiliary elements related to the current orbit
  323.      * @return the integration limits in L
  324.      */
  325.     protected abstract <T extends CalculusFieldElement<T>> T[] getLLimits(FieldSpacecraftState<T> state,
  326.             FieldAuxiliaryElements<T> auxiliaryElements);

  327.     /**
  328.      * Computes the mean equinoctial elements rates da<sub>i</sub> / dt.
  329.      *
  330.      * @param state      current state
  331.      * @param gauss      Gauss quadrature
  332.      * @param low        lower bound of the integral interval
  333.      * @param high       upper bound of the integral interval
  334.      * @param context    container for attributes
  335.      * @param parameters values of the force model parameters
  336.      * at state date (1 values for each parameters)
  337.      * @return the mean element rates
  338.      */
  339.     protected double[] getMeanElementRate(final SpacecraftState state, final GaussQuadrature gauss, final double low,
  340.             final double high, final AbstractGaussianContributionContext context, final double[] parameters) {

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

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

  344.         // Constant multiplier for integral
  345.         final double coef = 1. / (2. * FastMath.PI * auxiliaryElements.getB());
  346.         // Corrects mean element rates
  347.         for (int i = 0; i < 6; i++) {
  348.             meanElementRate[i] *= coef;
  349.         }
  350.         return meanElementRate;
  351.     }

  352.     /**
  353.      * Computes the mean equinoctial elements rates da<sub>i</sub> / dt.
  354.      *
  355.      * @param <T>        type of the elements
  356.      * @param state      current state
  357.      * @param gauss      Gauss quadrature
  358.      * @param low        lower bound of the integral interval
  359.      * @param high       upper bound of the integral interval
  360.      * @param context    container for attributes
  361.      * @param parameters values of the force model parameters(1 values for each parameters)
  362.      * @return the mean element rates
  363.      */
  364.     protected <T extends CalculusFieldElement<T>> T[] getMeanElementRate(final FieldSpacecraftState<T> state,
  365.             final GaussQuadrature gauss, final T low, final T high,
  366.             final FieldAbstractGaussianContributionContext<T> context, final T[] parameters) {

  367.         // Field
  368.         final Field<T> field = context.getA().getField();

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

  371.         final T[] meanElementRate = gauss.integrate(new FieldIntegrableFunction<>(state, true, 0, parameters, field),
  372.                 low, high, field);
  373.         // Constant multiplier for integral
  374.         final T coef = auxiliaryElements.getB().multiply(low.getPi()).multiply(2.).reciprocal();
  375.         // Corrects mean element rates
  376.         for (int i = 0; i < 6; i++) {
  377.             meanElementRate[i] = meanElementRate[i].multiply(coef);
  378.         }
  379.         return meanElementRate;
  380.     }

  381.     /**
  382.      * Estimates the weighted magnitude of the difference between 2 sets of
  383.      * equinoctial elements rates.
  384.      *
  385.      * @param meanRef reference rates
  386.      * @param meanCur current rates
  387.      * @param context container for attributes
  388.      * @return estimated magnitude of weighted differences
  389.      */
  390.     private double getRatesDiff(final double[] meanRef, final double[] meanCur,
  391.             final AbstractGaussianContributionContext context) {

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

  394.         double maxDiff = FastMath.abs(meanRef[0] - meanCur[0]) / auxiliaryElements.getSma();
  395.         // Corrects mean element rates
  396.         for (int i = 1; i < meanRef.length; i++) {
  397.             maxDiff = FastMath.max(maxDiff, FastMath.abs(meanRef[i] - meanCur[i]));
  398.         }
  399.         return maxDiff;
  400.     }

  401.     /**
  402.      * Estimates the weighted magnitude of the difference between 2 sets of
  403.      * equinoctial elements rates.
  404.      *
  405.      * @param <T>     type of the elements
  406.      * @param meanRef reference rates
  407.      * @param meanCur current rates
  408.      * @param context container for attributes
  409.      * @return estimated magnitude of weighted differences
  410.      */
  411.     private <T extends CalculusFieldElement<T>> T getRatesDiff(final T[] meanRef, final T[] meanCur,
  412.             final FieldAbstractGaussianContributionContext<T> context) {

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

  415.         T maxDiff = FastMath.abs(meanRef[0].subtract(meanCur[0])).divide(auxiliaryElements.getSma());
  416.         ;
  417.         // Corrects mean element rates
  418.         for (int i = 1; i < meanRef.length; i++) {
  419.             maxDiff = FastMath.max(maxDiff, FastMath.abs(meanRef[i].subtract(meanCur[i])));
  420.         }
  421.         return maxDiff;
  422.     }

  423.     /** {@inheritDoc} */
  424.     @Override
  425.     public void registerAttitudeProvider(final AttitudeProvider provider) {
  426.         this.attitudeProvider = provider;
  427.     }

  428.     /** {@inheritDoc} */
  429.     @Override
  430.     public void updateShortPeriodTerms(final double[] parameters, final SpacecraftState... meanStates) {

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

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

  435.             // Container of attributes
  436.             // Extract the proper parameters valid for the corresponding meanState date from the input array
  437.             final double[] extractedParameters = this.extractParameters(parameters, auxiliaryElements.getDate());
  438.             final AbstractGaussianContributionContext context = initializeStep(auxiliaryElements, extractedParameters);

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

  441.             // Generate the Cij and Sij coefficients
  442.             final FourierCjSjCoefficients fourierCjSj = new FourierCjSjCoefficients(meanState, JMAX, auxiliaryElements,
  443.                                                                                     extractedParameters);

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

  446.             gaussianSPCoefs.computeCoefficients(meanState, slot, fourierCjSj, uijvij, context.getMeanMotion(),
  447.                     auxiliaryElements.getSma());

  448.         }

  449.     }

  450.     /** {@inheritDoc} */
  451.     @Override
  452.     @SuppressWarnings("unchecked")
  453.     public <T extends CalculusFieldElement<T>> void updateShortPeriodTerms(final T[] parameters,
  454.             final FieldSpacecraftState<T>... meanStates) {

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

  457.         final FieldGaussianShortPeriodicCoefficients<T> fgspc = (FieldGaussianShortPeriodicCoefficients<T>) gaussianFieldSPCoefs
  458.                 .get(field);
  459.         final FieldSlot<T> slot = fgspc.createSlot(meanStates);
  460.         for (final FieldSpacecraftState<T> meanState : meanStates) {

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

  463.             // Container of attributes
  464.             // Extract the proper parameters valid for the corresponding meanState date from the input array
  465.             final T[] extractedParameters = this.extractParameters(parameters, auxiliaryElements.getDate());
  466.             final FieldAbstractGaussianContributionContext<T> context = initializeStep(auxiliaryElements, extractedParameters);

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

  469.             // Generate the Cij and Sij coefficients
  470.             final FieldFourierCjSjCoefficients<T> fourierCjSj = new FieldFourierCjSjCoefficients<>(meanState, JMAX,
  471.                     auxiliaryElements, extractedParameters, field);

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

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

  477.         }

  478.     }

  479.     /**
  480.      * Compute the auxiliary quantities ρ<sub>j</sub> and σ<sub>j</sub>.
  481.      * <p>
  482.      * The expressions used are equations 2.5.3-(4) from the Danielson paper. <br/>
  483.      * ρ<sub>j</sub> = (1+jB)(-b)<sup>j</sup>C<sub>j</sub>(k, h) <br/>
  484.      * σ<sub>j</sub> = (1+jB)(-b)<sup>j</sup>S<sub>j</sub>(k, h) <br/>
  485.      * </p>
  486.      * @param date              current date
  487.      * @param auxiliaryElements auxiliary elements related to the current orbit
  488.      * @return computed coefficients
  489.      */
  490.     private double[][] computeRhoSigmaCoefficients(final AbsoluteDate date, final AuxiliaryElements auxiliaryElements) {
  491.         final double[][] currentRhoSigmaj = new double[2][3 * JMAX + 1];
  492.         final CjSjCoefficient cjsjKH = new CjSjCoefficient(auxiliaryElements.getK(), auxiliaryElements.getH());
  493.         final double b = 1. / (1 + auxiliaryElements.getB());

  494.         // (-b)<sup>j</sup>
  495.         double mbtj = 1;

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

  497.             // Compute current rho and sigma;
  498.             mbtj *= -b;
  499.             final double coef = (1 + j * auxiliaryElements.getB()) * mbtj;
  500.             currentRhoSigmaj[0][j] = coef * cjsjKH.getCj(j);
  501.             currentRhoSigmaj[1][j] = coef * cjsjKH.getSj(j);
  502.         }
  503.         return currentRhoSigmaj;
  504.     }

  505.     /**
  506.      * Compute the auxiliary quantities ρ<sub>j</sub> and σ<sub>j</sub>.
  507.      * <p>
  508.      * The expressions used are equations 2.5.3-(4) from the Danielson paper. <br/>
  509.      * ρ<sub>j</sub> = (1+jB)(-b)<sup>j</sup>C<sub>j</sub>(k, h) <br/>
  510.      * σ<sub>j</sub> = (1+jB)(-b)<sup>j</sup>S<sub>j</sub>(k, h) <br/>
  511.      * </p>
  512.      * @param <T>     type of the elements
  513.      * @param date    current date
  514.      * @param context container for attributes
  515.      * @param field   field used by default
  516.      * @return computed coefficients
  517.      */
  518.     private <T extends CalculusFieldElement<T>> T[][] computeRhoSigmaCoefficients(final FieldAbsoluteDate<T> date,
  519.             final FieldAbstractGaussianContributionContext<T> context, final Field<T> field) {
  520.         // zero
  521.         final T zero = field.getZero();

  522.         final FieldAuxiliaryElements<T> auxiliaryElements = context.getFieldAuxiliaryElements();
  523.         final T[][] currentRhoSigmaj = MathArrays.buildArray(field, 2, 3 * JMAX + 1);
  524.         final FieldCjSjCoefficient<T> cjsjKH = new FieldCjSjCoefficient<>(auxiliaryElements.getK(),
  525.                 auxiliaryElements.getH(), field);
  526.         final T b = auxiliaryElements.getB().add(1.).reciprocal();

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

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

  530.             // Compute current rho and sigma;
  531.             mbtj = mbtj.multiply(b.negate());
  532.             final T coef = mbtj.multiply(auxiliaryElements.getB().multiply(j).add(1.));
  533.             currentRhoSigmaj[0][j] = coef.multiply(cjsjKH.getCj(j));
  534.             currentRhoSigmaj[1][j] = coef.multiply(cjsjKH.getSj(j));
  535.         }
  536.         return currentRhoSigmaj;
  537.     }

  538.     /**
  539.      * Internal class for numerical quadrature.
  540.      * <p>
  541.      * This class is a rewrite of {@link IntegrableFunction} for field elements
  542.      * </p>
  543.      * @param <T> type of the field elements
  544.      */
  545.     protected class FieldIntegrableFunction<T extends CalculusFieldElement<T>>
  546.             implements CalculusFieldUnivariateVectorFunction<T> {

  547.         /** Current state. */
  548.         private final FieldSpacecraftState<T> state;

  549.         /**
  550.          * Signal that this class is used to compute the values required by the mean
  551.          * element variations or by the short periodic element variations.
  552.          */
  553.         private final boolean meanMode;

  554.         /**
  555.          * The j index.
  556.          * <p>
  557.          * Used only for short periodic variation. Ignored for mean elements variation.
  558.          * </p>
  559.          */
  560.         private final int j;

  561.         /** Container for attributes. */
  562.         private final FieldAbstractGaussianContributionContext<T> context;

  563.         /** Auxiliary Elements. */
  564.         private final FieldAuxiliaryElements<T> auxiliaryElements;

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

  567.         /**
  568.          * Build a new instance with a new field.
  569.          * @param state      current state information: date, kinematics, attitude
  570.          * @param meanMode   if true return the value associated to the mean elements
  571.          *                   variation, if false return the values associated to the
  572.          *                   short periodic elements variation
  573.          * @param j          the j index. used only for short periodic variation.
  574.          *                   Ignored for mean elements variation.
  575.          * @param parameters values of the force model parameters (only 1 values
  576.          *                   for each parameters corresponding to state date) obtained by
  577.          *                   calling the extract parameter method {@link #extractParameters(double[], AbsoluteDate)}
  578.          *                   to selected the right value for state date or by getting the parameters for a specific date
  579.          * @param field      field utilized by default
  580.          */
  581.         public FieldIntegrableFunction(final FieldSpacecraftState<T> state, final boolean meanMode, final int j,
  582.                 final T[] parameters, final Field<T> field) {

  583.             this.meanMode = meanMode;
  584.             this.j = j;
  585.             this.parameters = parameters.clone();
  586.             this.auxiliaryElements = new FieldAuxiliaryElements<>(state.getOrbit(), I);
  587.             this.context = new FieldAbstractGaussianContributionContext<>(auxiliaryElements, this.parameters);
  588.             // remove derivatives from state
  589.             final T[] stateVector = MathArrays.buildArray(field, 6);
  590.             OrbitType.EQUINOCTIAL.mapOrbitToArray(state.getOrbit(), PositionAngleType.TRUE, stateVector, null);
  591.             final FieldOrbit<T> fixedOrbit = OrbitType.EQUINOCTIAL.mapArrayToOrbit(stateVector, null,
  592.                     PositionAngleType.TRUE, state.getDate(), context.getMu(), state.getFrame());
  593.             this.state = new FieldSpacecraftState<>(fixedOrbit, state.getAttitude(), state.getMass());
  594.         }

  595.         /** {@inheritDoc} */
  596.         @Override
  597.         public T[] value(final T x) {

  598.             // Parameters for array building
  599.             final Field<T> field = auxiliaryElements.getDate().getField();
  600.             final int dimension = 6;

  601.             // Compute the time difference from the true longitude difference
  602.             final T shiftedLm = trueToMean(x);
  603.             final T dLm = shiftedLm.subtract(auxiliaryElements.getLM());
  604.             final T dt = dLm.divide(context.getMeanMotion());

  605.             final FieldSinCos<T> scL = FastMath.sinCos(x);
  606.             final T cosL = scL.cos();
  607.             final T sinL = scL.sin();
  608.             final T roa  = auxiliaryElements.getB().multiply(auxiliaryElements.getB()).divide(auxiliaryElements.getH().multiply(sinL).add(auxiliaryElements.getK().multiply(cosL)).add(1.));
  609.             final T roa2 = roa.multiply(roa);
  610.             final T r = auxiliaryElements.getSma().multiply(roa);
  611.             final T X = r.multiply(cosL);
  612.             final T Y = r.multiply(sinL);
  613.             final T naob = context.getMeanMotion().multiply(auxiliaryElements.getSma())
  614.                     .divide(auxiliaryElements.getB());
  615.             final T Xdot = naob.multiply(auxiliaryElements.getH().add(sinL)).negate();
  616.             final T Ydot = naob.multiply(auxiliaryElements.getK().add(cosL));
  617.             final FieldVector3D<T> vel = new FieldVector3D<>(Xdot, auxiliaryElements.getVectorF(), Ydot,
  618.                     auxiliaryElements.getVectorG());

  619.             // Compute acceleration
  620.             FieldVector3D<T> acc = FieldVector3D.getZero(field);

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

  623.             // Recompose an orbit with time held fixed to be compliant with DSST theory
  624.             final FieldOrbit<T> recomposedOrbit = new FieldEquinoctialOrbit<>(shiftedOrbit.getA(),
  625.                     shiftedOrbit.getEquinoctialEx(), shiftedOrbit.getEquinoctialEy(), shiftedOrbit.getHx(),
  626.                     shiftedOrbit.getHy(), shiftedOrbit.getLv(), PositionAngleType.TRUE, shiftedOrbit.getFrame(),
  627.                     state.getDate(), context.getMu());

  628.             // Get the corresponding attitude
  629.             final FieldAttitude<T> recomposedAttitude = attitudeProvider.getAttitude(recomposedOrbit,
  630.                     recomposedOrbit.getDate(), recomposedOrbit.getFrame());

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

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

  635.             // Compute the derivatives of the elements by the speed
  636.             final T[] deriv = MathArrays.buildArray(field, dimension);
  637.             // da/dv
  638.             deriv[0] = getAoV(vel).dotProduct(acc);
  639.             // dex/dv
  640.             deriv[1] = getKoV(X, Y, Xdot, Ydot).dotProduct(acc);
  641.             // dey/dv
  642.             deriv[2] = getHoV(X, Y, Xdot, Ydot).dotProduct(acc);
  643.             // dhx/dv
  644.             deriv[3] = getQoV(X).dotProduct(acc);
  645.             // dhy/dv
  646.             deriv[4] = getPoV(Y).dotProduct(acc);
  647.             // dλ/dv
  648.             deriv[5] = getLoV(X, Y, Xdot, Ydot).dotProduct(acc);

  649.             // Compute mean elements rates
  650.             T[] val = null;
  651.             if (meanMode) {
  652.                 val = MathArrays.buildArray(field, dimension);
  653.                 for (int i = 0; i < 6; i++) {
  654.                     // da<sub>i</sub>/dt
  655.                     val[i] = deriv[i].multiply(roa2);
  656.                 }
  657.             } else {
  658.                 val = MathArrays.buildArray(field, dimension * 2);
  659.                 //Compute cos(j*L) and sin(j*L);
  660.                 final FieldSinCos<T> scjL = FastMath.sinCos(x.multiply(j));
  661.                 final T cosjL = j == 1 ? cosL : scjL.cos();
  662.                 final T sinjL = j == 1 ? sinL : scjL.sin();

  663.                 for (int i = 0; i < 6; i++) {
  664.                     // da<sub>i</sub>/dv * cos(jL)
  665.                     val[i] = deriv[i].multiply(cosjL);
  666.                     // da<sub>i</sub>/dv * sin(jL)
  667.                     val[i + 6] = deriv[i].multiply(sinjL);
  668.                 }
  669.             }

  670.             return val;
  671.         }

  672.         /**
  673.          * Converts true longitude to mean longitude.
  674.          * @param x True longitude
  675.          * @return Eccentric longitude
  676.          */
  677.         private T trueToMean(final T x) {
  678.             return eccentricToMean(trueToEccentric(x));
  679.         }

  680.         /**
  681.          * Converts true longitude to eccentric longitude.
  682.          * @param lv True longitude
  683.          * @return Eccentric longitude
  684.          */
  685.         private T trueToEccentric (final T lv) {
  686.             final FieldSinCos<T> sclV = FastMath.sinCos(lv);
  687.             final T cosLv   = sclV.cos();
  688.             final T sinLv   = sclV.sin();
  689.             final T num     = auxiliaryElements.getH().multiply(cosLv).subtract(auxiliaryElements.getK().multiply(sinLv));
  690.             final T den     = auxiliaryElements.getB().add(auxiliaryElements.getK().multiply(cosLv)).add(auxiliaryElements.getH().multiply(sinLv)).add(1.);
  691.             return FastMath.atan(num.divide(den)).multiply(2.).add(lv);
  692.         }

  693.         /**
  694.          * Converts eccentric longitude to mean longitude.
  695.          * @param le Eccentric longitude
  696.          * @return Mean longitude
  697.          */
  698.         private T eccentricToMean (final T le) {
  699.             final FieldSinCos<T> scle = FastMath.sinCos(le);
  700.             return le.subtract(auxiliaryElements.getK().multiply(scle.sin())).add(auxiliaryElements.getH().multiply(scle.cos()));
  701.         }

  702.         /**
  703.          * Compute δa/δv.
  704.          * @param vel satellite velocity
  705.          * @return δa/δv
  706.          */
  707.         private FieldVector3D<T> getAoV(final FieldVector3D<T> vel) {
  708.             return new FieldVector3D<>(context.getTon2a(), vel);
  709.         }

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

  731.         /**
  732.          * Compute δk/δv.
  733.          * @param X    satellite position component along f, equinoctial reference frame
  734.          *             1st vector
  735.          * @param Y    satellite position component along g, equinoctial reference frame
  736.          *             2nd vector
  737.          * @param Xdot satellite velocity component along f, equinoctial reference frame
  738.          *             1st vector
  739.          * @param Ydot satellite velocity component along g, equinoctial reference frame
  740.          *             2nd vector
  741.          * @return δk/δv
  742.          */
  743.         private FieldVector3D<T> getKoV(final T X, final T Y, final T Xdot, final T Ydot) {
  744.             final T kf = Y.multiply(Ydot).multiply(context.getOoMU());
  745.             final T kg = (X.multiply(Ydot).multiply(2.).subtract(Xdot.multiply(Y))).multiply(context.getOoMU());
  746.             final T kw = auxiliaryElements.getH().multiply(
  747.                     auxiliaryElements.getQ().multiply(Y).multiply(I).subtract(auxiliaryElements.getP().multiply(X)))
  748.                     .multiply(context.getOOAB());
  749.             return new FieldVector3D<>(kf.negate(), auxiliaryElements.getVectorF(), kg, auxiliaryElements.getVectorG(),
  750.                     kw.negate(), auxiliaryElements.getVectorW());
  751.         }

  752.         /**
  753.          * Compute δp/δv.
  754.          * @param Y satellite position component along g, equinoctial reference frame
  755.          *          2nd vector
  756.          * @return δp/δv
  757.          */
  758.         private FieldVector3D<T> getPoV(final T Y) {
  759.             return new FieldVector3D<>(context.getCo2AB().multiply(Y), auxiliaryElements.getVectorW());
  760.         }

  761.         /**
  762.          * Compute δq/δv.
  763.          * @param X satellite position component along f, equinoctial reference frame
  764.          *          1st vector
  765.          * @return δq/δv
  766.          */
  767.         private FieldVector3D<T> getQoV(final T X) {
  768.             return new FieldVector3D<>(context.getCo2AB().multiply(X).multiply(I), auxiliaryElements.getVectorW());
  769.         }

  770.         /**
  771.          * Compute δλ/δv.
  772.          * @param X    satellite position component along f, equinoctial reference frame
  773.          *             1st vector
  774.          * @param Y    satellite position component along g, equinoctial reference frame
  775.          *             2nd vector
  776.          * @param Xdot satellite velocity component along f, equinoctial reference frame
  777.          *             1st vector
  778.          * @param Ydot satellite velocity component along g, equinoctial reference frame
  779.          *             2nd vector
  780.          * @return δλ/δv
  781.          */
  782.         private FieldVector3D<T> getLoV(final T X, final T Y, final T Xdot, final T Ydot) {
  783.             final FieldVector3D<T> pos = new FieldVector3D<>(X, auxiliaryElements.getVectorF(), Y,
  784.                     auxiliaryElements.getVectorG());
  785.             final FieldVector3D<T> v2 = new FieldVector3D<>(auxiliaryElements.getK(), getHoV(X, Y, Xdot, Ydot),
  786.                     auxiliaryElements.getH().negate(), getKoV(X, Y, Xdot, Ydot));
  787.             return new FieldVector3D<>(context.getOOA().multiply(-2.), pos, context.getOoBpo(), v2,
  788.                     context.getOOA().multiply(auxiliaryElements.getQ().multiply(Y).multiply(I)
  789.                             .subtract(auxiliaryElements.getP().multiply(X))),
  790.                     auxiliaryElements.getVectorW());
  791.         }

  792.     }

  793.     /** Internal class for numerical quadrature. */
  794.     protected class IntegrableFunction implements UnivariateVectorFunction {

  795.         /** Current state. */
  796.         private final SpacecraftState state;

  797.         /**
  798.          * Signal that this class is used to compute the values required by the mean
  799.          * element variations or by the short periodic element variations.
  800.          */
  801.         private final boolean meanMode;

  802.         /**
  803.          * The j index.
  804.          * <p>
  805.          * Used only for short periodic variation. Ignored for mean elements variation.
  806.          * </p>
  807.          */
  808.         private final int j;

  809.         /** Container for attributes. */
  810.         private final AbstractGaussianContributionContext context;

  811.         /** Auxiliary Elements. */
  812.         private final AuxiliaryElements auxiliaryElements;

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

  815.         /**
  816.          * Build a new instance.
  817.          * @param state      current state information: date, kinematics, attitude
  818.          * @param meanMode   if true return the value associated to the mean elements
  819.          *                   variation, if false return the values associated to the
  820.          *                   short periodic elements variation
  821.          * @param j          the j index. used only for short periodic variation.
  822.          *                   Ignored for mean elements variation.
  823.          * @param parameters list of the estimated values for each driver at state date of the force model parameters
  824.          *                   only 1 value for each parameter
  825.          */
  826.         IntegrableFunction(final SpacecraftState state, final boolean meanMode, final int j,
  827.                 final double[] parameters) {

  828.             this.meanMode = meanMode;
  829.             this.j = j;
  830.             this.parameters = parameters.clone();
  831.             this.auxiliaryElements = new AuxiliaryElements(state.getOrbit(), I);
  832.             this.context = new AbstractGaussianContributionContext(auxiliaryElements, this.parameters);
  833.             // remove derivatives from state
  834.             final double[] stateVector = new double[6];
  835.             OrbitType.EQUINOCTIAL.mapOrbitToArray(state.getOrbit(), PositionAngleType.TRUE, stateVector, null);
  836.             final Orbit fixedOrbit = OrbitType.EQUINOCTIAL.mapArrayToOrbit(stateVector, null, PositionAngleType.TRUE,
  837.                     state.getDate(), context.getMu(), state.getFrame());
  838.             this.state = new SpacecraftState(fixedOrbit, state.getAttitude(), state.getMass());
  839.         }

  840.         /** {@inheritDoc} */
  841.         @Override
  842.         public double[] value(final double x) {

  843.             // Compute the time difference from the true longitude difference
  844.             final double shiftedLm = trueToMean(x);
  845.             final double dLm = shiftedLm - auxiliaryElements.getLM();
  846.             final double dt = dLm / context.getMeanMotion();

  847.             final SinCos scL  = FastMath.sinCos(x);
  848.             final double cosL = scL.cos();
  849.             final double sinL = scL.sin();
  850.             final double roa  = auxiliaryElements.getB() * auxiliaryElements.getB() / (1. + auxiliaryElements.getH() * sinL + auxiliaryElements.getK() * cosL);
  851.             final double roa2 = roa * roa;
  852.             final double r = auxiliaryElements.getSma() * roa;
  853.             final double X = r * cosL;
  854.             final double Y = r * sinL;
  855.             final double naob = context.getMeanMotion() * auxiliaryElements.getSma() / auxiliaryElements.getB();
  856.             final double Xdot = -naob * (auxiliaryElements.getH() + sinL);
  857.             final double Ydot = naob * (auxiliaryElements.getK() + cosL);
  858.             final Vector3D vel = new Vector3D(Xdot, auxiliaryElements.getVectorF(), Ydot,
  859.                     auxiliaryElements.getVectorG());

  860.             // Compute acceleration
  861.             Vector3D acc = Vector3D.ZERO;

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

  864.             // Recompose an orbit with time held fixed to be compliant with DSST theory
  865.             final Orbit recomposedOrbit = new EquinoctialOrbit(shiftedOrbit.getA(), shiftedOrbit.getEquinoctialEx(),
  866.                     shiftedOrbit.getEquinoctialEy(), shiftedOrbit.getHx(), shiftedOrbit.getHy(), shiftedOrbit.getLv(),
  867.                     PositionAngleType.TRUE, shiftedOrbit.getFrame(), state.getDate(), context.getMu());

  868.             // Get the corresponding attitude
  869.             final Attitude recomposedAttitude = attitudeProvider.getAttitude(recomposedOrbit, recomposedOrbit.getDate(),
  870.                     recomposedOrbit.getFrame());

  871.             // create shifted SpacecraftState with attitude at specified time
  872.             final SpacecraftState shiftedState = new SpacecraftState(recomposedOrbit, recomposedAttitude,
  873.                     state.getMass());

  874.             // here parameters is a list of all span values of each parameter driver
  875.             acc = contribution.acceleration(shiftedState, parameters);

  876.             // Compute the derivatives of the elements by the speed
  877.             final double[] deriv = new double[6];
  878.             // da/dv
  879.             deriv[0] = getAoV(vel).dotProduct(acc);
  880.             // dex/dv
  881.             deriv[1] = getKoV(X, Y, Xdot, Ydot).dotProduct(acc);
  882.             // dey/dv
  883.             deriv[2] = getHoV(X, Y, Xdot, Ydot).dotProduct(acc);
  884.             // dhx/dv
  885.             deriv[3] = getQoV(X).dotProduct(acc);
  886.             // dhy/dv
  887.             deriv[4] = getPoV(Y).dotProduct(acc);
  888.             // dλ/dv
  889.             deriv[5] = getLoV(X, Y, Xdot, Ydot).dotProduct(acc);

  890.             // Compute mean elements rates
  891.             double[] val = null;
  892.             if (meanMode) {
  893.                 val = new double[6];
  894.                 for (int i = 0; i < 6; i++) {
  895.                     // da<sub>i</sub>/dt
  896.                     val[i] = roa2 * deriv[i];
  897.                 }
  898.             } else {
  899.                 val = new double[12];
  900.                 //Compute cos(j*L) and sin(j*L);
  901.                 final SinCos scjL  = FastMath.sinCos(j * x);
  902.                 final double cosjL = j == 1 ? cosL : scjL.cos();
  903.                 final double sinjL = j == 1 ? sinL : scjL.sin();

  904.                 for (int i = 0; i < 6; i++) {
  905.                     // da<sub>i</sub>/dv * cos(jL)
  906.                     val[i] = cosjL * deriv[i];
  907.                     // da<sub>i</sub>/dv * sin(jL)
  908.                     val[i + 6] = sinjL * deriv[i];
  909.                 }
  910.             }
  911.             return val;
  912.         }

  913.         /**
  914.          * Converts true longitude to eccentric longitude.
  915.          * @param lv True longitude
  916.          * @return Eccentric longitude
  917.          */
  918.         private double trueToEccentric (final double lv) {
  919.             final SinCos scLv    = FastMath.sinCos(lv);
  920.             final double num     = auxiliaryElements.getH() * scLv.cos() - auxiliaryElements.getK() * scLv.sin();
  921.             final double den     = auxiliaryElements.getB() + 1. + auxiliaryElements.getK() * scLv.cos() + auxiliaryElements.getH() * scLv.sin();
  922.             return lv + 2. * FastMath.atan(num / den);
  923.         }

  924.         /**
  925.          * Converts eccentric longitude to mean longitude.
  926.          * @param le Eccentric longitude
  927.          * @return Mean longitude
  928.          */
  929.         private double eccentricToMean (final double le) {
  930.             final SinCos scLe = FastMath.sinCos(le);
  931.             return le - auxiliaryElements.getK() * scLe.sin() + auxiliaryElements.getH() * scLe.cos();
  932.         }

  933.         /**
  934.          * Converts true longitude to mean longitude.
  935.          * @param lv True longitude
  936.          * @return Eccentric longitude
  937.          */
  938.         private double trueToMean(final double lv) {
  939.             return eccentricToMean(trueToEccentric(lv));
  940.         }

  941.         /**
  942.          * Compute δa/δv.
  943.          * @param vel satellite velocity
  944.          * @return δa/δv
  945.          */
  946.         private Vector3D getAoV(final Vector3D vel) {
  947.             return new Vector3D(context.getTon2a(), vel);
  948.         }

  949.         /**
  950.          * Compute δh/δ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 δh/δv
  960.          */
  961.         private Vector3D getHoV(final double X, final double Y, final double Xdot, final double Ydot) {
  962.             final double kf = (2. * Xdot * Y - X * Ydot) * context.getOoMU();
  963.             final double kg = X * Xdot * context.getOoMU();
  964.             final double kw = auxiliaryElements.getK() *
  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 δk/δv.
  971.          * @param X    satellite position component along f, equinoctial reference frame
  972.          *             1st vector
  973.          * @param Y    satellite position component along g, equinoctial reference frame
  974.          *             2nd vector
  975.          * @param Xdot satellite velocity component along f, equinoctial reference frame
  976.          *             1st vector
  977.          * @param Ydot satellite velocity component along g, equinoctial reference frame
  978.          *             2nd vector
  979.          * @return δk/δv
  980.          */
  981.         private Vector3D getKoV(final double X, final double Y, final double Xdot, final double Ydot) {
  982.             final double kf = Y * Ydot * context.getOoMU();
  983.             final double kg = (2. * X * Ydot - Xdot * Y) * context.getOoMU();
  984.             final double kw = auxiliaryElements.getH() *
  985.                     (I * auxiliaryElements.getQ() * Y - auxiliaryElements.getP() * X) * context.getOOAB();
  986.             return new Vector3D(-kf, auxiliaryElements.getVectorF(), kg, auxiliaryElements.getVectorG(), -kw,
  987.                     auxiliaryElements.getVectorW());
  988.         }

  989.         /**
  990.          * Compute δp/δv.
  991.          * @param Y satellite position component along g, equinoctial reference frame
  992.          *          2nd vector
  993.          * @return δp/δv
  994.          */
  995.         private Vector3D getPoV(final double Y) {
  996.             return new Vector3D(context.getCo2AB() * Y, auxiliaryElements.getVectorW());
  997.         }

  998.         /**
  999.          * Compute δq/δv.
  1000.          * @param X satellite position component along f, equinoctial reference frame
  1001.          *          1st vector
  1002.          * @return δq/δv
  1003.          */
  1004.         private Vector3D getQoV(final double X) {
  1005.             return new Vector3D(I * context.getCo2AB() * X, auxiliaryElements.getVectorW());
  1006.         }

  1007.         /**
  1008.          * Compute δλ/δv.
  1009.          * @param X    satellite position component along f, equinoctial reference frame
  1010.          *             1st vector
  1011.          * @param Y    satellite position component along g, equinoctial reference frame
  1012.          *             2nd vector
  1013.          * @param Xdot satellite velocity component along f, equinoctial reference frame
  1014.          *             1st vector
  1015.          * @param Ydot satellite velocity component along g, equinoctial reference frame
  1016.          *             2nd vector
  1017.          * @return δλ/δv
  1018.          */
  1019.         private Vector3D getLoV(final double X, final double Y, final double Xdot, final double Ydot) {
  1020.             final Vector3D pos = new Vector3D(X, auxiliaryElements.getVectorF(), Y, auxiliaryElements.getVectorG());
  1021.             final Vector3D v2 = new Vector3D(auxiliaryElements.getK(), getHoV(X, Y, Xdot, Ydot),
  1022.                     -auxiliaryElements.getH(), getKoV(X, Y, Xdot, Ydot));
  1023.             return new Vector3D(-2. * context.getOOA(), pos, context.getOoBpo(), v2,
  1024.                     (I * auxiliaryElements.getQ() * Y - auxiliaryElements.getP() * X) * context.getOOA(),
  1025.                     auxiliaryElements.getVectorW());
  1026.         }

  1027.     }

  1028.     /**
  1029.      * Class used to {@link #integrate(UnivariateVectorFunction, double, double)
  1030.      * integrate} a {@link org.hipparchus.analysis.UnivariateVectorFunction
  1031.      * function} of the orbital elements using the Gaussian quadrature rule to get
  1032.      * the acceleration.
  1033.      */
  1034.     protected static class GaussQuadrature {

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

  1036.         /** Points for quadrature of order 12. */
  1037.         private static final double[] P_12 = { -0.98156063424671910000, -0.90411725637047490000,
  1038.             -0.76990267419430470000, -0.58731795428661740000, -0.36783149899818024000, -0.12523340851146890000,
  1039.             0.12523340851146890000, 0.36783149899818024000, 0.58731795428661740000, 0.76990267419430470000,
  1040.             0.90411725637047490000, 0.98156063424671910000 };

  1041.         /** Weights for quadrature of order 12. */
  1042.         private static final double[] W_12 = { 0.04717533638651220000, 0.10693932599531830000, 0.16007832854334633000,
  1043.             0.20316742672306584000, 0.23349253653835478000, 0.24914704581340286000, 0.24914704581340286000,
  1044.             0.23349253653835478000, 0.20316742672306584000, 0.16007832854334633000, 0.10693932599531830000,
  1045.             0.04717533638651220000 };

  1046.         /** Points for quadrature of order 16. */
  1047.         private static final double[] P_16 = { -0.98940093499164990000, -0.94457502307323260000,
  1048.             -0.86563120238783160000, -0.75540440835500310000, -0.61787624440264380000, -0.45801677765722737000,
  1049.             -0.28160355077925890000, -0.09501250983763745000, 0.09501250983763745000, 0.28160355077925890000,
  1050.             0.45801677765722737000, 0.61787624440264380000, 0.75540440835500310000, 0.86563120238783160000,
  1051.             0.94457502307323260000, 0.98940093499164990000 };

  1052.         /** Weights for quadrature of order 16. */
  1053.         private static final double[] W_16 = { 0.02715245941175405800, 0.06225352393864777000, 0.09515851168249283000,
  1054.             0.12462897125553388000, 0.14959598881657685000, 0.16915651939500256000, 0.18260341504492360000,
  1055.             0.18945061045506847000, 0.18945061045506847000, 0.18260341504492360000, 0.16915651939500256000,
  1056.             0.14959598881657685000, 0.12462897125553388000, 0.09515851168249283000, 0.06225352393864777000,
  1057.             0.02715245941175405800 };

  1058.         /** Points for quadrature of order 20. */
  1059.         private static final double[] P_20 = { -0.99312859918509490000, -0.96397192727791390000,
  1060.             -0.91223442825132600000, -0.83911697182221890000, -0.74633190646015080000, -0.63605368072651510000,
  1061.             -0.51086700195082700000, -0.37370608871541955000, -0.22778585114164507000, -0.07652652113349734000,
  1062.             0.07652652113349734000, 0.22778585114164507000, 0.37370608871541955000, 0.51086700195082700000,
  1063.             0.63605368072651510000, 0.74633190646015080000, 0.83911697182221890000, 0.91223442825132600000,
  1064.             0.96397192727791390000, 0.99312859918509490000 };

  1065.         /** Weights for quadrature of order 20. */
  1066.         private static final double[] W_20 = { 0.01761400713915226400, 0.04060142980038684000, 0.06267204833410904000,
  1067.             0.08327674157670477000, 0.10193011981724048000, 0.11819453196151844000, 0.13168863844917678000,
  1068.             0.14209610931838212000, 0.14917298647260380000, 0.15275338713072600000, 0.15275338713072600000,
  1069.             0.14917298647260380000, 0.14209610931838212000, 0.13168863844917678000, 0.11819453196151844000,
  1070.             0.10193011981724048000, 0.08327674157670477000, 0.06267204833410904000, 0.04060142980038684000,
  1071.             0.01761400713915226400 };

  1072.         /** Points for quadrature of order 24. */
  1073.         private static final double[] P_24 = { -0.99518721999702130000, -0.97472855597130950000,
  1074.             -0.93827455200273270000, -0.88641552700440100000, -0.82000198597390300000, -0.74012419157855440000,
  1075.             -0.64809365193697550000, -0.54542147138883950000, -0.43379350762604520000, -0.31504267969616340000,
  1076.             -0.19111886747361634000, -0.06405689286260563000, 0.06405689286260563000, 0.19111886747361634000,
  1077.             0.31504267969616340000, 0.43379350762604520000, 0.54542147138883950000, 0.64809365193697550000,
  1078.             0.74012419157855440000, 0.82000198597390300000, 0.88641552700440100000, 0.93827455200273270000,
  1079.             0.97472855597130950000, 0.99518721999702130000 };

  1080.         /** Weights for quadrature of order 24. */
  1081.         private static final double[] W_24 = { 0.01234122979998733500, 0.02853138862893380600, 0.04427743881741981000,
  1082.             0.05929858491543691500, 0.07334648141108027000, 0.08619016153195320000, 0.09761865210411391000,
  1083.             0.10744427011596558000, 0.11550566805372553000, 0.12167047292780335000, 0.12583745634682825000,
  1084.             0.12793819534675221000, 0.12793819534675221000, 0.12583745634682825000, 0.12167047292780335000,
  1085.             0.11550566805372553000, 0.10744427011596558000, 0.09761865210411391000, 0.08619016153195320000,
  1086.             0.07334648141108027000, 0.05929858491543691500, 0.04427743881741981000, 0.02853138862893380600,
  1087.             0.01234122979998733500 };

  1088.         /** Points for quadrature of order 32. */
  1089.         private static final double[] P_32 = { -0.99726386184948160000, -0.98561151154526840000,
  1090.             -0.96476225558750640000, -0.93490607593773970000, -0.89632115576605220000, -0.84936761373256990000,
  1091.             -0.79448379596794250000, -0.73218211874028970000, -0.66304426693021520000, -0.58771575724076230000,
  1092.             -0.50689990893222950000, -0.42135127613063540000, -0.33186860228212767000, -0.23928736225213710000,
  1093.             -0.14447196158279646000, -0.04830766568773831000, 0.04830766568773831000, 0.14447196158279646000,
  1094.             0.23928736225213710000, 0.33186860228212767000, 0.42135127613063540000, 0.50689990893222950000,
  1095.             0.58771575724076230000, 0.66304426693021520000, 0.73218211874028970000, 0.79448379596794250000,
  1096.             0.84936761373256990000, 0.89632115576605220000, 0.93490607593773970000, 0.96476225558750640000,
  1097.             0.98561151154526840000, 0.99726386184948160000 };

  1098.         /** Weights for quadrature of order 32. */
  1099.         private static final double[] W_32 = { 0.00701861000947013600, 0.01627439473090571200, 0.02539206530926214200,
  1100.             0.03427386291302141000, 0.04283589802222658600, 0.05099805926237621600, 0.05868409347853559000,
  1101.             0.06582222277636193000, 0.07234579410884862000, 0.07819389578707042000, 0.08331192422694673000,
  1102.             0.08765209300440380000, 0.09117387869576390000, 0.09384439908080441000, 0.09563872007927487000,
  1103.             0.09654008851472784000, 0.09654008851472784000, 0.09563872007927487000, 0.09384439908080441000,
  1104.             0.09117387869576390000, 0.08765209300440380000, 0.08331192422694673000, 0.07819389578707042000,
  1105.             0.07234579410884862000, 0.06582222277636193000, 0.05868409347853559000, 0.05099805926237621600,
  1106.             0.04283589802222658600, 0.03427386291302141000, 0.02539206530926214200, 0.01627439473090571200,
  1107.             0.00701861000947013600 };

  1108.         /** Points for quadrature of order 40. */
  1109.         private static final double[] P_40 = { -0.99823770971055930000, -0.99072623869945710000,
  1110.             -0.97725994998377420000, -0.95791681921379170000, -0.93281280827867660000, -0.90209880696887420000,
  1111.             -0.86595950321225960000, -0.82461223083331170000, -0.77830565142651940000, -0.72731825518992710000,
  1112.             -0.67195668461417960000, -0.61255388966798030000, -0.54946712509512820000, -0.48307580168617870000,
  1113.             -0.41377920437160500000, -0.34199409082575850000, -0.26815218500725370000, -0.19269758070137110000,
  1114.             -0.11608407067525522000, -0.03877241750605081600, 0.03877241750605081600, 0.11608407067525522000,
  1115.             0.19269758070137110000, 0.26815218500725370000, 0.34199409082575850000, 0.41377920437160500000,
  1116.             0.48307580168617870000, 0.54946712509512820000, 0.61255388966798030000, 0.67195668461417960000,
  1117.             0.72731825518992710000, 0.77830565142651940000, 0.82461223083331170000, 0.86595950321225960000,
  1118.             0.90209880696887420000, 0.93281280827867660000, 0.95791681921379170000, 0.97725994998377420000,
  1119.             0.99072623869945710000, 0.99823770971055930000 };

  1120.         /** Weights for quadrature of order 40. */
  1121.         private static final double[] W_40 = { 0.00452127709853309800, 0.01049828453115270400, 0.01642105838190797300,
  1122.             0.02224584919416689000, 0.02793700698002338000, 0.03346019528254786500, 0.03878216797447199000,
  1123.             0.04387090818567333000, 0.04869580763507221000, 0.05322784698393679000, 0.05743976909939157000,
  1124.             0.06130624249292891000, 0.06480401345660108000, 0.06791204581523394000, 0.07061164739128681000,
  1125.             0.07288658239580408000, 0.07472316905796833000, 0.07611036190062619000, 0.07703981816424793000,
  1126.             0.07750594797842482000, 0.07750594797842482000, 0.07703981816424793000, 0.07611036190062619000,
  1127.             0.07472316905796833000, 0.07288658239580408000, 0.07061164739128681000, 0.06791204581523394000,
  1128.             0.06480401345660108000, 0.06130624249292891000, 0.05743976909939157000, 0.05322784698393679000,
  1129.             0.04869580763507221000, 0.04387090818567333000, 0.03878216797447199000, 0.03346019528254786500,
  1130.             0.02793700698002338000, 0.02224584919416689000, 0.01642105838190797300, 0.01049828453115270400,
  1131.             0.00452127709853309800 };

  1132.         /** Points for quadrature of order 48. */
  1133.         private static final double[] P_48 = { -0.99877100725242610000, -0.99353017226635080000,
  1134.             -0.98412458372282700000, -0.97059159254624720000, -0.95298770316043080000, -0.93138669070655440000,
  1135.             -0.90587913671556960000, -0.87657202027424800000, -0.84358826162439350000, -0.80706620402944250000,
  1136.             -0.76715903251574020000, -0.72403413092381470000, -0.67787237963266400000, -0.62886739677651370000,
  1137.             -0.57722472608397270000, -0.52316097472223300000, -0.46690290475095840000, -0.40868648199071680000,
  1138.             -0.34875588629216070000, -0.28736248735545555000, -0.22476379039468908000, -0.16122235606889174000,
  1139.             -0.09700469920946270000, -0.03238017096286937000, 0.03238017096286937000, 0.09700469920946270000,
  1140.             0.16122235606889174000, 0.22476379039468908000, 0.28736248735545555000, 0.34875588629216070000,
  1141.             0.40868648199071680000, 0.46690290475095840000, 0.52316097472223300000, 0.57722472608397270000,
  1142.             0.62886739677651370000, 0.67787237963266400000, 0.72403413092381470000, 0.76715903251574020000,
  1143.             0.80706620402944250000, 0.84358826162439350000, 0.87657202027424800000, 0.90587913671556960000,
  1144.             0.93138669070655440000, 0.95298770316043080000, 0.97059159254624720000, 0.98412458372282700000,
  1145.             0.99353017226635080000, 0.99877100725242610000 };

  1146.         /** Weights for quadrature of order 48. */
  1147.         private static final double[] W_48 = { 0.00315334605230596250, 0.00732755390127620800, 0.01147723457923446900,
  1148.             0.01557931572294386600, 0.01961616045735556700, 0.02357076083932435600, 0.02742650970835688000,
  1149.             0.03116722783279807000, 0.03477722256477045000, 0.03824135106583080600, 0.04154508294346483000,
  1150.             0.04467456085669424000, 0.04761665849249054000, 0.05035903555385448000, 0.05289018948519365000,
  1151.             0.05519950369998416500, 0.05727729210040315000, 0.05911483969839566000, 0.06070443916589384000,
  1152.             0.06203942315989268000, 0.06311419228625403000, 0.06392423858464817000, 0.06446616443595010000,
  1153.             0.06473769681268386000, 0.06473769681268386000, 0.06446616443595010000, 0.06392423858464817000,
  1154.             0.06311419228625403000, 0.06203942315989268000, 0.06070443916589384000, 0.05911483969839566000,
  1155.             0.05727729210040315000, 0.05519950369998416500, 0.05289018948519365000, 0.05035903555385448000,
  1156.             0.04761665849249054000, 0.04467456085669424000, 0.04154508294346483000, 0.03824135106583080600,
  1157.             0.03477722256477045000, 0.03116722783279807000, 0.02742650970835688000, 0.02357076083932435600,
  1158.             0.01961616045735556700, 0.01557931572294386600, 0.01147723457923446900, 0.00732755390127620800,
  1159.             0.00315334605230596250 };

  1160.         /** Node points. */
  1161.         private final double[] nodePoints;

  1162.         /** Node weights. */
  1163.         private final double[] nodeWeights;

  1164.         /** Number of points. */
  1165.         private final int numberOfPoints;

  1166.         /**
  1167.          * Creates a Gauss integrator of the given order.
  1168.          *
  1169.          * @param numberOfPoints Order of the integration rule.
  1170.          */
  1171.         GaussQuadrature(final int numberOfPoints) {

  1172.             this.numberOfPoints = numberOfPoints;

  1173.             switch (numberOfPoints) {
  1174.                 case 12:
  1175.                     this.nodePoints = P_12.clone();
  1176.                     this.nodeWeights = W_12.clone();
  1177.                     break;
  1178.                 case 16:
  1179.                     this.nodePoints = P_16.clone();
  1180.                     this.nodeWeights = W_16.clone();
  1181.                     break;
  1182.                 case 20:
  1183.                     this.nodePoints = P_20.clone();
  1184.                     this.nodeWeights = W_20.clone();
  1185.                     break;
  1186.                 case 24:
  1187.                     this.nodePoints = P_24.clone();
  1188.                     this.nodeWeights = W_24.clone();
  1189.                     break;
  1190.                 case 32:
  1191.                     this.nodePoints = P_32.clone();
  1192.                     this.nodeWeights = W_32.clone();
  1193.                     break;
  1194.                 case 40:
  1195.                     this.nodePoints = P_40.clone();
  1196.                     this.nodeWeights = W_40.clone();
  1197.                     break;
  1198.                 case 48:
  1199.                 default:
  1200.                     this.nodePoints = P_48.clone();
  1201.                     this.nodeWeights = W_48.clone();
  1202.                     break;
  1203.             }

  1204.         }

  1205.         /**
  1206.          * Integrates a given function on the given interval.
  1207.          *
  1208.          * @param f          Function to integrate.
  1209.          * @param lowerBound Lower bound of the integration interval.
  1210.          * @param upperBound Upper bound of the integration interval.
  1211.          * @return the integral of the weighted function.
  1212.          */
  1213.         public double[] integrate(final UnivariateVectorFunction f, final double lowerBound, final double upperBound) {

  1214.             final double[] adaptedPoints = nodePoints.clone();
  1215.             final double[] adaptedWeights = nodeWeights.clone();
  1216.             transform(adaptedPoints, adaptedWeights, lowerBound, upperBound);
  1217.             return basicIntegrate(f, adaptedPoints, adaptedWeights);
  1218.         }

  1219.         /**
  1220.          * Integrates a given function on the given interval.
  1221.          *
  1222.          * @param <T>        the type of the field elements
  1223.          * @param f          Function to integrate.
  1224.          * @param lowerBound Lower bound of the integration interval.
  1225.          * @param upperBound Upper bound of the integration interval.
  1226.          * @param field      field utilized by default
  1227.          * @return the integral of the weighted function.
  1228.          */
  1229.         public <T extends CalculusFieldElement<T>> T[] integrate(final CalculusFieldUnivariateVectorFunction<T> f,
  1230.                 final T lowerBound, final T upperBound, final Field<T> field) {

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

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

  1234.             for (int i = 0; i < numberOfPoints; i++) {
  1235.                 adaptedPoints[i] = zero.add(nodePoints[i]);
  1236.                 adaptedWeights[i] = zero.add(nodeWeights[i]);
  1237.             }

  1238.             transform(adaptedPoints, adaptedWeights, lowerBound, upperBound);
  1239.             return basicIntegrate(f, adaptedPoints, adaptedWeights, field);
  1240.         }

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

  1262.         /**
  1263.          * Performs a change of variable so that the integration can be performed on an
  1264.          * arbitrary interval {@code [a, b]}.
  1265.          * <p>
  1266.          * It is assumed that the natural interval is {@code [-1, 1]}.
  1267.          * </p>
  1268.          * @param <T>     the type of the field elements
  1269.          * @param points  Points to adapt to the new interval.
  1270.          * @param weights Weights to adapt to the new interval.
  1271.          * @param a       Lower bound of the integration interval.
  1272.          * @param b       Lower bound of the integration interval
  1273.          */
  1274.         private <T extends CalculusFieldElement<T>> void transform(final T[] points, final T[] weights, final T a,
  1275.                 final T b) {
  1276.             // Scaling
  1277.             final T scale = (b.subtract(a)).divide(2.);
  1278.             final T shift = a.add(scale);
  1279.             for (int i = 0; i < points.length; i++) {
  1280.                 points[i] = scale.multiply(points[i]).add(shift);
  1281.                 weights[i] = scale.multiply(weights[i]);
  1282.             }
  1283.         }

  1284.         /**
  1285.          * Returns an estimate of the integral of {@code f(x) * w(x)}, where {@code w}
  1286.          * is a weight function that depends on the actual flavor of the Gauss
  1287.          * integration scheme.
  1288.          *
  1289.          * @param f       Function to integrate.
  1290.          * @param points  Nodes.
  1291.          * @param weights Nodes weights.
  1292.          * @return the integral of the weighted function.
  1293.          */
  1294.         private double[] basicIntegrate(final UnivariateVectorFunction f, final double[] points,
  1295.                 final double[] weights) {
  1296.             double x = points[0];
  1297.             double w = weights[0];
  1298.             double[] v = f.value(x);
  1299.             final double[] y = new double[v.length];
  1300.             for (int j = 0; j < v.length; j++) {
  1301.                 y[j] = w * v[j];
  1302.             }
  1303.             final double[] t = y.clone();
  1304.             final double[] c = new double[v.length];
  1305.             final double[] s = t.clone();
  1306.             for (int i = 1; i < points.length; i++) {
  1307.                 x = points[i];
  1308.                 w = weights[i];
  1309.                 v = f.value(x);
  1310.                 for (int j = 0; j < v.length; j++) {
  1311.                     y[j] = w * v[j] - c[j];
  1312.                     t[j] = s[j] + y[j];
  1313.                     c[j] = (t[j] - s[j]) - y[j];
  1314.                     s[j] = t[j];
  1315.                 }
  1316.             }
  1317.             return s;
  1318.         }

  1319.         /**
  1320.          * Returns an estimate of the integral of {@code f(x) * w(x)}, where {@code w}
  1321.          * is a weight function that depends on the actual flavor of the Gauss
  1322.          * integration scheme.
  1323.          *
  1324.          * @param <T>     the type of the field elements.
  1325.          * @param f       Function to integrate.
  1326.          * @param points  Nodes.
  1327.          * @param weights Nodes weight
  1328.          * @param field   field utilized by default
  1329.          * @return the integral of the weighted function.
  1330.          */
  1331.         private <T extends CalculusFieldElement<T>> T[] basicIntegrate(final CalculusFieldUnivariateVectorFunction<T> f,
  1332.                 final T[] points, final T[] weights, final Field<T> field) {

  1333.             T x = points[0];
  1334.             T w = weights[0];
  1335.             T[] v = f.value(x);

  1336.             final T[] y = MathArrays.buildArray(field, v.length);
  1337.             for (int j = 0; j < v.length; j++) {
  1338.                 y[j] = v[j].multiply(w);
  1339.             }
  1340.             final T[] t = y.clone();
  1341.             final T[] c = MathArrays.buildArray(field, v.length);
  1342.             ;
  1343.             final T[] s = t.clone();
  1344.             for (int i = 1; i < points.length; i++) {
  1345.                 x = points[i];
  1346.                 w = weights[i];
  1347.                 v = f.value(x);
  1348.                 for (int j = 0; j < v.length; j++) {
  1349.                     y[j] = v[j].multiply(w).subtract(c[j]);
  1350.                     t[j] = y[j].add(s[j]);
  1351.                     c[j] = (t[j].subtract(s[j])).subtract(y[j]);
  1352.                     s[j] = t[j];
  1353.                 }
  1354.             }
  1355.             return s;
  1356.         }

  1357.     }

  1358.     /**
  1359.      * Compute the C<sub>i</sub><sup>j</sup> and the S<sub>i</sub><sup>j</sup>
  1360.      * coefficients.
  1361.      * <p>
  1362.      * Those coefficients are given in Danielson paper by expression 4.4-(6)
  1363.      * </p>
  1364.      * @author Petre Bazavan
  1365.      * @author Lucian Barbulescu
  1366.      */
  1367.     protected class FourierCjSjCoefficients {

  1368.         /** Maximum possible value for j. */
  1369.         private final int jMax;

  1370.         /**
  1371.          * The C<sub>i</sub><sup>j</sup> coefficients.
  1372.          * <p>
  1373.          * the index i corresponds to the following elements: <br/>
  1374.          * - 0 for a <br>
  1375.          * - 1 for k <br>
  1376.          * - 2 for h <br>
  1377.          * - 3 for q <br>
  1378.          * - 4 for p <br>
  1379.          * - 5 for λ <br>
  1380.          * </p>
  1381.          */
  1382.         private final double[][] cCoef;

  1383.         /**
  1384.          * The C<sub>i</sub><sup>j</sup> coefficients.
  1385.          * <p>
  1386.          * the index i corresponds to the following elements: <br/>
  1387.          * - 0 for a <br>
  1388.          * - 1 for k <br>
  1389.          * - 2 for h <br>
  1390.          * - 3 for q <br>
  1391.          * - 4 for p <br>
  1392.          * - 5 for λ <br>
  1393.          * </p>
  1394.          */
  1395.         private final double[][] sCoef;

  1396.         /**
  1397.          * Standard constructor.
  1398.          * @param state             the current state
  1399.          * @param jMax              maximum value for j
  1400.          * @param auxiliaryElements auxiliary elements related to the current orbit
  1401.          * @param parameters        list of parameter values at state date for each driver
  1402.          * of the force model parameters (1 value per parameter)
  1403.          */
  1404.         FourierCjSjCoefficients(final SpacecraftState state, final int jMax, final AuxiliaryElements auxiliaryElements,
  1405.                 final double[] parameters) {

  1406.             // Initialise the fields
  1407.             this.jMax = jMax;

  1408.             // Allocate the arrays
  1409.             final int rows = jMax + 1;
  1410.             cCoef = new double[rows][6];
  1411.             sCoef = new double[rows][6];

  1412.             // Compute the coefficients
  1413.             computeCoefficients(state, auxiliaryElements, parameters);
  1414.         }

  1415.         /**
  1416.          * Compute the Fourrier coefficients.
  1417.          * <p>
  1418.          * Only the C<sub>i</sub><sup>j</sup> and S<sub>i</sub><sup>j</sup> coefficients
  1419.          * need to be computed as D<sub>i</sub><sup>m</sup> is always 0.
  1420.          * </p>
  1421.          * @param state             the current state
  1422.          * @param auxiliaryElements auxiliary elements related to the current orbit
  1423.          * @param parameters        list of parameter values at state date for each driver
  1424.          * of the force model parameters (1 value per parameter)
  1425.          */
  1426.         private void computeCoefficients(final SpacecraftState state, final AuxiliaryElements auxiliaryElements,
  1427.                 final double[] parameters) {

  1428.             // Computes the limits for the integral
  1429.             final double[] ll = getLLimits(state, auxiliaryElements);
  1430.             // Computes integrated mean element rates if Llow < Lhigh
  1431.             if (ll[0] < ll[1]) {
  1432.                 // Compute 1 / PI
  1433.                 final double ooPI = 1 / FastMath.PI;

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

  1438.                     // divide by PI and set the values for the coefficients
  1439.                     for (int i = 0; i < 6; i++) {
  1440.                         cCoef[j][i] = ooPI * curentCoefficients[i];
  1441.                         sCoef[j][i] = ooPI * curentCoefficients[i + 6];
  1442.                     }
  1443.                 }
  1444.             }
  1445.         }

  1446.         /**
  1447.          * Get the coefficient C<sub>i</sub><sup>j</sup>.
  1448.          * @param i i index - corresponds to the required variation
  1449.          * @param j j index
  1450.          * @return the coefficient C<sub>i</sub><sup>j</sup>
  1451.          */
  1452.         public double getCij(final int i, final int j) {
  1453.             return cCoef[j][i];
  1454.         }

  1455.         /**
  1456.          * Get the coefficient S<sub>i</sub><sup>j</sup>.
  1457.          * @param i i index - corresponds to the required variation
  1458.          * @param j j index
  1459.          * @return the coefficient S<sub>i</sub><sup>j</sup>
  1460.          */
  1461.         public double getSij(final int i, final int j) {
  1462.             return sCoef[j][i];
  1463.         }
  1464.     }

  1465.     /**
  1466.      * Compute the C<sub>i</sub><sup>j</sup> and the S<sub>i</sub><sup>j</sup>
  1467.      * coefficients with field elements.
  1468.      * <p>
  1469.      * Those coefficients are given in Danielson paper by expression 4.4-(6)
  1470.      * </p>
  1471.      * @author Petre Bazavan
  1472.      * @author Lucian Barbulescu
  1473.      * @param <T> type of the field elements
  1474.      */
  1475.     protected class FieldFourierCjSjCoefficients<T extends CalculusFieldElement<T>> {

  1476.         /** Maximum possible value for j. */
  1477.         private final int jMax;

  1478.         /**
  1479.          * The C<sub>i</sub><sup>j</sup> coefficients.
  1480.          * <p>
  1481.          * the index i corresponds to the following elements: <br/>
  1482.          * - 0 for a <br>
  1483.          * - 1 for k <br>
  1484.          * - 2 for h <br>
  1485.          * - 3 for q <br>
  1486.          * - 4 for p <br>
  1487.          * - 5 for λ <br>
  1488.          * </p>
  1489.          */
  1490.         private final T[][] cCoef;

  1491.         /**
  1492.          * The C<sub>i</sub><sup>j</sup> coefficients.
  1493.          * <p>
  1494.          * the index i corresponds to the following elements: <br/>
  1495.          * - 0 for a <br>
  1496.          * - 1 for k <br>
  1497.          * - 2 for h <br>
  1498.          * - 3 for q <br>
  1499.          * - 4 for p <br>
  1500.          * - 5 for λ <br>
  1501.          * </p>
  1502.          */
  1503.         private final T[][] sCoef;

  1504.         /**
  1505.          * Standard constructor.
  1506.          * @param state             the current state
  1507.          * @param jMax              maximum value for j
  1508.          * @param auxiliaryElements auxiliary elements related to the current orbit
  1509.          * @param parameters        values of the force model parameters
  1510.          * @param field             field used by default
  1511.          */
  1512.         FieldFourierCjSjCoefficients(final FieldSpacecraftState<T> state, final int jMax,
  1513.                 final FieldAuxiliaryElements<T> auxiliaryElements, final T[] parameters, final Field<T> field) {
  1514.             // Initialise the fields
  1515.             this.jMax = jMax;

  1516.             // Allocate the arrays
  1517.             final int rows = jMax + 1;
  1518.             cCoef = MathArrays.buildArray(field, rows, 6);
  1519.             sCoef = MathArrays.buildArray(field, rows, 6);

  1520.             // Compute the coefficients
  1521.             computeCoefficients(state, auxiliaryElements, parameters, field);
  1522.         }

  1523.         /**
  1524.          * Compute the Fourrier coefficients.
  1525.          * <p>
  1526.          * Only the C<sub>i</sub><sup>j</sup> and S<sub>i</sub><sup>j</sup> coefficients
  1527.          * need to be computed as D<sub>i</sub><sup>m</sup> is always 0.
  1528.          * </p>
  1529.          * @param state             the current state
  1530.          * @param auxiliaryElements auxiliary elements related to the current orbit
  1531.          * @param parameters        values of the force model parameters
  1532.          * @param field             field used by default
  1533.          */
  1534.         private void computeCoefficients(final FieldSpacecraftState<T> state,
  1535.                 final FieldAuxiliaryElements<T> auxiliaryElements, final T[] parameters, final Field<T> field) {
  1536.             // Zero
  1537.             final T zero = field.getZero();
  1538.             // Computes the limits for the integral
  1539.             final T[] ll = getLLimits(state, auxiliaryElements);
  1540.             // Computes integrated mean element rates if Llow < Lhigh
  1541.             if (ll[0].getReal() < ll[1].getReal()) {
  1542.                 // Compute 1 / PI
  1543.                 final T ooPI = zero.getPi().reciprocal();

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

  1548.                     // divide by PI and set the values for the coefficients
  1549.                     for (int i = 0; i < 6; i++) {
  1550.                         cCoef[j][i] = curentCoefficients[i].multiply(ooPI);
  1551.                         sCoef[j][i] = curentCoefficients[i + 6].multiply(ooPI);
  1552.                     }
  1553.                 }
  1554.             }
  1555.         }

  1556.         /**
  1557.          * Get the coefficient C<sub>i</sub><sup>j</sup>.
  1558.          * @param i i index - corresponds to the required variation
  1559.          * @param j j index
  1560.          * @return the coefficient C<sub>i</sub><sup>j</sup>
  1561.          */
  1562.         public T getCij(final int i, final int j) {
  1563.             return cCoef[j][i];
  1564.         }

  1565.         /**
  1566.          * Get the coefficient S<sub>i</sub><sup>j</sup>.
  1567.          * @param i i index - corresponds to the required variation
  1568.          * @param j j index
  1569.          * @return the coefficient S<sub>i</sub><sup>j</sup>
  1570.          */
  1571.         public T getSij(final int i, final int j) {
  1572.             return sCoef[j][i];
  1573.         }
  1574.     }

  1575.     /**
  1576.      * This class handles the short periodic coefficients described in Danielson
  1577.      * 2.5.3-26.
  1578.      *
  1579.      * <p>
  1580.      * The value of M is 0. Also, since the values of the Fourier coefficient
  1581.      * D<sub>i</sub><sup>m</sup> is 0 then the values of the coefficients
  1582.      * D<sub>i</sub><sup>m</sup> for m &gt; 2 are also 0.
  1583.      * </p>
  1584.      * @author Petre Bazavan
  1585.      * @author Lucian Barbulescu
  1586.      *
  1587.      */
  1588.     protected static class GaussianShortPeriodicCoefficients implements ShortPeriodTerms {

  1589.         /** Maximum value for j index. */
  1590.         private final int jMax;

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

  1593.         /** Prefix for coefficients keys. */
  1594.         private final String coefficientsKeyPrefix;

  1595.         /** All coefficients slots. */
  1596.         private final transient TimeSpanMap<Slot> slots;

  1597.         /**
  1598.          * Constructor.
  1599.          * @param coefficientsKeyPrefix prefix for coefficients keys
  1600.          * @param jMax                  maximum value for j index
  1601.          * @param interpolationPoints   number of points used in the interpolation
  1602.          *                              process
  1603.          * @param slots                 all coefficients slots
  1604.          */
  1605.         GaussianShortPeriodicCoefficients(final String coefficientsKeyPrefix, final int jMax,
  1606.                 final int interpolationPoints, final TimeSpanMap<Slot> slots) {
  1607.             // Initialize fields
  1608.             this.jMax = jMax;
  1609.             this.interpolationPoints = interpolationPoints;
  1610.             this.coefficientsKeyPrefix = coefficientsKeyPrefix;
  1611.             this.slots = slots;
  1612.         }

  1613.         /**
  1614.          * Get the slot valid for some date.
  1615.          * @param meanStates mean states defining the slot
  1616.          * @return slot valid at the specified date
  1617.          */
  1618.         public Slot createSlot(final SpacecraftState... meanStates) {
  1619.             final Slot slot = new Slot(jMax, interpolationPoints);
  1620.             final AbsoluteDate first = meanStates[0].getDate();
  1621.             final AbsoluteDate last = meanStates[meanStates.length - 1].getDate();
  1622.             final int compare = first.compareTo(last);
  1623.             if (compare < 0) {
  1624.                 slots.addValidAfter(slot, first, false);
  1625.             } else if (compare > 0) {
  1626.                 slots.addValidBefore(slot, first, false);
  1627.             } else {
  1628.                 // single date, valid for all time
  1629.                 slots.addValidAfter(slot, AbsoluteDate.PAST_INFINITY, false);
  1630.             }
  1631.             return slot;
  1632.         }

  1633.         /**
  1634.          * Compute the short periodic coefficients.
  1635.          *
  1636.          * @param state       current state information: date, kinematics, attitude
  1637.          * @param slot        coefficients slot
  1638.          * @param fourierCjSj Fourier coefficients
  1639.          * @param uijvij      U and V coefficients
  1640.          * @param n           Keplerian mean motion
  1641.          * @param a           semi major axis
  1642.          */
  1643.         private void computeCoefficients(final SpacecraftState state, final Slot slot,
  1644.                 final FourierCjSjCoefficients fourierCjSj, final UijVijCoefficients uijvij, final double n,
  1645.                 final double a) {

  1646.             // get the current date
  1647.             final AbsoluteDate date = state.getDate();

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

  1650.             // 1. / n
  1651.             final double oon = 1. / n;
  1652.             // 3. / (2 * a * n)
  1653.             final double to2an = 1.5 * oon / a;
  1654.             // 3. / (4 * a * n)
  1655.             final double to4an = to2an / 2;

  1656.             // Compute the coefficients for each element
  1657.             final int size = jMax + 1;
  1658.             final double[] di1 = new double[6];
  1659.             final double[] di2 = new double[6];
  1660.             final double[][] currentCij = new double[size][6];
  1661.             final double[][] currentSij = new double[size][6];
  1662.             for (int i = 0; i < 6; i++) {

  1663.                 // compute D<sub>i</sub>¹ and D<sub>i</sub>² (all others are 0)
  1664.                 di1[i] = -oon * fourierCjSj.getCij(i, 0);
  1665.                 if (i == 5) {
  1666.                     di1[i] += to2an * uijvij.getU1(0, 0);
  1667.                 }
  1668.                 di2[i] = 0.;
  1669.                 if (i == 5) {
  1670.                     di2[i] += -to4an * fourierCjSj.getCij(0, 0);
  1671.                 }

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

  1674.                 for (int j = 1; j <= jMax; j++) {
  1675.                     // compute the current C<sub>i</sub><sup>j</sup> and S<sub>i</sub><sup>j</sup>
  1676.                     currentCij[j][i] = oon * uijvij.getU1(j, i);
  1677.                     if (i == 5) {
  1678.                         currentCij[j][i] += -to2an * uijvij.getU2(j);
  1679.                     }
  1680.                     currentSij[j][i] = oon * uijvij.getV1(j, i);
  1681.                     if (i == 5) {
  1682.                         currentSij[j][i] += -to2an * uijvij.getV2(j);
  1683.                     }

  1684.                     // add the computed coefficients to C<sub>i</sub>⁰
  1685.                     currentCij[0][i] += -(currentCij[j][i] * uijvij.currentRhoSigmaj[0][j] +
  1686.                             currentSij[j][i] * uijvij.currentRhoSigmaj[1][j]);
  1687.                 }

  1688.             }

  1689.             // add the values to the interpolators
  1690.             slot.cij[0].addGridPoint(date, currentCij[0]);
  1691.             slot.dij[1].addGridPoint(date, di1);
  1692.             slot.dij[2].addGridPoint(date, di2);
  1693.             for (int j = 1; j <= jMax; j++) {
  1694.                 slot.cij[j].addGridPoint(date, currentCij[j]);
  1695.                 slot.sij[j].addGridPoint(date, currentSij[j]);
  1696.             }

  1697.         }

  1698.         /**
  1699.          * Compute the coefficient k₂⁰ by using the equation 2.5.3-(9a) from Danielson.
  1700.          * <p>
  1701.          * After inserting 2.5.3-(8) into 2.5.3-(9a) the result becomes:<br>
  1702.          * k₂⁰ = &Sigma;<sub>k=1</sub><sup>kMax</sup>[(2 / k²) * (σ<sub>k</sub>² +
  1703.          * ρ<sub>k</sub>²)]
  1704.          * </p>
  1705.          * @param kMax             max value fot k index
  1706.          * @param currentRhoSigmaj the current computed values for the ρ<sub>j</sub> and
  1707.          *                         σ<sub>j</sub> coefficients
  1708.          * @return the coefficient k₂⁰
  1709.          */
  1710.         private double computeK20(final int kMax, final double[][] currentRhoSigmaj) {
  1711.             double k20 = 0.;

  1712.             for (int kIndex = 1; kIndex <= kMax; kIndex++) {
  1713.                 // After inserting 2.5.3-(8) into 2.5.3-(9a) the result becomes:
  1714.                 // k₂⁰ = &Sigma;<sub>k=1</sub><sup>kMax</sup>[(2 / k²) * (σ<sub>k</sub>² +
  1715.                 // ρ<sub>k</sub>²)]
  1716.                 double currentTerm = currentRhoSigmaj[1][kIndex] * currentRhoSigmaj[1][kIndex] +
  1717.                         currentRhoSigmaj[0][kIndex] * currentRhoSigmaj[0][kIndex];

  1718.                 // multiply by 2 / k²
  1719.                 currentTerm *= 2. / (kIndex * kIndex);

  1720.                 // add the term to the result
  1721.                 k20 += currentTerm;
  1722.             }

  1723.             return k20;
  1724.         }

  1725.         /** {@inheritDoc} */
  1726.         @Override
  1727.         public double[] value(final Orbit meanOrbit) {

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

  1730.             // Get the True longitude L
  1731.             final double L = meanOrbit.getLv();

  1732.             // Compute the center (l - λ)
  1733.             final double center = L - meanOrbit.getLM();
  1734.             // Compute (l - λ)²
  1735.             final double center2 = center * center;

  1736.             // Initialize short periodic variations
  1737.             final double[] shortPeriodicVariation = slot.cij[0].value(meanOrbit.getDate());
  1738.             final double[] d1 = slot.dij[1].value(meanOrbit.getDate());
  1739.             final double[] d2 = slot.dij[2].value(meanOrbit.getDate());
  1740.             for (int i = 0; i < 6; i++) {
  1741.                 shortPeriodicVariation[i] += center * d1[i] + center2 * d2[i];
  1742.             }

  1743.             for (int j = 1; j <= JMAX; j++) {
  1744.                 final double[] c = slot.cij[j].value(meanOrbit.getDate());
  1745.                 final double[] s = slot.sij[j].value(meanOrbit.getDate());
  1746.                 final SinCos sc  = FastMath.sinCos(j * L);
  1747.                 final double cos = sc.cos();
  1748.                 final double sin = sc.sin();
  1749.                 for (int i = 0; i < 6; i++) {
  1750.                     // add corresponding term to the short periodic variation
  1751.                     shortPeriodicVariation[i] += c[i] * cos;
  1752.                     shortPeriodicVariation[i] += s[i] * sin;
  1753.                 }
  1754.             }

  1755.             return shortPeriodicVariation;

  1756.         }

  1757.         /** {@inheritDoc} */
  1758.         public String getCoefficientsKeyPrefix() {
  1759.             return coefficientsKeyPrefix;
  1760.         }

  1761.         /**
  1762.          * {@inheritDoc}
  1763.          * <p>
  1764.          * For Gaussian forces, there are JMAX cj coefficients, JMAX sj coefficients and
  1765.          * 3 dj coefficients. As JMAX = 12, this sums up to 27 coefficients. The j index
  1766.          * is the integer multiplier for the true longitude argument in the cj and sj
  1767.          * coefficients and to the degree in the polynomial dj coefficients.
  1768.          * </p>
  1769.          */
  1770.         @Override
  1771.         public Map<String, double[]> getCoefficients(final AbsoluteDate date, final Set<String> selected) {

  1772.             // select the coefficients slot
  1773.             final Slot slot = slots.get(date);

  1774.             final Map<String, double[]> coefficients = new HashMap<String, double[]>(2 * JMAX + 3);
  1775.             storeIfSelected(coefficients, selected, slot.cij[0].value(date), "d", 0);
  1776.             storeIfSelected(coefficients, selected, slot.dij[1].value(date), "d", 1);
  1777.             storeIfSelected(coefficients, selected, slot.dij[2].value(date), "d", 2);
  1778.             for (int j = 1; j <= JMAX; j++) {
  1779.                 storeIfSelected(coefficients, selected, slot.cij[j].value(date), "c", j);
  1780.                 storeIfSelected(coefficients, selected, slot.sij[j].value(date), "s", j);
  1781.             }

  1782.             return coefficients;

  1783.         }

  1784.         /**
  1785.          * Put a coefficient in a map if selected.
  1786.          * @param map      map to populate
  1787.          * @param selected set of coefficients that should be put in the map (empty set
  1788.          *                 means all coefficients are selected)
  1789.          * @param value    coefficient value
  1790.          * @param id       coefficient identifier
  1791.          * @param indices  list of coefficient indices
  1792.          */
  1793.         private void storeIfSelected(final Map<String, double[]> map, final Set<String> selected, final double[] value,
  1794.                 final String id, final int... indices) {
  1795.             final StringBuilder keyBuilder = new StringBuilder(getCoefficientsKeyPrefix());
  1796.             keyBuilder.append(id);
  1797.             for (int index : indices) {
  1798.                 keyBuilder.append('[').append(index).append(']');
  1799.             }
  1800.             final String key = keyBuilder.toString();
  1801.             if (selected.isEmpty() || selected.contains(key)) {
  1802.                 map.put(key, value);
  1803.             }
  1804.         }

  1805.     }

  1806.     /**
  1807.      * This class handles the short periodic coefficients described in Danielson
  1808.      * 2.5.3-26.
  1809.      *
  1810.      * <p>
  1811.      * The value of M is 0. Also, since the values of the Fourier coefficient
  1812.      * D<sub>i</sub><sup>m</sup> is 0 then the values of the coefficients
  1813.      * D<sub>i</sub><sup>m</sup> for m &gt; 2 are also 0.
  1814.      * </p>
  1815.      * @author Petre Bazavan
  1816.      * @author Lucian Barbulescu
  1817.      * @param <T> type of the field elements
  1818.      */
  1819.     protected static class FieldGaussianShortPeriodicCoefficients<T extends CalculusFieldElement<T>>
  1820.             implements FieldShortPeriodTerms<T> {

  1821.         /** Maximum value for j index. */
  1822.         private final int jMax;

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

  1825.         /** Prefix for coefficients keys. */
  1826.         private final String coefficientsKeyPrefix;

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

  1829.         /**
  1830.          * Constructor.
  1831.          * @param coefficientsKeyPrefix prefix for coefficients keys
  1832.          * @param jMax                  maximum value for j index
  1833.          * @param interpolationPoints   number of points used in the interpolation
  1834.          *                              process
  1835.          * @param slots                 all coefficients slots
  1836.          */
  1837.         FieldGaussianShortPeriodicCoefficients(final String coefficientsKeyPrefix, final int jMax,
  1838.                 final int interpolationPoints, final FieldTimeSpanMap<FieldSlot<T>, T> slots) {
  1839.             // Initialize fields
  1840.             this.jMax = jMax;
  1841.             this.interpolationPoints = interpolationPoints;
  1842.             this.coefficientsKeyPrefix = coefficientsKeyPrefix;
  1843.             this.slots = slots;
  1844.         }

  1845.         /**
  1846.          * Get the slot valid for some date.
  1847.          * @param meanStates mean states defining the slot
  1848.          * @return slot valid at the specified date
  1849.          */
  1850.         @SuppressWarnings("unchecked")
  1851.         public FieldSlot<T> createSlot(final FieldSpacecraftState<T>... meanStates) {
  1852.             final FieldSlot<T> slot = new FieldSlot<>(jMax, interpolationPoints);
  1853.             final FieldAbsoluteDate<T> first = meanStates[0].getDate();
  1854.             final FieldAbsoluteDate<T> last = meanStates[meanStates.length - 1].getDate();
  1855.             if (first.compareTo(last) <= 0) {
  1856.                 slots.addValidAfter(slot, first);
  1857.             } else {
  1858.                 slots.addValidBefore(slot, first);
  1859.             }
  1860.             return slot;
  1861.         }

  1862.         /**
  1863.          * Compute the short periodic coefficients.
  1864.          *
  1865.          * @param state       current state information: date, kinematics, attitude
  1866.          * @param slot        coefficients slot
  1867.          * @param fourierCjSj Fourier coefficients
  1868.          * @param uijvij      U and V coefficients
  1869.          * @param n           Keplerian mean motion
  1870.          * @param a           semi major axis
  1871.          * @param field       field used by default
  1872.          */
  1873.         private void computeCoefficients(final FieldSpacecraftState<T> state, final FieldSlot<T> slot,
  1874.                 final FieldFourierCjSjCoefficients<T> fourierCjSj, final FieldUijVijCoefficients<T> uijvij, final T n,
  1875.                 final T a, final Field<T> field) {

  1876.             // Zero
  1877.             final T zero = field.getZero();

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

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

  1882.             // 1. / n
  1883.             final T oon = n.reciprocal();
  1884.             // 3. / (2 * a * n)
  1885.             final T to2an = oon.multiply(1.5).divide(a);
  1886.             // 3. / (4 * a * n)
  1887.             final T to4an = to2an.divide(2.);

  1888.             // Compute the coefficients for each element
  1889.             final int size = jMax + 1;
  1890.             final T[] di1 = MathArrays.buildArray(field, 6);
  1891.             final T[] di2 = MathArrays.buildArray(field, 6);
  1892.             final T[][] currentCij = MathArrays.buildArray(field, size, 6);
  1893.             final T[][] currentSij = MathArrays.buildArray(field, size, 6);
  1894.             for (int i = 0; i < 6; i++) {

  1895.                 // compute D<sub>i</sub>¹ and D<sub>i</sub>² (all others are 0)
  1896.                 di1[i] = oon.negate().multiply(fourierCjSj.getCij(i, 0));
  1897.                 if (i == 5) {
  1898.                     di1[i] = di1[i].add(to2an.multiply(uijvij.getU1(0, 0)));
  1899.                 }
  1900.                 di2[i] = zero;
  1901.                 if (i == 5) {
  1902.                     di2[i] = di2[i].add(to4an.negate().multiply(fourierCjSj.getCij(0, 0)));
  1903.                 }

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

  1906.                 for (int j = 1; j <= jMax; j++) {
  1907.                     // compute the current C<sub>i</sub><sup>j</sup> and S<sub>i</sub><sup>j</sup>
  1908.                     currentCij[j][i] = oon.multiply(uijvij.getU1(j, i));
  1909.                     if (i == 5) {
  1910.                         currentCij[j][i] = currentCij[j][i].add(to2an.negate().multiply(uijvij.getU2(j)));
  1911.                     }
  1912.                     currentSij[j][i] = oon.multiply(uijvij.getV1(j, i));
  1913.                     if (i == 5) {
  1914.                         currentSij[j][i] = currentSij[j][i].add(to2an.negate().multiply(uijvij.getV2(j)));
  1915.                     }

  1916.                     // add the computed coefficients to C<sub>i</sub>⁰
  1917.                     currentCij[0][i] = currentCij[0][i].add(currentCij[j][i].multiply(uijvij.currentRhoSigmaj[0][j])
  1918.                             .add(currentSij[j][i].multiply(uijvij.currentRhoSigmaj[1][j])).negate());
  1919.                 }

  1920.             }

  1921.             // add the values to the interpolators
  1922.             slot.cij[0].addGridPoint(date, currentCij[0]);
  1923.             slot.dij[1].addGridPoint(date, di1);
  1924.             slot.dij[2].addGridPoint(date, di2);
  1925.             for (int j = 1; j <= jMax; j++) {
  1926.                 slot.cij[j].addGridPoint(date, currentCij[j]);
  1927.                 slot.sij[j].addGridPoint(date, currentSij[j]);
  1928.             }

  1929.         }

  1930.         /**
  1931.          * Compute the coefficient k₂⁰ by using the equation 2.5.3-(9a) from Danielson.
  1932.          * <p>
  1933.          * After inserting 2.5.3-(8) into 2.5.3-(9a) the result becomes:<br>
  1934.          * k₂⁰ = &Sigma;<sub>k=1</sub><sup>kMax</sup>[(2 / k²) * (σ<sub>k</sub>² +
  1935.          * ρ<sub>k</sub>²)]
  1936.          * </p>
  1937.          * @param kMax             max value fot k index
  1938.          * @param currentRhoSigmaj the current computed values for the ρ<sub>j</sub> and
  1939.          *                         σ<sub>j</sub> coefficients
  1940.          * @param field            field used by default
  1941.          * @return the coefficient k₂⁰
  1942.          */
  1943.         private T computeK20(final int kMax, final T[][] currentRhoSigmaj, final Field<T> field) {
  1944.             final T zero = field.getZero();
  1945.             T k20 = zero;

  1946.             for (int kIndex = 1; kIndex <= kMax; kIndex++) {
  1947.                 // After inserting 2.5.3-(8) into 2.5.3-(9a) the result becomes:
  1948.                 // k₂⁰ = &Sigma;<sub>k=1</sub><sup>kMax</sup>[(2 / k²) * (σ<sub>k</sub>² +
  1949.                 // ρ<sub>k</sub>²)]
  1950.                 T currentTerm = currentRhoSigmaj[1][kIndex].multiply(currentRhoSigmaj[1][kIndex])
  1951.                         .add(currentRhoSigmaj[0][kIndex].multiply(currentRhoSigmaj[0][kIndex]));

  1952.                 // multiply by 2 / k²
  1953.                 currentTerm = currentTerm.multiply(2. / (kIndex * kIndex));

  1954.                 // add the term to the result
  1955.                 k20 = k20.add(currentTerm);
  1956.             }

  1957.             return k20;
  1958.         }

  1959.         /** {@inheritDoc} */
  1960.         @Override
  1961.         public T[] value(final FieldOrbit<T> meanOrbit) {

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

  1964.             // Get the True longitude L
  1965.             final T L = meanOrbit.getLv();

  1966.             // Compute the center (l - λ)
  1967.             final T center = L.subtract(meanOrbit.getLM());
  1968.             // Compute (l - λ)²
  1969.             final T center2 = center.multiply(center);

  1970.             // Initialize short periodic variations
  1971.             final T[] shortPeriodicVariation = slot.cij[0].value(meanOrbit.getDate());
  1972.             final T[] d1 = slot.dij[1].value(meanOrbit.getDate());
  1973.             final T[] d2 = slot.dij[2].value(meanOrbit.getDate());
  1974.             for (int i = 0; i < 6; i++) {
  1975.                 shortPeriodicVariation[i] = shortPeriodicVariation[i]
  1976.                         .add(center.multiply(d1[i]).add(center2.multiply(d2[i])));
  1977.             }

  1978.             for (int j = 1; j <= JMAX; j++) {
  1979.                 final T[] c = slot.cij[j].value(meanOrbit.getDate());
  1980.                 final T[] s = slot.sij[j].value(meanOrbit.getDate());
  1981.                 final FieldSinCos<T> sc = FastMath.sinCos(L.multiply(j));
  1982.                 final T cos = sc.cos();
  1983.                 final T sin = sc.sin();
  1984.                 for (int i = 0; i < 6; i++) {
  1985.                     // add corresponding term to the short periodic variation
  1986.                     shortPeriodicVariation[i] = shortPeriodicVariation[i].add(c[i].multiply(cos));
  1987.                     shortPeriodicVariation[i] = shortPeriodicVariation[i].add(s[i].multiply(sin));
  1988.                 }
  1989.             }

  1990.             return shortPeriodicVariation;

  1991.         }

  1992.         /** {@inheritDoc} */
  1993.         public String getCoefficientsKeyPrefix() {
  1994.             return coefficientsKeyPrefix;
  1995.         }

  1996.         /**
  1997.          * {@inheritDoc}
  1998.          * <p>
  1999.          * For Gaussian forces, there are JMAX cj coefficients, JMAX sj coefficients and
  2000.          * 3 dj coefficients. As JMAX = 12, this sums up to 27 coefficients. The j index
  2001.          * is the integer multiplier for the true longitude argument in the cj and sj
  2002.          * coefficients and to the degree in the polynomial dj coefficients.
  2003.          * </p>
  2004.          */
  2005.         @Override
  2006.         public Map<String, T[]> getCoefficients(final FieldAbsoluteDate<T> date, final Set<String> selected) {

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

  2009.             final Map<String, T[]> coefficients = new HashMap<String, T[]>(2 * JMAX + 3);
  2010.             storeIfSelected(coefficients, selected, slot.cij[0].value(date), "d", 0);
  2011.             storeIfSelected(coefficients, selected, slot.dij[1].value(date), "d", 1);
  2012.             storeIfSelected(coefficients, selected, slot.dij[2].value(date), "d", 2);
  2013.             for (int j = 1; j <= JMAX; j++) {
  2014.                 storeIfSelected(coefficients, selected, slot.cij[j].value(date), "c", j);
  2015.                 storeIfSelected(coefficients, selected, slot.sij[j].value(date), "s", j);
  2016.             }

  2017.             return coefficients;

  2018.         }

  2019.         /**
  2020.          * Put a coefficient in a map if selected.
  2021.          * @param map      map to populate
  2022.          * @param selected set of coefficients that should be put in the map (empty set
  2023.          *                 means all coefficients are selected)
  2024.          * @param value    coefficient value
  2025.          * @param id       coefficient identifier
  2026.          * @param indices  list of coefficient indices
  2027.          */
  2028.         private void storeIfSelected(final Map<String, T[]> map, final Set<String> selected, final T[] value,
  2029.                 final String id, final int... indices) {
  2030.             final StringBuilder keyBuilder = new StringBuilder(getCoefficientsKeyPrefix());
  2031.             keyBuilder.append(id);
  2032.             for (int index : indices) {
  2033.                 keyBuilder.append('[').append(index).append(']');
  2034.             }
  2035.             final String key = keyBuilder.toString();
  2036.             if (selected.isEmpty() || selected.contains(key)) {
  2037.                 map.put(key, value);
  2038.             }
  2039.         }

  2040.     }

  2041.     /**
  2042.      * The U<sub>i</sub><sup>j</sup> and V<sub>i</sub><sup>j</sup> coefficients
  2043.      * described by equations 2.5.3-(21) and 2.5.3-(22) from Danielson.
  2044.      * <p>
  2045.      * The index i takes only the values 1 and 2<br>
  2046.      * For U only the index 0 for j is used.
  2047.      * </p>
  2048.      *
  2049.      * @author Petre Bazavan
  2050.      * @author Lucian Barbulescu
  2051.      */
  2052.     protected static class UijVijCoefficients {

  2053.         /**
  2054.          * The U₁<sup>j</sup> coefficients.
  2055.          * <p>
  2056.          * The first index identifies the Fourier coefficients used<br>
  2057.          * Those coefficients are computed for all Fourier C<sub>i</sub><sup>j</sup> and
  2058.          * S<sub>i</sub><sup>j</sup><br>
  2059.          * The only exception is when j = 0 when only the coefficient for fourier index
  2060.          * = 1 (i == 0) is needed.<br>
  2061.          * Also, for fourier index = 1 (i == 0), the coefficients up to 2 * jMax are
  2062.          * computed, because are required to compute the coefficients U₂<sup>j</sup>
  2063.          * </p>
  2064.          */
  2065.         private final double[][] u1ij;

  2066.         /**
  2067.          * The V₁<sup>j</sup> coefficients.
  2068.          * <p>
  2069.          * The first index identifies the Fourier coefficients used<br>
  2070.          * Those coefficients are computed for all Fourier C<sub>i</sub><sup>j</sup> and
  2071.          * S<sub>i</sub><sup>j</sup><br>
  2072.          * for fourier index = 1 (i == 0), the coefficients up to 2 * jMax are computed,
  2073.          * because are required to compute the coefficients V₂<sup>j</sup>
  2074.          * </p>
  2075.          */
  2076.         private final double[][] v1ij;

  2077.         /**
  2078.          * The U₂<sup>j</sup> coefficients.
  2079.          * <p>
  2080.          * Only the coefficients that use the Fourier index = 1 (i == 0) are computed as
  2081.          * they are the only ones required.
  2082.          * </p>
  2083.          */
  2084.         private final double[] u2ij;

  2085.         /**
  2086.          * The V₂<sup>j</sup> coefficients.
  2087.          * <p>
  2088.          * Only the coefficients that use the Fourier index = 1 (i == 0) are computed as
  2089.          * they are the only ones required.
  2090.          * </p>
  2091.          */
  2092.         private final double[] v2ij;

  2093.         /**
  2094.          * The current computed values for the ρ<sub>j</sub> and σ<sub>j</sub>
  2095.          * coefficients.
  2096.          */
  2097.         private final double[][] currentRhoSigmaj;

  2098.         /**
  2099.          * The C<sub>i</sub><sup>j</sup> and the S<sub>i</sub><sup>j</sup> Fourier
  2100.          * coefficients.
  2101.          */
  2102.         private final FourierCjSjCoefficients fourierCjSj;

  2103.         /** The maximum value for j index. */
  2104.         private final int jMax;

  2105.         /**
  2106.          * Constructor.
  2107.          * @param currentRhoSigmaj the current computed values for the ρ<sub>j</sub> and
  2108.          *                         σ<sub>j</sub> coefficients
  2109.          * @param fourierCjSj      the fourier coefficients C<sub>i</sub><sup>j</sup>
  2110.          *                         and the S<sub>i</sub><sup>j</sup>
  2111.          * @param jMax             maximum value for j index
  2112.          */
  2113.         UijVijCoefficients(final double[][] currentRhoSigmaj, final FourierCjSjCoefficients fourierCjSj,
  2114.                 final int jMax) {
  2115.             this.currentRhoSigmaj = currentRhoSigmaj;
  2116.             this.fourierCjSj = fourierCjSj;
  2117.             this.jMax = jMax;

  2118.             // initialize the internal arrays.
  2119.             this.u1ij = new double[6][2 * jMax + 1];
  2120.             this.v1ij = new double[6][2 * jMax + 1];
  2121.             this.u2ij = new double[jMax + 1];
  2122.             this.v2ij = new double[jMax + 1];

  2123.             // compute the coefficients
  2124.             computeU1V1Coefficients();
  2125.             computeU2V2Coefficients();
  2126.         }

  2127.         /** Build the U₁<sup>j</sup> and V₁<sup>j</sup> coefficients. */
  2128.         private void computeU1V1Coefficients() {
  2129.             // generate the U₁<sup>j</sup> and V₁<sup>j</sup> coefficients
  2130.             // for j >= 1
  2131.             // also the U₁⁰ for Fourier index = 1 (i == 0) coefficient will be computed
  2132.             u1ij[0][0] = 0;
  2133.             for (int j = 1; j <= jMax; j++) {
  2134.                 // compute 1 / j
  2135.                 final double ooj = 1. / j;

  2136.                 for (int i = 0; i < 6; i++) {
  2137.                     // j is aready between 1 and J
  2138.                     u1ij[i][j] = fourierCjSj.getSij(i, j);
  2139.                     v1ij[i][j] = fourierCjSj.getCij(i, j);

  2140.                     // 1 - δ<sub>1j</sub> is 1 for all j > 1
  2141.                     if (j > 1) {
  2142.                         // k starts with 1 because j-J is less than or equal to 0
  2143.                         for (int kIndex = 1; kIndex <= j - 1; kIndex++) {
  2144.                             // C<sub>i</sub><sup>j-k</sup> * σ<sub>k</sub> +
  2145.                             // S<sub>i</sub><sup>j-k</sup> * ρ<sub>k</sub>
  2146.                             u1ij[i][j] += fourierCjSj.getCij(i, j - kIndex) * currentRhoSigmaj[1][kIndex] +
  2147.                                     fourierCjSj.getSij(i, j - kIndex) * currentRhoSigmaj[0][kIndex];

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

  2154.                     // since j must be between 1 and J-1 and is already between 1 and J
  2155.                     // the following sum is skiped only for j = jMax
  2156.                     if (j != jMax) {
  2157.                         for (int kIndex = 1; kIndex <= jMax - j; kIndex++) {
  2158.                             // -C<sub>i</sub><sup>j+k</sup> * σ<sub>k</sub> +
  2159.                             // S<sub>i</sub><sup>j+k</sup> * ρ<sub>k</sub>
  2160.                             u1ij[i][j] += -fourierCjSj.getCij(i, j + kIndex) * currentRhoSigmaj[1][kIndex] +
  2161.                                     fourierCjSj.getSij(i, j + kIndex) * currentRhoSigmaj[0][kIndex];

  2162.                             // C<sub>i</sub><sup>j+k</sup> * ρ<sub>k</sub> +
  2163.                             // S<sub>i</sub><sup>j+k</sup> * σ<sub>k</sub>
  2164.                             v1ij[i][j] += fourierCjSj.getCij(i, j + kIndex) * currentRhoSigmaj[0][kIndex] +
  2165.                                     fourierCjSj.getSij(i, j + kIndex) * currentRhoSigmaj[1][kIndex];
  2166.                         }
  2167.                     }

  2168.                     for (int kIndex = 1; kIndex <= jMax; kIndex++) {
  2169.                         // C<sub>i</sub><sup>k</sup> * σ<sub>j+k</sub> -
  2170.                         // S<sub>i</sub><sup>k</sup> * ρ<sub>j+k</sub>
  2171.                         u1ij[i][j] += -fourierCjSj.getCij(i, kIndex) * currentRhoSigmaj[1][j + kIndex] -
  2172.                                 fourierCjSj.getSij(i, kIndex) * currentRhoSigmaj[0][j + kIndex];

  2173.                         // C<sub>i</sub><sup>k</sup> * ρ<sub>j+k</sub> +
  2174.                         // S<sub>i</sub><sup>k</sup> * σ<sub>j+k</sub>
  2175.                         v1ij[i][j] += fourierCjSj.getCij(i, kIndex) * currentRhoSigmaj[0][j + kIndex] +
  2176.                                 fourierCjSj.getSij(i, kIndex) * currentRhoSigmaj[1][j + kIndex];
  2177.                     }

  2178.                     // divide by 1 / j
  2179.                     u1ij[i][j] *= -ooj;
  2180.                     v1ij[i][j] *= ooj;

  2181.                     // if index = 1 (i == 0) add the computed terms to U₁⁰
  2182.                     if (i == 0) {
  2183.                         // - (U₁<sup>j</sup> * ρ<sub>j</sub> + V₁<sup>j</sup> * σ<sub>j</sub>
  2184.                         u1ij[0][0] += -u1ij[0][j] * currentRhoSigmaj[0][j] - v1ij[0][j] * currentRhoSigmaj[1][j];
  2185.                     }
  2186.                 }
  2187.             }

  2188.             // Terms with j > jMax are required only when computing the coefficients
  2189.             // U₂<sup>j</sup> and V₂<sup>j</sup>
  2190.             // and those coefficients are only required for Fourier index = 1 (i == 0).
  2191.             for (int j = jMax + 1; j <= 2 * jMax; j++) {
  2192.                 // compute 1 / j
  2193.                 final double ooj = 1. / j;
  2194.                 // the value of i is 0
  2195.                 u1ij[0][j] = 0.;
  2196.                 v1ij[0][j] = 0.;

  2197.                 // k starts from j-J as it is always greater than or equal to 1
  2198.                 for (int kIndex = j - jMax; kIndex <= j - 1; kIndex++) {
  2199.                     // C<sub>i</sub><sup>j-k</sup> * σ<sub>k</sub> +
  2200.                     // S<sub>i</sub><sup>j-k</sup> * ρ<sub>k</sub>
  2201.                     u1ij[0][j] += fourierCjSj.getCij(0, j - kIndex) * currentRhoSigmaj[1][kIndex] +
  2202.                             fourierCjSj.getSij(0, j - kIndex) * currentRhoSigmaj[0][kIndex];

  2203.                     // C<sub>i</sub><sup>j-k</sup> * ρ<sub>k</sub> -
  2204.                     // S<sub>i</sub><sup>j-k</sup> * σ<sub>k</sub>
  2205.                     v1ij[0][j] += fourierCjSj.getCij(0, j - kIndex) * currentRhoSigmaj[0][kIndex] -
  2206.                             fourierCjSj.getSij(0, j - kIndex) * currentRhoSigmaj[1][kIndex];
  2207.                 }
  2208.                 for (int kIndex = 1; kIndex <= jMax; kIndex++) {
  2209.                     // C<sub>i</sub><sup>k</sup> * σ<sub>j+k</sub> -
  2210.                     // S<sub>i</sub><sup>k</sup> * ρ<sub>j+k</sub>
  2211.                     u1ij[0][j] += -fourierCjSj.getCij(0, kIndex) * currentRhoSigmaj[1][j + kIndex] -
  2212.                             fourierCjSj.getSij(0, kIndex) * currentRhoSigmaj[0][j + kIndex];

  2213.                     // C<sub>i</sub><sup>k</sup> * ρ<sub>j+k</sub> +
  2214.                     // S<sub>i</sub><sup>k</sup> * σ<sub>j+k</sub>
  2215.                     v1ij[0][j] += fourierCjSj.getCij(0, kIndex) * currentRhoSigmaj[0][j + kIndex] +
  2216.                             fourierCjSj.getSij(0, kIndex) * currentRhoSigmaj[1][j + kIndex];
  2217.                 }

  2218.                 // divide by 1 / j
  2219.                 u1ij[0][j] *= -ooj;
  2220.                 v1ij[0][j] *= ooj;
  2221.             }
  2222.         }

  2223.         /**
  2224.          * Build the U₁<sup>j</sup> and V₁<sup>j</sup> coefficients.
  2225.          * <p>
  2226.          * Only the coefficients for Fourier index = 1 (i == 0) are required.
  2227.          * </p>
  2228.          */
  2229.         private void computeU2V2Coefficients() {
  2230.             for (int j = 1; j <= jMax; j++) {
  2231.                 // compute 1 / j
  2232.                 final double ooj = 1. / j;

  2233.                 // only the values for i == 0 are computed
  2234.                 u2ij[j] = v1ij[0][j];
  2235.                 v2ij[j] = u1ij[0][j];

  2236.                 // 1 - δ<sub>1j</sub> is 1 for all j > 1
  2237.                 if (j > 1) {
  2238.                     for (int l = 1; l <= j - 1; l++) {
  2239.                         // U₁<sup>j-l</sup> * σ<sub>l</sub> +
  2240.                         // V₁<sup>j-l</sup> * ρ<sub>l</sub>
  2241.                         u2ij[j] += u1ij[0][j - l] * currentRhoSigmaj[1][l] + v1ij[0][j - l] * currentRhoSigmaj[0][l];

  2242.                         // U₁<sup>j-l</sup> * ρ<sub>l</sub> -
  2243.                         // V₁<sup>j-l</sup> * σ<sub>l</sub>
  2244.                         v2ij[j] += u1ij[0][j - l] * currentRhoSigmaj[0][l] - v1ij[0][j - l] * currentRhoSigmaj[1][l];
  2245.                     }
  2246.                 }

  2247.                 for (int l = 1; l <= jMax; l++) {
  2248.                     // -U₁<sup>j+l</sup> * σ<sub>l</sub> +
  2249.                     // U₁<sup>l</sup> * σ<sub>j+l</sub> +
  2250.                     // V₁<sup>j+l</sup> * ρ<sub>l</sub> -
  2251.                     // V₁<sup>l</sup> * ρ<sub>j+l</sub>
  2252.                     u2ij[j] += -u1ij[0][j + l] * currentRhoSigmaj[1][l] + u1ij[0][l] * currentRhoSigmaj[1][j + l] +
  2253.                             v1ij[0][j + l] * currentRhoSigmaj[0][l] - v1ij[0][l] * currentRhoSigmaj[0][j + l];

  2254.                     // U₁<sup>j+l</sup> * ρ<sub>l</sub> +
  2255.                     // U₁<sup>l</sup> * ρ<sub>j+l</sub> +
  2256.                     // V₁<sup>j+l</sup> * σ<sub>l</sub> +
  2257.                     // V₁<sup>l</sup> * σ<sub>j+l</sub>
  2258.                     u2ij[j] += u1ij[0][j + l] * currentRhoSigmaj[0][l] + u1ij[0][l] * currentRhoSigmaj[0][j + l] +
  2259.                             v1ij[0][j + l] * currentRhoSigmaj[1][l] + v1ij[0][l] * currentRhoSigmaj[1][j + l];
  2260.                 }

  2261.                 // divide by 1 / j
  2262.                 u2ij[j] *= -ooj;
  2263.                 v2ij[j] *= ooj;
  2264.             }
  2265.         }

  2266.         /**
  2267.          * Get the coefficient U₁<sup>j</sup> for Fourier index i.
  2268.          *
  2269.          * @param j j index
  2270.          * @param i Fourier index (starts at 0)
  2271.          * @return the coefficient U₁<sup>j</sup> for the given Fourier index i
  2272.          */
  2273.         public double getU1(final int j, final int i) {
  2274.             return u1ij[i][j];
  2275.         }

  2276.         /**
  2277.          * Get the coefficient V₁<sup>j</sup> for Fourier index i.
  2278.          *
  2279.          * @param j j index
  2280.          * @param i Fourier index (starts at 0)
  2281.          * @return the coefficient V₁<sup>j</sup> for the given Fourier index i
  2282.          */
  2283.         public double getV1(final int j, final int i) {
  2284.             return v1ij[i][j];
  2285.         }

  2286.         /**
  2287.          * Get the coefficient U₂<sup>j</sup> for Fourier index = 1 (i == 0).
  2288.          *
  2289.          * @param j j index
  2290.          * @return the coefficient U₂<sup>j</sup> for Fourier index = 1 (i == 0)
  2291.          */
  2292.         public double getU2(final int j) {
  2293.             return u2ij[j];
  2294.         }

  2295.         /**
  2296.          * Get the coefficient V₂<sup>j</sup> for Fourier index = 1 (i == 0).
  2297.          *
  2298.          * @param j j index
  2299.          * @return the coefficient V₂<sup>j</sup> for Fourier index = 1 (i == 0)
  2300.          */
  2301.         public double getV2(final int j) {
  2302.             return v2ij[j];
  2303.         }
  2304.     }

  2305.     /**
  2306.      * The U<sub>i</sub><sup>j</sup> and V<sub>i</sub><sup>j</sup> coefficients
  2307.      * described by equations 2.5.3-(21) and 2.5.3-(22) from Danielson.
  2308.      * <p>
  2309.      * The index i takes only the values 1 and 2<br>
  2310.      * For U only the index 0 for j is used.
  2311.      * </p>
  2312.      *
  2313.      * @author Petre Bazavan
  2314.      * @author Lucian Barbulescu
  2315.      * @param <T> type of the field elements
  2316.      */
  2317.     protected static class FieldUijVijCoefficients<T extends CalculusFieldElement<T>> {

  2318.         /**
  2319.          * The U₁<sup>j</sup> coefficients.
  2320.          * <p>
  2321.          * The first index identifies the Fourier coefficients used<br>
  2322.          * Those coefficients are computed for all Fourier C<sub>i</sub><sup>j</sup> and
  2323.          * S<sub>i</sub><sup>j</sup><br>
  2324.          * The only exception is when j = 0 when only the coefficient for fourier index
  2325.          * = 1 (i == 0) is needed.<br>
  2326.          * Also, for fourier index = 1 (i == 0), the coefficients up to 2 * jMax are
  2327.          * computed, because are required to compute the coefficients U₂<sup>j</sup>
  2328.          * </p>
  2329.          */
  2330.         private final T[][] u1ij;

  2331.         /**
  2332.          * The V₁<sup>j</sup> coefficients.
  2333.          * <p>
  2334.          * The first index identifies the Fourier coefficients used<br>
  2335.          * Those coefficients are computed for all Fourier C<sub>i</sub><sup>j</sup> and
  2336.          * S<sub>i</sub><sup>j</sup><br>
  2337.          * for fourier index = 1 (i == 0), the coefficients up to 2 * jMax are computed,
  2338.          * because are required to compute the coefficients V₂<sup>j</sup>
  2339.          * </p>
  2340.          */
  2341.         private final T[][] v1ij;

  2342.         /**
  2343.          * The U₂<sup>j</sup> coefficients.
  2344.          * <p>
  2345.          * Only the coefficients that use the Fourier index = 1 (i == 0) are computed as
  2346.          * they are the only ones required.
  2347.          * </p>
  2348.          */
  2349.         private final T[] u2ij;

  2350.         /**
  2351.          * The V₂<sup>j</sup> coefficients.
  2352.          * <p>
  2353.          * Only the coefficients that use the Fourier index = 1 (i == 0) are computed as
  2354.          * they are the only ones required.
  2355.          * </p>
  2356.          */
  2357.         private final T[] v2ij;

  2358.         /**
  2359.          * The current computed values for the ρ<sub>j</sub> and σ<sub>j</sub>
  2360.          * coefficients.
  2361.          */
  2362.         private final T[][] currentRhoSigmaj;

  2363.         /**
  2364.          * The C<sub>i</sub><sup>j</sup> and the S<sub>i</sub><sup>j</sup> Fourier
  2365.          * coefficients.
  2366.          */
  2367.         private final FieldFourierCjSjCoefficients<T> fourierCjSj;

  2368.         /** The maximum value for j index. */
  2369.         private final int jMax;

  2370.         /**
  2371.          * Constructor.
  2372.          * @param currentRhoSigmaj the current computed values for the ρ<sub>j</sub> and
  2373.          *                         σ<sub>j</sub> coefficients
  2374.          * @param fourierCjSj      the fourier coefficients C<sub>i</sub><sup>j</sup>
  2375.          *                         and the S<sub>i</sub><sup>j</sup>
  2376.          * @param jMax             maximum value for j index
  2377.          * @param field            field used by default
  2378.          */
  2379.         FieldUijVijCoefficients(final T[][] currentRhoSigmaj, final FieldFourierCjSjCoefficients<T> fourierCjSj,
  2380.                 final int jMax, final Field<T> field) {
  2381.             this.currentRhoSigmaj = currentRhoSigmaj;
  2382.             this.fourierCjSj = fourierCjSj;
  2383.             this.jMax = jMax;

  2384.             // initialize the internal arrays.
  2385.             this.u1ij = MathArrays.buildArray(field, 6, 2 * jMax + 1);
  2386.             this.v1ij = MathArrays.buildArray(field, 6, 2 * jMax + 1);
  2387.             this.u2ij = MathArrays.buildArray(field, jMax + 1);
  2388.             this.v2ij = MathArrays.buildArray(field, jMax + 1);

  2389.             // compute the coefficients
  2390.             computeU1V1Coefficients(field);
  2391.             computeU2V2Coefficients(field);
  2392.         }

  2393.         /**
  2394.          * Build the U₁<sup>j</sup> and V₁<sup>j</sup> coefficients.
  2395.          * @param field field used by default
  2396.          */
  2397.         private void computeU1V1Coefficients(final Field<T> field) {
  2398.             // Zero
  2399.             final T zero = field.getZero();

  2400.             // generate the U₁<sup>j</sup> and V₁<sup>j</sup> coefficients
  2401.             // for j >= 1
  2402.             // also the U₁⁰ for Fourier index = 1 (i == 0) coefficient will be computed
  2403.             u1ij[0][0] = zero;
  2404.             for (int j = 1; j <= jMax; j++) {
  2405.                 // compute 1 / j
  2406.                 final double ooj = 1. / j;

  2407.                 for (int i = 0; i < 6; i++) {
  2408.                     // j is aready between 1 and J
  2409.                     u1ij[i][j] = fourierCjSj.getSij(i, j);
  2410.                     v1ij[i][j] = fourierCjSj.getCij(i, j);

  2411.                     // 1 - δ<sub>1j</sub> is 1 for all j > 1
  2412.                     if (j > 1) {
  2413.                         // k starts with 1 because j-J is less than or equal to 0
  2414.                         for (int kIndex = 1; kIndex <= j - 1; kIndex++) {
  2415.                             // C<sub>i</sub><sup>j-k</sup> * σ<sub>k</sub> +
  2416.                             // S<sub>i</sub><sup>j-k</sup> * ρ<sub>k</sub>
  2417.                             u1ij[i][j] = u1ij[i][j]
  2418.                                     .add(fourierCjSj.getCij(i, j - kIndex).multiply(currentRhoSigmaj[1][kIndex]).add(
  2419.                                             fourierCjSj.getSij(i, j - kIndex).multiply(currentRhoSigmaj[0][kIndex])));

  2420.                             // C<sub>i</sub><sup>j-k</sup> * ρ<sub>k</sub> -
  2421.                             // S<sub>i</sub><sup>j-k</sup> * σ<sub>k</sub>
  2422.                             v1ij[i][j] = v1ij[i][j].add(
  2423.                                     fourierCjSj.getCij(i, j - kIndex).multiply(currentRhoSigmaj[0][kIndex]).subtract(
  2424.                                             fourierCjSj.getSij(i, j - kIndex).multiply(currentRhoSigmaj[1][kIndex])));
  2425.                         }
  2426.                     }

  2427.                     // since j must be between 1 and J-1 and is already between 1 and J
  2428.                     // the following sum is skiped only for j = jMax
  2429.                     if (j != jMax) {
  2430.                         for (int kIndex = 1; kIndex <= jMax - j; kIndex++) {
  2431.                             // -C<sub>i</sub><sup>j+k</sup> * σ<sub>k</sub> +
  2432.                             // S<sub>i</sub><sup>j+k</sup> * ρ<sub>k</sub>
  2433.                             u1ij[i][j] = u1ij[i][j].add(fourierCjSj.getCij(i, j + kIndex).negate()
  2434.                                     .multiply(currentRhoSigmaj[1][kIndex])
  2435.                                     .add(fourierCjSj.getSij(i, j + kIndex).multiply(currentRhoSigmaj[0][kIndex])));

  2436.                             // C<sub>i</sub><sup>j+k</sup> * ρ<sub>k</sub> +
  2437.                             // S<sub>i</sub><sup>j+k</sup> * σ<sub>k</sub>
  2438.                             v1ij[i][j] = v1ij[i][j]
  2439.                                     .add(fourierCjSj.getCij(i, j + kIndex).multiply(currentRhoSigmaj[0][kIndex]).add(
  2440.                                             fourierCjSj.getSij(i, j + kIndex).multiply(currentRhoSigmaj[1][kIndex])));
  2441.                         }
  2442.                     }

  2443.                     for (int kIndex = 1; kIndex <= jMax; kIndex++) {
  2444.                         // C<sub>i</sub><sup>k</sup> * σ<sub>j+k</sub> -
  2445.                         // S<sub>i</sub><sup>k</sup> * ρ<sub>j+k</sub>
  2446.                         u1ij[i][j] = u1ij[i][j].add(fourierCjSj.getCij(i, kIndex).negate()
  2447.                                 .multiply(currentRhoSigmaj[1][j + kIndex])
  2448.                                 .subtract(fourierCjSj.getSij(i, kIndex).multiply(currentRhoSigmaj[0][j + kIndex])));

  2449.                         // C<sub>i</sub><sup>k</sup> * ρ<sub>j+k</sub> +
  2450.                         // S<sub>i</sub><sup>k</sup> * σ<sub>j+k</sub>
  2451.                         v1ij[i][j] = v1ij[i][j]
  2452.                                 .add(fourierCjSj.getCij(i, kIndex).multiply(currentRhoSigmaj[0][j + kIndex])
  2453.                                         .add(fourierCjSj.getSij(i, kIndex).multiply(currentRhoSigmaj[1][j + kIndex])));
  2454.                     }

  2455.                     // divide by 1 / j
  2456.                     u1ij[i][j] = u1ij[i][j].multiply(-ooj);
  2457.                     v1ij[i][j] = v1ij[i][j].multiply(ooj);

  2458.                     // if index = 1 (i == 0) add the computed terms to U₁⁰
  2459.                     if (i == 0) {
  2460.                         // - (U₁<sup>j</sup> * ρ<sub>j</sub> + V₁<sup>j</sup> * σ<sub>j</sub>
  2461.                         u1ij[0][0] = u1ij[0][0].add(u1ij[0][j].negate().multiply(currentRhoSigmaj[0][j])
  2462.                                 .subtract(v1ij[0][j].multiply(currentRhoSigmaj[1][j])));
  2463.                     }
  2464.                 }
  2465.             }

  2466.             // Terms with j > jMax are required only when computing the coefficients
  2467.             // U₂<sup>j</sup> and V₂<sup>j</sup>
  2468.             // and those coefficients are only required for Fourier index = 1 (i == 0).
  2469.             for (int j = jMax + 1; j <= 2 * jMax; j++) {
  2470.                 // compute 1 / j
  2471.                 final double ooj = 1. / j;
  2472.                 // the value of i is 0
  2473.                 u1ij[0][j] = zero;
  2474.                 v1ij[0][j] = zero;

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

  2481.                     // C<sub>i</sub><sup>j-k</sup> * ρ<sub>k</sub> -
  2482.                     // S<sub>i</sub><sup>j-k</sup> * σ<sub>k</sub>
  2483.                     v1ij[0][j] = v1ij[0][j].add(fourierCjSj.getCij(0, j - kIndex).multiply(currentRhoSigmaj[0][kIndex])
  2484.                             .subtract(fourierCjSj.getSij(0, j - kIndex).multiply(currentRhoSigmaj[1][kIndex])));
  2485.                 }
  2486.                 for (int kIndex = 1; kIndex <= jMax; kIndex++) {
  2487.                     // C<sub>i</sub><sup>k</sup> * σ<sub>j+k</sub> -
  2488.                     // S<sub>i</sub><sup>k</sup> * ρ<sub>j+k</sub>
  2489.                     u1ij[0][j] = u1ij[0][j]
  2490.                             .add(fourierCjSj.getCij(0, kIndex).negate().multiply(currentRhoSigmaj[1][j + kIndex])
  2491.                                     .subtract(fourierCjSj.getSij(0, kIndex).multiply(currentRhoSigmaj[0][j + kIndex])));

  2492.                     // C<sub>i</sub><sup>k</sup> * ρ<sub>j+k</sub> +
  2493.                     // S<sub>i</sub><sup>k</sup> * σ<sub>j+k</sub>
  2494.                     v1ij[0][j] = v1ij[0][j].add(fourierCjSj.getCij(0, kIndex).multiply(currentRhoSigmaj[0][j + kIndex])
  2495.                             .add(fourierCjSj.getSij(0, kIndex).multiply(currentRhoSigmaj[1][j + kIndex])));
  2496.                 }

  2497.                 // divide by 1 / j
  2498.                 u1ij[0][j] = u1ij[0][j].multiply(-ooj);
  2499.                 v1ij[0][j] = v1ij[0][j].multiply(ooj);
  2500.             }
  2501.         }

  2502.         /**
  2503.          * Build the U₁<sup>j</sup> and V₁<sup>j</sup> coefficients.
  2504.          * <p>
  2505.          * Only the coefficients for Fourier index = 1 (i == 0) are required.
  2506.          * </p>
  2507.          * @param field field used by default
  2508.          */
  2509.         private void computeU2V2Coefficients(final Field<T> field) {
  2510.             for (int j = 1; j <= jMax; j++) {
  2511.                 // compute 1 / j
  2512.                 final double ooj = 1. / j;

  2513.                 // only the values for i == 0 are computed
  2514.                 u2ij[j] = v1ij[0][j];
  2515.                 v2ij[j] = u1ij[0][j];

  2516.                 // 1 - δ<sub>1j</sub> is 1 for all j > 1
  2517.                 if (j > 1) {
  2518.                     for (int l = 1; l <= j - 1; l++) {
  2519.                         // U₁<sup>j-l</sup> * σ<sub>l</sub> +
  2520.                         // V₁<sup>j-l</sup> * ρ<sub>l</sub>
  2521.                         u2ij[j] = u2ij[j].add(u1ij[0][j - l].multiply(currentRhoSigmaj[1][l])
  2522.                                 .add(v1ij[0][j - l].multiply(currentRhoSigmaj[0][l])));

  2523.                         // U₁<sup>j-l</sup> * ρ<sub>l</sub> -
  2524.                         // V₁<sup>j-l</sup> * σ<sub>l</sub>
  2525.                         v2ij[j] = v2ij[j].add(u1ij[0][j - l].multiply(currentRhoSigmaj[0][l])
  2526.                                 .subtract(v1ij[0][j - l].multiply(currentRhoSigmaj[1][l])));
  2527.                     }
  2528.                 }

  2529.                 for (int l = 1; l <= jMax; l++) {
  2530.                     // -U₁<sup>j+l</sup> * σ<sub>l</sub> +
  2531.                     // U₁<sup>l</sup> * σ<sub>j+l</sub> +
  2532.                     // V₁<sup>j+l</sup> * ρ<sub>l</sub> -
  2533.                     // V₁<sup>l</sup> * ρ<sub>j+l</sub>
  2534.                     u2ij[j] = u2ij[j].add(u1ij[0][j + l].negate().multiply(currentRhoSigmaj[1][l])
  2535.                             .add(u1ij[0][l].multiply(currentRhoSigmaj[1][j + l]))
  2536.                             .add(v1ij[0][j + l].multiply(currentRhoSigmaj[0][l]))
  2537.                             .subtract(v1ij[0][l].multiply(currentRhoSigmaj[0][j + l])));

  2538.                     // U₁<sup>j+l</sup> * ρ<sub>l</sub> +
  2539.                     // U₁<sup>l</sup> * ρ<sub>j+l</sub> +
  2540.                     // V₁<sup>j+l</sup> * σ<sub>l</sub> +
  2541.                     // V₁<sup>l</sup> * σ<sub>j+l</sub>
  2542.                     u2ij[j] = u2ij[j].add(u1ij[0][j + l].multiply(currentRhoSigmaj[0][l])
  2543.                             .add(u1ij[0][l].multiply(currentRhoSigmaj[0][j + l]))
  2544.                             .add(v1ij[0][j + l].multiply(currentRhoSigmaj[1][l]))
  2545.                             .add(v1ij[0][l].multiply(currentRhoSigmaj[1][j + l])));
  2546.                 }

  2547.                 // divide by 1 / j
  2548.                 u2ij[j] = u2ij[j].multiply(-ooj);
  2549.                 v2ij[j] = v2ij[j].multiply(ooj);
  2550.             }
  2551.         }

  2552.         /**
  2553.          * Get the coefficient U₁<sup>j</sup> for Fourier index i.
  2554.          *
  2555.          * @param j j index
  2556.          * @param i Fourier index (starts at 0)
  2557.          * @return the coefficient U₁<sup>j</sup> for the given Fourier index i
  2558.          */
  2559.         public T getU1(final int j, final int i) {
  2560.             return u1ij[i][j];
  2561.         }

  2562.         /**
  2563.          * Get the coefficient V₁<sup>j</sup> for Fourier index i.
  2564.          *
  2565.          * @param j j index
  2566.          * @param i Fourier index (starts at 0)
  2567.          * @return the coefficient V₁<sup>j</sup> for the given Fourier index i
  2568.          */
  2569.         public T getV1(final int j, final int i) {
  2570.             return v1ij[i][j];
  2571.         }

  2572.         /**
  2573.          * Get the coefficient U₂<sup>j</sup> for Fourier index = 1 (i == 0).
  2574.          *
  2575.          * @param j j index
  2576.          * @return the coefficient U₂<sup>j</sup> for Fourier index = 1 (i == 0)
  2577.          */
  2578.         public T getU2(final int j) {
  2579.             return u2ij[j];
  2580.         }

  2581.         /**
  2582.          * Get the coefficient V₂<sup>j</sup> for Fourier index = 1 (i == 0).
  2583.          *
  2584.          * @param j j index
  2585.          * @return the coefficient V₂<sup>j</sup> for Fourier index = 1 (i == 0)
  2586.          */
  2587.         public T getV2(final int j) {
  2588.             return v2ij[j];
  2589.         }
  2590.     }

  2591.     /** Coefficients valid for one time slot. */
  2592.     protected static class Slot {

  2593.         /**
  2594.          * The coefficients D<sub>i</sub><sup>j</sup>.
  2595.          * <p>
  2596.          * Only for j = 1 and j = 2 the coefficients are not 0. <br>
  2597.          * i corresponds to the equinoctial element, as follows: - i=0 for a <br/>
  2598.          * - i=1 for k <br/>
  2599.          * - i=2 for h <br/>
  2600.          * - i=3 for q <br/>
  2601.          * - i=4 for p <br/>
  2602.          * - i=5 for λ <br/>
  2603.          * </p>
  2604.          */
  2605.         private final ShortPeriodicsInterpolatedCoefficient[] dij;

  2606.         /**
  2607.          * The coefficients C<sub>i</sub><sup>j</sup>.
  2608.          * <p>
  2609.          * The index order is cij[j][i] <br/>
  2610.          * i corresponds to the equinoctial element, as follows: <br/>
  2611.          * - i=0 for a <br/>
  2612.          * - i=1 for k <br/>
  2613.          * - i=2 for h <br/>
  2614.          * - i=3 for q <br/>
  2615.          * - i=4 for p <br/>
  2616.          * - i=5 for λ <br/>
  2617.          * </p>
  2618.          */
  2619.         private final ShortPeriodicsInterpolatedCoefficient[] cij;

  2620.         /**
  2621.          * The coefficients S<sub>i</sub><sup>j</sup>.
  2622.          * <p>
  2623.          * The index order is sij[j][i] <br/>
  2624.          * i corresponds to the equinoctial element, as follows: <br/>
  2625.          * - i=0 for a <br/>
  2626.          * - i=1 for k <br/>
  2627.          * - i=2 for h <br/>
  2628.          * - i=3 for q <br/>
  2629.          * - i=4 for p <br/>
  2630.          * - i=5 for λ <br/>
  2631.          * </p>
  2632.          */
  2633.         private final ShortPeriodicsInterpolatedCoefficient[] sij;

  2634.         /**
  2635.          * Simple constructor.
  2636.          * @param jMax                maximum value for j index
  2637.          * @param interpolationPoints number of points used in the interpolation process
  2638.          */
  2639.         Slot(final int jMax, final int interpolationPoints) {

  2640.             dij = new ShortPeriodicsInterpolatedCoefficient[3];
  2641.             cij = new ShortPeriodicsInterpolatedCoefficient[jMax + 1];
  2642.             sij = new ShortPeriodicsInterpolatedCoefficient[jMax + 1];

  2643.             // Initialize the C<sub>i</sub><sup>j</sup>, S<sub>i</sub><sup>j</sup> and
  2644.             // D<sub>i</sub><sup>j</sup> coefficients
  2645.             for (int j = 0; j <= jMax; j++) {
  2646.                 cij[j] = new ShortPeriodicsInterpolatedCoefficient(interpolationPoints);
  2647.                 if (j > 0) {
  2648.                     sij[j] = new ShortPeriodicsInterpolatedCoefficient(interpolationPoints);
  2649.                 }
  2650.                 // Initialize only the non-zero D<sub>i</sub><sup>j</sup> coefficients
  2651.                 if (j == 1 || j == 2) {
  2652.                     dij[j] = new ShortPeriodicsInterpolatedCoefficient(interpolationPoints);
  2653.                 }
  2654.             }

  2655.         }

  2656.     }

  2657.     /** Coefficients valid for one time slot.
  2658.      * @param <T> type of the field elements
  2659.      */
  2660.     protected static class FieldSlot<T extends CalculusFieldElement<T>> {

  2661.         /**
  2662.          * The coefficients D<sub>i</sub><sup>j</sup>.
  2663.          * <p>
  2664.          * Only for j = 1 and j = 2 the coefficients are not 0. <br>
  2665.          * i corresponds to the equinoctial element, as follows: - i=0 for a <br/>
  2666.          * - i=1 for k <br/>
  2667.          * - i=2 for h <br/>
  2668.          * - i=3 for q <br/>
  2669.          * - i=4 for p <br/>
  2670.          * - i=5 for λ <br/>
  2671.          * </p>
  2672.          */
  2673.         private final FieldShortPeriodicsInterpolatedCoefficient<T>[] dij;

  2674.         /**
  2675.          * The coefficients C<sub>i</sub><sup>j</sup>.
  2676.          * <p>
  2677.          * The index order is cij[j][i] <br/>
  2678.          * i corresponds to the equinoctial element, as follows: <br/>
  2679.          * - i=0 for a <br/>
  2680.          * - i=1 for k <br/>
  2681.          * - i=2 for h <br/>
  2682.          * - i=3 for q <br/>
  2683.          * - i=4 for p <br/>
  2684.          * - i=5 for λ <br/>
  2685.          * </p>
  2686.          */
  2687.         private final FieldShortPeriodicsInterpolatedCoefficient<T>[] cij;

  2688.         /**
  2689.          * The coefficients S<sub>i</sub><sup>j</sup>.
  2690.          * <p>
  2691.          * The index order is sij[j][i] <br/>
  2692.          * i corresponds to the equinoctial element, as follows: <br/>
  2693.          * - i=0 for a <br/>
  2694.          * - i=1 for k <br/>
  2695.          * - i=2 for h <br/>
  2696.          * - i=3 for q <br/>
  2697.          * - i=4 for p <br/>
  2698.          * - i=5 for λ <br/>
  2699.          * </p>
  2700.          */
  2701.         private final FieldShortPeriodicsInterpolatedCoefficient<T>[] sij;

  2702.         /**
  2703.          * Simple constructor.
  2704.          * @param jMax                maximum value for j index
  2705.          * @param interpolationPoints number of points used in the interpolation process
  2706.          */
  2707.         @SuppressWarnings("unchecked")
  2708.         FieldSlot(final int jMax, final int interpolationPoints) {

  2709.             dij = (FieldShortPeriodicsInterpolatedCoefficient<T>[]) Array
  2710.                     .newInstance(FieldShortPeriodicsInterpolatedCoefficient.class, 3);
  2711.             cij = (FieldShortPeriodicsInterpolatedCoefficient<T>[]) Array
  2712.                     .newInstance(FieldShortPeriodicsInterpolatedCoefficient.class, jMax + 1);
  2713.             sij = (FieldShortPeriodicsInterpolatedCoefficient<T>[]) Array
  2714.                     .newInstance(FieldShortPeriodicsInterpolatedCoefficient.class, jMax + 1);

  2715.             // Initialize the C<sub>i</sub><sup>j</sup>, S<sub>i</sub><sup>j</sup> and
  2716.             // D<sub>i</sub><sup>j</sup> coefficients
  2717.             for (int j = 0; j <= jMax; j++) {
  2718.                 cij[j] = new FieldShortPeriodicsInterpolatedCoefficient<>(interpolationPoints);
  2719.                 if (j > 0) {
  2720.                     sij[j] = new FieldShortPeriodicsInterpolatedCoefficient<>(interpolationPoints);
  2721.                 }
  2722.                 // Initialize only the non-zero D<sub>i</sub><sup>j</sup> coefficients
  2723.                 if (j == 1 || j == 2) {
  2724.                     dij[j] = new FieldShortPeriodicsInterpolatedCoefficient<>(interpolationPoints);
  2725.                 }
  2726.             }

  2727.         }

  2728.     }

  2729. }