FieldCartesianOrbit.java

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


  18. import java.util.Arrays;

  19. import org.hipparchus.CalculusFieldElement;
  20. import org.hipparchus.Field;
  21. import org.hipparchus.analysis.differentiation.FieldUnivariateDerivative2;
  22. import org.hipparchus.geometry.euclidean.threed.FieldVector3D;
  23. import org.hipparchus.geometry.euclidean.threed.Vector3D;
  24. import org.hipparchus.util.MathArrays;
  25. import org.orekit.frames.Frame;
  26. import org.orekit.time.AbsoluteDate;
  27. import org.orekit.time.FieldAbsoluteDate;
  28. import org.orekit.utils.FieldPVCoordinates;
  29. import org.orekit.utils.PVCoordinates;
  30. import org.orekit.utils.TimeStampedFieldPVCoordinates;


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

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

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

  66.     /** Indicator for non-Keplerian acceleration. */
  67.     private final transient boolean hasNonKeplerianAcceleration;

  68.     /** Underlying equinoctial orbit to which high-level methods are delegated. */
  69.     private transient FieldEquinoctialOrbit<T> 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(CalculusFieldElement)} and {@link #getPVCoordinates(FieldAbsoluteDate, 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 FieldCartesianOrbit(final TimeStampedFieldPVCoordinates<T> pvaCoordinates,
  85.                                final Frame frame, final T mu)
  86.         throws IllegalArgumentException {
  87.         super(pvaCoordinates, frame, mu);
  88.         hasNonKeplerianAcceleration = hasNonKeplerianAcceleration(pvaCoordinates, mu);
  89.         equinoctial = null;
  90.     }

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

  111.     /** Constructor from any kind of orbital parameters.
  112.      * @param op orbital parameters to copy
  113.      */
  114.     public FieldCartesianOrbit(final FieldOrbit<T> op) {
  115.         super(op.getPVCoordinates(), op.getFrame(), op.getMu());
  116.         hasNonKeplerianAcceleration = op.hasDerivatives();
  117.         if (op instanceof FieldEquinoctialOrbit) {
  118.             equinoctial = (FieldEquinoctialOrbit<T>) op;
  119.         } else if (op instanceof FieldCartesianOrbit) {
  120.             equinoctial = ((FieldCartesianOrbit<T>) op).equinoctial;
  121.         } else {
  122.             equinoctial = null;
  123.         }
  124.     }

  125.     /** Constructor from Field and CartesianOrbit.
  126.      * <p>Build a FieldCartesianOrbit from non-Field CartesianOrbit.</p>
  127.      * @param field CalculusField to base object on
  128.      * @param op non-field orbit with only "constant" terms
  129.      * @since 12.0
  130.      */
  131.     public FieldCartesianOrbit(final Field<T> field, final CartesianOrbit op) {
  132.         super(new TimeStampedFieldPVCoordinates<>(field, op.getPVCoordinates()), op.getFrame(),
  133.                 field.getZero().newInstance(op.getMu()));
  134.         hasNonKeplerianAcceleration = op.hasDerivatives();
  135.         if (op.isElliptical()) {
  136.             equinoctial = new FieldEquinoctialOrbit<>(field, new EquinoctialOrbit(op));
  137.         } else {
  138.             equinoctial = null;
  139.         }
  140.     }

  141.     /** Constructor from Field and Orbit.
  142.      * <p>Build a FieldCartesianOrbit from any non-Field Orbit.</p>
  143.      * @param field CalculusField to base object on
  144.      * @param op non-field orbit with only "constant" terms
  145.      * @since 12.0
  146.      */
  147.     public FieldCartesianOrbit(final Field<T> field, final Orbit op) {
  148.         this(field, new CartesianOrbit(op));
  149.     }

  150.     /** {@inheritDoc} */
  151.     public OrbitType getType() {
  152.         return OrbitType.CARTESIAN;
  153.     }

  154.     /** Lazy evaluation of equinoctial parameters. */
  155.     private void initEquinoctial() {
  156.         if (equinoctial == null) {
  157.             if (hasDerivatives()) {
  158.                 // getPVCoordinates includes accelerations that will be interpreted as derivatives
  159.                 equinoctial = new FieldEquinoctialOrbit<>(getPVCoordinates(), getFrame(), getDate(), getMu());
  160.             } else {
  161.                 // get rid of Keplerian acceleration so we don't assume
  162.                 // we have derivatives when in fact we don't have them
  163.                 final FieldPVCoordinates<T> pva = getPVCoordinates();
  164.                 equinoctial = new FieldEquinoctialOrbit<>(new FieldPVCoordinates<>(pva.getPosition(), pva.getVelocity()),
  165.                                                           getFrame(), getDate(), getMu());
  166.             }
  167.         }
  168.     }

  169.     /** Get the position/velocity with derivatives.
  170.      * @return position/velocity with derivatives
  171.      * @since 10.2
  172.      */
  173.     private FieldPVCoordinates<FieldUnivariateDerivative2<T>> getPVDerivatives() {
  174.         // PVA coordinates
  175.         final FieldPVCoordinates<T> pva = getPVCoordinates();
  176.         final FieldVector3D<T>      p   = pva.getPosition();
  177.         final FieldVector3D<T>      v   = pva.getVelocity();
  178.         final FieldVector3D<T>      a   = pva.getAcceleration();
  179.         // Field coordinates
  180.         final FieldVector3D<FieldUnivariateDerivative2<T>> pG = new FieldVector3D<>(new FieldUnivariateDerivative2<>(p.getX(), v.getX(), a.getX()),
  181.                                                              new FieldUnivariateDerivative2<>(p.getY(), v.getY(), a.getY()),
  182.                                                              new FieldUnivariateDerivative2<>(p.getZ(), v.getZ(), a.getZ()));
  183.         final FieldVector3D<FieldUnivariateDerivative2<T>> vG = new FieldVector3D<>(new FieldUnivariateDerivative2<>(v.getX(), a.getX(), getZero()),
  184.                                                              new FieldUnivariateDerivative2<>(v.getY(), a.getY(), getZero()),
  185.                                                              new FieldUnivariateDerivative2<>(v.getZ(), a.getZ(), getZero()));
  186.         return new FieldPVCoordinates<>(pG, vG);
  187.     }

  188.     /** {@inheritDoc} */
  189.     public T getA() {
  190.         // lazy evaluation of semi-major axis
  191.         final FieldPVCoordinates<T> pva = getPVCoordinates();
  192.         final T r  = pva.getPosition().getNorm();
  193.         final T V2 = pva.getVelocity().getNormSq();
  194.         return r.divide(r.negate().multiply(V2).divide(getMu()).add(2));
  195.     }

  196.     /** {@inheritDoc} */
  197.     public T getADot() {
  198.         if (hasDerivatives()) {
  199.             final FieldPVCoordinates<FieldUnivariateDerivative2<T>> pv = getPVDerivatives();
  200.             final FieldUnivariateDerivative2<T> r  = pv.getPosition().getNorm();
  201.             final FieldUnivariateDerivative2<T> V2 = pv.getVelocity().getNormSq();
  202.             final FieldUnivariateDerivative2<T> a  = r.divide(r.multiply(V2).divide(getMu()).subtract(2).negate());
  203.             return a.getDerivative(1);
  204.         } else {
  205.             return null;
  206.         }
  207.     }

  208.     /** {@inheritDoc} */
  209.     public T getE() {
  210.         final T muA = getA().multiply(getMu());
  211.         if (isElliptical()) {
  212.             // elliptic or circular orbit
  213.             final FieldPVCoordinates<T> pva = getPVCoordinates();
  214.             final FieldVector3D<T> pvP      = pva.getPosition();
  215.             final FieldVector3D<T> pvV      = pva.getVelocity();
  216.             final T rV2OnMu = pvP.getNorm().multiply(pvV.getNormSq()).divide(getMu());
  217.             final T eSE     = FieldVector3D.dotProduct(pvP, pvV).divide(muA.sqrt());
  218.             final T eCE     = rV2OnMu.subtract(1);
  219.             return (eCE.square().add(eSE.square())).sqrt();
  220.         } else {
  221.             // hyperbolic orbit
  222.             final FieldVector3D<T> pvM = getPVCoordinates().getMomentum();
  223.             return pvM.getNormSq().divide(muA).negate().add(1).sqrt();
  224.         }
  225.     }

  226.     /** {@inheritDoc} */
  227.     public T getEDot() {
  228.         if (hasDerivatives()) {
  229.             final FieldPVCoordinates<FieldUnivariateDerivative2<T>> pv = getPVDerivatives();
  230.             final FieldUnivariateDerivative2<T> r       = pv.getPosition().getNorm();
  231.             final FieldUnivariateDerivative2<T> V2      = pv.getVelocity().getNormSq();
  232.             final FieldUnivariateDerivative2<T> rV2OnMu = r.multiply(V2).divide(getMu());
  233.             final FieldUnivariateDerivative2<T> a       = r.divide(rV2OnMu.negate().add(2));
  234.             final FieldUnivariateDerivative2<T> eSE     = FieldVector3D.dotProduct(pv.getPosition(), pv.getVelocity()).divide(a.multiply(getMu()).sqrt());
  235.             final FieldUnivariateDerivative2<T> eCE     = rV2OnMu.subtract(1);
  236.             final FieldUnivariateDerivative2<T> e       = eCE.square().add(eSE.square()).sqrt();
  237.             return e.getDerivative(1);
  238.         } else {
  239.             return null;
  240.         }
  241.     }

  242.     /** {@inheritDoc} */
  243.     public T getI() {
  244.         return FieldVector3D.angle(new FieldVector3D<>(getZero(), getZero(), getOne()), getPVCoordinates().getMomentum());
  245.     }

  246.     /** {@inheritDoc} */
  247.     public T getIDot() {
  248.         if (hasDerivatives()) {
  249.             final FieldPVCoordinates<FieldUnivariateDerivative2<T>> pv = getPVDerivatives();
  250.             final FieldVector3D<FieldUnivariateDerivative2<T>> momentum =
  251.                             FieldVector3D.crossProduct(pv.getPosition(), pv.getVelocity());
  252.             final FieldUnivariateDerivative2<T> i = FieldVector3D.angle(Vector3D.PLUS_K, momentum);
  253.             return i.getDerivative(1);
  254.         } else {
  255.             return null;
  256.         }
  257.     }

  258.     /** {@inheritDoc} */
  259.     public T getEquinoctialEx() {
  260.         initEquinoctial();
  261.         return equinoctial.getEquinoctialEx();
  262.     }

  263.     /** {@inheritDoc} */
  264.     public T getEquinoctialExDot() {
  265.         initEquinoctial();
  266.         return equinoctial.getEquinoctialExDot();
  267.     }

  268.     /** {@inheritDoc} */
  269.     public T getEquinoctialEy() {
  270.         initEquinoctial();
  271.         return equinoctial.getEquinoctialEy();
  272.     }

  273.     /** {@inheritDoc} */
  274.     public T getEquinoctialEyDot() {
  275.         initEquinoctial();
  276.         return equinoctial.getEquinoctialEyDot();
  277.     }

  278.     /** {@inheritDoc} */
  279.     public T getHx() {
  280.         final FieldVector3D<T> w = getPVCoordinates().getMomentum().normalize();
  281.         // Check for equatorial retrograde orbit
  282.         final double x = w.getX().getReal();
  283.         final double y = w.getY().getReal();
  284.         final double z = w.getZ().getReal();
  285.         if ((x * x + y * y) == 0 && z < 0) {
  286.             return getZero().add(Double.NaN);
  287.         }
  288.         return w.getY().negate().divide(w.getZ().add(1));
  289.     }

  290.     /** {@inheritDoc} */
  291.     public T getHxDot() {
  292.         if (hasDerivatives()) {
  293.             final FieldPVCoordinates<FieldUnivariateDerivative2<T>> pv = getPVDerivatives();
  294.             final FieldVector3D<FieldUnivariateDerivative2<T>> w =
  295.                             FieldVector3D.crossProduct(pv.getPosition(), pv.getVelocity()).normalize();
  296.             // Check for equatorial retrograde orbit
  297.             final double x = w.getX().getValue().getReal();
  298.             final double y = w.getY().getValue().getReal();
  299.             final double z = w.getZ().getValue().getReal();
  300.             if ((x * x + y * y) == 0 && z < 0) {
  301.                 return getZero().add(Double.NaN);
  302.             }
  303.             final FieldUnivariateDerivative2<T> hx = w.getY().negate().divide(w.getZ().add(1));
  304.             return hx.getDerivative(1);
  305.         } else {
  306.             return null;
  307.         }
  308.     }

  309.     /** {@inheritDoc} */
  310.     public T getHy() {
  311.         final FieldVector3D<T> w = getPVCoordinates().getMomentum().normalize();
  312.         // Check for equatorial retrograde orbit
  313.         final double x = w.getX().getReal();
  314.         final double y = w.getY().getReal();
  315.         final double z = w.getZ().getReal();
  316.         if ((x * x + y * y) == 0 && z < 0) {
  317.             return getZero().add(Double.NaN);
  318.         }
  319.         return  w.getX().divide(w.getZ().add(1));
  320.     }

  321.     /** {@inheritDoc} */
  322.     public T getHyDot() {
  323.         if (hasDerivatives()) {
  324.             final FieldPVCoordinates<FieldUnivariateDerivative2<T>> pv = getPVDerivatives();
  325.             final FieldVector3D<FieldUnivariateDerivative2<T>> w =
  326.                             FieldVector3D.crossProduct(pv.getPosition(), pv.getVelocity()).normalize();
  327.             // Check for equatorial retrograde orbit
  328.             final double x = w.getX().getValue().getReal();
  329.             final double y = w.getY().getValue().getReal();
  330.             final double z = w.getZ().getValue().getReal();
  331.             if ((x * x + y * y) == 0 && z < 0) {
  332.                 return getZero().add(Double.NaN);
  333.             }
  334.             final FieldUnivariateDerivative2<T> hy = w.getX().divide(w.getZ().add(1));
  335.             return hy.getDerivative(1);
  336.         } else {
  337.             return null;
  338.         }
  339.     }

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

  345.     /** {@inheritDoc} */
  346.     public T getLvDot() {
  347.         initEquinoctial();
  348.         return equinoctial.getLvDot();
  349.     }

  350.     /** {@inheritDoc} */
  351.     public T getLE() {
  352.         initEquinoctial();
  353.         return equinoctial.getLE();
  354.     }

  355.     /** {@inheritDoc} */
  356.     public T getLEDot() {
  357.         initEquinoctial();
  358.         return equinoctial.getLEDot();
  359.     }

  360.     /** {@inheritDoc} */
  361.     public T getLM() {
  362.         initEquinoctial();
  363.         return equinoctial.getLM();
  364.     }

  365.     /** {@inheritDoc} */
  366.     public T getLMDot() {
  367.         initEquinoctial();
  368.         return equinoctial.getLMDot();
  369.     }

  370.     /** {@inheritDoc} */
  371.     public boolean hasDerivatives() {
  372.         return hasNonKeplerianAcceleration;
  373.     }

  374.     /** {@inheritDoc} */
  375.     protected FieldVector3D<T> initPosition() {
  376.         // nothing to do here, as the canonical elements are already the Cartesian ones
  377.         return getPVCoordinates().getPosition();
  378.     }

  379.     /** {@inheritDoc} */
  380.     protected TimeStampedFieldPVCoordinates<T> initPVCoordinates() {
  381.         // nothing to do here, as the canonical elements are already the Cartesian ones
  382.         return getPVCoordinates();
  383.     }

  384.     /** {@inheritDoc} */
  385.     public FieldCartesianOrbit<T> shiftedBy(final double dt) {
  386.         return shiftedBy(getZero().newInstance(dt));
  387.     }

  388.     /** {@inheritDoc} */
  389.     public FieldCartesianOrbit<T> shiftedBy(final T dt) {
  390.         final FieldPVCoordinates<T> shiftedPV = shiftPV(dt);
  391.         return new FieldCartesianOrbit<>(shiftedPV, getFrame(), getDate().shiftedBy(dt), getMu());
  392.     }

  393.     /** Compute shifted position and velocity.
  394.      * @param dt time shift
  395.      * @return shifted position and velocity
  396.      */
  397.     private FieldPVCoordinates<T> shiftPV(final T dt) {
  398.         final FieldPVCoordinates<T> pvCoordinates = getPVCoordinates();
  399.         final FieldPVCoordinates<T> shiftedPV = KeplerianMotionCartesianUtility.predictPositionVelocity(dt,
  400.                 pvCoordinates.getPosition(), pvCoordinates.getVelocity(), getMu());

  401.         if (hasNonKeplerianAcceleration) {

  402.             final FieldVector3D<T> pvP = getPosition();
  403.             final T r2 = pvP.getNormSq();
  404.             final T r = r2.sqrt();
  405.             // extract non-Keplerian part of the initial acceleration
  406.             final FieldVector3D<T> nonKeplerianAcceleration = new FieldVector3D<>(getOne(), getPVCoordinates().getAcceleration(),
  407.                                                                                   r.multiply(r2).reciprocal().multiply(getMu()), pvP);

  408.             // add the quadratic motion due to the non-Keplerian acceleration to the Keplerian motion
  409.             final FieldVector3D<T> shiftedP = shiftedPV.getPosition();
  410.             final FieldVector3D<T> shiftedV = shiftedPV.getVelocity();
  411.             final FieldVector3D<T> fixedP   = new FieldVector3D<>(getOne(), shiftedP,
  412.                                                                   dt.square().multiply(0.5), nonKeplerianAcceleration);
  413.             final T                fixedR2 = fixedP.getNormSq();
  414.             final T                fixedR  = fixedR2.sqrt();
  415.             final FieldVector3D<T> fixedV  = new FieldVector3D<>(getOne(), shiftedV,
  416.                                                                  dt, nonKeplerianAcceleration);
  417.             final FieldVector3D<T> fixedA  = new FieldVector3D<>(fixedR.multiply(fixedR2).reciprocal().multiply(getMu().negate()), shiftedP,
  418.                                                                  getOne(), nonKeplerianAcceleration);

  419.             return new FieldPVCoordinates<>(fixedP, fixedV, fixedA);

  420.         } else {
  421.             // don't include acceleration,
  422.             // so the shifted orbit is not considered to have derivatives
  423.             return shiftedPV;
  424.         }
  425.     }

  426.     /** Create a 6x6 identity matrix.
  427.      * @return 6x6 identity matrix
  428.      */
  429.     private T[][] create6x6Identity() {
  430.         // this is the fastest way to set the 6x6 identity matrix
  431.         final T[][] identity = MathArrays.buildArray(getField(), 6, 6);
  432.         for (int i = 0; i < 6; i++) {
  433.             Arrays.fill(identity[i], getZero());
  434.             identity[i][i] = getOne();
  435.         }
  436.         return identity;
  437.     }

  438.     @Override
  439.     protected T[][] computeJacobianMeanWrtCartesian() {
  440.         return create6x6Identity();
  441.     }

  442.     @Override
  443.     protected T[][] computeJacobianEccentricWrtCartesian() {
  444.         return create6x6Identity();
  445.     }

  446.     @Override
  447.     protected T[][] computeJacobianTrueWrtCartesian() {
  448.         return create6x6Identity();
  449.     }

  450.     /** {@inheritDoc} */
  451.     public void addKeplerContribution(final PositionAngleType type, final T gm,
  452.                                       final T[] pDot) {

  453.         final FieldPVCoordinates<T> pv = getPVCoordinates();

  454.         // position derivative is velocity
  455.         final FieldVector3D<T> velocity = pv.getVelocity();
  456.         pDot[0] = pDot[0].add(velocity.getX());
  457.         pDot[1] = pDot[1].add(velocity.getY());
  458.         pDot[2] = pDot[2].add(velocity.getZ());

  459.         // velocity derivative is Newtonian acceleration
  460.         final FieldVector3D<T> position = pv.getPosition();
  461.         final T r2         = position.getNormSq();
  462.         final T coeff      = r2.multiply(r2.sqrt()).reciprocal().negate().multiply(gm);
  463.         pDot[3] = pDot[3].add(coeff.multiply(position.getX()));
  464.         pDot[4] = pDot[4].add(coeff.multiply(position.getY()));
  465.         pDot[5] = pDot[5].add(coeff.multiply(position.getZ()));

  466.     }

  467.     /**  Returns a string representation of this Orbit object.
  468.      * @return a string representation of this object
  469.      */
  470.     public String toString() {
  471.         // use only the six defining elements, like the other Orbit.toString() methods
  472.         final String comma = ", ";
  473.         final PVCoordinates pv = getPVCoordinates().toPVCoordinates();
  474.         final Vector3D p = pv.getPosition();
  475.         final Vector3D v = pv.getVelocity();
  476.         return "Cartesian parameters: {P(" +
  477.                 p.getX() + comma +
  478.                 p.getY() + comma +
  479.                 p.getZ() + "), V(" +
  480.                 v.getX() + comma +
  481.                 v.getY() + comma +
  482.                 v.getZ() + ")}";
  483.     }

  484.     @Override
  485.     public CartesianOrbit toOrbit() {
  486.         final PVCoordinates pv = getPVCoordinates().toPVCoordinates();
  487.         final AbsoluteDate date = getPVCoordinates().getDate().toAbsoluteDate();
  488.         if (hasDerivatives()) {
  489.             // getPVCoordinates includes accelerations that will be interpreted as derivatives
  490.             return new CartesianOrbit(pv, getFrame(), date, getMu().getReal());
  491.         } else {
  492.             // get rid of Keplerian acceleration so we don't assume
  493.             // we have derivatives when in fact we don't have them
  494.             return new CartesianOrbit(new PVCoordinates(pv.getPosition(), pv.getVelocity()),
  495.                                       getFrame(), date, getMu().getReal());
  496.         }
  497.     }

  498. }