FieldEcksteinHechlerPropagator.java

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

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

  20. import org.hipparchus.Field;
  21. import org.hipparchus.CalculusFieldElement;
  22. import org.hipparchus.analysis.differentiation.FieldUnivariateDerivative2;
  23. import org.hipparchus.geometry.euclidean.threed.FieldVector3D;
  24. import org.hipparchus.util.FastMath;
  25. import org.hipparchus.util.FieldSinCos;
  26. import org.hipparchus.util.MathArrays;
  27. import org.hipparchus.util.MathUtils;
  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.FieldCartesianOrbit;
  35. import org.orekit.orbits.FieldCircularOrbit;
  36. import org.orekit.orbits.FieldOrbit;
  37. import org.orekit.orbits.OrbitType;
  38. import org.orekit.orbits.PositionAngle;
  39. import org.orekit.propagation.FieldSpacecraftState;
  40. import org.orekit.propagation.PropagationType;
  41. import org.orekit.time.FieldAbsoluteDate;
  42. import org.orekit.utils.FieldTimeSpanMap;
  43. import org.orekit.utils.ParameterDriver;
  44. import org.orekit.utils.TimeStampedFieldPVCoordinates;

  45. /** This class propagates a {@link org.orekit.propagation.FieldSpacecraftState}
  46.  *  using the analytical Eckstein-Hechler model.
  47.  * <p>The Eckstein-Hechler model is suited for near circular orbits
  48.  * (e &lt; 0.1, with poor accuracy between 0.005 and 0.1) and inclination
  49.  * neither equatorial (direct or retrograde) nor critical (direct or
  50.  * retrograde).</p>
  51.  * @see FieldOrbit
  52.  * @author Guylaine Prat
  53.  */
  54. public class FieldEcksteinHechlerPropagator<T extends CalculusFieldElement<T>> extends FieldAbstractAnalyticalPropagator<T> {

  55.     /** Initial Eckstein-Hechler model. */
  56.     private FieldEHModel<T> initialModel;

  57.     /** All models. */
  58.     private transient FieldTimeSpanMap<FieldEHModel<T>, T> models;

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

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

  63.     /** Un-normalized zonal coefficients. */
  64.     private double[] ck0;

  65.     /** Build a propagator from FieldOrbit and potential provider.
  66.      * <p>Mass and attitude provider are set to unspecified non-null arbitrary values.</p>
  67.      *
  68.      * <p>Using this constructor, an initial osculating orbit is considered.</p>
  69.      *
  70.      * @param initialOrbit initial FieldOrbit
  71.      * @param provider for un-normalized zonal coefficients
  72.      * @see #FieldEcksteinHechlerPropagator(FieldOrbit, AttitudeProvider,
  73.      * UnnormalizedSphericalHarmonicsProvider)
  74.      * @see #FieldEcksteinHechlerPropagator(FieldOrbit, UnnormalizedSphericalHarmonicsProvider,
  75.      * PropagationType)
  76.      */
  77.     public FieldEcksteinHechlerPropagator(final FieldOrbit<T> initialOrbit,
  78.                                           final UnnormalizedSphericalHarmonicsProvider provider) {
  79.         this(initialOrbit, InertialProvider.of(initialOrbit.getFrame()),
  80.              initialOrbit.getA().getField().getZero().add(DEFAULT_MASS), provider,
  81.              provider.onDate(initialOrbit.getDate().toAbsoluteDate()));
  82.     }

  83.     /**
  84.      * Private helper constructor.
  85.      * <p>Using this constructor, an initial osculating orbit is considered.</p>
  86.      * @param initialOrbit initial FieldOrbit
  87.      * @param attitude attitude provider
  88.      * @param mass spacecraft mass
  89.      * @param provider for un-normalized zonal coefficients
  90.      * @param harmonics {@code provider.onDate(initialOrbit.getDate())}
  91.      * @see #FieldEcksteinHechlerPropagator(FieldOrbit, AttitudeProvider, CalculusFieldElement,
  92.      * UnnormalizedSphericalHarmonicsProvider, UnnormalizedSphericalHarmonicsProvider.UnnormalizedSphericalHarmonics, PropagationType)
  93.      */
  94.     public FieldEcksteinHechlerPropagator(final FieldOrbit<T> initialOrbit,
  95.                                           final  AttitudeProvider attitude,
  96.                                           final T mass,
  97.                                           final UnnormalizedSphericalHarmonicsProvider provider,
  98.                                           final UnnormalizedSphericalHarmonicsProvider.UnnormalizedSphericalHarmonics harmonics) {
  99.         this(initialOrbit, attitude,  mass, provider.getAe(), initialOrbit.getA().getField().getZero().add(provider.getMu()),
  100.              harmonics.getUnnormalizedCnm(2, 0),
  101.              harmonics.getUnnormalizedCnm(3, 0),
  102.              harmonics.getUnnormalizedCnm(4, 0),
  103.              harmonics.getUnnormalizedCnm(5, 0),
  104.              harmonics.getUnnormalizedCnm(6, 0));
  105.     }

  106.     /** Build a propagator from FieldOrbit and potential.
  107.      * <p>Mass and attitude provider are set to unspecified non-null arbitrary values.</p>
  108.      * <p>The C<sub>n,0</sub> coefficients are the denormalized zonal coefficients, they
  109.      * are related to both the normalized coefficients
  110.      * <span style="text-decoration: overline">C</span><sub>n,0</sub>
  111.      *  and the J<sub>n</sub> one as follows:
  112.      * <p>
  113.      *   C<sub>n,0</sub> = [(2-δ<sub>0,m</sub>)(2n+1)(n-m)!/(n+m)!]<sup>½</sup>
  114.      *                      <span style="text-decoration: overline">C</span><sub>n,0</sub>
  115.      * <p>
  116.      *   C<sub>n,0</sub> = -J<sub>n</sub>
  117.      *
  118.      * <p>Using this constructor, an initial osculating orbit is considered.</p>
  119.      *
  120.      * @param initialOrbit initial FieldOrbit
  121.      * @param referenceRadius reference radius of the Earth for the potential model (m)
  122.      * @param mu central attraction coefficient (m³/s²)
  123.      * @param c20 un-normalized zonal coefficient (about -1.08e-3 for Earth)
  124.      * @param c30 un-normalized zonal coefficient (about +2.53e-6 for Earth)
  125.      * @param c40 un-normalized zonal coefficient (about +1.62e-6 for Earth)
  126.      * @param c50 un-normalized zonal coefficient (about +2.28e-7 for Earth)
  127.      * @param c60 un-normalized zonal coefficient (about -5.41e-7 for Earth)
  128.      * @see org.orekit.utils.Constants
  129.      * @see #FieldEcksteinHechlerPropagator(FieldOrbit, AttitudeProvider, double,
  130.      * CalculusFieldElement, double, double, double, double, double)
  131.      */
  132.     public FieldEcksteinHechlerPropagator(final FieldOrbit<T> initialOrbit,
  133.                                           final double referenceRadius, final T mu,
  134.                                           final double c20, final double c30, final double c40,
  135.                                           final double c50, final double c60) {
  136.         this(initialOrbit, InertialProvider.of(initialOrbit.getFrame()),
  137.              initialOrbit.getDate().getField().getZero().add(DEFAULT_MASS),
  138.              referenceRadius, mu, c20, c30, c40, c50, c60);
  139.     }

  140.     /** Build a propagator from FieldOrbit, mass and potential provider.
  141.      * <p>Attitude law is set to an unspecified non-null arbitrary value.</p>
  142.      *
  143.      * <p>Using this constructor, an initial osculating orbit is considered.</p>
  144.      *
  145.      * @param initialOrbit initial FieldOrbit
  146.      * @param mass spacecraft mass
  147.      * @param provider for un-normalized zonal coefficients
  148.      * @see #FieldEcksteinHechlerPropagator(FieldOrbit, AttitudeProvider,
  149.      * CalculusFieldElement, UnnormalizedSphericalHarmonicsProvider)
  150.      */
  151.     public FieldEcksteinHechlerPropagator(final FieldOrbit<T> initialOrbit, final T mass,
  152.                                           final UnnormalizedSphericalHarmonicsProvider provider) {
  153.         this(initialOrbit, InertialProvider.of(initialOrbit.getFrame()),
  154.                 mass, provider, provider.onDate(initialOrbit.getDate().toAbsoluteDate()));
  155.     }

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

  189.     /** Build a propagator from FieldOrbit, attitude provider and potential provider.
  190.      * <p>Mass is set to an unspecified non-null arbitrary value.</p>
  191.      * <p>Using this constructor, an initial osculating orbit is considered.</p>
  192.      * @param initialOrbit initial FieldOrbit
  193.      * @param attitudeProv attitude provider
  194.      * @param provider for un-normalized zonal coefficients
  195.      */
  196.     public FieldEcksteinHechlerPropagator(final FieldOrbit<T> initialOrbit,
  197.                                           final AttitudeProvider attitudeProv,
  198.                                           final UnnormalizedSphericalHarmonicsProvider provider) {
  199.         this(initialOrbit, attitudeProv, initialOrbit.getA().getField().getZero().add(DEFAULT_MASS), provider, provider.onDate(initialOrbit.getDate().toAbsoluteDate()));
  200.     }

  201.     /** Build a propagator from FieldOrbit, attitude provider and potential.
  202.      * <p>Mass is set to an unspecified non-null arbitrary value.</p>
  203.      * <p>The C<sub>n,0</sub> coefficients are the denormalized zonal coefficients, they
  204.      * are related to both the normalized coefficients
  205.      * <span style="text-decoration: overline">C</span><sub>n,0</sub>
  206.      *  and the J<sub>n</sub> one as follows:</p>
  207.      * <p>
  208.      *   C<sub>n,0</sub> = [(2-δ<sub>0,m</sub>)(2n+1)(n-m)!/(n+m)!]<sup>½</sup>
  209.      *                     <span style="text-decoration: overline">C</span><sub>n,0</sub>
  210.      * <p>
  211.      *   C<sub>n,0</sub> = -J<sub>n</sub>
  212.      *
  213.      * <p>Using this constructor, an initial osculating orbit is considered.</p>
  214.      *
  215.      * @param initialOrbit initial FieldOrbit
  216.      * @param attitudeProv attitude provider
  217.      * @param referenceRadius reference radius of the Earth for the potential model (m)
  218.      * @param mu central attraction coefficient (m³/s²)
  219.      * @param c20 un-normalized zonal coefficient (about -1.08e-3 for Earth)
  220.      * @param c30 un-normalized zonal coefficient (about +2.53e-6 for Earth)
  221.      * @param c40 un-normalized zonal coefficient (about +1.62e-6 for Earth)
  222.      * @param c50 un-normalized zonal coefficient (about +2.28e-7 for Earth)
  223.      * @param c60 un-normalized zonal coefficient (about -5.41e-7 for Earth)
  224.      */
  225.     public FieldEcksteinHechlerPropagator(final FieldOrbit<T> initialOrbit,
  226.                                           final AttitudeProvider attitudeProv,
  227.                                           final double referenceRadius, final T mu,
  228.                                           final double c20, final double c30, final double c40,
  229.                                           final double c50, final double c60) {
  230.         this(initialOrbit, attitudeProv, initialOrbit.getDate().getField().getZero().add(DEFAULT_MASS),
  231.              referenceRadius, mu, c20, c30, c40, c50, c60);
  232.     }

  233.     /** Build a propagator from FieldOrbit, attitude provider, mass and potential provider.
  234.      * <p>Using this constructor, an initial osculating orbit is considered.</p>
  235.      * @param initialOrbit initial FieldOrbit
  236.      * @param attitudeProv attitude provider
  237.      * @param mass spacecraft mass
  238.      * @param provider for un-normalized zonal coefficients
  239.      * @see #FieldEcksteinHechlerPropagator(FieldOrbit, AttitudeProvider, CalculusFieldElement,
  240.      * UnnormalizedSphericalHarmonicsProvider, PropagationType)
  241.      */
  242.     public FieldEcksteinHechlerPropagator(final FieldOrbit<T> initialOrbit,
  243.                                           final AttitudeProvider attitudeProv,
  244.                                           final T mass,
  245.                                           final UnnormalizedSphericalHarmonicsProvider provider) {
  246.         this(initialOrbit, attitudeProv, mass, provider, provider.onDate(initialOrbit.getDate().toAbsoluteDate()));
  247.     }

  248.     /** Build a propagator from FieldOrbit, attitude provider, mass and potential.
  249.      * <p>The C<sub>n,0</sub> coefficients are the denormalized zonal coefficients, they
  250.      * are related to both the normalized coefficients
  251.      * <span style="text-decoration: overline">C</span><sub>n,0</sub>
  252.      *  and the J<sub>n</sub> one as follows:</p>
  253.      * <p>
  254.      *   C<sub>n,0</sub> = [(2-δ<sub>0,m</sub>)(2n+1)(n-m)!/(n+m)!]<sup>½</sup>
  255.      *                      <span style="text-decoration: overline">C</span><sub>n,0</sub>
  256.      * <p>
  257.      *   C<sub>n,0</sub> = -J<sub>n</sub>
  258.      *
  259.      * <p>Using this constructor, an initial osculating orbit is considered.</p>
  260.      *
  261.      * @param initialOrbit initial FieldOrbit
  262.      * @param attitudeProv attitude provider
  263.      * @param mass spacecraft mass
  264.      * @param referenceRadius reference radius of the Earth for the potential model (m)
  265.      * @param mu central attraction coefficient (m³/s²)
  266.      * @param c20 un-normalized zonal coefficient (about -1.08e-3 for Earth)
  267.      * @param c30 un-normalized zonal coefficient (about +2.53e-6 for Earth)
  268.      * @param c40 un-normalized zonal coefficient (about +1.62e-6 for Earth)
  269.      * @param c50 un-normalized zonal coefficient (about +2.28e-7 for Earth)
  270.      * @param c60 un-normalized zonal coefficient (about -5.41e-7 for Earth)
  271.      * @see #FieldEcksteinHechlerPropagator(FieldOrbit, AttitudeProvider, CalculusFieldElement, double,
  272.      *                                      CalculusFieldElement, double, double, double, double, double, PropagationType)
  273.      */
  274.     public FieldEcksteinHechlerPropagator(final FieldOrbit<T> initialOrbit,
  275.                                           final AttitudeProvider attitudeProv,
  276.                                           final T mass,
  277.                                           final double referenceRadius, final T mu,
  278.                                           final double c20, final double c30, final double c40,
  279.                                           final double c50, final double c60) {
  280.         this(initialOrbit, attitudeProv, mass, referenceRadius, mu, c20, c30, c40, c50, c60, PropagationType.OSCULATING);
  281.     }

  282.     /** Build a propagator from orbit and potential provider.
  283.      * <p>Mass and attitude provider are set to unspecified non-null arbitrary values.</p>
  284.      *
  285.      * <p>Using this constructor, it is possible to define the initial orbit as
  286.      * a mean Eckstein-Hechler orbit or an osculating one.</p>
  287.      *
  288.      * @param initialOrbit initial orbit
  289.      * @param provider for un-normalized zonal coefficients
  290.      * @param initialType initial orbit type (mean Eckstein-Hechler orbit or osculating orbit)
  291.      * @since 10.2
  292.      */
  293.     public FieldEcksteinHechlerPropagator(final FieldOrbit<T> initialOrbit,
  294.                                           final UnnormalizedSphericalHarmonicsProvider provider,
  295.                                           final PropagationType initialType) {
  296.         this(initialOrbit, InertialProvider.of(initialOrbit.getFrame()),
  297.              initialOrbit.getA().getField().getZero().add(DEFAULT_MASS), provider,
  298.              provider.onDate(initialOrbit.getDate().toAbsoluteDate()), initialType);
  299.     }

  300.     /** Build a propagator from orbit, attitude provider, mass and potential provider.
  301.      * <p>Using this constructor, it is possible to define the initial orbit as
  302.      * a mean Eckstein-Hechler orbit or an osculating one.</p>
  303.      * @param initialOrbit initial orbit
  304.      * @param attitudeProv attitude provider
  305.      * @param mass spacecraft mass
  306.      * @param provider for un-normalized zonal coefficients
  307.      * @param initialType initial orbit type (mean Eckstein-Hechler orbit or osculating orbit)
  308.      * @since 10.2
  309.      */
  310.     public FieldEcksteinHechlerPropagator(final FieldOrbit<T> initialOrbit,
  311.                                           final AttitudeProvider attitudeProv,
  312.                                           final T mass,
  313.                                           final UnnormalizedSphericalHarmonicsProvider provider,
  314.                                           final PropagationType initialType) {
  315.         this(initialOrbit, attitudeProv, mass, provider, provider.onDate(initialOrbit.getDate().toAbsoluteDate()), initialType);
  316.     }

  317.     /**
  318.      * Private helper constructor.
  319.      * <p>Using this constructor, it is possible to define the initial orbit as
  320.      * a mean Eckstein-Hechler orbit or an osculating one.</p>
  321.      * @param initialOrbit initial orbit
  322.      * @param attitude attitude provider
  323.      * @param mass spacecraft mass
  324.      * @param provider for un-normalized zonal coefficients
  325.      * @param harmonics {@code provider.onDate(initialOrbit.getDate())}
  326.      * @param initialType initial orbit type (mean Eckstein-Hechler orbit or osculating orbit)
  327.      * @since 10.2
  328.      */
  329.     public FieldEcksteinHechlerPropagator(final FieldOrbit<T> initialOrbit,
  330.                                           final AttitudeProvider attitude,
  331.                                           final T mass,
  332.                                           final UnnormalizedSphericalHarmonicsProvider provider,
  333.                                           final UnnormalizedSphericalHarmonics harmonics,
  334.                                           final PropagationType initialType) {
  335.         this(initialOrbit, attitude, mass, provider.getAe(), initialOrbit.getA().getField().getZero().add(provider.getMu()),
  336.              harmonics.getUnnormalizedCnm(2, 0),
  337.              harmonics.getUnnormalizedCnm(3, 0),
  338.              harmonics.getUnnormalizedCnm(4, 0),
  339.              harmonics.getUnnormalizedCnm(5, 0),
  340.              harmonics.getUnnormalizedCnm(6, 0),
  341.              initialType);
  342.     }

  343.     /** Build a propagator from FieldOrbit, attitude provider, mass and potential.
  344.      * <p>The C<sub>n,0</sub> coefficients are the denormalized zonal coefficients, they
  345.      * are related to both the normalized coefficients
  346.      * <span style="text-decoration: overline">C</span><sub>n,0</sub>
  347.      *  and the J<sub>n</sub> one as follows:</p>
  348.      * <p>
  349.      *   C<sub>n,0</sub> = [(2-δ<sub>0,m</sub>)(2n+1)(n-m)!/(n+m)!]<sup>½</sup>
  350.      *                      <span style="text-decoration: overline">C</span><sub>n,0</sub>
  351.      * <p>
  352.      *   C<sub>n,0</sub> = -J<sub>n</sub>
  353.      *
  354.      * <p>Using this constructor, it is possible to define the initial orbit as
  355.      * a mean Eckstein-Hechler orbit or an osculating one.</p>
  356.      *
  357.      * @param initialOrbit initial FieldOrbit
  358.      * @param attitudeProv attitude provider
  359.      * @param mass spacecraft mass
  360.      * @param referenceRadius reference radius of the Earth for the potential model (m)
  361.      * @param mu central attraction coefficient (m³/s²)
  362.      * @param c20 un-normalized zonal coefficient (about -1.08e-3 for Earth)
  363.      * @param c30 un-normalized zonal coefficient (about +2.53e-6 for Earth)
  364.      * @param c40 un-normalized zonal coefficient (about +1.62e-6 for Earth)
  365.      * @param c50 un-normalized zonal coefficient (about +2.28e-7 for Earth)
  366.      * @param c60 un-normalized zonal coefficient (about -5.41e-7 for Earth)
  367.      * @param initialType initial orbit type (mean Eckstein-Hechler orbit or osculating orbit)
  368.      * @since 10.2
  369.      */
  370.     public FieldEcksteinHechlerPropagator(final FieldOrbit<T> initialOrbit,
  371.                                           final AttitudeProvider attitudeProv,
  372.                                           final T mass,
  373.                                           final double referenceRadius, final T mu,
  374.                                           final double c20, final double c30, final double c40,
  375.                                           final double c50, final double c60,
  376.                                           final PropagationType initialType) {

  377.         super(mass.getField(), attitudeProv);
  378.         try {

  379.             // store model coefficients
  380.             this.referenceRadius = referenceRadius;
  381.             this.mu  = mu;
  382.             this.ck0 = new double[] {
  383.                 0.0, 0.0, c20, c30, c40, c50, c60
  384.             };

  385.             // compute mean parameters if needed
  386.             // transform into circular adapted parameters used by the Eckstein-Hechler model
  387.             resetInitialState(new FieldSpacecraftState<>(initialOrbit,
  388.                                                          attitudeProv.getAttitude(initialOrbit,
  389.                                                                                   initialOrbit.getDate(),
  390.                                                                                   initialOrbit.getFrame()),
  391.                                                          mass),
  392.                               initialType);

  393.         } catch (OrekitException oe) {
  394.             throw new OrekitException(oe);
  395.         }
  396.     }

  397.     /** {@inheritDoc}
  398.      * <p>The new initial state to consider
  399.      * must be defined with an osculating orbit.</p>
  400.      * @see #resetInitialState(FieldSpacecraftState, PropagationType)
  401.      */
  402.     @Override
  403.     public void resetInitialState(final FieldSpacecraftState<T> state) {
  404.         resetInitialState(state, PropagationType.OSCULATING);
  405.     }

  406.     /** Reset the propagator initial state.
  407.      * @param state new initial state to consider
  408.      * @param stateType mean Eckstein-Hechler orbit or osculating orbit
  409.      * @since 10.2
  410.      */
  411.     public void resetInitialState(final FieldSpacecraftState<T> state, final PropagationType stateType) {
  412.         super.resetInitialState(state);
  413.         final FieldCircularOrbit<T> circular = (FieldCircularOrbit<T>) OrbitType.CIRCULAR.convertType(state.getOrbit());
  414.         this.initialModel = (stateType == PropagationType.MEAN) ?
  415.                              new FieldEHModel<>(circular, state.getMass(), referenceRadius, mu, ck0) :
  416.                              computeMeanParameters(circular, state.getMass());
  417.         this.models = new FieldTimeSpanMap<FieldEHModel<T>, T>(initialModel, state.getA().getField());
  418.     }

  419.     /** {@inheritDoc} */
  420.     @Override
  421.     protected void resetIntermediateState(final FieldSpacecraftState<T> state, final boolean forward) {
  422.         final FieldEHModel<T> newModel = computeMeanParameters((FieldCircularOrbit<T>) OrbitType.CIRCULAR.convertType(state.getOrbit()),
  423.                                                        state.getMass());
  424.         if (forward) {
  425.             models.addValidAfter(newModel, state.getDate());
  426.         } else {
  427.             models.addValidBefore(newModel, state.getDate());
  428.         }
  429.         stateChanged(state);
  430.     }

  431.     /** Compute mean parameters according to the Eckstein-Hechler analytical model.
  432.      * @param osculating osculating FieldOrbit
  433.      * @param mass constant mass
  434.      * @return Eckstein-Hechler mean model
  435.      */
  436.     private FieldEHModel<T> computeMeanParameters(final FieldCircularOrbit<T> osculating, final T mass) {

  437.         // sanity check
  438.         if (osculating.getA().getReal() < referenceRadius) {
  439.             throw new OrekitException(OrekitMessages.TRAJECTORY_INSIDE_BRILLOUIN_SPHERE,
  440.                                            osculating.getA());
  441.         }
  442.         final Field<T> field = mass.getField();
  443.         final T one = field.getOne();
  444.         final T zero = field.getZero();
  445.         // rough initialization of the mean parameters
  446.         FieldEHModel<T> current = new FieldEHModel<>(osculating, mass, referenceRadius, mu, ck0);
  447.         // threshold for each parameter
  448.         final T epsilon         = one.multiply(1.0e-13);
  449.         final T thresholdA      = epsilon.multiply(current.mean.getA().abs().add(1.0));
  450.         final T thresholdE      = epsilon.multiply(current.mean.getE().add(1.0));
  451.         final T thresholdAngles = epsilon.multiply(one.getPi());


  452.         int i = 0;
  453.         while (i++ < 100) {

  454.             // recompute the osculating parameters from the current mean parameters
  455.             final FieldUnivariateDerivative2<T>[] parameters = current.propagateParameters(current.mean.getDate());
  456.             // adapted parameters residuals
  457.             final T deltaA      = osculating.getA()         .subtract(parameters[0].getValue());
  458.             final T deltaEx     = osculating.getCircularEx().subtract(parameters[1].getValue());
  459.             final T deltaEy     = osculating.getCircularEy().subtract(parameters[2].getValue());
  460.             final T deltaI      = osculating.getI()         .subtract(parameters[3].getValue());
  461.             final T deltaRAAN   = MathUtils.normalizeAngle(osculating.getRightAscensionOfAscendingNode().subtract(
  462.                                                                 parameters[4].getValue()),
  463.                                                                 zero);
  464.             final T deltaAlphaM = MathUtils.normalizeAngle(osculating.getAlphaM().subtract(parameters[5].getValue()), zero);
  465.             // update mean parameters
  466.             current = new FieldEHModel<>(new FieldCircularOrbit<>(current.mean.getA().add(deltaA),
  467.                                                                   current.mean.getCircularEx().add( deltaEx),
  468.                                                                   current.mean.getCircularEy().add( deltaEy),
  469.                                                                   current.mean.getI()         .add( deltaI ),
  470.                                                                   current.mean.getRightAscensionOfAscendingNode().add(deltaRAAN),
  471.                                                                   current.mean.getAlphaM().add(deltaAlphaM),
  472.                                                                   PositionAngle.MEAN,
  473.                                                                   current.mean.getFrame(),
  474.                                                                   current.mean.getDate(), mu),
  475.                                   mass, referenceRadius, mu, ck0);
  476.             // check convergence
  477.             if (FastMath.abs(deltaA.getReal())      < thresholdA.getReal() &&
  478.                 FastMath.abs(deltaEx.getReal())     < thresholdE.getReal() &&
  479.                 FastMath.abs(deltaEy.getReal())     < thresholdE.getReal() &&
  480.                 FastMath.abs(deltaI.getReal())      < thresholdAngles.getReal() &&
  481.                 FastMath.abs(deltaRAAN.getReal())   < thresholdAngles.getReal() &&
  482.                 FastMath.abs(deltaAlphaM.getReal()) < thresholdAngles.getReal()) {
  483.                 return current;
  484.             }

  485.         }

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

  487.     }

  488.     /** {@inheritDoc} */
  489.     @Override
  490.     public FieldCartesianOrbit<T> propagateOrbit(final FieldAbsoluteDate<T> date, final T[] parameters) {
  491.         // compute Cartesian parameters, taking derivatives into account
  492.         // to make sure velocity and acceleration are consistent
  493.         final FieldEHModel<T> current = models.get(date);
  494.         return new FieldCartesianOrbit<>(toCartesian(date, current.propagateParameters(date)),
  495.                                          current.mean.getFrame(), mu);
  496.     }

  497.     /** Local class for Eckstein-Hechler model, with fixed mean parameters. */
  498.     private static class FieldEHModel<T extends CalculusFieldElement<T>> {

  499.         /** Mean FieldOrbit. */
  500.         private final FieldCircularOrbit<T> mean;

  501.         /** Constant mass. */
  502.         private final T mass;
  503.         // CHECKSTYLE: stop JavadocVariable check

  504.         // preprocessed values
  505.         private final T xnotDot;
  506.         private final T rdpom;
  507.         private final T rdpomp;
  508.         private final T eps1;
  509.         private final T eps2;
  510.         private final T xim;
  511.         private final T ommD;
  512.         private final T rdl;
  513.         private final T aMD;

  514.         private final T kh;
  515.         private final T kl;

  516.         private final T ax1;
  517.         private final T ay1;
  518.         private final T as1;
  519.         private final T ac2;
  520.         private final T axy3;
  521.         private final T as3;
  522.         private final T ac4;
  523.         private final T as5;
  524.         private final T ac6;

  525.         private final T ex1;
  526.         private final T exx2;
  527.         private final T exy2;
  528.         private final T ex3;
  529.         private final T ex4;

  530.         private final T ey1;
  531.         private final T eyx2;
  532.         private final T eyy2;
  533.         private final T ey3;
  534.         private final T ey4;

  535.         private final T rx1;
  536.         private final T ry1;
  537.         private final T r2;
  538.         private final T r3;
  539.         private final T rl;

  540.         private final T iy1;
  541.         private final T ix1;
  542.         private final T i2;
  543.         private final T i3;
  544.         private final T ih;

  545.         private final T lx1;
  546.         private final T ly1;
  547.         private final T l2;
  548.         private final T l3;
  549.         private final T ll;

  550.         // CHECKSTYLE: resume JavadocVariable check

  551.         /** Create a model for specified mean FieldOrbit.
  552.          * @param mean mean FieldOrbit
  553.          * @param mass constant mass
  554.          * @param referenceRadius reference radius of the central body attraction model (m)
  555.          * @param mu central attraction coefficient (m³/s²)
  556.          * @param ck0 un-normalized zonal coefficients
  557.          */
  558.         FieldEHModel(final FieldCircularOrbit<T> mean, final T mass,
  559.                      final double referenceRadius, final T mu, final double[] ck0) {

  560.             this.mean            = mean;
  561.             this.mass            = mass;
  562.             final T zero = mass.getField().getZero();
  563.             final T one  = mass.getField().getOne();
  564.             // preliminary processing
  565.             T q =  zero.add(referenceRadius).divide(mean.getA());
  566.             T ql = q.multiply(q);
  567.             final T g2 = ql.multiply(ck0[2]);
  568.             ql = ql.multiply(q);
  569.             final T g3 = ql.multiply(ck0[3]);
  570.             ql = ql.multiply(q);
  571.             final T g4 = ql.multiply(ck0[4]);
  572.             ql = ql.multiply(q);
  573.             final T g5 = ql.multiply(ck0[5]);
  574.             ql = ql.multiply(q);
  575.             final T g6 = ql.multiply(ck0[6]);

  576.             final FieldSinCos<T> sc = FastMath.sinCos(mean.getI());
  577.             final T cosI1 = sc.cos();
  578.             final T sinI1 = sc.sin();
  579.             final T sinI2 = sinI1.multiply(sinI1);
  580.             final T sinI4 = sinI2.multiply(sinI2);
  581.             final T sinI6 = sinI2.multiply(sinI4);

  582.             if (sinI2.getReal() < 1.0e-10) {
  583.                 throw new OrekitException(OrekitMessages.ALMOST_EQUATORIAL_ORBIT,
  584.                                                FastMath.toDegrees(mean.getI().getReal()));
  585.             }

  586.             if (FastMath.abs(sinI2.getReal() - 4.0 / 5.0) < 1.0e-3) {
  587.                 throw new OrekitException(OrekitMessages.ALMOST_CRITICALLY_INCLINED_ORBIT,
  588.                                                FastMath.toDegrees(mean.getI().getReal()));
  589.             }

  590.             if (mean.getE().getReal() > 0.1) {
  591.                 // if 0.005 < e < 0.1 no error is triggered, but accuracy is poor
  592.                 throw new OrekitException(OrekitMessages.TOO_LARGE_ECCENTRICITY_FOR_PROPAGATION_MODEL,
  593.                                                mean.getE());
  594.             }

  595.             xnotDot = mu.divide(mean.getA()).sqrt().divide(mean.getA());

  596.             rdpom = g2.multiply(-0.75).multiply(sinI2.multiply(-5.0).add(4.0));
  597.             rdpomp = g4.multiply(7.5).multiply(sinI2.multiply(-31.0 / 8.0).add(1.0).add( sinI4.multiply(49.0 / 16.0))).subtract(
  598.                     g6.multiply(13.125).multiply(one.subtract(sinI2.multiply(8.0)).add(sinI4.multiply(129.0 / 8.0)).subtract(sinI6.multiply(297.0 / 32.0)) ));


  599.             q = zero.add(3.0).divide(rdpom.multiply(32.0));
  600.             eps1 = q.multiply(g4).multiply(sinI2).multiply(sinI2.multiply(-35.0).add(30.0)).subtract(
  601.                    q.multiply(175.0).multiply(g6).multiply(sinI2).multiply(sinI2.multiply(-3.0).add(sinI4.multiply(2.0625)).add(1.0)));
  602.             q = sinI1.multiply(3.0).divide(rdpom.multiply(8.0));
  603.             eps2 = q.multiply(g3).multiply(sinI2.multiply(-5.0).add(4.0)).subtract(q.multiply(g5).multiply(sinI2.multiply(-35.0).add(sinI4.multiply(26.25)).add(10.0)));

  604.             xim = mean.getI();
  605.             ommD = cosI1.multiply(g2.multiply(1.50).subtract(g2.multiply(2.25).multiply(g2).multiply(sinI2.multiply(-19.0 / 6.0).add(2.5))).add(
  606.                             g4.multiply(0.9375).multiply(sinI2.multiply(7.0).subtract(4.0))).add(
  607.                             g6.multiply(3.28125).multiply(sinI2.multiply(-9.0).add(2.0).add(sinI4.multiply(8.25)))));

  608.             rdl = g2.multiply(-1.50).multiply(sinI2.multiply(-4.0).add(3.0)).add(1.0);
  609.             aMD = rdl.add(
  610.                     g2.multiply(2.25).multiply(g2.multiply(sinI2.multiply(-263.0 / 12.0 ).add(9.0).add(sinI4.multiply(341.0 / 24.0))))).add(
  611.                     g4.multiply(15.0 / 16.0).multiply(sinI2.multiply(-31.0).add(8.0).add(sinI4.multiply(24.5)))).add(
  612.                     g6.multiply(105.0 / 32.0).multiply(sinI2.multiply(25.0).add(-10.0 / 3.0).subtract(sinI4.multiply(48.75)).add(sinI6.multiply(27.5))));

  613.             final T qq   = g2.divide(rdl).multiply(-1.5);
  614.             final T qA   = g2.multiply(0.75).multiply(g2).multiply(sinI2);
  615.             final T qB   = g4.multiply(0.25).multiply(sinI2);
  616.             final T qC   = g6.multiply(105.0 / 16.0).multiply(sinI2);
  617.             final T qD   = g3.multiply(-0.75).multiply(sinI1);
  618.             final T qE   = g5.multiply(3.75).multiply(sinI1);
  619.             kh = zero.add(0.375).divide(rdpom);
  620.             kl = kh.divide(sinI1);

  621.             ax1 = qq.multiply(sinI2.multiply(-3.5).add(2.0));
  622.             ay1 = qq.multiply(sinI2.multiply(-2.5).add(2.0));
  623.             as1 = qD.multiply(sinI2.multiply(-5.0).add(4.0)).add(
  624.                   qE.multiply(sinI4.multiply(2.625).add(sinI2.multiply(-3.5)).add(1.0)));
  625.             ac2 = qq.multiply(sinI2).add(
  626.                   qA.multiply(7.0).multiply(sinI2.multiply(-3.0).add(2.0))).add(
  627.                   qB.multiply(sinI2.multiply(-17.5).add(15.0))).add(
  628.                   qC.multiply(sinI2.multiply(3.0).subtract(1.0).subtract(sinI4.multiply(33.0 / 16.0))));
  629.             axy3 = qq.multiply(3.5).multiply(sinI2);
  630.             as3 = qD.multiply(5.0 / 3.0).multiply(sinI2).add(
  631.                   qE.multiply(7.0 / 6.0).multiply(sinI2).multiply(sinI2.multiply(-1.125).add(1)));
  632.             ac4 = qA.multiply(sinI2).add(
  633.                   qB.multiply(4.375).multiply(sinI2)).add(
  634.                   qC.multiply(0.75).multiply(sinI4.multiply(1.1).subtract(sinI2)));

  635.             as5 = qE.multiply(21.0 / 80.0).multiply(sinI4);

  636.             ac6 = qC.multiply(-11.0 / 80.0).multiply(sinI4);

  637.             ex1 = qq.multiply(sinI2.multiply(-1.25).add(1.0));
  638.             exx2 = qq.multiply(0.5).multiply(sinI2.multiply(-5.0).add(3.0));
  639.             exy2 = qq.multiply(sinI2.multiply(-1.5).add(2.0));
  640.             ex3 = qq.multiply(7.0 / 12.0).multiply(sinI2);
  641.             ex4 = qq.multiply(17.0 / 8.0).multiply(sinI2);

  642.             ey1 = qq.multiply(sinI2.multiply(-1.75).add(1.0));
  643.             eyx2 = qq.multiply(sinI2.multiply(-3.0).add(1.0));
  644.             eyy2 = qq.multiply(sinI2.multiply(2.0).subtract(1.5));
  645.             ey3 = qq.multiply(7.0 / 12.0).multiply(sinI2);
  646.             ey4 = qq.multiply(17.0 / 8.0).multiply(sinI2);

  647.             q  = cosI1.multiply(qq).negate();
  648.             rx1 = q.multiply(3.5);
  649.             ry1 = q.multiply(-2.5);
  650.             r2 = q.multiply(-0.5);
  651.             r3 =  q.multiply(7.0 / 6.0);
  652.             rl = g3 .multiply( cosI1).multiply(sinI2.multiply(-15.0).add(4.0)).subtract(
  653.                  g5.multiply(2.5).multiply(cosI1).multiply(sinI2.multiply(-42.0).add(4.0).add(sinI4.multiply(52.5))));

  654.             q = qq.multiply(0.5).multiply(sinI1).multiply(cosI1);
  655.             iy1 =  q;
  656.             ix1 = q.negate();
  657.             i2 =  q;
  658.             i3 =  q.multiply(7.0 / 3.0);
  659.             ih = g3.negate().multiply(cosI1).multiply(sinI2.multiply(-5.0).add(4)).add(
  660.                  g5.multiply(2.5).multiply(cosI1).multiply(sinI2.multiply(-14.0).add(4.0).add(sinI4.multiply(10.5))));
  661.             lx1 = qq.multiply(sinI2.multiply(-77.0 / 8.0).add(7.0));
  662.             ly1 = qq.multiply(sinI2.multiply(55.0 / 8.0).subtract(7.50));
  663.             l2 = qq.multiply(sinI2.multiply(1.25).subtract(0.5));
  664.             l3 = qq.multiply(sinI2.multiply(77.0 / 24.0).subtract(7.0 / 6.0));
  665.             ll = g3.multiply(sinI2.multiply(53.0).subtract(4.0).add(sinI4.multiply(-57.5))).add(
  666.                  g5.multiply(2.5).multiply(sinI2.multiply(-96.0).add(4.0).add(sinI4.multiply(269.5).subtract(sinI6.multiply(183.75)))));

  667.         }

  668.         /** Extrapolate a FieldOrbit up to a specific target date.
  669.          * @param date target date for the FieldOrbit
  670.          * @return propagated parameters
  671.          */
  672.         public FieldUnivariateDerivative2<T>[] propagateParameters(final FieldAbsoluteDate<T> date) {
  673.             final Field<T> field = date.durationFrom(mean.getDate()).getField();
  674.             final T one = field.getOne();
  675.             final T zero = field.getZero();
  676.             // Keplerian evolution
  677.             final FieldUnivariateDerivative2<T> dt = new FieldUnivariateDerivative2<>(date.durationFrom(mean.getDate()), one, zero);
  678.             final FieldUnivariateDerivative2<T> xnot = dt.multiply(xnotDot);

  679.             // secular effects

  680.             // eccentricity
  681.             final FieldUnivariateDerivative2<T> x   = xnot.multiply(rdpom.add(rdpomp));
  682.             final FieldUnivariateDerivative2<T> cx  = x.cos();
  683.             final FieldUnivariateDerivative2<T> sx  = x.sin();
  684.             final FieldUnivariateDerivative2<T> exm = cx.multiply(mean.getCircularEx()).
  685.                                             add(sx.multiply(eps2.subtract(one.subtract(eps1).multiply(mean.getCircularEy()))));
  686.             final FieldUnivariateDerivative2<T> eym = sx.multiply(eps1.add(1.0).multiply(mean.getCircularEx())).
  687.                                             add(cx.multiply(mean.getCircularEy().subtract(eps2))).
  688.                                             add(eps2);
  689.             // no secular effect on inclination

  690.             // right ascension of ascending node
  691.             final FieldUnivariateDerivative2<T> omm =
  692.                            new FieldUnivariateDerivative2<>(MathUtils.normalizeAngle(mean.getRightAscensionOfAscendingNode().add(ommD.multiply(xnot.getValue())),
  693.                                                                                      one.getPi()),
  694.                                                             ommD.multiply(xnotDot),
  695.                                                             zero);
  696.             // latitude argument
  697.             final FieldUnivariateDerivative2<T> xlm =
  698.                             new FieldUnivariateDerivative2<>(MathUtils.normalizeAngle(mean.getAlphaM().add(aMD.multiply(xnot.getValue())),
  699.                                                                                       one.getPi()),
  700.                                                            aMD.multiply(xnotDot),
  701.                                                            zero);

  702.             // periodical terms
  703.             final FieldUnivariateDerivative2<T> cl1 = xlm.cos();
  704.             final FieldUnivariateDerivative2<T> sl1 = xlm.sin();
  705.             final FieldUnivariateDerivative2<T> cl2 = cl1.multiply(cl1).subtract(sl1.multiply(sl1));
  706.             final FieldUnivariateDerivative2<T> sl2 = cl1.multiply(sl1).add(sl1.multiply(cl1));
  707.             final FieldUnivariateDerivative2<T> cl3 = cl2.multiply(cl1).subtract(sl2.multiply(sl1));
  708.             final FieldUnivariateDerivative2<T> sl3 = cl2.multiply(sl1).add(sl2.multiply(cl1));
  709.             final FieldUnivariateDerivative2<T> cl4 = cl3.multiply(cl1).subtract(sl3.multiply(sl1));
  710.             final FieldUnivariateDerivative2<T> sl4 = cl3.multiply(sl1).add(sl3.multiply(cl1));
  711.             final FieldUnivariateDerivative2<T> cl5 = cl4.multiply(cl1).subtract(sl4.multiply(sl1));
  712.             final FieldUnivariateDerivative2<T> sl5 = cl4.multiply(sl1).add(sl4.multiply(cl1));
  713.             final FieldUnivariateDerivative2<T> cl6 = cl5.multiply(cl1).subtract(sl5.multiply(sl1));

  714.             final FieldUnivariateDerivative2<T> qh  = eym.subtract(eps2).multiply(kh);
  715.             final FieldUnivariateDerivative2<T> ql  = exm.multiply(kl);

  716.             final FieldUnivariateDerivative2<T> exmCl1 = exm.multiply(cl1);
  717.             final FieldUnivariateDerivative2<T> exmSl1 = exm.multiply(sl1);
  718.             final FieldUnivariateDerivative2<T> eymCl1 = eym.multiply(cl1);
  719.             final FieldUnivariateDerivative2<T> eymSl1 = eym.multiply(sl1);
  720.             final FieldUnivariateDerivative2<T> exmCl2 = exm.multiply(cl2);
  721.             final FieldUnivariateDerivative2<T> exmSl2 = exm.multiply(sl2);
  722.             final FieldUnivariateDerivative2<T> eymCl2 = eym.multiply(cl2);
  723.             final FieldUnivariateDerivative2<T> eymSl2 = eym.multiply(sl2);
  724.             final FieldUnivariateDerivative2<T> exmCl3 = exm.multiply(cl3);
  725.             final FieldUnivariateDerivative2<T> exmSl3 = exm.multiply(sl3);
  726.             final FieldUnivariateDerivative2<T> eymCl3 = eym.multiply(cl3);
  727.             final FieldUnivariateDerivative2<T> eymSl3 = eym.multiply(sl3);
  728.             final FieldUnivariateDerivative2<T> exmCl4 = exm.multiply(cl4);
  729.             final FieldUnivariateDerivative2<T> exmSl4 = exm.multiply(sl4);
  730.             final FieldUnivariateDerivative2<T> eymCl4 = eym.multiply(cl4);
  731.             final FieldUnivariateDerivative2<T> eymSl4 = eym.multiply(sl4);

  732.             // semi major axis
  733.             final FieldUnivariateDerivative2<T> rda = exmCl1.multiply(ax1).
  734.                                             add(eymSl1.multiply(ay1)).
  735.                                             add(sl1.multiply(as1)).
  736.                                             add(cl2.multiply(ac2)).
  737.                                             add(exmCl3.add(eymSl3).multiply(axy3)).
  738.                                             add(sl3.multiply(as3)).
  739.                                             add(cl4.multiply(ac4)).
  740.                                             add(sl5.multiply(as5)).
  741.                                             add(cl6.multiply(ac6));

  742.             // eccentricity
  743.             final FieldUnivariateDerivative2<T> rdex = cl1.multiply(ex1).
  744.                                              add(exmCl2.multiply(exx2)).
  745.                                              add(eymSl2.multiply(exy2)).
  746.                                              add(cl3.multiply(ex3)).
  747.                                              add(exmCl4.add(eymSl4).multiply(ex4));
  748.             final FieldUnivariateDerivative2<T> rdey = sl1.multiply(ey1).
  749.                                              add(exmSl2.multiply(eyx2)).
  750.                                              add(eymCl2.multiply(eyy2)).
  751.                                              add(sl3.multiply(ey3)).
  752.                                              add(exmSl4.subtract(eymCl4).multiply(ey4));

  753.             // ascending node
  754.             final FieldUnivariateDerivative2<T> rdom = exmSl1.multiply(rx1).
  755.                                              add(eymCl1.multiply(ry1)).
  756.                                              add(sl2.multiply(r2)).
  757.                                              add(eymCl3.subtract(exmSl3).multiply(r3)).
  758.                                              add(ql.multiply(rl));

  759.             // inclination
  760.             final FieldUnivariateDerivative2<T> rdxi = eymSl1.multiply(iy1).
  761.                                              add(exmCl1.multiply(ix1)).
  762.                                              add(cl2.multiply(i2)).
  763.                                              add(exmCl3.add(eymSl3).multiply(i3)).
  764.                                              add(qh.multiply(ih));

  765.             // latitude argument
  766.             final FieldUnivariateDerivative2<T> rdxl = exmSl1.multiply(lx1).
  767.                                              add(eymCl1.multiply(ly1)).
  768.                                              add(sl2.multiply(l2)).
  769.                                              add(exmSl3.subtract(eymCl3).multiply(l3)).
  770.                                              add(ql.multiply(ll));
  771.             // osculating parameters
  772.             final FieldUnivariateDerivative2<T>[] FTD = MathArrays.buildArray(rdxl.getField(), 6);

  773.             FTD[0] = rda.add(1.0).multiply(mean.getA());
  774.             FTD[1] = rdex.add(exm);
  775.             FTD[2] = rdey.add(eym);
  776.             FTD[3] = rdxi.add(xim);
  777.             FTD[4] = rdom.add(omm);
  778.             FTD[5] = rdxl.add(xlm);
  779.             return FTD;

  780.         }

  781.     }

  782.     /** Convert circular parameters <em>with derivatives</em> to Cartesian coordinates.
  783.      * @param date date of the parameters
  784.      * @param parameters circular parameters (a, ex, ey, i, raan, alphaM)
  785.      * @return Cartesian coordinates consistent with values and derivatives
  786.      */
  787.     private TimeStampedFieldPVCoordinates<T> toCartesian(final FieldAbsoluteDate<T> date, final FieldUnivariateDerivative2<T>[] parameters) {

  788.         // evaluate coordinates in the FieldOrbit canonical reference frame
  789.         final FieldUnivariateDerivative2<T> cosOmega = parameters[4].cos();
  790.         final FieldUnivariateDerivative2<T> sinOmega = parameters[4].sin();
  791.         final FieldUnivariateDerivative2<T> cosI     = parameters[3].cos();
  792.         final FieldUnivariateDerivative2<T> sinI     = parameters[3].sin();
  793.         final FieldUnivariateDerivative2<T> alphaE   = meanToEccentric(parameters[5], parameters[1], parameters[2]);
  794.         final FieldUnivariateDerivative2<T> cosAE    = alphaE.cos();
  795.         final FieldUnivariateDerivative2<T> sinAE    = alphaE.sin();
  796.         final FieldUnivariateDerivative2<T> ex2      = parameters[1].multiply(parameters[1]);
  797.         final FieldUnivariateDerivative2<T> ey2      = parameters[2].multiply(parameters[2]);
  798.         final FieldUnivariateDerivative2<T> exy      = parameters[1].multiply(parameters[2]);
  799.         final FieldUnivariateDerivative2<T> q        = ex2.add(ey2).subtract(1).negate().sqrt();
  800.         final FieldUnivariateDerivative2<T> beta     = q.add(1).reciprocal();
  801.         final FieldUnivariateDerivative2<T> bx2      = beta.multiply(ex2);
  802.         final FieldUnivariateDerivative2<T> by2      = beta.multiply(ey2);
  803.         final FieldUnivariateDerivative2<T> bxy      = beta.multiply(exy);
  804.         final FieldUnivariateDerivative2<T> u        = bxy.multiply(sinAE).subtract(parameters[1].add(by2.subtract(1).multiply(cosAE)));
  805.         final FieldUnivariateDerivative2<T> v        = bxy.multiply(cosAE).subtract(parameters[2].add(bx2.subtract(1).multiply(sinAE)));
  806.         final FieldUnivariateDerivative2<T> x        = parameters[0].multiply(u);
  807.         final FieldUnivariateDerivative2<T> y        = parameters[0].multiply(v);

  808.         // canonical FieldOrbit reference frame
  809.         final FieldVector3D<FieldUnivariateDerivative2<T>> p =
  810.                 new FieldVector3D<>(x.multiply(cosOmega).subtract(y.multiply(cosI.multiply(sinOmega))),
  811.                                     x.multiply(sinOmega).add(y.multiply(cosI.multiply(cosOmega))),
  812.                                     y.multiply(sinI));

  813.         // dispatch derivatives
  814.         final FieldVector3D<T> p0 = new FieldVector3D<>(p.getX().getValue(),
  815.                                                         p.getY().getValue(),
  816.                                                         p.getZ().getValue());
  817.         final FieldVector3D<T> p1 = new FieldVector3D<>(p.getX().getFirstDerivative(),
  818.                                                         p.getY().getFirstDerivative(),
  819.                                                         p.getZ().getFirstDerivative());
  820.         final FieldVector3D<T> p2 = new FieldVector3D<>(p.getX().getSecondDerivative(),
  821.                                                         p.getY().getSecondDerivative(),
  822.                                                         p.getZ().getSecondDerivative());
  823.         return new TimeStampedFieldPVCoordinates<>(date, p0, p1, p2);

  824.     }

  825.     /** Computes the eccentric latitude argument from the mean latitude argument.
  826.      * @param alphaM = M + Ω mean latitude argument (rad)
  827.      * @param ex e cos(Ω), first component of circular eccentricity vector
  828.      * @param ey e sin(Ω), second component of circular eccentricity vector
  829.      * @return the eccentric latitude argument.
  830.      */
  831.     private FieldUnivariateDerivative2<T> meanToEccentric(final FieldUnivariateDerivative2<T> alphaM,
  832.                                                 final FieldUnivariateDerivative2<T> ex,
  833.                                                 final FieldUnivariateDerivative2<T> ey) {
  834.         // Generalization of Kepler equation to circular parameters
  835.         // with alphaE = PA + E and
  836.         //      alphaM = PA + M = alphaE - ex.sin(alphaE) + ey.cos(alphaE)
  837.         FieldUnivariateDerivative2<T> alphaE        = alphaM;
  838.         FieldUnivariateDerivative2<T> shift         = alphaM.getField().getZero();
  839.         FieldUnivariateDerivative2<T> alphaEMalphaM = alphaM.getField().getZero();
  840.         FieldUnivariateDerivative2<T> cosAlphaE     = alphaE.cos();
  841.         FieldUnivariateDerivative2<T> sinAlphaE     = alphaE.sin();
  842.         int                 iter          = 0;
  843.         do {
  844.             final FieldUnivariateDerivative2<T> f2 = ex.multiply(sinAlphaE).subtract(ey.multiply(cosAlphaE));
  845.             final FieldUnivariateDerivative2<T> f1 = alphaM.getField().getOne().subtract(ex.multiply(cosAlphaE)).subtract(ey.multiply(sinAlphaE));
  846.             final FieldUnivariateDerivative2<T> f0 = alphaEMalphaM.subtract(f2);

  847.             final FieldUnivariateDerivative2<T> f12 = f1.multiply(2);
  848.             shift = f0.multiply(f12).divide(f1.multiply(f12).subtract(f0.multiply(f2)));

  849.             alphaEMalphaM  = alphaEMalphaM.subtract(shift);
  850.             alphaE         = alphaM.add(alphaEMalphaM);
  851.             cosAlphaE      = alphaE.cos();
  852.             sinAlphaE      = alphaE.sin();

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

  854.         return alphaE;

  855.     }

  856.     /** {@inheritDoc} */
  857.     @Override
  858.     protected T getMass(final FieldAbsoluteDate<T> date) {
  859.         return models.get(date).mass;
  860.     }

  861.     /** {@inheritDoc} */
  862.     @Override
  863.     protected List<ParameterDriver> getParametersDrivers() {
  864.         // Eckstein Hechler propagation model does not have parameter drivers.
  865.         return Collections.emptyList();
  866.     }

  867. }