FieldKeplerianOrbit.java

  1. /* Copyright 2002-2020 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.orbits;


  18. import java.util.List;
  19. import java.util.stream.Collectors;
  20. import java.util.stream.Stream;

  21. import org.hipparchus.RealFieldElement;
  22. import org.hipparchus.analysis.differentiation.FieldUnivariateDerivative1;
  23. import org.hipparchus.analysis.interpolation.FieldHermiteInterpolator;
  24. import org.hipparchus.exception.MathIllegalArgumentException;
  25. import org.hipparchus.geometry.euclidean.threed.FieldVector3D;
  26. import org.hipparchus.util.FastMath;
  27. import org.hipparchus.util.FieldSinCos;
  28. import org.hipparchus.util.MathArrays;
  29. import org.hipparchus.util.MathUtils;
  30. import org.hipparchus.util.Precision;
  31. import org.orekit.errors.OrekitException;
  32. import org.orekit.errors.OrekitIllegalArgumentException;
  33. import org.orekit.errors.OrekitInternalError;
  34. import org.orekit.errors.OrekitMessages;
  35. import org.orekit.frames.Frame;
  36. import org.orekit.time.FieldAbsoluteDate;
  37. import org.orekit.utils.FieldPVCoordinates;
  38. import org.orekit.utils.TimeStampedFieldPVCoordinates;


  39. /**
  40.  * This class handles traditional Keplerian orbital parameters.

  41.  * <p>
  42.  * The parameters used internally are the classical Keplerian elements:
  43.  *   <pre>
  44.  *     a
  45.  *     e
  46.  *     i
  47.  *     ω
  48.  *     Ω
  49.  *     v
  50.  *   </pre>
  51.  * where ω stands for the Perigee Argument, Ω stands for the
  52.  * Right Ascension of the Ascending Node and v stands for the true anomaly.
  53.  * <p>
  54.  * This class supports hyperbolic orbits, using the convention that semi major
  55.  * axis is negative for such orbits (and of course eccentricity is greater than 1).
  56.  * </p>
  57.  * <p>
  58.  * When orbit is either equatorial or circular, some Keplerian elements
  59.  * (more precisely ω and Ω) become ambiguous so this class should not
  60.  * be used for such orbits. For this reason, {@link EquinoctialOrbit equinoctial
  61.  * orbits} is the recommended way to represent orbits.
  62.  * </p>

  63.  * <p>
  64.  * The instance <code>KeplerianOrbit</code> is guaranteed to be immutable.
  65.  * </p>
  66.  * @see     Orbit
  67.  * @see    CircularOrbit
  68.  * @see    CartesianOrbit
  69.  * @see    EquinoctialOrbit
  70.  * @author Luc Maisonobe
  71.  * @author Guylaine Prat
  72.  * @author Fabien Maussion
  73.  * @author V&eacute;ronique Pommier-Maurussane
  74.  * @author Andrea Antolino
  75.  * @since 9.0
  76.  */
  77. public class FieldKeplerianOrbit<T extends RealFieldElement<T>> extends FieldOrbit<T> {

  78.     /** Name of the eccentricity parameter. */
  79.     private static final String ECCENTRICITY = "eccentricity";

  80.     /** First coefficient to compute Kepler equation solver starter. */
  81.     private static final double A;

  82.     /** Second coefficient to compute Kepler equation solver starter. */
  83.     private static final double B;

  84.     static {
  85.         final double k1 = 3 * FastMath.PI + 2;
  86.         final double k2 = FastMath.PI - 1;
  87.         final double k3 = 6 * FastMath.PI - 1;
  88.         A  = 3 * k2 * k2 / k1;
  89.         B  = k3 * k3 / (6 * k1);
  90.     }

  91.     /** Semi-major axis (m). */
  92.     private final T a;

  93.     /** Eccentricity. */
  94.     private final T e;

  95.     /** Inclination (rad). */
  96.     private final T i;

  97.     /** Perigee Argument (rad). */
  98.     private final T pa;

  99.     /** Right Ascension of Ascending Node (rad). */
  100.     private final T raan;

  101.     /** True anomaly (rad). */
  102.     private final T v;

  103.     /** Semi-major axis derivative (m/s). */
  104.     private final T aDot;

  105.     /** Eccentricity derivative. */
  106.     private final T eDot;

  107.     /** Inclination derivative (rad/s). */
  108.     private final T iDot;

  109.     /** Perigee Argument derivative (rad/s). */
  110.     private final T paDot;

  111.     /** Right Ascension of Ascending Node derivative (rad/s). */
  112.     private final T raanDot;

  113.     /** True anomaly derivative (rad/s). */
  114.     private final T vDot;

  115.     /** Partial Cartesian coordinates (position and velocity are valid, acceleration may be missing). */
  116.     private FieldPVCoordinates<T> partialPV;

  117.     /** Identity element. */
  118.     private final T one;

  119.     /** Zero element. */
  120.     private final T zero;

  121.     /** Third Canonical Vector. */
  122.     private final FieldVector3D<T> PLUS_K;

  123.     /** Creates a new instance.
  124.      * @param a  semi-major axis (m), negative for hyperbolic orbits
  125.      * @param e eccentricity (positive or equal to 0)
  126.      * @param i inclination (rad)
  127.      * @param pa perigee argument (ω, rad)
  128.      * @param raan right ascension of ascending node (Ω, rad)
  129.      * @param anomaly mean, eccentric or true anomaly (rad)
  130.      * @param type type of anomaly
  131.      * @param frame the frame in which the parameters are defined
  132.      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
  133.      * @param date date of the orbital parameters
  134.      * @param mu central attraction coefficient (m³/s²)
  135.      * @exception IllegalArgumentException if frame is not a {@link
  136.      * Frame#isPseudoInertial pseudo-inertial frame} or a and e don't match for hyperbolic orbits,
  137.      * or v is out of range for hyperbolic orbits
  138.      */
  139.     public FieldKeplerianOrbit(final T a, final T e, final T i,
  140.                                final T pa, final T raan,
  141.                                final T anomaly, final PositionAngle type,
  142.                                final Frame frame, final FieldAbsoluteDate<T> date, final T mu)
  143.         throws IllegalArgumentException {
  144.         this(a, e, i, pa, raan, anomaly,
  145.              null, null, null, null, null, null,
  146.              type, frame, date, mu);
  147.     }

  148.     /** Creates a new instance.
  149.      * @param a  semi-major axis (m), negative for hyperbolic orbits
  150.      * @param e eccentricity (positive or equal to 0)
  151.      * @param i inclination (rad)
  152.      * @param pa perigee argument (ω, rad)
  153.      * @param raan right ascension of ascending node (Ω, rad)
  154.      * @param anomaly mean, eccentric or true anomaly (rad)
  155.      * @param aDot  semi-major axis derivative, null if unknown (m/s)
  156.      * @param eDot eccentricity derivative, null if unknown
  157.      * @param iDot inclination derivative, null if unknown (rad/s)
  158.      * @param paDot perigee argument derivative, null if unknown (rad/s)
  159.      * @param raanDot right ascension of ascending node derivative, null if unknown (rad/s)
  160.      * @param anomalyDot mean, eccentric or true anomaly derivative, null if unknown (rad/s)
  161.      * @param type type of anomaly
  162.      * @param frame the frame in which the parameters are defined
  163.      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
  164.      * @param date date of the orbital parameters
  165.      * @param mu central attraction coefficient (m³/s²)
  166.      * @exception IllegalArgumentException if frame is not a {@link
  167.      * Frame#isPseudoInertial pseudo-inertial frame} or a and e don't match for hyperbolic orbits,
  168.      * or v is out of range for hyperbolic orbits
  169.      */
  170.     public FieldKeplerianOrbit(final T a, final T e, final T i,
  171.                                final T pa, final T raan, final T anomaly,
  172.                                final T aDot, final T eDot, final T iDot,
  173.                                final T paDot, final T raanDot, final T anomalyDot,
  174.                                final PositionAngle type,
  175.                                final Frame frame, final FieldAbsoluteDate<T> date, final T mu)
  176.         throws IllegalArgumentException {
  177.         super(frame, date, mu);
  178.         if (a.multiply(e.negate().add(1)).getReal() < 0) {
  179.             throw new OrekitIllegalArgumentException(OrekitMessages.ORBIT_A_E_MISMATCH_WITH_CONIC_TYPE, a.getReal(), e.getReal());
  180.         }

  181.         // Checking eccentricity range
  182.         checkParameterRangeInclusive(ECCENTRICITY, e.getReal(), 0.0, Double.POSITIVE_INFINITY);

  183.         this.a       =    a;
  184.         this.aDot    =    aDot;
  185.         this.e       =    e;
  186.         this.eDot    =    eDot;
  187.         this.i       =    i;
  188.         this.iDot    =    iDot;
  189.         this.pa      =    pa;
  190.         this.paDot   =    paDot;
  191.         this.raan    =    raan;
  192.         this.raanDot =    raanDot;

  193.         /** Identity element. */
  194.         this.one = a.getField().getOne();

  195.         /** Zero element. */
  196.         this.zero = a.getField().getZero();

  197.         /**Third canonical vector. */
  198.         this.PLUS_K = FieldVector3D.getPlusK(a.getField());

  199.         if (hasDerivatives()) {
  200.             final FieldUnivariateDerivative1<T> eUD        = new FieldUnivariateDerivative1<>(e, eDot);
  201.             final FieldUnivariateDerivative1<T> anomalyUD  = new FieldUnivariateDerivative1<>(anomaly,  anomalyDot);
  202.             final FieldUnivariateDerivative1<T> vUD;
  203.             switch (type) {
  204.                 case MEAN :
  205.                     vUD = (a.getReal() < 0) ?
  206.                           hyperbolicEccentricToTrue(meanToHyperbolicEccentric(anomalyUD, eUD), eUD) :
  207.                           ellipticEccentricToTrue(meanToEllipticEccentric(anomalyUD, eUD), eUD);
  208.                     break;
  209.                 case ECCENTRIC :
  210.                     vUD = (a.getReal() < 0) ?
  211.                           hyperbolicEccentricToTrue(anomalyUD, eUD) :
  212.                           ellipticEccentricToTrue(anomalyUD, eUD);
  213.                     break;
  214.                 case TRUE :
  215.                     vUD = anomalyUD;
  216.                     break;
  217.                 default : // this should never happen
  218.                     throw new OrekitInternalError(null);
  219.             }
  220.             this.v    = vUD.getValue();
  221.             this.vDot = vUD.getDerivative(1);
  222.         } else {
  223.             switch (type) {
  224.                 case MEAN :

  225.                     this.v = (a.getReal() < 0) ?
  226.                              hyperbolicEccentricToTrue(meanToHyperbolicEccentric(anomaly, e), e) :
  227.                              ellipticEccentricToTrue(meanToEllipticEccentric(anomaly, e), e);

  228.                     break;
  229.                 case ECCENTRIC :
  230.                     this.v = (a.getReal() < 0) ?
  231.                              hyperbolicEccentricToTrue(anomaly, e) :
  232.                              ellipticEccentricToTrue(anomaly, e);

  233.                     break;
  234.                 case TRUE :
  235.                     this.v = anomaly;
  236.                     break;
  237.                 default : // this should never happen
  238.                     throw new OrekitInternalError(null);
  239.             }
  240.             this.vDot = null;
  241.         }

  242.         // check true anomaly range
  243.         if (e.multiply(v.cos()).add(1).getReal() <= 0) {
  244.             final double vMax = e.reciprocal().negate().acos().getReal();
  245.             throw new OrekitIllegalArgumentException(OrekitMessages.ORBIT_ANOMALY_OUT_OF_HYPERBOLIC_RANGE,
  246.                                                      v.getReal(), e.getReal(), -vMax, vMax);
  247.         }

  248.         this.partialPV = null;

  249.     }

  250.     /** Constructor from Cartesian parameters.
  251.      *
  252.      * <p> The acceleration provided in {@code FieldPVCoordinates} is accessible using
  253.      * {@link #getPVCoordinates()} and {@link #getPVCoordinates(Frame)}. All other methods
  254.      * use {@code mu} and the position to compute the acceleration, including
  255.      * {@link #shiftedBy(RealFieldElement)} and {@link #getPVCoordinates(FieldAbsoluteDate, Frame)}.
  256.      *
  257.      * @param pvCoordinates the PVCoordinates of the satellite
  258.      * @param frame the frame in which are defined the {@link FieldPVCoordinates}
  259.      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
  260.      * @param mu central attraction coefficient (m³/s²)
  261.      * @exception IllegalArgumentException if frame is not a {@link
  262.      * Frame#isPseudoInertial pseudo-inertial frame}
  263.      */
  264.     public FieldKeplerianOrbit(final TimeStampedFieldPVCoordinates<T> pvCoordinates,
  265.                                final Frame frame, final T mu)
  266.         throws IllegalArgumentException {
  267.         this(pvCoordinates, frame, mu, hasNonKeplerianAcceleration(pvCoordinates, mu));
  268.     }

  269.     /** Constructor from Cartesian parameters.
  270.      *
  271.      * <p> The acceleration provided in {@code FieldPVCoordinates} is accessible using
  272.      * {@link #getPVCoordinates()} and {@link #getPVCoordinates(Frame)}. All other methods
  273.      * use {@code mu} and the position to compute the acceleration, including
  274.      * {@link #shiftedBy(RealFieldElement)} and {@link #getPVCoordinates(FieldAbsoluteDate, Frame)}.
  275.      *
  276.      * @param pvCoordinates the PVCoordinates of the satellite
  277.      * @param frame the frame in which are defined the {@link FieldPVCoordinates}
  278.      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
  279.      * @param mu central attraction coefficient (m³/s²)
  280.      * @param reliableAcceleration if true, the acceleration is considered to be reliable
  281.      * @exception IllegalArgumentException if frame is not a {@link
  282.      * Frame#isPseudoInertial pseudo-inertial frame}
  283.      */
  284.     private FieldKeplerianOrbit(final TimeStampedFieldPVCoordinates<T> pvCoordinates,
  285.                                 final Frame frame, final T mu,
  286.                                 final boolean reliableAcceleration)
  287.         throws IllegalArgumentException {

  288.         super(pvCoordinates, frame, mu);

  289.         // identity element
  290.         this.one = pvCoordinates.getPosition().getX().getField().getOne();

  291.         // zero element
  292.         this.zero = one.getField().getZero();

  293.         // third canonical vector
  294.         this.PLUS_K = FieldVector3D.getPlusK(one.getField());

  295.         // compute inclination
  296.         final FieldVector3D<T> momentum = pvCoordinates.getMomentum();
  297.         final T m2 = momentum.getNormSq();

  298.         i = FieldVector3D.angle(momentum, PLUS_K);
  299.         // compute right ascension of ascending node
  300.         raan = FieldVector3D.crossProduct(PLUS_K, momentum).getAlpha();
  301.         // preliminary computations for parameters depending on orbit shape (elliptic or hyperbolic)
  302.         final FieldVector3D<T> pvP     = pvCoordinates.getPosition();
  303.         final FieldVector3D<T> pvV     = pvCoordinates.getVelocity();
  304.         final FieldVector3D<T> pvA     = pvCoordinates.getAcceleration();

  305.         final T   r2      = pvP.getNormSq();
  306.         final T   r       = r2.sqrt();
  307.         final T   V2      = pvV.getNormSq();
  308.         final T   rV2OnMu = r.multiply(V2).divide(mu);

  309.         // compute semi-major axis (will be negative for hyperbolic orbits)
  310.         a = r.divide(rV2OnMu.negate().add(2.0));
  311.         final T muA = a.multiply(mu);

  312.         // compute true anomaly
  313.         if (a.getReal() > 0) {
  314.             // elliptic or circular orbit
  315.             final T eSE = FieldVector3D.dotProduct(pvP, pvV).divide(muA.sqrt());
  316.             final T eCE = rV2OnMu.subtract(1);
  317.             e = (eSE.multiply(eSE).add(eCE.multiply(eCE))).sqrt();
  318.             v = ellipticEccentricToTrue(eSE.atan2(eCE), e);
  319.         } else {
  320.             // hyperbolic orbit
  321.             final T eSH = FieldVector3D.dotProduct(pvP, pvV).divide(muA.negate().sqrt());
  322.             final T eCH = rV2OnMu.subtract(1);
  323.             e = (m2.negate().divide(muA).add(1)).sqrt();
  324.             v = hyperbolicEccentricToTrue((eCH.add(eSH)).divide(eCH.subtract(eSH)).log().divide(2), e);
  325.         }

  326.         // Checking eccentricity range
  327.         checkParameterRangeInclusive(ECCENTRICITY, e.getReal(), 0.0, Double.POSITIVE_INFINITY);

  328.         // compute perigee argument
  329.         final FieldVector3D<T> node = new FieldVector3D<>(raan, zero);
  330.         final T px = FieldVector3D.dotProduct(pvP, node);
  331.         final T py = FieldVector3D.dotProduct(pvP, FieldVector3D.crossProduct(momentum, node)).divide(m2.sqrt());
  332.         pa = py.atan2(px).subtract(v);

  333.         partialPV = pvCoordinates;

  334.         if (reliableAcceleration) {
  335.             // we have a relevant acceleration, we can compute derivatives

  336.             final T[][] jacobian = MathArrays.buildArray(a.getField(), 6, 6);
  337.             getJacobianWrtCartesian(PositionAngle.MEAN, jacobian);

  338.             final FieldVector3D<T> keplerianAcceleration    = new FieldVector3D<>(r.multiply(r2).reciprocal().multiply(mu.negate()), pvP);
  339.             final FieldVector3D<T> nonKeplerianAcceleration = pvA.subtract(keplerianAcceleration);
  340.             final T   aX                       = nonKeplerianAcceleration.getX();
  341.             final T   aY                       = nonKeplerianAcceleration.getY();
  342.             final T   aZ                       = nonKeplerianAcceleration.getZ();
  343.             aDot    = jacobian[0][3].multiply(aX).add(jacobian[0][4].multiply(aY)).add(jacobian[0][5].multiply(aZ));
  344.             eDot    = jacobian[1][3].multiply(aX).add(jacobian[1][4].multiply(aY)).add(jacobian[1][5].multiply(aZ));
  345.             iDot    = jacobian[2][3].multiply(aX).add(jacobian[2][4].multiply(aY)).add(jacobian[2][5].multiply(aZ));
  346.             paDot   = jacobian[3][3].multiply(aX).add(jacobian[3][4].multiply(aY)).add(jacobian[3][5].multiply(aZ));
  347.             raanDot = jacobian[4][3].multiply(aX).add(jacobian[4][4].multiply(aY)).add(jacobian[4][5].multiply(aZ));

  348.             // in order to compute true anomaly derivative, we must compute
  349.             // mean anomaly derivative including Keplerian motion and convert to true anomaly
  350.             final T MDot = getKeplerianMeanMotion().
  351.                            add(jacobian[5][3].multiply(aX)).add(jacobian[5][4].multiply(aY)).add(jacobian[5][5].multiply(aZ));
  352.             final FieldUnivariateDerivative1<T> eUD = new FieldUnivariateDerivative1<>(e, eDot);
  353.             final FieldUnivariateDerivative1<T> MUD = new FieldUnivariateDerivative1<>(getMeanAnomaly(), MDot);
  354.             final FieldUnivariateDerivative1<T> vUD = (a.getReal() < 0) ?
  355.                                             FieldKeplerianOrbit.hyperbolicEccentricToTrue(FieldKeplerianOrbit.meanToHyperbolicEccentric(MUD, eUD), eUD) :
  356.                                             FieldKeplerianOrbit.ellipticEccentricToTrue(FieldKeplerianOrbit.meanToEllipticEccentric(MUD, eUD), eUD);
  357.             vDot = vUD.getDerivative(1);

  358.         } else {
  359.             // acceleration is either almost zero or NaN,
  360.             // we assume acceleration was not known
  361.             // we don't set up derivatives
  362.             aDot    = null;
  363.             eDot    = null;
  364.             iDot    = null;
  365.             paDot   = null;
  366.             raanDot = null;
  367.             vDot    = null;
  368.         }

  369.     }

  370.     /** Constructor from Cartesian parameters.
  371.      *
  372.      * <p> The acceleration provided in {@code FieldPVCoordinates} is accessible using
  373.      * {@link #getPVCoordinates()} and {@link #getPVCoordinates(Frame)}. All other methods
  374.      * use {@code mu} and the position to compute the acceleration, including
  375.      * {@link #shiftedBy(RealFieldElement)} and {@link #getPVCoordinates(FieldAbsoluteDate, Frame)}.
  376.      *
  377.      * @param FieldPVCoordinates the PVCoordinates of the satellite
  378.      * @param frame the frame in which are defined the {@link FieldPVCoordinates}
  379.      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
  380.      * @param date date of the orbital parameters
  381.      * @param mu central attraction coefficient (m³/s²)
  382.      * @exception IllegalArgumentException if frame is not a {@link
  383.      * Frame#isPseudoInertial pseudo-inertial frame}
  384.      */
  385.     public FieldKeplerianOrbit(final FieldPVCoordinates<T> FieldPVCoordinates,
  386.                                final Frame frame, final FieldAbsoluteDate<T> date, final T mu)
  387.         throws IllegalArgumentException {
  388.         this(new TimeStampedFieldPVCoordinates<>(date, FieldPVCoordinates), frame, mu);
  389.     }

  390.     /** Constructor from any kind of orbital parameters.
  391.      * @param op orbital parameters to copy
  392.      */
  393.     public FieldKeplerianOrbit(final FieldOrbit<T> op) {
  394.         this(op.getPVCoordinates(), op.getFrame(), op.getMu(), op.hasDerivatives());
  395.     }

  396.     /** {@inheritDoc} */
  397.     public OrbitType getType() {
  398.         return OrbitType.KEPLERIAN;
  399.     }

  400.     /** {@inheritDoc} */
  401.     public T getA() {
  402.         return a;
  403.     }

  404.     /** {@inheritDoc} */
  405.     public T getADot() {
  406.         return aDot;
  407.     }

  408.     /** {@inheritDoc} */
  409.     public T getE() {
  410.         return e;
  411.     }

  412.     /** {@inheritDoc} */
  413.     public T getEDot() {
  414.         return eDot;
  415.     }

  416.     /** {@inheritDoc} */
  417.     public T getI() {
  418.         return i;
  419.     }

  420.     /** {@inheritDoc} */
  421.     public T getIDot() {
  422.         return iDot;
  423.     }

  424.     /** Get the perigee argument.
  425.      * @return perigee argument (rad)
  426.      */
  427.     public T getPerigeeArgument() {
  428.         return pa;
  429.     }

  430.     /** Get the perigee argument derivative.
  431.      * <p>
  432.      * If the orbit was created without derivatives, the value returned is null.
  433.      * </p>
  434.      * @return perigee argument derivative (rad/s)
  435.      */
  436.     public T getPerigeeArgumentDot() {
  437.         return paDot;
  438.     }

  439.     /** Get the right ascension of the ascending node.
  440.      * @return right ascension of the ascending node (rad)
  441.      */
  442.     public T getRightAscensionOfAscendingNode() {
  443.         return raan;
  444.     }

  445.     /** Get the right ascension of the ascending node derivative.
  446.      * <p>
  447.      * If the orbit was created without derivatives, the value returned is null.
  448.      * </p>
  449.      * @return right ascension of the ascending node derivative (rad/s)
  450.      */
  451.     public T getRightAscensionOfAscendingNodeDot() {
  452.         return raanDot;
  453.     }

  454.     /** Get the true anomaly.
  455.      * @return true anomaly (rad)
  456.      */
  457.     public T getTrueAnomaly() {
  458.         return v;
  459.     }

  460.     /** Get the true anomaly derivative.
  461.      * <p>
  462.      * If the orbit was created without derivatives, the value returned is null.
  463.      * </p>
  464.      * @return true anomaly derivative (rad/s)
  465.      */
  466.     public T getTrueAnomalyDot() {
  467.         return vDot;
  468.     }

  469.     /** Get the eccentric anomaly.
  470.      * @return eccentric anomaly (rad)
  471.      */
  472.     public T getEccentricAnomaly() {
  473.         return (a.getReal() < 0) ? trueToHyperbolicEccentric(v, e) : trueToEllipticEccentric(v, e);
  474.     }

  475.     /** Get the eccentric anomaly derivative.
  476.      * <p>
  477.      * If the orbit was created without derivatives, the value returned is null.
  478.      * </p>
  479.      * @return eccentric anomaly derivative (rad/s)
  480.      */
  481.     public T getEccentricAnomalyDot() {

  482.         if (!hasDerivatives()) {
  483.             return null;
  484.         }

  485.         final FieldUnivariateDerivative1<T> eUD = new FieldUnivariateDerivative1<>(e, eDot);
  486.         final FieldUnivariateDerivative1<T> vUD = new FieldUnivariateDerivative1<>(v, vDot);
  487.         final FieldUnivariateDerivative1<T> EUD = (a.getReal() < 0) ?
  488.                                                 trueToHyperbolicEccentric(vUD, eUD) :
  489.                                                 trueToEllipticEccentric(vUD, eUD);
  490.         return EUD.getDerivative(1);

  491.     }

  492.     /** Get the mean anomaly.
  493.      * @return mean anomaly (rad)
  494.      */
  495.     public T getMeanAnomaly() {
  496.         return (a.getReal() < 0) ?
  497.                hyperbolicEccentricToMean(trueToHyperbolicEccentric(v, e), e) :
  498.                ellipticEccentricToMean(trueToEllipticEccentric(v, e), e);
  499.     }

  500.     /** Get the mean anomaly derivative.
  501.      * <p>
  502.      * If the orbit was created without derivatives, the value returned is null.
  503.      * </p>
  504.      * @return mean anomaly derivative (rad/s)
  505.      */
  506.     public T getMeanAnomalyDot() {

  507.         if (!hasDerivatives()) {
  508.             return null;
  509.         }

  510.         final FieldUnivariateDerivative1<T> eUD = new FieldUnivariateDerivative1<>(e, eDot);
  511.         final FieldUnivariateDerivative1<T> vUD = new FieldUnivariateDerivative1<>(v, vDot);
  512.         final FieldUnivariateDerivative1<T> MUD = (a.getReal() < 0) ?
  513.                                                 hyperbolicEccentricToMean(trueToHyperbolicEccentric(vUD, eUD), eUD) :
  514.                                                 ellipticEccentricToMean(trueToEllipticEccentric(vUD, eUD), eUD);
  515.         return MUD.getDerivative(1);

  516.     }

  517.     /** Get the anomaly.
  518.      * @param type type of the angle
  519.      * @return anomaly (rad)
  520.      */
  521.     public T getAnomaly(final PositionAngle type) {
  522.         return (type == PositionAngle.MEAN) ? getMeanAnomaly() :
  523.                                               ((type == PositionAngle.ECCENTRIC) ? getEccentricAnomaly() :
  524.                                                                                    getTrueAnomaly());
  525.     }

  526.     /** Get the anomaly derivative.
  527.      * <p>
  528.      * If the orbit was created without derivatives, the value returned is null.
  529.      * </p>
  530.      * @param type type of the angle
  531.      * @return anomaly derivative (rad/s)
  532.      */
  533.     public T getAnomalyDot(final PositionAngle type) {
  534.         return (type == PositionAngle.MEAN) ? getMeanAnomalyDot() :
  535.                                               ((type == PositionAngle.ECCENTRIC) ? getEccentricAnomalyDot() :
  536.                                                                                    getTrueAnomalyDot());
  537.     }

  538.     /** {@inheritDoc} */
  539.     @Override
  540.     public boolean hasDerivatives() {
  541.         return aDot != null;
  542.     }

  543.     /** Computes the true anomaly from the elliptic eccentric anomaly.
  544.      * @param E eccentric anomaly (rad)
  545.      * @param e eccentricity
  546.      * @param <T> type of the field elements
  547.      * @return v the true anomaly
  548.      */
  549.     public static <T extends RealFieldElement<T>> T ellipticEccentricToTrue(final T E, final T e) {
  550.         final T beta = e.divide(e.multiply(e).negate().add(1).sqrt().add(1));
  551.         final FieldSinCos<T> scE = FastMath.sinCos(E);
  552.         return E.add(beta.multiply(scE.sin()).divide(beta.multiply(scE.cos()).subtract(1).negate()).atan().multiply(2));
  553.     }

  554.     /** Computes the elliptic eccentric anomaly from the true anomaly.
  555.      * @param v true anomaly (rad)
  556.      * @param e eccentricity
  557.      * @param <T> type of the field elements
  558.      * @return E the elliptic eccentric anomaly
  559.      */
  560.     public static <T extends RealFieldElement<T>> T trueToEllipticEccentric(final T v, final T e) {
  561.         final T beta = e.divide(e.multiply(e).negate().add(1).sqrt().add(1));
  562.         final FieldSinCos<T> scv = FastMath.sinCos(v);
  563.         return v.subtract((beta.multiply(scv.sin()).divide(beta.multiply(scv.cos()).add(1))).atan().multiply(2));
  564.     }

  565.     /** Computes the true anomaly from the hyperbolic eccentric anomaly.
  566.      * @param H hyperbolic eccentric anomaly (rad)
  567.      * @param e eccentricity
  568.      * @param <T> type of the field elements
  569.      * @return v the true anomaly
  570.      */
  571.     public static <T extends RealFieldElement<T>> T hyperbolicEccentricToTrue(final T H, final T e) {
  572.         final T s    = e.add(1).divide(e.subtract(1)).sqrt();
  573.         final T tanH = H.multiply(0.5).tanh();
  574.         return s.multiply(tanH).atan().multiply(2);
  575.     }

  576.     /** Computes the hyperbolic eccentric anomaly from the true anomaly.
  577.      * @param v true anomaly (rad)
  578.      * @param e eccentricity
  579.      * @param <T> type of the field elements
  580.      * @return H the hyperbolic eccentric anomaly
  581.      */
  582.     public static <T extends RealFieldElement<T>> T trueToHyperbolicEccentric(final T v, final T e) {
  583.         final FieldSinCos<T> scv = FastMath.sinCos(v);
  584.         final T sinhH = e.multiply(e).subtract(1).sqrt().multiply(scv.sin()).divide(e.multiply(scv.cos()).add(1));
  585.         return sinhH.asinh();
  586.     }

  587.     /** Computes the mean anomaly from the hyperbolic eccentric anomaly.
  588.      * @param H hyperbolic eccentric anomaly (rad)
  589.      * @param e eccentricity
  590.      * @param <T> type of the field elements
  591.      * @return M the mean anomaly
  592.      */
  593.     public static <T extends RealFieldElement<T>> T hyperbolicEccentricToMean(final T H, final T e) {
  594.         return e.multiply(H.sinh()).subtract(H);
  595.     }

  596.     /** Computes the elliptic eccentric anomaly from the mean anomaly.
  597.      * <p>
  598.      * The algorithm used here for solving Kepler equation has been published
  599.      * in: "Procedures for  solving Kepler's Equation", A. W. Odell and
  600.      * R. H. Gooding, Celestial Mechanics 38 (1986) 307-334
  601.      * </p>
  602.      * @param M mean anomaly (rad)
  603.      * @param e eccentricity
  604.      * @param <T> type of the field elements
  605.      * @return E the eccentric anomaly
  606.      */
  607.     public static <T extends RealFieldElement<T>> T meanToEllipticEccentric(final T M, final T e) {
  608.         // reduce M to [-PI PI) interval
  609.         final T reducedM = MathUtils.normalizeAngle(M, M.getField().getZero());

  610.         // compute start value according to A. W. Odell and R. H. Gooding S12 starter
  611.         T E;
  612.         if (reducedM.abs().getReal() < 1.0 / 6.0) {
  613.             if (FastMath.abs(reducedM.getReal()) < Precision.SAFE_MIN) {
  614.                 // this is an Orekit change to the S12 starter, mainly used when T is some kind of derivative structure.
  615.                 // If reducedM is 0.0, the derivative of cbrt is infinite which induces NaN appearing later in
  616.                 // the computation. As in this case E and M are almost equal, we initialize E with reducedM
  617.                 E = reducedM;
  618.             } else {
  619.                 // this is the standard S12 starter
  620.                 E = reducedM.add(e.multiply( (reducedM.multiply(6).cbrt()).subtract(reducedM)));
  621.             }
  622.         } else {
  623.             if (reducedM.getReal() < 0) {
  624.                 final T w = reducedM.add(FastMath.PI);
  625.                 E = reducedM.add(e.multiply(w.multiply(A).divide(w.negate().add(B)).subtract(FastMath.PI).subtract(reducedM)));
  626.             } else {
  627.                 final T w = reducedM.negate().add(FastMath.PI);
  628.                 E = reducedM.add(e.multiply(w.multiply(A).divide(w.negate().add(B)).negate().subtract(reducedM).add(FastMath.PI)));
  629.             }
  630.         }

  631.         final T e1 = e.negate().add(1);
  632.         final boolean noCancellationRisk = (e1.getReal() + E.getReal() * E.getReal() / 6) >= 0.1;

  633.         // perform two iterations, each consisting of one Halley step and one Newton-Raphson step
  634.         for (int j = 0; j < 2; ++j) {

  635.             final T f;
  636.             T fd;
  637.             final FieldSinCos<T> scE = FastMath.sinCos(E);
  638.             final T fdd  = e.multiply(scE.sin());
  639.             final T fddd = e.multiply(scE.cos());

  640.             if (noCancellationRisk) {

  641.                 f  = (E.subtract(fdd)).subtract(reducedM);
  642.                 fd = fddd.negate().add(1);
  643.             } else {


  644.                 f  = eMeSinE(E, e).subtract(reducedM);
  645.                 final T s = E.multiply(0.5).sin();
  646.                 fd = e1.add(e.multiply(s).multiply(s).multiply(2));
  647.             }
  648.             final T dee = f.multiply(fd).divide(f.multiply(fdd).multiply(0.5).subtract(fd.multiply(fd)));

  649.             // update eccentric anomaly, using expressions that limit underflow problems
  650.             final T w = fd.add(dee.multiply(fdd.add(dee.multiply(fddd.divide(3)))).multiply(0.5));
  651.             fd = fd.add(dee.multiply(fdd.add(dee.multiply(fddd).multiply(0.5))));
  652.             E = E.subtract(f.subtract(dee.multiply(fd.subtract(w))).divide(fd));

  653.         }

  654.         // expand the result back to original range
  655.         E = E.add(M).subtract(reducedM);
  656.         return E;
  657.     }

  658.     /** Accurate computation of E - e sin(E).
  659.      * <p>
  660.      * This method is used when E is close to 0 and e close to 1,
  661.      * i.e. near the perigee of almost parabolic orbits
  662.      * </p>
  663.      * @param E eccentric anomaly
  664.      * @param e eccentricity
  665.      * @param <T> Type of the field elements
  666.      * @return E - e sin(E)
  667.      */
  668.     private static <T extends RealFieldElement<T>> T eMeSinE(final T E, final T e) {

  669.         T x = (e.negate().add(1)).multiply(E.sin());
  670.         final T mE2 = E.negate().multiply(E);
  671.         T term = E;
  672.         double d    = 0;
  673.         // the inequality test below IS intentional and should NOT be replaced by a check with a small tolerance
  674.         for (T x0 = E.getField().getZero().add(Double.NaN); !Double.valueOf(x.getReal()).equals(Double.valueOf(x0.getReal()));) {
  675.             d += 2;
  676.             term = term.multiply(mE2.divide(d * (d + 1)));
  677.             x0 = x;
  678.             x = x.subtract(term);
  679.         }
  680.         return x;
  681.     }

  682.     /** Computes the hyperbolic eccentric anomaly from the mean anomaly.
  683.      * <p>
  684.      * The algorithm used here for solving hyperbolic Kepler equation is
  685.      * Danby's iterative method (3rd order) with Vallado's initial guess.
  686.      * </p>
  687.      * @param M mean anomaly (rad)
  688.      * @param e eccentricity
  689.      * @param <T> Type of the field elements
  690.      * @return H the hyperbolic eccentric anomaly
  691.      */
  692.     public static <T extends RealFieldElement<T>> T meanToHyperbolicEccentric(final T M, final T e) {

  693.         // Resolution of hyperbolic Kepler equation for Keplerian parameters

  694.         // Initial guess
  695.         T H;
  696.         if (e.getReal() < 1.6) {
  697.             if ((-FastMath.PI < M.getReal() && M.getReal() < 0.) || M.getReal() > FastMath.PI) {
  698.                 H = M.subtract(e);
  699.             } else {
  700.                 H = M.add(e);
  701.             }
  702.         } else {
  703.             if (e.getReal() < 3.6 && M.abs().getReal() > FastMath.PI) {
  704.                 H = M.subtract(e.copySign(M));
  705.             } else {
  706.                 H = M.divide(e.subtract(1));
  707.             }
  708.         }

  709.         // Iterative computation
  710.         int iter = 0;
  711.         do {
  712.             final T f3  = e.multiply(H.cosh());
  713.             final T f2  = e.multiply(H.sinh());
  714.             final T f1  = f3.subtract(1);
  715.             final T f0  = f2.subtract(H).subtract(M);
  716.             final T f12 = f1.multiply(2);
  717.             final T d   = f0.divide(f12);
  718.             final T fdf = f1.subtract(d.multiply(f2));
  719.             final T ds  = f0.divide(fdf);

  720.             final T shift = f0.divide(fdf.add(ds.multiply(ds.multiply(f3.divide(6)))));

  721.             H = H.subtract(shift);

  722.             if ((shift.abs().getReal()) <= 1.0e-12) {
  723.                 return H;
  724.             }

  725.         } while (++iter < 50);

  726.         throw new MathIllegalArgumentException(OrekitMessages.UNABLE_TO_COMPUTE_HYPERBOLIC_ECCENTRIC_ANOMALY,
  727.                                                iter);
  728.     }

  729.     /** Computes the mean anomaly from the elliptic eccentric anomaly.
  730.      * @param E eccentric anomaly (rad)
  731.      * @param e eccentricity
  732.      * @param <T> type of the field elements
  733.      * @return M the mean anomaly
  734.      */
  735.     public static <T extends RealFieldElement<T>> T ellipticEccentricToMean(final T E, final T e) {
  736.         return E.subtract(e.multiply(E.sin()));
  737.     }

  738.     /** {@inheritDoc} */
  739.     public T getEquinoctialEx() {
  740.         return e.multiply(pa.add(raan).cos());
  741.     }

  742.     /** {@inheritDoc} */
  743.     public T getEquinoctialExDot() {

  744.         if (!hasDerivatives()) {
  745.             return null;
  746.         }

  747.         final FieldUnivariateDerivative1<T> eUD    = new FieldUnivariateDerivative1<>(e,    eDot);
  748.         final FieldUnivariateDerivative1<T> paUD   = new FieldUnivariateDerivative1<>(pa,   paDot);
  749.         final FieldUnivariateDerivative1<T> raanUD = new FieldUnivariateDerivative1<>(raan, raanDot);
  750.         return eUD.multiply(paUD.add(raanUD).cos()).getDerivative(1);

  751.     }

  752.     /** {@inheritDoc} */
  753.     public T getEquinoctialEy() {
  754.         return  e.multiply((pa.add(raan)).sin());
  755.     }

  756.     /** {@inheritDoc} */
  757.     public T getEquinoctialEyDot() {

  758.         if (!hasDerivatives()) {
  759.             return null;
  760.         }

  761.         final FieldUnivariateDerivative1<T> eUD    = new FieldUnivariateDerivative1<>(e,    eDot);
  762.         final FieldUnivariateDerivative1<T> paUD   = new FieldUnivariateDerivative1<>(pa,   paDot);
  763.         final FieldUnivariateDerivative1<T> raanUD = new FieldUnivariateDerivative1<>(raan, raanDot);
  764.         return eUD.multiply(paUD.add(raanUD).sin()).getDerivative(1);

  765.     }

  766.     /** {@inheritDoc} */
  767.     public T getHx() {
  768.         // Check for equatorial retrograde orbit
  769.         if (FastMath.abs(i.getReal() - FastMath.PI) < 1.0e-10) {
  770.             return this.zero.add(Double.NaN);
  771.         }
  772.         return  raan.cos().multiply(i.divide(2).tan());
  773.     }

  774.     /** {@inheritDoc} */
  775.     public T getHxDot() {

  776.         if (!hasDerivatives()) {
  777.             return null;
  778.         }

  779.         // Check for equatorial retrograde orbit
  780.         if (FastMath.abs(i.getReal() - FastMath.PI) < 1.0e-10) {
  781.             return this.zero.add(Double.NaN);
  782.         }

  783.         final FieldUnivariateDerivative1<T> iUD    = new FieldUnivariateDerivative1<>(i,    iDot);
  784.         final FieldUnivariateDerivative1<T> raanUD = new FieldUnivariateDerivative1<>(raan, raanDot);
  785.         return raanUD.cos().multiply(iUD.multiply(0.5).tan()).getDerivative(1);

  786.     }

  787.     /** {@inheritDoc} */
  788.     public T getHy() {
  789.         // Check for equatorial retrograde orbit
  790.         if (FastMath.abs(i.getReal() - FastMath.PI) < 1.0e-10) {
  791.             return this.zero.add(Double.NaN);
  792.         }
  793.         return  raan.sin().multiply(i.divide(2).tan());
  794.     }

  795.     /** {@inheritDoc} */
  796.     public T getHyDot() {

  797.         if (!hasDerivatives()) {
  798.             return null;
  799.         }

  800.         // Check for equatorial retrograde orbit
  801.         if (FastMath.abs(i.getReal() - FastMath.PI) < 1.0e-10) {
  802.             return this.zero.add(Double.NaN);
  803.         }

  804.         final FieldUnivariateDerivative1<T> iUD    = new FieldUnivariateDerivative1<>(i,    iDot);
  805.         final FieldUnivariateDerivative1<T> raanUD = new FieldUnivariateDerivative1<>(raan, raanDot);
  806.         return raanUD.sin().multiply(iUD.multiply(0.5).tan()).getDerivative(1);

  807.     }

  808.     /** {@inheritDoc} */
  809.     public T getLv() {
  810.         return pa.add(raan).add(v);
  811.     }

  812.     /** {@inheritDoc} */
  813.     public T getLvDot() {
  814.         return hasDerivatives() ?
  815.                paDot.add(raanDot).add(vDot) :
  816.                null;
  817.     }

  818.     /** {@inheritDoc} */
  819.     public T getLE() {
  820.         return pa.add(raan).add(getEccentricAnomaly());
  821.     }

  822.     /** {@inheritDoc} */
  823.     public T getLEDot() {
  824.         return hasDerivatives() ?
  825.                paDot.add(raanDot).add(getEccentricAnomalyDot()) :
  826.                null;
  827.     }

  828.     /** {@inheritDoc} */
  829.     public T getLM() {
  830.         return pa.add(raan).add(getMeanAnomaly());
  831.     }

  832.     /** {@inheritDoc} */
  833.     public T getLMDot() {
  834.         return hasDerivatives() ?
  835.                paDot.add(raanDot).add(getMeanAnomalyDot()) :
  836.                null;
  837.     }

  838.     /** Compute position and velocity but not acceleration.
  839.      */
  840.     private void computePVWithoutA() {

  841.         if (partialPV != null) {
  842.             // already computed
  843.             return;
  844.         }

  845.         // preliminary variables
  846.         final FieldSinCos<T> scRaan = FastMath.sinCos(raan);
  847.         final FieldSinCos<T> scPa   = FastMath.sinCos(pa);
  848.         final FieldSinCos<T> scI    = FastMath.sinCos(i);
  849.         final T cosRaan = scRaan.cos();
  850.         final T sinRaan = scRaan.sin();
  851.         final T cosPa   = scPa.cos();
  852.         final T sinPa   = scPa.sin();
  853.         final T cosI    = scI.cos();
  854.         final T sinI    = scI.sin();
  855.         final T crcp    = cosRaan.multiply(cosPa);
  856.         final T crsp    = cosRaan.multiply(sinPa);
  857.         final T srcp    = sinRaan.multiply(cosPa);
  858.         final T srsp    = sinRaan.multiply(sinPa);

  859.         // reference axes defining the orbital plane
  860.         final FieldVector3D<T> p = new FieldVector3D<>(crcp.subtract(cosI.multiply(srsp)),  srcp.add(cosI.multiply(crsp)), sinI.multiply(sinPa));
  861.         final FieldVector3D<T> q = new FieldVector3D<>(crsp.add(cosI.multiply(srcp)).negate(), cosI.multiply(crcp).subtract(srsp), sinI.multiply(cosPa));

  862.         if (a.getReal() > 0) {

  863.             // elliptical case

  864.             // elliptic eccentric anomaly
  865.             final T uME2             = e.negate().add(1).multiply(e.add(1));
  866.             final T s1Me2            = uME2.sqrt();
  867.             final FieldSinCos<T> scE = FastMath.sinCos(getEccentricAnomaly());
  868.             final T cosE             = scE.cos();
  869.             final T sinE             = scE.sin();

  870.             // coordinates of position and velocity in the orbital plane
  871.             final T x      = a.multiply(cosE.subtract(e));
  872.             final T y      = a.multiply(sinE).multiply(s1Me2);
  873.             final T factor = FastMath.sqrt(getMu().divide(a)).divide(e.negate().multiply(cosE).add(1));
  874.             final T xDot   = sinE.negate().multiply(factor);
  875.             final T yDot   = cosE.multiply(s1Me2).multiply(factor);

  876.             final FieldVector3D<T> position = new FieldVector3D<>(x, p, y, q);
  877.             final FieldVector3D<T> velocity = new FieldVector3D<>(xDot, p, yDot, q);
  878.             partialPV = new FieldPVCoordinates<>(position, velocity);

  879.         } else {

  880.             // hyperbolic case

  881.             // compute position and velocity factors
  882.             final FieldSinCos<T> scV = FastMath.sinCos(v);
  883.             final T sinV             = scV.sin();
  884.             final T cosV             = scV.cos();
  885.             final T f                = a.multiply(e.multiply(e).negate().add(1));
  886.             final T posFactor        = f.divide(e.multiply(cosV).add(1));
  887.             final T velFactor        = FastMath.sqrt(getMu().divide(f));

  888.             final FieldVector3D<T> position     = new FieldVector3D<>(posFactor.multiply(cosV), p, posFactor.multiply(sinV), q);
  889.             final FieldVector3D<T> velocity     = new FieldVector3D<>(velFactor.multiply(sinV).negate(), p, velFactor.multiply(e.add(cosV)), q);
  890.             partialPV = new FieldPVCoordinates<>(position, velocity);

  891.         }

  892.     }

  893.     /** Compute non-Keplerian part of the acceleration from first time derivatives.
  894.      * <p>
  895.      * This method should be called only when {@link #hasDerivatives()} returns true.
  896.      * </p>
  897.      * @return non-Keplerian part of the acceleration
  898.      */
  899.     private FieldVector3D<T> nonKeplerianAcceleration() {

  900.         final T[][] dCdP = MathArrays.buildArray(a.getField(), 6, 6);
  901.         getJacobianWrtParameters(PositionAngle.MEAN, dCdP);

  902.         final T nonKeplerianMeanMotion = getMeanAnomalyDot().subtract(getKeplerianMeanMotion());
  903.         final T nonKeplerianAx =     dCdP[3][0].multiply(aDot).
  904.                                  add(dCdP[3][1].multiply(eDot)).
  905.                                  add(dCdP[3][2].multiply(iDot)).
  906.                                  add(dCdP[3][3].multiply(paDot)).
  907.                                  add(dCdP[3][4].multiply(raanDot)).
  908.                                  add(dCdP[3][5].multiply(nonKeplerianMeanMotion));
  909.         final T nonKeplerianAy =     dCdP[4][0].multiply(aDot).
  910.                                  add(dCdP[4][1].multiply(eDot)).
  911.                                  add(dCdP[4][2].multiply(iDot)).
  912.                                  add(dCdP[4][3].multiply(paDot)).
  913.                                  add(dCdP[4][4].multiply(raanDot)).
  914.                                  add(dCdP[4][5].multiply(nonKeplerianMeanMotion));
  915.         final T nonKeplerianAz =     dCdP[5][0].multiply(aDot).
  916.                                  add(dCdP[5][1].multiply(eDot)).
  917.                                  add(dCdP[5][2].multiply(iDot)).
  918.                                  add(dCdP[5][3].multiply(paDot)).
  919.                                  add(dCdP[5][4].multiply(raanDot)).
  920.                                  add(dCdP[5][5].multiply(nonKeplerianMeanMotion));

  921.         return new FieldVector3D<>(nonKeplerianAx, nonKeplerianAy, nonKeplerianAz);

  922.     }

  923.     /** {@inheritDoc} */
  924.     protected TimeStampedFieldPVCoordinates<T> initPVCoordinates() {

  925.         // position and velocity
  926.         computePVWithoutA();

  927.         // acceleration
  928.         final T r2 = partialPV.getPosition().getNormSq();
  929.         final FieldVector3D<T> keplerianAcceleration = new FieldVector3D<>(r2.multiply(FastMath.sqrt(r2)).reciprocal().multiply(getMu().negate()),
  930.                                                                            partialPV.getPosition());
  931.         final FieldVector3D<T> acceleration = hasDerivatives() ?
  932.                                               keplerianAcceleration.add(nonKeplerianAcceleration()) :
  933.                                               keplerianAcceleration;

  934.         return new TimeStampedFieldPVCoordinates<>(getDate(), partialPV.getPosition(), partialPV.getVelocity(), acceleration);

  935.     }

  936.     /** {@inheritDoc} */
  937.     public FieldKeplerianOrbit<T> shiftedBy(final double dt) {
  938.         return shiftedBy(getDate().getField().getZero().add(dt));
  939.     }

  940.     /** {@inheritDoc} */
  941.     public FieldKeplerianOrbit<T> shiftedBy(final T dt) {

  942.         // use Keplerian-only motion
  943.         final FieldKeplerianOrbit<T> keplerianShifted = new FieldKeplerianOrbit<>(a, e, i, pa, raan,
  944.                                                                                   getKeplerianMeanMotion().multiply(dt).add(getMeanAnomaly()),
  945.                                                                                   PositionAngle.MEAN, getFrame(), getDate().shiftedBy(dt), getMu());

  946.         if (hasDerivatives()) {

  947.             // extract non-Keplerian acceleration from first time derivatives
  948.             final FieldVector3D<T> nonKeplerianAcceleration = nonKeplerianAcceleration();

  949.             // add quadratic effect of non-Keplerian acceleration to Keplerian-only shift
  950.             keplerianShifted.computePVWithoutA();
  951.             final FieldVector3D<T> fixedP   = new FieldVector3D<>(one, keplerianShifted.partialPV.getPosition(),
  952.                                                                   dt.multiply(dt).multiply(0.5), nonKeplerianAcceleration);
  953.             final T   fixedR2 = fixedP.getNormSq();
  954.             final T   fixedR  = fixedR2.sqrt();
  955.             final FieldVector3D<T> fixedV  = new FieldVector3D<>(one, keplerianShifted.partialPV.getVelocity(),
  956.                                                                  dt, nonKeplerianAcceleration);
  957.             final FieldVector3D<T> fixedA  = new FieldVector3D<>(fixedR2.multiply(fixedR).reciprocal().multiply(getMu().negate()),
  958.                                                                  keplerianShifted.partialPV.getPosition(),
  959.                                                                  one, nonKeplerianAcceleration);

  960.             // build a new orbit, taking non-Keplerian acceleration into account
  961.             return new FieldKeplerianOrbit<>(new TimeStampedFieldPVCoordinates<>(keplerianShifted.getDate(),
  962.                                                                                  fixedP, fixedV, fixedA),
  963.                                              keplerianShifted.getFrame(), keplerianShifted.getMu());

  964.         } else {
  965.             // Keplerian-only motion is all we can do
  966.             return keplerianShifted;
  967.         }

  968.     }

  969.     /** {@inheritDoc}
  970.      * <p>
  971.      * The interpolated instance is created by polynomial Hermite interpolation
  972.      * on Keplerian elements, without derivatives (which means the interpolation
  973.      * falls back to Lagrange interpolation only).
  974.      * </p>
  975.      * <p>
  976.      * As this implementation of interpolation is polynomial, it should be used only
  977.      * with small samples (about 10-20 points) in order to avoid <a
  978.      * href="http://en.wikipedia.org/wiki/Runge%27s_phenomenon">Runge's phenomenon</a>
  979.      * and numerical problems (including NaN appearing).
  980.      * </p>
  981.      * <p>
  982.      * If orbit interpolation on large samples is needed, using the {@link
  983.      * org.orekit.propagation.analytical.Ephemeris} class is a better way than using this
  984.      * low-level interpolation. The Ephemeris class automatically handles selection of
  985.      * a neighboring sub-sample with a predefined number of point from a large global sample
  986.      * in a thread-safe way.
  987.      * </p>
  988.      */
  989.     public FieldKeplerianOrbit<T> interpolate(final FieldAbsoluteDate<T> date, final Stream<FieldOrbit<T>> sample) {

  990.         // first pass to check if derivatives are available throughout the sample
  991.         final List<FieldOrbit<T>> list = sample.collect(Collectors.toList());
  992.         boolean useDerivatives = true;
  993.         for (final FieldOrbit<T> orbit : list) {
  994.             useDerivatives = useDerivatives && orbit.hasDerivatives();
  995.         }

  996.         // set up an interpolator
  997.         final FieldHermiteInterpolator<T> interpolator = new FieldHermiteInterpolator<>();

  998.         // second pass to feed interpolator
  999.         FieldAbsoluteDate<T> previousDate = null;
  1000.         T                    previousPA   = zero.add(Double.NaN);
  1001.         T                    previousRAAN = zero.add(Double.NaN);
  1002.         T                    previousM    = zero.add(Double.NaN);
  1003.         for (final FieldOrbit<T> orbit : list) {
  1004.             final FieldKeplerianOrbit<T> kep = (FieldKeplerianOrbit<T>) OrbitType.KEPLERIAN.convertType(orbit);
  1005.             final T continuousPA;
  1006.             final T continuousRAAN;
  1007.             final T continuousM;
  1008.             if (previousDate == null) {
  1009.                 continuousPA   = kep.getPerigeeArgument();
  1010.                 continuousRAAN = kep.getRightAscensionOfAscendingNode();
  1011.                 continuousM    = kep.getMeanAnomaly();
  1012.             } else {
  1013.                 final T dt      = kep.getDate().durationFrom(previousDate);
  1014.                 final T keplerM = previousM.add(kep.getKeplerianMeanMotion().multiply(dt));
  1015.                 continuousPA   = MathUtils.normalizeAngle(kep.getPerigeeArgument(), previousPA);
  1016.                 continuousRAAN = MathUtils.normalizeAngle(kep.getRightAscensionOfAscendingNode(), previousRAAN);
  1017.                 continuousM    = MathUtils.normalizeAngle(kep.getMeanAnomaly(), keplerM);
  1018.             }
  1019.             previousDate = kep.getDate();
  1020.             previousPA   = continuousPA;
  1021.             previousRAAN = continuousRAAN;
  1022.             previousM    = continuousM;
  1023.             final T[] toAdd = MathArrays.buildArray(getA().getField(), 6);
  1024.             toAdd[0] = kep.getA();
  1025.             toAdd[1] = kep.getE();
  1026.             toAdd[2] = kep.getI();
  1027.             toAdd[3] = continuousPA;
  1028.             toAdd[4] = continuousRAAN;
  1029.             toAdd[5] = continuousM;
  1030.             if (useDerivatives) {
  1031.                 final T[] toAddDot = MathArrays.buildArray(one.getField(), 6);
  1032.                 toAddDot[0] = kep.getADot();
  1033.                 toAddDot[1] = kep.getEDot();
  1034.                 toAddDot[2] = kep.getIDot();
  1035.                 toAddDot[3] = kep.getPerigeeArgumentDot();
  1036.                 toAddDot[4] = kep.getRightAscensionOfAscendingNodeDot();
  1037.                 toAddDot[5] = kep.getMeanAnomalyDot();
  1038.                 interpolator.addSamplePoint(kep.getDate().durationFrom(date),
  1039.                                             toAdd, toAddDot);
  1040.             } else {
  1041.                 interpolator.addSamplePoint(this.zero.add(kep.getDate().durationFrom(date)),
  1042.                                             toAdd);
  1043.             }
  1044.         }

  1045.         // interpolate
  1046.         final T[][] interpolated = interpolator.derivatives(zero, 1);

  1047.         // build a new interpolated instance
  1048.         return new FieldKeplerianOrbit<>(interpolated[0][0], interpolated[0][1], interpolated[0][2],
  1049.                                          interpolated[0][3], interpolated[0][4], interpolated[0][5],
  1050.                                          interpolated[1][0], interpolated[1][1], interpolated[1][2],
  1051.                                          interpolated[1][3], interpolated[1][4], interpolated[1][5],
  1052.                                          PositionAngle.MEAN, getFrame(), date, getMu());

  1053.     }

  1054.     /** {@inheritDoc} */
  1055.     protected T[][] computeJacobianMeanWrtCartesian() {
  1056.         if (a.getReal() > 0) {
  1057.             return computeJacobianMeanWrtCartesianElliptical();
  1058.         } else {
  1059.             return computeJacobianMeanWrtCartesianHyperbolic();
  1060.         }
  1061.     }

  1062.     /** Compute the Jacobian of the orbital parameters with respect to the Cartesian parameters.
  1063.      * <p>
  1064.      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
  1065.      * respect to Cartesian coordinate j (x for j=0, y for j=1, z for j=2, xDot for j=3,
  1066.      * yDot for j=4, zDot for j=5).
  1067.      * </p>
  1068.      * @return 6x6 Jacobian matrix
  1069.      */
  1070.     private T[][] computeJacobianMeanWrtCartesianElliptical() {

  1071.         final T[][] jacobian = MathArrays.buildArray(getA().getField(), 6, 6);

  1072.         // compute various intermediate parameters
  1073.         computePVWithoutA();
  1074.         final FieldVector3D<T> position = partialPV.getPosition();
  1075.         final FieldVector3D<T> velocity = partialPV.getVelocity();
  1076.         final FieldVector3D<T> momentum = partialPV.getMomentum();
  1077.         final T v2         = velocity.getNormSq();
  1078.         final T r2         = position.getNormSq();
  1079.         final T r          = r2.sqrt();
  1080.         final T r3         = r.multiply(r2);

  1081.         final T px         = position.getX();
  1082.         final T py         = position.getY();
  1083.         final T pz         = position.getZ();
  1084.         final T vx         = velocity.getX();
  1085.         final T vy         = velocity.getY();
  1086.         final T vz         = velocity.getZ();
  1087.         final T mx         = momentum.getX();
  1088.         final T my         = momentum.getY();
  1089.         final T mz         = momentum.getZ();

  1090.         final T mu         = getMu();
  1091.         final T sqrtMuA    = FastMath.sqrt(a.multiply(mu));
  1092.         final T sqrtAoMu   = FastMath.sqrt(a.divide(mu));
  1093.         final T a2         = a.multiply(a);
  1094.         final T twoA       = a.multiply(2);
  1095.         final T rOnA       = r.divide(a);

  1096.         final T oMe2       = e.multiply(e).negate().add(1);
  1097.         final T epsilon    = oMe2.sqrt();
  1098.         final T sqrtRec    = epsilon.reciprocal();

  1099.         final FieldSinCos<T> scI  = FastMath.sinCos(i);
  1100.         final FieldSinCos<T> scPA = FastMath.sinCos(pa);
  1101.         final T cosI       = scI.cos();
  1102.         final T sinI       = scI.sin();
  1103.         final T cosPA      = scPA.cos();
  1104.         final T sinPA      = scPA.sin();

  1105.         final T pv         = FieldVector3D.dotProduct(position, velocity);
  1106.         final T cosE       = a.subtract(r).divide(a.multiply(e));
  1107.         final T sinE       = pv.divide(e.multiply(sqrtMuA));

  1108.         // da
  1109.         final FieldVector3D<T> vectorAR = new FieldVector3D<>(a2.multiply(2).divide(r3), position);
  1110.         final FieldVector3D<T> vectorARDot = velocity.scalarMultiply(a2.multiply(mu.divide(2.).reciprocal()));
  1111.         fillHalfRow(this.one, vectorAR,    jacobian[0], 0);
  1112.         fillHalfRow(this.one, vectorARDot, jacobian[0], 3);

  1113.         // de
  1114.         final T factorER3 = pv.divide(twoA);
  1115.         final FieldVector3D<T> vectorER   = new FieldVector3D<>(cosE.multiply(v2).divide(r.multiply(mu)), position,
  1116.                                                                 sinE.divide(sqrtMuA), velocity,
  1117.                                                                 factorER3.negate().multiply(sinE).divide(sqrtMuA), vectorAR);
  1118.         final FieldVector3D<T> vectorERDot = new FieldVector3D<>(sinE.divide(sqrtMuA), position,
  1119.                                                                  cosE.multiply(mu.divide(2.).reciprocal()).multiply(r), velocity,
  1120.                                                                  factorER3.negate().multiply(sinE).divide(sqrtMuA), vectorARDot);
  1121.         fillHalfRow(this.one, vectorER,    jacobian[1], 0);
  1122.         fillHalfRow(this.one, vectorERDot, jacobian[1], 3);

  1123.         // dE / dr (Eccentric anomaly)
  1124.         final T coefE = cosE.divide(e.multiply(sqrtMuA));
  1125.         final FieldVector3D<T>  vectorEAnR =
  1126.             new FieldVector3D<>(sinE.negate().multiply(v2).divide(e.multiply(r).multiply(mu)), position, coefE, velocity,
  1127.                                 factorER3.negate().multiply(coefE), vectorAR);

  1128.         // dE / drDot
  1129.         final FieldVector3D<T>  vectorEAnRDot =
  1130.             new FieldVector3D<>(sinE.multiply(-2).multiply(r).divide(e.multiply(mu)), velocity, coefE, position,
  1131.                                 factorER3.negate().multiply(coefE), vectorARDot);

  1132.         // precomputing some more factors
  1133.         final T s1 = sinE.negate().multiply(pz).divide(r).subtract(cosE.multiply(vz).multiply(sqrtAoMu));
  1134.         final T s2 = cosE.negate().multiply(pz).divide(r3);
  1135.         final T s3 = sinE.multiply(vz).divide(sqrtMuA.multiply(-2));
  1136.         final T t1 = sqrtRec.multiply(cosE.multiply(pz).divide(r).subtract(sinE.multiply(vz).multiply(sqrtAoMu)));
  1137.         final T t2 = sqrtRec.multiply(sinE.negate().multiply(pz).divide(r3));
  1138.         final T t3 = sqrtRec.multiply(cosE.subtract(e)).multiply(vz).divide(sqrtMuA.multiply(2));
  1139.         final T t4 = sqrtRec.multiply(e.multiply(sinI).multiply(cosPA).multiply(sqrtRec).subtract(vz.multiply(sqrtAoMu)));
  1140.         final FieldVector3D<T> s = new FieldVector3D<>(cosE.divide(r), this.PLUS_K,
  1141.                                                        s1,       vectorEAnR,
  1142.                                                        s2,       position,
  1143.                                                        s3,       vectorAR);
  1144.         final FieldVector3D<T> sDot = new FieldVector3D<>(sinE.negate().multiply(sqrtAoMu), this.PLUS_K,
  1145.                                                           s1,               vectorEAnRDot,
  1146.                                                           s3,               vectorARDot);
  1147.         final FieldVector3D<T> t =
  1148.             new FieldVector3D<>(sqrtRec.multiply(sinE).divide(r), this.PLUS_K).add(new FieldVector3D<>(t1, vectorEAnR,
  1149.                                                                                                        t2, position,
  1150.                                                                                                        t3, vectorAR,
  1151.                                                                                                        t4, vectorER));
  1152.         final FieldVector3D<T> tDot = new FieldVector3D<>(sqrtRec.multiply(cosE.subtract(e)).multiply(sqrtAoMu), this.PLUS_K,
  1153.                                                           t1,                                                    vectorEAnRDot,
  1154.                                                           t3,                                                    vectorARDot,
  1155.                                                           t4,                                                    vectorERDot);

  1156.         // di
  1157.         final T factorI1 = sinI.negate().multiply(sqrtRec).divide(sqrtMuA);
  1158.         final T i1 =  factorI1;
  1159.         final T i2 =  factorI1.negate().multiply(mz).divide(twoA);
  1160.         final T i3 =  factorI1.multiply(mz).multiply(e).divide(oMe2);
  1161.         final T i4 = cosI.multiply(sinPA);
  1162.         final T i5 = cosI.multiply(cosPA);
  1163.         fillHalfRow(i1, new FieldVector3D<>(vy, vx.negate(), this.zero), i2, vectorAR, i3, vectorER, i4, s, i5, t,
  1164.                     jacobian[2], 0);
  1165.         fillHalfRow(i1, new FieldVector3D<>(py.negate(), px, this.zero), i2, vectorARDot, i3, vectorERDot, i4, sDot, i5, tDot,
  1166.                     jacobian[2], 3);

  1167.         // dpa
  1168.         fillHalfRow(cosPA.divide(sinI), s,    sinPA.negate().divide(sinI), t,    jacobian[3], 0);
  1169.         fillHalfRow(cosPA.divide(sinI), sDot, sinPA.negate().divide(sinI), tDot, jacobian[3], 3);

  1170.         // dRaan
  1171.         final T factorRaanR = (a.multiply(mu).multiply(oMe2).multiply(sinI).multiply(sinI)).reciprocal();
  1172.         fillHalfRow( factorRaanR.negate().multiply(my), new FieldVector3D<>(zero, vz, vy.negate()),
  1173.                      factorRaanR.multiply(mx), new FieldVector3D<>(vz.negate(), zero,  vx),
  1174.                      jacobian[4], 0);
  1175.         fillHalfRow(factorRaanR.negate().multiply(my), new FieldVector3D<>(zero, pz.negate(),  py),
  1176.                      factorRaanR.multiply(mx), new FieldVector3D<>(pz, zero, px.negate()),
  1177.                      jacobian[4], 3);

  1178.         // dM
  1179.         fillHalfRow(rOnA, vectorEAnR,    sinE.negate(), vectorER,    jacobian[5], 0);
  1180.         fillHalfRow(rOnA, vectorEAnRDot, sinE.negate(), vectorERDot, jacobian[5], 3);

  1181.         return jacobian;

  1182.     }

  1183.     /** Compute the Jacobian of the orbital parameters with respect to the Cartesian parameters.
  1184.      * <p>
  1185.      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
  1186.      * respect to Cartesian coordinate j (x for j=0, y for j=1, z for j=2, xDot for j=3,
  1187.      * yDot for j=4, zDot for j=5).
  1188.      * </p>
  1189.      * @return 6x6 Jacobian matrix
  1190.      */
  1191.     private T[][] computeJacobianMeanWrtCartesianHyperbolic() {

  1192.         final T[][] jacobian = MathArrays.buildArray(getA().getField(), 6, 6);

  1193.         // compute various intermediate parameters
  1194.         computePVWithoutA();
  1195.         final FieldVector3D<T> position = partialPV.getPosition();
  1196.         final FieldVector3D<T> velocity = partialPV.getVelocity();
  1197.         final FieldVector3D<T> momentum = partialPV.getMomentum();
  1198.         final T r2         = position.getNormSq();
  1199.         final T r          = r2.sqrt();
  1200.         final T r3         = r.multiply(r2);

  1201.         final T x          = position.getX();
  1202.         final T y          = position.getY();
  1203.         final T z          = position.getZ();
  1204.         final T vx         = velocity.getX();
  1205.         final T vy         = velocity.getY();
  1206.         final T vz         = velocity.getZ();
  1207.         final T mx         = momentum.getX();
  1208.         final T my         = momentum.getY();
  1209.         final T mz         = momentum.getZ();

  1210.         final T mu         = getMu();
  1211.         final T absA       = a.negate();
  1212.         final T sqrtMuA    = absA.multiply(mu).sqrt();
  1213.         final T a2         = a.multiply(a);
  1214.         final T rOa        = r.divide(absA);

  1215.         final FieldSinCos<T> scI = FastMath.sinCos(i);
  1216.         final T cosI       = scI.cos();
  1217.         final T sinI       = scI.sin();

  1218.         final T pv         = FieldVector3D.dotProduct(position, velocity);

  1219.         // da
  1220.         final FieldVector3D<T> vectorAR = new FieldVector3D<>(a2.multiply(-2).divide(r3), position);
  1221.         final FieldVector3D<T> vectorARDot = velocity.scalarMultiply(a2.multiply(-2).divide(mu));
  1222.         fillHalfRow(this.one.negate(), vectorAR,    jacobian[0], 0);
  1223.         fillHalfRow(this.one.negate(), vectorARDot, jacobian[0], 3);

  1224.         // differentials of the momentum
  1225.         final T m      = momentum.getNorm();
  1226.         final T oOm    = m.reciprocal();
  1227.         final FieldVector3D<T> dcXP = new FieldVector3D<>(this.zero, vz, vy.negate());
  1228.         final FieldVector3D<T> dcYP = new FieldVector3D<>(vz.negate(),   this.zero,  vx);
  1229.         final FieldVector3D<T> dcZP = new FieldVector3D<>( vy, vx.negate(),   this.zero);
  1230.         final FieldVector3D<T> dcXV = new FieldVector3D<>(  this.zero,  z.negate(),   y);
  1231.         final FieldVector3D<T> dcYV = new FieldVector3D<>(  z,   this.zero,  x.negate());
  1232.         final FieldVector3D<T> dcZV = new FieldVector3D<>( y.negate(),   x,   this.zero);
  1233.         final FieldVector3D<T> dCP  = new FieldVector3D<>(mx.multiply(oOm), dcXP, my.multiply(oOm), dcYP, mz.multiply(oOm), dcZP);
  1234.         final FieldVector3D<T> dCV  = new FieldVector3D<>(mx.multiply(oOm), dcXV, my.multiply(oOm), dcYV, mz.multiply(oOm), dcZV);

  1235.         // dp
  1236.         final T mOMu   = m.divide(mu);
  1237.         final FieldVector3D<T> dpP  = new FieldVector3D<>(mOMu.multiply(2), dCP);
  1238.         final FieldVector3D<T> dpV  = new FieldVector3D<>(mOMu.multiply(2), dCV);

  1239.         // de
  1240.         final T p      = m.multiply(mOMu);
  1241.         final T moO2ae = absA.multiply(2).multiply(e).reciprocal();
  1242.         final T m2OaMu = p.negate().divide(absA);
  1243.         fillHalfRow(moO2ae, dpP, m2OaMu.multiply(moO2ae), vectorAR,    jacobian[1], 0);
  1244.         fillHalfRow(moO2ae, dpV, m2OaMu.multiply(moO2ae), vectorARDot, jacobian[1], 3);

  1245.         // di
  1246.         final T cI1 = m.multiply(sinI).reciprocal();
  1247.         final T cI2 = cosI.multiply(cI1);
  1248.         fillHalfRow(cI2, dCP, cI1.negate(), dcZP, jacobian[2], 0);
  1249.         fillHalfRow(cI2, dCV, cI1.negate(), dcZV, jacobian[2], 3);


  1250.         // dPA
  1251.         final T cP1     =  y.multiply(oOm);
  1252.         final T cP2     =  x.negate().multiply(oOm);
  1253.         final T cP3     =  mx.multiply(cP1).add(my.multiply(cP2)).negate();
  1254.         final T cP4     =  cP3.multiply(oOm);
  1255.         final T cP5     =  r2.multiply(sinI).multiply(sinI).negate().reciprocal();
  1256.         final T cP6     = z.multiply(cP5);
  1257.         final T cP7     = cP3.multiply(cP5);
  1258.         final FieldVector3D<T> dacP  = new FieldVector3D<>(cP1, dcXP, cP2, dcYP, cP4, dCP, oOm, new FieldVector3D<>(my.negate(), mx, this.zero));
  1259.         final FieldVector3D<T> dacV  = new FieldVector3D<>(cP1, dcXV, cP2, dcYV, cP4, dCV);
  1260.         final FieldVector3D<T> dpoP  = new FieldVector3D<>(cP6, dacP, cP7, this.PLUS_K);
  1261.         final FieldVector3D<T> dpoV  = new FieldVector3D<>(cP6, dacV);

  1262.         final T re2     = r2.multiply(e).multiply(e);
  1263.         final T recOre2 = p.subtract(r).divide(re2);
  1264.         final T resOre2 = pv.multiply(mOMu).divide(re2);
  1265.         final FieldVector3D<T> dreP  = new FieldVector3D<>(mOMu, velocity, pv.divide(mu), dCP);
  1266.         final FieldVector3D<T> dreV  = new FieldVector3D<>(mOMu, position, pv.divide(mu), dCV);
  1267.         final FieldVector3D<T> davP  = new FieldVector3D<>(resOre2.negate(), dpP, recOre2, dreP, resOre2.divide(r), position);
  1268.         final FieldVector3D<T> davV  = new FieldVector3D<>(resOre2.negate(), dpV, recOre2, dreV);
  1269.         fillHalfRow(this.one, dpoP, this.one.negate(), davP, jacobian[3], 0);
  1270.         fillHalfRow(this.one, dpoV, this.one.negate(), davV, jacobian[3], 3);

  1271.         // dRAAN
  1272.         final T cO0 = cI1.multiply(cI1);
  1273.         final T cO1 =  mx.multiply(cO0);
  1274.         final T cO2 =  my.negate().multiply(cO0);
  1275.         fillHalfRow(cO1, dcYP, cO2, dcXP, jacobian[4], 0);
  1276.         fillHalfRow(cO1, dcYV, cO2, dcXV, jacobian[4], 3);

  1277.         // dM
  1278.         final T s2a    = pv.divide(absA.multiply(2));
  1279.         final T oObux  = m.multiply(m).add(absA.multiply(mu)).sqrt().reciprocal();
  1280.         final T scasbu = pv.multiply(oObux);
  1281.         final FieldVector3D<T> dauP = new FieldVector3D<>(sqrtMuA.reciprocal(), velocity, s2a.negate().divide(sqrtMuA), vectorAR);
  1282.         final FieldVector3D<T> dauV = new FieldVector3D<>(sqrtMuA.reciprocal(), position, s2a.negate().divide(sqrtMuA), vectorARDot);
  1283.         final FieldVector3D<T> dbuP = new FieldVector3D<>(oObux.multiply(mu.divide(2.)), vectorAR,    m.multiply(oObux), dCP);
  1284.         final FieldVector3D<T> dbuV = new FieldVector3D<>(oObux.multiply(mu.divide(2.)), vectorARDot, m.multiply(oObux), dCV);
  1285.         final FieldVector3D<T> dcuP = new FieldVector3D<>(oObux, velocity, scasbu.negate().multiply(oObux), dbuP);
  1286.         final FieldVector3D<T> dcuV = new FieldVector3D<>(oObux, position, scasbu.negate().multiply(oObux), dbuV);
  1287.         fillHalfRow(this.one, dauP, e.negate().divide(rOa.add(1)), dcuP, jacobian[5], 0);
  1288.         fillHalfRow(this.one, dauV, e.negate().divide(rOa.add(1)), dcuV, jacobian[5], 3);

  1289.         return jacobian;

  1290.     }

  1291.     /** {@inheritDoc} */
  1292.     protected T[][] computeJacobianEccentricWrtCartesian() {
  1293.         if (a.getReal() > 0) {
  1294.             return computeJacobianEccentricWrtCartesianElliptical();
  1295.         } else {
  1296.             return computeJacobianEccentricWrtCartesianHyperbolic();
  1297.         }
  1298.     }

  1299.     /** Compute the Jacobian of the orbital parameters with respect to the Cartesian parameters.
  1300.      * <p>
  1301.      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
  1302.      * respect to Cartesian coordinate j (x for j=0, y for j=1, z for j=2, xDot for j=3,
  1303.      * yDot for j=4, zDot for j=5).
  1304.      * </p>
  1305.      * @return 6x6 Jacobian matrix
  1306.      */
  1307.     private T[][] computeJacobianEccentricWrtCartesianElliptical() {

  1308.         // start by computing the Jacobian with mean angle
  1309.         final T[][] jacobian = computeJacobianMeanWrtCartesianElliptical();

  1310.         // Differentiating the Kepler equation M = E - e sin E leads to:
  1311.         // dM = (1 - e cos E) dE - sin E de
  1312.         // which is inverted and rewritten as:
  1313.         // dE = a/r dM + sin E a/r de
  1314.         final FieldSinCos<T> scE = FastMath.sinCos(getEccentricAnomaly());
  1315.         final T aOr              = e.negate().multiply(scE.cos()).add(1).reciprocal();

  1316.         // update anomaly row
  1317.         final T[] eRow           = jacobian[1];
  1318.         final T[] anomalyRow     = jacobian[5];
  1319.         for (int j = 0; j < anomalyRow.length; ++j) {
  1320.             anomalyRow[j] = aOr.multiply(anomalyRow[j].add(scE.sin().multiply(eRow[j])));
  1321.         }

  1322.         return jacobian;

  1323.     }

  1324.     /** Compute the Jacobian of the orbital parameters with respect to the Cartesian parameters.
  1325.      * <p>
  1326.      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
  1327.      * respect to Cartesian coordinate j (x for j=0, y for j=1, z for j=2, xDot for j=3,
  1328.      * yDot for j=4, zDot for j=5).
  1329.      * </p>
  1330.      * @return 6x6 Jacobian matrix
  1331.      */
  1332.     private T[][] computeJacobianEccentricWrtCartesianHyperbolic() {

  1333.         // start by computing the Jacobian with mean angle
  1334.         final T[][] jacobian = computeJacobianMeanWrtCartesianHyperbolic();

  1335.         // Differentiating the Kepler equation M = e sinh H - H leads to:
  1336.         // dM = (e cosh H - 1) dH + sinh H de
  1337.         // which is inverted and rewritten as:
  1338.         // dH = 1 / (e cosh H - 1) dM - sinh H / (e cosh H - 1) de
  1339.         final T H      = getEccentricAnomaly();
  1340.         final T coshH  = H.cosh();
  1341.         final T sinhH  = H.sinh();
  1342.         final T absaOr = e.multiply(coshH).subtract(1).reciprocal();
  1343.         // update anomaly row
  1344.         final T[] eRow       = jacobian[1];
  1345.         final T[] anomalyRow = jacobian[5];

  1346.         for (int j = 0; j < anomalyRow.length; ++j) {
  1347.             anomalyRow[j] = absaOr.multiply(anomalyRow[j].subtract(sinhH.multiply(eRow[j])));

  1348.         }

  1349.         return jacobian;

  1350.     }

  1351.     /** {@inheritDoc} */
  1352.     protected T[][] computeJacobianTrueWrtCartesian() {
  1353.         if (a.getReal() > 0) {
  1354.             return computeJacobianTrueWrtCartesianElliptical();
  1355.         } else {
  1356.             return computeJacobianTrueWrtCartesianHyperbolic();
  1357.         }
  1358.     }

  1359.     /** Compute the Jacobian of the orbital parameters with respect to the Cartesian parameters.
  1360.      * <p>
  1361.      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
  1362.      * respect to Cartesian coordinate j (x for j=0, y for j=1, z for j=2, xDot for j=3,
  1363.      * yDot for j=4, zDot for j=5).
  1364.      * </p>
  1365.      * @return 6x6 Jacobian matrix
  1366.      */
  1367.     private T[][] computeJacobianTrueWrtCartesianElliptical() {

  1368.         // start by computing the Jacobian with eccentric angle
  1369.         final T[][] jacobian = computeJacobianEccentricWrtCartesianElliptical();
  1370.         // Differentiating the eccentric anomaly equation sin E = sqrt(1-e^2) sin v / (1 + e cos v)
  1371.         // and using cos E = (e + cos v) / (1 + e cos v) to get rid of cos E leads to:
  1372.         // dE = [sqrt (1 - e^2) / (1 + e cos v)] dv - [sin E / (1 - e^2)] de
  1373.         // which is inverted and rewritten as:
  1374.         // dv = sqrt (1 - e^2) a/r dE + [sin E / sqrt (1 - e^2)] a/r de
  1375.         final T e2               = e.multiply(e);
  1376.         final T oMe2             = e2.negate().add(1);
  1377.         final T epsilon          = oMe2.sqrt();
  1378.         final FieldSinCos<T> scE = FastMath.sinCos(getEccentricAnomaly());
  1379.         final T aOr              = e.multiply(scE.cos()).negate().add(1).reciprocal();
  1380.         final T aFactor          = epsilon.multiply(aOr);
  1381.         final T eFactor          = scE.sin().multiply(aOr).divide(epsilon);

  1382.         // update anomaly row
  1383.         final T[] eRow           = jacobian[1];
  1384.         final T[] anomalyRow     = jacobian[5];
  1385.         for (int j = 0; j < anomalyRow.length; ++j) {
  1386.             anomalyRow[j] = aFactor.multiply(anomalyRow[j]).add(eFactor.multiply(eRow[j]));
  1387.         }
  1388.         return jacobian;

  1389.     }

  1390.     /** Compute the Jacobian of the orbital parameters with respect to the Cartesian parameters.
  1391.      * <p>
  1392.      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
  1393.      * respect to Cartesian coordinate j (x for j=0, y for j=1, z for j=2, xDot for j=3,
  1394.      * yDot for j=4, zDot for j=5).
  1395.      * </p>
  1396.      * @return 6x6 Jacobian matrix
  1397.      */
  1398.     private T[][] computeJacobianTrueWrtCartesianHyperbolic() {

  1399.         // start by computing the Jacobian with eccentric angle
  1400.         final T[][] jacobian = computeJacobianEccentricWrtCartesianHyperbolic();

  1401.         // Differentiating the eccentric anomaly equation sinh H = sqrt(e^2-1) sin v / (1 + e cos v)
  1402.         // and using cosh H = (e + cos v) / (1 + e cos v) to get rid of cosh H leads to:
  1403.         // dH = [sqrt (e^2 - 1) / (1 + e cos v)] dv + [sinh H / (e^2 - 1)] de
  1404.         // which is inverted and rewritten as:
  1405.         // dv = sqrt (1 - e^2) a/r dH - [sinh H / sqrt (e^2 - 1)] a/r de
  1406.         final T e2       = e.multiply(e);
  1407.         final T e2Mo     = e2.subtract(1);
  1408.         final T epsilon  = e2Mo.sqrt();
  1409.         final T H        = getEccentricAnomaly();
  1410.         final T coshH    = H.cosh();
  1411.         final T sinhH    = H.sinh();
  1412.         final T aOr      = e.multiply(coshH).subtract(1).reciprocal();
  1413.         final T aFactor  = epsilon.multiply(aOr);
  1414.         final T eFactor  = sinhH .multiply(aOr).divide(epsilon);

  1415.         // update anomaly row
  1416.         final T[] eRow           = jacobian[1];
  1417.         final T[] anomalyRow     = jacobian[5];
  1418.         for (int j = 0; j < anomalyRow.length; ++j) {
  1419.             anomalyRow[j] = aFactor.multiply(anomalyRow[j]).subtract(eFactor.multiply(eRow[j]));
  1420.         }

  1421.         return jacobian;

  1422.     }

  1423.     /** {@inheritDoc} */
  1424.     public void addKeplerContribution(final PositionAngle type, final T gm,
  1425.                                       final T[] pDot) {
  1426.         final T oMe2;
  1427.         final T ksi;
  1428.         final T absA = a.abs();
  1429.         final T n    = absA.reciprocal().multiply(gm).sqrt().divide(absA);
  1430.         switch (type) {
  1431.             case MEAN :
  1432.                 pDot[5] = pDot[5].add(n);
  1433.                 break;
  1434.             case ECCENTRIC :
  1435.                 oMe2 = e.multiply(e).negate().add(1).abs();
  1436.                 ksi  = e.multiply(v.cos()).add(1);
  1437.                 pDot[5] = pDot[5].add( n.multiply(ksi).divide(oMe2));
  1438.                 break;
  1439.             case TRUE :
  1440.                 oMe2 = e.multiply(e).negate().add(1).abs();
  1441.                 ksi  = e.multiply(v.cos()).add(1);
  1442.                 pDot[5] = pDot[5].add(n.multiply(ksi).multiply(ksi).divide(oMe2.multiply(oMe2.sqrt())));
  1443.                 break;
  1444.             default :
  1445.                 throw new OrekitInternalError(null);
  1446.         }
  1447.     }

  1448.     /**  Returns a string representation of this Keplerian parameters object.
  1449.      * @return a string representation of this object
  1450.      */
  1451.     public String toString() {
  1452.         return new StringBuffer().append("Keplerian parameters: ").append('{').
  1453.                                   append("a: ").append(a.getReal()).
  1454.                                   append("; e: ").append(e.getReal()).
  1455.                                   append("; i: ").append(FastMath.toDegrees(i.getReal())).
  1456.                                   append("; pa: ").append(FastMath.toDegrees(pa.getReal())).
  1457.                                   append("; raan: ").append(FastMath.toDegrees(raan.getReal())).
  1458.                                   append("; v: ").append(FastMath.toDegrees(v.getReal())).
  1459.                                   append(";}").toString();
  1460.     }

  1461.     /** Check if the given parameter is within an acceptable range.
  1462.      * The bounds are inclusive: an exception is raised when either of those conditions are met:
  1463.      * <ul>
  1464.      *     <li>The parameter is strictly greater than upperBound</li>
  1465.      *     <li>The parameter is strictly lower than lowerBound</li>
  1466.      * </ul>
  1467.      * <p>
  1468.      * In either of these cases, an OrekitException is raised.
  1469.      * </p>
  1470.      * @param parameterName name of the parameter
  1471.      * @param parameter value of the parameter
  1472.      * @param lowerBound lower bound of the acceptable range (inclusive)
  1473.      * @param upperBound upper bound of the acceptable range (inclusive)
  1474.      */
  1475.     private void checkParameterRangeInclusive(final String parameterName, final double parameter,
  1476.                                               final double lowerBound, final double upperBound) {
  1477.         if ((parameter < lowerBound) || (parameter > upperBound)) {
  1478.             throw new OrekitException(OrekitMessages.INVALID_PARAMETER_RANGE, parameterName,
  1479.                                       parameter, lowerBound, upperBound);
  1480.         }
  1481.     }

  1482.     /**
  1483.      * Normalize an angle in a 2&pi; wide interval around a center value.
  1484.      * <p>This method has three main uses:</p>
  1485.      * <ul>
  1486.      *   <li>normalize an angle between 0 and 2&pi;:<br>
  1487.      *       {@code a = MathUtils.normalizeAngle(a, FastMath.PI);}</li>
  1488.      *   <li>normalize an angle between -&pi; and +&pi;<br>
  1489.      *       {@code a = MathUtils.normalizeAngle(a, 0.0);}</li>
  1490.      *   <li>compute the angle between two defining angular positions:<br>
  1491.      *       {@code angle = MathUtils.normalizeAngle(end, start) - start;}</li>
  1492.      * </ul>
  1493.      * <p>Note that due to numerical accuracy and since &pi; cannot be represented
  1494.      * exactly, the result interval is <em>closed</em>, it cannot be half-closed
  1495.      * as would be more satisfactory in a purely mathematical view.</p>
  1496.      * @param a angle to normalize
  1497.      * @param center center of the desired 2&pi; interval for the result
  1498.      * @param <T> the type of the field elements
  1499.      * @return a-2k&pi; with integer k and center-&pi; &lt;= a-2k&pi; &lt;= center+&pi;
  1500.      * @deprecated replaced by {@link MathUtils#normalizeAngle(RealFieldElement, RealFieldElement)}
  1501.      */
  1502.     @Deprecated
  1503.     public static <T extends RealFieldElement<T>> T normalizeAngle(final T a, final T center) {
  1504.         return a.subtract(2 * FastMath.PI * FastMath.floor((a.getReal() + FastMath.PI - center.getReal()) / (2 * FastMath.PI)));
  1505.     }

  1506.     @Override
  1507.     public KeplerianOrbit toOrbit() {
  1508.         if (hasDerivatives()) {
  1509.             return new KeplerianOrbit(a.getReal(), e.getReal(), i.getReal(),
  1510.                                       pa.getReal(), raan.getReal(), v.getReal(),
  1511.                                       aDot.getReal(), eDot.getReal(), iDot.getReal(),
  1512.                                       paDot.getReal(), raanDot.getReal(), vDot.getReal(),
  1513.                                       PositionAngle.TRUE,
  1514.                                       getFrame(), getDate().toAbsoluteDate(), getMu().getReal());
  1515.         } else {
  1516.             return new KeplerianOrbit(a.getReal(), e.getReal(), i.getReal(),
  1517.                                       pa.getReal(), raan.getReal(), v.getReal(),
  1518.                                       PositionAngle.TRUE,
  1519.                                       getFrame(), getDate().toAbsoluteDate(), getMu().getReal());
  1520.         }
  1521.     }


  1522. }