Orbit.java

  1. /* Copyright 2002-2018 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 org.hipparchus.geometry.euclidean.threed.Vector3D;
  20. import org.hipparchus.linear.DecompositionSolver;
  21. import org.hipparchus.linear.MatrixUtils;
  22. import org.hipparchus.linear.QRDecomposition;
  23. import org.hipparchus.linear.RealMatrix;
  24. import org.hipparchus.util.FastMath;
  25. import org.hipparchus.util.MathArrays;
  26. import org.orekit.errors.OrekitIllegalArgumentException;
  27. import org.orekit.errors.OrekitException;
  28. import org.orekit.errors.OrekitInternalError;
  29. import org.orekit.errors.OrekitMessages;
  30. import org.orekit.frames.Frame;
  31. import org.orekit.frames.Transform;
  32. import org.orekit.time.AbsoluteDate;
  33. import org.orekit.time.TimeInterpolable;
  34. import org.orekit.time.TimeShiftable;
  35. import org.orekit.time.TimeStamped;
  36. import org.orekit.utils.PVCoordinates;
  37. import org.orekit.utils.PVCoordinatesProvider;
  38. import org.orekit.utils.TimeStampedPVCoordinates;

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

  41.  * <p>
  42.  * For user convenience, both the Cartesian and the equinoctial elements
  43.  * are provided by this class, regardless of the canonical representation
  44.  * implemented in the derived class (which may be classical Keplerian
  45.  * elements for example).
  46.  * </p>
  47.  * <p>
  48.  * The parameters are defined in a frame specified by the user. It is important
  49.  * to make sure this frame is consistent: it probably is inertial and centered
  50.  * on the central body. This information is used for example by some
  51.  * force models.
  52.  * </p>
  53.  * <p>
  54.  * Instance of this class are guaranteed to be immutable.
  55.  * </p>
  56.  * @author Luc Maisonobe
  57.  * @author Guylaine Prat
  58.  * @author Fabien Maussion
  59.  * @author V&eacute;ronique Pommier-Maurussane
  60.  */
  61. public abstract class Orbit
  62.     implements TimeStamped, TimeShiftable<Orbit>, TimeInterpolable<Orbit>,
  63.                Serializable, PVCoordinatesProvider {

  64.     /** Serializable UID. */
  65.     private static final long serialVersionUID = 438733454597999578L;

  66.     /** Frame in which are defined the orbital parameters. */
  67.     private final Frame frame;

  68.     /** Date of the orbital parameters. */
  69.     private final AbsoluteDate date;

  70.     /** Value of mu used to compute position and velocity (m³/s²). */
  71.     private final double mu;

  72.     /** Computed PVCoordinates. */
  73.     private transient TimeStampedPVCoordinates pvCoordinates;

  74.     /** Jacobian of the orbital parameters with mean angle with respect to the Cartesian coordinates. */
  75.     private transient double[][] jacobianMeanWrtCartesian;

  76.     /** Jacobian of the Cartesian coordinates with respect to the orbital parameters with mean angle. */
  77.     private transient double[][] jacobianWrtParametersMean;

  78.     /** Jacobian of the orbital parameters with eccentric angle with respect to the Cartesian coordinates. */
  79.     private transient double[][] jacobianEccentricWrtCartesian;

  80.     /** Jacobian of the Cartesian coordinates with respect to the orbital parameters with eccentric angle. */
  81.     private transient double[][] jacobianWrtParametersEccentric;

  82.     /** Jacobian of the orbital parameters with true angle with respect to the Cartesian coordinates. */
  83.     private transient double[][] jacobianTrueWrtCartesian;

  84.     /** Jacobian of the Cartesian coordinates with respect to the orbital parameters with true angle. */
  85.     private transient double[][] jacobianWrtParametersTrue;

  86.     /** Default constructor.
  87.      * Build a new instance with arbitrary default elements.
  88.      * @param frame the frame in which the parameters are defined
  89.      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
  90.      * @param date date of the orbital parameters
  91.      * @param mu central attraction coefficient (m^3/s^2)
  92.      * @exception IllegalArgumentException if frame is not a {@link
  93.      * Frame#isPseudoInertial pseudo-inertial frame}
  94.      */
  95.     protected Orbit(final Frame frame, final AbsoluteDate date, final double mu)
  96.         throws IllegalArgumentException {
  97.         ensurePseudoInertialFrame(frame);
  98.         this.date                      = date;
  99.         this.mu                        = mu;
  100.         this.pvCoordinates             = null;
  101.         this.frame                     = frame;
  102.         jacobianMeanWrtCartesian       = null;
  103.         jacobianWrtParametersMean      = null;
  104.         jacobianEccentricWrtCartesian  = null;
  105.         jacobianWrtParametersEccentric = null;
  106.         jacobianTrueWrtCartesian       = null;
  107.         jacobianWrtParametersTrue      = null;
  108.     }

  109.     /** Set the orbit from Cartesian parameters.
  110.      *
  111.      * <p> The acceleration provided in {@code pvCoordinates} is accessible using
  112.      * {@link #getPVCoordinates()} and {@link #getPVCoordinates(Frame)}. All other methods
  113.      * use {@code mu} and the position to compute the acceleration, including
  114.      * {@link #shiftedBy(double)} and {@link #getPVCoordinates(AbsoluteDate, Frame)}.
  115.      *
  116.      * @param pvCoordinates the position and velocity in the inertial frame
  117.      * @param frame the frame in which the {@link TimeStampedPVCoordinates} are defined
  118.      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
  119.      * @param mu central attraction coefficient (m^3/s^2)
  120.      * @exception IllegalArgumentException if frame is not a {@link
  121.      * Frame#isPseudoInertial pseudo-inertial frame}
  122.      */
  123.     protected Orbit(final TimeStampedPVCoordinates pvCoordinates, final Frame frame, final double mu)
  124.         throws IllegalArgumentException {
  125.         ensurePseudoInertialFrame(frame);
  126.         this.date = pvCoordinates.getDate();
  127.         this.mu = mu;
  128.         if (pvCoordinates.getAcceleration().getNormSq() == 0) {
  129.             // the acceleration was not provided,
  130.             // compute it from Newtonian attraction
  131.             final double r2 = pvCoordinates.getPosition().getNormSq();
  132.             final double r3 = r2 * FastMath.sqrt(r2);
  133.             this.pvCoordinates = new TimeStampedPVCoordinates(pvCoordinates.getDate(),
  134.                                                               pvCoordinates.getPosition(),
  135.                                                               pvCoordinates.getVelocity(),
  136.                                                               new Vector3D(-mu / r3, pvCoordinates.getPosition()));
  137.         } else {
  138.             this.pvCoordinates = pvCoordinates;
  139.         }
  140.         this.frame = frame;
  141.     }

  142.     /** Check if Cartesian coordinates include non-Keplerian acceleration.
  143.      * @param pva Cartesian coordinates
  144.      * @param mu central attraction coefficient
  145.      * @return true if Cartesian coordinates include non-Keplerian acceleration
  146.      */
  147.     protected static boolean hasNonKeplerianAcceleration(final PVCoordinates pva, final double mu) {

  148.         final Vector3D a = pva.getAcceleration();
  149.         if (a == null) {
  150.             return false;
  151.         }

  152.         final Vector3D p = pva.getPosition();
  153.         final double r2 = p.getNormSq();
  154.         final double r  = FastMath.sqrt(r2);
  155.         final Vector3D keplerianAcceleration = new Vector3D(-mu / (r * r2), p);
  156.         if (a.getNorm() > 1.0e-9 * keplerianAcceleration.getNorm()) {
  157.             // we have a relevant acceleration, we can compute derivatives
  158.             return true;
  159.         } else {
  160.             // the provided acceleration is either too small to be reliable (probably even 0), or NaN
  161.             return false;
  162.         }

  163.     }

  164.     /** Get the orbit type.
  165.      * @return orbit type
  166.      */
  167.     public abstract OrbitType getType();

  168.     /** Ensure the defining frame is a pseudo-inertial frame.
  169.      * @param frame frame to check
  170.      * @exception IllegalArgumentException if frame is not a {@link
  171.      * Frame#isPseudoInertial pseudo-inertial frame}
  172.      */
  173.     private static void ensurePseudoInertialFrame(final Frame frame)
  174.         throws IllegalArgumentException {
  175.         if (!frame.isPseudoInertial()) {
  176.             throw new OrekitIllegalArgumentException(OrekitMessages.NON_PSEUDO_INERTIAL_FRAME,
  177.                                                      frame.getName());
  178.         }
  179.     }

  180.     /** Get the frame in which the orbital parameters are defined.
  181.      * @return frame in which the orbital parameters are defined
  182.      */
  183.     public Frame getFrame() {
  184.         return frame;
  185.     }

  186.     /** Get the semi-major axis.
  187.      * <p>Note that the semi-major axis is considered negative for hyperbolic orbits.</p>
  188.      * @return semi-major axis (m)
  189.      */
  190.     public abstract double getA();

  191.     /** Get the semi-major axis derivative.
  192.      * <p>Note that the semi-major axis is considered negative for hyperbolic orbits.</p>
  193.      * <p>
  194.      * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
  195.      * </p>
  196.      * @return semi-major axis  derivative (m/s)
  197.      * @see #hasDerivatives()
  198.      * @since 9.0
  199.      */
  200.     public abstract double getADot();

  201.     /** Get the first component of the equinoctial eccentricity vector derivative.
  202.      * @return first component of the equinoctial eccentricity vector derivative
  203.      */
  204.     public abstract double getEquinoctialEx();

  205.     /** Get the first component of the equinoctial eccentricity vector.
  206.      * <p>
  207.      * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
  208.      * </p>
  209.      * @return first component of the equinoctial eccentricity vector
  210.      * @see #hasDerivatives()
  211.      * @since 9.0
  212.      */
  213.     public abstract double getEquinoctialExDot();

  214.     /** Get the second component of the equinoctial eccentricity vector derivative.
  215.      * @return second component of the equinoctial eccentricity vector derivative
  216.      */
  217.     public abstract double getEquinoctialEy();

  218.     /** Get the second component of the equinoctial eccentricity vector.
  219.      * <p>
  220.      * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
  221.      * </p>
  222.      * @return second component of the equinoctial eccentricity vector
  223.      * @see #hasDerivatives()
  224.      * @since 9.0
  225.      */
  226.     public abstract double getEquinoctialEyDot();

  227.     /** Get the first component of the inclination vector.
  228.      * @return first component of the inclination vector
  229.      */
  230.     public abstract double getHx();

  231.     /** Get the first component of the inclination vector derivative.
  232.      * <p>
  233.      * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
  234.      * </p>
  235.      * @return first component of the inclination vector derivative
  236.      * @see #hasDerivatives()
  237.      * @since 9.0
  238.      */
  239.     public abstract double getHxDot();

  240.     /** Get the second component of the inclination vector.
  241.      * @return second component of the inclination vector
  242.      */
  243.     public abstract double getHy();

  244.     /** Get the second component of the inclination vector derivative.
  245.      * <p>
  246.      * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
  247.      * </p>
  248.      * @return second component of the inclination vector derivative
  249.      * @see #hasDerivatives()
  250.      * @since 9.0
  251.      */
  252.     public abstract double getHyDot();

  253.     /** Get the eccentric longitude argument.
  254.      * @return E + ω + Ω eccentric longitude argument (rad)
  255.      */
  256.     public abstract double getLE();

  257.     /** Get the eccentric longitude argument derivative.
  258.      * <p>
  259.      * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
  260.      * </p>
  261.      * @return d(E + ω + Ω)/dt eccentric longitude argument derivative (rad/s)
  262.      * @see #hasDerivatives()
  263.      * @since 9.0
  264.      */
  265.     public abstract double getLEDot();

  266.     /** Get the true longitude argument.
  267.      * @return v + ω + Ω true longitude argument (rad)
  268.      */
  269.     public abstract double getLv();

  270.     /** Get the true longitude argument derivative.
  271.      * <p>
  272.      * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
  273.      * </p>
  274.      * @return d(v + ω + Ω)/dt true longitude argument derivative (rad/s)
  275.      * @see #hasDerivatives()
  276.      * @since 9.0
  277.      */
  278.     public abstract double getLvDot();

  279.     /** Get the mean longitude argument.
  280.      * @return M + ω + Ω mean longitude argument (rad)
  281.      */
  282.     public abstract double getLM();

  283.     /** Get the mean longitude argument derivative.
  284.      * <p>
  285.      * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
  286.      * </p>
  287.      * @return d(M + ω + Ω)/dt mean longitude argument derivative (rad/s)
  288.      * @see #hasDerivatives()
  289.      * @since 9.0
  290.      */
  291.     public abstract double getLMDot();

  292.     // Additional orbital elements

  293.     /** Get the eccentricity.
  294.      * @return eccentricity
  295.      */
  296.     public abstract double getE();

  297.     /** Get the eccentricity derivative.
  298.      * <p>
  299.      * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
  300.      * </p>
  301.      * @return eccentricity derivative
  302.      * @see #hasDerivatives()
  303.      * @since 9.0
  304.      */
  305.     public abstract double getEDot();

  306.     /** Get the inclination.
  307.      * @return inclination (rad)
  308.      */
  309.     public abstract double getI();

  310.     /** Get the inclination derivative.
  311.      * <p>
  312.      * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
  313.      * </p>
  314.      * @return inclination derivative (rad/s)
  315.      * @see #hasDerivatives()
  316.      * @since 9.0
  317.      */
  318.     public abstract double getIDot();

  319.     /** Check if orbit includes derivatives.
  320.      * @return true if orbit includes derivatives
  321.      * @see #getADot()
  322.      * @see #getEquinoctialExDot()
  323.      * @see #getEquinoctialEyDot()
  324.      * @see #getHxDot()
  325.      * @see #getHyDot()
  326.      * @see #getLEDot()
  327.      * @see #getLvDot()
  328.      * @see #getLMDot()
  329.      * @see #getEDot()
  330.      * @see #getIDot()
  331.      * @since 9.0
  332.      */
  333.     public boolean hasDerivatives() {
  334.         return !Double.isNaN(getADot());
  335.     }

  336.     /** Get the central acceleration constant.
  337.      * @return central acceleration constant
  338.      */
  339.     public double getMu() {
  340.         return mu;
  341.     }

  342.     /** Get the Keplerian period.
  343.      * <p>The Keplerian period is computed directly from semi major axis
  344.      * and central acceleration constant.</p>
  345.      * @return Keplerian period in seconds, or positive infinity for hyperbolic orbits
  346.      */
  347.     public double getKeplerianPeriod() {
  348.         final double a = getA();
  349.         return (a < 0) ? Double.POSITIVE_INFINITY : 2.0 * FastMath.PI * a * FastMath.sqrt(a / mu);
  350.     }

  351.     /** Get the Keplerian mean motion.
  352.      * <p>The Keplerian mean motion is computed directly from semi major axis
  353.      * and central acceleration constant.</p>
  354.      * @return Keplerian mean motion in radians per second
  355.      */
  356.     public double getKeplerianMeanMotion() {
  357.         final double absA = FastMath.abs(getA());
  358.         return FastMath.sqrt(mu / absA) / absA;
  359.     }

  360.     /** Get the date of orbital parameters.
  361.      * @return date of the orbital parameters
  362.      */
  363.     public AbsoluteDate getDate() {
  364.         return date;
  365.     }

  366.     /** Get the {@link TimeStampedPVCoordinates} in a specified frame.
  367.      * @param outputFrame frame in which the position/velocity coordinates shall be computed
  368.      * @return pvCoordinates in the specified output frame
  369.      * @exception OrekitException if transformation between frames cannot be computed
  370.      * @see #getPVCoordinates()
  371.      */
  372.     public TimeStampedPVCoordinates getPVCoordinates(final Frame outputFrame)
  373.         throws OrekitException {
  374.         if (pvCoordinates == null) {
  375.             pvCoordinates = initPVCoordinates();
  376.         }

  377.         // If output frame requested is the same as definition frame,
  378.         // PV coordinates are returned directly
  379.         if (outputFrame == frame) {
  380.             return pvCoordinates;
  381.         }

  382.         // Else, PV coordinates are transformed to output frame
  383.         final Transform t = frame.getTransformTo(outputFrame, date);
  384.         return t.transformPVCoordinates(pvCoordinates);
  385.     }

  386.     /** {@inheritDoc} */
  387.     public TimeStampedPVCoordinates getPVCoordinates(final AbsoluteDate otherDate, final Frame otherFrame)
  388.         throws OrekitException {
  389.         return shiftedBy(otherDate.durationFrom(getDate())).getPVCoordinates(otherFrame);
  390.     }


  391.     /** Get the {@link TimeStampedPVCoordinates} in definition frame.
  392.      * @return pvCoordinates in the definition frame
  393.      * @see #getPVCoordinates(Frame)
  394.      */
  395.     public TimeStampedPVCoordinates getPVCoordinates() {
  396.         if (pvCoordinates == null) {
  397.             pvCoordinates = initPVCoordinates();
  398.         }
  399.         return pvCoordinates;
  400.     }

  401.     /** Compute the position/velocity coordinates from the canonical parameters.
  402.      * @return computed position/velocity coordinates
  403.      */
  404.     protected abstract TimeStampedPVCoordinates initPVCoordinates();

  405.     /** Get a time-shifted orbit.
  406.      * <p>
  407.      * The orbit can be slightly shifted to close dates. The shifting model is a
  408.      * Keplerian one if no derivatives are available in the orbit, or Keplerian
  409.      * plus quadratic effect of the non-Keplerian acceleration if derivatives are
  410.      * available. Shifting is <em>not</em> intended as a replacement for proper
  411.      * orbit propagation but should be sufficient for small time shifts or coarse
  412.      * accuracy.
  413.      * </p>
  414.      * @param dt time shift in seconds
  415.      * @return a new orbit, shifted with respect to the instance (which is immutable)
  416.      */
  417.     public abstract Orbit shiftedBy(double dt);

  418.     /** Compute the Jacobian of the orbital parameters with respect to the Cartesian parameters.
  419.      * <p>
  420.      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
  421.      * respect to Cartesian coordinate j. This means each row corresponds to one orbital parameter
  422.      * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.
  423.      * </p>
  424.      * @param type type of the position angle to use
  425.      * @param jacobian placeholder 6x6 (or larger) matrix to be filled with the Jacobian, if matrix
  426.      * is larger than 6x6, only the 6x6 upper left corner will be modified
  427.      */
  428.     public void getJacobianWrtCartesian(final PositionAngle type, final double[][] jacobian) {

  429.         final double[][] cachedJacobian;
  430.         synchronized (this) {
  431.             switch (type) {
  432.                 case MEAN :
  433.                     if (jacobianMeanWrtCartesian == null) {
  434.                         // first call, we need to compute the Jacobian and cache it
  435.                         jacobianMeanWrtCartesian = computeJacobianMeanWrtCartesian();
  436.                     }
  437.                     cachedJacobian = jacobianMeanWrtCartesian;
  438.                     break;
  439.                 case ECCENTRIC :
  440.                     if (jacobianEccentricWrtCartesian == null) {
  441.                         // first call, we need to compute the Jacobian and cache it
  442.                         jacobianEccentricWrtCartesian = computeJacobianEccentricWrtCartesian();
  443.                     }
  444.                     cachedJacobian = jacobianEccentricWrtCartesian;
  445.                     break;
  446.                 case TRUE :
  447.                     if (jacobianTrueWrtCartesian == null) {
  448.                         // first call, we need to compute the Jacobian and cache it
  449.                         jacobianTrueWrtCartesian = computeJacobianTrueWrtCartesian();
  450.                     }
  451.                     cachedJacobian = jacobianTrueWrtCartesian;
  452.                     break;
  453.                 default :
  454.                     throw new OrekitInternalError(null);
  455.             }
  456.         }

  457.         // fill the user provided array
  458.         for (int i = 0; i < cachedJacobian.length; ++i) {
  459.             System.arraycopy(cachedJacobian[i], 0, jacobian[i], 0, cachedJacobian[i].length);
  460.         }

  461.     }

  462.     /** Compute the Jacobian of the Cartesian parameters with respect to the orbital parameters.
  463.      * <p>
  464.      * Element {@code jacobian[i][j]} is the derivative of Cartesian coordinate i of the orbit with
  465.      * respect to orbital parameter j. This means each row corresponds to one Cartesian coordinate
  466.      * x, y, z, xdot, ydot, zdot whereas columns 0 to 5 correspond to the orbital parameters.
  467.      * </p>
  468.      * @param type type of the position angle to use
  469.      * @param jacobian placeholder 6x6 (or larger) matrix to be filled with the Jacobian, if matrix
  470.      * is larger than 6x6, only the 6x6 upper left corner will be modified
  471.      */
  472.     public void getJacobianWrtParameters(final PositionAngle type, final double[][] jacobian) {

  473.         final double[][] cachedJacobian;
  474.         synchronized (this) {
  475.             switch (type) {
  476.                 case MEAN :
  477.                     if (jacobianWrtParametersMean == null) {
  478.                         // first call, we need to compute the Jacobian and cache it
  479.                         jacobianWrtParametersMean = createInverseJacobian(type);
  480.                     }
  481.                     cachedJacobian = jacobianWrtParametersMean;
  482.                     break;
  483.                 case ECCENTRIC :
  484.                     if (jacobianWrtParametersEccentric == null) {
  485.                         // first call, we need to compute the Jacobian and cache it
  486.                         jacobianWrtParametersEccentric = createInverseJacobian(type);
  487.                     }
  488.                     cachedJacobian = jacobianWrtParametersEccentric;
  489.                     break;
  490.                 case TRUE :
  491.                     if (jacobianWrtParametersTrue == null) {
  492.                         // first call, we need to compute the Jacobian and cache it
  493.                         jacobianWrtParametersTrue = createInverseJacobian(type);
  494.                     }
  495.                     cachedJacobian = jacobianWrtParametersTrue;
  496.                     break;
  497.                 default :
  498.                     throw new OrekitInternalError(null);
  499.             }
  500.         }

  501.         // fill the user-provided array
  502.         for (int i = 0; i < cachedJacobian.length; ++i) {
  503.             System.arraycopy(cachedJacobian[i], 0, jacobian[i], 0, cachedJacobian[i].length);
  504.         }

  505.     }

  506.     /** Create an inverse Jacobian.
  507.      * @param type type of the position angle to use
  508.      * @return inverse Jacobian
  509.      */
  510.     private double[][] createInverseJacobian(final PositionAngle type) {

  511.         // get the direct Jacobian
  512.         final double[][] directJacobian = new double[6][6];
  513.         getJacobianWrtCartesian(type, directJacobian);

  514.         // invert the direct Jacobian
  515.         final RealMatrix matrix = MatrixUtils.createRealMatrix(directJacobian);
  516.         final DecompositionSolver solver = new QRDecomposition(matrix).getSolver();
  517.         return solver.getInverse().getData();

  518.     }

  519.     /** Compute the Jacobian of the orbital parameters with mean angle with respect to the Cartesian parameters.
  520.      * <p>
  521.      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
  522.      * respect to Cartesian coordinate j. This means each row correspond to one orbital parameter
  523.      * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.
  524.      * </p>
  525.      * @return 6x6 Jacobian matrix
  526.      * @see #computeJacobianEccentricWrtCartesian()
  527.      * @see #computeJacobianTrueWrtCartesian()
  528.      */
  529.     protected abstract double[][] computeJacobianMeanWrtCartesian();

  530.     /** Compute the Jacobian of the orbital parameters with eccentric angle with respect to the Cartesian parameters.
  531.      * <p>
  532.      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
  533.      * respect to Cartesian coordinate j. This means each row correspond to one orbital parameter
  534.      * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.
  535.      * </p>
  536.      * @return 6x6 Jacobian matrix
  537.      * @see #computeJacobianMeanWrtCartesian()
  538.      * @see #computeJacobianTrueWrtCartesian()
  539.      */
  540.     protected abstract double[][] computeJacobianEccentricWrtCartesian();

  541.     /** Compute the Jacobian of the orbital parameters with true angle with respect to the Cartesian parameters.
  542.      * <p>
  543.      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
  544.      * respect to Cartesian coordinate j. This means each row correspond to one orbital parameter
  545.      * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.
  546.      * </p>
  547.      * @return 6x6 Jacobian matrix
  548.      * @see #computeJacobianMeanWrtCartesian()
  549.      * @see #computeJacobianEccentricWrtCartesian()
  550.      */
  551.     protected abstract double[][] computeJacobianTrueWrtCartesian();

  552.     /** Add the contribution of the Keplerian motion to parameters derivatives
  553.      * <p>
  554.      * This method is used by integration-based propagators to evaluate the part of Keplerian
  555.      * motion to evolution of the orbital state.
  556.      * </p>
  557.      * @param type type of the position angle in the state
  558.      * @param gm attraction coefficient to use
  559.      * @param pDot array containing orbital state derivatives to update (the Keplerian
  560.      * part must be <em>added</em> to the array components, as the array may already
  561.      * contain some non-zero elements corresponding to non-Keplerian parts)
  562.      */
  563.     public abstract void addKeplerContribution(PositionAngle type, double gm, double[] pDot);

  564.         /** Fill a Jacobian half row with a single vector.
  565.      * @param a coefficient of the vector
  566.      * @param v vector
  567.      * @param row Jacobian matrix row
  568.      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
  569.      */
  570.     protected static void fillHalfRow(final double a, final Vector3D v, final double[] row, final int j) {
  571.         row[j]     = a * v.getX();
  572.         row[j + 1] = a * v.getY();
  573.         row[j + 2] = a * v.getZ();
  574.     }

  575.     /** Fill a Jacobian half row with a linear combination of vectors.
  576.      * @param a1 coefficient of the first vector
  577.      * @param v1 first vector
  578.      * @param a2 coefficient of the second vector
  579.      * @param v2 second vector
  580.      * @param row Jacobian matrix row
  581.      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
  582.      */
  583.     protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2,
  584.                                       final double[] row, final int j) {
  585.         row[j]     = MathArrays.linearCombination(a1, v1.getX(), a2, v2.getX());
  586.         row[j + 1] = MathArrays.linearCombination(a1, v1.getY(), a2, v2.getY());
  587.         row[j + 2] = MathArrays.linearCombination(a1, v1.getZ(), a2, v2.getZ());
  588.     }

  589.     /** Fill a Jacobian half row with a linear combination of vectors.
  590.      * @param a1 coefficient of the first vector
  591.      * @param v1 first vector
  592.      * @param a2 coefficient of the second vector
  593.      * @param v2 second vector
  594.      * @param a3 coefficient of the third vector
  595.      * @param v3 third vector
  596.      * @param row Jacobian matrix row
  597.      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
  598.      */
  599.     protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2,
  600.                                       final double a3, final Vector3D v3,
  601.                                       final double[] row, final int j) {
  602.         row[j]     = MathArrays.linearCombination(a1, v1.getX(), a2, v2.getX(), a3, v3.getX());
  603.         row[j + 1] = MathArrays.linearCombination(a1, v1.getY(), a2, v2.getY(), a3, v3.getY());
  604.         row[j + 2] = MathArrays.linearCombination(a1, v1.getZ(), a2, v2.getZ(), a3, v3.getZ());
  605.     }

  606.     /** Fill a Jacobian half row with a linear combination of vectors.
  607.      * @param a1 coefficient of the first vector
  608.      * @param v1 first vector
  609.      * @param a2 coefficient of the second vector
  610.      * @param v2 second vector
  611.      * @param a3 coefficient of the third vector
  612.      * @param v3 third vector
  613.      * @param a4 coefficient of the fourth vector
  614.      * @param v4 fourth vector
  615.      * @param row Jacobian matrix row
  616.      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
  617.      */
  618.     protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2,
  619.                                       final double a3, final Vector3D v3, final double a4, final Vector3D v4,
  620.                                       final double[] row, final int j) {
  621.         row[j]     = MathArrays.linearCombination(a1, v1.getX(), a2, v2.getX(), a3, v3.getX(), a4, v4.getX());
  622.         row[j + 1] = MathArrays.linearCombination(a1, v1.getY(), a2, v2.getY(), a3, v3.getY(), a4, v4.getY());
  623.         row[j + 2] = MathArrays.linearCombination(a1, v1.getZ(), a2, v2.getZ(), a3, v3.getZ(), a4, v4.getZ());
  624.     }

  625.     /** Fill a Jacobian half row with a linear combination of vectors.
  626.      * @param a1 coefficient of the first vector
  627.      * @param v1 first vector
  628.      * @param a2 coefficient of the second vector
  629.      * @param v2 second vector
  630.      * @param a3 coefficient of the third vector
  631.      * @param v3 third vector
  632.      * @param a4 coefficient of the fourth vector
  633.      * @param v4 fourth vector
  634.      * @param a5 coefficient of the fifth vector
  635.      * @param v5 fifth vector
  636.      * @param row Jacobian matrix row
  637.      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
  638.      */
  639.     protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2,
  640.                                       final double a3, final Vector3D v3, final double a4, final Vector3D v4,
  641.                                       final double a5, final Vector3D v5,
  642.                                       final double[] row, final int j) {
  643.         final double[] a = new double[] {
  644.             a1, a2, a3, a4, a5
  645.         };
  646.         row[j]     = MathArrays.linearCombination(a, new double[] {
  647.             v1.getX(), v2.getX(), v3.getX(), v4.getX(), v5.getX()
  648.         });
  649.         row[j + 1] = MathArrays.linearCombination(a, new double[] {
  650.             v1.getY(), v2.getY(), v3.getY(), v4.getY(), v5.getY()
  651.         });
  652.         row[j + 2] = MathArrays.linearCombination(a, new double[] {
  653.             v1.getZ(), v2.getZ(), v3.getZ(), v4.getZ(), v5.getZ()
  654.         });
  655.     }

  656.     /** Fill a Jacobian half row with a linear combination of vectors.
  657.      * @param a1 coefficient of the first vector
  658.      * @param v1 first vector
  659.      * @param a2 coefficient of the second vector
  660.      * @param v2 second vector
  661.      * @param a3 coefficient of the third vector
  662.      * @param v3 third vector
  663.      * @param a4 coefficient of the fourth vector
  664.      * @param v4 fourth vector
  665.      * @param a5 coefficient of the fifth vector
  666.      * @param v5 fifth vector
  667.      * @param a6 coefficient of the sixth vector
  668.      * @param v6 sixth vector
  669.      * @param row Jacobian matrix row
  670.      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
  671.      */
  672.     protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2,
  673.                                       final double a3, final Vector3D v3, final double a4, final Vector3D v4,
  674.                                       final double a5, final Vector3D v5, final double a6, final Vector3D v6,
  675.                                       final double[] row, final int j) {
  676.         final double[] a = new double[] {
  677.             a1, a2, a3, a4, a5, a6
  678.         };
  679.         row[j]     = MathArrays.linearCombination(a, new double[] {
  680.             v1.getX(), v2.getX(), v3.getX(), v4.getX(), v5.getX(), v6.getX()
  681.         });
  682.         row[j + 1] = MathArrays.linearCombination(a, new double[] {
  683.             v1.getY(), v2.getY(), v3.getY(), v4.getY(), v5.getY(), v6.getY()
  684.         });
  685.         row[j + 2] = MathArrays.linearCombination(a, new double[] {
  686.             v1.getZ(), v2.getZ(), v3.getZ(), v4.getZ(), v5.getZ(), v6.getZ()
  687.         });
  688.     }

  689. }