BrouwerLyddanePropagator.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.ArrayList;
  19. import java.util.Collections;
  20. import java.util.List;

  21. import org.hipparchus.analysis.differentiation.UnivariateDerivative2;
  22. import org.hipparchus.linear.RealMatrix;
  23. import org.hipparchus.util.CombinatoricsUtils;
  24. import org.hipparchus.util.FastMath;
  25. import org.hipparchus.util.MathUtils;
  26. import org.hipparchus.util.Precision;
  27. import org.hipparchus.util.SinCos;
  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.KeplerianOrbit;
  35. import org.orekit.orbits.Orbit;
  36. import org.orekit.orbits.OrbitType;
  37. import org.orekit.orbits.PositionAngle;
  38. import org.orekit.propagation.AbstractMatricesHarvester;
  39. import org.orekit.propagation.MatricesHarvester;
  40. import org.orekit.propagation.PropagationType;
  41. import org.orekit.propagation.SpacecraftState;
  42. import org.orekit.propagation.analytical.tle.TLE;
  43. import org.orekit.time.AbsoluteDate;
  44. import org.orekit.utils.DoubleArrayDictionary;
  45. import org.orekit.utils.ParameterDriver;
  46. import org.orekit.utils.TimeSpanMap;

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

  88.     /** Parameter name for M2 coefficient. */
  89.     public static final String M2_NAME = "M2";

  90.     /** Default value for M2 coefficient. */
  91.     public static final double M2 = 0.0;

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

  94.     /** Default value for maxIterations. */
  95.     private static final int MAX_ITERATIONS_DEFAULT = 200;

  96.     /** Parameters scaling factor.
  97.      * <p>
  98.      * We use a power of 2 to avoid numeric noise introduction
  99.      * in the multiplications/divisions sequences.
  100.      * </p>
  101.      */
  102.     private static final double SCALE = FastMath.scalb(1.0, -32);

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

  105.     /** Initial Brouwer-Lyddane model. */
  106.     private BLModel initialModel;

  107.     /** All models. */
  108.     private transient TimeSpanMap<BLModel> models;

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

  111.     /** Central attraction coefficient (m³/s²). */
  112.     private double mu;

  113.     /** Un-normalized zonal coefficients. */
  114.     private double[] ck0;

  115.     /** Empirical coefficient used in the drag modeling. */
  116.     private final ParameterDriver M2Driver;

  117.     /** Build a propagator from orbit and potential provider.
  118.      * <p>Mass and attitude provider are set to unspecified non-null arbitrary values.</p>
  119.      *
  120.      * <p>Using this constructor, an initial osculating orbit is considered.</p>
  121.      *
  122.      * @param initialOrbit initial orbit
  123.      * @param provider for un-normalized zonal coefficients
  124.      * @param M2 value of empirical drag coefficient in rad/s².
  125.      *        If equal to {@link #M2} drag is not computed
  126.      * @see #BrouwerLyddanePropagator(Orbit, AttitudeProvider, UnnormalizedSphericalHarmonicsProvider, double)
  127.      * @see #BrouwerLyddanePropagator(Orbit, UnnormalizedSphericalHarmonicsProvider, PropagationType, double)
  128.      */
  129.     public BrouwerLyddanePropagator(final Orbit initialOrbit,
  130.                                     final UnnormalizedSphericalHarmonicsProvider provider,
  131.                                     final double M2) {
  132.         this(initialOrbit, InertialProvider.of(initialOrbit.getFrame()),
  133.              DEFAULT_MASS, provider, provider.onDate(initialOrbit.getDate()), M2);
  134.     }

  135.     /**
  136.      * Private helper constructor.
  137.      * <p>Using this constructor, an initial osculating orbit is considered.</p>
  138.      * @param initialOrbit initial orbit
  139.      * @param attitude attitude provider
  140.      * @param mass spacecraft mass
  141.      * @param provider for un-normalized zonal coefficients
  142.      * @param harmonics {@code provider.onDate(initialOrbit.getDate())}
  143.      * @param M2 value of empirical drag coefficient in rad/s².
  144.      *        If equal to {@link #M2} drag is not computed
  145.      * @see #BrouwerLyddanePropagator(Orbit, AttitudeProvider, double,
  146.      *                                 UnnormalizedSphericalHarmonicsProvider,
  147.      *                                 UnnormalizedSphericalHarmonicsProvider.UnnormalizedSphericalHarmonics,
  148.      *                                 PropagationType, double)
  149.      */
  150.     public BrouwerLyddanePropagator(final Orbit initialOrbit,
  151.                                     final AttitudeProvider attitude,
  152.                                     final double mass,
  153.                                     final UnnormalizedSphericalHarmonicsProvider provider,
  154.                                     final UnnormalizedSphericalHarmonics harmonics,
  155.                                     final double M2) {
  156.         this(initialOrbit, attitude, mass, provider.getAe(), provider.getMu(),
  157.              harmonics.getUnnormalizedCnm(2, 0),
  158.              harmonics.getUnnormalizedCnm(3, 0),
  159.              harmonics.getUnnormalizedCnm(4, 0),
  160.              harmonics.getUnnormalizedCnm(5, 0),
  161.              M2);
  162.     }

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

  197.     /** Build a propagator from orbit, mass and potential provider.
  198.      * <p>Attitude law is set to an unspecified non-null arbitrary value.</p>
  199.      *
  200.      * <p>Using this constructor, an initial osculating orbit is considered.</p>
  201.      *
  202.      * @param initialOrbit initial orbit
  203.      * @param mass spacecraft mass
  204.      * @param provider for un-normalized zonal coefficients
  205.      * @param M2 value of empirical drag coefficient in rad/s².
  206.      *        If equal to {@link #M2} drag is not computed
  207.      * @see #BrouwerLyddanePropagator(Orbit, AttitudeProvider, double, UnnormalizedSphericalHarmonicsProvider, double)
  208.      */
  209.     public BrouwerLyddanePropagator(final Orbit initialOrbit, final double mass,
  210.                                     final UnnormalizedSphericalHarmonicsProvider provider,
  211.                                     final double M2) {
  212.         this(initialOrbit, InertialProvider.of(initialOrbit.getFrame()),
  213.              mass, provider, provider.onDate(initialOrbit.getDate()), M2);
  214.     }

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

  249.     /** Build a propagator from orbit, attitude provider and potential provider.
  250.      * <p>Mass is set to an unspecified non-null arbitrary value.</p>
  251.      * <p>Using this constructor, an initial osculating orbit is considered.</p>
  252.      * @param initialOrbit initial orbit
  253.      * @param attitudeProv attitude provider
  254.      * @param provider for un-normalized zonal coefficients
  255.      * @param M2 value of empirical drag coefficient in rad/s².
  256.      *        If equal to {@link #M2} drag is not computed
  257.      */
  258.     public BrouwerLyddanePropagator(final Orbit initialOrbit,
  259.                                     final AttitudeProvider attitudeProv,
  260.                                     final UnnormalizedSphericalHarmonicsProvider provider,
  261.                                     final double M2) {
  262.         this(initialOrbit, attitudeProv, DEFAULT_MASS, provider, provider.onDate(initialOrbit.getDate()), M2);
  263.     }

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

  296.     /** Build a propagator from orbit, attitude provider, mass and potential provider.
  297.      * <p>Using this constructor, an initial osculating orbit is considered.</p>
  298.      * @param initialOrbit initial orbit
  299.      * @param attitudeProv attitude provider
  300.      * @param mass spacecraft mass
  301.      * @param provider for un-normalized zonal coefficients
  302.      * @param M2 value of empirical drag coefficient in rad/s².
  303.      *        If equal to {@link #M2} drag is not computed
  304.      * @see #BrouwerLyddanePropagator(Orbit, AttitudeProvider, double,
  305.      *                                UnnormalizedSphericalHarmonicsProvider, PropagationType, double)
  306.      */
  307.     public BrouwerLyddanePropagator(final Orbit initialOrbit,
  308.                                     final AttitudeProvider attitudeProv,
  309.                                     final double mass,
  310.                                     final UnnormalizedSphericalHarmonicsProvider provider,
  311.                                     final double M2) {
  312.         this(initialOrbit, attitudeProv, mass, provider, provider.onDate(initialOrbit.getDate()), M2);
  313.     }

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


  350.     /** Build a propagator from orbit and potential provider.
  351.      * <p>Mass and attitude provider are set to unspecified non-null arbitrary values.</p>
  352.      *
  353.      * <p>Using this constructor, it is possible to define the initial orbit as
  354.      * a mean Brouwer-Lyddane orbit or an osculating one.</p>
  355.      *
  356.      * @param initialOrbit initial orbit
  357.      * @param provider for un-normalized zonal coefficients
  358.      * @param initialType initial orbit type (mean Brouwer-Lyddane orbit or osculating orbit)
  359.      * @param M2 value of empirical drag coefficient in rad/s².
  360.      *        If equal to {@link #M2} drag is not computed
  361.      */
  362.     public BrouwerLyddanePropagator(final Orbit initialOrbit,
  363.                                     final UnnormalizedSphericalHarmonicsProvider provider,
  364.                                     final PropagationType initialType, final double M2) {
  365.         this(initialOrbit, InertialProvider.of(initialOrbit.getFrame()),
  366.              DEFAULT_MASS, provider, provider.onDate(initialOrbit.getDate()), initialType, M2);
  367.     }

  368.     /** Build a propagator from orbit, attitude provider, mass and potential provider.
  369.      * <p>Using this constructor, it is possible to define the initial orbit as
  370.      * a mean Brouwer-Lyddane orbit or an osculating one.</p>
  371.      * @param initialOrbit initial orbit
  372.      * @param attitudeProv attitude provider
  373.      * @param mass spacecraft mass
  374.      * @param provider for un-normalized zonal coefficients
  375.      * @param initialType initial orbit type (mean Brouwer-Lyddane orbit or osculating orbit)
  376.      * @param M2 value of empirical drag coefficient in rad/s².
  377.      *        If equal to {@link #M2} drag is not computed
  378.      */
  379.     public BrouwerLyddanePropagator(final Orbit initialOrbit,
  380.                                     final AttitudeProvider attitudeProv,
  381.                                     final double mass,
  382.                                     final UnnormalizedSphericalHarmonicsProvider provider,
  383.                                     final PropagationType initialType, final double M2) {
  384.         this(initialOrbit, attitudeProv, mass, provider, provider.onDate(initialOrbit.getDate()), initialType, M2);
  385.     }

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

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

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

  487.         super(attitudeProv);

  488.         // store model coefficients
  489.         this.referenceRadius = referenceRadius;
  490.         this.mu  = mu;
  491.         this.ck0 = new double[] {0.0, 0.0, c20, c30, c40, c50};

  492.         // initialize M2 driver
  493.         this.M2Driver = new ParameterDriver(M2_NAME, M2, SCALE,
  494.                                             Double.NEGATIVE_INFINITY,
  495.                                             Double.POSITIVE_INFINITY);

  496.         // compute mean parameters if needed
  497.         resetInitialState(new SpacecraftState(initialOrbit,
  498.                                               attitudeProv.getAttitude(initialOrbit,
  499.                                                                        initialOrbit.getDate(),
  500.                                                                        initialOrbit.getFrame()),
  501.                                               mass),
  502.                           initialType, epsilon, maxIterations);

  503.     }

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

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

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

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

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

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

  642.     /** {@inheritDoc} */
  643.     protected void resetIntermediateState(final SpacecraftState state, final boolean forward) {
  644.         resetIntermediateState(state, forward, EPSILON_DEFAULT, MAX_ITERATIONS_DEFAULT);
  645.     }

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

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

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

  680.         // rough initialization of the mean parameters
  681.         BLModel current = new BLModel(osculating, mass, referenceRadius, mu, ck0);

  682.         // threshold for each parameter
  683.         final double thresholdA      = epsilon * (1 + FastMath.abs(current.mean.getA()));
  684.         final double thresholdE      = epsilon * (1 + current.mean.getE());
  685.         final double thresholdAngles = epsilon * FastMath.PI;

  686.         int i = 0;
  687.         while (i++ < maxIterations) {

  688.             // recompute the osculating parameters from the current mean parameters
  689.             final KeplerianOrbit parameters = current.propagateParameters(current.mean.getDate());

  690.             // adapted parameters residuals
  691.             final double deltaA     = osculating.getA() - parameters.getA();
  692.             final double deltaE     = osculating.getE() - parameters.getE();
  693.             final double deltaI     = osculating.getI() - parameters.getI();
  694.             final double deltaOmega = MathUtils.normalizeAngle(osculating.getPerigeeArgument() -
  695.                                                                parameters.getPerigeeArgument(),
  696.                                                                0.0);
  697.             final double deltaRAAN  = MathUtils.normalizeAngle(osculating.getRightAscensionOfAscendingNode() -
  698.                                                                parameters.getRightAscensionOfAscendingNode(),
  699.                                                                0.0);
  700.             final double deltaAnom = MathUtils.normalizeAngle(osculating.getMeanAnomaly() -
  701.                                                               parameters.getMeanAnomaly(),
  702.                                                               0.0);


  703.             // update mean parameters
  704.             current = new BLModel(new KeplerianOrbit(current.mean.getA() + deltaA,
  705.                                                      FastMath.max(current.mean.getE() + deltaE, 0.0),
  706.                                                      current.mean.getI() + deltaI,
  707.                                                      current.mean.getPerigeeArgument() + deltaOmega,
  708.                                                      current.mean.getRightAscensionOfAscendingNode() + deltaRAAN,
  709.                                                      current.mean.getMeanAnomaly() + deltaAnom,
  710.                                                      PositionAngle.MEAN,
  711.                                                      current.mean.getFrame(),
  712.                                                      current.mean.getDate(), mu),
  713.                                   mass, referenceRadius, mu, ck0);
  714.             // check convergence
  715.             if (FastMath.abs(deltaA)     < thresholdA &&
  716.                 FastMath.abs(deltaE)     < thresholdE &&
  717.                 FastMath.abs(deltaI)     < thresholdAngles &&
  718.                 FastMath.abs(deltaOmega) < thresholdAngles &&
  719.                 FastMath.abs(deltaRAAN)  < thresholdAngles &&
  720.                 FastMath.abs(deltaAnom)  < thresholdAngles) {
  721.                 return current;
  722.             }
  723.         }
  724.         throw new OrekitException(OrekitMessages.UNABLE_TO_COMPUTE_BROUWER_LYDDANE_MEAN_PARAMETERS, i);
  725.     }


  726.     /** {@inheritDoc} */
  727.     public KeplerianOrbit propagateOrbit(final AbsoluteDate date) {
  728.         // compute keplerian parameters, taking derivatives into account
  729.         final BLModel current = models.get(date);
  730.         return current.propagateParameters(date);
  731.     }

  732.     /**
  733.      * Get the value of the M2 drag parameter.
  734.      * @return the value of the M2 drag parameter
  735.      */
  736.     public double getM2() {
  737.         return M2Driver.getValue();
  738.     }

  739.     /**
  740.      * Get the central attraction coefficient μ.
  741.      * @return mu central attraction coefficient (m³/s²)
  742.      */
  743.     public double getMu() {
  744.         return mu;
  745.     }

  746.     /**
  747.      * Get the un-normalized zonal coefficients.
  748.      * @return the un-normalized zonal coefficients
  749.      */
  750.     public double[] getCk0() {
  751.         return ck0.clone();
  752.     }

  753.     /**
  754.      * Get the reference radius of the central body attraction model.
  755.      * @return the reference radius in meters
  756.      */
  757.     public double getReferenceRadius() {
  758.         return referenceRadius;
  759.     }

  760.     /**
  761.      * Get the parameters driver for propagation model.
  762.      * @return drivers for propagation model
  763.      */
  764.     public List<ParameterDriver> getParametersDrivers() {
  765.         return Collections.singletonList(M2Driver);
  766.     }

  767.     /** {@inheritDoc} */
  768.     @Override
  769.     protected AbstractMatricesHarvester createHarvester(final String stmName, final RealMatrix initialStm,
  770.                                                         final DoubleArrayDictionary initialJacobianColumns) {
  771.         // Create the harvester
  772.         final BrouwerLyddaneHarvester harvester = new BrouwerLyddaneHarvester(this, stmName, initialStm, initialJacobianColumns);
  773.         // Update the list of additional state provider
  774.         addAdditionalStateProvider(harvester);
  775.         // Return the configured harvester
  776.         return harvester;
  777.     }

  778.     /**
  779.      * Get the names of the parameters in the matrix returned by {@link MatricesHarvester#getParametersJacobian}.
  780.      * @return names of the parameters (i.e. columns) of the Jacobian matrix
  781.      */
  782.     protected List<String> getJacobiansColumnsNames() {
  783.         final List<String> columnsNames = new ArrayList<>();
  784.         for (final ParameterDriver driver : getParametersDrivers()) {
  785.             if (driver.isSelected() && !columnsNames.contains(driver.getName())) {
  786.                 columnsNames.add(driver.getName());
  787.             }
  788.         }
  789.         Collections.sort(columnsNames);
  790.         return columnsNames;
  791.     }

  792.     /** Local class for Brouwer-Lyddane model. */
  793.     private class BLModel {

  794.         /** Mean orbit. */
  795.         private final KeplerianOrbit mean;

  796.         /** Constant mass. */
  797.         private final double mass;

  798.         // CHECKSTYLE: stop JavadocVariable check

  799.         // preprocessed values

  800.         // Constant for secular terms l'', g'', h''
  801.         // l standing for true anomaly, g for perigee argument and h for raan
  802.         private final double xnotDot;
  803.         private final double n;
  804.         private final double lt;
  805.         private final double gt;
  806.         private final double ht;


  807.         // Long period terms
  808.         private final double dei3sg;
  809.         private final double de2sg;
  810.         private final double deisg;
  811.         private final double de;


  812.         private final double dlgs2g;
  813.         private final double dlgc3g;
  814.         private final double dlgcg;


  815.         private final double dh2sgcg;
  816.         private final double dhsgcg;
  817.         private final double dhcg;


  818.         // Short period terms
  819.         private final double aC;
  820.         private final double aCbis;
  821.         private final double ac2g2f;


  822.         private final double eC;
  823.         private final double ecf;
  824.         private final double e2cf;
  825.         private final double e3cf;
  826.         private final double ec2f2g;
  827.         private final double ecfc2f2g;
  828.         private final double e2cfc2f2g;
  829.         private final double e3cfc2f2g;
  830.         private final double ec2gf;
  831.         private final double ec2g3f;


  832.         private final double ide;
  833.         private final double isfs2f2g;
  834.         private final double icfc2f2g;
  835.         private final double ic2f2g;


  836.         private final double glf;
  837.         private final double gll;
  838.         private final double glsf;
  839.         private final double glosf;
  840.         private final double gls2f2g;
  841.         private final double gls2gf;
  842.         private final double glos2gf;
  843.         private final double gls2g3f;
  844.         private final double glos2g3f;


  845.         private final double hf;
  846.         private final double hl;
  847.         private final double hsf;
  848.         private final double hcfs2g2f;
  849.         private final double hs2g2f;
  850.         private final double hsfc2g2f;


  851.         private final double edls2g;
  852.         private final double edlcg;
  853.         private final double edlc3g;
  854.         private final double edlsf;
  855.         private final double edls2gf;
  856.         private final double edls2g3f;

  857.         // Drag terms
  858.         private final double aRate;
  859.         private final double eRate;

  860.         // CHECKSTYLE: resume JavadocVariable check

  861.         /** Create a model for specified mean orbit.
  862.          * @param mean mean orbit
  863.          * @param mass constant mass
  864.          * @param referenceRadius reference radius of the central body attraction model (m)
  865.          * @param mu central attraction coefficient (m³/s²)
  866.          * @param ck0 un-normalized zonal coefficients
  867.          */
  868.         BLModel(final KeplerianOrbit mean, final double mass,
  869.                 final double referenceRadius, final double mu, final double[] ck0) {

  870.             this.mean = mean;
  871.             this.mass = mass;

  872.             final double app = mean.getA();
  873.             xnotDot = FastMath.sqrt(mu / app) / app;

  874.             // preliminary processing
  875.             final double q = referenceRadius / app;
  876.             double ql = q * q;
  877.             final double y2 = -0.5 * ck0[2] * ql;

  878.             n = FastMath.sqrt(1 - mean.getE() * mean.getE());
  879.             final double n2 = n * n;
  880.             final double n3 = n2 * n;
  881.             final double n4 = n2 * n2;
  882.             final double n6 = n4 * n2;
  883.             final double n8 = n4 * n4;
  884.             final double n10 = n8 * n2;

  885.             final double yp2 = y2 / n4;
  886.             ql *= q;
  887.             final double yp3 = ck0[3] * ql / n6;
  888.             ql *= q;
  889.             final double yp4 = 0.375 * ck0[4] * ql / n8;
  890.             ql *= q;
  891.             final double yp5 = ck0[5] * ql / n10;


  892.             final SinCos sc    = FastMath.sinCos(mean.getI());
  893.             final double sinI1 = sc.sin();
  894.             final double sinI2 = sinI1 * sinI1;
  895.             final double cosI1 = sc.cos();
  896.             final double cosI2 = cosI1 * cosI1;
  897.             final double cosI3 = cosI2 * cosI1;
  898.             final double cosI4 = cosI2 * cosI2;
  899.             final double cosI6 = cosI4 * cosI2;
  900.             final double C5c2 = 1.0 / T2(cosI1);
  901.             final double C3c2 = 3.0 * cosI2 - 1.0;

  902.             final double epp = mean.getE();
  903.             final double epp2 = epp * epp;
  904.             final double epp3 = epp2 * epp;
  905.             final double epp4 = epp2 * epp2;

  906.             if (epp >= 1) {
  907.                 // Only for elliptical (e < 1) orbits
  908.                 throw new OrekitException(OrekitMessages.TOO_LARGE_ECCENTRICITY_FOR_PROPAGATION_MODEL,
  909.                                           mean.getE());
  910.             }

  911.             // secular multiplicative
  912.             lt = 1 +
  913.                     1.5 * yp2 * n * C3c2 +
  914.                     0.09375 * yp2 * yp2 * n * (-15.0 + 16.0 * n + 25.0 * n2 + (30.0 - 96.0 * n - 90.0 * n2) * cosI2 + (105.0 + 144.0 * n + 25.0 * n2) * cosI4) +
  915.                     0.9375 * yp4 * n * epp2 * (3.0 - 30.0 * cosI2 + 35.0 * cosI4);

  916.             gt = -1.5 * yp2 * C5c2 +
  917.                     0.09375 * yp2 * yp2 * (-35.0 + 24.0 * n + 25.0 * n2 + (90.0 - 192.0 * n - 126.0 * n2) * cosI2 + (385.0 + 360.0 * n + 45.0 * n2) * cosI4) +
  918.                     0.3125 * yp4 * (21.0 - 9.0 * n2 + (-270.0 + 126.0 * n2) * cosI2 + (385.0 - 189.0 * n2) * cosI4 );

  919.             ht = -3.0 * yp2 * cosI1 +
  920.                     0.375 * yp2 * yp2 * ((-5.0 + 12.0 * n + 9.0 * n2) * cosI1 + (-35.0 - 36.0 * n - 5.0 * n2) * cosI3) +
  921.                     1.25 * yp4 * (5.0 - 3.0 * n2) * cosI1 * (3.0 - 7.0 * cosI2);

  922.             final double cA = 1.0 - 11.0 * cosI2 - 40.0 * cosI4 / C5c2;
  923.             final double cB = 1.0 - 3.0 * cosI2 - 8.0 * cosI4 / C5c2;
  924.             final double cC = 1.0 - 9.0 * cosI2 - 24.0 * cosI4 / C5c2;
  925.             final double cD = 1.0 - 5.0 * cosI2 - 16.0 * cosI4 / C5c2;

  926.             final double qyp2_4 = 3.0 * yp2 * yp2 * cA -
  927.                                   10.0 * yp4 * cB;
  928.             final double qyp52 = epp3 * cosI1 * (0.5 * cD / sinI1 +
  929.                                  sinI1 * (5.0 + 32.0 * cosI2 / C5c2 + 80.0 * cosI4 / C5c2 / C5c2));
  930.             final double qyp22 = 2.0 + epp2 - 11.0 * (2.0 + 3.0 * epp2) * cosI2 -
  931.                                  40.0 * (2.0 + 5.0 * epp2) * cosI4 / C5c2 -
  932.                                  400.0 * epp2 * cosI6 / C5c2 / C5c2;
  933.             final double qyp42 = ( qyp22 + 4.0 * (2.0 + epp2 - (2.0 + 3.0 * epp2) * cosI2) ) / 5.0;
  934.             final double qyp52bis = epp * cosI1 * sinI1 *
  935.                                     (4.0 + 3.0 * epp2) *
  936.                                     (3.0 + 16.0 * cosI2 / C5c2 + 40.0 * cosI4 / C5c2 / C5c2);

  937.             // long periodic multiplicative
  938.             dei3sg =  35.0 / 96.0 * yp5 / yp2 * epp2 * n2 * cD * sinI1;
  939.             de2sg = -1.0 / 12.0 * epp * n2 / yp2 * qyp2_4;
  940.             deisg = ( -35.0 / 128.0 * yp5 / yp2 * epp2 * n2 * cD +
  941.                     1.0 / 4.0 * n2 / yp2 * (yp3 + 5.0 / 16.0 * yp5 * (4.0 + 3.0 * epp2) * cC)) * sinI1;
  942.             de = epp2 * n2 / 24.0 / yp2 * qyp2_4;

  943.             final double qyp52quotient = epp * (-32.0 + 81.0 * epp4) / (4.0 + 3.0 * epp2 + n * (4.0 + 9.0 * epp2));
  944.             dlgs2g = 1.0 / 48.0 / yp2 * (-3.0 * yp2 * yp2 * qyp22 +
  945.                      10.0 * yp4 * qyp42 ) +
  946.                      n3 / yp2 * qyp2_4 / 24.0;
  947.             dlgc3g = 35.0 / 384.0 * yp5 / yp2 * n3 * epp * cD * sinI1 +
  948.                      35.0 / 1152.0 * yp5 / yp2 * (2.0 * qyp52 * cosI1 - epp * cD * sinI1 * (3.0 + 2.0 * epp2));
  949.             dlgcg = -yp3 * epp * cosI2 / ( 4.0 * yp2 * sinI1) +
  950.                     0.078125 * yp5 / yp2 * (-epp * cosI2 / sinI1 * (4.0 + 3.0 * epp2) + epp2 * sinI1 * (26.0 + 9.0 * epp2)) * cC -
  951.                     0.46875 * yp5 / yp2 * qyp52bis * cosI1 +
  952.                     0.25 * yp3 / yp2 * sinI1 * epp / (1.0 + n3) * (3.0 - epp2 * (3.0 - epp2)) +
  953.                     0.078125 * yp5 / yp2 * n2 * cC * qyp52quotient * sinI1;


  954.             final double qyp24 = 3.0 * yp2 * yp2 * (11.0 + 80.0 * cosI2 / sinI1 + 200.0 * cosI4 / sinI2) -
  955.                                  10.0 * yp4 * (3.0 + 16.0 * cosI2 / sinI1 + 40.0 * cosI4 / sinI2);
  956.             dh2sgcg = 35.0 / 144.0 * yp5 / yp2 * qyp52;
  957.             dhsgcg = -epp2 * cosI1 / (12.0 * yp2) * qyp24;
  958.             dhcg = -35.0 / 576.0 * yp5 / yp2 * qyp52 +
  959.                    epp * cosI1 / (4.0 * yp2 * sinI1) * (yp3 + 0.3125 * yp5 * (4.0 + 3.0 * epp2) * cC) +
  960.                    1.875 / (4.0 * yp2) * yp5 * qyp52bis;

  961.             // short periodic multiplicative
  962.             aC = -yp2 * C3c2 * app / n3;
  963.             aCbis = y2 * app * C3c2;
  964.             ac2g2f = y2 * app * 3.0 * sinI2;

  965.             double qe = 0.5 * n2 * y2 * C3c2 / n6;
  966.             eC = qe * epp / (1.0 + n3) * (3.0 - epp2 * (3.0 - epp2));
  967.             ecf = 3.0 * qe;
  968.             e2cf = 3.0 * epp * qe;
  969.             e3cf = epp2 * qe;
  970.             qe = 0.5 * n2 * y2 * 3.0 * (1.0 - cosI2) / n6;
  971.             ec2f2g = qe * epp;
  972.             ecfc2f2g = 3.0 * qe;
  973.             e2cfc2f2g = 3.0 * epp * qe;
  974.             e3cfc2f2g = epp2 * qe;
  975.             qe = -0.5 * yp2 * n2 * (1.0 - cosI2);
  976.             ec2gf = 3.0 * qe;
  977.             ec2g3f = qe;

  978.             double qi = epp * yp2 * cosI1 * sinI1;
  979.             ide = -epp * cosI1 / (n2 * sinI1);
  980.             isfs2f2g = qi;
  981.             icfc2f2g = 2.0 * qi;
  982.             qi = yp2 * cosI1 * sinI1;
  983.             ic2f2g = 1.5 * qi;

  984.             double qgl1 = 0.25 * yp2;
  985.             double qgl2 = 0.25 * yp2 * epp * n2 / (1.0 + n);
  986.             glf = qgl1 * -6.0 * C5c2;
  987.             gll = qgl1 * 6.0 * C5c2;
  988.             glsf = qgl1 * -6.0 * C5c2 * epp +
  989.                    qgl2 * 2.0 * C3c2;
  990.             glosf = qgl2 * 2.0 * C3c2;
  991.             qgl1 = qgl1 * (3.0 - 5.0 * cosI2);
  992.             qgl2 = qgl2 * 3.0 * (1.0 - cosI2);
  993.             gls2f2g = 3.0 * qgl1;
  994.             gls2gf = 3.0 * epp * qgl1 +
  995.                      qgl2;
  996.             glos2gf = -1.0 * qgl2;
  997.             gls2g3f = qgl1 * epp +
  998.                       1.0 / 3.0 * qgl2;
  999.             glos2g3f = qgl2;

  1000.             final double qh = 3.0 * yp2 * cosI1;
  1001.             hf = -qh;
  1002.             hl = qh;
  1003.             hsf = -epp * qh;
  1004.             hcfs2g2f = 2.0 * epp * yp2 * cosI1;
  1005.             hs2g2f = 1.5 * yp2 * cosI1;
  1006.             hsfc2g2f = -epp * yp2 * cosI1;

  1007.             final double qedl = -0.25 * yp2 * n3;
  1008.             edls2g = 1.0 / 24.0 * epp * n3 / yp2 * qyp2_4;
  1009.             edlcg = -0.25 * yp3 / yp2 * n3 * sinI1 -
  1010.                     0.078125 * yp5 / yp2 * n3 * sinI1 * (4.0 + 9.0 * epp2) * cC;
  1011.             edlc3g = 35.0 / 384.0 * yp5 / yp2 * n3 * epp2 * cD * sinI1;
  1012.             edlsf = 2.0 * qedl * C3c2;
  1013.             edls2gf = 3.0 * qedl * (1.0 - cosI2);
  1014.             edls2g3f = 1.0 / 3.0 * qedl;

  1015.             // secular rates of the mean semi-major axis and eccentricity
  1016.             // Eq. 2.41 and Eq. 2.45 of Phipps' 1992 thesis
  1017.             aRate = -4.0 * app / (3.0 * xnotDot);
  1018.             eRate = -4.0 * epp * n * n / (3.0 * xnotDot);

  1019.         }

  1020.         /**
  1021.          * Accurate computation of E - e sin(E).
  1022.          *
  1023.          * @param E eccentric anomaly
  1024.          * @return E - e sin(E)
  1025.          */
  1026.         private UnivariateDerivative2 eMeSinE(final UnivariateDerivative2 E) {
  1027.             UnivariateDerivative2 x = E.sin().multiply(1 - mean.getE());
  1028.             final UnivariateDerivative2 mE2 = E.negate().multiply(E);
  1029.             UnivariateDerivative2 term = E;
  1030.             UnivariateDerivative2 d    = E.getField().getZero();
  1031.             // the inequality test below IS intentional and should NOT be replaced by a check with a small tolerance
  1032.             for (UnivariateDerivative2 x0 = d.add(Double.NaN); !Double.valueOf(x.getValue()).equals(Double.valueOf(x0.getValue()));) {
  1033.                 d = d.add(2);
  1034.                 term = term.multiply(mE2.divide(d.multiply(d.add(1))));
  1035.                 x0 = x;
  1036.                 x = x.subtract(term);
  1037.             }
  1038.             return x;
  1039.         }

  1040.         /**
  1041.          * Gets eccentric anomaly from mean anomaly.
  1042.          * <p>The algorithm used to solve the Kepler equation has been published in:
  1043.          * "Procedures for  solving Kepler's Equation", A. W. Odell and R. H. Gooding,
  1044.          * Celestial Mechanics 38 (1986) 307-334</p>
  1045.          * <p>It has been copied from the OREKIT library (KeplerianOrbit class).</p>
  1046.          *
  1047.          * @param mk the mean anomaly (rad)
  1048.          * @return the eccentric anomaly (rad)
  1049.          */
  1050.         private UnivariateDerivative2 getEccentricAnomaly(final UnivariateDerivative2 mk) {
  1051.             final double k1 = 3 * FastMath.PI + 2;
  1052.             final double k2 = FastMath.PI - 1;
  1053.             final double k3 = 6 * FastMath.PI - 1;
  1054.             final double A  = 3.0 * k2 * k2 / k1;
  1055.             final double B  = k3 * k3 / (6.0 * k1);
  1056.             // reduce M to [-PI PI] interval
  1057.             final UnivariateDerivative2 reducedM = new UnivariateDerivative2(MathUtils.normalizeAngle(mk.getValue(), 0.0),
  1058.                                                                              mk.getFirstDerivative(),
  1059.                                                                              mk.getSecondDerivative());

  1060.             // compute start value according to A. W. Odell and R. H. Gooding S12 starter
  1061.             UnivariateDerivative2 ek;
  1062.             if (FastMath.abs(reducedM.getValue()) < 1.0 / 6.0) {
  1063.                 if (FastMath.abs(reducedM.getValue()) < Precision.SAFE_MIN) {
  1064.                     // this is an Orekit change to the S12 starter.
  1065.                     // If reducedM is 0.0, the derivative of cbrt is infinite which induces NaN appearing later in
  1066.                     // the computation. As in this case E and M are almost equal, we initialize ek with reducedM
  1067.                     ek = reducedM;
  1068.                 } else {
  1069.                     // this is the standard S12 starter
  1070.                     ek = reducedM.add(reducedM.multiply(6).cbrt().subtract(reducedM).multiply(mean.getE()));
  1071.                 }
  1072.             } else {
  1073.                 if (reducedM.getValue() < 0) {
  1074.                     final UnivariateDerivative2 w = reducedM.add(FastMath.PI);
  1075.                     ek = reducedM.add(w.multiply(-A).divide(w.subtract(B)).subtract(FastMath.PI).subtract(reducedM).multiply(mean.getE()));
  1076.                 } else {
  1077.                     final UnivariateDerivative2 minusW = reducedM.subtract(FastMath.PI);
  1078.                     ek = reducedM.add(minusW.multiply(A).divide(minusW.add(B)).add(FastMath.PI).subtract(reducedM).multiply(mean.getE()));
  1079.                 }
  1080.             }

  1081.             final double e1 = 1 - mean.getE();
  1082.             final boolean noCancellationRisk = (e1 + ek.getValue() * ek.getValue() / 6) >= 0.1;

  1083.             // perform two iterations, each consisting of one Halley step and one Newton-Raphson step
  1084.             for (int j = 0; j < 2; ++j) {
  1085.                 final UnivariateDerivative2 f;
  1086.                 UnivariateDerivative2 fd;
  1087.                 final UnivariateDerivative2 fdd  = ek.sin().multiply(mean.getE());
  1088.                 final UnivariateDerivative2 fddd = ek.cos().multiply(mean.getE());
  1089.                 if (noCancellationRisk) {
  1090.                     f  = ek.subtract(fdd).subtract(reducedM);
  1091.                     fd = fddd.subtract(1).negate();
  1092.                 } else {
  1093.                     f  = eMeSinE(ek).subtract(reducedM);
  1094.                     final UnivariateDerivative2 s = ek.multiply(0.5).sin();
  1095.                     fd = s.multiply(s).multiply(2 * mean.getE()).add(e1);
  1096.                 }
  1097.                 final UnivariateDerivative2 dee = f.multiply(fd).divide(f.multiply(0.5).multiply(fdd).subtract(fd.multiply(fd)));

  1098.                 // update eccentric anomaly, using expressions that limit underflow problems
  1099.                 final UnivariateDerivative2 w = fd.add(dee.multiply(0.5).multiply(fdd.add(dee.multiply(fdd).divide(3))));
  1100.                 fd = fd.add(dee.multiply(fdd.add(dee.multiply(0.5).multiply(fdd))));
  1101.                 ek = ek.subtract(f.subtract(dee.multiply(fd.subtract(w))).divide(fd));
  1102.             }

  1103.             // expand the result back to original range
  1104.             ek = ek.add(mk.getValue() - reducedM.getValue());

  1105.             // Returns the eccentric anomaly
  1106.             return ek;
  1107.         }

  1108.         /**
  1109.          * This method is used in Brouwer-Lyddane model to avoid singularity at the
  1110.          * critical inclination (i = 63.4°).
  1111.          * <p>
  1112.          * This method, based on Warren Phipps's 1992 thesis (Eq. 2.47 and 2.48),
  1113.          * approximate the factor (1.0 - 5.0 * cos²(inc))^-1 (causing the singularity)
  1114.          * by a function, named T2 in the thesis.
  1115.          * </p>
  1116.          * @param cosInc cosine of the mean inclination
  1117.          * @return an approximation of (1.0 - 5.0 * cos²(inc))^-1 term
  1118.          */
  1119.         private double T2(final double cosInc) {

  1120.             // X = (1.0 - 5.0 * cos²(inc))
  1121.             final double x  = 1.0 - 5.0 * cosInc * cosInc;
  1122.             final double x2 = x * x;

  1123.             // Eq. 2.48
  1124.             double sum = 0.0;
  1125.             for (int i = 0; i <= 12; i++) {
  1126.                 final double sign = i % 2 == 0 ? +1.0 : -1.0;
  1127.                 sum += sign * FastMath.pow(BETA, i) * FastMath.pow(x2, i) / CombinatoricsUtils.factorialDouble(i + 1);
  1128.             }

  1129.             // Right term of equation 2.47
  1130.             double product = 1.0;
  1131.             for (int i = 0; i <= 10; i++) {
  1132.                 product *= 1 + FastMath.exp(FastMath.scalb(-1.0, i) * BETA * x2);
  1133.             }

  1134.             // Return (Eq. 2.47)
  1135.             return BETA * x * sum * product;

  1136.         }

  1137.         /** Extrapolate an orbit up to a specific target date.
  1138.          * @param date target date for the orbit
  1139.          * @return propagated parameters
  1140.          */
  1141.         public KeplerianOrbit propagateParameters(final AbsoluteDate date) {

  1142.             // Empirical drag coefficient M2
  1143.             final double m2 = M2Driver.getValue();

  1144.             // Keplerian evolution
  1145.             final UnivariateDerivative2 dt = new UnivariateDerivative2(date.durationFrom(mean.getDate()), 1.0, 0.0);
  1146.             final UnivariateDerivative2 xnot = dt.multiply(xnotDot);

  1147.             //____________________________________
  1148.             // secular effects

  1149.             // mean mean anomaly (with drag Eq. 2.38 of Phipps' 1992 thesis)
  1150.             final UnivariateDerivative2 dtM2  = dt.multiply(m2);
  1151.             final UnivariateDerivative2 dt2M2 = dt.multiply(dtM2);
  1152.             final UnivariateDerivative2 lpp = new UnivariateDerivative2(MathUtils.normalizeAngle(mean.getMeanAnomaly() + lt * xnot.getValue() + dt2M2.getValue(), 0),
  1153.                                                                       lt * xnotDot + 2.0 * dtM2.getValue(),
  1154.                                                                       2.0 * m2);
  1155.             // mean argument of perigee
  1156.             final UnivariateDerivative2 gpp = new UnivariateDerivative2(MathUtils.normalizeAngle(mean.getPerigeeArgument() + gt * xnot.getValue(), 0),
  1157.                                                                       gt * xnotDot,
  1158.                                                                       0.0);
  1159.             // mean longitude of ascending node
  1160.             final UnivariateDerivative2 hpp = new UnivariateDerivative2(MathUtils.normalizeAngle(mean.getRightAscensionOfAscendingNode() + ht * xnot.getValue(), 0),
  1161.                                                                       ht * xnotDot,
  1162.                                                                       0.0);

  1163.             // ________________________________________________
  1164.             // secular rates of the mean semi-major axis and eccentricity

  1165.             // semi-major axis
  1166.             final UnivariateDerivative2 appDrag = dt.multiply(aRate * m2);

  1167.             // eccentricity
  1168.             final UnivariateDerivative2 eppDrag = dt.multiply(eRate * m2);

  1169.             //____________________________________
  1170.             // Long periodical terms
  1171.             final UnivariateDerivative2 cg1 = gpp.cos();
  1172.             final UnivariateDerivative2 sg1 = gpp.sin();
  1173.             final UnivariateDerivative2 c2g = cg1.multiply(cg1).subtract(sg1.multiply(sg1));
  1174.             final UnivariateDerivative2 s2g = cg1.multiply(sg1).add(sg1.multiply(cg1));
  1175.             final UnivariateDerivative2 c3g = c2g.multiply(cg1).subtract(s2g.multiply(sg1));
  1176.             final UnivariateDerivative2 sg2 = sg1.multiply(sg1);
  1177.             final UnivariateDerivative2 sg3 = sg1.multiply(sg2);


  1178.             // de eccentricity
  1179.             final UnivariateDerivative2 d1e = sg3.multiply(dei3sg).
  1180.                                               add(sg1.multiply(deisg)).
  1181.                                               add(sg2.multiply(de2sg)).
  1182.                                               add(de);

  1183.             // l' + g'
  1184.             final UnivariateDerivative2 lp_p_gp = s2g.multiply(dlgs2g).
  1185.                                                add(c3g.multiply(dlgc3g)).
  1186.                                                add(cg1.multiply(dlgcg)).
  1187.                                                add(lpp).
  1188.                                                add(gpp);

  1189.             // h'
  1190.             final UnivariateDerivative2 hp = sg2.multiply(cg1).multiply(dh2sgcg).
  1191.                                                add(sg1.multiply(cg1).multiply(dhsgcg)).
  1192.                                                add(cg1.multiply(dhcg)).
  1193.                                                add(hpp);

  1194.             // Short periodical terms
  1195.             // eccentric anomaly
  1196.             final UnivariateDerivative2 Ep = getEccentricAnomaly(lpp);
  1197.             final UnivariateDerivative2 cf1 = (Ep.cos().subtract(mean.getE())).divide(Ep.cos().multiply(-mean.getE()).add(1.0));
  1198.             final UnivariateDerivative2 sf1 = (Ep.sin().multiply(n)).divide(Ep.cos().multiply(-mean.getE()).add(1.0));
  1199.             final UnivariateDerivative2 f = FastMath.atan2(sf1, cf1);

  1200.             final UnivariateDerivative2 c2f = cf1.multiply(cf1).subtract(sf1.multiply(sf1));
  1201.             final UnivariateDerivative2 s2f = cf1.multiply(sf1).add(sf1.multiply(cf1));
  1202.             final UnivariateDerivative2 c3f = c2f.multiply(cf1).subtract(s2f.multiply(sf1));
  1203.             final UnivariateDerivative2 s3f = c2f.multiply(sf1).add(s2f.multiply(cf1));
  1204.             final UnivariateDerivative2 cf2 = cf1.multiply(cf1);
  1205.             final UnivariateDerivative2 cf3 = cf1.multiply(cf2);

  1206.             final UnivariateDerivative2 c2g1f = cf1.multiply(c2g).subtract(sf1.multiply(s2g));
  1207.             final UnivariateDerivative2 c2g2f = c2f.multiply(c2g).subtract(s2f.multiply(s2g));
  1208.             final UnivariateDerivative2 c2g3f = c3f.multiply(c2g).subtract(s3f.multiply(s2g));
  1209.             final UnivariateDerivative2 s2g1f = cf1.multiply(s2g).add(c2g.multiply(sf1));
  1210.             final UnivariateDerivative2 s2g2f = c2f.multiply(s2g).add(c2g.multiply(s2f));
  1211.             final UnivariateDerivative2 s2g3f = c3f.multiply(s2g).add(c2g.multiply(s3f));

  1212.             final UnivariateDerivative2 eE = (Ep.cos().multiply(-mean.getE()).add(1.0)).reciprocal();
  1213.             final UnivariateDerivative2 eE3 = eE.multiply(eE).multiply(eE);
  1214.             final UnivariateDerivative2 sigma = eE.multiply(n * n).multiply(eE).add(eE);

  1215.             // Semi-major axis (with drag Eq. 2.41 of Phipps' 1992 thesis)
  1216.             final UnivariateDerivative2 a = eE3.multiply(aCbis).add(appDrag.add(mean.getA())).
  1217.                                             add(aC).
  1218.                                             add(eE3.multiply(c2g2f).multiply(ac2g2f));

  1219.             // Eccentricity (with drag Eq. 2.45 of Phipps' 1992 thesis)
  1220.             final UnivariateDerivative2 e = d1e.add(eppDrag.add(mean.getE())).
  1221.                                             add(eC).
  1222.                                             add(cf1.multiply(ecf)).
  1223.                                             add(cf2.multiply(e2cf)).
  1224.                                             add(cf3.multiply(e3cf)).
  1225.                                             add(c2g2f.multiply(ec2f2g)).
  1226.                                             add(c2g2f.multiply(cf1).multiply(ecfc2f2g)).
  1227.                                             add(c2g2f.multiply(cf2).multiply(e2cfc2f2g)).
  1228.                                             add(c2g2f.multiply(cf3).multiply(e3cfc2f2g)).
  1229.                                             add(c2g1f.multiply(ec2gf)).
  1230.                                             add(c2g3f.multiply(ec2g3f));

  1231.             // Inclination
  1232.             final UnivariateDerivative2 i = d1e.multiply(ide).
  1233.                                             add(mean.getI()).
  1234.                                             add(sf1.multiply(s2g2f.multiply(isfs2f2g))).
  1235.                                             add(cf1.multiply(c2g2f.multiply(icfc2f2g))).
  1236.                                             add(c2g2f.multiply(ic2f2g));

  1237.             // Argument of perigee + True anomaly
  1238.             final UnivariateDerivative2 g_p_l = lp_p_gp.add(f.multiply(glf)).
  1239.                                                 add(lpp.multiply(gll)).
  1240.                                                 add(sf1.multiply(glsf)).
  1241.                                                 add(sigma.multiply(sf1).multiply(glosf)).
  1242.                                                 add(s2g2f.multiply(gls2f2g)).
  1243.                                                 add(s2g1f.multiply(gls2gf)).
  1244.                                                 add(sigma.multiply(s2g1f).multiply(glos2gf)).
  1245.                                                 add(s2g3f.multiply(gls2g3f)).
  1246.                                                 add(sigma.multiply(s2g3f).multiply(glos2g3f));


  1247.             // Longitude of ascending node
  1248.             final UnivariateDerivative2 h = hp.add(f.multiply(hf)).
  1249.                                             add(lpp.multiply(hl)).
  1250.                                             add(sf1.multiply(hsf)).
  1251.                                             add(cf1.multiply(s2g2f).multiply(hcfs2g2f)).
  1252.                                             add(s2g2f.multiply(hs2g2f)).
  1253.                                             add(c2g2f.multiply(sf1).multiply(hsfc2g2f));

  1254.             final UnivariateDerivative2 edl = s2g.multiply(edls2g).
  1255.                                             add(cg1.multiply(edlcg)).
  1256.                                             add(c3g.multiply(edlc3g)).
  1257.                                             add(sf1.multiply(edlsf)).
  1258.                                             add(s2g1f.multiply(edls2gf)).
  1259.                                             add(s2g3f.multiply(edls2g3f)).
  1260.                                             add(sf1.multiply(sigma).multiply(edlsf)).
  1261.                                             add(s2g1f.multiply(sigma).multiply(-edls2gf)).
  1262.                                             add(s2g3f.multiply(sigma).multiply(3.0 * edls2g3f));

  1263.             final UnivariateDerivative2 A = e.multiply(lpp.cos()).subtract(edl.multiply(lpp.sin()));
  1264.             final UnivariateDerivative2 B = e.multiply(lpp.sin()).add(edl.multiply(lpp.cos()));

  1265.             // True anomaly
  1266.             final UnivariateDerivative2 l = FastMath.atan2(B, A);

  1267.             // Argument of perigee
  1268.             final UnivariateDerivative2 g = g_p_l.subtract(l);

  1269.             // Return a Keplerian orbit
  1270.             return new KeplerianOrbit(a.getValue(), e.getValue(), i.getValue(),
  1271.                                       g.getValue(), h.getValue(), l.getValue(),
  1272.                                       a.getFirstDerivative(), e.getFirstDerivative(), i.getFirstDerivative(),
  1273.                                       g.getFirstDerivative(), h.getFirstDerivative(), l.getFirstDerivative(),
  1274.                                       PositionAngle.MEAN, mean.getFrame(), date, mu);

  1275.         }

  1276.     }

  1277.     /** {@inheritDoc} */
  1278.     protected double getMass(final AbsoluteDate date) {
  1279.         return models.get(date).mass;
  1280.     }

  1281. }