EcksteinHechlerPropagator.java

  1. /* Copyright 2002-2024 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.io.Serializable;

  19. import org.hipparchus.analysis.differentiation.UnivariateDerivative2;
  20. import org.hipparchus.geometry.euclidean.threed.FieldVector3D;
  21. import org.hipparchus.geometry.euclidean.threed.Vector3D;
  22. import org.hipparchus.linear.RealMatrix;
  23. import org.hipparchus.util.FastMath;
  24. import org.hipparchus.util.MathUtils;
  25. import org.hipparchus.util.SinCos;
  26. import org.orekit.attitudes.AttitudeProvider;
  27. import org.orekit.attitudes.FrameAlignedProvider;
  28. import org.orekit.errors.OrekitException;
  29. import org.orekit.errors.OrekitMessages;
  30. import org.orekit.forces.gravity.potential.UnnormalizedSphericalHarmonicsProvider;
  31. import org.orekit.forces.gravity.potential.UnnormalizedSphericalHarmonicsProvider.UnnormalizedSphericalHarmonics;
  32. import org.orekit.orbits.CartesianOrbit;
  33. import org.orekit.orbits.CircularOrbit;
  34. import org.orekit.orbits.Orbit;
  35. import org.orekit.orbits.OrbitType;
  36. import org.orekit.orbits.PositionAngleType;
  37. import org.orekit.propagation.AbstractMatricesHarvester;
  38. import org.orekit.propagation.PropagationType;
  39. import org.orekit.propagation.SpacecraftState;
  40. import org.orekit.time.AbsoluteDate;
  41. import org.orekit.utils.DoubleArrayDictionary;
  42. import org.orekit.utils.TimeSpanMap;
  43. import org.orekit.utils.TimeStampedPVCoordinates;

  44. /** This class propagates a {@link org.orekit.propagation.SpacecraftState}
  45.  *  using the analytical Eckstein-Hechler model.
  46.  * <p>The Eckstein-Hechler model is suited for near circular orbits
  47.  * (e &lt; 0.1, with poor accuracy between 0.005 and 0.1) and inclination
  48.  * neither equatorial (direct or retrograde) nor critical (direct or
  49.  * retrograde).</p>
  50.  * <p>
  51.  * Note that before version 7.0, there was a large inconsistency in the generated
  52.  * orbits, and it was fixed as of version 7.0 of Orekit, with a visible side effect.
  53.  * The problems is that if the circular parameters produced by the Eckstein-Hechler
  54.  * model are used to build an orbit considered to be osculating, the velocity deduced
  55.  * from this orbit was <em>inconsistent with the position evolution</em>! The reason is
  56.  * that the model includes non-Keplerian effects but it does not include a corresponding
  57.  * circular/Cartesian conversion. As a consequence, all subsequent computation involving
  58.  * velocity were wrong. This includes attitude modes like yaw compensation and Doppler
  59.  * effect. As this effect was considered serious enough and as accurate velocities were
  60.  * considered important, the propagator now generates {@link CartesianOrbit Cartesian
  61.  * orbits} which are built in a special way to ensure consistency throughout propagation.
  62.  * A side effect is that if circular parameters are rebuilt by user from these propagated
  63.  * Cartesian orbit, the circular parameters will generally <em>not</em> match the initial
  64.  * orbit (differences in semi-major axis can exceed 120 m). The position however <em>will</em>
  65.  * match to sub-micrometer level, and this position will be identical to the positions
  66.  * that were generated by previous versions (in other words, the internals of the models
  67.  * have not been changed, only the output parameters have been changed). The correctness
  68.  * of the initialization has been assessed and is good, as it allows the subsequent orbit
  69.  * to remain close to a numerical reference orbit.
  70.  * </p>
  71.  * <p>
  72.  * If users need a more definitive initialization of an Eckstein-Hechler propagator, they
  73.  * should consider using a {@link org.orekit.propagation.conversion.PropagatorConverter
  74.  * propagator converter} to initialize their Eckstein-Hechler propagator using a complete
  75.  * sample instead of just a single initial orbit.
  76.  * </p>
  77.  * @see Orbit
  78.  * @author Guylaine Prat
  79.  */
  80. public class EcksteinHechlerPropagator extends AbstractAnalyticalPropagator {

  81.     /** Default convergence threshold for mean parameters conversion. */
  82.     private static final double EPSILON_DEFAULT = 1.0e-11;

  83.     /** Default value for maxIterations. */
  84.     private static final int MAX_ITERATIONS_DEFAULT = 100;

  85.     /** Initial Eckstein-Hechler model. */
  86.     private EHModel initialModel;

  87.     /** All models. */
  88.     private transient TimeSpanMap<EHModel> models;

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

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

  93.     /** Un-normalized zonal coefficients. */
  94.     private double[] ck0;

  95.     /** Build a propagator from orbit and potential provider.
  96.      * <p>Mass and attitude provider are set to unspecified non-null arbitrary values.</p>
  97.      *
  98.      * <p>Using this constructor, an initial osculating orbit is considered.</p>
  99.      *
  100.      * @param initialOrbit initial orbit
  101.      * @param provider for un-normalized zonal coefficients
  102.      * @see #EcksteinHechlerPropagator(Orbit, AttitudeProvider,
  103.      * UnnormalizedSphericalHarmonicsProvider)
  104.      * @see #EcksteinHechlerPropagator(Orbit, UnnormalizedSphericalHarmonicsProvider,
  105.      *                                 PropagationType)
  106.      */
  107.     public EcksteinHechlerPropagator(final Orbit initialOrbit,
  108.                                      final UnnormalizedSphericalHarmonicsProvider provider) {
  109.         this(initialOrbit, FrameAlignedProvider.of(initialOrbit.getFrame()),
  110.              DEFAULT_MASS, provider, provider.onDate(initialOrbit.getDate()));
  111.     }

  112.     /**
  113.      * Private helper constructor.
  114.      * <p>Using this constructor, an initial osculating orbit is considered.</p>
  115.      * @param initialOrbit initial orbit
  116.      * @param attitude attitude provider
  117.      * @param mass spacecraft mass
  118.      * @param provider for un-normalized zonal coefficients
  119.      * @param harmonics {@code provider.onDate(initialOrbit.getDate())}
  120.      * @see #EcksteinHechlerPropagator(Orbit, AttitudeProvider, double,
  121.      *                                 UnnormalizedSphericalHarmonicsProvider,
  122.      *                                 UnnormalizedSphericalHarmonicsProvider.UnnormalizedSphericalHarmonics,
  123.      *                                 PropagationType)
  124.      */
  125.     public EcksteinHechlerPropagator(final Orbit initialOrbit,
  126.                                      final AttitudeProvider attitude,
  127.                                      final double mass,
  128.                                      final UnnormalizedSphericalHarmonicsProvider provider,
  129.                                      final UnnormalizedSphericalHarmonics harmonics) {
  130.         this(initialOrbit, attitude, mass, provider.getAe(), provider.getMu(),
  131.              harmonics.getUnnormalizedCnm(2, 0),
  132.              harmonics.getUnnormalizedCnm(3, 0),
  133.              harmonics.getUnnormalizedCnm(4, 0),
  134.              harmonics.getUnnormalizedCnm(5, 0),
  135.              harmonics.getUnnormalizedCnm(6, 0));
  136.     }

  137.     /** Build a propagator from orbit and potential.
  138.      * <p>Mass and attitude provider are set to unspecified non-null arbitrary values.</p>
  139.      * <p>The C<sub>n,0</sub> coefficients are the denormalized zonal coefficients, they
  140.      * are related to both the normalized coefficients
  141.      * <span style="text-decoration: overline">C</span><sub>n,0</sub>
  142.      *  and the J<sub>n</sub> one as follows:</p>
  143.      *
  144.      * <p> C<sub>n,0</sub> = [(2-δ<sub>0,m</sub>)(2n+1)(n-m)!/(n+m)!]<sup>½</sup>
  145.      * <span style="text-decoration: overline">C</span><sub>n,0</sub>
  146.      *
  147.      * <p> C<sub>n,0</sub> = -J<sub>n</sub>
  148.      *
  149.      * <p>Using this constructor, an initial osculating orbit is considered.</p>
  150.      *
  151.      * @param initialOrbit initial orbit
  152.      * @param referenceRadius reference radius of the Earth for the potential model (m)
  153.      * @param mu central attraction coefficient (m³/s²)
  154.      * @param c20 un-normalized zonal coefficient (about -1.08e-3 for Earth)
  155.      * @param c30 un-normalized zonal coefficient (about +2.53e-6 for Earth)
  156.      * @param c40 un-normalized zonal coefficient (about +1.62e-6 for Earth)
  157.      * @param c50 un-normalized zonal coefficient (about +2.28e-7 for Earth)
  158.      * @param c60 un-normalized zonal coefficient (about -5.41e-7 for Earth)
  159.      * @see org.orekit.utils.Constants
  160.      * @see #EcksteinHechlerPropagator(Orbit, AttitudeProvider, double, double, double,
  161.      * double, double, double, double, double)
  162.      */
  163.     public EcksteinHechlerPropagator(final Orbit initialOrbit,
  164.                                      final double referenceRadius, final double mu,
  165.                                      final double c20, final double c30, final double c40,
  166.                                      final double c50, final double c60) {
  167.         this(initialOrbit, FrameAlignedProvider.of(initialOrbit.getFrame()),
  168.              DEFAULT_MASS, referenceRadius, mu, c20, c30, c40, c50, c60);
  169.     }

  170.     /** Build a propagator from orbit, mass and potential provider.
  171.      * <p>Attitude law is set to an unspecified non-null arbitrary value.</p>
  172.      *
  173.      * <p>Using this constructor, an initial osculating orbit is considered.</p>
  174.      *
  175.      * @param initialOrbit initial orbit
  176.      * @param mass spacecraft mass
  177.      * @param provider for un-normalized zonal coefficients
  178.      * @see #EcksteinHechlerPropagator(Orbit, AttitudeProvider, double,
  179.      * UnnormalizedSphericalHarmonicsProvider)
  180.      */
  181.     public EcksteinHechlerPropagator(final Orbit initialOrbit, final double mass,
  182.                                      final UnnormalizedSphericalHarmonicsProvider provider) {
  183.         this(initialOrbit, FrameAlignedProvider.of(initialOrbit.getFrame()),
  184.              mass, provider, provider.onDate(initialOrbit.getDate()));
  185.     }

  186.     /** Build a propagator from orbit, mass and potential.
  187.      * <p>Attitude law is set to an unspecified non-null arbitrary value.</p>
  188.      * <p>The C<sub>n,0</sub> coefficients are the denormalized zonal coefficients, they
  189.      * are related to both the normalized coefficients
  190.      * <span style="text-decoration: overline">C</span><sub>n,0</sub>
  191.      *  and the J<sub>n</sub> one as follows:</p>
  192.      *
  193.      * <p> C<sub>n,0</sub> = [(2-δ<sub>0,m</sub>)(2n+1)(n-m)!/(n+m)!]<sup>½</sup>
  194.      * <span style="text-decoration: overline">C</span><sub>n,0</sub>
  195.      *
  196.      * <p> C<sub>n,0</sub> = -J<sub>n</sub>
  197.      *
  198.      * <p>Using this constructor, an initial osculating orbit is considered.</p>
  199.      *
  200.      * @param initialOrbit initial orbit
  201.      * @param mass spacecraft mass
  202.      * @param referenceRadius reference radius of the Earth for the potential model (m)
  203.      * @param mu central attraction coefficient (m³/s²)
  204.      * @param c20 un-normalized zonal coefficient (about -1.08e-3 for Earth)
  205.      * @param c30 un-normalized zonal coefficient (about +2.53e-6 for Earth)
  206.      * @param c40 un-normalized zonal coefficient (about +1.62e-6 for Earth)
  207.      * @param c50 un-normalized zonal coefficient (about +2.28e-7 for Earth)
  208.      * @param c60 un-normalized zonal coefficient (about -5.41e-7 for Earth)
  209.      * @see #EcksteinHechlerPropagator(Orbit, AttitudeProvider, double, double, double,
  210.      * double, double, double, double, double)
  211.      */
  212.     public EcksteinHechlerPropagator(final Orbit initialOrbit, final double mass,
  213.                                      final double referenceRadius, final double mu,
  214.                                      final double c20, final double c30, final double c40,
  215.                                      final double c50, final double c60) {
  216.         this(initialOrbit, FrameAlignedProvider.of(initialOrbit.getFrame()),
  217.              mass, referenceRadius, mu, c20, c30, c40, c50, c60);
  218.     }

  219.     /** Build a propagator from orbit, attitude provider and potential provider.
  220.      * <p>Mass is set to an unspecified non-null arbitrary value.</p>
  221.      * <p>Using this constructor, an initial osculating orbit is considered.</p>
  222.      * @param initialOrbit initial orbit
  223.      * @param attitudeProv attitude provider
  224.      * @param provider for un-normalized zonal coefficients
  225.      */
  226.     public EcksteinHechlerPropagator(final Orbit initialOrbit,
  227.                                      final AttitudeProvider attitudeProv,
  228.                                      final UnnormalizedSphericalHarmonicsProvider provider) {
  229.         this(initialOrbit, attitudeProv, DEFAULT_MASS, provider, provider.onDate(initialOrbit.getDate()));
  230.     }

  231.     /** Build a propagator from orbit, attitude provider and potential.
  232.      * <p>Mass is set to an unspecified non-null arbitrary value.</p>
  233.      * <p>The C<sub>n,0</sub> coefficients are the denormalized zonal coefficients, they
  234.      * are related to both the normalized coefficients
  235.      * <span style="text-decoration: overline">C</span><sub>n,0</sub>
  236.      *  and the J<sub>n</sub> one as follows:</p>
  237.      *
  238.      * <p> C<sub>n,0</sub> = [(2-δ<sub>0,m</sub>)(2n+1)(n-m)!/(n+m)!]<sup>½</sup>
  239.      * <span style="text-decoration: overline">C</span><sub>n,0</sub>
  240.      *
  241.      * <p> C<sub>n,0</sub> = -J<sub>n</sub>
  242.      *
  243.      * <p>Using this constructor, an initial osculating orbit is considered.</p>
  244.      *
  245.      * @param initialOrbit initial orbit
  246.      * @param attitudeProv attitude provider
  247.      * @param referenceRadius reference radius of the Earth for the potential model (m)
  248.      * @param mu central attraction coefficient (m³/s²)
  249.      * @param c20 un-normalized zonal coefficient (about -1.08e-3 for Earth)
  250.      * @param c30 un-normalized zonal coefficient (about +2.53e-6 for Earth)
  251.      * @param c40 un-normalized zonal coefficient (about +1.62e-6 for Earth)
  252.      * @param c50 un-normalized zonal coefficient (about +2.28e-7 for Earth)
  253.      * @param c60 un-normalized zonal coefficient (about -5.41e-7 for Earth)
  254.      */
  255.     public EcksteinHechlerPropagator(final Orbit initialOrbit,
  256.                                      final AttitudeProvider attitudeProv,
  257.                                      final double referenceRadius, final double mu,
  258.                                      final double c20, final double c30, final double c40,
  259.                                      final double c50, final double c60) {
  260.         this(initialOrbit, attitudeProv, DEFAULT_MASS, referenceRadius, mu, c20, c30, c40, c50, c60);
  261.     }

  262.     /** Build a propagator from orbit, attitude provider, mass and potential provider.
  263.      * <p>Using this constructor, an initial osculating orbit is considered.</p>
  264.      * @param initialOrbit initial orbit
  265.      * @param attitudeProv attitude provider
  266.      * @param mass spacecraft mass
  267.      * @param provider for un-normalized zonal coefficients
  268.      * @see #EcksteinHechlerPropagator(Orbit, AttitudeProvider, double,
  269.      *                                 UnnormalizedSphericalHarmonicsProvider, PropagationType)
  270.      */
  271.     public EcksteinHechlerPropagator(final Orbit initialOrbit,
  272.                                      final AttitudeProvider attitudeProv,
  273.                                      final double mass,
  274.                                      final UnnormalizedSphericalHarmonicsProvider provider) {
  275.         this(initialOrbit, attitudeProv, mass, provider, provider.onDate(initialOrbit.getDate()));
  276.     }

  277.     /** Build a propagator from orbit, attitude provider, mass and potential.
  278.      * <p>The C<sub>n,0</sub> coefficients are the denormalized zonal coefficients, they
  279.      * are related to both the normalized coefficients
  280.      * <span style="text-decoration: overline">C</span><sub>n,0</sub>
  281.      *  and the J<sub>n</sub> one as follows:</p>
  282.      *
  283.      * <p> C<sub>n,0</sub> = [(2-δ<sub>0,m</sub>)(2n+1)(n-m)!/(n+m)!]<sup>½</sup>
  284.      * <span style="text-decoration: overline">C</span><sub>n,0</sub>
  285.      *
  286.      * <p> C<sub>n,0</sub> = -J<sub>n</sub>
  287.      *
  288.      * <p>Using this constructor, an initial osculating orbit is considered.</p>
  289.      *
  290.      * @param initialOrbit initial orbit
  291.      * @param attitudeProv attitude provider
  292.      * @param mass spacecraft mass
  293.      * @param referenceRadius reference radius of the Earth for the potential model (m)
  294.      * @param mu central attraction coefficient (m³/s²)
  295.      * @param c20 un-normalized zonal coefficient (about -1.08e-3 for Earth)
  296.      * @param c30 un-normalized zonal coefficient (about +2.53e-6 for Earth)
  297.      * @param c40 un-normalized zonal coefficient (about +1.62e-6 for Earth)
  298.      * @param c50 un-normalized zonal coefficient (about +2.28e-7 for Earth)
  299.      * @param c60 un-normalized zonal coefficient (about -5.41e-7 for Earth)
  300.      * @see #EcksteinHechlerPropagator(Orbit, AttitudeProvider, double, double, double,
  301.      *                                 double, double, double, double, double, PropagationType)
  302.      */
  303.     public EcksteinHechlerPropagator(final Orbit initialOrbit,
  304.                                      final AttitudeProvider attitudeProv,
  305.                                      final double mass,
  306.                                      final double referenceRadius, final double mu,
  307.                                      final double c20, final double c30, final double c40,
  308.                                      final double c50, final double c60) {
  309.         this(initialOrbit, attitudeProv, mass, referenceRadius, mu, c20, c30, c40, c50, c60,
  310.              PropagationType.OSCULATING);
  311.     }


  312.     /** Build a propagator from orbit and potential provider.
  313.      * <p>Mass and attitude provider are set to unspecified non-null arbitrary values.</p>
  314.      *
  315.      * <p>Using this constructor, it is possible to define the initial orbit as
  316.      * a mean Eckstein-Hechler orbit or an osculating one.</p>
  317.      *
  318.      * @param initialOrbit initial orbit
  319.      * @param provider for un-normalized zonal coefficients
  320.      * @param initialType initial orbit type (mean Eckstein-Hechler orbit or osculating orbit)
  321.      * @since 10.2
  322.      */
  323.     public EcksteinHechlerPropagator(final Orbit initialOrbit,
  324.                                      final UnnormalizedSphericalHarmonicsProvider provider,
  325.                                      final PropagationType initialType) {
  326.         this(initialOrbit, FrameAlignedProvider.of(initialOrbit.getFrame()),
  327.              DEFAULT_MASS, provider, provider.onDate(initialOrbit.getDate()), initialType);
  328.     }

  329.     /** Build a propagator from orbit, attitude provider, mass and potential provider.
  330.      * <p>Using this constructor, it is possible to define the initial orbit as
  331.      * a mean Eckstein-Hechler orbit or an osculating one.</p>
  332.      * @param initialOrbit initial orbit
  333.      * @param attitudeProv attitude provider
  334.      * @param mass spacecraft mass
  335.      * @param provider for un-normalized zonal coefficients
  336.      * @param initialType initial orbit type (mean Eckstein-Hechler orbit or osculating orbit)
  337.      * @since 10.2
  338.      */
  339.     public EcksteinHechlerPropagator(final Orbit initialOrbit,
  340.                                      final AttitudeProvider attitudeProv,
  341.                                      final double mass,
  342.                                      final UnnormalizedSphericalHarmonicsProvider provider,
  343.                                      final PropagationType initialType) {
  344.         this(initialOrbit, attitudeProv, mass, provider, provider.onDate(initialOrbit.getDate()), initialType);
  345.     }

  346.     /**
  347.      * Private helper constructor.
  348.      * <p>Using this constructor, it is possible to define the initial orbit as
  349.      * a mean Eckstein-Hechler orbit or an osculating one.</p>
  350.      * @param initialOrbit initial orbit
  351.      * @param attitude attitude provider
  352.      * @param mass spacecraft mass
  353.      * @param provider for un-normalized zonal coefficients
  354.      * @param harmonics {@code provider.onDate(initialOrbit.getDate())}
  355.      * @param initialType initial orbit type (mean Eckstein-Hechler orbit or osculating orbit)
  356.      * @since 10.2
  357.      */
  358.     public EcksteinHechlerPropagator(final Orbit initialOrbit,
  359.                                      final AttitudeProvider attitude,
  360.                                      final double mass,
  361.                                      final UnnormalizedSphericalHarmonicsProvider provider,
  362.                                      final UnnormalizedSphericalHarmonics harmonics,
  363.                                      final PropagationType initialType) {
  364.         this(initialOrbit, attitude, mass, provider.getAe(), provider.getMu(),
  365.              harmonics.getUnnormalizedCnm(2, 0),
  366.              harmonics.getUnnormalizedCnm(3, 0),
  367.              harmonics.getUnnormalizedCnm(4, 0),
  368.              harmonics.getUnnormalizedCnm(5, 0),
  369.              harmonics.getUnnormalizedCnm(6, 0),
  370.              initialType);
  371.     }

  372.     /** Build a propagator from orbit, attitude provider, mass and potential.
  373.      * <p>The C<sub>n,0</sub> coefficients are the denormalized zonal coefficients, they
  374.      * are related to both the normalized coefficients
  375.      * <span style="text-decoration: overline">C</span><sub>n,0</sub>
  376.      *  and the J<sub>n</sub> one as follows:</p>
  377.      *
  378.      * <p> C<sub>n,0</sub> = [(2-δ<sub>0,m</sub>)(2n+1)(n-m)!/(n+m)!]<sup>½</sup>
  379.      * <span style="text-decoration: overline">C</span><sub>n,0</sub>
  380.      *
  381.      * <p> C<sub>n,0</sub> = -J<sub>n</sub>
  382.      *
  383.      * <p>Using this constructor, it is possible to define the initial orbit as
  384.      * a mean Eckstein-Hechler orbit or an osculating one.</p>
  385.      *
  386.      * @param initialOrbit initial orbit
  387.      * @param attitudeProv attitude provider
  388.      * @param mass spacecraft mass
  389.      * @param referenceRadius reference radius of the Earth for the potential model (m)
  390.      * @param mu central attraction coefficient (m³/s²)
  391.      * @param c20 un-normalized zonal coefficient (about -1.08e-3 for Earth)
  392.      * @param c30 un-normalized zonal coefficient (about +2.53e-6 for Earth)
  393.      * @param c40 un-normalized zonal coefficient (about +1.62e-6 for Earth)
  394.      * @param c50 un-normalized zonal coefficient (about +2.28e-7 for Earth)
  395.      * @param c60 un-normalized zonal coefficient (about -5.41e-7 for Earth)
  396.      * @param initialType initial orbit type (mean Eckstein-Hechler orbit or osculating orbit)
  397.      * @since 10.2
  398.      */
  399.     public EcksteinHechlerPropagator(final Orbit initialOrbit,
  400.                                      final AttitudeProvider attitudeProv,
  401.                                      final double mass,
  402.                                      final double referenceRadius, final double mu,
  403.                                      final double c20, final double c30, final double c40,
  404.                                      final double c50, final double c60,
  405.                                      final PropagationType initialType) {
  406.         this(initialOrbit, attitudeProv, mass, referenceRadius, mu,
  407.              c20, c30, c40, c50, c60, initialType, EPSILON_DEFAULT, MAX_ITERATIONS_DEFAULT);
  408.     }

  409.     /** Build a propagator from orbit, attitude provider, mass and potential.
  410.      * <p>The C<sub>n,0</sub> coefficients are the denormalized zonal coefficients, they
  411.      * are related to both the normalized coefficients
  412.      * <span style="text-decoration: overline">C</span><sub>n,0</sub>
  413.      *  and the J<sub>n</sub> one as follows:</p>
  414.      *
  415.      * <p> C<sub>n,0</sub> = [(2-δ<sub>0,m</sub>)(2n+1)(n-m)!/(n+m)!]<sup>½</sup>
  416.      * <span style="text-decoration: overline">C</span><sub>n,0</sub>
  417.      *
  418.      * <p> C<sub>n,0</sub> = -J<sub>n</sub>
  419.      *
  420.      * <p>Using this constructor, it is possible to define the initial orbit as
  421.      * a mean Eckstein-Hechler orbit or an osculating one.</p>
  422.      *
  423.      * @param initialOrbit initial orbit
  424.      * @param attitudeProv attitude provider
  425.      * @param mass spacecraft mass
  426.      * @param referenceRadius reference radius of the Earth for the potential model (m)
  427.      * @param mu central attraction coefficient (m³/s²)
  428.      * @param c20 un-normalized zonal coefficient (about -1.08e-3 for Earth)
  429.      * @param c30 un-normalized zonal coefficient (about +2.53e-6 for Earth)
  430.      * @param c40 un-normalized zonal coefficient (about +1.62e-6 for Earth)
  431.      * @param c50 un-normalized zonal coefficient (about +2.28e-7 for Earth)
  432.      * @param c60 un-normalized zonal coefficient (about -5.41e-7 for Earth)
  433.      * @param epsilon convergence threshold for mean parameters conversion
  434.      * @param maxIterations maximum iterations for mean parameters conversion
  435.      * @param initialType initial orbit type (mean Eckstein-Hechler orbit or osculating orbit)
  436.      * @since 11.2
  437.      */
  438.     public EcksteinHechlerPropagator(final Orbit initialOrbit,
  439.                                      final AttitudeProvider attitudeProv,
  440.                                      final double mass,
  441.                                      final double referenceRadius, final double mu,
  442.                                      final double c20, final double c30, final double c40,
  443.                                      final double c50, final double c60,
  444.                                      final PropagationType initialType,
  445.                                      final double epsilon, final int maxIterations) {
  446.         super(attitudeProv);

  447.         // store model coefficients
  448.         this.referenceRadius = referenceRadius;
  449.         this.mu  = mu;
  450.         this.ck0 = new double[] {
  451.             0.0, 0.0, c20, c30, c40, c50, c60
  452.         };

  453.         // compute mean parameters if needed
  454.         // transform into circular adapted parameters used by the Eckstein-Hechler model
  455.         resetInitialState(new SpacecraftState(initialOrbit,
  456.                                               attitudeProv.getAttitude(initialOrbit,
  457.                                                                        initialOrbit.getDate(),
  458.                                                                        initialOrbit.getFrame()),
  459.                                               mass),
  460.                           initialType, epsilon, maxIterations);

  461.     }

  462.     /** Conversion from osculating to mean orbit.
  463.      * <p>
  464.      * Compute mean orbit <b>in a Eckstein-Hechler sense</b>, corresponding to the
  465.      * osculating SpacecraftState in input.
  466.      * </p>
  467.      * <p>
  468.      * Since the osculating orbit is obtained with the computation of
  469.      * short-periodic variation, the resulting output will depend on
  470.      * the gravity field parameterized in input.
  471.      * </p>
  472.      * <p>
  473.      * The computation is done through a fixed-point iteration process.
  474.      * </p>
  475.      * @param osculating osculating orbit to convert
  476.      * @param provider for un-normalized zonal coefficients
  477.      * @param harmonics {@code provider.onDate(osculating.getDate())}
  478.      * @return mean orbit in a Eckstein-Hechler sense
  479.      * @since 11.2
  480.      */
  481.     public static CircularOrbit computeMeanOrbit(final Orbit osculating,
  482.                                                  final UnnormalizedSphericalHarmonicsProvider provider,
  483.                                                  final UnnormalizedSphericalHarmonics harmonics) {
  484.         return computeMeanOrbit(osculating, provider, harmonics, EPSILON_DEFAULT, MAX_ITERATIONS_DEFAULT);
  485.     }

  486.     /** Conversion from osculating to mean orbit.
  487.      * <p>
  488.      * Compute mean orbit <b>in a Eckstein-Hechler sense</b>, corresponding to the
  489.      * osculating SpacecraftState in input.
  490.      * </p>
  491.      * <p>
  492.      * Since the osculating orbit is obtained with the computation of
  493.      * short-periodic variation, the resulting output will depend on
  494.      * the gravity field parameterized in input.
  495.      * </p>
  496.      * <p>
  497.      * The computation is done through a fixed-point iteration process.
  498.      * </p>
  499.      * @param osculating osculating orbit to convert
  500.      * @param provider for un-normalized zonal coefficients
  501.      * @param harmonics {@code provider.onDate(osculating.getDate())}
  502.      * @param epsilon convergence threshold for mean parameters conversion
  503.      * @param maxIterations maximum iterations for mean parameters conversion
  504.      * @return mean orbit in a Eckstein-Hechler sense
  505.      * @since 11.2
  506.      */
  507.     public static CircularOrbit computeMeanOrbit(final Orbit osculating,
  508.                                                  final UnnormalizedSphericalHarmonicsProvider provider,
  509.                                                  final UnnormalizedSphericalHarmonics harmonics,
  510.                                                  final double epsilon, final int maxIterations) {
  511.         return computeMeanOrbit(osculating,
  512.                                 provider.getAe(), provider.getMu(),
  513.                                 harmonics.getUnnormalizedCnm(2, 0),
  514.                                 harmonics.getUnnormalizedCnm(3, 0),
  515.                                 harmonics.getUnnormalizedCnm(4, 0),
  516.                                 harmonics.getUnnormalizedCnm(5, 0),
  517.                                 harmonics.getUnnormalizedCnm(6, 0),
  518.                                 epsilon, maxIterations);
  519.     }

  520.     /** Conversion from osculating to mean orbit.
  521.      * <p>
  522.      * Compute mean orbit <b>in a Eckstein-Hechler sense</b>, corresponding to the
  523.      * osculating SpacecraftState in input.
  524.      * </p>
  525.      * <p>
  526.      * Since the osculating orbit is obtained with the computation of
  527.      * short-periodic variation, the resulting output will depend on
  528.      * the gravity field parameterized in input.
  529.      * </p>
  530.      * <p>
  531.      * The computation is done through a fixed-point iteration process.
  532.      * </p>
  533.      * @param osculating osculating orbit to convert
  534.      * @param referenceRadius reference radius of the Earth for the potential model (m)
  535.      * @param mu central attraction coefficient (m³/s²)
  536.      * @param c20 un-normalized zonal coefficient (about -1.08e-3 for Earth)
  537.      * @param c30 un-normalized zonal coefficient (about +2.53e-6 for Earth)
  538.      * @param c40 un-normalized zonal coefficient (about +1.62e-6 for Earth)
  539.      * @param c50 un-normalized zonal coefficient (about +2.28e-7 for Earth)
  540.      * @param c60 un-normalized zonal coefficient (about -5.41e-7 for Earth)
  541.      * @param epsilon convergence threshold for mean parameters conversion
  542.      * @param maxIterations maximum iterations for mean parameters conversion
  543.      * @return mean orbit in a Eckstein-Hechler sense
  544.      * @since 11.2
  545.      */
  546.     public static CircularOrbit computeMeanOrbit(final Orbit osculating,
  547.                                                  final double referenceRadius, final double mu,
  548.                                                  final double c20, final double c30, final double c40,
  549.                                                  final double c50, final double c60,
  550.                                                  final double epsilon, final int maxIterations) {
  551.         final EcksteinHechlerPropagator propagator =
  552.                         new EcksteinHechlerPropagator(osculating,
  553.                                                       FrameAlignedProvider.of(osculating.getFrame()),
  554.                                                       DEFAULT_MASS,
  555.                                                       referenceRadius, mu, c20, c30, c40, c50, c60,
  556.                                                       PropagationType.OSCULATING, epsilon, maxIterations);
  557.         return propagator.initialModel.mean;
  558.     }

  559.     /** {@inheritDoc}
  560.      * <p>The new initial state to consider
  561.      * must be defined with an osculating orbit.</p>
  562.      * @see #resetInitialState(SpacecraftState, PropagationType)
  563.      */
  564.     public void resetInitialState(final SpacecraftState state) {
  565.         resetInitialState(state, PropagationType.OSCULATING);
  566.     }

  567.     /** Reset the propagator initial state.
  568.      * @param state new initial state to consider
  569.      * @param stateType mean Eckstein-Hechler orbit or osculating orbit
  570.      * @since 10.2
  571.      */
  572.     public void resetInitialState(final SpacecraftState state, final PropagationType stateType) {
  573.         resetInitialState(state, stateType, EPSILON_DEFAULT, MAX_ITERATIONS_DEFAULT);
  574.     }

  575.     /** Reset the propagator initial state.
  576.      * @param state new initial state to consider
  577.      * @param stateType mean Eckstein-Hechler orbit or osculating orbit
  578.      * @param epsilon convergence threshold for mean parameters conversion
  579.      * @param maxIterations maximum iterations for mean parameters conversion
  580.      * @since 11.2
  581.      */
  582.     public void resetInitialState(final SpacecraftState state, final PropagationType stateType,
  583.                                   final double epsilon, final int maxIterations) {
  584.         super.resetInitialState(state);
  585.         final CircularOrbit circular = (CircularOrbit) OrbitType.CIRCULAR.convertType(state.getOrbit());
  586.         this.initialModel = (stateType == PropagationType.MEAN) ?
  587.                              new EHModel(circular, state.getMass(), referenceRadius, mu, ck0) :
  588.                              computeMeanParameters(circular, state.getMass(), epsilon, maxIterations);
  589.         this.models = new TimeSpanMap<EHModel>(initialModel);
  590.     }

  591.     /** {@inheritDoc} */
  592.     protected void resetIntermediateState(final SpacecraftState state, final boolean forward) {
  593.         resetIntermediateState(state, forward, EPSILON_DEFAULT, MAX_ITERATIONS_DEFAULT);
  594.     }

  595.     /** Reset an intermediate state.
  596.      * @param state new intermediate state to consider
  597.      * @param forward if true, the intermediate state is valid for
  598.      * propagations after itself
  599.      * @param epsilon convergence threshold for mean parameters conversion
  600.      * @param maxIterations maximum iterations for mean parameters conversion
  601.      * @since 11.2
  602.      */
  603.     protected void resetIntermediateState(final SpacecraftState state, final boolean forward,
  604.                                           final double epsilon, final int maxIterations) {
  605.         final EHModel newModel = computeMeanParameters((CircularOrbit) OrbitType.CIRCULAR.convertType(state.getOrbit()),
  606.                                                        state.getMass(), epsilon, maxIterations);
  607.         if (forward) {
  608.             models.addValidAfter(newModel, state.getDate(), false);
  609.         } else {
  610.             models.addValidBefore(newModel, state.getDate(), false);
  611.         }
  612.         stateChanged(state);
  613.     }


  614.     /** Compute mean parameters according to the Eckstein-Hechler analytical model.
  615.      * @param osculating osculating orbit
  616.      * @param mass constant mass
  617.      * @param epsilon convergence threshold for mean parameters conversion
  618.      * @param maxIterations maximum iterations for mean parameters conversion
  619.      * @return Eckstein-Hechler mean model
  620.      */
  621.     private EHModel computeMeanParameters(final CircularOrbit osculating, final double mass,
  622.                                           final double epsilon, final int maxIterations) {

  623.         // sanity check
  624.         if (osculating.getA() < referenceRadius) {
  625.             throw new OrekitException(OrekitMessages.TRAJECTORY_INSIDE_BRILLOUIN_SPHERE,
  626.                                       osculating.getA());
  627.         }

  628.         // rough initialization of the mean parameters
  629.         EHModel current = new EHModel(osculating, mass, referenceRadius, mu, ck0);

  630.         // threshold for each parameter
  631.         final double thresholdA      = epsilon * (1 + FastMath.abs(current.mean.getA()));
  632.         final double thresholdE      = epsilon * (1 + current.mean.getE());
  633.         final double thresholdAngles = epsilon * MathUtils.TWO_PI;

  634.         int i = 0;
  635.         while (i++ < maxIterations) {

  636.             // recompute the osculating parameters from the current mean parameters
  637.             final UnivariateDerivative2[] parameters = current.propagateParameters(current.mean.getDate());

  638.             // adapted parameters residuals
  639.             final double deltaA      = osculating.getA()          - parameters[0].getValue();
  640.             final double deltaEx     = osculating.getCircularEx() - parameters[1].getValue();
  641.             final double deltaEy     = osculating.getCircularEy() - parameters[2].getValue();
  642.             final double deltaI      = osculating.getI()          - parameters[3].getValue();
  643.             final double deltaRAAN   = MathUtils.normalizeAngle(osculating.getRightAscensionOfAscendingNode() -
  644.                                                                 parameters[4].getValue(),
  645.                                                                 0.0);
  646.             final double deltaAlphaM = MathUtils.normalizeAngle(osculating.getAlphaM() - parameters[5].getValue(), 0.0);

  647.             // update mean parameters
  648.             current = new EHModel(new CircularOrbit(current.mean.getA()          + deltaA,
  649.                                                     current.mean.getCircularEx() + deltaEx,
  650.                                                     current.mean.getCircularEy() + deltaEy,
  651.                                                     current.mean.getI()          + deltaI,
  652.                                                     current.mean.getRightAscensionOfAscendingNode() + deltaRAAN,
  653.                                                     current.mean.getAlphaM()     + deltaAlphaM,
  654.                                                     PositionAngleType.MEAN,
  655.                                                     current.mean.getFrame(),
  656.                                                     current.mean.getDate(), mu),
  657.                                   mass, referenceRadius, mu, ck0);

  658.             // check convergence
  659.             if (FastMath.abs(deltaA)      < thresholdA &&
  660.                 FastMath.abs(deltaEx)     < thresholdE &&
  661.                 FastMath.abs(deltaEy)     < thresholdE &&
  662.                 FastMath.abs(deltaI)      < thresholdAngles &&
  663.                 FastMath.abs(deltaRAAN)   < thresholdAngles &&
  664.                 FastMath.abs(deltaAlphaM) < thresholdAngles) {
  665.                 return current;
  666.             }

  667.         }

  668.         throw new OrekitException(OrekitMessages.UNABLE_TO_COMPUTE_ECKSTEIN_HECHLER_MEAN_PARAMETERS, i);

  669.     }

  670.     /** {@inheritDoc} */
  671.     public CartesianOrbit propagateOrbit(final AbsoluteDate date) {
  672.         // compute Cartesian parameters, taking derivatives into account
  673.         // to make sure velocity and acceleration are consistent
  674.         final EHModel current = models.get(date);
  675.         return new CartesianOrbit(toCartesian(date, current.propagateParameters(date)),
  676.                                   current.mean.getFrame(), mu);
  677.     }

  678.     /**
  679.      * Get the central attraction coefficient μ.
  680.      * @return mu central attraction coefficient (m³/s²)
  681.      * @since 11.1
  682.      */
  683.     public double getMu() {
  684.         return mu;
  685.     }

  686.     /**
  687.      * Get the un-normalized zonal coefficients.
  688.      * @return the un-normalized zonal coefficients
  689.      * @since 11.1
  690.      */
  691.     public double[] getCk0() {
  692.         return ck0.clone();
  693.     }

  694.     /**
  695.      * Get the reference radius of the central body attraction model.
  696.      * @return the reference radius in meters
  697.      * @since 11.1
  698.      */
  699.     public double getReferenceRadius() {
  700.         return referenceRadius;
  701.     }

  702.     /** {@inheritDoc} */
  703.     @Override
  704.     protected AbstractMatricesHarvester createHarvester(final String stmName, final RealMatrix initialStm,
  705.                                                         final DoubleArrayDictionary initialJacobianColumns) {
  706.         // Create the harvester
  707.         final EcksteinHechlerHarvester harvester = new EcksteinHechlerHarvester(this, stmName, initialStm, initialJacobianColumns);
  708.         // Update the list of additional state provider
  709.         addAdditionalStateProvider(harvester);
  710.         // Return the configured harvester
  711.         return harvester;
  712.     }

  713.     /** Local class for Eckstein-Hechler model, with fixed mean parameters. */
  714.     private static class EHModel implements Serializable {

  715.         /** Serializable UID. */
  716.         private static final long serialVersionUID = 20160115L;

  717.         /** Mean orbit. */
  718.         private final CircularOrbit mean;

  719.         /** Constant mass. */
  720.         private final double mass;

  721.         // CHECKSTYLE: stop JavadocVariable check

  722.         // preprocessed values
  723.         private final double xnotDot;
  724.         private final double rdpom;
  725.         private final double rdpomp;
  726.         private final double eps1;
  727.         private final double eps2;
  728.         private final double xim;
  729.         private final double ommD;
  730.         private final double rdl;
  731.         private final double aMD;

  732.         private final double kh;
  733.         private final double kl;

  734.         private final double ax1;
  735.         private final double ay1;
  736.         private final double as1;
  737.         private final double ac2;
  738.         private final double axy3;
  739.         private final double as3;
  740.         private final double ac4;
  741.         private final double as5;
  742.         private final double ac6;

  743.         private final double ex1;
  744.         private final double exx2;
  745.         private final double exy2;
  746.         private final double ex3;
  747.         private final double ex4;

  748.         private final double ey1;
  749.         private final double eyx2;
  750.         private final double eyy2;
  751.         private final double ey3;
  752.         private final double ey4;

  753.         private final double rx1;
  754.         private final double ry1;
  755.         private final double r2;
  756.         private final double r3;
  757.         private final double rl;

  758.         private final double iy1;
  759.         private final double ix1;
  760.         private final double i2;
  761.         private final double i3;
  762.         private final double ih;

  763.         private final double lx1;
  764.         private final double ly1;
  765.         private final double l2;
  766.         private final double l3;
  767.         private final double ll;

  768.         // CHECKSTYLE: resume JavadocVariable check

  769.         /** Create a model for specified mean orbit.
  770.          * @param mean mean orbit
  771.          * @param mass constant mass
  772.          * @param referenceRadius reference radius of the central body attraction model (m)
  773.          * @param mu central attraction coefficient (m³/s²)
  774.          * @param ck0 un-normalized zonal coefficients
  775.          */
  776.         EHModel(final CircularOrbit mean, final double mass,
  777.                 final double referenceRadius, final double mu, final double[] ck0) {

  778.             this.mean            = mean;
  779.             this.mass            = mass;

  780.             // preliminary processing
  781.             double q = referenceRadius / mean.getA();
  782.             double ql = q * q;
  783.             final double g2 = ck0[2] * ql;
  784.             ql *= q;
  785.             final double g3 = ck0[3] * ql;
  786.             ql *= q;
  787.             final double g4 = ck0[4] * ql;
  788.             ql *= q;
  789.             final double g5 = ck0[5] * ql;
  790.             ql *= q;
  791.             final double g6 = ck0[6] * ql;

  792.             final SinCos sc    = FastMath.sinCos(mean.getI());
  793.             final double cosI1 = sc.cos();
  794.             final double sinI1 = sc.sin();
  795.             final double sinI2 = sinI1 * sinI1;
  796.             final double sinI4 = sinI2 * sinI2;
  797.             final double sinI6 = sinI2 * sinI4;

  798.             if (sinI2 < 1.0e-10) {
  799.                 throw new OrekitException(OrekitMessages.ALMOST_EQUATORIAL_ORBIT,
  800.                                           FastMath.toDegrees(mean.getI()));
  801.             }

  802.             if (FastMath.abs(sinI2 - 4.0 / 5.0) < 1.0e-3) {
  803.                 throw new OrekitException(OrekitMessages.ALMOST_CRITICALLY_INCLINED_ORBIT,
  804.                                           FastMath.toDegrees(mean.getI()));
  805.             }

  806.             if (mean.getE() > 0.1) {
  807.                 // if 0.005 < e < 0.1 no error is triggered, but accuracy is poor
  808.                 throw new OrekitException(OrekitMessages.TOO_LARGE_ECCENTRICITY_FOR_PROPAGATION_MODEL,
  809.                                           mean.getE());
  810.             }

  811.             xnotDot = FastMath.sqrt(mu / mean.getA()) / mean.getA();

  812.             rdpom = -0.75 * g2 * (4.0 - 5.0 * sinI2);
  813.             rdpomp = 7.5 * g4 * (1.0 - 31.0 / 8.0 * sinI2 + 49.0 / 16.0 * sinI4) -
  814.                     13.125 * g6 * (1.0 - 8.0 * sinI2 + 129.0 / 8.0 * sinI4 - 297.0 / 32.0 * sinI6);

  815.             q = 3.0 / (32.0 * rdpom);
  816.             eps1 = q * g4 * sinI2 * (30.0 - 35.0 * sinI2) -
  817.                     175.0 * q * g6 * sinI2 * (1.0 - 3.0 * sinI2 + 2.0625 * sinI4);
  818.             q = 3.0 * sinI1 / (8.0 * rdpom);
  819.             eps2 = q * g3 * (4.0 - 5.0 * sinI2) - q * g5 * (10.0 - 35.0 * sinI2 + 26.25 * sinI4);

  820.             xim = mean.getI();
  821.             ommD = cosI1 * (1.50    * g2 - 2.25 * g2 * g2 * (2.5 - 19.0 / 6.0 * sinI2) +
  822.                             0.9375  * g4 * (7.0 * sinI2 - 4.0) +
  823.                             3.28125 * g6 * (2.0 - 9.0 * sinI2 + 8.25 * sinI4));

  824.             rdl = 1.0 - 1.50 * g2 * (3.0 - 4.0 * sinI2);
  825.             aMD = rdl +
  826.                     2.25 * g2 * g2 * (9.0 - 263.0 / 12.0 * sinI2 + 341.0 / 24.0 * sinI4) +
  827.                     15.0 / 16.0 * g4 * (8.0 - 31.0 * sinI2 + 24.5 * sinI4) +
  828.                     105.0 / 32.0 * g6 * (-10.0 / 3.0 + 25.0 * sinI2 - 48.75 * sinI4 + 27.5 * sinI6);

  829.             final double qq = -1.5 * g2 / rdl;
  830.             final double qA   = 0.75 * g2 * g2 * sinI2;
  831.             final double qB   = 0.25 * g4 * sinI2;
  832.             final double qC   = 105.0 / 16.0 * g6 * sinI2;
  833.             final double qD   = -0.75 * g3 * sinI1;
  834.             final double qE   = 3.75 * g5 * sinI1;
  835.             kh = 0.375 / rdpom;
  836.             kl = kh / sinI1;

  837.             ax1 = qq * (2.0 - 3.5 * sinI2);
  838.             ay1 = qq * (2.0 - 2.5 * sinI2);
  839.             as1 = qD * (4.0 - 5.0 * sinI2) +
  840.                   qE * (2.625 * sinI4 - 3.5 * sinI2 + 1.0);
  841.             ac2 = qq * sinI2 +
  842.                   qA * 7.0 * (2.0 - 3.0 * sinI2) +
  843.                   qB * (15.0 - 17.5 * sinI2) +
  844.                   qC * (3.0 * sinI2 - 1.0 - 33.0 / 16.0 * sinI4);
  845.             axy3 = qq * 3.5 * sinI2;
  846.             as3 = qD * 5.0 / 3.0 * sinI2 +
  847.                   qE * 7.0 / 6.0 * sinI2 * (1.0 - 1.125 * sinI2);
  848.             ac4 = qA * sinI2 +
  849.                   qB * 4.375 * sinI2 +
  850.                   qC * 0.75 * (1.1 * sinI4 - sinI2);

  851.             as5 = qE * 21.0 / 80.0 * sinI4;

  852.             ac6 = qC * -11.0 / 80.0 * sinI4;

  853.             ex1 = qq * (1.0 - 1.25 * sinI2);
  854.             exx2 = qq * 0.5 * (3.0 - 5.0 * sinI2);
  855.             exy2 = qq * (2.0 - 1.5 * sinI2);
  856.             ex3 = qq * 7.0 / 12.0 * sinI2;
  857.             ex4 = qq * 17.0 / 8.0 * sinI2;

  858.             ey1 = qq * (1.0 - 1.75 * sinI2);
  859.             eyx2 = qq * (1.0 - 3.0 * sinI2);
  860.             eyy2 = qq * (2.0 * sinI2 - 1.5);
  861.             ey3 = qq * 7.0 / 12.0 * sinI2;
  862.             ey4 = qq * 17.0 / 8.0 * sinI2;

  863.             q  = -qq * cosI1;
  864.             rx1 =  3.5 * q;
  865.             ry1 = -2.5 * q;
  866.             r2 = -0.5 * q;
  867.             r3 =  7.0 / 6.0 * q;
  868.             rl = g3 * cosI1 * (4.0 - 15.0 * sinI2) -
  869.                  2.5 * g5 * cosI1 * (4.0 - 42.0 * sinI2 + 52.5 * sinI4);

  870.             q = 0.5 * qq * sinI1 * cosI1;
  871.             iy1 =  q;
  872.             ix1 = -q;
  873.             i2 =  q;
  874.             i3 =  q * 7.0 / 3.0;
  875.             ih = -g3 * cosI1 * (4.0 - 5.0 * sinI2) +
  876.                  2.5 * g5 * cosI1 * (4.0 - 14.0 * sinI2 + 10.5 * sinI4);

  877.             lx1 = qq * (7.0 - 77.0 / 8.0 * sinI2);
  878.             ly1 = qq * (55.0 / 8.0 * sinI2 - 7.50);
  879.             l2 = qq * (1.25 * sinI2 - 0.5);
  880.             l3 = qq * (77.0 / 24.0 * sinI2 - 7.0 / 6.0);
  881.             ll = g3 * (53.0 * sinI2 - 4.0 - 57.5 * sinI4) +
  882.                  2.5 * g5 * (4.0 - 96.0 * sinI2 + 269.5 * sinI4 - 183.75 * sinI6);

  883.         }

  884.         /** Extrapolate an orbit up to a specific target date.
  885.          * @param date target date for the orbit
  886.          * @return propagated parameters
  887.          */
  888.         public UnivariateDerivative2[] propagateParameters(final AbsoluteDate date) {

  889.             // Keplerian evolution
  890.             final UnivariateDerivative2 dt = new UnivariateDerivative2(date.durationFrom(mean.getDate()), 1.0, 0.0);
  891.             final UnivariateDerivative2 xnot = dt.multiply(xnotDot);

  892.             // secular effects

  893.             // eccentricity
  894.             final UnivariateDerivative2 x   = xnot.multiply(rdpom + rdpomp);
  895.             final UnivariateDerivative2 cx  = x.cos();
  896.             final UnivariateDerivative2 sx  = x.sin();
  897.             final UnivariateDerivative2 exm = cx.multiply(mean.getCircularEx()).
  898.                                               add(sx.multiply(eps2 - (1.0 - eps1) * mean.getCircularEy()));
  899.             final UnivariateDerivative2 eym = sx.multiply((1.0 + eps1) * mean.getCircularEx()).
  900.                                               add(cx.multiply(mean.getCircularEy() - eps2)).
  901.                                               add(eps2);

  902.             // no secular effect on inclination

  903.             // right ascension of ascending node
  904.             final UnivariateDerivative2 omm = new UnivariateDerivative2(MathUtils.normalizeAngle(mean.getRightAscensionOfAscendingNode() + ommD * xnot.getValue(),
  905.                                                                                                FastMath.PI),
  906.                                                                       ommD * xnotDot,
  907.                                                                       0.0);

  908.             // latitude argument
  909.             final UnivariateDerivative2 xlm = new UnivariateDerivative2(MathUtils.normalizeAngle(mean.getAlphaM() + aMD * xnot.getValue(),
  910.                                                                                                  FastMath.PI),
  911.                                                                         aMD * xnotDot,
  912.                                                                         0.0);

  913.             // periodical terms
  914.             final UnivariateDerivative2 cl1 = xlm.cos();
  915.             final UnivariateDerivative2 sl1 = xlm.sin();
  916.             final UnivariateDerivative2 cl2 = cl1.multiply(cl1).subtract(sl1.multiply(sl1));
  917.             final UnivariateDerivative2 sl2 = cl1.multiply(sl1).add(sl1.multiply(cl1));
  918.             final UnivariateDerivative2 cl3 = cl2.multiply(cl1).subtract(sl2.multiply(sl1));
  919.             final UnivariateDerivative2 sl3 = cl2.multiply(sl1).add(sl2.multiply(cl1));
  920.             final UnivariateDerivative2 cl4 = cl3.multiply(cl1).subtract(sl3.multiply(sl1));
  921.             final UnivariateDerivative2 sl4 = cl3.multiply(sl1).add(sl3.multiply(cl1));
  922.             final UnivariateDerivative2 cl5 = cl4.multiply(cl1).subtract(sl4.multiply(sl1));
  923.             final UnivariateDerivative2 sl5 = cl4.multiply(sl1).add(sl4.multiply(cl1));
  924.             final UnivariateDerivative2 cl6 = cl5.multiply(cl1).subtract(sl5.multiply(sl1));

  925.             final UnivariateDerivative2 qh  = eym.subtract(eps2).multiply(kh);
  926.             final UnivariateDerivative2 ql  = exm.multiply(kl);

  927.             final UnivariateDerivative2 exmCl1 = exm.multiply(cl1);
  928.             final UnivariateDerivative2 exmSl1 = exm.multiply(sl1);
  929.             final UnivariateDerivative2 eymCl1 = eym.multiply(cl1);
  930.             final UnivariateDerivative2 eymSl1 = eym.multiply(sl1);
  931.             final UnivariateDerivative2 exmCl2 = exm.multiply(cl2);
  932.             final UnivariateDerivative2 exmSl2 = exm.multiply(sl2);
  933.             final UnivariateDerivative2 eymCl2 = eym.multiply(cl2);
  934.             final UnivariateDerivative2 eymSl2 = eym.multiply(sl2);
  935.             final UnivariateDerivative2 exmCl3 = exm.multiply(cl3);
  936.             final UnivariateDerivative2 exmSl3 = exm.multiply(sl3);
  937.             final UnivariateDerivative2 eymCl3 = eym.multiply(cl3);
  938.             final UnivariateDerivative2 eymSl3 = eym.multiply(sl3);
  939.             final UnivariateDerivative2 exmCl4 = exm.multiply(cl4);
  940.             final UnivariateDerivative2 exmSl4 = exm.multiply(sl4);
  941.             final UnivariateDerivative2 eymCl4 = eym.multiply(cl4);
  942.             final UnivariateDerivative2 eymSl4 = eym.multiply(sl4);

  943.             // semi major axis
  944.             final UnivariateDerivative2 rda = exmCl1.multiply(ax1).
  945.                                               add(eymSl1.multiply(ay1)).
  946.                                               add(sl1.multiply(as1)).
  947.                                               add(cl2.multiply(ac2)).
  948.                                               add(exmCl3.add(eymSl3).multiply(axy3)).
  949.                                               add(sl3.multiply(as3)).
  950.                                               add(cl4.multiply(ac4)).
  951.                                               add(sl5.multiply(as5)).
  952.                                               add(cl6.multiply(ac6));

  953.             // eccentricity
  954.             final UnivariateDerivative2 rdex = cl1.multiply(ex1).
  955.                                                add(exmCl2.multiply(exx2)).
  956.                                                add(eymSl2.multiply(exy2)).
  957.                                                add(cl3.multiply(ex3)).
  958.                                                add(exmCl4.add(eymSl4).multiply(ex4));
  959.             final UnivariateDerivative2 rdey = sl1.multiply(ey1).
  960.                                                add(exmSl2.multiply(eyx2)).
  961.                                                add(eymCl2.multiply(eyy2)).
  962.                                                add(sl3.multiply(ey3)).
  963.                                                add(exmSl4.subtract(eymCl4).multiply(ey4));

  964.             // ascending node
  965.             final UnivariateDerivative2 rdom = exmSl1.multiply(rx1).
  966.                                                add(eymCl1.multiply(ry1)).
  967.                                                add(sl2.multiply(r2)).
  968.                                                add(eymCl3.subtract(exmSl3).multiply(r3)).
  969.                                                add(ql.multiply(rl));

  970.             // inclination
  971.             final UnivariateDerivative2 rdxi = eymSl1.multiply(iy1).
  972.                                                add(exmCl1.multiply(ix1)).
  973.                                                add(cl2.multiply(i2)).
  974.                                                add(exmCl3.add(eymSl3).multiply(i3)).
  975.                                                add(qh.multiply(ih));

  976.             // latitude argument
  977.             final UnivariateDerivative2 rdxl = exmSl1.multiply(lx1).
  978.                                                add(eymCl1.multiply(ly1)).
  979.                                                add(sl2.multiply(l2)).
  980.                                                add(exmSl3.subtract(eymCl3).multiply(l3)).
  981.                                                add(ql.multiply(ll));

  982.             // osculating parameters
  983.             return new UnivariateDerivative2[] {
  984.                 rda.add(1.0).multiply(mean.getA()),
  985.                 rdex.add(exm),
  986.                 rdey.add(eym),
  987.                 rdxi.add(xim),
  988.                 rdom.add(omm),
  989.                 rdxl.add(xlm)
  990.             };

  991.         }

  992.     }

  993.     /** Convert circular parameters <em>with derivatives</em> to Cartesian coordinates.
  994.      * @param date date of the orbital parameters
  995.      * @param parameters circular parameters (a, ex, ey, i, raan, alphaM)
  996.      * @return Cartesian coordinates consistent with values and derivatives
  997.      */
  998.     private TimeStampedPVCoordinates toCartesian(final AbsoluteDate date, final UnivariateDerivative2[] parameters) {

  999.         // evaluate coordinates in the orbit canonical reference frame
  1000.         final UnivariateDerivative2 cosOmega = parameters[4].cos();
  1001.         final UnivariateDerivative2 sinOmega = parameters[4].sin();
  1002.         final UnivariateDerivative2 cosI     = parameters[3].cos();
  1003.         final UnivariateDerivative2 sinI     = parameters[3].sin();
  1004.         final UnivariateDerivative2 alphaE   = meanToEccentric(parameters[5], parameters[1], parameters[2]);
  1005.         final UnivariateDerivative2 cosAE    = alphaE.cos();
  1006.         final UnivariateDerivative2 sinAE    = alphaE.sin();
  1007.         final UnivariateDerivative2 ex2      = parameters[1].square();
  1008.         final UnivariateDerivative2 ey2      = parameters[2].square();
  1009.         final UnivariateDerivative2 exy      = parameters[1].multiply(parameters[2]);
  1010.         final UnivariateDerivative2 q        = ex2.add(ey2).subtract(1).negate().sqrt();
  1011.         final UnivariateDerivative2 beta     = q.add(1).reciprocal();
  1012.         final UnivariateDerivative2 bx2      = beta.multiply(ex2);
  1013.         final UnivariateDerivative2 by2      = beta.multiply(ey2);
  1014.         final UnivariateDerivative2 bxy      = beta.multiply(exy);
  1015.         final UnivariateDerivative2 u        = bxy.multiply(sinAE).subtract(parameters[1].add(by2.subtract(1).multiply(cosAE)));
  1016.         final UnivariateDerivative2 v        = bxy.multiply(cosAE).subtract(parameters[2].add(bx2.subtract(1).multiply(sinAE)));
  1017.         final UnivariateDerivative2 x        = parameters[0].multiply(u);
  1018.         final UnivariateDerivative2 y        = parameters[0].multiply(v);

  1019.         // canonical orbit reference frame
  1020.         final FieldVector3D<UnivariateDerivative2> p =
  1021.                 new FieldVector3D<>(x.multiply(cosOmega).subtract(y.multiply(cosI.multiply(sinOmega))),
  1022.                                     x.multiply(sinOmega).add(y.multiply(cosI.multiply(cosOmega))),
  1023.                                     y.multiply(sinI));

  1024.         // dispatch derivatives
  1025.         final Vector3D p0 = new Vector3D(p.getX().getValue(),
  1026.                                          p.getY().getValue(),
  1027.                                          p.getZ().getValue());
  1028.         final Vector3D p1 = new Vector3D(p.getX().getFirstDerivative(),
  1029.                                          p.getY().getFirstDerivative(),
  1030.                                          p.getZ().getFirstDerivative());
  1031.         final Vector3D p2 = new Vector3D(p.getX().getSecondDerivative(),
  1032.                                          p.getY().getSecondDerivative(),
  1033.                                          p.getZ().getSecondDerivative());
  1034.         return new TimeStampedPVCoordinates(date, p0, p1, p2);

  1035.     }

  1036.     /** Computes the eccentric latitude argument from the mean latitude argument.
  1037.      * @param alphaM = M + Ω mean latitude argument (rad)
  1038.      * @param ex e cos(Ω), first component of circular eccentricity vector
  1039.      * @param ey e sin(Ω), second component of circular eccentricity vector
  1040.      * @return the eccentric latitude argument.
  1041.      */
  1042.     private UnivariateDerivative2 meanToEccentric(final UnivariateDerivative2 alphaM,
  1043.                                                   final UnivariateDerivative2 ex,
  1044.                                                   final UnivariateDerivative2 ey) {
  1045.         // Generalization of Kepler equation to circular parameters
  1046.         // with alphaE = PA + E and
  1047.         //      alphaM = PA + M = alphaE - ex.sin(alphaE) + ey.cos(alphaE)
  1048.         UnivariateDerivative2 alphaE        = alphaM;
  1049.         UnivariateDerivative2 shift         = alphaM.getField().getZero();
  1050.         UnivariateDerivative2 alphaEMalphaM = alphaM.getField().getZero();
  1051.         UnivariateDerivative2 cosAlphaE     = alphaE.cos();
  1052.         UnivariateDerivative2 sinAlphaE     = alphaE.sin();
  1053.         int                 iter          = 0;
  1054.         do {
  1055.             final UnivariateDerivative2 f2 = ex.multiply(sinAlphaE).subtract(ey.multiply(cosAlphaE));
  1056.             final UnivariateDerivative2 f1 = alphaM.getField().getOne().subtract(ex.multiply(cosAlphaE)).subtract(ey.multiply(sinAlphaE));
  1057.             final UnivariateDerivative2 f0 = alphaEMalphaM.subtract(f2);

  1058.             final UnivariateDerivative2 f12 = f1.multiply(2);
  1059.             shift = f0.multiply(f12).divide(f1.multiply(f12).subtract(f0.multiply(f2)));

  1060.             alphaEMalphaM  = alphaEMalphaM.subtract(shift);
  1061.             alphaE         = alphaM.add(alphaEMalphaM);
  1062.             cosAlphaE      = alphaE.cos();
  1063.             sinAlphaE      = alphaE.sin();

  1064.         } while (++iter < 50 && FastMath.abs(shift.getValue()) > 1.0e-12);

  1065.         return alphaE;

  1066.     }

  1067.     /** {@inheritDoc} */
  1068.     protected double getMass(final AbsoluteDate date) {
  1069.         return models.get(date).mass;
  1070.     }

  1071. }