FieldKeplerianOrbit.java

  1. /* Copyright 2002-2025 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.lang.reflect.Array;

  19. import org.hipparchus.CalculusFieldElement;
  20. import org.hipparchus.Field;
  21. import org.hipparchus.analysis.differentiation.FieldUnivariateDerivative1;
  22. import org.hipparchus.geometry.euclidean.threed.FieldVector3D;
  23. import org.hipparchus.util.FastMath;
  24. import org.hipparchus.util.FieldSinCos;
  25. import org.hipparchus.util.MathArrays;
  26. import org.orekit.errors.OrekitException;
  27. import org.orekit.errors.OrekitIllegalArgumentException;
  28. import org.orekit.errors.OrekitInternalError;
  29. import org.orekit.errors.OrekitMessages;
  30. import org.orekit.frames.FieldKinematicTransform;
  31. import org.orekit.frames.Frame;
  32. import org.orekit.time.FieldAbsoluteDate;
  33. import org.orekit.utils.FieldPVCoordinates;
  34. import org.orekit.utils.TimeStampedFieldPVCoordinates;


  35. /**
  36.  * This class handles traditional Keplerian orbital parameters.

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

  59.  * <p>
  60.  * The instance <code>KeplerianOrbit</code> is guaranteed to be immutable.
  61.  * </p>
  62.  * @see     Orbit
  63.  * @see    CircularOrbit
  64.  * @see    CartesianOrbit
  65.  * @see    EquinoctialOrbit
  66.  * @author Luc Maisonobe
  67.  * @author Guylaine Prat
  68.  * @author Fabien Maussion
  69.  * @author V&eacute;ronique Pommier-Maurussane
  70.  * @author Andrea Antolino
  71.  * @since 9.0
  72.  * @param <T> type of the field elements
  73.  */
  74. public class FieldKeplerianOrbit<T extends CalculusFieldElement<T>> extends FieldOrbit<T>
  75.         implements PositionAngleBased<FieldKeplerianOrbit<T>> {

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

  78.     /** Semi-major axis (m). */
  79.     private final T a;

  80.     /** Eccentricity. */
  81.     private final T e;

  82.     /** Inclination (rad). */
  83.     private final T i;

  84.     /** Perigee Argument (rad). */
  85.     private final T pa;

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

  88.     /** Cached anomaly (rad). */
  89.     private final T cachedAnomaly;

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

  92.     /** Eccentricity derivative. */
  93.     private final T eDot;

  94.     /** Inclination derivative (rad/s). */
  95.     private final T iDot;

  96.     /** Perigee Argument derivative (rad/s). */
  97.     private final T paDot;

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

  100.     /** Derivative of cached anomaly (rad/s). */
  101.     private final T cachedAnomalyDot;

  102.     /** Cached type of position angle. */
  103.     private final PositionAngleType cachedPositionAngleType;

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

  106.     /** PThird Canonical Vector. */
  107.     private final FieldVector3D<T> PLUS_K;

  108.     /** Creates a new instance.
  109.      * @param a  semi-major axis (m), negative for hyperbolic orbits
  110.      * @param e eccentricity (positive or equal to 0)
  111.      * @param i inclination (rad)
  112.      * @param pa perigee argument (ω, rad)
  113.      * @param raan right ascension of ascending node (Ω, rad)
  114.      * @param anomaly mean, eccentric or true anomaly (rad)
  115.      * @param type type of anomaly
  116.      * @param cachedPositionAngleType type of cached anomaly
  117.      * @param frame the frame in which the parameters are defined
  118.      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
  119.      * @param date date of the orbital parameters
  120.      * @param mu central attraction coefficient (m³/s²)
  121.      * @exception IllegalArgumentException if frame is not a {@link
  122.      * Frame#isPseudoInertial pseudo-inertial frame} or a and e don't match for hyperbolic orbits,
  123.      * or v is out of range for hyperbolic orbits
  124.      * @since 12.1
  125.      */
  126.     public FieldKeplerianOrbit(final T a, final T e, final T i,
  127.                                final T pa, final T raan,
  128.                                final T anomaly, final PositionAngleType type,
  129.                                final PositionAngleType cachedPositionAngleType,
  130.                                final Frame frame, final FieldAbsoluteDate<T> date, final T mu)
  131.         throws IllegalArgumentException {
  132.         this(a, e, i, pa, raan, anomaly,
  133.              a.getField().getZero(), a.getField().getZero(), a.getField().getZero(), a.getField().getZero(), a.getField().getZero(),
  134.              computeKeplerianAnomalyDot(type, a, e, mu, anomaly, type), type, cachedPositionAngleType, frame, date, mu);
  135.     }

  136.     /** Creates a new instance.
  137.      * @param a  semi-major axis (m), negative for hyperbolic orbits
  138.      * @param e eccentricity (positive or equal to 0)
  139.      * @param i inclination (rad)
  140.      * @param pa perigee argument (ω, rad)
  141.      * @param raan right ascension of ascending node (Ω, rad)
  142.      * @param anomaly mean, eccentric or true anomaly (rad)
  143.      * @param type type of anomaly
  144.      * @param frame the frame in which the parameters are defined
  145.      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
  146.      * @param date date of the orbital parameters
  147.      * @param mu central attraction coefficient (m³/s²)
  148.      * @exception IllegalArgumentException if frame is not a {@link
  149.      * Frame#isPseudoInertial pseudo-inertial frame} or a and e don't match for hyperbolic orbits,
  150.      * or v is out of range for hyperbolic orbits
  151.      * @since 12.1
  152.      */
  153.     public FieldKeplerianOrbit(final T a, final T e, final T i,
  154.                                final T pa, final T raan,
  155.                                final T anomaly, final PositionAngleType type,
  156.                                final Frame frame, final FieldAbsoluteDate<T> date, final T mu)
  157.             throws IllegalArgumentException {
  158.         this(a, e, i, pa, raan, anomaly, type, type, frame, date, mu);
  159.     }

  160.     /** Creates a new instance.
  161.      * @param a  semi-major axis (m), negative for hyperbolic orbits
  162.      * @param e eccentricity (positive or equal to 0)
  163.      * @param i inclination (rad)
  164.      * @param pa perigee argument (ω, rad)
  165.      * @param raan right ascension of ascending node (Ω, rad)
  166.      * @param anomaly mean, eccentric or true anomaly (rad)
  167.      * @param aDot  semi-major axis derivative, null if unknown (m/s)
  168.      * @param eDot eccentricity derivative, null if unknown
  169.      * @param iDot inclination derivative, null if unknown (rad/s)
  170.      * @param paDot perigee argument derivative, null if unknown (rad/s)
  171.      * @param raanDot right ascension of ascending node derivative, null if unknown (rad/s)
  172.      * @param anomalyDot mean, eccentric or true anomaly derivative, null if unknown (rad/s)
  173.      * @param type type of anomaly
  174.      * @param cachedPositionAngleType type of cached anomaly
  175.      * @param frame the frame in which the parameters are defined
  176.      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
  177.      * @param date date of the orbital parameters
  178.      * @param mu central attraction coefficient (m³/s²)
  179.      * @exception IllegalArgumentException if frame is not a {@link
  180.      * Frame#isPseudoInertial pseudo-inertial frame} or a and e don't match for hyperbolic orbits,
  181.      * or v is out of range for hyperbolic orbits
  182.      * @since 12.1
  183.      */
  184.     public FieldKeplerianOrbit(final T a, final T e, final T i,
  185.                                final T pa, final T raan, final T anomaly,
  186.                                final T aDot, final T eDot, final T iDot,
  187.                                final T paDot, final T raanDot, final T anomalyDot,
  188.                                final PositionAngleType type, final PositionAngleType cachedPositionAngleType,
  189.                                final Frame frame, final FieldAbsoluteDate<T> date, final T mu)
  190.         throws IllegalArgumentException {
  191.         super(frame, date, mu);
  192.         this.cachedPositionAngleType = cachedPositionAngleType;

  193.         if (a.multiply(e.negate().add(1)).getReal() < 0) {
  194.             throw new OrekitIllegalArgumentException(OrekitMessages.ORBIT_A_E_MISMATCH_WITH_CONIC_TYPE, a.getReal(), e.getReal());
  195.         }

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

  198.         this.a       =    a;
  199.         this.aDot    =    aDot;
  200.         this.e       =    e;
  201.         this.eDot    =    eDot;
  202.         this.i       =    i;
  203.         this.iDot    =    iDot;
  204.         this.pa      =    pa;
  205.         this.paDot   =    paDot;
  206.         this.raan    =    raan;
  207.         this.raanDot =    raanDot;

  208.         this.PLUS_K = FieldVector3D.getPlusK(a.getField());  // third canonical vector

  209.         final FieldUnivariateDerivative1<T> cachedAnomalyUD = initializeCachedAnomaly(anomaly, anomalyDot, type);
  210.         this.cachedAnomaly = cachedAnomalyUD.getValue();
  211.         this.cachedAnomalyDot = cachedAnomalyUD.getFirstDerivative();

  212.         // check true anomaly range
  213.         if (!isElliptical()) {
  214.             final T trueAnomaly = getTrueAnomaly();
  215.             if (e.multiply(trueAnomaly.cos()).add(1).getReal() <= 0) {
  216.                 final double vMax = e.reciprocal().negate().acos().getReal();
  217.                 throw new OrekitIllegalArgumentException(OrekitMessages.ORBIT_ANOMALY_OUT_OF_HYPERBOLIC_RANGE,
  218.                         trueAnomaly.getReal(), e.getReal(), -vMax, vMax);
  219.             }
  220.         }

  221.         this.partialPV = null;

  222.     }

  223.     /** Creates a new instance.
  224.      * @param a  semi-major axis (m), negative for hyperbolic orbits
  225.      * @param e eccentricity (positive or equal to 0)
  226.      * @param i inclination (rad)
  227.      * @param pa perigee argument (ω, rad)
  228.      * @param raan right ascension of ascending node (Ω, rad)
  229.      * @param anomaly mean, eccentric or true anomaly (rad)
  230.      * @param aDot  semi-major axis derivative, null if unknown (m/s)
  231.      * @param eDot eccentricity derivative, null if unknown
  232.      * @param iDot inclination derivative, null if unknown (rad/s)
  233.      * @param paDot perigee argument derivative, null if unknown (rad/s)
  234.      * @param raanDot right ascension of ascending node derivative, null if unknown (rad/s)
  235.      * @param anomalyDot mean, eccentric or true anomaly derivative, null if unknown (rad/s)
  236.      * @param type type of anomaly
  237.      * @param frame the frame in which the parameters are defined
  238.      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
  239.      * @param date date of the orbital parameters
  240.      * @param mu central attraction coefficient (m³/s²)
  241.      * @exception IllegalArgumentException if frame is not a {@link
  242.      * Frame#isPseudoInertial pseudo-inertial frame} or a and e don't match for hyperbolic orbits,
  243.      * or v is out of range for hyperbolic orbits
  244.      */
  245.     public FieldKeplerianOrbit(final T a, final T e, final T i,
  246.                                final T pa, final T raan, final T anomaly,
  247.                                final T aDot, final T eDot, final T iDot,
  248.                                final T paDot, final T raanDot, final T anomalyDot,
  249.                                final PositionAngleType type,
  250.                                final Frame frame, final FieldAbsoluteDate<T> date, final T mu)
  251.             throws IllegalArgumentException {
  252.         this(a, e, i, pa, raan, anomaly, aDot, eDot, iDot, paDot, raanDot, anomalyDot,
  253.                 type, type, frame, date, mu);
  254.     }

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

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

  293.         super(pvCoordinates, frame, mu);

  294.         // third canonical vector
  295.         this.PLUS_K = FieldVector3D.getPlusK(getOne().getField());

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

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

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

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

  313.         // compute true anomaly
  314.         if (isElliptical()) {
  315.             // elliptic or circular orbit
  316.             final T eSE = FieldVector3D.dotProduct(pvP, pvV).divide(muA.sqrt());
  317.             final T eCE = rV2OnMu.subtract(1);
  318.             e = (eSE.multiply(eSE).add(eCE.multiply(eCE))).sqrt();
  319.             this.cachedPositionAngleType = PositionAngleType.ECCENTRIC;
  320.             cachedAnomaly = eSE.atan2(eCE);
  321.         } else {
  322.             // hyperbolic orbit
  323.             final T eSH = FieldVector3D.dotProduct(pvP, pvV).divide(muA.negate().sqrt());
  324.             final T eCH = rV2OnMu.subtract(1);
  325.             e = (m2.negate().divide(muA).add(1)).sqrt();
  326.             this.cachedPositionAngleType = PositionAngleType.TRUE;
  327.             cachedAnomaly = FieldKeplerianAnomalyUtility.hyperbolicEccentricToTrue(e, (eCH.add(eSH)).divide(eCH.subtract(eSH)).log().divide(2));
  328.         }

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

  331.         // compute perigee argument
  332.         final FieldVector3D<T> node = new FieldVector3D<>(raan, getZero());
  333.         final T px = FieldVector3D.dotProduct(pvP, node);
  334.         final T py = FieldVector3D.dotProduct(pvP, FieldVector3D.crossProduct(momentum, node)).divide(m2.sqrt());
  335.         pa = py.atan2(px).subtract(getTrueAnomaly());

  336.         partialPV = pvCoordinates;

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

  339.             final T[][] jacobian = MathArrays.buildArray(a.getField(), 6, 6);
  340.             getJacobianWrtCartesian(PositionAngleType.MEAN, jacobian);

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

  351.             // in order to compute true anomaly derivative, we must compute
  352.             // mean anomaly derivative including Keplerian motion and convert to true anomaly
  353.             final T MDot = getKeplerianMeanMotion().
  354.                            add(jacobian[5][3].multiply(aX)).add(jacobian[5][4].multiply(aY)).add(jacobian[5][5].multiply(aZ));
  355.             final FieldUnivariateDerivative1<T> eUD = new FieldUnivariateDerivative1<>(e, eDot);
  356.             final FieldUnivariateDerivative1<T> MUD = new FieldUnivariateDerivative1<>(getMeanAnomaly(), MDot);
  357.             if (cachedPositionAngleType == PositionAngleType.ECCENTRIC) {
  358.                 final FieldUnivariateDerivative1<T> EUD = (a.getReal() < 0) ?
  359.                          FieldKeplerianAnomalyUtility.hyperbolicMeanToEccentric(eUD, MUD) :
  360.                          FieldKeplerianAnomalyUtility.ellipticMeanToEccentric(eUD, MUD);
  361.                 cachedAnomalyDot = EUD.getFirstDerivative();
  362.             } else { // TRUE
  363.                 final FieldUnivariateDerivative1<T> vUD = (a.getReal() < 0) ?
  364.                          FieldKeplerianAnomalyUtility.hyperbolicMeanToTrue(eUD, MUD) :
  365.                          FieldKeplerianAnomalyUtility.ellipticMeanToTrue(eUD, MUD);
  366.                 cachedAnomalyDot = vUD.getFirstDerivative();
  367.             }

  368.         } else {
  369.             // acceleration is either almost zero or NaN,
  370.             // we assume acceleration was not known
  371.             // we don't set up derivatives
  372.             aDot    = getZero();
  373.             eDot    = getZero();
  374.             iDot    = getZero();
  375.             paDot   = getZero();
  376.             raanDot = getZero();
  377.             cachedAnomalyDot    = computeKeplerianAnomalyDot(cachedPositionAngleType, a, e, mu, cachedAnomaly, cachedPositionAngleType);
  378.         }

  379.     }

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

  400.     /** Constructor from any kind of orbital parameters.
  401.      * @param op orbital parameters to copy
  402.      */
  403.     public FieldKeplerianOrbit(final FieldOrbit<T> op) {
  404.         this(op.getPVCoordinates(), op.getFrame(), op.getMu());
  405.     }

  406.     /** Constructor from Field and KeplerianOrbit.
  407.      * <p>Build a FieldKeplerianOrbit from non-Field KeplerianOrbit.</p>
  408.      * @param field CalculusField to base object on
  409.      * @param op non-field orbit with only "constant" terms
  410.      * @since 12.0
  411.      */
  412.     public FieldKeplerianOrbit(final Field<T> field, final KeplerianOrbit op) {
  413.         this(field.getZero().newInstance(op.getA()), field.getZero().newInstance(op.getE()), field.getZero().newInstance(op.getI()),
  414.                 field.getZero().newInstance(op.getPerigeeArgument()), field.getZero().newInstance(op.getRightAscensionOfAscendingNode()),
  415.                 field.getZero().newInstance(op.getAnomaly(op.getCachedPositionAngleType())),
  416.                 field.getZero().newInstance(op.getADot()),
  417.                 field.getZero().newInstance(op.getEDot()),
  418.                 field.getZero().newInstance(op.getIDot()),
  419.                 field.getZero().newInstance(op.getPerigeeArgumentDot()),
  420.                 field.getZero().newInstance(op.getRightAscensionOfAscendingNodeDot()),
  421.                 field.getZero().newInstance(op.getAnomalyDot(op.getCachedPositionAngleType())),
  422.                 op.getCachedPositionAngleType(), op.getFrame(),
  423.                 new FieldAbsoluteDate<>(field, op.getDate()), field.getZero().newInstance(op.getMu()));
  424.     }

  425.     /** Constructor from Field and Orbit.
  426.      * <p>Build a FieldKeplerianOrbit from any non-Field Orbit.</p>
  427.      * @param field CalculusField to base object on
  428.      * @param op non-field orbit with only "constant" terms
  429.      * @since 12.0
  430.      */
  431.     public FieldKeplerianOrbit(final Field<T> field, final Orbit op) {
  432.         this(field, (KeplerianOrbit) OrbitType.KEPLERIAN.convertType(op));
  433.     }

  434.     /** {@inheritDoc} */
  435.     @Override
  436.     public OrbitType getType() {
  437.         return OrbitType.KEPLERIAN;
  438.     }

  439.     /** {@inheritDoc} */
  440.     @Override
  441.     public T getA() {
  442.         return a;
  443.     }

  444.     /** {@inheritDoc} */
  445.     @Override
  446.     public T getADot() {
  447.         return aDot;
  448.     }

  449.     /** {@inheritDoc} */
  450.     @Override
  451.     public T getE() {
  452.         return e;
  453.     }

  454.     /** {@inheritDoc} */
  455.     @Override
  456.     public T getEDot() {
  457.         return eDot;
  458.     }

  459.     /** {@inheritDoc} */
  460.     @Override
  461.     public T getI() {
  462.         return i;
  463.     }

  464.     /** {@inheritDoc} */
  465.     @Override
  466.     public T getIDot() {
  467.         return iDot;
  468.     }

  469.     /** Get the perigee argument.
  470.      * @return perigee argument (rad)
  471.      */
  472.     public T getPerigeeArgument() {
  473.         return pa;
  474.     }

  475.     /** Get the perigee argument derivative.
  476.      * @return perigee argument derivative (rad/s)
  477.      */
  478.     public T getPerigeeArgumentDot() {
  479.         return paDot;
  480.     }

  481.     /** Get the right ascension of the ascending node.
  482.      * @return right ascension of the ascending node (rad)
  483.      */
  484.     public T getRightAscensionOfAscendingNode() {
  485.         return raan;
  486.     }

  487.     /** Get the right ascension of the ascending node derivative.
  488.      * @return right ascension of the ascending node derivative (rad/s)
  489.      */
  490.     public T getRightAscensionOfAscendingNodeDot() {
  491.         return raanDot;
  492.     }

  493.     /** Get the true anomaly.
  494.      * @return true anomaly (rad)
  495.      */
  496.     public T getTrueAnomaly() {
  497.         switch (cachedPositionAngleType) {
  498.             case MEAN: return (a.getReal() < 0) ? FieldKeplerianAnomalyUtility.hyperbolicMeanToTrue(e, cachedAnomaly) :
  499.                     FieldKeplerianAnomalyUtility.ellipticMeanToTrue(e, cachedAnomaly);

  500.             case TRUE: return cachedAnomaly;

  501.             case ECCENTRIC: return (a.getReal() < 0) ? FieldKeplerianAnomalyUtility.hyperbolicEccentricToTrue(e, cachedAnomaly) :
  502.                     FieldKeplerianAnomalyUtility.ellipticEccentricToTrue(e, cachedAnomaly);

  503.             default: throw new OrekitInternalError(null);
  504.         }
  505.     }

  506.     /** Get the true anomaly derivative.
  507.      * @return true anomaly derivative (rad/s)
  508.      */
  509.     public T getTrueAnomalyDot() {
  510.         switch (cachedPositionAngleType) {
  511.             case MEAN:
  512.                 final FieldUnivariateDerivative1<T> eUD = new FieldUnivariateDerivative1<>(e, eDot);
  513.                 final FieldUnivariateDerivative1<T> MUD = new FieldUnivariateDerivative1<>(cachedAnomaly, cachedAnomalyDot);
  514.                 final FieldUnivariateDerivative1<T> vUD = (a.getReal() < 0) ?
  515.                          FieldKeplerianAnomalyUtility.hyperbolicMeanToTrue(eUD, MUD) :
  516.                          FieldKeplerianAnomalyUtility.ellipticMeanToTrue(eUD, MUD);
  517.                 return vUD.getFirstDerivative();

  518.             case TRUE:
  519.                 return cachedAnomalyDot;

  520.             case ECCENTRIC:
  521.                 final FieldUnivariateDerivative1<T> eUD2 = new FieldUnivariateDerivative1<>(e, eDot);
  522.                 final FieldUnivariateDerivative1<T> EUD = new FieldUnivariateDerivative1<>(cachedAnomaly, cachedAnomalyDot);
  523.                 final FieldUnivariateDerivative1<T> vUD2 = (a.getReal() < 0) ?
  524.                          FieldKeplerianAnomalyUtility.hyperbolicEccentricToTrue(eUD2, EUD) :
  525.                          FieldKeplerianAnomalyUtility.ellipticEccentricToTrue(eUD2, EUD);
  526.                 return vUD2.getFirstDerivative();

  527.             default:
  528.                 throw new OrekitInternalError(null);
  529.         }
  530.     }

  531.     /** Get the eccentric anomaly.
  532.      * @return eccentric anomaly (rad)
  533.      */
  534.     public T getEccentricAnomaly() {
  535.         switch (cachedPositionAngleType) {
  536.             case MEAN:
  537.                 return (a.getReal() < 0) ? FieldKeplerianAnomalyUtility.hyperbolicMeanToEccentric(e, cachedAnomaly) :
  538.                     FieldKeplerianAnomalyUtility.ellipticMeanToEccentric(e, cachedAnomaly);

  539.             case ECCENTRIC:
  540.                 return cachedAnomaly;

  541.             case TRUE:
  542.                 return (a.getReal() < 0) ? FieldKeplerianAnomalyUtility.hyperbolicTrueToEccentric(e, cachedAnomaly) :
  543.                     FieldKeplerianAnomalyUtility.ellipticTrueToEccentric(e, cachedAnomaly);

  544.             default:
  545.                 throw new OrekitInternalError(null);
  546.         }
  547.     }

  548.     /** Get the eccentric anomaly derivative.
  549.      * @return eccentric anomaly derivative (rad/s)
  550.      */
  551.     public T getEccentricAnomalyDot() {
  552.         switch (cachedPositionAngleType) {
  553.             case ECCENTRIC:
  554.                 return cachedAnomalyDot;

  555.             case TRUE:
  556.                 final FieldUnivariateDerivative1<T> eUD = new FieldUnivariateDerivative1<>(e, eDot);
  557.                 final FieldUnivariateDerivative1<T> vUD = new FieldUnivariateDerivative1<>(cachedAnomaly, cachedAnomalyDot);
  558.                 final FieldUnivariateDerivative1<T> EUD = (a.getReal() < 0) ?
  559.                          FieldKeplerianAnomalyUtility.hyperbolicTrueToEccentric(eUD, vUD) :
  560.                          FieldKeplerianAnomalyUtility.ellipticTrueToEccentric(eUD, vUD);
  561.                 return EUD.getFirstDerivative();

  562.             case MEAN:
  563.                 final FieldUnivariateDerivative1<T> eUD2 = new FieldUnivariateDerivative1<>(e, eDot);
  564.                 final FieldUnivariateDerivative1<T> MUD = new FieldUnivariateDerivative1<>(cachedAnomaly, cachedAnomalyDot);
  565.                 final FieldUnivariateDerivative1<T> EUD2 = (a.getReal() < 0) ?
  566.                          FieldKeplerianAnomalyUtility.hyperbolicMeanToEccentric(eUD2, MUD) :
  567.                          FieldKeplerianAnomalyUtility.ellipticMeanToEccentric(eUD2, MUD);
  568.                 return EUD2.getFirstDerivative();

  569.             default:
  570.                 throw new OrekitInternalError(null);
  571.         }
  572.     }

  573.     /** Get the mean anomaly.
  574.      * @return mean anomaly (rad)
  575.      */
  576.     public T getMeanAnomaly() {
  577.         switch (cachedPositionAngleType) {
  578.             case ECCENTRIC: return (a.getReal() < 0) ? FieldKeplerianAnomalyUtility.hyperbolicEccentricToMean(e, cachedAnomaly) :
  579.                     FieldKeplerianAnomalyUtility.ellipticEccentricToMean(e, cachedAnomaly);

  580.             case MEAN: return cachedAnomaly;

  581.             case TRUE: return (a.getReal() < 0) ? FieldKeplerianAnomalyUtility.hyperbolicTrueToMean(e, cachedAnomaly) :
  582.                     FieldKeplerianAnomalyUtility.ellipticTrueToMean(e, cachedAnomaly);

  583.             default: throw new OrekitInternalError(null);
  584.         }
  585.     }

  586.     /** Get the mean anomaly derivative.

  587.      * @return mean anomaly derivative (rad/s)
  588.      */
  589.     public T getMeanAnomalyDot() {
  590.         switch (cachedPositionAngleType) {
  591.             case MEAN:
  592.                 return cachedAnomalyDot;

  593.             case ECCENTRIC:
  594.                 final FieldUnivariateDerivative1<T> eUD = new FieldUnivariateDerivative1<>(e, eDot);
  595.                 final FieldUnivariateDerivative1<T> EUD = new FieldUnivariateDerivative1<>(cachedAnomaly, cachedAnomalyDot);
  596.                 final FieldUnivariateDerivative1<T> MUD = (a.getReal() < 0) ?
  597.                          FieldKeplerianAnomalyUtility.hyperbolicEccentricToMean(eUD, EUD) :
  598.                          FieldKeplerianAnomalyUtility.ellipticEccentricToMean(eUD, EUD);
  599.                 return MUD.getFirstDerivative();

  600.             case TRUE:
  601.                 final FieldUnivariateDerivative1<T> eUD2 = new FieldUnivariateDerivative1<>(e, eDot);
  602.                 final FieldUnivariateDerivative1<T> vUD = new FieldUnivariateDerivative1<>(cachedAnomaly, cachedAnomalyDot);
  603.                 final FieldUnivariateDerivative1<T> MUD2 = (a.getReal() < 0) ?
  604.                          FieldKeplerianAnomalyUtility.hyperbolicTrueToMean(eUD2, vUD) :
  605.                          FieldKeplerianAnomalyUtility.ellipticTrueToMean(eUD2, vUD);
  606.                 return MUD2.getFirstDerivative();

  607.             default:
  608.                 throw new OrekitInternalError(null);
  609.         }
  610.     }

  611.     /** Get the anomaly.
  612.      * @param type type of the angle
  613.      * @return anomaly (rad)
  614.      */
  615.     public T getAnomaly(final PositionAngleType type) {
  616.         return (type == PositionAngleType.MEAN) ? getMeanAnomaly() :
  617.                                               ((type == PositionAngleType.ECCENTRIC) ? getEccentricAnomaly() :
  618.                                                                                    getTrueAnomaly());
  619.     }

  620.     /** Get the anomaly derivative.
  621.      * @param type type of the angle
  622.      * @return anomaly derivative (rad/s)
  623.      */
  624.     public T getAnomalyDot(final PositionAngleType type) {
  625.         return (type == PositionAngleType.MEAN) ? getMeanAnomalyDot() :
  626.                                               ((type == PositionAngleType.ECCENTRIC) ? getEccentricAnomalyDot() :
  627.                                                                                    getTrueAnomalyDot());
  628.     }

  629.     /** {@inheritDoc} */
  630.     @Override
  631.     public boolean hasNonKeplerianAcceleration() {
  632.         return aDot.getReal() != 0. || eDot.getReal() != 0 || iDot.getReal() != 0. || raanDot.getReal() != 0. || paDot.getReal() != 0. ||
  633.                 FastMath.abs(cachedAnomalyDot.subtract(computeKeplerianAnomalyDot(cachedPositionAngleType, a, e, getMu(), cachedAnomaly, cachedPositionAngleType)).getReal()) > TOLERANCE_POSITION_ANGLE_RATE;
  634.     }

  635.     /** {@inheritDoc} */
  636.     @Override
  637.     public T getEquinoctialEx() {
  638.         return e.multiply(pa.add(raan).cos());
  639.     }

  640.     /** {@inheritDoc} */
  641.     @Override
  642.     public T getEquinoctialExDot() {
  643.         if (!hasNonKeplerianRates()) {
  644.             return getZero();
  645.         }
  646.         final FieldUnivariateDerivative1<T> eUD    = new FieldUnivariateDerivative1<>(e,    eDot);
  647.         final FieldUnivariateDerivative1<T> paUD   = new FieldUnivariateDerivative1<>(pa,   paDot);
  648.         final FieldUnivariateDerivative1<T> raanUD = new FieldUnivariateDerivative1<>(raan, raanDot);
  649.         return eUD.multiply(paUD.add(raanUD).cos()).getFirstDerivative();

  650.     }

  651.     /** {@inheritDoc} */
  652.     @Override
  653.     public T getEquinoctialEy() {
  654.         return e.multiply((pa.add(raan)).sin());
  655.     }

  656.     /** {@inheritDoc} */
  657.     @Override
  658.     public T getEquinoctialEyDot() {
  659.         if (!hasNonKeplerianRates()) {
  660.             return getZero();
  661.         }
  662.         final FieldUnivariateDerivative1<T> eUD    = new FieldUnivariateDerivative1<>(e,    eDot);
  663.         final FieldUnivariateDerivative1<T> paUD   = new FieldUnivariateDerivative1<>(pa,   paDot);
  664.         final FieldUnivariateDerivative1<T> raanUD = new FieldUnivariateDerivative1<>(raan, raanDot);
  665.         return eUD.multiply(paUD.add(raanUD).sin()).getFirstDerivative();

  666.     }

  667.     /** {@inheritDoc} */
  668.     @Override
  669.     public T getHx() {
  670.         // Check for equatorial retrograde orbit
  671.         if (FastMath.abs(i.subtract(i.getPi()).getReal()) < 1.0e-10) {
  672.             return getZero().add(Double.NaN);
  673.         }
  674.         return raan.cos().multiply(i.divide(2).tan());
  675.     }

  676.     /** {@inheritDoc} */
  677.     @Override
  678.     public T getHxDot() {

  679.         // Check for equatorial retrograde orbit
  680.         if (FastMath.abs(i.subtract(i.getPi()).getReal()) < 1.0e-10) {
  681.             return getZero().add(Double.NaN);
  682.         }
  683.         if (!hasNonKeplerianRates()) {
  684.             return getZero();
  685.         }

  686.         final FieldUnivariateDerivative1<T> iUD    = new FieldUnivariateDerivative1<>(i,    iDot);
  687.         final FieldUnivariateDerivative1<T> raanUD = new FieldUnivariateDerivative1<>(raan, raanDot);
  688.         return raanUD.cos().multiply(iUD.multiply(0.5).tan()).getFirstDerivative();

  689.     }

  690.     /** {@inheritDoc} */
  691.     @Override
  692.     public T getHy() {
  693.         // Check for equatorial retrograde orbit
  694.         if (FastMath.abs(i.subtract(i.getPi()).getReal()) < 1.0e-10) {
  695.             return getZero().add(Double.NaN);
  696.         }
  697.         return raan.sin().multiply(i.divide(2).tan());
  698.     }

  699.     /** {@inheritDoc} */
  700.     @Override
  701.     public T getHyDot() {

  702.         // Check for equatorial retrograde orbit
  703.         if (FastMath.abs(i.subtract(i.getPi()).getReal()) < 1.0e-10) {
  704.             return getZero().add(Double.NaN);
  705.         }
  706.         if (!hasNonKeplerianRates()) {
  707.             return getZero();
  708.         }

  709.         final FieldUnivariateDerivative1<T> iUD    = new FieldUnivariateDerivative1<>(i,    iDot);
  710.         final FieldUnivariateDerivative1<T> raanUD = new FieldUnivariateDerivative1<>(raan, raanDot);
  711.         return raanUD.sin().multiply(iUD.multiply(0.5).tan()).getFirstDerivative();

  712.     }

  713.     /** {@inheritDoc} */
  714.     @Override
  715.     public T getLv() {
  716.         return pa.add(raan).add(getTrueAnomaly());
  717.     }

  718.     /** {@inheritDoc} */
  719.     @Override
  720.     public T getLvDot() {
  721.         return paDot.add(raanDot).add(getTrueAnomalyDot());
  722.     }

  723.     /** {@inheritDoc} */
  724.     @Override
  725.     public T getLE() {
  726.         return pa.add(raan).add(getEccentricAnomaly());
  727.     }

  728.     /** {@inheritDoc} */
  729.     @Override
  730.     public T getLEDot() {
  731.         return paDot.add(raanDot).add(getEccentricAnomalyDot());
  732.     }

  733.     /** {@inheritDoc} */
  734.     @Override
  735.     public T getLM() {
  736.         return pa.add(raan).add(getMeanAnomaly());
  737.     }

  738.     /** {@inheritDoc} */
  739.     @Override
  740.     public T getLMDot() {
  741.         return paDot.add(raanDot).add(getMeanAnomalyDot());
  742.     }

  743.     /** Initialize cached anomaly with rate.
  744.      * @param anomaly input anomaly
  745.      * @param anomalyDot rate of input anomaly
  746.      * @param inputType position angle type passed as input
  747.      * @return anomaly to cache with rate
  748.      * @since 12.1
  749.      */
  750.     private FieldUnivariateDerivative1<T> initializeCachedAnomaly(final T anomaly, final T anomalyDot,
  751.                                                                   final PositionAngleType inputType) {
  752.         if (cachedPositionAngleType == inputType) {
  753.             return new FieldUnivariateDerivative1<>(anomaly, anomalyDot);

  754.         } else {
  755.             final FieldUnivariateDerivative1<T> eUD = new FieldUnivariateDerivative1<>(e, eDot);
  756.             final FieldUnivariateDerivative1<T> anomalyUD = new FieldUnivariateDerivative1<>(anomaly, anomalyDot);

  757.             if (a.getReal() < 0) {
  758.                 switch (cachedPositionAngleType) {
  759.                     case MEAN:
  760.                         if (inputType == PositionAngleType.ECCENTRIC) {
  761.                             return FieldKeplerianAnomalyUtility.hyperbolicEccentricToMean(eUD, anomalyUD);
  762.                         } else {
  763.                             return FieldKeplerianAnomalyUtility.hyperbolicTrueToMean(eUD, anomalyUD);
  764.                         }

  765.                     case ECCENTRIC:
  766.                         if (inputType == PositionAngleType.MEAN) {
  767.                             return FieldKeplerianAnomalyUtility.hyperbolicMeanToEccentric(eUD, anomalyUD);
  768.                         } else {
  769.                             return FieldKeplerianAnomalyUtility.hyperbolicTrueToEccentric(eUD, anomalyUD);
  770.                         }

  771.                     case TRUE:
  772.                         if (inputType == PositionAngleType.MEAN) {
  773.                             return FieldKeplerianAnomalyUtility.hyperbolicMeanToTrue(eUD, anomalyUD);
  774.                         } else {
  775.                             return FieldKeplerianAnomalyUtility.hyperbolicEccentricToTrue(eUD, anomalyUD);
  776.                         }

  777.                     default:
  778.                         break;
  779.                 }

  780.             } else {
  781.                 switch (cachedPositionAngleType) {
  782.                     case MEAN:
  783.                         if (inputType == PositionAngleType.ECCENTRIC) {
  784.                             return FieldKeplerianAnomalyUtility.ellipticEccentricToMean(eUD, anomalyUD);
  785.                         } else {
  786.                             return FieldKeplerianAnomalyUtility.ellipticTrueToMean(eUD, anomalyUD);
  787.                         }

  788.                     case ECCENTRIC:
  789.                         if (inputType == PositionAngleType.MEAN) {
  790.                             return FieldKeplerianAnomalyUtility.ellipticMeanToEccentric(eUD, anomalyUD);
  791.                         } else {
  792.                             return FieldKeplerianAnomalyUtility.ellipticTrueToEccentric(eUD, anomalyUD);
  793.                         }

  794.                     case TRUE:
  795.                         if (inputType == PositionAngleType.MEAN) {
  796.                             return FieldKeplerianAnomalyUtility.ellipticMeanToTrue(eUD, anomalyUD);
  797.                         } else {
  798.                             return FieldKeplerianAnomalyUtility.ellipticEccentricToTrue(eUD, anomalyUD);
  799.                         }

  800.                     default:
  801.                         break;
  802.                 }

  803.             }
  804.             throw new OrekitInternalError(null);
  805.         }

  806.     }

  807.     /** Compute reference axes.
  808.      * @return reference axes
  809.      * @since 12.0
  810.      */
  811.     private FieldVector3D<T>[] referenceAxes() {

  812.         // preliminary variables
  813.         final FieldSinCos<T> scRaan = FastMath.sinCos(raan);
  814.         final FieldSinCos<T> scPa   = FastMath.sinCos(pa);
  815.         final FieldSinCos<T> scI    = FastMath.sinCos(i);
  816.         final T cosRaan = scRaan.cos();
  817.         final T sinRaan = scRaan.sin();
  818.         final T cosPa   = scPa.cos();
  819.         final T sinPa   = scPa.sin();
  820.         final T cosI    = scI.cos();
  821.         final T sinI    = scI.sin();
  822.         final T crcp    = cosRaan.multiply(cosPa);
  823.         final T crsp    = cosRaan.multiply(sinPa);
  824.         final T srcp    = sinRaan.multiply(cosPa);
  825.         final T srsp    = sinRaan.multiply(sinPa);

  826.         // reference axes defining the orbital plane
  827.         @SuppressWarnings("unchecked")
  828.         final FieldVector3D<T>[] axes = (FieldVector3D<T>[]) Array.newInstance(FieldVector3D.class, 2);
  829.         axes[0] = new FieldVector3D<>(crcp.subtract(cosI.multiply(srsp)),  srcp.add(cosI.multiply(crsp)), sinI.multiply(sinPa));
  830.         axes[1] = new FieldVector3D<>(crsp.add(cosI.multiply(srcp)).negate(), cosI.multiply(crcp).subtract(srsp), sinI.multiply(cosPa));

  831.         return axes;

  832.     }

  833.     /** Compute position and velocity but not acceleration.
  834.      */
  835.     private void computePVWithoutA() {

  836.         if (partialPV != null) {
  837.             // already computed
  838.             return;
  839.         }

  840.         final FieldVector3D<T>[] axes = referenceAxes();

  841.         if (isElliptical()) {

  842.             // elliptical case

  843.             // elliptic eccentric anomaly
  844.             final T uME2             = e.negate().add(1).multiply(e.add(1));
  845.             final T s1Me2            = uME2.sqrt();
  846.             final FieldSinCos<T> scE = FastMath.sinCos(getEccentricAnomaly());
  847.             final T cosE             = scE.cos();
  848.             final T sinE             = scE.sin();

  849.             // coordinates of position and velocity in the orbital plane
  850.             final T x      = a.multiply(cosE.subtract(e));
  851.             final T y      = a.multiply(sinE).multiply(s1Me2);
  852.             final T factor = FastMath.sqrt(getMu().divide(a)).divide(e.negate().multiply(cosE).add(1));
  853.             final T xDot   = sinE.negate().multiply(factor);
  854.             final T yDot   = cosE.multiply(s1Me2).multiply(factor);

  855.             final FieldVector3D<T> position = new FieldVector3D<>(x, axes[0], y, axes[1]);
  856.             final FieldVector3D<T> velocity = new FieldVector3D<>(xDot, axes[0], yDot, axes[1]);
  857.             partialPV = new FieldPVCoordinates<>(position, velocity);

  858.         } else {

  859.             // hyperbolic case

  860.             // compute position and velocity factors
  861.             final FieldSinCos<T> scV = FastMath.sinCos(getTrueAnomaly());
  862.             final T sinV             = scV.sin();
  863.             final T cosV             = scV.cos();
  864.             final T f                = a.multiply(e.square().negate().add(1));
  865.             final T posFactor        = f.divide(e.multiply(cosV).add(1));
  866.             final T velFactor        = FastMath.sqrt(getMu().divide(f));

  867.             final FieldVector3D<T> position     = new FieldVector3D<>(posFactor.multiply(cosV), axes[0], posFactor.multiply(sinV), axes[1]);
  868.             final FieldVector3D<T> velocity     = new FieldVector3D<>(velFactor.multiply(sinV).negate(), axes[0], velFactor.multiply(e.add(cosV)), axes[1]);
  869.             partialPV = new FieldPVCoordinates<>(position, velocity);

  870.         }

  871.     }

  872.     /** Compute non-Keplerian part of the acceleration from first time derivatives.
  873.      * @return non-Keplerian part of the acceleration
  874.      */
  875.     private FieldVector3D<T> nonKeplerianAcceleration() {

  876.         final T[][] dCdP = MathArrays.buildArray(a.getField(), 6, 6);
  877.         getJacobianWrtParameters(PositionAngleType.MEAN, dCdP);

  878.         final T nonKeplerianMeanMotion = getMeanAnomalyDot().subtract(getKeplerianMeanMotion());
  879.         final T nonKeplerianAx =     dCdP[3][0].multiply(aDot).
  880.                                  add(dCdP[3][1].multiply(eDot)).
  881.                                  add(dCdP[3][2].multiply(iDot)).
  882.                                  add(dCdP[3][3].multiply(paDot)).
  883.                                  add(dCdP[3][4].multiply(raanDot)).
  884.                                  add(dCdP[3][5].multiply(nonKeplerianMeanMotion));
  885.         final T nonKeplerianAy =     dCdP[4][0].multiply(aDot).
  886.                                  add(dCdP[4][1].multiply(eDot)).
  887.                                  add(dCdP[4][2].multiply(iDot)).
  888.                                  add(dCdP[4][3].multiply(paDot)).
  889.                                  add(dCdP[4][4].multiply(raanDot)).
  890.                                  add(dCdP[4][5].multiply(nonKeplerianMeanMotion));
  891.         final T nonKeplerianAz =     dCdP[5][0].multiply(aDot).
  892.                                  add(dCdP[5][1].multiply(eDot)).
  893.                                  add(dCdP[5][2].multiply(iDot)).
  894.                                  add(dCdP[5][3].multiply(paDot)).
  895.                                  add(dCdP[5][4].multiply(raanDot)).
  896.                                  add(dCdP[5][5].multiply(nonKeplerianMeanMotion));

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

  898.     }

  899.     /** {@inheritDoc} */
  900.     @Override
  901.     protected FieldVector3D<T> initPosition() {
  902.         final FieldVector3D<T>[] axes = referenceAxes();

  903.         if (isElliptical()) {

  904.             // elliptical case

  905.             // elliptic eccentric anomaly
  906.             final T uME2             = e.negate().add(1).multiply(e.add(1));
  907.             final T s1Me2            = uME2.sqrt();
  908.             final FieldSinCos<T> scE = FastMath.sinCos(getEccentricAnomaly());
  909.             final T cosE             = scE.cos();
  910.             final T sinE             = scE.sin();

  911.             return new FieldVector3D<>(a.multiply(cosE.subtract(e)), axes[0], a.multiply(sinE).multiply(s1Me2), axes[1]);

  912.         } else {

  913.             // hyperbolic case

  914.             // compute position and velocity factors
  915.             final FieldSinCos<T> scV = FastMath.sinCos(getTrueAnomaly());
  916.             final T sinV             = scV.sin();
  917.             final T cosV             = scV.cos();
  918.             final T f                = a.multiply(e.square().negate().add(1));
  919.             final T posFactor        = f.divide(e.multiply(cosV).add(1));

  920.             return new FieldVector3D<>(posFactor.multiply(cosV), axes[0], posFactor.multiply(sinV), axes[1]);

  921.         }

  922.     }

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

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

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

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

  936.     }

  937.     /** {@inheritDoc} */
  938.     @Override
  939.     public FieldKeplerianOrbit<T> inFrame(final Frame inertialFrame) {
  940.         final FieldPVCoordinates<T> fieldPVCoordinates;
  941.         if (hasNonKeplerianAcceleration()) {
  942.             fieldPVCoordinates = getPVCoordinates(inertialFrame);
  943.         } else {
  944.             final FieldKinematicTransform<T> transform = getFrame().getKinematicTransformTo(inertialFrame, getDate());
  945.             fieldPVCoordinates = transform.transformOnlyPV(getPVCoordinates());
  946.         }
  947.         final FieldKeplerianOrbit<T> fieldOrbit = new FieldKeplerianOrbit<>(fieldPVCoordinates, inertialFrame, getDate(), getMu());
  948.         if (fieldOrbit.getCachedPositionAngleType() == getCachedPositionAngleType()) {
  949.             return fieldOrbit;
  950.         } else {
  951.             return fieldOrbit.withCachedPositionAngleType(getCachedPositionAngleType());
  952.         }
  953.     }

  954.     /** {@inheritDoc} */
  955.     @Override
  956.     public FieldKeplerianOrbit<T> withCachedPositionAngleType(final PositionAngleType positionAngleType) {
  957.         return new FieldKeplerianOrbit<>(a, e, i, pa, raan, getAnomaly(positionAngleType), aDot, eDot, iDot, paDot, raanDot,
  958.                 getAnomalyDot(positionAngleType), positionAngleType, getFrame(), getDate(), getMu());
  959.     }

  960.     /** {@inheritDoc} */
  961.     @Override
  962.     public FieldKeplerianOrbit<T> shiftedBy(final double dt) {
  963.         return shiftedBy(getZero().newInstance(dt));
  964.     }

  965.     /** {@inheritDoc} */
  966.     @Override
  967.     public FieldKeplerianOrbit<T> shiftedBy(final T dt) {

  968.         // use Keplerian-only motion
  969.         final FieldKeplerianOrbit<T> keplerianShifted = new FieldKeplerianOrbit<>(a, e, i, pa, raan,
  970.                                                                                   getKeplerianMeanMotion().multiply(dt).add(getMeanAnomaly()),
  971.                                                                                   PositionAngleType.MEAN, cachedPositionAngleType, getFrame(), getDate().shiftedBy(dt), getMu());

  972.         if (hasNonKeplerianRates()) {

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

  975.             // add quadratic effect of non-Keplerian acceleration to Keplerian-only shift
  976.             keplerianShifted.computePVWithoutA();
  977.             final FieldVector3D<T> fixedP   = new FieldVector3D<>(getOne(), keplerianShifted.partialPV.getPosition(),
  978.                                                                   dt.square().multiply(0.5), nonKeplerianAcceleration);
  979.             final T   fixedR2 = fixedP.getNormSq();
  980.             final T   fixedR  = fixedR2.sqrt();
  981.             final FieldVector3D<T> fixedV  = new FieldVector3D<>(getOne(), keplerianShifted.partialPV.getVelocity(),
  982.                                                                  dt, nonKeplerianAcceleration);
  983.             final FieldVector3D<T> fixedA  = new FieldVector3D<>(fixedR2.multiply(fixedR).reciprocal().multiply(getMu().negate()),
  984.                                                                  keplerianShifted.partialPV.getPosition(),
  985.                                                                  getOne(), nonKeplerianAcceleration);

  986.             // build a new orbit, taking non-Keplerian acceleration into account
  987.             return new FieldKeplerianOrbit<>(new TimeStampedFieldPVCoordinates<>(keplerianShifted.getDate(),
  988.                                                                                  fixedP, fixedV, fixedA),
  989.                                              keplerianShifted.getFrame(), keplerianShifted.getMu());

  990.         } else {
  991.             // Keplerian-only motion is all we can do
  992.             return keplerianShifted;
  993.         }

  994.     }

  995.     /** {@inheritDoc} */
  996.     @Override
  997.     protected T[][] computeJacobianMeanWrtCartesian() {
  998.         if (isElliptical()) {
  999.             return computeJacobianMeanWrtCartesianElliptical();
  1000.         } else {
  1001.             return computeJacobianMeanWrtCartesianHyperbolic();
  1002.         }
  1003.     }

  1004.     /** Compute the Jacobian of the orbital parameters with respect to the Cartesian parameters.
  1005.      * <p>
  1006.      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
  1007.      * respect to Cartesian coordinate j (x for j=0, y for j=1, z for j=2, xDot for j=3,
  1008.      * yDot for j=4, zDot for j=5).
  1009.      * </p>
  1010.      * @return 6x6 Jacobian matrix
  1011.      */
  1012.     private T[][] computeJacobianMeanWrtCartesianElliptical() {

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

  1014.         // compute various intermediate parameters
  1015.         computePVWithoutA();
  1016.         final FieldVector3D<T> position = partialPV.getPosition();
  1017.         final FieldVector3D<T> velocity = partialPV.getVelocity();
  1018.         final FieldVector3D<T> momentum = partialPV.getMomentum();
  1019.         final T v2         = velocity.getNormSq();
  1020.         final T r2         = position.getNormSq();
  1021.         final T r          = r2.sqrt();
  1022.         final T r3         = r.multiply(r2);

  1023.         final T px         = position.getX();
  1024.         final T py         = position.getY();
  1025.         final T pz         = position.getZ();
  1026.         final T vx         = velocity.getX();
  1027.         final T vy         = velocity.getY();
  1028.         final T vz         = velocity.getZ();
  1029.         final T mx         = momentum.getX();
  1030.         final T my         = momentum.getY();
  1031.         final T mz         = momentum.getZ();

  1032.         final T mu         = getMu();
  1033.         final T sqrtMuA    = FastMath.sqrt(a.multiply(mu));
  1034.         final T sqrtAoMu   = FastMath.sqrt(a.divide(mu));
  1035.         final T a2         = a.square();
  1036.         final T twoA       = a.multiply(2);
  1037.         final T rOnA       = r.divide(a);

  1038.         final T oMe2       = e.square().negate().add(1);
  1039.         final T epsilon    = oMe2.sqrt();
  1040.         final T sqrtRec    = epsilon.reciprocal();

  1041.         final FieldSinCos<T> scI  = FastMath.sinCos(i);
  1042.         final FieldSinCos<T> scPA = FastMath.sinCos(pa);
  1043.         final T cosI       = scI.cos();
  1044.         final T sinI       = scI.sin();
  1045.         final T cosPA      = scPA.cos();
  1046.         final T sinPA      = scPA.sin();

  1047.         final T pv         = FieldVector3D.dotProduct(position, velocity);
  1048.         final T cosE       = a.subtract(r).divide(a.multiply(e));
  1049.         final T sinE       = pv.divide(e.multiply(sqrtMuA));

  1050.         // da
  1051.         final FieldVector3D<T> vectorAR = new FieldVector3D<>(a2.multiply(2).divide(r3), position);
  1052.         final FieldVector3D<T> vectorARDot = velocity.scalarMultiply(a2.multiply(mu.divide(2.).reciprocal()));
  1053.         fillHalfRow(getOne(), vectorAR,    jacobian[0], 0);
  1054.         fillHalfRow(getOne(), vectorARDot, jacobian[0], 3);

  1055.         // de
  1056.         final T factorER3 = pv.divide(twoA);
  1057.         final FieldVector3D<T> vectorER   = new FieldVector3D<>(cosE.multiply(v2).divide(r.multiply(mu)), position,
  1058.                                                                 sinE.divide(sqrtMuA), velocity,
  1059.                                                                 factorER3.negate().multiply(sinE).divide(sqrtMuA), vectorAR);
  1060.         final FieldVector3D<T> vectorERDot = new FieldVector3D<>(sinE.divide(sqrtMuA), position,
  1061.                                                                  cosE.multiply(mu.divide(2.).reciprocal()).multiply(r), velocity,
  1062.                                                                  factorER3.negate().multiply(sinE).divide(sqrtMuA), vectorARDot);
  1063.         fillHalfRow(getOne(), vectorER,    jacobian[1], 0);
  1064.         fillHalfRow(getOne(), vectorERDot, jacobian[1], 3);

  1065.         // dE / dr (Eccentric anomaly)
  1066.         final T coefE = cosE.divide(e.multiply(sqrtMuA));
  1067.         final FieldVector3D<T>  vectorEAnR =
  1068.             new FieldVector3D<>(sinE.negate().multiply(v2).divide(e.multiply(r).multiply(mu)), position, coefE, velocity,
  1069.                                 factorER3.negate().multiply(coefE), vectorAR);

  1070.         // dE / drDot
  1071.         final FieldVector3D<T>  vectorEAnRDot =
  1072.             new FieldVector3D<>(sinE.multiply(-2).multiply(r).divide(e.multiply(mu)), velocity, coefE, position,
  1073.                                 factorER3.negate().multiply(coefE), vectorARDot);

  1074.         // precomputing some more factors
  1075.         final T s1 = sinE.negate().multiply(pz).divide(r).subtract(cosE.multiply(vz).multiply(sqrtAoMu));
  1076.         final T s2 = cosE.negate().multiply(pz).divide(r3);
  1077.         final T s3 = sinE.multiply(vz).divide(sqrtMuA.multiply(-2));
  1078.         final T t1 = sqrtRec.multiply(cosE.multiply(pz).divide(r).subtract(sinE.multiply(vz).multiply(sqrtAoMu)));
  1079.         final T t2 = sqrtRec.multiply(sinE.negate().multiply(pz).divide(r3));
  1080.         final T t3 = sqrtRec.multiply(cosE.subtract(e)).multiply(vz).divide(sqrtMuA.multiply(2));
  1081.         final T t4 = sqrtRec.multiply(e.multiply(sinI).multiply(cosPA).multiply(sqrtRec).subtract(vz.multiply(sqrtAoMu)));
  1082.         final FieldVector3D<T> s = new FieldVector3D<>(cosE.divide(r), this.PLUS_K,
  1083.                                                        s1,       vectorEAnR,
  1084.                                                        s2,       position,
  1085.                                                        s3,       vectorAR);
  1086.         final FieldVector3D<T> sDot = new FieldVector3D<>(sinE.negate().multiply(sqrtAoMu), this.PLUS_K,
  1087.                                                           s1,               vectorEAnRDot,
  1088.                                                           s3,               vectorARDot);
  1089.         final FieldVector3D<T> t =
  1090.             new FieldVector3D<>(sqrtRec.multiply(sinE).divide(r), this.PLUS_K).add(new FieldVector3D<>(t1, vectorEAnR,
  1091.                                                                                                        t2, position,
  1092.                                                                                                        t3, vectorAR,
  1093.                                                                                                        t4, vectorER));
  1094.         final FieldVector3D<T> tDot = new FieldVector3D<>(sqrtRec.multiply(cosE.subtract(e)).multiply(sqrtAoMu), this.PLUS_K,
  1095.                                                           t1,                                                    vectorEAnRDot,
  1096.                                                           t3,                                                    vectorARDot,
  1097.                                                           t4,                                                    vectorERDot);

  1098.         // di
  1099.         final T factorI1 = sinI.negate().multiply(sqrtRec).divide(sqrtMuA);
  1100.         final T i1 =  factorI1;
  1101.         final T i2 =  factorI1.negate().multiply(mz).divide(twoA);
  1102.         final T i3 =  factorI1.multiply(mz).multiply(e).divide(oMe2);
  1103.         final T i4 = cosI.multiply(sinPA);
  1104.         final T i5 = cosI.multiply(cosPA);
  1105.         fillHalfRow(i1, new FieldVector3D<>(vy, vx.negate(), getZero()), i2, vectorAR, i3, vectorER, i4, s, i5, t,
  1106.                     jacobian[2], 0);
  1107.         fillHalfRow(i1, new FieldVector3D<>(py.negate(), px, getZero()), i2, vectorARDot, i3, vectorERDot, i4, sDot, i5, tDot,
  1108.                     jacobian[2], 3);

  1109.         // dpa
  1110.         fillHalfRow(cosPA.divide(sinI), s,    sinPA.negate().divide(sinI), t,    jacobian[3], 0);
  1111.         fillHalfRow(cosPA.divide(sinI), sDot, sinPA.negate().divide(sinI), tDot, jacobian[3], 3);

  1112.         // dRaan
  1113.         final T factorRaanR = (a.multiply(mu).multiply(oMe2).multiply(sinI).multiply(sinI)).reciprocal();
  1114.         fillHalfRow( factorRaanR.negate().multiply(my), new FieldVector3D<>(getZero(), vz, vy.negate()),
  1115.                      factorRaanR.multiply(mx), new FieldVector3D<>(vz.negate(), getZero(),  vx),
  1116.                      jacobian[4], 0);
  1117.         fillHalfRow(factorRaanR.negate().multiply(my), new FieldVector3D<>(getZero(), pz.negate(),  py),
  1118.                      factorRaanR.multiply(mx), new FieldVector3D<>(pz, getZero(), px.negate()),
  1119.                      jacobian[4], 3);

  1120.         // dM
  1121.         fillHalfRow(rOnA, vectorEAnR,    sinE.negate(), vectorER,    jacobian[5], 0);
  1122.         fillHalfRow(rOnA, vectorEAnRDot, sinE.negate(), vectorERDot, jacobian[5], 3);

  1123.         return jacobian;

  1124.     }

  1125.     /** Compute the Jacobian of the orbital parameters with respect to the Cartesian parameters.
  1126.      * <p>
  1127.      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
  1128.      * respect to Cartesian coordinate j (x for j=0, y for j=1, z for j=2, xDot for j=3,
  1129.      * yDot for j=4, zDot for j=5).
  1130.      * </p>
  1131.      * @return 6x6 Jacobian matrix
  1132.      */
  1133.     private T[][] computeJacobianMeanWrtCartesianHyperbolic() {

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

  1135.         // compute various intermediate parameters
  1136.         computePVWithoutA();
  1137.         final FieldVector3D<T> position = partialPV.getPosition();
  1138.         final FieldVector3D<T> velocity = partialPV.getVelocity();
  1139.         final FieldVector3D<T> momentum = partialPV.getMomentum();
  1140.         final T r2         = position.getNormSq();
  1141.         final T r          = r2.sqrt();
  1142.         final T r3         = r.multiply(r2);

  1143.         final T x          = position.getX();
  1144.         final T y          = position.getY();
  1145.         final T z          = position.getZ();
  1146.         final T vx         = velocity.getX();
  1147.         final T vy         = velocity.getY();
  1148.         final T vz         = velocity.getZ();
  1149.         final T mx         = momentum.getX();
  1150.         final T my         = momentum.getY();
  1151.         final T mz         = momentum.getZ();

  1152.         final T mu         = getMu();
  1153.         final T absA       = a.negate();
  1154.         final T sqrtMuA    = absA.multiply(mu).sqrt();
  1155.         final T a2         = a.square();
  1156.         final T rOa        = r.divide(absA);

  1157.         final FieldSinCos<T> scI = FastMath.sinCos(i);
  1158.         final T cosI       = scI.cos();
  1159.         final T sinI       = scI.sin();

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

  1161.         // da
  1162.         final FieldVector3D<T> vectorAR = new FieldVector3D<>(a2.multiply(-2).divide(r3), position);
  1163.         final FieldVector3D<T> vectorARDot = velocity.scalarMultiply(a2.multiply(-2).divide(mu));
  1164.         fillHalfRow(getOne().negate(), vectorAR,    jacobian[0], 0);
  1165.         fillHalfRow(getOne().negate(), vectorARDot, jacobian[0], 3);

  1166.         // differentials of the momentum
  1167.         final T m      = momentum.getNorm();
  1168.         final T oOm    = m.reciprocal();
  1169.         final FieldVector3D<T> dcXP = new FieldVector3D<>(getZero(), vz, vy.negate());
  1170.         final FieldVector3D<T> dcYP = new FieldVector3D<>(vz.negate(),   getZero(),  vx);
  1171.         final FieldVector3D<T> dcZP = new FieldVector3D<>( vy, vx.negate(),   getZero());
  1172.         final FieldVector3D<T> dcXV = new FieldVector3D<>(  getZero(),  z.negate(),   y);
  1173.         final FieldVector3D<T> dcYV = new FieldVector3D<>(  z,   getZero(),  x.negate());
  1174.         final FieldVector3D<T> dcZV = new FieldVector3D<>( y.negate(),   x,   getZero());
  1175.         final FieldVector3D<T> dCP  = new FieldVector3D<>(mx.multiply(oOm), dcXP, my.multiply(oOm), dcYP, mz.multiply(oOm), dcZP);
  1176.         final FieldVector3D<T> dCV  = new FieldVector3D<>(mx.multiply(oOm), dcXV, my.multiply(oOm), dcYV, mz.multiply(oOm), dcZV);

  1177.         // dp
  1178.         final T mOMu   = m.divide(mu);
  1179.         final FieldVector3D<T> dpP  = new FieldVector3D<>(mOMu.multiply(2), dCP);
  1180.         final FieldVector3D<T> dpV  = new FieldVector3D<>(mOMu.multiply(2), dCV);

  1181.         // de
  1182.         final T p      = m.multiply(mOMu);
  1183.         final T moO2ae = absA.multiply(2).multiply(e).reciprocal();
  1184.         final T m2OaMu = p.negate().divide(absA);
  1185.         fillHalfRow(moO2ae, dpP, m2OaMu.multiply(moO2ae), vectorAR,    jacobian[1], 0);
  1186.         fillHalfRow(moO2ae, dpV, m2OaMu.multiply(moO2ae), vectorARDot, jacobian[1], 3);

  1187.         // di
  1188.         final T cI1 = m.multiply(sinI).reciprocal();
  1189.         final T cI2 = cosI.multiply(cI1);
  1190.         fillHalfRow(cI2, dCP, cI1.negate(), dcZP, jacobian[2], 0);
  1191.         fillHalfRow(cI2, dCV, cI1.negate(), dcZV, jacobian[2], 3);


  1192.         // dPA
  1193.         final T cP1     =  y.multiply(oOm);
  1194.         final T cP2     =  x.negate().multiply(oOm);
  1195.         final T cP3     =  mx.multiply(cP1).add(my.multiply(cP2)).negate();
  1196.         final T cP4     =  cP3.multiply(oOm);
  1197.         final T cP5     =  r2.multiply(sinI).multiply(sinI).negate().reciprocal();
  1198.         final T cP6     = z.multiply(cP5);
  1199.         final T cP7     = cP3.multiply(cP5);
  1200.         final FieldVector3D<T> dacP  = new FieldVector3D<>(cP1, dcXP, cP2, dcYP, cP4, dCP, oOm, new FieldVector3D<>(my.negate(), mx, getZero()));
  1201.         final FieldVector3D<T> dacV  = new FieldVector3D<>(cP1, dcXV, cP2, dcYV, cP4, dCV);
  1202.         final FieldVector3D<T> dpoP  = new FieldVector3D<>(cP6, dacP, cP7, this.PLUS_K);
  1203.         final FieldVector3D<T> dpoV  = new FieldVector3D<>(cP6, dacV);

  1204.         final T re2     = r2.multiply(e.square());
  1205.         final T recOre2 = p.subtract(r).divide(re2);
  1206.         final T resOre2 = pv.multiply(mOMu).divide(re2);
  1207.         final FieldVector3D<T> dreP  = new FieldVector3D<>(mOMu, velocity, pv.divide(mu), dCP);
  1208.         final FieldVector3D<T> dreV  = new FieldVector3D<>(mOMu, position, pv.divide(mu), dCV);
  1209.         final FieldVector3D<T> davP  = new FieldVector3D<>(resOre2.negate(), dpP, recOre2, dreP, resOre2.divide(r), position);
  1210.         final FieldVector3D<T> davV  = new FieldVector3D<>(resOre2.negate(), dpV, recOre2, dreV);
  1211.         fillHalfRow(getOne(), dpoP, getOne().negate(), davP, jacobian[3], 0);
  1212.         fillHalfRow(getOne(), dpoV, getOne().negate(), davV, jacobian[3], 3);

  1213.         // dRAAN
  1214.         final T cO0 = cI1.square();
  1215.         final T cO1 =  mx.multiply(cO0);
  1216.         final T cO2 =  my.negate().multiply(cO0);
  1217.         fillHalfRow(cO1, dcYP, cO2, dcXP, jacobian[4], 0);
  1218.         fillHalfRow(cO1, dcYV, cO2, dcXV, jacobian[4], 3);

  1219.         // dM
  1220.         final T s2a    = pv.divide(absA.multiply(2));
  1221.         final T oObux  = m.square().add(absA.multiply(mu)).sqrt().reciprocal();
  1222.         final T scasbu = pv.multiply(oObux);
  1223.         final FieldVector3D<T> dauP = new FieldVector3D<>(sqrtMuA.reciprocal(), velocity, s2a.negate().divide(sqrtMuA), vectorAR);
  1224.         final FieldVector3D<T> dauV = new FieldVector3D<>(sqrtMuA.reciprocal(), position, s2a.negate().divide(sqrtMuA), vectorARDot);
  1225.         final FieldVector3D<T> dbuP = new FieldVector3D<>(oObux.multiply(mu.divide(2.)), vectorAR,    m.multiply(oObux), dCP);
  1226.         final FieldVector3D<T> dbuV = new FieldVector3D<>(oObux.multiply(mu.divide(2.)), vectorARDot, m.multiply(oObux), dCV);
  1227.         final FieldVector3D<T> dcuP = new FieldVector3D<>(oObux, velocity, scasbu.negate().multiply(oObux), dbuP);
  1228.         final FieldVector3D<T> dcuV = new FieldVector3D<>(oObux, position, scasbu.negate().multiply(oObux), dbuV);
  1229.         fillHalfRow(getOne(), dauP, e.negate().divide(rOa.add(1)), dcuP, jacobian[5], 0);
  1230.         fillHalfRow(getOne(), dauV, e.negate().divide(rOa.add(1)), dcuV, jacobian[5], 3);

  1231.         return jacobian;

  1232.     }

  1233.     /** {@inheritDoc} */
  1234.     @Override
  1235.     protected T[][] computeJacobianEccentricWrtCartesian() {
  1236.         if (isElliptical()) {
  1237.             return computeJacobianEccentricWrtCartesianElliptical();
  1238.         } else {
  1239.             return computeJacobianEccentricWrtCartesianHyperbolic();
  1240.         }
  1241.     }

  1242.     /** Compute the Jacobian of the orbital parameters with respect to the Cartesian parameters.
  1243.      * <p>
  1244.      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
  1245.      * respect to Cartesian coordinate j (x for j=0, y for j=1, z for j=2, xDot for j=3,
  1246.      * yDot for j=4, zDot for j=5).
  1247.      * </p>
  1248.      * @return 6x6 Jacobian matrix
  1249.      */
  1250.     private T[][] computeJacobianEccentricWrtCartesianElliptical() {

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

  1253.         // Differentiating the Kepler equation M = E - e sin E leads to:
  1254.         // dM = (1 - e cos E) dE - sin E de
  1255.         // which is inverted and rewritten as:
  1256.         // dE = a/r dM + sin E a/r de
  1257.         final FieldSinCos<T> scE = FastMath.sinCos(getEccentricAnomaly());
  1258.         final T aOr              = e.negate().multiply(scE.cos()).add(1).reciprocal();

  1259.         // update anomaly row
  1260.         final T[] eRow           = jacobian[1];
  1261.         final T[] anomalyRow     = jacobian[5];
  1262.         for (int j = 0; j < anomalyRow.length; ++j) {
  1263.             anomalyRow[j] = aOr.multiply(anomalyRow[j].add(scE.sin().multiply(eRow[j])));
  1264.         }

  1265.         return jacobian;

  1266.     }

  1267.     /** Compute the Jacobian of the orbital parameters with respect to the Cartesian parameters.
  1268.      * <p>
  1269.      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
  1270.      * respect to Cartesian coordinate j (x for j=0, y for j=1, z for j=2, xDot for j=3,
  1271.      * yDot for j=4, zDot for j=5).
  1272.      * </p>
  1273.      * @return 6x6 Jacobian matrix
  1274.      */
  1275.     private T[][] computeJacobianEccentricWrtCartesianHyperbolic() {

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

  1278.         // Differentiating the Kepler equation M = e sinh H - H leads to:
  1279.         // dM = (e cosh H - 1) dH + sinh H de
  1280.         // which is inverted and rewritten as:
  1281.         // dH = 1 / (e cosh H - 1) dM - sinh H / (e cosh H - 1) de
  1282.         final T H      = getEccentricAnomaly();
  1283.         final T coshH  = H.cosh();
  1284.         final T sinhH  = H.sinh();
  1285.         final T absaOr = e.multiply(coshH).subtract(1).reciprocal();
  1286.         // update anomaly row
  1287.         final T[] eRow       = jacobian[1];
  1288.         final T[] anomalyRow = jacobian[5];

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

  1291.         }

  1292.         return jacobian;

  1293.     }

  1294.     /** {@inheritDoc} */
  1295.     @Override
  1296.     protected T[][] computeJacobianTrueWrtCartesian() {
  1297.         if (isElliptical()) {
  1298.             return computeJacobianTrueWrtCartesianElliptical();
  1299.         } else {
  1300.             return computeJacobianTrueWrtCartesianHyperbolic();
  1301.         }
  1302.     }

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

  1312.         // start by computing the Jacobian with eccentric angle
  1313.         final T[][] jacobian = computeJacobianEccentricWrtCartesianElliptical();
  1314.         // Differentiating the eccentric anomaly equation sin E = sqrt(1-e^2) sin v / (1 + e cos v)
  1315.         // and using cos E = (e + cos v) / (1 + e cos v) to get rid of cos E leads to:
  1316.         // dE = [sqrt (1 - e^2) / (1 + e cos v)] dv - [sin E / (1 - e^2)] de
  1317.         // which is inverted and rewritten as:
  1318.         // dv = sqrt (1 - e^2) a/r dE + [sin E / sqrt (1 - e^2)] a/r de
  1319.         final T e2               = e.square();
  1320.         final T oMe2             = e2.negate().add(1);
  1321.         final T epsilon          = oMe2.sqrt();
  1322.         final FieldSinCos<T> scE = FastMath.sinCos(getEccentricAnomaly());
  1323.         final T aOr              = e.multiply(scE.cos()).negate().add(1).reciprocal();
  1324.         final T aFactor          = epsilon.multiply(aOr);
  1325.         final T eFactor          = scE.sin().multiply(aOr).divide(epsilon);

  1326.         // update anomaly row
  1327.         final T[] eRow           = jacobian[1];
  1328.         final T[] anomalyRow     = jacobian[5];
  1329.         for (int j = 0; j < anomalyRow.length; ++j) {
  1330.             anomalyRow[j] = aFactor.multiply(anomalyRow[j]).add(eFactor.multiply(eRow[j]));
  1331.         }
  1332.         return jacobian;

  1333.     }

  1334.     /** Compute the Jacobian of the orbital parameters with respect to the Cartesian parameters.
  1335.      * <p>
  1336.      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
  1337.      * respect to Cartesian coordinate j (x for j=0, y for j=1, z for j=2, xDot for j=3,
  1338.      * yDot for j=4, zDot for j=5).
  1339.      * </p>
  1340.      * @return 6x6 Jacobian matrix
  1341.      */
  1342.     private T[][] computeJacobianTrueWrtCartesianHyperbolic() {

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

  1345.         // Differentiating the eccentric anomaly equation sinh H = sqrt(e^2-1) sin v / (1 + e cos v)
  1346.         // and using cosh H = (e + cos v) / (1 + e cos v) to get rid of cosh H leads to:
  1347.         // dH = [sqrt (e^2 - 1) / (1 + e cos v)] dv + [sinh H / (e^2 - 1)] de
  1348.         // which is inverted and rewritten as:
  1349.         // dv = sqrt (1 - e^2) a/r dH - [sinh H / sqrt (e^2 - 1)] a/r de
  1350.         final T e2       = e.square();
  1351.         final T e2Mo     = e2.subtract(1);
  1352.         final T epsilon  = e2Mo.sqrt();
  1353.         final T H        = getEccentricAnomaly();
  1354.         final T coshH    = H.cosh();
  1355.         final T sinhH    = H.sinh();
  1356.         final T aOr      = e.multiply(coshH).subtract(1).reciprocal();
  1357.         final T aFactor  = epsilon.multiply(aOr);
  1358.         final T eFactor  = sinhH.multiply(aOr).divide(epsilon);

  1359.         // update anomaly row
  1360.         final T[] eRow           = jacobian[1];
  1361.         final T[] anomalyRow     = jacobian[5];
  1362.         for (int j = 0; j < anomalyRow.length; ++j) {
  1363.             anomalyRow[j] = aFactor.multiply(anomalyRow[j]).subtract(eFactor.multiply(eRow[j]));
  1364.         }

  1365.         return jacobian;

  1366.     }

  1367.     /** {@inheritDoc} */
  1368.     @Override
  1369.     public void addKeplerContribution(final PositionAngleType type, final T gm,
  1370.                                       final T[] pDot) {
  1371.         pDot[5] = pDot[5].add(computeKeplerianAnomalyDot(type, a, e, gm, cachedAnomaly, cachedPositionAngleType));
  1372.     }

  1373.     /**
  1374.      * Compute rate of argument of latitude.
  1375.      * @param type position angle type of rate
  1376.      * @param a semi major axis
  1377.      * @param e eccentricity
  1378.      * @param mu mu
  1379.      * @param anomaly anomaly
  1380.      * @param cachedType position angle type of passed anomaly
  1381.      * @param <T> field type
  1382.      * @return first-order time derivative for anomaly
  1383.      * @since 12.2
  1384.      */
  1385.     private static <T extends CalculusFieldElement<T>> T computeKeplerianAnomalyDot(final PositionAngleType type, final T a, final T e,
  1386.                                                                                     final T mu, final T anomaly, final PositionAngleType cachedType) {
  1387.         final T absA = a.abs();
  1388.         final T n    = absA.reciprocal().multiply(mu).sqrt().divide(absA);
  1389.         if (type == PositionAngleType.MEAN) {
  1390.             return n;
  1391.         }
  1392.         final T ksi;
  1393.         final T oMe2;
  1394.         final T trueAnomaly = FieldKeplerianAnomalyUtility.convertAnomaly(cachedType, anomaly, e, PositionAngleType.TRUE);
  1395.         if (type == PositionAngleType.ECCENTRIC) {
  1396.             oMe2 = e.square().negate().add(1).abs();
  1397.             ksi = e.multiply(trueAnomaly.cos()).add(1);
  1398.             return n.multiply(ksi).divide(oMe2);
  1399.         } else {
  1400.             oMe2 = e.square().negate().add(1).abs();
  1401.             ksi  = e.multiply(trueAnomaly.cos()).add(1);
  1402.             return n.multiply(ksi).multiply(ksi).divide(oMe2.multiply(oMe2.sqrt()));
  1403.         }
  1404.     }

  1405.     /**  Returns a string representation of this Keplerian parameters object.
  1406.      * @return a string representation of this object
  1407.      */
  1408.     public String toString() {
  1409.         return new StringBuilder().append("Keplerian parameters: ").append('{').
  1410.                                   append("a: ").append(a.getReal()).
  1411.                                   append("; e: ").append(e.getReal()).
  1412.                                   append("; i: ").append(FastMath.toDegrees(i.getReal())).
  1413.                                   append("; pa: ").append(FastMath.toDegrees(pa.getReal())).
  1414.                                   append("; raan: ").append(FastMath.toDegrees(raan.getReal())).
  1415.                                   append("; v: ").append(FastMath.toDegrees(getTrueAnomaly().getReal())).
  1416.                                   append(";}").toString();
  1417.     }

  1418.     /** {@inheritDoc} */
  1419.     @Override
  1420.     public PositionAngleType getCachedPositionAngleType() {
  1421.         return cachedPositionAngleType;
  1422.     }

  1423.     /** {@inheritDoc} */
  1424.     @Override
  1425.     public boolean hasNonKeplerianRates() {
  1426.         return hasNonKeplerianAcceleration();
  1427.     }

  1428.     /** {@inheritDoc} */
  1429.     @Override
  1430.     public FieldKeplerianOrbit<T> withKeplerianRates() {
  1431.         return new FieldKeplerianOrbit<>(getA(), getE(), getI(), getPerigeeArgument(), getRightAscensionOfAscendingNode(),
  1432.                 cachedAnomaly, cachedPositionAngleType, getFrame(), getDate(), getMu());
  1433.     }

  1434.     /** Check if the given parameter is within an acceptable range.
  1435.      * The bounds are inclusive: an exception is raised when either of those conditions are met:
  1436.      * <ul>
  1437.      *     <li>The parameter is strictly greater than upperBound</li>
  1438.      *     <li>The parameter is strictly lower than lowerBound</li>
  1439.      * </ul>
  1440.      * <p>
  1441.      * In either of these cases, an OrekitException is raised.
  1442.      * </p>
  1443.      * @param parameterName name of the parameter
  1444.      * @param parameter value of the parameter
  1445.      * @param lowerBound lower bound of the acceptable range (inclusive)
  1446.      * @param upperBound upper bound of the acceptable range (inclusive)
  1447.      */
  1448.     private void checkParameterRangeInclusive(final String parameterName, final double parameter,
  1449.                                               final double lowerBound, final double upperBound) {
  1450.         if (parameter < lowerBound || parameter > upperBound) {
  1451.             throw new OrekitException(OrekitMessages.INVALID_PARAMETER_RANGE, parameterName,
  1452.                                       parameter, lowerBound, upperBound);
  1453.         }
  1454.     }

  1455.     /** {@inheritDoc} */
  1456.     @Override
  1457.     public KeplerianOrbit toOrbit() {
  1458.         final double cachedPositionAngle = cachedAnomaly.getReal();
  1459.         return new KeplerianOrbit(a.getReal(), e.getReal(), i.getReal(),
  1460.                                   pa.getReal(), raan.getReal(), cachedPositionAngle,
  1461.                                   aDot.getReal(), eDot.getReal(), iDot.getReal(),
  1462.                                   paDot.getReal(), raanDot.getReal(),
  1463.                                   cachedAnomalyDot.getReal(),
  1464.                                   cachedPositionAngleType, cachedPositionAngleType,
  1465.                                   getFrame(), getDate().toAbsoluteDate(), getMu().getReal());
  1466.     }


  1467. }