Orbit.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 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.OrekitInternalError;
  28. import org.orekit.errors.OrekitMessages;
  29. import org.orekit.frames.Frame;
  30. import org.orekit.frames.Transform;
  31. import org.orekit.time.AbsoluteDate;
  32. import org.orekit.time.TimeInterpolable;
  33. import org.orekit.time.TimeShiftable;
  34. import org.orekit.time.TimeStamped;
  35. import org.orekit.utils.PVCoordinates;
  36. import org.orekit.utils.PVCoordinatesProvider;
  37. import org.orekit.utils.TimeStampedPVCoordinates;

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

  147.         final Vector3D p = pva.getPosition();
  148.         final double r2 = p.getNormSq();
  149.         final double r  = FastMath.sqrt(r2);
  150.         final Vector3D keplerianAcceleration = new Vector3D(-mu / (r * r2), p);

  151.         // Check if acceleration is null or relatively close to 0 compared to the keplerain acceleration
  152.         final Vector3D a = pva.getAcceleration();
  153.         if (a == null || a.getNorm() < 1e-9 * keplerianAcceleration.getNorm()) {
  154.             return false;
  155.         }

  156.         final Vector3D nonKeplerianAcceleration = a.subtract(keplerianAcceleration);

  157.         if ( nonKeplerianAcceleration.getNorm() > 1e-9 * keplerianAcceleration.getNorm()) {
  158.             // we have a relevant acceleration, we can compute derivatives
  159.             return true;
  160.         } else {
  161.             // the provided acceleration is either too small to be reliable (probably even 0), or NaN
  162.             return false;
  163.         }
  164.     }

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

  293.     // Additional orbital elements

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

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

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

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

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

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

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

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

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

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

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

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

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


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

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

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

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

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

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

  459.     }

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

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

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

  503.     }

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

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

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

  516.     }

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

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

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

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

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

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

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

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

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

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

  687. }