CartesianOrbit.java

  1. /* Copyright 2002-2016 CS Systèmes d'Information
  2.  * Licensed to CS Systèmes d'Information (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.ArrayList;
  20. import java.util.Collection;
  21. import java.util.List;

  22. import org.apache.commons.math3.exception.ConvergenceException;
  23. import org.apache.commons.math3.geometry.euclidean.threed.Rotation;
  24. import org.apache.commons.math3.geometry.euclidean.threed.RotationConvention;
  25. import org.apache.commons.math3.geometry.euclidean.threed.Vector3D;
  26. import org.apache.commons.math3.util.FastMath;
  27. import org.orekit.errors.OrekitMessages;
  28. import org.orekit.frames.Frame;
  29. import org.orekit.time.AbsoluteDate;
  30. import org.orekit.utils.CartesianDerivativesFilter;
  31. import org.orekit.utils.PVCoordinates;
  32. import org.orekit.utils.TimeStampedPVCoordinates;


  33. /** This class holds cartesian orbital parameters.

  34.  * <p>
  35.  * The parameters used internally are the cartesian coordinates:
  36.  *   <ul>
  37.  *     <li>x</li>
  38.  *     <li>y</li>
  39.  *     <li>z</li>
  40.  *     <li>xDot</li>
  41.  *     <li>yDot</li>
  42.  *     <li>zDot</li>
  43.  *   </ul>
  44.  * contained in {@link PVCoordinates}.
  45.  * </p>

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

  66.     /** Serializable UID. */
  67.     private static final long serialVersionUID = 20140723L;

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

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

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

  110.     /** Constructor from any kind of orbital parameters.
  111.      * @param op orbital parameters to copy
  112.      */
  113.     public CartesianOrbit(final Orbit op) {
  114.         super(op.getPVCoordinates(), op.getFrame(), op.getMu());
  115.         if (op instanceof EquinoctialOrbit) {
  116.             equinoctial = (EquinoctialOrbit) op;
  117.         } else if (op instanceof CartesianOrbit) {
  118.             equinoctial = ((CartesianOrbit) op).equinoctial;
  119.         } else {
  120.             equinoctial = null;
  121.         }
  122.     }

  123.     /** {@inheritDoc} */
  124.     public OrbitType getType() {
  125.         return OrbitType.CARTESIAN;
  126.     }

  127.     /** Lazy evaluation of equinoctial parameters. */
  128.     private void initEquinoctial() {
  129.         if (equinoctial == null) {
  130.             equinoctial = new EquinoctialOrbit(getPVCoordinates(), getFrame(), getDate(), getMu());
  131.         }
  132.     }

  133.     /** {@inheritDoc} */
  134.     public double getA() {
  135.         // lazy evaluation of semi-major axis
  136.         final double r  = getPVCoordinates().getPosition().getNorm();
  137.         final double V2 = getPVCoordinates().getVelocity().getNormSq();
  138.         return r / (2 - r * V2 / getMu());
  139.     }

  140.     /** {@inheritDoc} */
  141.     public double getE() {
  142.         final Vector3D pvP   = getPVCoordinates().getPosition();
  143.         final Vector3D pvV   = getPVCoordinates().getVelocity();
  144.         final double rV2OnMu = pvP.getNorm() * pvV.getNormSq() / getMu();
  145.         final double eSE     = Vector3D.dotProduct(pvP, pvV) / FastMath.sqrt(getMu() * getA());
  146.         final double eCE     = rV2OnMu - 1;
  147.         return FastMath.sqrt(eCE * eCE + eSE * eSE);
  148.     }

  149.     /** {@inheritDoc} */
  150.     public double getI() {
  151.         return Vector3D.angle(Vector3D.PLUS_K, getPVCoordinates().getMomentum());
  152.     }

  153.     /** {@inheritDoc} */
  154.     public double getEquinoctialEx() {
  155.         initEquinoctial();
  156.         return equinoctial.getEquinoctialEx();
  157.     }

  158.     /** {@inheritDoc} */
  159.     public double getEquinoctialEy() {
  160.         initEquinoctial();
  161.         return equinoctial.getEquinoctialEy();
  162.     }

  163.     /** {@inheritDoc} */
  164.     public double getHx() {
  165.         final Vector3D w = getPVCoordinates().getMomentum().normalize();
  166.         // Check for equatorial retrograde orbit
  167.         if (((w.getX() * w.getX() + w.getY() * w.getY()) == 0) && w.getZ() < 0) {
  168.             return Double.NaN;
  169.         }
  170.         return -w.getY() / (1 + w.getZ());
  171.     }

  172.     /** {@inheritDoc} */
  173.     public double getHy() {
  174.         final Vector3D w = getPVCoordinates().getMomentum().normalize();
  175.         // Check for equatorial retrograde orbit
  176.         if (((w.getX() * w.getX() + w.getY() * w.getY()) == 0) && w.getZ() < 0) {
  177.             return Double.NaN;
  178.         }
  179.         return  w.getX() / (1 + w.getZ());
  180.     }

  181.     /** {@inheritDoc} */
  182.     public double getLv() {
  183.         initEquinoctial();
  184.         return equinoctial.getLv();
  185.     }

  186.     /** {@inheritDoc} */
  187.     public double getLE() {
  188.         initEquinoctial();
  189.         return equinoctial.getLE();
  190.     }

  191.     /** {@inheritDoc} */
  192.     public double getLM() {
  193.         initEquinoctial();
  194.         return equinoctial.getLM();
  195.     }

  196.     /** {@inheritDoc} */
  197.     protected TimeStampedPVCoordinates initPVCoordinates() {
  198.         // nothing to do here, as the canonical elements are already the cartesian ones
  199.         return getPVCoordinates();
  200.     }

  201.     /** {@inheritDoc} */
  202.     public CartesianOrbit shiftedBy(final double dt) {
  203.         final PVCoordinates shiftedPV = (getA() < 0) ? shiftPVHyperbolic(dt) : shiftPVElliptic(dt);
  204.         return new CartesianOrbit(shiftedPV, getFrame(), getDate().shiftedBy(dt), getMu());
  205.     }

  206.     /** {@inheritDoc}
  207.      * <p>
  208.      * The interpolated instance is created by polynomial Hermite interpolation
  209.      * ensuring velocity remains the exact derivative of position.
  210.      * </p>
  211.      * <p>
  212.      * As this implementation of interpolation is polynomial, it should be used only
  213.      * with small samples (about 10-20 points) in order to avoid <a
  214.      * href="http://en.wikipedia.org/wiki/Runge%27s_phenomenon">Runge's phenomenon</a>
  215.      * and numerical problems (including NaN appearing).
  216.      * </p>
  217.      * <p>
  218.      * If orbit interpolation on large samples is needed, using the {@link
  219.      * org.orekit.propagation.analytical.Ephemeris} class is a better way than using this
  220.      * low-level interpolation. The Ephemeris class automatically handles selection of
  221.      * a neighboring sub-sample with a predefined number of point from a large global sample
  222.      * in a thread-safe way.
  223.      * </p>
  224.      */
  225.     public CartesianOrbit interpolate(final AbsoluteDate date, final Collection<Orbit> sample) {
  226.         final List<TimeStampedPVCoordinates> datedPV = new ArrayList<TimeStampedPVCoordinates>(sample.size());
  227.         for (final Orbit o : sample) {
  228.             datedPV.add(new TimeStampedPVCoordinates(o.getDate(),
  229.                                                      o.getPVCoordinates().getPosition(),
  230.                                                      o.getPVCoordinates().getVelocity(),
  231.                                                      o.getPVCoordinates().getAcceleration()));
  232.         }
  233.         final TimeStampedPVCoordinates interpolated =
  234.                 TimeStampedPVCoordinates.interpolate(date, CartesianDerivativesFilter.USE_PVA, datedPV);
  235.         return new CartesianOrbit(interpolated, getFrame(), date, getMu());
  236.     }

  237.     /** Compute shifted position and velocity in elliptic case.
  238.      * @param dt time shift
  239.      * @return shifted position and velocity
  240.      */
  241.     private PVCoordinates shiftPVElliptic(final double dt) {

  242.         // preliminary computation
  243.         final Vector3D pvP   = getPVCoordinates().getPosition();
  244.         final Vector3D pvV   = getPVCoordinates().getVelocity();
  245.         final double r       = pvP.getNorm();
  246.         final double rV2OnMu = r * pvV.getNormSq() / getMu();
  247.         final double a       = getA();
  248.         final double eSE     = Vector3D.dotProduct(pvP, pvV) / FastMath.sqrt(getMu() * a);
  249.         final double eCE     = rV2OnMu - 1;
  250.         final double e2      = eCE * eCE + eSE * eSE;

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

  255.         // the following equations rely on the specific choice of u explained above,
  256.         // some coefficients that vanish to 0 in this case have already been removed here
  257.         final double ex      = (eCE - e2) * a / r;
  258.         final double ey      = -FastMath.sqrt(1 - e2) * eSE * a / r;
  259.         final double beta    = 1 / (1 + FastMath.sqrt(1 - e2));
  260.         final double thetaE0 = FastMath.atan2(ey + eSE * beta * ex, r / a + ex - eSE * beta * ey);
  261.         final double thetaM0 = thetaE0 - ex * FastMath.sin(thetaE0) + ey * FastMath.cos(thetaE0);

  262.         // compute in-plane shifted eccentric argument
  263.         final double thetaM1 = thetaM0 + getKeplerianMeanMotion() * dt;
  264.         final double thetaE1 = meanToEccentric(thetaM1, ex, ey);
  265.         final double cTE     = FastMath.cos(thetaE1);
  266.         final double sTE     = FastMath.sin(thetaE1);

  267.         // compute shifted in-plane cartesian coordinates
  268.         final double exey   = ex * ey;
  269.         final double exCeyS = ex * cTE + ey * sTE;
  270.         final double x      = a * ((1 - beta * ey * ey) * cTE + beta * exey * sTE - ex);
  271.         final double y      = a * ((1 - beta * ex * ex) * sTE + beta * exey * cTE - ey);
  272.         final double factor = FastMath.sqrt(getMu() / a) / (1 - exCeyS);
  273.         final double xDot   = factor * (-sTE + beta * ey * exCeyS);
  274.         final double yDot   = factor * ( cTE - beta * ex * exCeyS);

  275.         final Vector3D shiftedP = new Vector3D(x, u, y, v);
  276.         final double   r2       = x * x + y * y;
  277.         final Vector3D shiftedV = new Vector3D(xDot, u, yDot, v);
  278.         final Vector3D shiftedA = new Vector3D(-getMu() / (r2 * FastMath.sqrt(r2)), shiftedP);
  279.         return new PVCoordinates(shiftedP, shiftedV, shiftedA);

  280.     }

  281.     /** Compute shifted position and velocity in hyperbolic case.
  282.      * @param dt time shift
  283.      * @return shifted position and velocity
  284.      */
  285.     private PVCoordinates shiftPVHyperbolic(final double dt) {

  286.         final PVCoordinates pv = getPVCoordinates();
  287.         final Vector3D pvP   = pv.getPosition();
  288.         final Vector3D pvV   = pv.getVelocity();
  289.         final Vector3D pvM   = pv.getMomentum();
  290.         final double r       = pvP.getNorm();
  291.         final double rV2OnMu = r * pvV.getNormSq() / getMu();
  292.         final double a       = getA();
  293.         final double muA     = getMu() * a;
  294.         final double e       = FastMath.sqrt(1 - Vector3D.dotProduct(pvM, pvM) / muA);
  295.         final double sqrt    = FastMath.sqrt((e + 1) / (e - 1));

  296.         // compute mean anomaly
  297.         final double eSH     = Vector3D.dotProduct(pvP, pvV) / FastMath.sqrt(-muA);
  298.         final double eCH     = rV2OnMu - 1;
  299.         final double H0      = FastMath.log((eCH + eSH) / (eCH - eSH)) / 2;
  300.         final double M0      = e * FastMath.sinh(H0) - H0;

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

  305.         // compute shifted eccentric anomaly
  306.         final double M1      = M0 + getKeplerianMeanMotion() * dt;
  307.         final double H1      = meanToHyperbolicEccentric(M1, e);

  308.         // compute shifted in-plane cartesian coordinates
  309.         final double cH     = FastMath.cosh(H1);
  310.         final double sH     = FastMath.sinh(H1);
  311.         final double sE2m1  = FastMath.sqrt((e - 1) * (e + 1));

  312.         // coordinates of position and velocity in the orbital plane
  313.         final double x      = a * (cH - e);
  314.         final double y      = -a * sE2m1 * sH;
  315.         final double factor = FastMath.sqrt(getMu() / -a) / (e * cH - 1);
  316.         final double xDot   = -factor * sH;
  317.         final double yDot   =  factor * sE2m1 * cH;

  318.         final Vector3D shiftedP = new Vector3D(x, p, y, q);
  319.         final double   r2       = x * x + y * y;
  320.         final Vector3D shiftedV = new Vector3D(xDot, p, yDot, q);
  321.         final Vector3D shiftedA = new Vector3D(-getMu() / (r2 * FastMath.sqrt(r2)), shiftedP);
  322.         return new PVCoordinates(shiftedP, shiftedV, shiftedA);

  323.     }

  324.     /** Computes the eccentric in-plane argument from the mean in-plane argument.
  325.      * @param thetaM = mean in-plane argument (rad)
  326.      * @param ex first component of eccentricity vector
  327.      * @param ey second component of eccentricity vector
  328.      * @return the eccentric in-plane argument.
  329.      */
  330.     private double meanToEccentric(final double thetaM, final double ex, final double ey) {
  331.         // Generalization of Kepler equation to in-plane parameters
  332.         // with thetaE = eta + E and
  333.         //      thetaM = eta + M = thetaE - ex.sin(thetaE) + ey.cos(thetaE)
  334.         // and eta being counted from an arbitrary reference in the orbital plane
  335.         double thetaE        = thetaM;
  336.         double thetaEMthetaM = 0.0;
  337.         int    iter          = 0;
  338.         do {
  339.             final double cosThetaE = FastMath.cos(thetaE);
  340.             final double sinThetaE = FastMath.sin(thetaE);

  341.             final double f2 = ex * sinThetaE - ey * cosThetaE;
  342.             final double f1 = 1.0 - ex * cosThetaE - ey * sinThetaE;
  343.             final double f0 = thetaEMthetaM - f2;

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

  346.             thetaEMthetaM -= shift;
  347.             thetaE         = thetaM + thetaEMthetaM;

  348.             if (FastMath.abs(shift) <= 1.0e-12) {
  349.                 return thetaE;
  350.             }

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

  352.         throw new ConvergenceException();

  353.     }

  354.     /** Computes the hyperbolic eccentric anomaly from the mean anomaly.
  355.      * <p>
  356.      * The algorithm used here for solving hyperbolic Kepler equation is
  357.      * Danby's iterative method (3rd order) with Vallado's initial guess.
  358.      * </p>
  359.      * @param M mean anomaly (rad)
  360.      * @param ecc eccentricity
  361.      * @return the hyperbolic eccentric anomaly
  362.      */
  363.     private double meanToHyperbolicEccentric(final double M, final double ecc) {

  364.         // Resolution of hyperbolic Kepler equation for keplerian parameters

  365.         // Initial guess
  366.         double H;
  367.         if (ecc < 1.6) {
  368.             if ((-FastMath.PI < M && M < 0.) || M > FastMath.PI) {
  369.                 H = M - ecc;
  370.             } else {
  371.                 H = M + ecc;
  372.             }
  373.         } else {
  374.             if (ecc < 3.6 && FastMath.abs(M) > FastMath.PI) {
  375.                 H = M - FastMath.copySign(ecc, M);
  376.             } else {
  377.                 H = M / (ecc - 1.);
  378.             }
  379.         }

  380.         // Iterative computation
  381.         int iter = 0;
  382.         do {
  383.             final double f3  = ecc * FastMath.cosh(H);
  384.             final double f2  = ecc * FastMath.sinh(H);
  385.             final double f1  = f3 - 1.;
  386.             final double f0  = f2 - H - M;
  387.             final double f12 = 2. * f1;
  388.             final double d   = f0 / f12;
  389.             final double fdf = f1 - d * f2;
  390.             final double ds  = f0 / fdf;

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

  392.             H -= shift;

  393.             if (FastMath.abs(shift) <= 1.0e-12) {
  394.                 return H;
  395.             }

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

  397.         throw new ConvergenceException(OrekitMessages.UNABLE_TO_COMPUTE_HYPERBOLIC_ECCENTRIC_ANOMALY,
  398.                                        iter);
  399.     }

  400.     @Override
  401.     public void getJacobianWrtCartesian(final PositionAngle type, final double[][] jacobian) {
  402.         // this is the fastest way to set the 6x6 identity matrix
  403.         for (int i = 0; i < 6; i++) {
  404.             for (int j = 0; j < 6; j++) {
  405.                 jacobian[i][j] = 0;
  406.             }
  407.             jacobian[i][i] = 1;
  408.         }
  409.     }

  410.     @Override
  411.     protected double[][] computeJacobianMeanWrtCartesian() {
  412.         // not used
  413.         return null;
  414.     }
  415.     @Override
  416.     protected double[][] computeJacobianEccentricWrtCartesian() {
  417.         // not used
  418.         return null;
  419.     }

  420.     @Override
  421.     protected double[][] computeJacobianTrueWrtCartesian() {
  422.         // not used
  423.         return null;
  424.     }

  425.     /** {@inheritDoc} */
  426.     public void addKeplerContribution(final PositionAngle type, final double gm,
  427.                                       final double[] pDot) {

  428.         final PVCoordinates pv = getPVCoordinates();

  429.         // position derivative is velocity
  430.         final Vector3D velocity = pv.getVelocity();
  431.         pDot[0] += velocity.getX();
  432.         pDot[1] += velocity.getY();
  433.         pDot[2] += velocity.getZ();

  434.         // velocity derivative is Newtonian acceleration
  435.         final Vector3D position = pv.getPosition();
  436.         final double r2         = position.getNormSq();
  437.         final double coeff      = -gm / (r2 * FastMath.sqrt(r2));
  438.         pDot[3] += coeff * position.getX();
  439.         pDot[4] += coeff * position.getY();
  440.         pDot[5] += coeff * position.getZ();

  441.     }

  442.     /**  Returns a string representation of this Orbit object.
  443.      * @return a string representation of this object
  444.      */
  445.     public String toString() {
  446.         return "cartesian parameters: " + getPVCoordinates().toString();
  447.     }

  448.     /** Replace the instance with a data transfer object for serialization.
  449.      * <p>
  450.      * This intermediate class serializes all needed parameters,
  451.      * including position-velocity which are <em>not</em> serialized by parent
  452.      * {@link Orbit} class.
  453.      * </p>
  454.      * @return data transfer object that will be serialized
  455.      */
  456.     private Object writeReplace() {
  457.         return new DTO(this);
  458.     }

  459.     /** Internal class used only for serialization. */
  460.     private static class DTO implements Serializable {

  461.         /** Serializable UID. */
  462.         private static final long serialVersionUID = 20140723L;

  463.         /** Double values. */
  464.         private double[] d;

  465.         /** Frame in which are defined the orbital parameters. */
  466.         private final Frame frame;

  467.         /** Simple constructor.
  468.          * @param orbit instance to serialize
  469.          */
  470.         private DTO(final CartesianOrbit orbit) {

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

  472.             // decompose date
  473.             final double epoch  = FastMath.floor(pv.getDate().durationFrom(AbsoluteDate.J2000_EPOCH));
  474.             final double offset = pv.getDate().durationFrom(AbsoluteDate.J2000_EPOCH.shiftedBy(epoch));

  475.             this.d = new double[] {
  476.                 epoch, offset, orbit.getMu(),
  477.                 pv.getPosition().getX(),     pv.getPosition().getY(),     pv.getPosition().getZ(),
  478.                 pv.getVelocity().getX(),     pv.getVelocity().getY(),     pv.getVelocity().getZ(),
  479.                 pv.getAcceleration().getX(), pv.getAcceleration().getY(), pv.getAcceleration().getZ(),
  480.             };

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

  482.         }

  483.         /** Replace the deserialized data transfer object with a {@link CartesianOrbit}.
  484.          * @return replacement {@link CartesianOrbit}
  485.          */
  486.         private Object readResolve() {
  487.             return new CartesianOrbit(new TimeStampedPVCoordinates(AbsoluteDate.J2000_EPOCH.shiftedBy(d[0]).shiftedBy(d[1]),
  488.                                                                    new Vector3D(d[3], d[ 4], d[ 5]),
  489.                                                                    new Vector3D(d[6], d[ 7], d[ 8]),
  490.                                                                    new Vector3D(d[9], d[10], d[11])),
  491.                                       frame, d[2]);
  492.         }

  493.     }

  494. }