FieldOrbit.java

  1. /* Copyright 2002-2019 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 org.hipparchus.RealFieldElement;
  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.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.FieldAbsoluteDate;
  32. import org.orekit.time.FieldTimeInterpolable;
  33. import org.orekit.time.FieldTimeShiftable;
  34. import org.orekit.time.FieldTimeStamped;
  35. import org.orekit.utils.FieldPVCoordinates;
  36. import org.orekit.utils.FieldPVCoordinatesProvider;
  37. import org.orekit.utils.TimeStampedFieldPVCoordinates;
  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.  * @since 9.0
  61.  */
  62. public abstract class FieldOrbit<T extends RealFieldElement<T>>
  63.     implements FieldPVCoordinatesProvider<T>, FieldTimeStamped<T>, FieldTimeShiftable<FieldOrbit<T>, T>, FieldTimeInterpolable<FieldOrbit<T>, T> {

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

  162.     }

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

  276.     // Additional orbital elements

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

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

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

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

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

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

  317.     public double getMu() {
  318.         return mu;
  319.     }

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

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

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

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

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

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

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


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

  373.         }
  374.         return pvCoordinates;
  375.     }

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

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

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

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

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

  434.     }

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

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

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

  478.     }

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

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

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

  491.     }

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

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

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

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

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

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

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

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



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

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

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

  629. //        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()));
  630. //        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()));
  631. //        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()));
  632.     }

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


  656. //        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()));
  657. //        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()));
  658. //        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()));
  659.     }


  660. }