FieldOrbit.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 org.hipparchus.CalculusFieldElement;
  19. import org.hipparchus.geometry.euclidean.threed.FieldVector3D;
  20. import org.hipparchus.linear.FieldDecompositionSolver;
  21. import org.hipparchus.linear.FieldLUDecomposition;
  22. import org.hipparchus.linear.FieldMatrix;
  23. import org.hipparchus.linear.MatrixUtils;
  24. import org.hipparchus.util.MathArrays;
  25. import org.orekit.errors.OrekitIllegalArgumentException;
  26. import org.orekit.errors.OrekitInternalError;
  27. import org.orekit.errors.OrekitMessages;
  28. import org.orekit.frames.Frame;
  29. import org.orekit.frames.Transform;
  30. import org.orekit.time.FieldAbsoluteDate;
  31. import org.orekit.time.FieldTimeInterpolable;
  32. import org.orekit.time.FieldTimeShiftable;
  33. import org.orekit.time.FieldTimeStamped;
  34. import org.orekit.utils.FieldPVCoordinates;
  35. import org.orekit.utils.FieldPVCoordinatesProvider;
  36. import org.orekit.utils.TimeStampedFieldPVCoordinates;
  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.  * @since 9.0
  60.  */
  61. public abstract class FieldOrbit<T extends CalculusFieldElement<T>>
  62.     implements FieldPVCoordinatesProvider<T>, FieldTimeStamped<T>, FieldTimeShiftable<FieldOrbit<T>, T>, FieldTimeInterpolable<FieldOrbit<T>, T> {

  63.     /** Frame in which are defined the orbital parameters. */
  64.     private final Frame frame;

  65.     /** Date of the orbital parameters. */
  66.     private final FieldAbsoluteDate<T> date;

  67.     /** Value of mu used to compute position and velocity (m³/s²). */
  68.     private final T mu;

  69.     /** Computed PVCoordinates. */
  70.     private transient TimeStampedFieldPVCoordinates<T> pvCoordinates;

  71.     /** Jacobian of the orbital parameters with mean angle with respect to the Cartesian coordinates. */
  72.     private transient T[][] jacobianMeanWrtCartesian;

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

  75.     /** Jacobian of the orbital parameters with eccentric angle with respect to the Cartesian coordinates. */
  76.     private transient T[][] jacobianEccentricWrtCartesian;

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

  79.     /** Jacobian of the orbital parameters with true angle with respect to the Cartesian coordinates. */
  80.     private transient T[][] jacobianTrueWrtCartesian;

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

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

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

  139.     /** Check if Cartesian coordinates include non-Keplerian acceleration.
  140.      * @param pva Cartesian coordinates
  141.      * @param mu central attraction coefficient
  142.      * @param <T> type of the field elements
  143.      * @return true if Cartesian coordinates include non-Keplerian acceleration
  144.      */
  145.     protected static <T extends CalculusFieldElement<T>> boolean hasNonKeplerianAcceleration(final FieldPVCoordinates<T> pva, final T mu) {

  146.         final FieldVector3D<T> a = pva.getAcceleration();
  147.         if (a == null) {
  148.             return false;
  149.         }

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

  161.     }

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

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

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

  184.     /**Transforms the FieldOrbit instance into an Orbit instance.
  185.      * @return Orbit instance with same properties
  186.      */
  187.     public abstract Orbit toOrbit();

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

  193.     /** Get the semi-major axis derivative.
  194.      * <p>Note that the semi-major axis is considered negative for hyperbolic orbits.</p>
  195.      * <p>
  196.      * If the orbit was created without derivatives, the value returned is null.
  197.      * </p>
  198.      * @return semi-major axis  derivative (m/s)
  199.      */
  200.     public abstract T getADot();

  201.     /** Get the first component of the equinoctial eccentricity vector.
  202.      * @return first component of the equinoctial eccentricity vector
  203.      */
  204.     public abstract T 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 null.
  208.      * </p>
  209.      * @return first component of the equinoctial eccentricity vector
  210.      */
  211.     public abstract T getEquinoctialExDot();

  212.     /** Get the second component of the equinoctial eccentricity vector.
  213.      * @return second component of the equinoctial eccentricity vector
  214.      */
  215.     public abstract T getEquinoctialEy();

  216.     /** Get the second component of the equinoctial eccentricity vector.
  217.      * <p>
  218.      * If the orbit was created without derivatives, the value returned is null.
  219.      * </p>
  220.      * @return second component of the equinoctial eccentricity vector
  221.      */
  222.     public abstract T getEquinoctialEyDot();

  223.     /** Get the first component of the inclination vector.
  224.      * @return first component of the inclination vector
  225.      */
  226.     public abstract T getHx();

  227.     /** Get the first component of the inclination vector derivative.
  228.      * <p>
  229.      * If the orbit was created without derivatives, the value returned is null.
  230.      * </p>
  231.      * @return first component of the inclination vector derivative
  232.      */
  233.     public abstract T getHxDot();

  234.     /** Get the second component of the inclination vector.
  235.      * @return second component of the inclination vector
  236.      */
  237.     public abstract T getHy();

  238.     /** Get the second component of the inclination vector derivative.
  239.      * @return second component of the inclination vector derivative
  240.      */
  241.     public abstract T getHyDot();

  242.     /** Get the eccentric longitude argument.
  243.      * @return E + ω + Ω eccentric longitude argument (rad)
  244.      */
  245.     public abstract T getLE();

  246.     /** Get the eccentric longitude argument derivative.
  247.      * <p>
  248.      * If the orbit was created without derivatives, the value returned is null.
  249.      * </p>
  250.      * @return d(E + ω + Ω)/dt eccentric longitude argument derivative (rad/s)
  251.      */
  252.     public abstract T getLEDot();

  253.     /** Get the true longitude argument.
  254.      * @return v + ω + Ω true longitude argument (rad)
  255.      */
  256.     public abstract T getLv();

  257.     /** Get the true longitude argument derivative.
  258.      * <p>
  259.      * If the orbit was created without derivatives, the value returned is null.
  260.      * </p>
  261.      * @return d(v + ω + Ω)/dt true longitude argument derivative (rad/s)
  262.      */
  263.     public abstract T getLvDot();

  264.     /** Get the mean longitude argument.
  265.      * @return M + ω + Ω mean longitude argument (rad)
  266.      */
  267.     public abstract T getLM();

  268.     /** Get the mean longitude argument derivative.
  269.      * <p>
  270.      * If the orbit was created without derivatives, the value returned is null.
  271.      * </p>
  272.      * @return d(M + ω + Ω)/dt mean longitude argument derivative (rad/s)
  273.      */
  274.     public abstract T getLMDot();

  275.     // Additional orbital elements

  276.     /** Get the eccentricity.
  277.      * @return eccentricity
  278.      */
  279.     public abstract T getE();

  280.     /** Get the eccentricity derivative.
  281.      * <p>
  282.      * If the orbit was created without derivatives, the value returned is null.
  283.      * </p>
  284.      * @return eccentricity derivative
  285.      */
  286.     public abstract T getEDot();

  287.     /** Get the inclination.
  288.      * <p>
  289.      * If the orbit was created without derivatives, the value returned is null.
  290.      * </p>
  291.      * @return inclination (rad)
  292.      */
  293.     public abstract T getI();

  294.     /** Get the inclination derivative.
  295.      * @return inclination derivative (rad/s)
  296.      */
  297.     public abstract T getIDot();

  298.     /** Get the central acceleration constant.
  299.      * @return central acceleration constant
  300.      */

  301.     /** Check if orbit includes derivatives.
  302.      * @return true if orbit includes derivatives
  303.      * @see #getADot()
  304.      * @see #getEquinoctialExDot()
  305.      * @see #getEquinoctialEyDot()
  306.      * @see #getHxDot()
  307.      * @see #getHyDot()
  308.      * @see #getLEDot()
  309.      * @see #getLvDot()
  310.      * @see #getLMDot()
  311.      * @see #getEDot()
  312.      * @see #getIDot()
  313.      * @since 9.0
  314.      */
  315.     public abstract boolean hasDerivatives();

  316.     /** Get the central attraction coefficient used for position and velocity conversions (m³/s²).
  317.      * @return central attraction coefficient used for position and velocity conversions (m³/s²)
  318.      */
  319.     public T getMu() {
  320.         return mu;
  321.     }

  322.     /** Get the Keplerian period.
  323.      * <p>The Keplerian period is computed directly from semi major axis
  324.      * and central acceleration constant.</p>
  325.      * @return Keplerian period in seconds, or positive infinity for hyperbolic orbits
  326.      */
  327.     public T getKeplerianPeriod() {
  328.         final T a = getA();
  329.         return (a.getReal() < 0) ? getA().getField().getZero().add(Double.POSITIVE_INFINITY) : a.multiply(a.getPi().multiply(2.0)).multiply(a.divide(mu).sqrt());
  330.     }

  331.     /** Get the Keplerian mean motion.
  332.      * <p>The Keplerian mean motion is computed directly from semi major axis
  333.      * and central acceleration constant.</p>
  334.      * @return Keplerian mean motion in radians per second
  335.      */
  336.     public T getKeplerianMeanMotion() {
  337.         final T absA = getA().abs();
  338.         return absA.reciprocal().multiply(mu).sqrt().divide(absA);
  339.     }

  340.     /** Get the date of orbital parameters.
  341.      * @return date of the orbital parameters
  342.      */
  343.     public FieldAbsoluteDate<T> getDate() {
  344.         return date;
  345.     }

  346.     /** Get the {@link TimeStampedPVCoordinates} in a specified frame.
  347.      * @param outputFrame frame in which the position/velocity coordinates shall be computed
  348.      * @return FieldPVCoordinates in the specified output frame
  349.           * @see #getPVCoordinates()
  350.      */
  351.     public TimeStampedFieldPVCoordinates<T> getPVCoordinates(final Frame outputFrame) {
  352.         if (pvCoordinates == null) {
  353.             pvCoordinates = initPVCoordinates();
  354.         }

  355.         // If output frame requested is the same as definition frame,
  356.         // PV coordinates are returned directly
  357.         if (outputFrame == frame) {
  358.             return pvCoordinates;
  359.         }

  360.         // Else, PV coordinates are transformed to output frame
  361.         final Transform t = frame.getTransformTo(outputFrame, date.toAbsoluteDate()); //TODO CHECK THIS
  362.         return t.transformPVCoordinates(pvCoordinates);
  363.     }

  364.     /** {@inheritDoc} */
  365.     public TimeStampedFieldPVCoordinates<T> getPVCoordinates(final FieldAbsoluteDate<T> otherDate, final Frame otherFrame) {
  366.         return shiftedBy(otherDate.durationFrom(getDate())).getPVCoordinates(otherFrame);
  367.     }


  368.     /** Get the {@link TimeStampedPVCoordinates} in definition frame.
  369.      * @return FieldPVCoordinates in the definition frame
  370.      * @see #getPVCoordinates(Frame)
  371.      */
  372.     public TimeStampedFieldPVCoordinates<T> getPVCoordinates() {
  373.         if (pvCoordinates == null) {
  374.             pvCoordinates = initPVCoordinates();

  375.         }
  376.         return pvCoordinates;
  377.     }

  378.     /** Compute the position/velocity coordinates from the canonical parameters.
  379.      * @return computed position/velocity coordinates
  380.      */
  381.     protected abstract TimeStampedFieldPVCoordinates<T> initPVCoordinates();

  382.     /** Get a time-shifted orbit.
  383.      * <p>
  384.      * The orbit can be slightly shifted to close dates. This shift is based on
  385.      * a simple Keplerian model. It is <em>not</em> intended as a replacement
  386.      * for proper orbit and attitude propagation but should be sufficient for
  387.      * small time shifts or coarse accuracy.
  388.      * </p>
  389.      * @param dt time shift in seconds
  390.      * @return a new orbit, shifted with respect to the instance (which is immutable)
  391.      */
  392.     public abstract FieldOrbit<T> shiftedBy(T dt);

  393.     /** Compute the Jacobian of the orbital parameters with respect to the Cartesian parameters.
  394.      * <p>
  395.      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
  396.      * respect to Cartesian coordinate j. This means each row corresponds to one orbital parameter
  397.      * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.
  398.      * </p>
  399.      * @param type type of the position angle to use
  400.      * @param jacobian placeholder 6x6 (or larger) matrix to be filled with the Jacobian, if matrix
  401.      * is larger than 6x6, only the 6x6 upper left corner will be modified
  402.      */
  403.     public void getJacobianWrtCartesian(final PositionAngle type, final T[][] jacobian) {

  404.         final T[][] cachedJacobian;
  405.         synchronized (this) {
  406.             switch (type) {
  407.                 case MEAN :
  408.                     if (jacobianMeanWrtCartesian == null) {
  409.                         // first call, we need to compute the Jacobian and cache it
  410.                         jacobianMeanWrtCartesian = computeJacobianMeanWrtCartesian();
  411.                     }
  412.                     cachedJacobian = jacobianMeanWrtCartesian;
  413.                     break;
  414.                 case ECCENTRIC :
  415.                     if (jacobianEccentricWrtCartesian == null) {
  416.                         // first call, we need to compute the Jacobian and cache it
  417.                         jacobianEccentricWrtCartesian = computeJacobianEccentricWrtCartesian();
  418.                     }
  419.                     cachedJacobian = jacobianEccentricWrtCartesian;
  420.                     break;
  421.                 case TRUE :
  422.                     if (jacobianTrueWrtCartesian == null) {
  423.                         // first call, we need to compute the Jacobian and cache it
  424.                         jacobianTrueWrtCartesian = computeJacobianTrueWrtCartesian();
  425.                     }
  426.                     cachedJacobian = jacobianTrueWrtCartesian;
  427.                     break;
  428.                 default :
  429.                     throw new OrekitInternalError(null);
  430.             }
  431.         }

  432.         // fill the user provided array
  433.         for (int i = 0; i < cachedJacobian.length; ++i) {
  434.             System.arraycopy(cachedJacobian[i], 0, jacobian[i], 0, cachedJacobian[i].length);
  435.         }

  436.     }

  437.     /** Compute the Jacobian of the Cartesian parameters with respect to the orbital parameters.
  438.      * <p>
  439.      * Element {@code jacobian[i][j]} is the derivative of Cartesian coordinate i of the orbit with
  440.      * respect to orbital parameter j. This means each row corresponds to one Cartesian coordinate
  441.      * x, y, z, xdot, ydot, zdot whereas columns 0 to 5 correspond to the orbital parameters.
  442.      * </p>
  443.      * @param type type of the position angle to use
  444.      * @param jacobian placeholder 6x6 (or larger) matrix to be filled with the Jacobian, if matrix
  445.      * is larger than 6x6, only the 6x6 upper left corner will be modified
  446.      */
  447.     public void getJacobianWrtParameters(final PositionAngle type, final T[][] jacobian) {

  448.         final T[][] cachedJacobian;
  449.         synchronized (this) {
  450.             switch (type) {
  451.                 case MEAN :
  452.                     if (jacobianWrtParametersMean == null) {
  453.                         // first call, we need to compute the Jacobian and cache it
  454.                         jacobianWrtParametersMean = createInverseJacobian(type);
  455.                     }
  456.                     cachedJacobian = jacobianWrtParametersMean;
  457.                     break;
  458.                 case ECCENTRIC :
  459.                     if (jacobianWrtParametersEccentric == null) {
  460.                         // first call, we need to compute the Jacobian and cache it
  461.                         jacobianWrtParametersEccentric = createInverseJacobian(type);
  462.                     }
  463.                     cachedJacobian = jacobianWrtParametersEccentric;
  464.                     break;
  465.                 case TRUE :
  466.                     if (jacobianWrtParametersTrue == null) {
  467.                         // first call, we need to compute the Jacobian and cache it
  468.                         jacobianWrtParametersTrue = createInverseJacobian(type);
  469.                     }
  470.                     cachedJacobian = jacobianWrtParametersTrue;
  471.                     break;
  472.                 default :
  473.                     throw new OrekitInternalError(null);
  474.             }
  475.         }

  476.         // fill the user-provided array
  477.         for (int i = 0; i < cachedJacobian.length; ++i) {
  478.             System.arraycopy(cachedJacobian[i], 0, jacobian[i], 0, cachedJacobian[i].length);
  479.         }

  480.     }

  481.     /** Create an inverse Jacobian.
  482.      * @param type type of the position angle to use
  483.      * @return inverse Jacobian
  484.      */
  485.     private T[][] createInverseJacobian(final PositionAngle type) {

  486.         // get the direct Jacobian
  487.         final T[][] directJacobian = MathArrays.buildArray(getA().getField(), 6, 6);
  488.         getJacobianWrtCartesian(type, directJacobian);

  489.         // invert the direct Jacobian
  490.         final FieldMatrix<T> matrix = MatrixUtils.createFieldMatrix(directJacobian);
  491.         final FieldDecompositionSolver<T> solver = new FieldLUDecomposition<>(matrix).getSolver();
  492.         return solver.getInverse().getData();

  493.     }

  494.     /** Compute the Jacobian of the orbital parameters with mean angle with respect to the Cartesian parameters.
  495.      * <p>
  496.      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
  497.      * respect to Cartesian coordinate j. This means each row correspond to one orbital parameter
  498.      * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.
  499.      * </p>
  500.      * @return 6x6 Jacobian matrix
  501.      * @see #computeJacobianEccentricWrtCartesian()
  502.      * @see #computeJacobianTrueWrtCartesian()
  503.      */
  504.     protected abstract T[][] computeJacobianMeanWrtCartesian();

  505.     /** Compute the Jacobian of the orbital parameters with eccentric angle with respect to the Cartesian parameters.
  506.      * <p>
  507.      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
  508.      * respect to Cartesian coordinate j. This means each row correspond to one orbital parameter
  509.      * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.
  510.      * </p>
  511.      * @return 6x6 Jacobian matrix
  512.      * @see #computeJacobianMeanWrtCartesian()
  513.      * @see #computeJacobianTrueWrtCartesian()
  514.      */
  515.     protected abstract T[][] computeJacobianEccentricWrtCartesian();

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

  527.     /** Add the contribution of the Keplerian motion to parameters derivatives
  528.      * <p>
  529.      * This method is used by integration-based propagators to evaluate the part of Keplerian
  530.      * motion to evolution of the orbital state.
  531.      * </p>
  532.      * @param type type of the position angle in the state
  533.      * @param gm attraction coefficient to use
  534.      * @param pDot array containing orbital state derivatives to update (the Keplerian
  535.      * part must be <em>added</em> to the array components, as the array may already
  536.      * contain some non-zero elements corresponding to non-Keplerian parts)
  537.      */
  538.     public abstract void addKeplerContribution(PositionAngle type, T gm, T[] pDot);

  539.         /** Fill a Jacobian half row with a single vector.
  540.      * @param a coefficient of the vector
  541.      * @param v vector
  542.      * @param row Jacobian matrix row
  543.      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
  544.      */

  545.     protected void fillHalfRow(final T a, final FieldVector3D<T> v, final T[] row, final int j) {
  546.         row[j]     = a.multiply(v.getX());
  547.         row[j + 1] = a.multiply(v.getY());
  548.         row[j + 2] = a.multiply(v.getZ());
  549.     }

  550.     /** Fill a Jacobian half row with a linear combination of vectors.
  551.      * @param a1 coefficient of the first vector
  552.      * @param v1 first vector
  553.      * @param a2 coefficient of the second vector
  554.      * @param v2 second vector
  555.      * @param row Jacobian matrix row
  556.      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
  557.      */
  558.     protected void fillHalfRow(final T a1, final FieldVector3D<T> v1, final T a2, final FieldVector3D<T> v2,
  559.                                       final T[] row, final int j) {
  560.         row[j]     = a1.linearCombination(a1, v1.getX(), a2, v2.getX());
  561.         row[j + 1] = a1.linearCombination(a1, v1.getY(), a2, v2.getY());
  562.         row[j + 2] = a1.linearCombination(a1, v1.getZ(), a2, v2.getZ());
  563. //        row[j]     = a1.multiply(v1.getX()).add(a2.multiply(v2.getX()));
  564. //        row[j + 1] = a1.multiply(v1.getY()).add(a2.multiply(v2.getY()));
  565. //        row[j + 2] = a1.multiply(v1.getZ()).add(a2.multiply(v2.getZ()));
  566.     }

  567.     /** Fill a Jacobian half row with a linear combination of vectors.
  568.      * @param a1 coefficient of the first vector
  569.      * @param v1 first vector
  570.      * @param a2 coefficient of the second vector
  571.      * @param v2 second vector
  572.      * @param a3 coefficient of the third vector
  573.      * @param v3 third vector
  574.      * @param row Jacobian matrix row
  575.      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
  576.      */
  577.     protected void fillHalfRow(final T a1, final FieldVector3D<T> v1,
  578.                                final T a2, final FieldVector3D<T> v2,
  579.                                final T a3, final FieldVector3D<T> v3,
  580.                                       final T[] row, final int j) {
  581.         row[j]     = a1.linearCombination(a1, v1.getX(), a2, v2.getX(), a3, v3.getX());
  582.         row[j + 1] = a1.linearCombination(a1, v1.getY(), a2, v2.getY(), a3, v3.getY());
  583.         row[j + 2] = a1.linearCombination(a1, v1.getZ(), a2, v2.getZ(), a3, v3.getZ());



  584. //        row[j]     = a1.multiply(v1.getX()).add(a2.multiply(v2.getX())).add(a3.multiply(v3.getX()));
  585. //        row[j + 1] = a1.multiply(v1.getY()).add(a2.multiply(v2.getY())).add(a3.multiply(v3.getY()));
  586. //        row[j + 2] = a1.multiply(v1.getZ()).add(a2.multiply(v2.getZ())).add(a3.multiply(v3.getZ()));
  587.     }

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

  610.     /** Fill a Jacobian half row with a linear combination of vectors.
  611.      * @param a1 coefficient of the first vector
  612.      * @param v1 first vector
  613.      * @param a2 coefficient of the second vector
  614.      * @param v2 second vector
  615.      * @param a3 coefficient of the third vector
  616.      * @param v3 third vector
  617.      * @param a4 coefficient of the fourth vector
  618.      * @param v4 fourth vector
  619.      * @param a5 coefficient of the fifth vector
  620.      * @param v5 fifth vector
  621.      * @param row Jacobian matrix row
  622.      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
  623.      */
  624.     protected void fillHalfRow(final T a1, final FieldVector3D<T> v1, final T a2, final FieldVector3D<T> v2,
  625.                                       final T a3, final FieldVector3D<T> v3, final T a4, final FieldVector3D<T> v4,
  626.                                       final T a5, final FieldVector3D<T> v5,
  627.                                       final T[] row, final int j) {
  628.         row[j]     = a1.linearCombination(a1, v1.getX(), a2, v2.getX(), a3, v3.getX(), a4, v4.getX()).add(a5.multiply(v5.getX()));
  629.         row[j + 1] = a1.linearCombination(a1, v1.getY(), a2, v2.getY(), a3, v3.getY(), a4, v4.getY()).add(a5.multiply(v5.getY()));
  630.         row[j + 2] = a1.linearCombination(a1, v1.getZ(), a2, v2.getZ(), a3, v3.getZ(), a4, v4.getZ()).add(a5.multiply(v5.getZ()));

  631. //        row[j]     = a1.multiply(v1.getX()).add(a2.multiply(v2.getX())).add(a3.multiply(v3.getX())).add(a4.multiply(v4.getX())).add(a5.multiply(v5.getX()));
  632. //        row[j + 1] = a1.multiply(v1.getY()).add(a2.multiply(v2.getY())).add(a3.multiply(v3.getY())).add(a4.multiply(v4.getY())).add(a5.multiply(v5.getY()));
  633. //        row[j + 2] = a1.multiply(v1.getZ()).add(a2.multiply(v2.getZ())).add(a3.multiply(v3.getZ())).add(a4.multiply(v4.getZ())).add(a5.multiply(v5.getZ()));
  634.     }

  635.     /** Fill a Jacobian half row with a linear combination of vectors.
  636.      * @param a1 coefficient of the first vector
  637.      * @param v1 first vector
  638.      * @param a2 coefficient of the second vector
  639.      * @param v2 second vector
  640.      * @param a3 coefficient of the third vector
  641.      * @param v3 third vector
  642.      * @param a4 coefficient of the fourth vector
  643.      * @param v4 fourth vector
  644.      * @param a5 coefficient of the fifth vector
  645.      * @param v5 fifth vector
  646.      * @param a6 coefficient of the sixth vector
  647.      * @param v6 sixth vector
  648.      * @param row Jacobian matrix row
  649.      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
  650.      */
  651.     protected void fillHalfRow(final T a1, final FieldVector3D<T> v1, final T a2, final FieldVector3D<T> v2,
  652.                                       final T a3, final FieldVector3D<T> v3, final T a4, final FieldVector3D<T> v4,
  653.                                       final T a5, final FieldVector3D<T> v5, final T a6, final FieldVector3D<T> v6,
  654.                                       final T[] row, final int j) {
  655.         row[j]     = a1.linearCombination(a1, v1.getX(), a2, v2.getX(), a3, v3.getX(), a4, v4.getX()).add(a1.linearCombination(a5, v5.getX(), a6, v6.getX()));
  656.         row[j + 1] = a1.linearCombination(a1, v1.getY(), a2, v2.getY(), a3, v3.getY(), a4, v4.getY()).add(a1.linearCombination(a5, v5.getY(), a6, v6.getY()));
  657.         row[j + 2] = a1.linearCombination(a1, v1.getZ(), a2, v2.getZ(), a3, v3.getZ(), a4, v4.getZ()).add(a1.linearCombination(a5, v5.getZ(), a6, v6.getZ()));


  658. //        row[j]     = a1.multiply(v1.getX()).add(a2.multiply(v2.getX())).add(a3.multiply(v3.getX())).add(a4.multiply(v4.getX())).add(a5.multiply(v5.getX())).add(a6.multiply(v6.getX()));
  659. //        row[j + 1] = a1.multiply(v1.getY()).add(a2.multiply(v2.getY())).add(a3.multiply(v3.getY())).add(a4.multiply(v4.getY())).add(a5.multiply(v5.getY())).add(a6.multiply(v6.getY()));
  660. //        row[j + 2] = a1.multiply(v1.getZ()).add(a2.multiply(v2.getZ())).add(a3.multiply(v3.getZ())).add(a4.multiply(v4.getZ())).add(a5.multiply(v5.getZ())).add(a6.multiply(v6.getZ()));
  661.     }


  662. }