FieldOrbit.java

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



  18. import org.hipparchus.CalculusFieldElement;
  19. import org.hipparchus.Field;
  20. import org.hipparchus.analysis.differentiation.Derivative;
  21. import org.hipparchus.geometry.euclidean.threed.FieldVector3D;
  22. import org.hipparchus.linear.FieldDecompositionSolver;
  23. import org.hipparchus.linear.FieldLUDecomposition;
  24. import org.hipparchus.linear.FieldMatrix;
  25. import org.hipparchus.linear.MatrixUtils;
  26. import org.hipparchus.util.MathArrays;
  27. import org.orekit.errors.OrekitIllegalArgumentException;
  28. import org.orekit.errors.OrekitInternalError;
  29. import org.orekit.errors.OrekitMessages;
  30. import org.orekit.frames.FieldStaticTransform;
  31. import org.orekit.frames.FieldTransform;
  32. import org.orekit.frames.Frame;
  33. import org.orekit.time.FieldAbsoluteDate;
  34. import org.orekit.time.FieldTimeShiftable;
  35. import org.orekit.time.FieldTimeStamped;
  36. import org.orekit.utils.FieldPVCoordinates;
  37. import org.orekit.utils.FieldPVCoordinatesProvider;
  38. import org.orekit.utils.TimeStampedFieldPVCoordinates;
  39. import org.orekit.utils.TimeStampedPVCoordinates;

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

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

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

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

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

  73.     /** Field-valued one.
  74.      * @since 12.0
  75.      * */
  76.     private final T one;

  77.     /** Field-valued zero.
  78.      * @since 12.0
  79.      * */
  80.     private final T zero;

  81.     /** Field used by this class.
  82.      * @since 12.0
  83.      */
  84.     private final Field<T> field;

  85.     /** Computed position.
  86.      * @since 12.0
  87.      */
  88.     private transient FieldVector3D<T> position;

  89.     /** Computed PVCoordinates. */
  90.     private transient TimeStampedFieldPVCoordinates<T> pvCoordinates;

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

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

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

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

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

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

  103.     /** Default constructor.
  104.      * Build a new instance with arbitrary default elements.
  105.      * @param frame the frame in which the parameters are defined
  106.      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
  107.      * @param date date of the orbital parameters
  108.      * @param mu central attraction coefficient (m^3/s^2)
  109.      * @exception IllegalArgumentException if frame is not a {@link
  110.      * Frame#isPseudoInertial pseudo-inertial frame}
  111.      */
  112.     protected FieldOrbit(final Frame frame, final FieldAbsoluteDate<T> date, final T mu)
  113.         throws IllegalArgumentException {
  114.         ensurePseudoInertialFrame(frame);
  115.         this.field = mu.getField();
  116.         this.zero = this.field.getZero();
  117.         this.one = this.field.getOne();
  118.         this.date                      = date;
  119.         this.mu                        = mu;
  120.         this.pvCoordinates             = null;
  121.         this.frame                     = frame;
  122.         jacobianMeanWrtCartesian       = null;
  123.         jacobianWrtParametersMean      = null;
  124.         jacobianEccentricWrtCartesian  = null;
  125.         jacobianWrtParametersEccentric = null;
  126.         jacobianTrueWrtCartesian       = null;
  127.         jacobianWrtParametersTrue      = null;
  128.     }

  129.     /** Set the orbit from Cartesian parameters.
  130.      *
  131.      * <p> The acceleration provided in {@code FieldPVCoordinates} is accessible using
  132.      * {@link #getPVCoordinates()} and {@link #getPVCoordinates(Frame)}. All other methods
  133.      * use {@code mu} and the position to compute the acceleration, including
  134.      * {@link #shiftedBy(CalculusFieldElement)} and {@link #getPVCoordinates(FieldAbsoluteDate, Frame)}.
  135.      *
  136.      * @param fieldPVCoordinates the position and velocity in the inertial frame
  137.      * @param frame the frame in which the {@link TimeStampedPVCoordinates} are defined
  138.      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
  139.      * @param mu central attraction coefficient (m^3/s^2)
  140.      * @exception IllegalArgumentException if frame is not a {@link
  141.      * Frame#isPseudoInertial pseudo-inertial frame}
  142.      */
  143.     protected FieldOrbit(final TimeStampedFieldPVCoordinates<T> fieldPVCoordinates, final Frame frame, final T mu)
  144.         throws IllegalArgumentException {
  145.         ensurePseudoInertialFrame(frame);
  146.         this.field = mu.getField();
  147.         this.zero = this.field.getZero();
  148.         this.one = this.field.getOne();
  149.         this.date = fieldPVCoordinates.getDate();
  150.         this.mu = mu;
  151.         if (fieldPVCoordinates.getAcceleration().getNormSq().getReal() == 0.0) {
  152.             // the acceleration was not provided,
  153.             // compute it from Newtonian attraction
  154.             final T r2 = fieldPVCoordinates.getPosition().getNormSq();
  155.             final T r3 = r2.multiply(r2.sqrt());
  156.             this.pvCoordinates = new TimeStampedFieldPVCoordinates<>(fieldPVCoordinates.getDate(),
  157.                                                                      fieldPVCoordinates.getPosition(),
  158.                                                                      fieldPVCoordinates.getVelocity(),
  159.                                                                      new FieldVector3D<>(r3.reciprocal().multiply(mu.negate()), fieldPVCoordinates.getPosition()));
  160.         } else {
  161.             this.pvCoordinates = fieldPVCoordinates;
  162.         }
  163.         this.frame = frame;
  164.     }

  165.     /** Check if Cartesian coordinates include non-Keplerian acceleration.
  166.      * @param pva Cartesian coordinates
  167.      * @param mu central attraction coefficient
  168.      * @param <T> type of the field elements
  169.      * @return true if Cartesian coordinates include non-Keplerian acceleration
  170.      */
  171.     protected static <T extends CalculusFieldElement<T>> boolean hasNonKeplerianAcceleration(final FieldPVCoordinates<T> pva, final T mu) {

  172.         if (mu.getField().getZero() instanceof Derivative<?>) {
  173.             return Orbit.hasNonKeplerianAcceleration(pva.toPVCoordinates(), mu.getReal()); // for performance

  174.         } else {
  175.             final FieldVector3D<T> a = pva.getAcceleration();
  176.             if (a == null) {
  177.                 return false;
  178.             }

  179.             final FieldVector3D<T> p = pva.getPosition();
  180.             final T r2 = p.getNormSq();

  181.             // Check if acceleration is relatively close to 0 compared to the Keplerian acceleration
  182.             final double tolerance = mu.getReal() * 1e-9;
  183.             final FieldVector3D<T> aTimesR2 = a.scalarMultiply(r2);
  184.             if (aTimesR2.getNorm().getReal() < tolerance) {
  185.                 return false;
  186.             }

  187.             if ((aTimesR2.add(p.normalize().scalarMultiply(mu))).getNorm().getReal() > tolerance) {
  188.                 // we have a relevant acceleration, we can compute derivatives
  189.                 return true;
  190.             } else {
  191.                 // the provided acceleration is either too small to be reliable (probably even 0), or NaN
  192.                 return false;
  193.             }
  194.         }

  195.     }

  196.     /** Returns true if and only if the orbit is elliptical i.e. has a non-negative semi-major axis.
  197.      * @return true if getA() is strictly greater than 0
  198.      * @since 12.0
  199.      */
  200.     public boolean isElliptical() {
  201.         return getA().getReal() > 0.;
  202.     }

  203.     /** Get the orbit type.
  204.      * @return orbit type
  205.      */
  206.     public abstract OrbitType getType();

  207.     /** Ensure the defining frame is a pseudo-inertial frame.
  208.      * @param frame frame to check
  209.      * @exception IllegalArgumentException if frame is not a {@link
  210.      * Frame#isPseudoInertial pseudo-inertial frame}
  211.      */
  212.     private static void ensurePseudoInertialFrame(final Frame frame)
  213.         throws IllegalArgumentException {
  214.         if (!frame.isPseudoInertial()) {
  215.             throw new OrekitIllegalArgumentException(OrekitMessages.NON_PSEUDO_INERTIAL_FRAME,
  216.                                                      frame.getName());
  217.         }
  218.     }

  219.     /** Get the frame in which the orbital parameters are defined.
  220.      * @return frame in which the orbital parameters are defined
  221.      */
  222.     public Frame getFrame() {
  223.         return frame;
  224.     }

  225.     /**Transforms the FieldOrbit instance into an Orbit instance.
  226.      * @return Orbit instance with same properties
  227.      */
  228.     public abstract Orbit toOrbit();

  229.     /** Get the semi-major axis.
  230.      * <p>Note that the semi-major axis is considered negative for hyperbolic orbits.</p>
  231.      * @return semi-major axis (m)
  232.      */
  233.     public abstract T getA();

  234.     /** Get the semi-major axis derivative.
  235.      * <p>Note that the semi-major axis is considered negative for hyperbolic orbits.</p>
  236.      * <p>
  237.      * If the orbit was created without derivatives, the value returned is null.
  238.      * </p>
  239.      * @return semi-major axis  derivative (m/s)
  240.      */
  241.     public abstract T getADot();

  242.     /** Get the first component of the equinoctial eccentricity vector.
  243.      * @return first component of the equinoctial eccentricity vector
  244.      */
  245.     public abstract T getEquinoctialEx();

  246.     /** Get the first component of the equinoctial eccentricity vector derivative.
  247.      * <p>
  248.      * If the orbit was created without derivatives, the value returned is null.
  249.      * </p>
  250.      * @return first component of the equinoctial eccentricity vector derivative
  251.      */
  252.     public abstract T getEquinoctialExDot();

  253.     /** Get the second component of the equinoctial eccentricity vector.
  254.      * @return second component of the equinoctial eccentricity vector
  255.      */
  256.     public abstract T getEquinoctialEy();

  257.     /** Get the second component of the equinoctial eccentricity vector derivative.
  258.      * <p>
  259.      * If the orbit was created without derivatives, the value returned is null.
  260.      * </p>
  261.      * @return second component of the equinoctial eccentricity vector derivative
  262.      */
  263.     public abstract T getEquinoctialEyDot();

  264.     /** Get the first component of the inclination vector.
  265.      * @return first component of the inclination vector
  266.      */
  267.     public abstract T getHx();

  268.     /** Get the first component of the inclination vector derivative.
  269.      * <p>
  270.      * If the orbit was created without derivatives, the value returned is null.
  271.      * </p>
  272.      * @return first component of the inclination vector derivative
  273.      */
  274.     public abstract T getHxDot();

  275.     /** Get the second component of the inclination vector.
  276.      * @return second component of the inclination vector
  277.      */
  278.     public abstract T getHy();

  279.     /** Get the second component of the inclination vector derivative.
  280.      * @return second component of the inclination vector derivative
  281.      */
  282.     public abstract T getHyDot();

  283.     /** Get the eccentric longitude argument.
  284.      * @return E + ω + Ω eccentric longitude argument (rad)
  285.      */
  286.     public abstract T getLE();

  287.     /** Get the eccentric longitude argument derivative.
  288.      * <p>
  289.      * If the orbit was created without derivatives, the value returned is null.
  290.      * </p>
  291.      * @return d(E + ω + Ω)/dt eccentric longitude argument derivative (rad/s)
  292.      */
  293.     public abstract T getLEDot();

  294.     /** Get the true longitude argument.
  295.      * @return v + ω + Ω true longitude argument (rad)
  296.      */
  297.     public abstract T getLv();

  298.     /** Get the true longitude argument derivative.
  299.      * <p>
  300.      * If the orbit was created without derivatives, the value returned is null.
  301.      * </p>
  302.      * @return d(v + ω + Ω)/dt true longitude argument derivative (rad/s)
  303.      */
  304.     public abstract T getLvDot();

  305.     /** Get the mean longitude argument.
  306.      * @return M + ω + Ω mean longitude argument (rad)
  307.      */
  308.     public abstract T getLM();

  309.     /** Get the mean longitude argument derivative.
  310.      * <p>
  311.      * If the orbit was created without derivatives, the value returned is null.
  312.      * </p>
  313.      * @return d(M + ω + Ω)/dt mean longitude argument derivative (rad/s)
  314.      */
  315.     public abstract T getLMDot();

  316.     // Additional orbital elements

  317.     /** Get the eccentricity.
  318.      * @return eccentricity
  319.      */
  320.     public abstract T getE();

  321.     /** Get the eccentricity derivative.
  322.      * <p>
  323.      * If the orbit was created without derivatives, the value returned is null.
  324.      * </p>
  325.      * @return eccentricity derivative
  326.      */
  327.     public abstract T getEDot();

  328.     /** Get the inclination.
  329.      * <p>
  330.      * If the orbit was created without derivatives, the value returned is null.
  331.      * </p>
  332.      * @return inclination (rad)
  333.      */
  334.     public abstract T getI();

  335.     /** Get the inclination derivative.
  336.      * @return inclination derivative (rad/s)
  337.      */
  338.     public abstract T getIDot();

  339.     /** Check if orbit includes derivatives.
  340.      * @return true if orbit includes derivatives
  341.      * @see #getADot()
  342.      * @see #getEquinoctialExDot()
  343.      * @see #getEquinoctialEyDot()
  344.      * @see #getHxDot()
  345.      * @see #getHyDot()
  346.      * @see #getLEDot()
  347.      * @see #getLvDot()
  348.      * @see #getLMDot()
  349.      * @see #getEDot()
  350.      * @see #getIDot()
  351.      * @since 9.0
  352.      */
  353.     public abstract boolean hasDerivatives();

  354.     /** Get the central attraction coefficient used for position and velocity conversions (m³/s²).
  355.      * @return central attraction coefficient used for position and velocity conversions (m³/s²)
  356.      */
  357.     public T getMu() {
  358.         return mu;
  359.     }

  360.     /** Get the Keplerian period.
  361.      * <p>The Keplerian period is computed directly from semi major axis
  362.      * and central acceleration constant.</p>
  363.      * @return Keplerian period in seconds, or positive infinity for hyperbolic orbits
  364.      */
  365.     public T getKeplerianPeriod() {
  366.         final T a = getA();
  367.         return isElliptical() ?
  368.                a.multiply(a.getPi().multiply(2.0)).multiply(a.divide(mu).sqrt()) :
  369.                zero.add(Double.POSITIVE_INFINITY);
  370.     }

  371.     /** Get the Keplerian mean motion.
  372.      * <p>The Keplerian mean motion is computed directly from semi major axis
  373.      * and central acceleration constant.</p>
  374.      * @return Keplerian mean motion in radians per second
  375.      */
  376.     public T getKeplerianMeanMotion() {
  377.         final T absA = getA().abs();
  378.         return absA.reciprocal().multiply(mu).sqrt().divide(absA);
  379.     }

  380.     /** Get the derivative of the mean anomaly with respect to the semi major axis.
  381.      * @return derivative of the mean anomaly with respect to the semi major axis
  382.      */
  383.     public T getMeanAnomalyDotWrtA() {
  384.         return getKeplerianMeanMotion().divide(getA()).multiply(-1.5);
  385.     }

  386.     /** Get the date of orbital parameters.
  387.      * @return date of the orbital parameters
  388.      */
  389.     public FieldAbsoluteDate<T> getDate() {
  390.         return date;
  391.     }

  392.     /** Get the {@link TimeStampedPVCoordinates} in a specified frame.
  393.      * @param outputFrame frame in which the position/velocity coordinates shall be computed
  394.      * @return FieldPVCoordinates in the specified output frame
  395.           * @see #getPVCoordinates()
  396.      */
  397.     public TimeStampedFieldPVCoordinates<T> getPVCoordinates(final Frame outputFrame) {
  398.         if (pvCoordinates == null) {
  399.             pvCoordinates = initPVCoordinates();
  400.         }

  401.         // If output frame requested is the same as definition frame,
  402.         // PV coordinates are returned directly
  403.         if (outputFrame == frame) {
  404.             return pvCoordinates;
  405.         }

  406.         // Else, PV coordinates are transformed to output frame
  407.         final FieldTransform<T> t = frame.getTransformTo(outputFrame, date);
  408.         return t.transformPVCoordinates(pvCoordinates);
  409.     }

  410.     /** {@inheritDoc} */
  411.     public TimeStampedFieldPVCoordinates<T> getPVCoordinates(final FieldAbsoluteDate<T> otherDate, final Frame otherFrame) {
  412.         return shiftedBy(otherDate.durationFrom(getDate())).getPVCoordinates(otherFrame);
  413.     }

  414.     /** {@inheritDoc} */
  415.     @Override
  416.     public FieldVector3D<T> getPosition(final FieldAbsoluteDate<T> otherDate, final Frame otherFrame) {
  417.         return shiftedBy(otherDate.durationFrom(getDate())).getPosition(otherFrame);
  418.     }

  419.     /** Get the position in a specified frame.
  420.      * @param outputFrame frame in which the position coordinates shall be computed
  421.      * @return position in the specified output frame
  422.      * @see #getPosition()
  423.      * @since 12.0
  424.      */
  425.     public FieldVector3D<T> getPosition(final Frame outputFrame) {
  426.         if (position == null) {
  427.             position = initPosition();
  428.         }

  429.         // If output frame requested is the same as definition frame,
  430.         // Position vector is returned directly
  431.         if (outputFrame == frame) {
  432.             return position;
  433.         }

  434.         // Else, position vector is transformed to output frame
  435.         final FieldStaticTransform<T> t = frame.getStaticTransformTo(outputFrame, date);
  436.         return t.transformPosition(position);

  437.     }

  438.     /** Get the position in definition frame.
  439.      * @return position in the definition frame
  440.      * @see #getPVCoordinates()
  441.      * @since 12.0
  442.      */
  443.     public FieldVector3D<T> getPosition() {
  444.         if (position == null) {
  445.             position = initPosition();
  446.         }
  447.         return position;
  448.     }

  449.     /** Get the {@link TimeStampedPVCoordinates} in definition frame.
  450.      * @return FieldPVCoordinates in the definition frame
  451.      * @see #getPVCoordinates(Frame)
  452.      */
  453.     public TimeStampedFieldPVCoordinates<T> getPVCoordinates() {
  454.         if (pvCoordinates == null) {
  455.             pvCoordinates = initPVCoordinates();
  456.             position      = pvCoordinates.getPosition();
  457.         }
  458.         return pvCoordinates;
  459.     }

  460.     /** Getter for Field-valued one.
  461.      * @return one
  462.      * */
  463.     protected T getOne() {
  464.         return one;
  465.     }

  466.     /** Getter for Field-valued zero.
  467.      * @return zero
  468.      * */
  469.     protected T getZero() {
  470.         return zero;
  471.     }

  472.     /** Getter for Field.
  473.      * @return CalculusField
  474.      * */
  475.     protected Field<T> getField() {
  476.         return field;
  477.     }

  478.     /** Compute the position coordinates from the canonical parameters.
  479.      * @return computed position coordinates
  480.      * @since 12.0
  481.      */
  482.     protected abstract FieldVector3D<T> initPosition();

  483.     /** Compute the position/velocity coordinates from the canonical parameters.
  484.      * @return computed position/velocity coordinates
  485.      */
  486.     protected abstract TimeStampedFieldPVCoordinates<T> initPVCoordinates();

  487.     /** Get a time-shifted orbit.
  488.      * <p>
  489.      * The orbit can be slightly shifted to close dates. This shift is based on
  490.      * a simple Keplerian model. It is <em>not</em> intended as a replacement
  491.      * for proper orbit and attitude propagation but should be sufficient for
  492.      * small time shifts or coarse accuracy.
  493.      * </p>
  494.      * @param dt time shift in seconds
  495.      * @return a new orbit, shifted with respect to the instance (which is immutable)
  496.      */
  497.     public abstract FieldOrbit<T> shiftedBy(T dt);

  498.     /** Compute the Jacobian of the orbital parameters with respect to the Cartesian parameters.
  499.      * <p>
  500.      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
  501.      * respect to Cartesian coordinate j. This means each row corresponds to one orbital parameter
  502.      * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.
  503.      * </p>
  504.      * @param type type of the position angle to use
  505.      * @param jacobian placeholder 6x6 (or larger) matrix to be filled with the Jacobian, if matrix
  506.      * is larger than 6x6, only the 6x6 upper left corner will be modified
  507.      */
  508.     public void getJacobianWrtCartesian(final PositionAngleType type, final T[][] jacobian) {

  509.         final T[][] cachedJacobian;
  510.         synchronized (this) {
  511.             switch (type) {
  512.                 case MEAN :
  513.                     if (jacobianMeanWrtCartesian == null) {
  514.                         // first call, we need to compute the Jacobian and cache it
  515.                         jacobianMeanWrtCartesian = computeJacobianMeanWrtCartesian();
  516.                     }
  517.                     cachedJacobian = jacobianMeanWrtCartesian;
  518.                     break;
  519.                 case ECCENTRIC :
  520.                     if (jacobianEccentricWrtCartesian == null) {
  521.                         // first call, we need to compute the Jacobian and cache it
  522.                         jacobianEccentricWrtCartesian = computeJacobianEccentricWrtCartesian();
  523.                     }
  524.                     cachedJacobian = jacobianEccentricWrtCartesian;
  525.                     break;
  526.                 case TRUE :
  527.                     if (jacobianTrueWrtCartesian == null) {
  528.                         // first call, we need to compute the Jacobian and cache it
  529.                         jacobianTrueWrtCartesian = computeJacobianTrueWrtCartesian();
  530.                     }
  531.                     cachedJacobian = jacobianTrueWrtCartesian;
  532.                     break;
  533.                 default :
  534.                     throw new OrekitInternalError(null);
  535.             }
  536.         }

  537.         // fill the user provided array
  538.         for (int i = 0; i < cachedJacobian.length; ++i) {
  539.             System.arraycopy(cachedJacobian[i], 0, jacobian[i], 0, cachedJacobian[i].length);
  540.         }

  541.     }

  542.     /** Compute the Jacobian of the Cartesian parameters with respect to the orbital parameters.
  543.      * <p>
  544.      * Element {@code jacobian[i][j]} is the derivative of Cartesian coordinate i of the orbit with
  545.      * respect to orbital parameter j. This means each row corresponds to one Cartesian coordinate
  546.      * x, y, z, xdot, ydot, zdot whereas columns 0 to 5 correspond to the orbital parameters.
  547.      * </p>
  548.      * @param type type of the position angle to use
  549.      * @param jacobian placeholder 6x6 (or larger) matrix to be filled with the Jacobian, if matrix
  550.      * is larger than 6x6, only the 6x6 upper left corner will be modified
  551.      */
  552.     public void getJacobianWrtParameters(final PositionAngleType type, final T[][] jacobian) {

  553.         final T[][] cachedJacobian;
  554.         synchronized (this) {
  555.             switch (type) {
  556.                 case MEAN :
  557.                     if (jacobianWrtParametersMean == null) {
  558.                         // first call, we need to compute the Jacobian and cache it
  559.                         jacobianWrtParametersMean = createInverseJacobian(type);
  560.                     }
  561.                     cachedJacobian = jacobianWrtParametersMean;
  562.                     break;
  563.                 case ECCENTRIC :
  564.                     if (jacobianWrtParametersEccentric == null) {
  565.                         // first call, we need to compute the Jacobian and cache it
  566.                         jacobianWrtParametersEccentric = createInverseJacobian(type);
  567.                     }
  568.                     cachedJacobian = jacobianWrtParametersEccentric;
  569.                     break;
  570.                 case TRUE :
  571.                     if (jacobianWrtParametersTrue == null) {
  572.                         // first call, we need to compute the Jacobian and cache it
  573.                         jacobianWrtParametersTrue = createInverseJacobian(type);
  574.                     }
  575.                     cachedJacobian = jacobianWrtParametersTrue;
  576.                     break;
  577.                 default :
  578.                     throw new OrekitInternalError(null);
  579.             }
  580.         }

  581.         // fill the user-provided array
  582.         for (int i = 0; i < cachedJacobian.length; ++i) {
  583.             System.arraycopy(cachedJacobian[i], 0, jacobian[i], 0, cachedJacobian[i].length);
  584.         }

  585.     }

  586.     /** Create an inverse Jacobian.
  587.      * @param type type of the position angle to use
  588.      * @return inverse Jacobian
  589.      */
  590.     private T[][] createInverseJacobian(final PositionAngleType type) {

  591.         // get the direct Jacobian
  592.         final T[][] directJacobian = MathArrays.buildArray(getA().getField(), 6, 6);
  593.         getJacobianWrtCartesian(type, directJacobian);

  594.         // invert the direct Jacobian
  595.         final FieldMatrix<T> matrix = MatrixUtils.createFieldMatrix(directJacobian);
  596.         final FieldDecompositionSolver<T> solver = new FieldLUDecomposition<>(matrix).getSolver();
  597.         return solver.getInverse().getData();

  598.     }

  599.     /** Compute the Jacobian of the orbital parameters with mean angle with respect to the Cartesian parameters.
  600.      * <p>
  601.      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
  602.      * respect to Cartesian coordinate j. This means each row correspond to one orbital parameter
  603.      * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.
  604.      * </p>
  605.      * @return 6x6 Jacobian matrix
  606.      * @see #computeJacobianEccentricWrtCartesian()
  607.      * @see #computeJacobianTrueWrtCartesian()
  608.      */
  609.     protected abstract T[][] computeJacobianMeanWrtCartesian();

  610.     /** Compute the Jacobian of the orbital parameters with eccentric angle with respect to the Cartesian parameters.
  611.      * <p>
  612.      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
  613.      * respect to Cartesian coordinate j. This means each row correspond to one orbital parameter
  614.      * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.
  615.      * </p>
  616.      * @return 6x6 Jacobian matrix
  617.      * @see #computeJacobianMeanWrtCartesian()
  618.      * @see #computeJacobianTrueWrtCartesian()
  619.      */
  620.     protected abstract T[][] computeJacobianEccentricWrtCartesian();

  621.     /** Compute the Jacobian of the orbital parameters with true angle with respect to the Cartesian parameters.
  622.      * <p>
  623.      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
  624.      * respect to Cartesian coordinate j. This means each row correspond to one orbital parameter
  625.      * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.
  626.      * </p>
  627.      * @return 6x6 Jacobian matrix
  628.      * @see #computeJacobianMeanWrtCartesian()
  629.      * @see #computeJacobianEccentricWrtCartesian()
  630.      */
  631.     protected abstract T[][] computeJacobianTrueWrtCartesian();

  632.     /** Add the contribution of the Keplerian motion to parameters derivatives
  633.      * <p>
  634.      * This method is used by integration-based propagators to evaluate the part of Keplerian
  635.      * motion to evolution of the orbital state.
  636.      * </p>
  637.      * @param type type of the position angle in the state
  638.      * @param gm attraction coefficient to use
  639.      * @param pDot array containing orbital state derivatives to update (the Keplerian
  640.      * part must be <em>added</em> to the array components, as the array may already
  641.      * contain some non-zero elements corresponding to non-Keplerian parts)
  642.      */
  643.     public abstract void addKeplerContribution(PositionAngleType type, T gm, T[] pDot);

  644.         /** Fill a Jacobian half row with a single vector.
  645.      * @param a coefficient of the vector
  646.      * @param v vector
  647.      * @param row Jacobian matrix row
  648.      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
  649.      */

  650.     protected void fillHalfRow(final T a, final FieldVector3D<T> v, final T[] row, final int j) {
  651.         row[j]     = a.multiply(v.getX());
  652.         row[j + 1] = a.multiply(v.getY());
  653.         row[j + 2] = a.multiply(v.getZ());
  654.     }

  655.     /** Fill a Jacobian half row with a linear combination of vectors.
  656.      * @param a1 coefficient of the first vector
  657.      * @param v1 first vector
  658.      * @param a2 coefficient of the second vector
  659.      * @param v2 second vector
  660.      * @param row Jacobian matrix row
  661.      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
  662.      */
  663.     protected void fillHalfRow(final T a1, final FieldVector3D<T> v1, final T a2, final FieldVector3D<T> v2,
  664.                                       final T[] row, final int j) {
  665.         row[j]     = a1.linearCombination(a1, v1.getX(), a2, v2.getX());
  666.         row[j + 1] = a1.linearCombination(a1, v1.getY(), a2, v2.getY());
  667.         row[j + 2] = a1.linearCombination(a1, v1.getZ(), a2, v2.getZ());

  668.     }

  669.     /** Fill a Jacobian half row with a linear combination of vectors.
  670.      * @param a1 coefficient of the first vector
  671.      * @param v1 first vector
  672.      * @param a2 coefficient of the second vector
  673.      * @param v2 second vector
  674.      * @param a3 coefficient of the third vector
  675.      * @param v3 third vector
  676.      * @param row Jacobian matrix row
  677.      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
  678.      */
  679.     protected void fillHalfRow(final T a1, final FieldVector3D<T> v1,
  680.                                final T a2, final FieldVector3D<T> v2,
  681.                                final T a3, final FieldVector3D<T> v3,
  682.                                       final T[] row, final int j) {
  683.         row[j]     = a1.linearCombination(a1, v1.getX(), a2, v2.getX(), a3, v3.getX());
  684.         row[j + 1] = a1.linearCombination(a1, v1.getY(), a2, v2.getY(), a3, v3.getY());
  685.         row[j + 2] = a1.linearCombination(a1, v1.getZ(), a2, v2.getZ(), a3, v3.getZ());

  686.     }

  687.     /** Fill a Jacobian half row with a linear combination of vectors.
  688.      * @param a1 coefficient of the first vector
  689.      * @param v1 first vector
  690.      * @param a2 coefficient of the second vector
  691.      * @param v2 second vector
  692.      * @param a3 coefficient of the third vector
  693.      * @param v3 third vector
  694.      * @param a4 coefficient of the fourth vector
  695.      * @param v4 fourth vector
  696.      * @param row Jacobian matrix row
  697.      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
  698.      */
  699.     protected void fillHalfRow(final T a1, final FieldVector3D<T> v1, final T a2, final FieldVector3D<T> v2,
  700.                                       final T a3, final FieldVector3D<T> v3, final T a4, final FieldVector3D<T> v4,
  701.                                       final T[] row, final int j) {
  702.         row[j]     = a1.linearCombination(a1, v1.getX(), a2, v2.getX(), a3, v3.getX(), a4, v4.getX());
  703.         row[j + 1] = a1.linearCombination(a1, v1.getY(), a2, v2.getY(), a3, v3.getY(), a4, v4.getY());
  704.         row[j + 2] = a1.linearCombination(a1, v1.getZ(), a2, v2.getZ(), a3, v3.getZ(), a4, v4.getZ());

  705.     }

  706.     /** Fill a Jacobian half row with a linear combination of vectors.
  707.      * @param a1 coefficient of the first vector
  708.      * @param v1 first vector
  709.      * @param a2 coefficient of the second vector
  710.      * @param v2 second vector
  711.      * @param a3 coefficient of the third vector
  712.      * @param v3 third vector
  713.      * @param a4 coefficient of the fourth vector
  714.      * @param v4 fourth vector
  715.      * @param a5 coefficient of the fifth vector
  716.      * @param v5 fifth vector
  717.      * @param row Jacobian matrix row
  718.      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
  719.      */
  720.     protected void fillHalfRow(final T a1, final FieldVector3D<T> v1, final T a2, final FieldVector3D<T> v2,
  721.                                       final T a3, final FieldVector3D<T> v3, final T a4, final FieldVector3D<T> v4,
  722.                                       final T a5, final FieldVector3D<T> v5,
  723.                                       final T[] row, final int j) {
  724.         row[j]     = a1.linearCombination(a1, v1.getX(), a2, v2.getX(), a3, v3.getX(), a4, v4.getX()).add(a5.multiply(v5.getX()));
  725.         row[j + 1] = a1.linearCombination(a1, v1.getY(), a2, v2.getY(), a3, v3.getY(), a4, v4.getY()).add(a5.multiply(v5.getY()));
  726.         row[j + 2] = a1.linearCombination(a1, v1.getZ(), a2, v2.getZ(), a3, v3.getZ(), a4, v4.getZ()).add(a5.multiply(v5.getZ()));

  727.     }

  728.     /** Fill a Jacobian half row with a linear combination of vectors.
  729.      * @param a1 coefficient of the first vector
  730.      * @param v1 first vector
  731.      * @param a2 coefficient of the second vector
  732.      * @param v2 second vector
  733.      * @param a3 coefficient of the third vector
  734.      * @param v3 third vector
  735.      * @param a4 coefficient of the fourth vector
  736.      * @param v4 fourth vector
  737.      * @param a5 coefficient of the fifth vector
  738.      * @param v5 fifth vector
  739.      * @param a6 coefficient of the sixth vector
  740.      * @param v6 sixth vector
  741.      * @param row Jacobian matrix row
  742.      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
  743.      */
  744.     protected void fillHalfRow(final T a1, final FieldVector3D<T> v1, final T a2, final FieldVector3D<T> v2,
  745.                                       final T a3, final FieldVector3D<T> v3, final T a4, final FieldVector3D<T> v4,
  746.                                       final T a5, final FieldVector3D<T> v5, final T a6, final FieldVector3D<T> v6,
  747.                                       final T[] row, final int j) {
  748.         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()));
  749.         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()));
  750.         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()));

  751.     }


  752. }