CartesianOrbit.java

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

  18. import java.io.Serializable;
  19. import java.util.stream.Stream;

  20. import org.hipparchus.analysis.differentiation.UnivariateDerivative2;
  21. import org.hipparchus.exception.LocalizedCoreFormats;
  22. import org.hipparchus.exception.MathIllegalStateException;
  23. import org.hipparchus.geometry.euclidean.threed.FieldVector3D;
  24. import org.hipparchus.geometry.euclidean.threed.Rotation;
  25. import org.hipparchus.geometry.euclidean.threed.RotationConvention;
  26. import org.hipparchus.geometry.euclidean.threed.Vector3D;
  27. import org.hipparchus.util.FastMath;
  28. import org.hipparchus.util.SinCos;
  29. import org.orekit.annotation.DefaultDataContext;
  30. import org.orekit.data.DataContext;
  31. import org.orekit.errors.OrekitMessages;
  32. import org.orekit.frames.Frame;
  33. import org.orekit.time.AbsoluteDate;
  34. import org.orekit.utils.CartesianDerivativesFilter;
  35. import org.orekit.utils.FieldPVCoordinates;
  36. import org.orekit.utils.PVCoordinates;
  37. import org.orekit.utils.TimeStampedPVCoordinates;


  38. /** This class holds Cartesian orbital parameters.

  39.  * <p>
  40.  * The parameters used internally are the Cartesian coordinates:
  41.  *   <ul>
  42.  *     <li>x</li>
  43.  *     <li>y</li>
  44.  *     <li>z</li>
  45.  *     <li>xDot</li>
  46.  *     <li>yDot</li>
  47.  *     <li>zDot</li>
  48.  *   </ul>
  49.  * contained in {@link PVCoordinates}.
  50.  *

  51.  * <p>
  52.  * Note that the implementation of this class delegates all non-Cartesian related
  53.  * computations ({@link #getA()}, {@link #getEquinoctialEx()}, ...) to an underlying
  54.  * instance of the {@link EquinoctialOrbit} class. This implies that using this class
  55.  * only for analytical computations which are always based on non-Cartesian
  56.  * parameters is perfectly possible but somewhat sub-optimal.
  57.  * </p>
  58.  * <p>
  59.  * The instance <code>CartesianOrbit</code> is guaranteed to be immutable.
  60.  * </p>
  61.  * @see    Orbit
  62.  * @see    KeplerianOrbit
  63.  * @see    CircularOrbit
  64.  * @see    EquinoctialOrbit
  65.  * @author Luc Maisonobe
  66.  * @author Guylaine Prat
  67.  * @author Fabien Maussion
  68.  * @author V&eacute;ronique Pommier-Maurussane
  69.  */
  70. public class CartesianOrbit extends Orbit {

  71.     /** Serializable UID. */
  72.     private static final long serialVersionUID = 20170414L;

  73.     /** Indicator for non-Keplerian derivatives. */
  74.     private final transient boolean hasNonKeplerianAcceleration;

  75.     /** Underlying equinoctial orbit to which high-level methods are delegated. */
  76.     private transient EquinoctialOrbit equinoctial;

  77.     /** Constructor from Cartesian parameters.
  78.      *
  79.      * <p> The acceleration provided in {@code pvCoordinates} is accessible using
  80.      * {@link #getPVCoordinates()} and {@link #getPVCoordinates(Frame)}. All other methods
  81.      * use {@code mu} and the position to compute the acceleration, including
  82.      * {@link #shiftedBy(double)} and {@link #getPVCoordinates(AbsoluteDate, Frame)}.
  83.      *
  84.      * @param pvaCoordinates the position, velocity and acceleration of the satellite.
  85.      * @param frame the frame in which the {@link PVCoordinates} are defined
  86.      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
  87.      * @param mu central attraction coefficient (m³/s²)
  88.      * @exception IllegalArgumentException if frame is not a {@link
  89.      * Frame#isPseudoInertial pseudo-inertial frame}
  90.      */
  91.     public CartesianOrbit(final TimeStampedPVCoordinates pvaCoordinates,
  92.                           final Frame frame, final double mu)
  93.         throws IllegalArgumentException {
  94.         super(pvaCoordinates, frame, mu);
  95.         hasNonKeplerianAcceleration = hasNonKeplerianAcceleration(pvaCoordinates, mu);
  96.         equinoctial = null;
  97.     }

  98.     /** Constructor from Cartesian parameters.
  99.      *
  100.      * <p> The acceleration provided in {@code pvCoordinates} is accessible using
  101.      * {@link #getPVCoordinates()} and {@link #getPVCoordinates(Frame)}. All other methods
  102.      * use {@code mu} and the position to compute the acceleration, including
  103.      * {@link #shiftedBy(double)} and {@link #getPVCoordinates(AbsoluteDate, Frame)}.
  104.      *
  105.      * @param pvaCoordinates the position and velocity of the satellite.
  106.      * @param frame the frame in which the {@link PVCoordinates} are defined
  107.      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
  108.      * @param date date of the orbital parameters
  109.      * @param mu central attraction coefficient (m³/s²)
  110.      * @exception IllegalArgumentException if frame is not a {@link
  111.      * Frame#isPseudoInertial pseudo-inertial frame}
  112.      */
  113.     public CartesianOrbit(final PVCoordinates pvaCoordinates, final Frame frame,
  114.                           final AbsoluteDate date, final double mu)
  115.         throws IllegalArgumentException {
  116.         this(new TimeStampedPVCoordinates(date, pvaCoordinates), frame, mu);
  117.     }

  118.     /** Constructor from any kind of orbital parameters.
  119.      * @param op orbital parameters to copy
  120.      */
  121.     public CartesianOrbit(final Orbit op) {
  122.         super(op.getPVCoordinates(), op.getFrame(), op.getMu());
  123.         hasNonKeplerianAcceleration = op.hasDerivatives();
  124.         if (op instanceof EquinoctialOrbit) {
  125.             equinoctial = (EquinoctialOrbit) op;
  126.         } else if (op instanceof CartesianOrbit) {
  127.             equinoctial = ((CartesianOrbit) op).equinoctial;
  128.         } else {
  129.             equinoctial = null;
  130.         }
  131.     }

  132.     /** {@inheritDoc} */
  133.     public OrbitType getType() {
  134.         return OrbitType.CARTESIAN;
  135.     }

  136.     /** Lazy evaluation of equinoctial parameters. */
  137.     private void initEquinoctial() {
  138.         if (equinoctial == null) {
  139.             if (hasDerivatives()) {
  140.                 // getPVCoordinates includes accelerations that will be interpreted as derivatives
  141.                 equinoctial = new EquinoctialOrbit(getPVCoordinates(), getFrame(), getDate(), getMu());
  142.             } else {
  143.                 // get rid of Keplerian acceleration so we don't assume
  144.                 // we have derivatives when in fact we don't have them
  145.                 equinoctial = new EquinoctialOrbit(new PVCoordinates(getPVCoordinates().getPosition(),
  146.                                                                      getPVCoordinates().getVelocity()),
  147.                                                    getFrame(), getDate(), getMu());
  148.             }
  149.         }
  150.     }

  151.     /** Get the position/velocity with derivatives.
  152.      * @return position/velocity with derivatives
  153.      * @since 10.2
  154.      */
  155.     private FieldPVCoordinates<UnivariateDerivative2> getPVDerivatives() {
  156.         // PVA coordinates
  157.         final PVCoordinates pva = getPVCoordinates();
  158.         final Vector3D      p   = pva.getPosition();
  159.         final Vector3D      v   = pva.getVelocity();
  160.         final Vector3D      a   = pva.getAcceleration();
  161.         // Field coordinates
  162.         final FieldVector3D<UnivariateDerivative2> pG = new FieldVector3D<>(new UnivariateDerivative2(p.getX(), v.getX(), a.getX()),
  163.                                                                new UnivariateDerivative2(p.getY(), v.getY(), a.getY()),
  164.                                                                new UnivariateDerivative2(p.getZ(), v.getZ(), a.getZ()));
  165.         final FieldVector3D<UnivariateDerivative2> vG = new FieldVector3D<>(new UnivariateDerivative2(v.getX(), a.getX(), 0.0),
  166.                                                                new UnivariateDerivative2(v.getY(), a.getY(), 0.0),
  167.                                                                new UnivariateDerivative2(v.getZ(), a.getZ(), 0.0));
  168.         return new FieldPVCoordinates<>(pG, vG);
  169.     }

  170.     /** {@inheritDoc} */
  171.     public double getA() {
  172.         final double r  = getPVCoordinates().getPosition().getNorm();
  173.         final double V2 = getPVCoordinates().getVelocity().getNormSq();
  174.         return r / (2 - r * V2 / getMu());
  175.     }

  176.     /** {@inheritDoc} */
  177.     public double getADot() {
  178.         if (hasDerivatives()) {
  179.             final FieldPVCoordinates<UnivariateDerivative2> pv = getPVDerivatives();
  180.             final UnivariateDerivative2 r  = pv.getPosition().getNorm();
  181.             final UnivariateDerivative2 V2 = pv.getVelocity().getNormSq();
  182.             final UnivariateDerivative2 a  = r.divide(r.multiply(V2).divide(getMu()).subtract(2).negate());
  183.             return a.getDerivative(1);
  184.         } else {
  185.             return Double.NaN;
  186.         }
  187.     }

  188.     /** {@inheritDoc} */
  189.     public double getE() {
  190.         final double muA = getMu() * getA();
  191.         if (muA > 0) {
  192.             // elliptic or circular orbit
  193.             final Vector3D pvP   = getPVCoordinates().getPosition();
  194.             final Vector3D pvV   = getPVCoordinates().getVelocity();
  195.             final double rV2OnMu = pvP.getNorm() * pvV.getNormSq() / getMu();
  196.             final double eSE     = Vector3D.dotProduct(pvP, pvV) / FastMath.sqrt(muA);
  197.             final double eCE     = rV2OnMu - 1;
  198.             return FastMath.sqrt(eCE * eCE + eSE * eSE);
  199.         } else {
  200.             // hyperbolic orbit
  201.             final Vector3D pvM = getPVCoordinates().getMomentum();
  202.             return FastMath.sqrt(1 - pvM.getNormSq() / muA);
  203.         }
  204.     }

  205.     /** {@inheritDoc} */
  206.     public double getEDot() {
  207.         if (hasDerivatives()) {
  208.             final FieldPVCoordinates<UnivariateDerivative2> pv = getPVDerivatives();
  209.             final FieldVector3D<UnivariateDerivative2> pvP   = pv.getPosition();
  210.             final FieldVector3D<UnivariateDerivative2> pvV   = pv.getVelocity();
  211.             final UnivariateDerivative2 r       = pvP.getNorm();
  212.             final UnivariateDerivative2 V2      = pvV.getNormSq();
  213.             final UnivariateDerivative2 rV2OnMu = r.multiply(V2).divide(getMu());
  214.             final UnivariateDerivative2 a       = r.divide(rV2OnMu.negate().add(2));
  215.             final UnivariateDerivative2 eSE     = FieldVector3D.dotProduct(pvP, pvV).divide(a.multiply(getMu()).sqrt());
  216.             final UnivariateDerivative2 eCE     = rV2OnMu.subtract(1);
  217.             final UnivariateDerivative2 e       = eCE.multiply(eCE).add(eSE.multiply(eSE)).sqrt();
  218.             return e.getDerivative(1);
  219.         } else {
  220.             return Double.NaN;
  221.         }
  222.     }

  223.     /** {@inheritDoc} */
  224.     public double getI() {
  225.         return Vector3D.angle(Vector3D.PLUS_K, getPVCoordinates().getMomentum());
  226.     }

  227.     /** {@inheritDoc} */
  228.     public double getIDot() {
  229.         if (hasDerivatives()) {
  230.             final FieldPVCoordinates<UnivariateDerivative2> pv = getPVDerivatives();
  231.             final FieldVector3D<UnivariateDerivative2> momentum =
  232.                             FieldVector3D.crossProduct(pv.getPosition(), pv.getVelocity());
  233.             final UnivariateDerivative2 i = FieldVector3D.angle(Vector3D.PLUS_K, momentum);
  234.             return i.getDerivative(1);
  235.         } else {
  236.             return Double.NaN;
  237.         }
  238.     }

  239.     /** {@inheritDoc} */
  240.     public double getEquinoctialEx() {
  241.         initEquinoctial();
  242.         return equinoctial.getEquinoctialEx();
  243.     }

  244.     /** {@inheritDoc} */
  245.     public double getEquinoctialExDot() {
  246.         initEquinoctial();
  247.         return equinoctial.getEquinoctialExDot();
  248.     }

  249.     /** {@inheritDoc} */
  250.     public double getEquinoctialEy() {
  251.         initEquinoctial();
  252.         return equinoctial.getEquinoctialEy();
  253.     }

  254.     /** {@inheritDoc} */
  255.     public double getEquinoctialEyDot() {
  256.         initEquinoctial();
  257.         return equinoctial.getEquinoctialEyDot();
  258.     }

  259.     /** {@inheritDoc} */
  260.     public double getHx() {
  261.         final Vector3D w = getPVCoordinates().getMomentum().normalize();
  262.         // Check for equatorial retrograde orbit
  263.         if ((w.getX() * w.getX() + w.getY() * w.getY()) == 0 && w.getZ() < 0) {
  264.             return Double.NaN;
  265.         }
  266.         return -w.getY() / (1 + w.getZ());
  267.     }

  268.     /** {@inheritDoc} */
  269.     public double getHxDot() {
  270.         if (hasDerivatives()) {
  271.             final FieldPVCoordinates<UnivariateDerivative2> pv = getPVDerivatives();
  272.             final FieldVector3D<UnivariateDerivative2> w =
  273.                             FieldVector3D.crossProduct(pv.getPosition(), pv.getVelocity()).normalize();
  274.             // Check for equatorial retrograde orbit
  275.             final double x = w.getX().getValue();
  276.             final double y = w.getY().getValue();
  277.             final double z = w.getZ().getValue();
  278.             if ((x * x + y * y) == 0 && z < 0) {
  279.                 return Double.NaN;
  280.             }
  281.             final UnivariateDerivative2 hx = w.getY().negate().divide(w.getZ().add(1));
  282.             return hx.getDerivative(1);
  283.         } else {
  284.             return Double.NaN;
  285.         }
  286.     }

  287.     /** {@inheritDoc} */
  288.     public double getHy() {
  289.         final Vector3D w = getPVCoordinates().getMomentum().normalize();
  290.         // Check for equatorial retrograde orbit
  291.         if ((w.getX() * w.getX() + w.getY() * w.getY()) == 0 && w.getZ() < 0) {
  292.             return Double.NaN;
  293.         }
  294.         return  w.getX() / (1 + w.getZ());
  295.     }

  296.     /** {@inheritDoc} */
  297.     public double getHyDot() {
  298.         if (hasDerivatives()) {
  299.             final FieldPVCoordinates<UnivariateDerivative2> pv = getPVDerivatives();
  300.             final FieldVector3D<UnivariateDerivative2> w =
  301.                             FieldVector3D.crossProduct(pv.getPosition(), pv.getVelocity()).normalize();
  302.             // Check for equatorial retrograde orbit
  303.             final double x = w.getX().getValue();
  304.             final double y = w.getY().getValue();
  305.             final double z = w.getZ().getValue();
  306.             if ((x * x + y * y) == 0 && z < 0) {
  307.                 return Double.NaN;
  308.             }
  309.             final UnivariateDerivative2 hy = w.getX().divide(w.getZ().add(1));
  310.             return hy.getDerivative(1);
  311.         } else {
  312.             return Double.NaN;
  313.         }
  314.     }

  315.     /** {@inheritDoc} */
  316.     public double getLv() {
  317.         initEquinoctial();
  318.         return equinoctial.getLv();
  319.     }

  320.     /** {@inheritDoc} */
  321.     public double getLvDot() {
  322.         initEquinoctial();
  323.         return equinoctial.getLvDot();
  324.     }

  325.     /** {@inheritDoc} */
  326.     public double getLE() {
  327.         initEquinoctial();
  328.         return equinoctial.getLE();
  329.     }

  330.     /** {@inheritDoc} */
  331.     public double getLEDot() {
  332.         initEquinoctial();
  333.         return equinoctial.getLEDot();
  334.     }

  335.     /** {@inheritDoc} */
  336.     public double getLM() {
  337.         initEquinoctial();
  338.         return equinoctial.getLM();
  339.     }

  340.     /** {@inheritDoc} */
  341.     public double getLMDot() {
  342.         initEquinoctial();
  343.         return equinoctial.getLMDot();
  344.     }

  345.     /** {@inheritDoc} */
  346.     public boolean hasDerivatives() {
  347.         return hasNonKeplerianAcceleration;
  348.     }

  349.     /** {@inheritDoc} */
  350.     protected TimeStampedPVCoordinates initPVCoordinates() {
  351.         // nothing to do here, as the canonical elements are already the Cartesian ones
  352.         return getPVCoordinates();
  353.     }

  354.     /** {@inheritDoc} */
  355.     public CartesianOrbit shiftedBy(final double dt) {
  356.         final PVCoordinates shiftedPV = (getA() < 0) ? shiftPVHyperbolic(dt) : shiftPVElliptic(dt);
  357.         return new CartesianOrbit(shiftedPV, getFrame(), getDate().shiftedBy(dt), getMu());
  358.     }

  359.     /** {@inheritDoc}
  360.      * <p>
  361.      * The interpolated instance is created by polynomial Hermite interpolation
  362.      * ensuring velocity remains the exact derivative of position.
  363.      * </p>
  364.      * <p>
  365.      * As this implementation of interpolation is polynomial, it should be used only
  366.      * with small samples (about 10-20 points) in order to avoid <a
  367.      * href="http://en.wikipedia.org/wiki/Runge%27s_phenomenon">Runge's phenomenon</a>
  368.      * and numerical problems (including NaN appearing).
  369.      * </p>
  370.      * <p>
  371.      * If orbit interpolation on large samples is needed, using the {@link
  372.      * org.orekit.propagation.analytical.Ephemeris} class is a better way than using this
  373.      * low-level interpolation. The Ephemeris class automatically handles selection of
  374.      * a neighboring sub-sample with a predefined number of point from a large global sample
  375.      * in a thread-safe way.
  376.      * </p>
  377.      */
  378.     public CartesianOrbit interpolate(final AbsoluteDate date, final Stream<Orbit> sample) {
  379.         final TimeStampedPVCoordinates interpolated =
  380.                 TimeStampedPVCoordinates.interpolate(date, CartesianDerivativesFilter.USE_PVA,
  381.                                                      sample.map(orbit -> orbit.getPVCoordinates()));
  382.         return new CartesianOrbit(interpolated, getFrame(), date, getMu());
  383.     }

  384.     /** Compute shifted position and velocity in elliptic case.
  385.      * @param dt time shift
  386.      * @return shifted position and velocity
  387.      */
  388.     private PVCoordinates shiftPVElliptic(final double dt) {

  389.         // preliminary computation
  390.         final Vector3D pvP   = getPVCoordinates().getPosition();
  391.         final Vector3D pvV   = getPVCoordinates().getVelocity();
  392.         final double r2      = pvP.getNormSq();
  393.         final double r       = FastMath.sqrt(r2);
  394.         final double rV2OnMu = r * pvV.getNormSq() / getMu();
  395.         final double a       = getA();
  396.         final double eSE     = Vector3D.dotProduct(pvP, pvV) / FastMath.sqrt(getMu() * a);
  397.         final double eCE     = rV2OnMu - 1;
  398.         final double e2      = eCE * eCE + eSE * eSE;

  399.         // we can use any arbitrary reference 2D frame in the orbital plane
  400.         // in order to simplify some equations below, we use the current position as the u axis
  401.         final Vector3D u     = pvP.normalize();
  402.         final Vector3D v     = Vector3D.crossProduct(getPVCoordinates().getMomentum(), u).normalize();

  403.         // the following equations rely on the specific choice of u explained above,
  404.         // some coefficients that vanish to 0 in this case have already been removed here
  405.         final double ex        = (eCE - e2) * a / r;
  406.         final double ey        = -FastMath.sqrt(1 - e2) * eSE * a / r;
  407.         final double beta      = 1 / (1 + FastMath.sqrt(1 - e2));
  408.         final double thetaE0   = FastMath.atan2(ey + eSE * beta * ex, r / a + ex - eSE * beta * ey);
  409.         final SinCos scThetaE0 = FastMath.sinCos(thetaE0);
  410.         final double thetaM0   = thetaE0 - ex * scThetaE0.sin() + ey * scThetaE0.cos();

  411.         // compute in-plane shifted eccentric argument
  412.         final double thetaM1 = thetaM0 + getKeplerianMeanMotion() * dt;
  413.         final double thetaE1 = meanToEccentric(thetaM1, ex, ey);
  414.         final SinCos scTE    = FastMath.sinCos(thetaE1);
  415.         final double cTE     = scTE.cos();
  416.         final double sTE     = scTE.sin();

  417.         // compute shifted in-plane Cartesian coordinates
  418.         final double exey   = ex * ey;
  419.         final double exCeyS = ex * cTE + ey * sTE;
  420.         final double x      = a * ((1 - beta * ey * ey) * cTE + beta * exey * sTE - ex);
  421.         final double y      = a * ((1 - beta * ex * ex) * sTE + beta * exey * cTE - ey);
  422.         final double factor = FastMath.sqrt(getMu() / a) / (1 - exCeyS);
  423.         final double xDot   = factor * (-sTE + beta * ey * exCeyS);
  424.         final double yDot   = factor * ( cTE - beta * ex * exCeyS);

  425.         final Vector3D shiftedP = new Vector3D(x, u, y, v);
  426.         final Vector3D shiftedV = new Vector3D(xDot, u, yDot, v);
  427.         if (hasNonKeplerianAcceleration) {

  428.             // extract non-Keplerian part of the initial acceleration
  429.             final Vector3D nonKeplerianAcceleration = new Vector3D(1, getPVCoordinates().getAcceleration(),
  430.                                                                    getMu() / (r2 * r), pvP);

  431.             // add the quadratic motion due to the non-Keplerian acceleration to the Keplerian motion
  432.             final Vector3D fixedP   = new Vector3D(1, shiftedP,
  433.                                                    0.5 * dt * dt, nonKeplerianAcceleration);
  434.             final double   fixedR2 = fixedP.getNormSq();
  435.             final double   fixedR  = FastMath.sqrt(fixedR2);
  436.             final Vector3D fixedV  = new Vector3D(1, shiftedV,
  437.                                                   dt, nonKeplerianAcceleration);
  438.             final Vector3D fixedA  = new Vector3D(-getMu() / (fixedR2 * fixedR), shiftedP,
  439.                                                   1, nonKeplerianAcceleration);

  440.             return new PVCoordinates(fixedP, fixedV, fixedA);

  441.         } else {
  442.             // don't include acceleration,
  443.             // so the shifted orbit is not considered to have derivatives
  444.             return new PVCoordinates(shiftedP, shiftedV);
  445.         }

  446.     }

  447.     /** Compute shifted position and velocity in hyperbolic case.
  448.      * @param dt time shift
  449.      * @return shifted position and velocity
  450.      */
  451.     private PVCoordinates shiftPVHyperbolic(final double dt) {

  452.         final PVCoordinates pv = getPVCoordinates();
  453.         final Vector3D pvP   = pv.getPosition();
  454.         final Vector3D pvV   = pv.getVelocity();
  455.         final Vector3D pvM   = pv.getMomentum();
  456.         final double r2      = pvP.getNormSq();
  457.         final double r       = FastMath.sqrt(r2);
  458.         final double rV2OnMu = r * pvV.getNormSq() / getMu();
  459.         final double a       = getA();
  460.         final double muA     = getMu() * a;
  461.         final double e       = FastMath.sqrt(1 - Vector3D.dotProduct(pvM, pvM) / muA);
  462.         final double sqrt    = FastMath.sqrt((e + 1) / (e - 1));

  463.         // compute mean anomaly
  464.         final double eSH     = Vector3D.dotProduct(pvP, pvV) / FastMath.sqrt(-muA);
  465.         final double eCH     = rV2OnMu - 1;
  466.         final double H0      = FastMath.log((eCH + eSH) / (eCH - eSH)) / 2;
  467.         final double M0      = e * FastMath.sinh(H0) - H0;

  468.         // find canonical 2D frame with p pointing to perigee
  469.         final double v0      = 2 * FastMath.atan(sqrt * FastMath.tanh(H0 / 2));
  470.         final Vector3D p     = new Rotation(pvM, v0, RotationConvention.FRAME_TRANSFORM).applyTo(pvP).normalize();
  471.         final Vector3D q     = Vector3D.crossProduct(pvM, p).normalize();

  472.         // compute shifted eccentric anomaly
  473.         final double M1      = M0 + getKeplerianMeanMotion() * dt;
  474.         final double H1      = meanToHyperbolicEccentric(M1, e);

  475.         // compute shifted in-plane Cartesian coordinates
  476.         final double cH     = FastMath.cosh(H1);
  477.         final double sH     = FastMath.sinh(H1);
  478.         final double sE2m1  = FastMath.sqrt((e - 1) * (e + 1));

  479.         // coordinates of position and velocity in the orbital plane
  480.         final double x      = a * (cH - e);
  481.         final double y      = -a * sE2m1 * sH;
  482.         final double factor = FastMath.sqrt(getMu() / -a) / (e * cH - 1);
  483.         final double xDot   = -factor * sH;
  484.         final double yDot   =  factor * sE2m1 * cH;

  485.         final Vector3D shiftedP = new Vector3D(x, p, y, q);
  486.         final Vector3D shiftedV = new Vector3D(xDot, p, yDot, q);
  487.         if (hasNonKeplerianAcceleration) {

  488.             // extract non-Keplerian part of the initial acceleration
  489.             final Vector3D nonKeplerianAcceleration = new Vector3D(1, getPVCoordinates().getAcceleration(),
  490.                                                                    getMu() / (r2 * r), pvP);

  491.             // add the quadratic motion due to the non-Keplerian acceleration to the Keplerian motion
  492.             final Vector3D fixedP   = new Vector3D(1, shiftedP,
  493.                                                    0.5 * dt * dt, nonKeplerianAcceleration);
  494.             final double   fixedR2 = fixedP.getNormSq();
  495.             final double   fixedR  = FastMath.sqrt(fixedR2);
  496.             final Vector3D fixedV  = new Vector3D(1, shiftedV,
  497.                                                   dt, nonKeplerianAcceleration);
  498.             final Vector3D fixedA  = new Vector3D(-getMu() / (fixedR2 * fixedR), shiftedP,
  499.                                                   1, nonKeplerianAcceleration);

  500.             return new PVCoordinates(fixedP, fixedV, fixedA);

  501.         } else {
  502.             // don't include acceleration,
  503.             // so the shifted orbit is not considered to have derivatives
  504.             return new PVCoordinates(shiftedP, shiftedV);
  505.         }

  506.     }

  507.     /** Computes the eccentric in-plane argument from the mean in-plane argument.
  508.      * @param thetaM = mean in-plane argument (rad)
  509.      * @param ex first component of eccentricity vector
  510.      * @param ey second component of eccentricity vector
  511.      * @return the eccentric in-plane argument.
  512.      */
  513.     private double meanToEccentric(final double thetaM, final double ex, final double ey) {
  514.         // Generalization of Kepler equation to in-plane parameters
  515.         // with thetaE = eta + E and
  516.         //      thetaM = eta + M = thetaE - ex.sin(thetaE) + ey.cos(thetaE)
  517.         // and eta being counted from an arbitrary reference in the orbital plane
  518.         double thetaE        = thetaM;
  519.         double thetaEMthetaM = 0.0;
  520.         int    iter          = 0;
  521.         do {
  522.             final SinCos scThetaE = FastMath.sinCos(thetaE);

  523.             final double f2 = ex * scThetaE.sin() - ey * scThetaE.cos();
  524.             final double f1 = 1.0 - ex * scThetaE.cos() - ey * scThetaE.sin();
  525.             final double f0 = thetaEMthetaM - f2;

  526.             final double f12 = 2.0 * f1;
  527.             final double shift = f0 * f12 / (f1 * f12 - f0 * f2);

  528.             thetaEMthetaM -= shift;
  529.             thetaE         = thetaM + thetaEMthetaM;

  530.             if (FastMath.abs(shift) <= 1.0e-12) {
  531.                 return thetaE;
  532.             }

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

  534.         throw new MathIllegalStateException(LocalizedCoreFormats.CONVERGENCE_FAILED);

  535.     }

  536.     /** Computes the hyperbolic eccentric anomaly from the mean anomaly.
  537.      * <p>
  538.      * The algorithm used here for solving hyperbolic Kepler equation is
  539.      * Danby's iterative method (3rd order) with Vallado's initial guess.
  540.      * </p>
  541.      * @param M mean anomaly (rad)
  542.      * @param ecc eccentricity
  543.      * @return the hyperbolic eccentric anomaly
  544.      */
  545.     private double meanToHyperbolicEccentric(final double M, final double ecc) {

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

  547.         // Initial guess
  548.         double H;
  549.         if (ecc < 1.6) {
  550.             if (-FastMath.PI < M && M < 0. || M > FastMath.PI) {
  551.                 H = M - ecc;
  552.             } else {
  553.                 H = M + ecc;
  554.             }
  555.         } else {
  556.             if (ecc < 3.6 && FastMath.abs(M) > FastMath.PI) {
  557.                 H = M - FastMath.copySign(ecc, M);
  558.             } else {
  559.                 H = M / (ecc - 1.);
  560.             }
  561.         }

  562.         // Iterative computation
  563.         int iter = 0;
  564.         do {
  565.             final double f3  = ecc * FastMath.cosh(H);
  566.             final double f2  = ecc * FastMath.sinh(H);
  567.             final double f1  = f3 - 1.;
  568.             final double f0  = f2 - H - M;
  569.             final double f12 = 2. * f1;
  570.             final double d   = f0 / f12;
  571.             final double fdf = f1 - d * f2;
  572.             final double ds  = f0 / fdf;

  573.             final double shift = f0 / (fdf + ds * ds * f3 / 6.);

  574.             H -= shift;

  575.             if (FastMath.abs(shift) <= 1.0e-12) {
  576.                 return H;
  577.             }

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

  579.         throw new MathIllegalStateException(OrekitMessages.UNABLE_TO_COMPUTE_HYPERBOLIC_ECCENTRIC_ANOMALY,
  580.                                             iter);
  581.     }

  582.     /** Create a 6x6 identity matrix.
  583.      * @return 6x6 identity matrix
  584.      */
  585.     private double[][] create6x6Identity() {
  586.         // this is the fastest way to set the 6x6 identity matrix
  587.         final double[][] identity = new double[6][6];
  588.         for (int i = 0; i < 6; i++) {
  589.             identity[i][i] = 1.0;
  590.         }
  591.         return identity;
  592.     }

  593.     @Override
  594.     protected double[][] computeJacobianMeanWrtCartesian() {
  595.         return create6x6Identity();
  596.     }

  597.     @Override
  598.     protected double[][] computeJacobianEccentricWrtCartesian() {
  599.         return create6x6Identity();
  600.     }

  601.     @Override
  602.     protected double[][] computeJacobianTrueWrtCartesian() {
  603.         return create6x6Identity();
  604.     }

  605.     /** {@inheritDoc} */
  606.     public void addKeplerContribution(final PositionAngle type, final double gm,
  607.                                       final double[] pDot) {

  608.         final PVCoordinates pv = getPVCoordinates();

  609.         // position derivative is velocity
  610.         final Vector3D velocity = pv.getVelocity();
  611.         pDot[0] += velocity.getX();
  612.         pDot[1] += velocity.getY();
  613.         pDot[2] += velocity.getZ();

  614.         // velocity derivative is Newtonian acceleration
  615.         final Vector3D position = pv.getPosition();
  616.         final double r2         = position.getNormSq();
  617.         final double coeff      = -gm / (r2 * FastMath.sqrt(r2));
  618.         pDot[3] += coeff * position.getX();
  619.         pDot[4] += coeff * position.getY();
  620.         pDot[5] += coeff * position.getZ();

  621.     }

  622.     /**  Returns a string representation of this Orbit object.
  623.      * @return a string representation of this object
  624.      */
  625.     public String toString() {
  626.         // use only the six defining elements, like the other Orbit.toString() methods
  627.         final String comma = ", ";
  628.         final PVCoordinates pv = getPVCoordinates();
  629.         final Vector3D p = pv.getPosition();
  630.         final Vector3D v = pv.getVelocity();
  631.         return "Cartesian parameters: {P(" +
  632.                 p.getX() + comma +
  633.                 p.getY() + comma +
  634.                 p.getZ() + "), V(" +
  635.                 v.getX() + comma +
  636.                 v.getY() + comma +
  637.                 v.getZ() + ")}";
  638.     }

  639.     /** Replace the instance with a data transfer object for serialization.
  640.      * <p>
  641.      * This intermediate class serializes all needed parameters,
  642.      * including position-velocity which are <em>not</em> serialized by parent
  643.      * {@link Orbit} class.
  644.      * </p>
  645.      * @return data transfer object that will be serialized
  646.      */
  647.     @DefaultDataContext
  648.     private Object writeReplace() {
  649.         return new DTO(this);
  650.     }

  651.     /** Internal class used only for serialization. */
  652.     @DefaultDataContext
  653.     private static class DTO implements Serializable {

  654.         /** Serializable UID. */
  655.         private static final long serialVersionUID = 20170414L;

  656.         /** Double values. */
  657.         private double[] d;

  658.         /** Frame in which are defined the orbital parameters. */
  659.         private final Frame frame;

  660.         /** Simple constructor.
  661.          * @param orbit instance to serialize
  662.          */
  663.         private DTO(final CartesianOrbit orbit) {

  664.             final TimeStampedPVCoordinates pv = orbit.getPVCoordinates();

  665.             // decompose date
  666.             final AbsoluteDate j2000Epoch =
  667.                     DataContext.getDefault().getTimeScales().getJ2000Epoch();
  668.             final double epoch  = FastMath.floor(pv.getDate().durationFrom(j2000Epoch));
  669.             final double offset = pv.getDate().durationFrom(j2000Epoch.shiftedBy(epoch));

  670.             if (orbit.hasDerivatives()) {
  671.                 this.d = new double[] {
  672.                     epoch, offset, orbit.getMu(),
  673.                     pv.getPosition().getX(),     pv.getPosition().getY(),     pv.getPosition().getZ(),
  674.                     pv.getVelocity().getX(),     pv.getVelocity().getY(),     pv.getVelocity().getZ(),
  675.                     pv.getAcceleration().getX(), pv.getAcceleration().getY(), pv.getAcceleration().getZ()
  676.                 };
  677.             } else {
  678.                 this.d = new double[] {
  679.                     epoch, offset, orbit.getMu(),
  680.                     pv.getPosition().getX(),     pv.getPosition().getY(),     pv.getPosition().getZ(),
  681.                     pv.getVelocity().getX(),     pv.getVelocity().getY(),     pv.getVelocity().getZ()
  682.                 };
  683.             }

  684.             this.frame = orbit.getFrame();

  685.         }

  686.         /** Replace the deserialized data transfer object with a {@link CartesianOrbit}.
  687.          * @return replacement {@link CartesianOrbit}
  688.          */
  689.         private Object readResolve() {
  690.             final AbsoluteDate j2000Epoch =
  691.                     DataContext.getDefault().getTimeScales().getJ2000Epoch();
  692.             if (d.length >= 12) {
  693.                 // we have derivatives
  694.                 return new CartesianOrbit(new TimeStampedPVCoordinates(j2000Epoch.shiftedBy(d[0]).shiftedBy(d[1]),
  695.                                                                        new Vector3D(d[3], d[ 4], d[ 5]),
  696.                                                                        new Vector3D(d[6], d[ 7], d[ 8]),
  697.                                                                        new Vector3D(d[9], d[10], d[11])),
  698.                                           frame, d[2]);
  699.             } else {
  700.                 // we don't have derivatives
  701.                 return new CartesianOrbit(new TimeStampedPVCoordinates(j2000Epoch.shiftedBy(d[0]).shiftedBy(d[1]),
  702.                                                                        new Vector3D(d[3], d[ 4], d[ 5]),
  703.                                                                        new Vector3D(d[6], d[ 7], d[ 8])),
  704.                                           frame, d[2]);
  705.             }
  706.         }

  707.     }

  708. }