FieldBrouwerLyddanePropagator.java

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

  18. import java.util.Collections;
  19. import java.util.List;

  20. import org.hipparchus.CalculusFieldElement;
  21. import org.hipparchus.Field;
  22. import org.hipparchus.analysis.differentiation.FieldUnivariateDerivative2;
  23. import org.hipparchus.util.CombinatoricsUtils;
  24. import org.hipparchus.util.FastMath;
  25. import org.hipparchus.util.FieldSinCos;
  26. import org.hipparchus.util.MathUtils;
  27. import org.hipparchus.util.Precision;
  28. import org.orekit.attitudes.AttitudeProvider;
  29. import org.orekit.attitudes.InertialProvider;
  30. import org.orekit.errors.OrekitException;
  31. import org.orekit.errors.OrekitMessages;
  32. import org.orekit.forces.gravity.potential.UnnormalizedSphericalHarmonicsProvider;
  33. import org.orekit.forces.gravity.potential.UnnormalizedSphericalHarmonicsProvider.UnnormalizedSphericalHarmonics;
  34. import org.orekit.orbits.FieldKeplerianOrbit;
  35. import org.orekit.orbits.FieldOrbit;
  36. import org.orekit.orbits.OrbitType;
  37. import org.orekit.orbits.PositionAngle;
  38. import org.orekit.propagation.FieldSpacecraftState;
  39. import org.orekit.propagation.PropagationType;
  40. import org.orekit.propagation.analytical.tle.FieldTLE;
  41. import org.orekit.time.FieldAbsoluteDate;
  42. import org.orekit.utils.FieldTimeSpanMap;
  43. import org.orekit.utils.ParameterDriver;

  44. /** This class propagates a {@link org.orekit.propagation.FieldSpacecraftState}
  45.  *  using the analytical Brouwer-Lyddane model (from J2 to J5 zonal harmonics).
  46.  * <p>
  47.  * At the opposite of the {@link FieldEcksteinHechlerPropagator}, the Brouwer-Lyddane model is
  48.  * suited for elliptical orbits, there is no problem having a rather small eccentricity or inclination
  49.  * (Lyddane helped to solve this issue with the Brouwer model). Singularity for the critical
  50.  * inclination i = 63.4° is avoided using the method developed in Warren Phipps' 1992 thesis.
  51.  * <p>
  52.  * By default, Brouwer-Lyddane model considers only the perturbations due to zonal harmonics.
  53.  * However, for low Earth orbits, the magnitude of the perturbative acceleration due to
  54.  * atmospheric drag can be significant. Warren Phipps' 1992 thesis considered the atmospheric
  55.  * drag by time derivatives of the <i>mean</i> mean anomaly using the catch-all coefficient
  56.  * {@link #M2Driver}.
  57.  *
  58.  * Usually, M2 is adjusted during an orbit determination process and it represents the
  59.  * combination of all unmodeled secular along-track effects (i.e. not just the atmospheric drag).
  60.  * The behavior of M2 is close to the {@link FieldTLE#getBStar()} parameter for the TLE.
  61.  *
  62.  * If the value of M2 is equal to {@link BrouwerLyddanePropagator#M2 0.0}, the along-track secular
  63.  * effects are not considered in the dynamical model. Typical values for M2 are not known.
  64.  * It depends on the orbit type. However, the value of M2 must be very small (e.g. between 1.0e-14 and 1.0e-15).
  65.  * The unit of M2 is rad/s².
  66.  *
  67.  * The along-track effects, represented by the secular rates of the mean semi-major axis
  68.  * and eccentricity, are computed following Eq. 2.38, 2.41, and 2.45 of Warren Phipps' thesis.
  69.  *
  70.  * @see "Brouwer, Dirk. Solution of the problem of artificial satellite theory without drag.
  71.  *       YALE UNIV NEW HAVEN CT NEW HAVEN United States, 1959."
  72.  *
  73.  * @see "Lyddane, R. H. Small eccentricities or inclinations in the Brouwer theory of the
  74.  *       artificial satellite. The Astronomical Journal 68 (1963): 555."
  75.  *
  76.  * @see "Phipps Jr, Warren E. Parallelization of the Navy Space Surveillance Center
  77.  *       (NAVSPASUR) Satellite Model. NAVAL POSTGRADUATE SCHOOL MONTEREY CA, 1992."
  78.  *
  79.  * @author Melina Vanel
  80.  * @author Bryan Cazabonne
  81.  * @since 11.1
  82.  */
  83. public class FieldBrouwerLyddanePropagator<T extends CalculusFieldElement<T>> extends FieldAbstractAnalyticalPropagator<T>  {

  84.     /** Default convergence threshold for mean parameters conversion. */
  85.     private static final double EPSILON_DEFAULT = 1.0e-13;

  86.     /** Default value for maxIterations. */
  87.     private static final int MAX_ITERATIONS_DEFAULT = 200;

  88.     /** Parameters scaling factor.
  89.      * <p>
  90.      * We use a power of 2 to avoid numeric noise introduction
  91.      * in the multiplications/divisions sequences.
  92.      * </p>
  93.      */
  94.     private static final double SCALE = FastMath.scalb(1.0, -20);

  95.     /** Beta constant used by T2 function. */
  96.     private static final double BETA = FastMath.scalb(100, -11);

  97.     /** Initial Brouwer-Lyddane model. */
  98.     private FieldBLModel<T> initialModel;

  99.     /** All models. */
  100.     private transient FieldTimeSpanMap<FieldBLModel<T>, T> models;

  101.     /** Reference radius of the central body attraction model (m). */
  102.     private double referenceRadius;

  103.     /** Central attraction coefficient (m³/s²). */
  104.     private T mu;

  105.     /** Un-normalized zonal coefficients. */
  106.     private double[] ck0;

  107.     /** Empirical coefficient used in the drag modeling. */
  108.     private final ParameterDriver M2Driver;

  109.     /** Build a propagator from orbit and potential provider.
  110.      * <p>Mass and attitude provider are set to unspecified non-null arbitrary values.</p>
  111.      *
  112.      * <p>Using this constructor, an initial osculating orbit is considered.</p>
  113.      *
  114.      * @param initialOrbit initial orbit
  115.      * @param provider for un-normalized zonal coefficients
  116.      * @param M2 value of empirical drag coefficient in rad/s².
  117.      *        If equal to {@link BrouwerLyddanePropagator#M2} drag is not computed
  118.      * @see #FieldBrouwerLyddanePropagator(FieldOrbit, UnnormalizedSphericalHarmonicsProvider, PropagationType, double)
  119.      */
  120.     public FieldBrouwerLyddanePropagator(final FieldOrbit<T> initialOrbit,
  121.                                          final UnnormalizedSphericalHarmonicsProvider provider,
  122.                                          final double M2) {
  123.         this(initialOrbit, InertialProvider.of(initialOrbit.getFrame()),
  124.              initialOrbit.getMu().newInstance(DEFAULT_MASS), provider,
  125.              provider.onDate(initialOrbit.getDate().toAbsoluteDate()), M2);
  126.     }

  127.     /**
  128.      * Private helper constructor.
  129.      * <p>Using this constructor, an initial osculating orbit is considered.</p>
  130.      * @param initialOrbit initial orbit
  131.      * @param attitude attitude provider
  132.      * @param mass spacecraft mass
  133.      * @param provider for un-normalized zonal coefficients
  134.      * @param harmonics {@code provider.onDate(initialOrbit.getDate())}
  135.      * @param M2 value of empirical drag coefficient in rad/s².
  136.      *        If equal to {@link BrouwerLyddanePropagator#M2} drag is not computed
  137.      * @see #FieldBrouwerLyddanePropagator(FieldOrbit, AttitudeProvider, CalculusFieldElement,
  138.      * UnnormalizedSphericalHarmonicsProvider, UnnormalizedSphericalHarmonicsProvider.UnnormalizedSphericalHarmonics, PropagationType, double)
  139.      */
  140.     public FieldBrouwerLyddanePropagator(final FieldOrbit<T> initialOrbit,
  141.                                          final AttitudeProvider attitude,
  142.                                          final T mass,
  143.                                          final UnnormalizedSphericalHarmonicsProvider provider,
  144.                                          final UnnormalizedSphericalHarmonics harmonics,
  145.                                          final double M2) {
  146.         this(initialOrbit, attitude,  mass, provider.getAe(), initialOrbit.getMu().newInstance(provider.getMu()),
  147.              harmonics.getUnnormalizedCnm(2, 0),
  148.              harmonics.getUnnormalizedCnm(3, 0),
  149.              harmonics.getUnnormalizedCnm(4, 0),
  150.              harmonics.getUnnormalizedCnm(5, 0),
  151.              M2);
  152.     }

  153.     /** Build a propagator from orbit and potential.
  154.      * <p>Mass and attitude provider are set to unspecified non-null arbitrary values.</p>
  155.      * <p>The C<sub>n,0</sub> coefficients are the denormalized zonal coefficients, they
  156.      * are related to both the normalized coefficients
  157.      * <span style="text-decoration: overline">C</span><sub>n,0</sub>
  158.      *  and the J<sub>n</sub> one as follows:</p>
  159.      *
  160.      * <p> C<sub>n,0</sub> = [(2-δ<sub>0,m</sub>)(2n+1)(n-m)!/(n+m)!]<sup>½</sup>
  161.      * <span style="text-decoration: overline">C</span><sub>n,0</sub>
  162.      *
  163.      * <p> C<sub>n,0</sub> = -J<sub>n</sub>
  164.      *
  165.      * <p>Using this constructor, an initial osculating orbit is considered.</p>
  166.      *
  167.      * @param initialOrbit initial orbit
  168.      * @param referenceRadius reference radius of the Earth for the potential model (m)
  169.      * @param mu central attraction coefficient (m³/s²)
  170.      * @param c20 un-normalized zonal coefficient (about -1.08e-3 for Earth)
  171.      * @param c30 un-normalized zonal coefficient (about +2.53e-6 for Earth)
  172.      * @param c40 un-normalized zonal coefficient (about +1.62e-6 for Earth)
  173.      * @param c50 un-normalized zonal coefficient (about +2.28e-7 for Earth)
  174.      * @param M2 value of empirical drag coefficient in rad/s².
  175.      *        If equal to {@link BrouwerLyddanePropagator#M2} drag is not computed
  176.      * @see org.orekit.utils.Constants
  177.      * @see #FieldBrouwerLyddanePropagator(FieldOrbit, AttitudeProvider, double, CalculusFieldElement, double, double, double, double, double)
  178.      */
  179.     public FieldBrouwerLyddanePropagator(final FieldOrbit<T> initialOrbit,
  180.                                          final double referenceRadius, final T mu,
  181.                                          final double c20, final double c30, final double c40,
  182.                                          final double c50, final double M2) {
  183.         this(initialOrbit, InertialProvider.of(initialOrbit.getFrame()),
  184.              initialOrbit.getMu().newInstance(DEFAULT_MASS),
  185.              referenceRadius, mu, c20, c30, c40, c50, M2);
  186.     }

  187.     /** Build a propagator from orbit, mass and potential provider.
  188.      * <p>Attitude law is set to an unspecified non-null arbitrary value.</p>
  189.      *
  190.      * <p>Using this constructor, an initial osculating orbit is considered.</p>
  191.      *
  192.      * @param initialOrbit initial orbit
  193.      * @param mass spacecraft mass
  194.      * @param provider for un-normalized zonal coefficients
  195.      * @param M2 value of empirical drag coefficient in rad/s².
  196.      *        If equal to {@link BrouwerLyddanePropagator#M2} drag is not computed
  197.      * @see #FieldBrouwerLyddanePropagator(FieldOrbit, AttitudeProvider, CalculusFieldElement, UnnormalizedSphericalHarmonicsProvider, double)
  198.      */
  199.     public FieldBrouwerLyddanePropagator(final FieldOrbit<T> initialOrbit, final T mass,
  200.                                          final UnnormalizedSphericalHarmonicsProvider provider,
  201.                                          final double M2) {
  202.         this(initialOrbit, InertialProvider.of(initialOrbit.getFrame()),
  203.              mass, provider, provider.onDate(initialOrbit.getDate().toAbsoluteDate()), M2);
  204.     }

  205.     /** Build a propagator from orbit, mass and potential.
  206.      * <p>Attitude law is set to an unspecified non-null arbitrary value.</p>
  207.      * <p>The C<sub>n,0</sub> coefficients are the denormalized zonal coefficients, they
  208.      * are related to both the normalized coefficients
  209.      * <span style="text-decoration: overline">C</span><sub>n,0</sub>
  210.      *  and the J<sub>n</sub> one as follows:</p>
  211.      *
  212.      * <p> C<sub>n,0</sub> = [(2-δ<sub>0,m</sub>)(2n+1)(n-m)!/(n+m)!]<sup>½</sup>
  213.      * <span style="text-decoration: overline">C</span><sub>n,0</sub>
  214.      *
  215.      * <p> C<sub>n,0</sub> = -J<sub>n</sub>
  216.      *
  217.      * <p>Using this constructor, an initial osculating orbit is considered.</p>
  218.      *
  219.      * @param initialOrbit initial orbit
  220.      * @param mass spacecraft mass
  221.      * @param referenceRadius reference radius of the Earth for the potential model (m)
  222.      * @param mu central attraction coefficient (m³/s²)
  223.      * @param c20 un-normalized zonal coefficient (about -1.08e-3 for Earth)
  224.      * @param c30 un-normalized zonal coefficient (about +2.53e-6 for Earth)
  225.      * @param c40 un-normalized zonal coefficient (about +1.62e-6 for Earth)
  226.      * @param c50 un-normalized zonal coefficient (about +2.28e-7 for Earth)
  227.      * @param M2 value of empirical drag coefficient in rad/s².
  228.      *        If equal to {@link BrouwerLyddanePropagator#M2} drag is not computed
  229.      * @see #FieldBrouwerLyddanePropagator(FieldOrbit, AttitudeProvider, CalculusFieldElement, double, CalculusFieldElement, double, double, double, double, double)
  230.      */
  231.     public FieldBrouwerLyddanePropagator(final FieldOrbit<T> initialOrbit, final T mass,
  232.                                          final double referenceRadius, final T mu,
  233.                                          final double c20, final double c30, final double c40,
  234.                                          final double c50, final double M2) {
  235.         this(initialOrbit, InertialProvider.of(initialOrbit.getFrame()),
  236.              mass, referenceRadius, mu, c20, c30, c40, c50, M2);
  237.     }

  238.     /** Build a propagator from orbit, attitude provider and potential provider.
  239.      * <p>Mass is set to an unspecified non-null arbitrary value.</p>
  240.      * <p>Using this constructor, an initial osculating orbit is considered.</p>
  241.      * @param initialOrbit initial orbit
  242.      * @param attitudeProv attitude provider
  243.      * @param provider for un-normalized zonal coefficients
  244.      * @param M2 value of empirical drag coefficient in rad/s².
  245.      *        If equal to {@link BrouwerLyddanePropagator#M2} drag is not computed
  246.      */
  247.     public FieldBrouwerLyddanePropagator(final FieldOrbit<T> initialOrbit,
  248.                                          final AttitudeProvider attitudeProv,
  249.                                          final UnnormalizedSphericalHarmonicsProvider provider,
  250.                                          final double M2) {
  251.         this(initialOrbit, attitudeProv, initialOrbit.getMu().newInstance(DEFAULT_MASS), provider,
  252.              provider.onDate(initialOrbit.getDate().toAbsoluteDate()), M2);
  253.     }

  254.     /** Build a propagator from orbit, attitude provider and potential.
  255.      * <p>Mass is set to an unspecified non-null arbitrary value.</p>
  256.      * <p>The C<sub>n,0</sub> coefficients are the denormalized zonal coefficients, they
  257.      * are related to both the normalized coefficients
  258.      * <span style="text-decoration: overline">C</span><sub>n,0</sub>
  259.      *  and the J<sub>n</sub> one as follows:</p>
  260.      *
  261.      * <p> C<sub>n,0</sub> = [(2-δ<sub>0,m</sub>)(2n+1)(n-m)!/(n+m)!]<sup>½</sup>
  262.      * <span style="text-decoration: overline">C</span><sub>n,0</sub>
  263.      *
  264.      * <p> C<sub>n,0</sub> = -J<sub>n</sub>
  265.      *
  266.      * <p>Using this constructor, an initial osculating orbit is considered.</p>
  267.      *
  268.      * @param initialOrbit initial orbit
  269.      * @param attitudeProv attitude provider
  270.      * @param referenceRadius reference radius of the Earth for the potential model (m)
  271.      * @param mu central attraction coefficient (m³/s²)
  272.      * @param c20 un-normalized zonal coefficient (about -1.08e-3 for Earth)
  273.      * @param c30 un-normalized zonal coefficient (about +2.53e-6 for Earth)
  274.      * @param c40 un-normalized zonal coefficient (about +1.62e-6 for Earth)
  275.      * @param c50 un-normalized zonal coefficient (about +2.28e-7 for Earth)
  276.      * @param M2 value of empirical drag coefficient in rad/s².
  277.      *        If equal to {@link BrouwerLyddanePropagator#M2} drag is not computed
  278.      */
  279.     public FieldBrouwerLyddanePropagator(final FieldOrbit<T> initialOrbit,
  280.                                          final AttitudeProvider attitudeProv,
  281.                                          final double referenceRadius, final T mu,
  282.                                          final double c20, final double c30, final double c40,
  283.                                          final double c50, final double M2) {
  284.         this(initialOrbit, attitudeProv, initialOrbit.getMu().newInstance(DEFAULT_MASS),
  285.              referenceRadius, mu, c20, c30, c40, c50, M2);
  286.     }

  287.     /** Build a propagator from orbit, attitude provider, mass and potential provider.
  288.      * <p>Using this constructor, an initial osculating orbit is considered.</p>
  289.      * @param initialOrbit initial orbit
  290.      * @param attitudeProv attitude provider
  291.      * @param mass spacecraft mass
  292.      * @param provider for un-normalized zonal coefficients
  293.      * @param M2 value of empirical drag coefficient in rad/s².
  294.      *        If equal to {@link BrouwerLyddanePropagator#M2} drag is not computed
  295.      * @see #FieldBrouwerLyddanePropagator(FieldOrbit, AttitudeProvider, CalculusFieldElement, UnnormalizedSphericalHarmonicsProvider, PropagationType, double)
  296.      */
  297.     public FieldBrouwerLyddanePropagator(final FieldOrbit<T> initialOrbit,
  298.                                          final AttitudeProvider attitudeProv,
  299.                                          final T mass,
  300.                                          final UnnormalizedSphericalHarmonicsProvider provider,
  301.                                          final double M2) {
  302.         this(initialOrbit, attitudeProv, mass, provider, provider.onDate(initialOrbit.getDate().toAbsoluteDate()), M2);
  303.     }

  304.     /** Build a propagator from orbit, attitude provider, mass and potential.
  305.      * <p>The C<sub>n,0</sub> coefficients are the denormalized zonal coefficients, they
  306.      * are related to both the normalized coefficients
  307.      * <span style="text-decoration: overline">C</span><sub>n,0</sub>
  308.      *  and the J<sub>n</sub> one as follows:</p>
  309.      *
  310.      * <p> C<sub>n,0</sub> = [(2-δ<sub>0,m</sub>)(2n+1)(n-m)!/(n+m)!]<sup>½</sup>
  311.      * <span style="text-decoration: overline">C</span><sub>n,0</sub>
  312.      *
  313.      * <p> C<sub>n,0</sub> = -J<sub>n</sub>
  314.      *
  315.      * <p>Using this constructor, an initial osculating orbit is considered.</p>
  316.      *
  317.      * @param initialOrbit initial orbit
  318.      * @param attitudeProv attitude provider
  319.      * @param mass spacecraft mass
  320.      * @param referenceRadius reference radius of the Earth for the potential model (m)
  321.      * @param mu central attraction coefficient (m³/s²)
  322.      * @param c20 un-normalized zonal coefficient (about -1.08e-3 for Earth)
  323.      * @param c30 un-normalized zonal coefficient (about +2.53e-6 for Earth)
  324.      * @param c40 un-normalized zonal coefficient (about +1.62e-6 for Earth)
  325.      * @param c50 un-normalized zonal coefficient (about +2.28e-7 for Earth)
  326.      * @param M2 value of empirical drag coefficient in rad/s².
  327.      *        If equal to {@link BrouwerLyddanePropagator#M2} drag is not computed
  328.      * @see #FieldBrouwerLyddanePropagator(FieldOrbit, AttitudeProvider, CalculusFieldElement, double, CalculusFieldElement, double, double, double, double, PropagationType, double)
  329.      */
  330.     public FieldBrouwerLyddanePropagator(final FieldOrbit<T> initialOrbit,
  331.                                          final AttitudeProvider attitudeProv,
  332.                                          final T mass,
  333.                                          final double referenceRadius, final T mu,
  334.                                          final double c20, final double c30, final double c40,
  335.                                          final double c50, final double M2) {
  336.         this(initialOrbit, attitudeProv, mass, referenceRadius, mu, c20, c30, c40, c50, PropagationType.OSCULATING, M2);
  337.     }


  338.     /** Build a propagator from orbit and potential provider.
  339.      * <p>Mass and attitude provider are set to unspecified non-null arbitrary values.</p>
  340.      *
  341.      * <p>Using this constructor, it is possible to define the initial orbit as
  342.      * a mean Brouwer-Lyddane orbit or an osculating one.</p>
  343.      *
  344.      * @param initialOrbit initial orbit
  345.      * @param provider for un-normalized zonal coefficients
  346.      * @param initialType initial orbit type (mean Brouwer-Lyddane orbit or osculating orbit)
  347.      * @param M2 value of empirical drag coefficient in rad/s².
  348.      *        If equal to {@link BrouwerLyddanePropagator#M2} drag is not computed
  349.      */
  350.     public FieldBrouwerLyddanePropagator(final FieldOrbit<T> initialOrbit,
  351.                                          final UnnormalizedSphericalHarmonicsProvider provider,
  352.                                          final PropagationType initialType,
  353.                                          final double M2) {
  354.         this(initialOrbit, InertialProvider.of(initialOrbit.getFrame()),
  355.              initialOrbit.getMu().newInstance(DEFAULT_MASS), provider,
  356.              provider.onDate(initialOrbit.getDate().toAbsoluteDate()), initialType, M2);
  357.     }

  358.     /** Build a propagator from orbit, attitude provider, mass and potential provider.
  359.      * <p>Using this constructor, it is possible to define the initial orbit as
  360.      * a mean Brouwer-Lyddane orbit or an osculating one.</p>
  361.      * @param initialOrbit initial orbit
  362.      * @param attitudeProv attitude provider
  363.      * @param mass spacecraft mass
  364.      * @param provider for un-normalized zonal coefficients
  365.      * @param initialType initial orbit type (mean Brouwer-Lyddane orbit or osculating orbit)
  366.      * @param M2 value of empirical drag coefficient in rad/s².
  367.      *        If equal to {@link BrouwerLyddanePropagator#M2} drag is not computed
  368.      */
  369.     public FieldBrouwerLyddanePropagator(final FieldOrbit<T> initialOrbit,
  370.                                          final AttitudeProvider attitudeProv,
  371.                                          final T mass,
  372.                                          final UnnormalizedSphericalHarmonicsProvider provider,
  373.                                          final PropagationType initialType,
  374.                                          final double M2) {
  375.         this(initialOrbit, attitudeProv, mass, provider, provider.onDate(initialOrbit.getDate().toAbsoluteDate()), initialType, M2);
  376.     }

  377.     /**
  378.      * Private helper constructor.
  379.      * <p>Using this constructor, it is possible to define the initial orbit as
  380.      * a mean Brouwer-Lyddane orbit or an osculating one.</p>
  381.      * @param initialOrbit initial orbit
  382.      * @param attitude attitude provider
  383.      * @param mass spacecraft mass
  384.      * @param provider for un-normalized zonal coefficients
  385.      * @param harmonics {@code provider.onDate(initialOrbit.getDate())}
  386.      * @param initialType initial orbit type (mean Brouwer-Lyddane orbit or osculating orbit)
  387.      * @param M2 value of empirical drag coefficient in rad/s².
  388.      *        If equal to {@link BrouwerLyddanePropagator#M2} drag is not computed
  389.      */
  390.     public FieldBrouwerLyddanePropagator(final FieldOrbit<T> initialOrbit,
  391.                                          final AttitudeProvider attitude,
  392.                                          final T mass,
  393.                                          final UnnormalizedSphericalHarmonicsProvider provider,
  394.                                          final UnnormalizedSphericalHarmonics harmonics,
  395.                                          final PropagationType initialType,
  396.                                          final double M2) {
  397.         this(initialOrbit, attitude, mass, provider.getAe(), initialOrbit.getMu().newInstance(provider.getMu()),
  398.              harmonics.getUnnormalizedCnm(2, 0),
  399.              harmonics.getUnnormalizedCnm(3, 0),
  400.              harmonics.getUnnormalizedCnm(4, 0),
  401.              harmonics.getUnnormalizedCnm(5, 0),
  402.              initialType, M2);
  403.     }

  404.     /** Build a propagator from orbit, attitude provider, mass and potential.
  405.      * <p>The C<sub>n,0</sub> coefficients are the denormalized zonal coefficients, they
  406.      * are related to both the normalized coefficients
  407.      * <span style="text-decoration: overline">C</span><sub>n,0</sub>
  408.      *  and the J<sub>n</sub> one as follows:</p>
  409.      *
  410.      * <p> C<sub>n,0</sub> = [(2-δ<sub>0,m</sub>)(2n+1)(n-m)!/(n+m)!]<sup>½</sup>
  411.      * <span style="text-decoration: overline">C</span><sub>n,0</sub>
  412.      *
  413.      * <p> C<sub>n,0</sub> = -J<sub>n</sub>
  414.      *
  415.      * <p>Using this constructor, it is possible to define the initial orbit as
  416.      * a mean Brouwer-Lyddane orbit or an osculating one.</p>
  417.      *
  418.      * @param initialOrbit initial orbit
  419.      * @param attitudeProv attitude provider
  420.      * @param mass spacecraft mass
  421.      * @param referenceRadius reference radius of the Earth for the potential model (m)
  422.      * @param mu central attraction coefficient (m³/s²)
  423.      * @param c20 un-normalized zonal coefficient (about -1.08e-3 for Earth)
  424.      * @param c30 un-normalized zonal coefficient (about +2.53e-6 for Earth)
  425.      * @param c40 un-normalized zonal coefficient (about +1.62e-6 for Earth)
  426.      * @param c50 un-normalized zonal coefficient (about +2.28e-7 for Earth)
  427.      * @param initialType initial orbit type (mean Brouwer-Lyddane orbit or osculating orbit)
  428.      * @param M2 value of empirical drag coefficient in rad/s².
  429.      *        If equal to {@link BrouwerLyddanePropagator#M2} drag is not computed
  430.      */
  431.     public FieldBrouwerLyddanePropagator(final FieldOrbit<T> initialOrbit,
  432.                                          final AttitudeProvider attitudeProv,
  433.                                          final T mass,
  434.                                          final double referenceRadius, final T mu,
  435.                                          final double c20, final double c30, final double c40,
  436.                                          final double c50,
  437.                                          final PropagationType initialType,
  438.                                          final double M2) {
  439.         this(initialOrbit, attitudeProv, mass, referenceRadius, mu,
  440.              c20, c30, c40, c50, initialType, M2, EPSILON_DEFAULT, MAX_ITERATIONS_DEFAULT);
  441.     }

  442.     /** Build a propagator from orbit, attitude provider, mass and potential.
  443.      * <p>The C<sub>n,0</sub> coefficients are the denormalized zonal coefficients, they
  444.      * are related to both the normalized coefficients
  445.      * <span style="text-decoration: overline">C</span><sub>n,0</sub>
  446.      *  and the J<sub>n</sub> one as follows:</p>
  447.      *
  448.      * <p> C<sub>n,0</sub> = [(2-δ<sub>0,m</sub>)(2n+1)(n-m)!/(n+m)!]<sup>½</sup>
  449.      * <span style="text-decoration: overline">C</span><sub>n,0</sub>
  450.      *
  451.      * <p> C<sub>n,0</sub> = -J<sub>n</sub>
  452.      *
  453.      * <p>Using this constructor, it is possible to define the initial orbit as
  454.      * a mean Brouwer-Lyddane orbit or an osculating one.</p>
  455.      *
  456.      * @param initialOrbit initial orbit
  457.      * @param attitudeProv attitude provider
  458.      * @param mass spacecraft mass
  459.      * @param referenceRadius reference radius of the Earth for the potential model (m)
  460.      * @param mu central attraction coefficient (m³/s²)
  461.      * @param c20 un-normalized zonal coefficient (about -1.08e-3 for Earth)
  462.      * @param c30 un-normalized zonal coefficient (about +2.53e-6 for Earth)
  463.      * @param c40 un-normalized zonal coefficient (about +1.62e-6 for Earth)
  464.      * @param c50 un-normalized zonal coefficient (about +2.28e-7 for Earth)
  465.      * @param initialType initial orbit type (mean Brouwer-Lyddane orbit or osculating orbit)
  466.      * @param M2 value of empirical drag coefficient in rad/s².
  467.      *        If equal to {@link BrouwerLyddanePropagator#M2} drag is not computed
  468.      * @param epsilon convergence threshold for mean parameters conversion
  469.      * @param maxIterations maximum iterations for mean parameters conversion
  470.      * @since 11.2
  471.      */
  472.     public FieldBrouwerLyddanePropagator(final FieldOrbit<T> initialOrbit,
  473.                                          final AttitudeProvider attitudeProv,
  474.                                          final T mass,
  475.                                          final double referenceRadius, final T mu,
  476.                                          final double c20, final double c30, final double c40,
  477.                                          final double c50,
  478.                                          final PropagationType initialType,
  479.                                          final double M2, final double epsilon, final int maxIterations) {

  480.         super(mass.getField(), attitudeProv);

  481.         // store model coefficients
  482.         this.referenceRadius = referenceRadius;
  483.         this.mu  = mu;
  484.         this.ck0 = new double[] {0.0, 0.0, c20, c30, c40, c50};

  485.         // initialize M2 driver
  486.         this.M2Driver = new ParameterDriver(BrouwerLyddanePropagator.M2_NAME, M2, SCALE,
  487.                                             Double.NEGATIVE_INFINITY,
  488.                                             Double.POSITIVE_INFINITY);

  489.         // compute mean parameters if needed
  490.         resetInitialState(new FieldSpacecraftState<>(initialOrbit,
  491.                                                  attitudeProv.getAttitude(initialOrbit,
  492.                                                                           initialOrbit.getDate(),
  493.                                                                           initialOrbit.getFrame()),
  494.                                                  mass),
  495.                                                  initialType, epsilon, maxIterations);

  496.     }

  497.     /** Conversion from osculating to mean orbit.
  498.      * <p>
  499.      * Compute mean orbit <b>in a Brouwer-Lyddane sense</b>, corresponding to the
  500.      * osculating SpacecraftState in input.
  501.      * </p>
  502.      * <p>
  503.      * Since the osculating orbit is obtained with the computation of
  504.      * short-periodic variation, the resulting output will depend on
  505.      * both the gravity field parameterized in input and the
  506.      * atmospheric drag represented by the {@code m2} parameter.
  507.      * </p>
  508.      * <p>
  509.      * The computation is done through a fixed-point iteration process.
  510.      * </p>
  511.      * @param <T> type of the filed elements
  512.      * @param osculating osculating orbit to convert
  513.      * @param provider for un-normalized zonal coefficients
  514.      * @param harmonics {@code provider.onDate(osculating.getDate())}
  515.      * @param M2Value value of empirical drag coefficient in rad/s².
  516.      *        If equal to {@code BrouwerLyddanePropagator.M2} drag is not considered
  517.      * @return mean orbit in a Brouwer-Lyddane sense
  518.      * @since 11.2
  519.      */
  520.     public static <T extends CalculusFieldElement<T>> FieldKeplerianOrbit<T> computeMeanOrbit(final FieldOrbit<T> osculating,
  521.                                                                                               final UnnormalizedSphericalHarmonicsProvider provider,
  522.                                                                                               final UnnormalizedSphericalHarmonics harmonics,
  523.                                                                                               final double M2Value) {
  524.         return computeMeanOrbit(osculating, provider, harmonics, M2Value, EPSILON_DEFAULT, MAX_ITERATIONS_DEFAULT);
  525.     }

  526.     /** Conversion from osculating to mean orbit.
  527.      * <p>
  528.      * Compute mean orbit <b>in a Brouwer-Lyddane sense</b>, corresponding to the
  529.      * osculating SpacecraftState in input.
  530.      * </p>
  531.      * <p>
  532.      * Since the osculating orbit is obtained with the computation of
  533.      * short-periodic variation, the resulting output will depend on
  534.      * both the gravity field parameterized in input and the
  535.      * atmospheric drag represented by the {@code m2} parameter.
  536.      * </p>
  537.      * <p>
  538.      * The computation is done through a fixed-point iteration process.
  539.      * </p>
  540.      * @param <T> type of the filed elements
  541.      * @param osculating osculating orbit to convert
  542.      * @param provider for un-normalized zonal coefficients
  543.      * @param harmonics {@code provider.onDate(osculating.getDate())}
  544.      * @param M2Value value of empirical drag coefficient in rad/s².
  545.      *        If equal to {@code BrouwerLyddanePropagator.M2} drag is not considered
  546.      * @param epsilon convergence threshold for mean parameters conversion
  547.      * @param maxIterations maximum iterations for mean parameters conversion
  548.      * @return mean orbit in a Brouwer-Lyddane sense
  549.      * @since 11.2
  550.      */
  551.     public static <T extends CalculusFieldElement<T>> FieldKeplerianOrbit<T> computeMeanOrbit(final FieldOrbit<T> osculating,
  552.                                                                                               final UnnormalizedSphericalHarmonicsProvider provider,
  553.                                                                                               final UnnormalizedSphericalHarmonics harmonics,
  554.                                                                                               final double M2Value,
  555.                                                                                               final double epsilon, final int maxIterations) {
  556.         return computeMeanOrbit(osculating,
  557.                                 provider.getAe(), provider.getMu(),
  558.                                 harmonics.getUnnormalizedCnm(2, 0),
  559.                                 harmonics.getUnnormalizedCnm(3, 0),
  560.                                 harmonics.getUnnormalizedCnm(4, 0),
  561.                                 harmonics.getUnnormalizedCnm(5, 0),
  562.                                 M2Value, epsilon, maxIterations);
  563.     }

  564.     /** Conversion from osculating to mean orbit.
  565.      * <p>
  566.      * Compute mean orbit <b>in a Brouwer-Lyddane sense</b>, corresponding to the
  567.      * osculating SpacecraftState in input.
  568.      * </p>
  569.      * <p>
  570.      * Since the osculating orbit is obtained with the computation of
  571.      * short-periodic variation, the resulting output will depend on
  572.      * both the gravity field parameterized in input and the
  573.      * atmospheric drag represented by the {@code m2} parameter.
  574.      * </p>
  575.      * <p>
  576.      * The computation is done through a fixed-point iteration process.
  577.      * </p>
  578.      * @param <T> type of the filed elements
  579.      * @param osculating osculating orbit to convert
  580.      * @param referenceRadius reference radius of the Earth for the potential model (m)
  581.      * @param mu central attraction coefficient (m³/s²)
  582.      * @param c20 un-normalized zonal coefficient (about -1.08e-3 for Earth)
  583.      * @param c30 un-normalized zonal coefficient (about +2.53e-6 for Earth)
  584.      * @param c40 un-normalized zonal coefficient (about +1.62e-6 for Earth)
  585.      * @param c50 un-normalized zonal coefficient (about +2.28e-7 for Earth)
  586.      * @param M2Value value of empirical drag coefficient in rad/s².
  587.      *        If equal to {@code BrouwerLyddanePropagator.M2} drag is not considered
  588.      * @param epsilon convergence threshold for mean parameters conversion
  589.      * @param maxIterations maximum iterations for mean parameters conversion
  590.      * @return mean orbit in a Brouwer-Lyddane sense
  591.      * @since 11.2
  592.      */
  593.     public static <T extends CalculusFieldElement<T>> FieldKeplerianOrbit<T> computeMeanOrbit(final FieldOrbit<T> osculating,
  594.                                                                                              final double referenceRadius, final double mu,
  595.                                                                                              final double c20, final double c30, final double c40,
  596.                                                                                              final double c50, final double M2Value,
  597.                                                                                              final double epsilon, final int maxIterations) {
  598.         final FieldBrouwerLyddanePropagator<T> propagator =
  599.                         new FieldBrouwerLyddanePropagator<>(osculating,
  600.                                                             InertialProvider.of(osculating.getFrame()),
  601.                                                             osculating.getMu().newInstance(DEFAULT_MASS),
  602.                                                             referenceRadius, osculating.getMu().newInstance(mu),
  603.                                                             c20, c30, c40, c50,
  604.                                                             PropagationType.OSCULATING, M2Value,
  605.                                                             epsilon, maxIterations);
  606.         return propagator.initialModel.mean;
  607.     }

  608.     /** {@inheritDoc}
  609.      * <p>The new initial state to consider
  610.      * must be defined with an osculating orbit.</p>
  611.      * @see #resetInitialState(FieldSpacecraftState, PropagationType)
  612.      */
  613.     @Override
  614.     public void resetInitialState(final FieldSpacecraftState<T> state) {
  615.         resetInitialState(state, PropagationType.OSCULATING);
  616.     }

  617.     /** Reset the propagator initial state.
  618.      * @param state new initial state to consider
  619.      * @param stateType mean Brouwer-Lyddane orbit or osculating orbit
  620.      */
  621.     public void resetInitialState(final FieldSpacecraftState<T> state, final PropagationType stateType) {
  622.         resetInitialState(state, stateType, EPSILON_DEFAULT, MAX_ITERATIONS_DEFAULT);
  623.     }

  624.     /** Reset the propagator initial state.
  625.      * @param state new initial state to consider
  626.      * @param stateType mean Brouwer-Lyddane orbit or osculating orbit
  627.      * @param epsilon convergence threshold for mean parameters conversion
  628.      * @param maxIterations maximum iterations for mean parameters conversion
  629.      * @since 11.2
  630.      */
  631.     public void resetInitialState(final FieldSpacecraftState<T> state, final PropagationType stateType,
  632.                                   final double epsilon, final int maxIterations) {
  633.         super.resetInitialState(state);
  634.         final FieldKeplerianOrbit<T> keplerian = (FieldKeplerianOrbit<T>) OrbitType.KEPLERIAN.convertType(state.getOrbit());
  635.         this.initialModel = (stateType == PropagationType.MEAN) ?
  636.                              new FieldBLModel<>(keplerian, state.getMass(), referenceRadius, mu, ck0) :
  637.                              computeMeanParameters(keplerian, state.getMass(), epsilon, maxIterations);
  638.         this.models = new FieldTimeSpanMap<FieldBLModel<T>, T>(initialModel, state.getA().getField());
  639.     }

  640.     /** {@inheritDoc} */
  641.     @Override
  642.     protected void resetIntermediateState(final FieldSpacecraftState<T> state, final boolean forward) {
  643.         resetIntermediateState(state, forward, EPSILON_DEFAULT, MAX_ITERATIONS_DEFAULT);
  644.     }

  645.     /** Reset an intermediate state.
  646.      * @param state new intermediate state to consider
  647.      * @param forward if true, the intermediate state is valid for
  648.      * propagations after itself
  649.      * @param epsilon convergence threshold for mean parameters conversion
  650.      * @param maxIterations maximum iterations for mean parameters conversion
  651.      * @since 11.2
  652.      */
  653.     protected void resetIntermediateState(final FieldSpacecraftState<T> state, final boolean forward,
  654.                                           final double epsilon, final int maxIterations) {
  655.         final FieldBLModel<T> newModel = computeMeanParameters((FieldKeplerianOrbit<T>) OrbitType.KEPLERIAN.convertType(state.getOrbit()),
  656.                                                                state.getMass(), epsilon, maxIterations);
  657.         if (forward) {
  658.             models.addValidAfter(newModel, state.getDate());
  659.         } else {
  660.             models.addValidBefore(newModel, state.getDate());
  661.         }
  662.         stateChanged(state);
  663.     }

  664.     /** Compute mean parameters according to the Brouwer-Lyddane analytical model computation
  665.      * in order to do the propagation.
  666.      * @param osculating osculating orbit
  667.      * @param mass constant mass
  668.      * @param epsilon convergence threshold for mean parameters conversion
  669.      * @param maxIterations maximum iterations for mean parameters conversion
  670.      * @return Brouwer-Lyddane mean model
  671.      */
  672.     private FieldBLModel<T> computeMeanParameters(final FieldKeplerianOrbit<T> osculating, final T mass,
  673.                                                   final double epsilon, final int maxIterations) {

  674.         // sanity check
  675.         if (osculating.getA().getReal() < referenceRadius) {
  676.             throw new OrekitException(OrekitMessages.TRAJECTORY_INSIDE_BRILLOUIN_SPHERE,
  677.                                            osculating.getA());
  678.         }

  679.         final Field<T> field = mass.getField();
  680.         final T one = field.getOne();
  681.         final T zero = field.getZero();
  682.         // rough initialization of the mean parameters
  683.         FieldBLModel<T> current = new FieldBLModel<T>(osculating, mass, referenceRadius, mu, ck0);

  684.         // threshold for each parameter
  685.         final T thresholdA      = current.mean.getA().abs().add(1.0).multiply(epsilon);
  686.         final T thresholdE      = current.mean.getE().add(1.0).multiply(epsilon);
  687.         final T thresholdAngles = one.getPi().multiply(epsilon);

  688.         int i = 0;
  689.         while (i++ < maxIterations) {

  690.             // recompute the osculating parameters from the current mean parameters
  691.             final FieldKeplerianOrbit<T> parameters = current.propagateParameters(current.mean.getDate(), getParameters(mass.getField()));

  692.             // adapted parameters residuals
  693.             final T deltaA     = osculating.getA()  .subtract(parameters.getA());
  694.             final T deltaE     = osculating.getE()  .subtract(parameters.getE());
  695.             final T deltaI     = osculating.getI()  .subtract(parameters.getI());
  696.             final T deltaOmega = MathUtils.normalizeAngle(osculating.getPerigeeArgument().subtract(parameters.getPerigeeArgument()), zero);
  697.             final T deltaRAAN  = MathUtils.normalizeAngle(osculating.getRightAscensionOfAscendingNode().subtract(parameters.getRightAscensionOfAscendingNode()), zero);
  698.             final T deltaAnom  = MathUtils.normalizeAngle(osculating.getMeanAnomaly().subtract(parameters.getMeanAnomaly()), zero);


  699.             // update mean parameters
  700.             current = new FieldBLModel<T>(new FieldKeplerianOrbit<T>(current.mean.getA()            .add(deltaA),
  701.                                                      FastMath.max(current.mean.getE().add(deltaE), zero),
  702.                                                      current.mean.getI()                            .add(deltaI),
  703.                                                      current.mean.getPerigeeArgument()              .add(deltaOmega),
  704.                                                      current.mean.getRightAscensionOfAscendingNode().add(deltaRAAN),
  705.                                                      current.mean.getMeanAnomaly()                  .add(deltaAnom),
  706.                                                      PositionAngle.MEAN,
  707.                                                      current.mean.getFrame(),
  708.                                                      current.mean.getDate(), mu),
  709.                                   mass, referenceRadius, mu, ck0);
  710.             // check convergence
  711.             if (FastMath.abs(deltaA.getReal())     < thresholdA.getReal() &&
  712.                 FastMath.abs(deltaE.getReal())     < thresholdE.getReal() &&
  713.                 FastMath.abs(deltaI.getReal())     < thresholdAngles.getReal() &&
  714.                 FastMath.abs(deltaOmega.getReal()) < thresholdAngles.getReal() &&
  715.                 FastMath.abs(deltaRAAN.getReal())  < thresholdAngles.getReal() &&
  716.                 FastMath.abs(deltaAnom.getReal())  < thresholdAngles.getReal()) {
  717.                 return current;
  718.             }
  719.         }
  720.         throw new OrekitException(OrekitMessages.UNABLE_TO_COMPUTE_BROUWER_LYDDANE_MEAN_PARAMETERS, i);
  721.     }


  722.     /** {@inheritDoc} */
  723.     public FieldKeplerianOrbit<T> propagateOrbit(final FieldAbsoluteDate<T> date, final T[] parameters) {
  724.         // compute Cartesian parameters, taking derivatives into account
  725.         final FieldBLModel<T> current = models.get(date);
  726.         return current.propagateParameters(date, parameters);
  727.     }

  728.     /**
  729.      * Get the value of the M2 drag parameter.
  730.      * @return the value of the M2 drag parameter
  731.      */
  732.     public double getM2() {
  733.         return M2Driver.getValue();
  734.     }

  735.     /** Local class for Brouwer-Lyddane model. */
  736.     private static class FieldBLModel<T extends CalculusFieldElement<T>> {

  737.         /** Mean orbit. */
  738.         private final FieldKeplerianOrbit<T> mean;

  739.         /** Constant mass. */
  740.         private final T mass;

  741.         /** Central attraction coefficient. */
  742.         private final T mu;

  743.         // CHECKSTYLE: stop JavadocVariable check

  744.         // preprocessed values

  745.         // Constant for secular terms l'', g'', h''
  746.         // l standing for true anomaly, g for perigee argument and h for raan
  747.         private final T xnotDot;
  748.         private final T n;
  749.         private final T lt;
  750.         private final T gt;
  751.         private final T ht;


  752.         // Long periodT
  753.         private final T dei3sg;
  754.         private final T de2sg;
  755.         private final T deisg;
  756.         private final T de;


  757.         private final T dlgs2g;
  758.         private final T dlgc3g;
  759.         private final T dlgcg;


  760.         private final T dh2sgcg;
  761.         private final T dhsgcg;
  762.         private final T dhcg;


  763.         // Short perioTs
  764.         private final T aC;
  765.         private final T aCbis;
  766.         private final T ac2g2f;


  767.         private final T eC;
  768.         private final T ecf;
  769.         private final T e2cf;
  770.         private final T e3cf;
  771.         private final T ec2f2g;
  772.         private final T ecfc2f2g;
  773.         private final T e2cfc2f2g;
  774.         private final T e3cfc2f2g;
  775.         private final T ec2gf;
  776.         private final T ec2g3f;


  777.         private final T ide;
  778.         private final T isfs2f2g;
  779.         private final T icfc2f2g;
  780.         private final T ic2f2g;


  781.         private final T glf;
  782.         private final T gll;
  783.         private final T glsf;
  784.         private final T glosf;
  785.         private final T gls2f2g;
  786.         private final T gls2gf;
  787.         private final T glos2gf;
  788.         private final T gls2g3f;
  789.         private final T glos2g3f;


  790.         private final T hf;
  791.         private final T hl;
  792.         private final T hsf;
  793.         private final T hcfs2g2f;
  794.         private final T hs2g2f;
  795.         private final T hsfc2g2f;


  796.         private final T edls2g;
  797.         private final T edlcg;
  798.         private final T edlc3g;
  799.         private final T edlsf;
  800.         private final T edls2gf;
  801.         private final T edls2g3f;

  802.         // Drag terms
  803.         private final T aRate;
  804.         private final T eRate;

  805.         // CHECKSTYLE: resume JavadocVariable check

  806.         /** Create a model for specified mean orbit.
  807.          * @param mean mean Fieldorbit
  808.          * @param mass constant mass
  809.          * @param referenceRadius reference radius of the central body attraction model (m)
  810.          * @param mu central attraction coefficient (m³/s²)
  811.          * @param ck0 un-normalized zonal coefficients
  812.          */
  813.         FieldBLModel(final FieldKeplerianOrbit<T> mean, final T mass,
  814.                 final double referenceRadius, final T mu, final double[] ck0) {

  815.             this.mean = mean;
  816.             this.mass = mass;
  817.             this.mu   = mu;
  818.             final T one  = mass.getField().getOne();

  819.             final T app = mean.getA();
  820.             xnotDot = mu.divide(app).sqrt().divide(app);

  821.             // preliminary processing
  822.             final T q = app.divide(referenceRadius).reciprocal();
  823.             T ql = q.multiply(q);
  824.             final T y2 = ql.multiply(-0.5 * ck0[2]);

  825.             n = ((mean.getE().multiply(mean.getE()).negate()).add(1.0)).sqrt();
  826.             final T n2 = n.multiply(n);
  827.             final T n3 = n2.multiply(n);
  828.             final T n4 = n2.multiply(n2);
  829.             final T n6 = n4.multiply(n2);
  830.             final T n8 = n4.multiply(n4);
  831.             final T n10 = n8.multiply(n2);

  832.             final T yp2 = y2.divide(n4);
  833.             ql = ql.multiply(q);
  834.             final T yp3 = ql.multiply(ck0[3]).divide(n6);
  835.             ql = ql.multiply(q);
  836.             final T yp4 = ql.multiply(0.375 * ck0[4]).divide(n8);
  837.             ql = ql.multiply(q);
  838.             final T yp5 = ql.multiply(ck0[5]).divide(n10);


  839.             final FieldSinCos<T> sc = FastMath.sinCos(mean.getI());
  840.             final T sinI1 = sc.sin();
  841.             final T sinI2 = sinI1.multiply(sinI1);
  842.             final T cosI1 = sc.cos();
  843.             final T cosI2 = cosI1.multiply(cosI1);
  844.             final T cosI3 = cosI2.multiply(cosI1);
  845.             final T cosI4 = cosI2.multiply(cosI2);
  846.             final T cosI6 = cosI4.multiply(cosI2);
  847.             final T C5c2 = T2(cosI1).reciprocal();
  848.             final T C3c2 = cosI2.multiply(3.0).subtract(1.0);

  849.             final T epp = mean.getE();
  850.             final T epp2 = epp.multiply(epp);
  851.             final T epp3 = epp2.multiply(epp);
  852.             final T epp4 = epp2.multiply(epp2);

  853.             if (epp.getReal() >= 1) {
  854.                 // Only for elliptical (e < 1) orbits
  855.                 throw new OrekitException(OrekitMessages.TOO_LARGE_ECCENTRICITY_FOR_PROPAGATION_MODEL,
  856.                                           mean.getE().getReal());
  857.             }

  858.             // secular multiplicative
  859.             lt = one.add(yp2.multiply(n).multiply(C3c2).multiply(1.5)).
  860.                      add(yp2.multiply(0.09375).multiply(yp2).multiply(n).multiply(n2.multiply(25.0).add(n.multiply(16.0)).add(-15.0).
  861.                              add(cosI2.multiply(n2.multiply(-90.0).add(n.multiply(-96.0)).add(30.0))).
  862.                              add(cosI4.multiply(n2.multiply(25.0).add(n.multiply(144.0)).add(105.0))))).
  863.                      add(yp4.multiply(0.9375).multiply(n).multiply(epp2).multiply(cosI4.multiply(35.0).add(cosI2.multiply(-30.0)).add(3.0)));

  864.             gt = yp2.multiply(-1.5).multiply(C5c2).
  865.                      add(yp2.multiply(0.09375).multiply(yp2).multiply(n2.multiply(25.0).add(n.multiply(24.0)).add(-35.0).
  866.                             add(cosI2.multiply(n2.multiply(-126.0).add(n.multiply(-192.0)).add(90.0))).
  867.                             add(cosI4.multiply(n2.multiply(45.0).add(n.multiply(360.0)).add(385.0))))).
  868.                      add(yp4.multiply(0.3125).multiply(n2.multiply(-9.0).add(21.0).
  869.                             add(cosI2.multiply(n2.multiply(126.0).add(-270.0))).
  870.                             add(cosI4.multiply(n2.multiply(-189.0).add(385.0)))));

  871.             ht = yp2.multiply(-3.0).multiply(cosI1).
  872.                      add(yp2.multiply(0.375).multiply(yp2).multiply(cosI1.multiply(n2.multiply(9.0).add(n.multiply(12.0)).add(-5.0)).
  873.                                                                     add(cosI3.multiply(n2.multiply(-5.0).add(n.multiply(-36.0)).add(-35.0))))).
  874.                      add(yp4.multiply(1.25).multiply(cosI1).multiply(n2.multiply(-3.0).add(5.0)).multiply(cosI2.multiply(-7.0).add(3.0)));

  875.             final T cA = one.subtract(cosI2.multiply(11.0)).subtract(cosI4.multiply(40.0).divide(C5c2));
  876.             final T cB = one.subtract(cosI2.multiply(3.0)).subtract(cosI4.multiply(8.0).divide(C5c2));
  877.             final T cC = one.subtract(cosI2.multiply(9.0)).subtract(cosI4.multiply(24.0).divide(C5c2));
  878.             final T cD = one.subtract(cosI2.multiply(5.0)).subtract(cosI4.multiply(16.0).divide(C5c2));

  879.             final T qyp2_4 = yp2.multiply(3.0).multiply(yp2).multiply(cA).
  880.                              subtract(yp4.multiply(10.0).multiply(cB));
  881.             final T qyp52 = cosI1.multiply(epp3).multiply(cD.multiply(0.5).divide(sinI1).
  882.                                                           add(sinI1.multiply(cosI4.divide(C5c2).divide(C5c2).multiply(80.0).
  883.                                                                              add(cosI2.divide(C5c2).multiply(32.0).
  884.                                                                              add(5.0)))));
  885.             final T qyp22 = epp2.add(2.0).subtract(epp2.multiply(3.0).add(2.0).multiply(11.0).multiply(cosI2)).
  886.                             subtract(epp2.multiply(5.0).add(2.0).multiply(40.0).multiply(cosI4.divide(C5c2))).
  887.                             subtract(epp2.multiply(400.0).multiply(cosI6).divide(C5c2).divide(C5c2));
  888.             final T qyp42 = one.divide(5.0).multiply(qyp22.
  889.                                                      add(one.multiply(4.0).multiply(epp2.
  890.                                                                                     add(2.0).
  891.                                                                                     subtract(cosI2.multiply(epp2.multiply(3.0).add(2.0))))));
  892.             final T qyp52bis = cosI1.multiply(sinI1).multiply(epp).multiply(epp2.multiply(3.0).add(4.0)).
  893.                                                                    multiply(cosI4.divide(C5c2).divide(C5c2).multiply(40.0).
  894.                                                                             add(cosI2.divide(C5c2).multiply(16.0)).
  895.                                                                             add(3.0));
  896.            // long periodic multiplicative
  897.             dei3sg =  yp5.divide(yp2).multiply(35.0 / 96.0).multiply(epp2).multiply(n2).multiply(cD).multiply(sinI1);
  898.             de2sg = qyp2_4.divide(yp2).multiply(epp).multiply(n2).multiply(-1.0 / 12.0);
  899.             deisg = sinI1.multiply(yp5.multiply(-35.0 / 128.0).divide(yp2).multiply(epp2).multiply(n2).multiply(cD).
  900.                             add(n2.multiply(0.25).divide(yp2).multiply(yp3.add(yp5.multiply(5.0 / 16.0).multiply(epp2.multiply(3.0).add(4.0)).multiply(cC)))));
  901.             de = epp2.multiply(n2).multiply(qyp2_4).divide(24.0).divide(yp2);

  902.             final T qyp52quotient = epp.multiply(epp4.multiply(81.0).add(-32.0)).divide(n.multiply(epp2.multiply(9.0).add(4.0)).add(epp2.multiply(3.0)).add(4.0));
  903.             dlgs2g = yp2.multiply(48.0).reciprocal().multiply(yp2.multiply(-3.0).multiply(yp2).multiply(qyp22).
  904.                             add(yp4.multiply(10.0).multiply(qyp42))).
  905.                             add(n3.divide(yp2).multiply(qyp2_4).divide(24.0));
  906.             dlgc3g =  yp5.multiply(35.0 / 384.0).divide(yp2).multiply(n3).multiply(epp).multiply(cD).multiply(sinI1).
  907.                             add(yp5.multiply(35.0 / 1152.0).divide(yp2).multiply(qyp52.multiply(2.0).multiply(cosI1).subtract(epp.multiply(cD).multiply(sinI1).multiply(epp2.multiply(2.0).add(3.0)))));
  908.             dlgcg = yp3.negate().multiply(cosI2).multiply(epp).divide(yp2.multiply(sinI1).multiply(4.0)).
  909.                     add(yp5.divide(yp2).multiply(0.078125).multiply(cC).multiply(cosI2.divide(sinI1).multiply(epp.negate()).multiply(epp2.multiply(3.0).add(4.0)).
  910.                                                                     add(sinI1.multiply(epp2).multiply(epp2.multiply(9.0).add(26.0)))).
  911.                     subtract(yp5.divide(yp2).multiply(0.46875).multiply(qyp52bis).multiply(cosI1)).
  912.                     add(yp3.divide(yp2).multiply(0.25).multiply(sinI1).multiply(epp).divide(n3.add(1.0)).multiply(epp2.multiply(epp2.subtract(3.0)).add(3.0))).
  913.                     add(yp5.divide(yp2).multiply(0.078125).multiply(n2).multiply(cC).multiply(qyp52quotient).multiply(sinI1)));
  914.             final T qyp24 = yp2.multiply(yp2).multiply(3.0).multiply(cosI4.divide(sinI2).multiply(200.0).add(cosI2.divide(sinI1).multiply(80.0)).add(11.0)).
  915.                             subtract(yp4.multiply(10.0).multiply(cosI4.divide(sinI2).multiply(40.0).add(cosI2.divide(sinI1).multiply(16.0)).add(3.0)));
  916.             dh2sgcg = yp5.divide(yp2).multiply(qyp52).multiply(35.0 / 144.0);
  917.             dhsgcg = qyp24.multiply(cosI1).multiply(epp2.negate()).divide(yp2.multiply(12.0));
  918.             dhcg = yp5.divide(yp2).multiply(qyp52).multiply(-35.0 / 576.0).
  919.                    add(cosI1.multiply(epp).divide(yp2.multiply(sinI1).multiply(4.0)).multiply(yp3.add(yp5.multiply(0.3125).multiply(cC).multiply(epp2.multiply(3.0).add(4.0))))).
  920.                    add(yp5.multiply(qyp52bis).multiply(1.875).divide(yp2.multiply(4.0)));

  921.             // short periodic multiplicative
  922.             aC = yp2.negate().multiply(C3c2).multiply(app).divide(n3);
  923.             aCbis = y2.multiply(app).multiply(C3c2);
  924.             ac2g2f = y2.multiply(app).multiply(sinI2).multiply(3.0);

  925.             T qe = y2.multiply(C3c2).multiply(0.5).multiply(n2).divide(n6);
  926.             eC = qe.multiply(epp).divide(n3.add(1.0)).multiply(epp2.multiply(epp2.subtract(3.0)).add(3.0));
  927.             ecf = qe.multiply(3.0);
  928.             e2cf = qe.multiply(3.0).multiply(epp);
  929.             e3cf = qe.multiply(epp2);
  930.             qe = y2.multiply(0.5).multiply(n2).multiply(3.0).multiply(cosI2.negate().add(1.0)).divide(n6);
  931.             ec2f2g = qe.multiply(epp);
  932.             ecfc2f2g = qe.multiply(3.0);
  933.             e2cfc2f2g = qe.multiply(3.0).multiply(epp);
  934.             e3cfc2f2g = qe.multiply(epp2);
  935.             qe = yp2.multiply(-0.5).multiply(n2).multiply(cosI2.negate().add(1.0));
  936.             ec2gf = qe.multiply(3.0);
  937.             ec2g3f = qe;

  938.             T qi = yp2.multiply(epp).multiply(cosI1).multiply(sinI1);
  939.             ide = cosI1.multiply(epp.negate()).divide(sinI1.multiply(n2));
  940.             isfs2f2g = qi;
  941.             icfc2f2g = qi.multiply(2.0);
  942.             qi = yp2.multiply(cosI1).multiply(sinI1);
  943.             ic2f2g = qi.multiply(1.5);

  944.             T qgl1 = yp2.multiply(0.25);
  945.             T qgl2 =  yp2.multiply(epp).multiply(n2).multiply(0.25).divide(n.add(1.0));
  946.             glf = qgl1.multiply(C5c2).multiply(-6.0);
  947.             gll = qgl1.multiply(C5c2).multiply(6.0);
  948.             glsf = qgl1.multiply(C5c2).multiply(-6.0).multiply(epp).
  949.                    add(qgl2.multiply(C3c2).multiply(2.0));
  950.             glosf = qgl2.multiply(C3c2).multiply(2.0);
  951.             qgl1 = qgl1.multiply(cosI2.multiply(-5.0).add(3.0));
  952.             qgl2 = qgl2.multiply(3.0).multiply(cosI2.negate().add(1.0));
  953.             gls2f2g = qgl1.multiply(3.0);
  954.             gls2gf = qgl1.multiply(3.0).multiply(epp).
  955.                      add(qgl2);
  956.             glos2gf = qgl2.negate();
  957.             gls2g3f = qgl1.multiply(epp).
  958.                       add(qgl2.multiply(1.0 / 3.0));
  959.             glos2g3f = qgl2;

  960.             final T qh = yp2.multiply(cosI1).multiply(3.0);
  961.             hf = qh.negate();
  962.             hl = qh;
  963.             hsf = qh.multiply(epp).negate();
  964.             hcfs2g2f = yp2.multiply(cosI1).multiply(epp).multiply(2.0);
  965.             hs2g2f = yp2.multiply(cosI1).multiply(1.5);
  966.             hsfc2g2f = yp2.multiply(cosI1).multiply(epp).negate();

  967.             final T qedl = yp2.multiply(n3).multiply(-0.25);
  968.             edls2g = qyp2_4.multiply(1.0 / 24.0).multiply(epp).multiply(n3).divide(yp2);
  969.             edlcg = yp3.divide(yp2).multiply(-0.25).multiply(n3).multiply(sinI1).
  970.                     subtract(yp5.divide(yp2).multiply(0.078125).multiply(n3).multiply(sinI1).multiply(cC).multiply(epp2.multiply(9.0).add(4.0)));
  971.             edlc3g = yp5.divide(yp2).multiply(n3).multiply(epp2).multiply(cD).multiply(sinI1).multiply(35.0 / 384.0);
  972.             edlsf = qedl.multiply(C3c2).multiply(2.0);
  973.             edls2gf = qedl.multiply(3.0).multiply(cosI2.negate().add(1.0));
  974.             edls2g3f = qedl.multiply(1.0 / 3.0);

  975.             // secular rates of the mean semi-major axis and eccentricity
  976.             // Eq. 2.41 and Eq. 2.45 of Phipps' 1992 thesis
  977.             aRate = app.multiply(-4.0).divide(xnotDot.multiply(3.0));
  978.             eRate = epp.multiply(n).multiply(n).multiply(-4.0).divide(xnotDot.multiply(3.0));

  979.         }

  980.         /**
  981.          * Accurate computation of E - e sin(E).
  982.          *
  983.          * @param E eccentric anomaly
  984.          * @return E - e sin(E)
  985.          */
  986.         private FieldUnivariateDerivative2<T> eMeSinE(final FieldUnivariateDerivative2<T> E) {
  987.             FieldUnivariateDerivative2<T> x = E.sin().multiply(mean.getE().negate().add(1.0));
  988.             final FieldUnivariateDerivative2<T> mE2 = E.negate().multiply(E);
  989.             FieldUnivariateDerivative2<T> term = E;
  990.             FieldUnivariateDerivative2<T> d    = E.getField().getZero();
  991.             // the inequality test below IS intentional and should NOT be replaced by a check with a small tolerance
  992.             for (FieldUnivariateDerivative2<T> x0 = d.add(Double.NaN); !Double.valueOf(x.getValue().getReal()).equals(Double.valueOf(x0.getValue().getReal()));) {
  993.                 d = d.add(2);
  994.                 term = term.multiply(mE2.divide(d.multiply(d.add(1))));
  995.                 x0 = x;
  996.                 x = x.subtract(term);
  997.             }
  998.             return x;
  999.         }

  1000.         /**
  1001.          * Gets eccentric anomaly from mean anomaly.
  1002.          * <p>The algorithm used to solve the Kepler equation has been published in:
  1003.          * "Procedures for  solving Kepler's Equation", A. W. Odell and R. H. Gooding,
  1004.          * Celestial Mechanics 38 (1986) 307-334</p>
  1005.          * <p>It has been copied from the OREKIT library (KeplerianOrbit class).</p>
  1006.          *
  1007.          * @param mk the mean anomaly (rad)
  1008.          * @return the eccentric anomaly (rad)
  1009.          */
  1010.         private FieldUnivariateDerivative2<T> getEccentricAnomaly(final FieldUnivariateDerivative2<T> mk) {
  1011.             final double k1 = 3 * FastMath.PI + 2;
  1012.             final double k2 = FastMath.PI - 1;
  1013.             final double k3 = 6 * FastMath.PI - 1;
  1014.             final double A  = 3.0 * k2 * k2 / k1;
  1015.             final double B  = k3 * k3 / (6.0 * k1);

  1016.             final T zero = mean.getE().getField().getZero();

  1017.             // reduce M to [-PI PI] interval
  1018.             final FieldUnivariateDerivative2<T> reducedM = new FieldUnivariateDerivative2<T>(MathUtils.normalizeAngle(mk.getValue(), zero ),
  1019.                                                                              mk.getFirstDerivative(),
  1020.                                                                              mk.getSecondDerivative());

  1021.             // compute start value according to A. W. Odell and R. H. Gooding S12 starter
  1022.             FieldUnivariateDerivative2<T> ek;
  1023.             if (reducedM.getValue().abs().getReal() < 1.0 / 6.0) {
  1024.                 if (reducedM.getValue().abs().getReal() < Precision.SAFE_MIN) {
  1025.                     // this is an Orekit change to the S12 starter.
  1026.                     // If reducedM is 0.0, the derivative of cbrt is infinite which induces NaN appearing later in
  1027.                     // the computation. As in this case E and M are almost equal, we initialize ek with reducedM
  1028.                     ek = reducedM;
  1029.                 } else {
  1030.                     // this is the standard S12 starter
  1031.                     ek = reducedM.add(reducedM.multiply(6).cbrt().subtract(reducedM).multiply(mean.getE()));
  1032.                 }
  1033.             } else {
  1034.                 if (reducedM.getValue().getReal() < 0) {
  1035.                     final FieldUnivariateDerivative2<T> w = reducedM.add(FastMath.PI);
  1036.                     ek = reducedM.add(w.multiply(-A).divide(w.subtract(B)).subtract(FastMath.PI).subtract(reducedM).multiply(mean.getE()));
  1037.                 } else {
  1038.                     final FieldUnivariateDerivative2<T> minusW = reducedM.subtract(FastMath.PI);
  1039.                     ek = reducedM.add(minusW.multiply(A).divide(minusW.add(B)).add(FastMath.PI).subtract(reducedM).multiply(mean.getE()));
  1040.                 }
  1041.             }

  1042.             final T e1 = mean.getE().negate().add(1.0);
  1043.             final boolean noCancellationRisk = (e1.add(ek.getValue().multiply(ek.getValue())).getReal() / 6) >= 0.1;

  1044.             // perform two iterations, each consisting of one Halley step and one Newton-Raphson step
  1045.             for (int j = 0; j < 2; ++j) {
  1046.                 final FieldUnivariateDerivative2<T> f;
  1047.                 FieldUnivariateDerivative2<T> fd;
  1048.                 final FieldUnivariateDerivative2<T> fdd  = ek.sin().multiply(mean.getE());
  1049.                 final FieldUnivariateDerivative2<T> fddd = ek.cos().multiply(mean.getE());
  1050.                 if (noCancellationRisk) {
  1051.                     f  = ek.subtract(fdd).subtract(reducedM);
  1052.                     fd = fddd.subtract(1).negate();
  1053.                 } else {
  1054.                     f  = eMeSinE(ek).subtract(reducedM);
  1055.                     final FieldUnivariateDerivative2<T> s = ek.multiply(0.5).sin();
  1056.                     fd = s.multiply(s).multiply(mean.getE().multiply(2.0)).add(e1);
  1057.                 }
  1058.                 final FieldUnivariateDerivative2<T> dee = f.multiply(fd).divide(f.multiply(0.5).multiply(fdd).subtract(fd.multiply(fd)));

  1059.                 // update eccentric anomaly, using expressions that limit underflow problems
  1060.                 final FieldUnivariateDerivative2<T> w = fd.add(dee.multiply(0.5).multiply(fdd.add(dee.multiply(fdd).divide(3))));
  1061.                 fd = fd.add(dee.multiply(fdd.add(dee.multiply(0.5).multiply(fdd))));
  1062.                 ek = ek.subtract(f.subtract(dee.multiply(fd.subtract(w))).divide(fd));
  1063.             }

  1064.             // expand the result back to original range
  1065.             ek = ek.add(mk.getValue().subtract(reducedM.getValue()));

  1066.             // Returns the eccentric anomaly
  1067.             return ek;
  1068.         }

  1069.         /**
  1070.          * This method is used in Brouwer-Lyddane model to avoid singularity at the
  1071.          * critical inclination (i = 63.4°).
  1072.          * <p>
  1073.          * This method, based on Warren Phipps's 1992 thesis (Eq. 2.47 and 2.48),
  1074.          * approximate the factor (1.0 - 5.0 * cos²(inc))^-1 (causing the singularity)
  1075.          * by a function, named T2 in the thesis.
  1076.          * </p>
  1077.          * @param cosInc cosine of the mean inclination
  1078.          * @return an approximation of (1.0 - 5.0 * cos²(inc))^-1 term
  1079.          */
  1080.         private T T2(final T cosInc) {

  1081.             // X = (1.0 - 5.0 * cos²(inc))
  1082.             final T x  = cosInc.multiply(cosInc).multiply(-5.0).add(1.0);
  1083.             final T x2 = x.multiply(x);

  1084.             // Eq. 2.48
  1085.             T sum = x.getField().getZero();
  1086.             for (int i = 0; i <= 12; i++) {
  1087.                 final double sign = i % 2 == 0 ? +1.0 : -1.0;
  1088.                 sum = sum.add(FastMath.pow(x2, i).multiply(FastMath.pow(BETA, i)).multiply(sign).divide(CombinatoricsUtils.factorialDouble(i + 1)));
  1089.             }

  1090.             // Right term of equation 2.47
  1091.             T product = x.getField().getOne();
  1092.             for (int i = 0; i <= 10; i++) {
  1093.                 product = product.multiply(FastMath.exp(x2.multiply(BETA).multiply(FastMath.scalb(-1.0, i))).add(1.0));
  1094.             }

  1095.             // Return (Eq. 2.47)
  1096.             return x.multiply(BETA).multiply(sum).multiply(product);

  1097.         }

  1098.         /** Extrapolate an orbit up to a specific target date.
  1099.          * @param date target date for the orbit
  1100.          * @param parameters model parameters
  1101.          * @return propagated parameters
  1102.          */
  1103.         public FieldKeplerianOrbit<T> propagateParameters(final FieldAbsoluteDate<T> date, final T[] parameters) {

  1104.             // Field
  1105.             final Field<T> field = date.getField();
  1106.             final T one  = field.getOne();
  1107.             final T zero = field.getZero();

  1108.             // Empirical drag coefficient M2
  1109.             final T m2 = parameters[0];

  1110.             // Keplerian evolution
  1111.             final FieldUnivariateDerivative2<T> dt = new FieldUnivariateDerivative2<>(date.durationFrom(mean.getDate()), one, zero);
  1112.             final FieldUnivariateDerivative2<T> xnot = dt.multiply(xnotDot);

  1113.             //____________________________________
  1114.             // secular effects

  1115.             // mean mean anomaly
  1116.             final FieldUnivariateDerivative2<T> dtM2  = dt.multiply(m2);
  1117.             final FieldUnivariateDerivative2<T> dt2M2 = dt.multiply(dtM2);
  1118.             final FieldUnivariateDerivative2<T> lpp = new FieldUnivariateDerivative2<T>(MathUtils.normalizeAngle(mean.getMeanAnomaly().add(lt.multiply(xnot.getValue())).add(dt2M2.getValue()), zero),
  1119.                                                                                         lt.multiply(xnotDot).add(dtM2.multiply(2.0).getValue()),
  1120.                                                                                         m2.multiply(2.0));
  1121.             // mean argument of perigee
  1122.             final FieldUnivariateDerivative2<T> gpp = new FieldUnivariateDerivative2<T>(MathUtils.normalizeAngle(mean.getPerigeeArgument().add(gt.multiply(xnot.getValue())), zero),
  1123.                                                                                         gt.multiply(xnotDot),
  1124.                                                                                         zero);
  1125.             // mean longitude of ascending node
  1126.             final FieldUnivariateDerivative2<T> hpp = new FieldUnivariateDerivative2<T>(MathUtils.normalizeAngle(mean.getRightAscensionOfAscendingNode().add(ht.multiply(xnot.getValue())), zero),
  1127.                                                                                         ht.multiply(xnotDot),
  1128.                                                                                         zero);

  1129.             // ________________________________________________
  1130.             // secular rates of the mean semi-major axis and eccentricity

  1131.             // semi-major axis
  1132.             final FieldUnivariateDerivative2<T> appDrag = dt.multiply(aRate.multiply(m2));

  1133.             // eccentricity
  1134.             final FieldUnivariateDerivative2<T> eppDrag = dt.multiply(eRate.multiply(m2));

  1135.             //____________________________________
  1136.             // Long periodical terms
  1137.             final FieldUnivariateDerivative2<T> cg1 = gpp.cos();
  1138.             final FieldUnivariateDerivative2<T> sg1 = gpp.sin();
  1139.             final FieldUnivariateDerivative2<T> c2g = cg1.multiply(cg1).subtract(sg1.multiply(sg1));
  1140.             final FieldUnivariateDerivative2<T> s2g = cg1.multiply(sg1).add(sg1.multiply(cg1));
  1141.             final FieldUnivariateDerivative2<T> c3g = c2g.multiply(cg1).subtract(s2g.multiply(sg1));
  1142.             final FieldUnivariateDerivative2<T> sg2 = sg1.multiply(sg1);
  1143.             final FieldUnivariateDerivative2<T> sg3 = sg1.multiply(sg2);


  1144.             // de eccentricity
  1145.             final FieldUnivariateDerivative2<T> d1e = sg3.multiply(dei3sg).
  1146.                                                add(sg1.multiply(deisg)).
  1147.                                                add(sg2.multiply(de2sg)).
  1148.                                                add(de);

  1149.             // l' + g'
  1150.             final FieldUnivariateDerivative2<T> lp_p_gp = s2g.multiply(dlgs2g).
  1151.                                                add(c3g.multiply(dlgc3g)).
  1152.                                                add(cg1.multiply(dlgcg)).
  1153.                                                add(lpp).
  1154.                                                add(gpp);

  1155.             // h'
  1156.             final FieldUnivariateDerivative2<T> hp = sg2.multiply(cg1).multiply(dh2sgcg).
  1157.                                                add(sg1.multiply(cg1).multiply(dhsgcg)).
  1158.                                                add(cg1.multiply(dhcg)).
  1159.                                                add(hpp);

  1160.             // Short periodical terms
  1161.             // eccentric anomaly
  1162.             final FieldUnivariateDerivative2<T> Ep = getEccentricAnomaly(lpp);
  1163.             final FieldUnivariateDerivative2<T> cf1 = (Ep.cos().subtract(mean.getE())).divide(Ep.cos().multiply(mean.getE().negate()).add(1.0));
  1164.             final FieldUnivariateDerivative2<T> sf1 = (Ep.sin().multiply(n)).divide(Ep.cos().multiply(mean.getE().negate()).add(1.0));
  1165.             final FieldUnivariateDerivative2<T> f = FastMath.atan2(sf1, cf1);

  1166.             final FieldUnivariateDerivative2<T> c2f = cf1.multiply(cf1).subtract(sf1.multiply(sf1));
  1167.             final FieldUnivariateDerivative2<T> s2f = cf1.multiply(sf1).add(sf1.multiply(cf1));
  1168.             final FieldUnivariateDerivative2<T> c3f = c2f.multiply(cf1).subtract(s2f.multiply(sf1));
  1169.             final FieldUnivariateDerivative2<T> s3f = c2f.multiply(sf1).add(s2f.multiply(cf1));
  1170.             final FieldUnivariateDerivative2<T> cf2 = cf1.multiply(cf1);
  1171.             final FieldUnivariateDerivative2<T> cf3 = cf1.multiply(cf2);

  1172.             final FieldUnivariateDerivative2<T> c2g1f = cf1.multiply(c2g).subtract(sf1.multiply(s2g));
  1173.             final FieldUnivariateDerivative2<T> c2g2f = c2f.multiply(c2g).subtract(s2f.multiply(s2g));
  1174.             final FieldUnivariateDerivative2<T> c2g3f = c3f.multiply(c2g).subtract(s3f.multiply(s2g));
  1175.             final FieldUnivariateDerivative2<T> s2g1f = cf1.multiply(s2g).add(c2g.multiply(sf1));
  1176.             final FieldUnivariateDerivative2<T> s2g2f = c2f.multiply(s2g).add(c2g.multiply(s2f));
  1177.             final FieldUnivariateDerivative2<T> s2g3f = c3f.multiply(s2g).add(c2g.multiply(s3f));

  1178.             final FieldUnivariateDerivative2<T> eE = (Ep.cos().multiply(mean.getE().negate()).add(1.0)).reciprocal();
  1179.             final FieldUnivariateDerivative2<T> eE3 = eE.multiply(eE).multiply(eE);
  1180.             final FieldUnivariateDerivative2<T> sigma = eE.multiply(n.multiply(n)).multiply(eE).add(eE);

  1181.             // Semi-major axis
  1182.             final FieldUnivariateDerivative2<T> a = eE3.multiply(aCbis).add(appDrag.add(mean.getA())).
  1183.                                             add(aC).
  1184.                                             add(eE3.multiply(c2g2f).multiply(ac2g2f));

  1185.             // Eccentricity
  1186.             final FieldUnivariateDerivative2<T> e = d1e.add(eppDrag.add(mean.getE())).
  1187.                                             add(eC).
  1188.                                             add(cf1.multiply(ecf)).
  1189.                                             add(cf2.multiply(e2cf)).
  1190.                                             add(cf3.multiply(e3cf)).
  1191.                                             add(c2g2f.multiply(ec2f2g)).
  1192.                                             add(c2g2f.multiply(cf1).multiply(ecfc2f2g)).
  1193.                                             add(c2g2f.multiply(cf2).multiply(e2cfc2f2g)).
  1194.                                             add(c2g2f.multiply(cf3).multiply(e3cfc2f2g)).
  1195.                                             add(c2g1f.multiply(ec2gf)).
  1196.                                             add(c2g3f.multiply(ec2g3f));

  1197.             // Inclination
  1198.             final FieldUnivariateDerivative2<T> i = d1e.multiply(ide).
  1199.                                             add(mean.getI()).
  1200.                                             add(sf1.multiply(s2g2f.multiply(isfs2f2g))).
  1201.                                             add(cf1.multiply(c2g2f.multiply(icfc2f2g))).
  1202.                                             add(c2g2f.multiply(ic2f2g));

  1203.             // Argument of perigee + True anomaly
  1204.             final FieldUnivariateDerivative2<T> g_p_l = lp_p_gp.add(f.multiply(glf)).
  1205.                                              add(lpp.multiply(gll)).
  1206.                                              add(sf1.multiply(glsf)).
  1207.                                              add(sigma.multiply(sf1).multiply(glosf)).
  1208.                                              add(s2g2f.multiply(gls2f2g)).
  1209.                                              add(s2g1f.multiply(gls2gf)).
  1210.                                              add(sigma.multiply(s2g1f).multiply(glos2gf)).
  1211.                                              add(s2g3f.multiply(gls2g3f)).
  1212.                                              add(sigma.multiply(s2g3f).multiply(glos2g3f));


  1213.             // Longitude of ascending node
  1214.             final FieldUnivariateDerivative2<T> h = hp.add(f.multiply(hf)).
  1215.                                             add(lpp.multiply(hl)).
  1216.                                             add(sf1.multiply(hsf)).
  1217.                                             add(cf1.multiply(s2g2f).multiply(hcfs2g2f)).
  1218.                                             add(s2g2f.multiply(hs2g2f)).
  1219.                                             add(c2g2f.multiply(sf1).multiply(hsfc2g2f));

  1220.             final FieldUnivariateDerivative2<T> edl = s2g.multiply(edls2g).
  1221.                                             add(cg1.multiply(edlcg)).
  1222.                                             add(c3g.multiply(edlc3g)).
  1223.                                             add(sf1.multiply(edlsf)).
  1224.                                             add(s2g1f.multiply(edls2gf)).
  1225.                                             add(s2g3f.multiply(edls2g3f)).
  1226.                                             add(sf1.multiply(sigma).multiply(edlsf)).
  1227.                                             add(s2g1f.multiply(sigma).multiply(edls2gf.negate())).
  1228.                                             add(s2g3f.multiply(sigma).multiply(edls2g3f.multiply(3.0)));

  1229.             final FieldUnivariateDerivative2<T> A = e.multiply(lpp.cos()).subtract(edl.multiply(lpp.sin()));
  1230.             final FieldUnivariateDerivative2<T> B = e.multiply(lpp.sin()).add(edl.multiply(lpp.cos()));

  1231.             // True anomaly
  1232.             final FieldUnivariateDerivative2<T> l = FastMath.atan2(B, A);

  1233.             // Argument of perigee
  1234.             final FieldUnivariateDerivative2<T> g = g_p_l.subtract(l);

  1235.             // Return a keplerian orbit
  1236.             return new FieldKeplerianOrbit<>(a.getValue(), e.getValue(), i.getValue(),
  1237.                                              g.getValue(), h.getValue(), l.getValue(),
  1238.                                              a.getFirstDerivative(), e.getFirstDerivative(), i.getFirstDerivative(),
  1239.                                              g.getFirstDerivative(), h.getFirstDerivative(), l.getFirstDerivative(),
  1240.                                              PositionAngle.MEAN, mean.getFrame(), date, this.mu);

  1241.         }
  1242.     }

  1243.     /** {@inheritDoc} */
  1244.     @Override
  1245.     protected T getMass(final FieldAbsoluteDate<T> date) {
  1246.         return models.get(date).mass;
  1247.     }

  1248.     /** {@inheritDoc} */
  1249.     @Override
  1250.     protected List<ParameterDriver> getParametersDrivers() {
  1251.         return Collections.singletonList(M2Driver);
  1252.     }

  1253. }