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  
19  import java.io.Serializable;
20  
21  import org.hipparchus.geometry.euclidean.threed.Vector3D;
22  import org.hipparchus.linear.DecompositionSolver;
23  import org.hipparchus.linear.MatrixUtils;
24  import org.hipparchus.linear.QRDecomposition;
25  import org.hipparchus.linear.RealMatrix;
26  import org.hipparchus.util.FastMath;
27  import org.hipparchus.util.MathArrays;
28  import org.orekit.errors.OrekitIllegalArgumentException;
29  import org.orekit.errors.OrekitInternalError;
30  import org.orekit.errors.OrekitMessages;
31  import org.orekit.frames.Frame;
32  import org.orekit.frames.StaticTransform;
33  import org.orekit.frames.Transform;
34  import org.orekit.time.AbsoluteDate;
35  import org.orekit.time.TimeOffset;
36  import org.orekit.time.TimeShiftable;
37  import org.orekit.time.TimeStamped;
38  import org.orekit.utils.PVCoordinates;
39  import org.orekit.utils.PVCoordinatesProvider;
40  import org.orekit.utils.TimeStampedPVCoordinates;
41  
42  /**
43   * This class handles orbital parameters.
44  
45   * <p>
46   * For user convenience, both the Cartesian and the equinoctial elements
47   * are provided by this class, regardless of the canonical representation
48   * implemented in the derived class (which may be classical Keplerian
49   * elements for example).
50   * </p>
51   * <p>
52   * The parameters are defined in a frame specified by the user. It is important
53   * to make sure this frame is consistent: it probably is inertial and centered
54   * on the central body. This information is used for example by some
55   * force models.
56   * </p>
57   * <p>
58   * Instance of this class are guaranteed to be immutable.
59   * </p>
60   * @author Luc Maisonobe
61   * @author Guylaine Prat
62   * @author Fabien Maussion
63   * @author V&eacute;ronique Pommier-Maurussane
64   */
65  public abstract class Orbit
66      implements TimeStamped, TimeShiftable<Orbit>, Serializable, PVCoordinatesProvider {
67  
68      /** Absolute tolerance when checking if the rate of the position angle is Keplerian or not. */
69      protected static final double TOLERANCE_POSITION_ANGLE_RATE = 1e-15;
70  
71      /** Serializable UID. */
72      private static final long serialVersionUID = 438733454597999578L;
73  
74      /** Frame in which are defined the orbital parameters. */
75      private final Frame frame;
76  
77      /** Date of the orbital parameters. */
78      private final AbsoluteDate date;
79  
80      /** Value of mu used to compute position and velocity (m³/s²). */
81      private final double mu;
82  
83      /** Computed position.
84       * @since 12.0
85       */
86      private transient Vector3D position;
87  
88      /** Computed PVCoordinates. */
89      private transient TimeStampedPVCoordinates pvCoordinates;
90  
91      /** Jacobian of the orbital parameters with mean angle with respect to the Cartesian coordinates. */
92      private transient double[][] jacobianMeanWrtCartesian;
93  
94      /** Jacobian of the Cartesian coordinates with respect to the orbital parameters with mean angle. */
95      private transient double[][] jacobianWrtParametersMean;
96  
97      /** Jacobian of the orbital parameters with eccentric angle with respect to the Cartesian coordinates. */
98      private transient double[][] jacobianEccentricWrtCartesian;
99  
100     /** Jacobian of the Cartesian coordinates with respect to the orbital parameters with eccentric angle. */
101     private transient double[][] jacobianWrtParametersEccentric;
102 
103     /** Jacobian of the orbital parameters with true angle with respect to the Cartesian coordinates. */
104     private transient double[][] jacobianTrueWrtCartesian;
105 
106     /** Jacobian of the Cartesian coordinates with respect to the orbital parameters with true angle. */
107     private transient double[][] jacobianWrtParametersTrue;
108 
109     /** Default constructor.
110      * Build a new instance with arbitrary default elements.
111      * @param frame the frame in which the parameters are defined
112      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
113      * @param date date of the orbital parameters
114      * @param mu central attraction coefficient (m^3/s^2)
115      * @exception IllegalArgumentException if frame is not a {@link
116      * Frame#isPseudoInertial pseudo-inertial frame}
117      */
118     protected Orbit(final Frame frame, final AbsoluteDate date, final double mu)
119         throws IllegalArgumentException {
120         ensurePseudoInertialFrame(frame);
121         this.date                      = date;
122         this.mu                        = mu;
123         this.pvCoordinates             = null;
124         this.frame                     = frame;
125         jacobianMeanWrtCartesian       = null;
126         jacobianWrtParametersMean      = null;
127         jacobianEccentricWrtCartesian  = null;
128         jacobianWrtParametersEccentric = null;
129         jacobianTrueWrtCartesian       = null;
130         jacobianWrtParametersTrue      = null;
131     }
132 
133     /** Set the orbit from Cartesian parameters.
134      *
135      * <p> The acceleration provided in {@code pvCoordinates} is accessible using
136      * {@link #getPVCoordinates()} and {@link #getPVCoordinates(Frame)}. All other methods
137      * use {@code mu} and the position to compute the acceleration, including
138      * {@link #shiftedBy(double)} and {@link #getPVCoordinates(AbsoluteDate, Frame)}.
139      *
140      * @param pvCoordinates the position and velocity in the inertial frame
141      * @param frame the frame in which the {@link TimeStampedPVCoordinates} are defined
142      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
143      * @param mu central attraction coefficient (m^3/s^2)
144      * @exception IllegalArgumentException if frame is not a {@link
145      * Frame#isPseudoInertial pseudo-inertial frame}
146      */
147     protected Orbit(final TimeStampedPVCoordinates pvCoordinates, final Frame frame, final double mu)
148         throws IllegalArgumentException {
149         ensurePseudoInertialFrame(frame);
150         this.date = pvCoordinates.getDate();
151         this.mu = mu;
152         if (pvCoordinates.getAcceleration().getNormSq() == 0) {
153             // the acceleration was not provided,
154             // compute it from Newtonian attraction
155             final double r2 = pvCoordinates.getPosition().getNormSq();
156             final double r3 = r2 * FastMath.sqrt(r2);
157             this.pvCoordinates = new TimeStampedPVCoordinates(pvCoordinates.getDate(),
158                                                               pvCoordinates.getPosition(),
159                                                               pvCoordinates.getVelocity(),
160                                                               new Vector3D(-mu / r3, pvCoordinates.getPosition()));
161         } else {
162             this.pvCoordinates = pvCoordinates;
163         }
164         this.frame = frame;
165     }
166 
167     /** Check if Cartesian coordinates include non-Keplerian acceleration.
168      * @param pva Cartesian coordinates
169      * @param mu central attraction coefficient
170      * @return true if Cartesian coordinates include non-Keplerian acceleration
171      */
172     protected static boolean hasNonKeplerianAcceleration(final PVCoordinates pva, final double mu) {
173 
174         final Vector3D a = pva.getAcceleration();
175         if (a == null) {
176             return false;
177         }
178 
179         final Vector3D p = pva.getPosition();
180         final double r2 = p.getNormSq();
181 
182         // Check if acceleration is relatively close to 0 compared to the Keplerian acceleration
183         final double tolerance = mu * 1e-9;
184         final Vector3D aTimesR2 = a.scalarMultiply(r2);
185         if (aTimesR2.getNorm() < tolerance) {
186             return false;
187         }
188 
189         if ((aTimesR2.add(p.normalize().scalarMultiply(mu))).getNorm() > tolerance) {
190             // we have a relevant acceleration, we can compute derivatives
191             return true;
192         } else {
193             // the provided acceleration is either too small to be reliable (probably even 0), or NaN
194             return false;
195         }
196     }
197 
198     /** Returns true if and only if the orbit is elliptical i.e. has a non-negative semi-major axis.
199      * @return true if getA() is strictly greater than 0
200      * @since 12.0
201      */
202     public boolean isElliptical() {
203         return getA() > 0.;
204     }
205 
206     /** Get the orbit type.
207      * @return orbit type
208      */
209     public abstract OrbitType getType();
210 
211     /** Ensure the defining frame is a pseudo-inertial frame.
212      * @param frame frame to check
213      * @exception IllegalArgumentException if frame is not a {@link
214      * Frame#isPseudoInertial pseudo-inertial frame}
215      */
216     private static void ensurePseudoInertialFrame(final Frame frame)
217         throws IllegalArgumentException {
218         if (!frame.isPseudoInertial()) {
219             throw new OrekitIllegalArgumentException(OrekitMessages.NON_PSEUDO_INERTIAL_FRAME,
220                                                      frame.getName());
221         }
222     }
223 
224     /** Get the frame in which the orbital parameters are defined.
225      * @return frame in which the orbital parameters are defined
226      */
227     public Frame getFrame() {
228         return frame;
229     }
230 
231     /** Get the semi-major axis.
232      * <p>Note that the semi-major axis is considered negative for hyperbolic orbits.</p>
233      * @return semi-major axis (m)
234      */
235     public abstract double getA();
236 
237     /** Get the semi-major axis derivative.
238      * <p>Note that the semi-major axis is considered negative for hyperbolic orbits.</p>
239      * <p>
240      * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
241      * </p>
242      * @return semi-major axis  derivative (m/s)
243      * @since 9.0
244      */
245     public abstract double getADot();
246 
247     /** Get the first component of the equinoctial eccentricity vector.
248      * @return first component of the equinoctial eccentricity vector
249      */
250     public abstract double getEquinoctialEx();
251 
252     /** Get the first component of the equinoctial eccentricity vector derivative.
253      * <p>
254      * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
255      * </p>
256      * @return first component of the equinoctial eccentricity vector derivative
257      * @since 9.0
258      */
259     public abstract double getEquinoctialExDot();
260 
261     /** Get the second component of the equinoctial eccentricity vector.
262      * @return second component of the equinoctial eccentricity vector
263      */
264     public abstract double getEquinoctialEy();
265 
266     /** Get the second component of the equinoctial eccentricity vector derivative.
267      * <p>
268      * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
269      * </p>
270      * @return second component of the equinoctial eccentricity vector derivative
271      * @since 9.0
272      */
273     public abstract double getEquinoctialEyDot();
274 
275     /** Get the first component of the inclination vector.
276      * @return first component of the inclination vector
277      */
278     public abstract double getHx();
279 
280     /** Get the first component of the inclination vector derivative.
281      * <p>
282      * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
283      * </p>
284      * @return first component of the inclination vector derivative
285 
286      * @since 9.0
287      */
288     public abstract double getHxDot();
289 
290     /** Get the second component of the inclination vector.
291      * @return second component of the inclination vector
292      */
293     public abstract double getHy();
294 
295     /** Get the second component of the inclination vector derivative.
296      * <p>
297      * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
298      * </p>
299      * @return second component of the inclination vector derivative
300      * @since 9.0
301      */
302     public abstract double getHyDot();
303 
304     /** Get the eccentric longitude argument.
305      * @return E + ω + Ω eccentric longitude argument (rad)
306      */
307     public abstract double getLE();
308 
309     /** Get the eccentric longitude argument derivative.
310      * <p>
311      * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
312      * </p>
313      * @return d(E + ω + Ω)/dt eccentric longitude argument derivative (rad/s)
314      * @since 9.0
315      */
316     public abstract double getLEDot();
317 
318     /** Get the true longitude argument.
319      * @return v + ω + Ω true longitude argument (rad)
320      */
321     public abstract double getLv();
322 
323     /** Get the true longitude argument derivative.
324      * <p>
325      * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
326      * </p>
327      * @return d(v + ω + Ω)/dt true longitude argument derivative (rad/s)
328      * @since 9.0
329      */
330     public abstract double getLvDot();
331 
332     /** Get the mean longitude argument.
333      * @return M + ω + Ω mean longitude argument (rad)
334      */
335     public abstract double getLM();
336 
337     /** Get the mean longitude argument derivative.
338      * <p>
339      * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
340      * </p>
341      * @return d(M + ω + Ω)/dt mean longitude argument derivative (rad/s)
342      * @since 9.0
343      */
344     public abstract double getLMDot();
345 
346     // Additional orbital elements
347 
348     /** Get the eccentricity.
349      * @return eccentricity
350      */
351     public abstract double getE();
352 
353     /** Get the eccentricity derivative.
354      * <p>
355      * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
356      * </p>
357      * @return eccentricity derivative
358      * @since 9.0
359      */
360     public abstract double getEDot();
361 
362     /** Get the inclination.
363      * @return inclination (rad)
364      */
365     public abstract double getI();
366 
367     /** Get the inclination derivative.
368      * <p>
369      * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
370      * </p>
371      * @return inclination derivative (rad/s)
372      * @since 9.0
373      */
374     public abstract double getIDot();
375 
376     /** Check if orbit includes non-Keplerian rates.
377      * @return true if orbit includes non-Keplerian derivatives
378      * @see #getADot()
379      * @see #getEquinoctialExDot()
380      * @see #getEquinoctialEyDot()
381      * @see #getHxDot()
382      * @see #getHyDot()
383      * @see #getLEDot()
384      * @see #getLvDot()
385      * @see #getLMDot()
386      * @see #getEDot()
387      * @see #getIDot()
388      * @since 13.0
389      */
390     public boolean hasNonKeplerianAcceleration() {
391         return hasNonKeplerianAcceleration(getPVCoordinates(), getMu());
392     }
393 
394     /** Get the central acceleration constant.
395      * @return central acceleration constant
396      */
397     public double getMu() {
398         return mu;
399     }
400 
401     /** Get the Keplerian period.
402      * <p>The Keplerian period is computed directly from semi major axis
403      * and central acceleration constant.</p>
404      * @return Keplerian period in seconds, or positive infinity for hyperbolic orbits
405      */
406     public double getKeplerianPeriod() {
407         final double a = getA();
408         return isElliptical() ? 2.0 * FastMath.PI * a * FastMath.sqrt(a / mu) : Double.POSITIVE_INFINITY;
409     }
410 
411     /** Get the Keplerian mean motion.
412      * <p>The Keplerian mean motion is computed directly from semi major axis
413      * and central acceleration constant.</p>
414      * @return Keplerian mean motion in radians per second
415      */
416     public double getKeplerianMeanMotion() {
417         final double absA = FastMath.abs(getA());
418         return FastMath.sqrt(mu / absA) / absA;
419     }
420 
421     /** Get the derivative of the mean anomaly with respect to the semi major axis.
422      * @return derivative of the mean anomaly with respect to the semi major axis
423      */
424     public double getMeanAnomalyDotWrtA() {
425         return -1.5 * getKeplerianMeanMotion() / getA();
426     }
427 
428     /** Get the date of orbital parameters.
429      * @return date of the orbital parameters
430      */
431     public AbsoluteDate getDate() {
432         return date;
433     }
434 
435     /** Get the {@link TimeStampedPVCoordinates} in a specified frame.
436      * @param outputFrame frame in which the position/velocity coordinates shall be computed
437      * @return pvCoordinates in the specified output frame
438           * @see #getPVCoordinates()
439      */
440     public TimeStampedPVCoordinates getPVCoordinates(final Frame outputFrame) {
441         if (pvCoordinates == null) {
442             pvCoordinates = initPVCoordinates();
443         }
444 
445         // If output frame requested is the same as definition frame,
446         // PV coordinates are returned directly
447         if (outputFrame == frame) {
448             return pvCoordinates;
449         }
450 
451         // Else, PV coordinates are transformed to output frame
452         final Transform t = frame.getTransformTo(outputFrame, date);
453         return t.transformPVCoordinates(pvCoordinates);
454     }
455 
456     /** {@inheritDoc} */
457     public TimeStampedPVCoordinates getPVCoordinates(final AbsoluteDate otherDate, final Frame otherFrame) {
458         return shiftedBy(otherDate.durationFrom(getDate())).getPVCoordinates(otherFrame);
459     }
460 
461     /** {@inheritDoc} */
462     @Override
463     public Vector3D getPosition(final AbsoluteDate otherDate, final Frame otherFrame) {
464         return shiftedBy(otherDate.durationFrom(getDate())).getPosition(otherFrame);
465     }
466 
467     /** Get the position in a specified frame.
468      * @param outputFrame frame in which the position coordinates shall be computed
469      * @return position in the specified output frame
470      * @see #getPosition()
471      * @since 12.0
472      */
473     public Vector3D getPosition(final Frame outputFrame) {
474         if (position == null) {
475             position = initPosition();
476         }
477 
478         // If output frame requested is the same as definition frame,
479         // Position vector is returned directly
480         if (outputFrame == frame) {
481             return position;
482         }
483 
484         // Else, position vector is transformed to output frame
485         final StaticTransform t = frame.getStaticTransformTo(outputFrame, date);
486         return t.transformPosition(position);
487 
488     }
489 
490     /** Get the position in definition frame.
491      * @return position in the definition frame
492      * @see #getPVCoordinates()
493      * @since 12.0
494      */
495     public Vector3D getPosition() {
496         if (position == null) {
497             position = initPosition();
498         }
499         return position;
500     }
501 
502     /** Get the {@link TimeStampedPVCoordinates} in definition frame.
503      * @return pvCoordinates in the definition frame
504      * @see #getPVCoordinates(Frame)
505      */
506     public TimeStampedPVCoordinates getPVCoordinates() {
507         if (pvCoordinates == null) {
508             pvCoordinates = initPVCoordinates();
509             position      = pvCoordinates.getPosition();
510         }
511         return pvCoordinates;
512     }
513 
514     /** Compute the position coordinates from the canonical parameters.
515      * @return computed position coordinates
516      * @since 12.0
517      */
518     protected abstract Vector3D initPosition();
519 
520     /** Compute the position/velocity coordinates from the canonical parameters.
521      * @return computed position/velocity coordinates
522      */
523     protected abstract TimeStampedPVCoordinates initPVCoordinates();
524 
525     /**
526      * Create a new object representing the same physical orbital state, but attached to a different reference frame.
527      * If the new frame is not inertial, an exception will be thrown.
528      *
529      * @param inertialFrame reference frame of output orbit
530      * @return orbit with different frame
531      * @since 13.0
532      */
533     public abstract Orbit withFrame(Frame inertialFrame);
534 
535     /** Get a time-shifted orbit.
536      * <p>
537      * The orbit can be slightly shifted to close dates. The shifting model is a
538      * Keplerian one if no derivatives are available in the orbit, or Keplerian
539      * plus quadratic effect of the non-Keplerian acceleration if derivatives are
540      * available. Shifting is <em>not</em> intended as a replacement for proper
541      * orbit propagation but should be sufficient for small time shifts or coarse
542      * accuracy.
543      * </p>
544      * @param dt time shift in seconds
545      * @return a new orbit, shifted with respect to the instance (which is immutable)
546      */
547     @Override
548     public abstract Orbit shiftedBy(double dt);
549 
550     /** Get a time-shifted orbit.
551      * <p>
552      * The orbit can be slightly shifted to close dates. The shifting model is a
553      * Keplerian one if no derivatives are available in the orbit, or Keplerian
554      * plus quadratic effect of the non-Keplerian acceleration if derivatives are
555      * available. Shifting is <em>not</em> intended as a replacement for proper
556      * orbit propagation but should be sufficient for small time shifts or coarse
557      * accuracy.
558      * </p>
559      * @param dt time shift
560      * @return a new orbit, shifted with respect to the instance (which is immutable)
561      */
562     @Override
563     public abstract Orbit shiftedBy(TimeOffset dt);
564 
565     /** Compute the Jacobian of the orbital parameters with respect to the Cartesian parameters.
566      * <p>
567      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
568      * respect to Cartesian coordinate j. This means each row corresponds to one orbital parameter
569      * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.
570      * </p>
571      * @param type type of the position angle to use
572      * @param jacobian placeholder 6x6 (or larger) matrix to be filled with the Jacobian, if matrix
573      * is larger than 6x6, only the 6x6 upper left corner will be modified
574      */
575     public void getJacobianWrtCartesian(final PositionAngleType type, final double[][] jacobian) {
576 
577         final double[][] cachedJacobian;
578         synchronized (this) {
579             switch (type) {
580                 case MEAN :
581                     if (jacobianMeanWrtCartesian == null) {
582                         // first call, we need to compute the Jacobian and cache it
583                         jacobianMeanWrtCartesian = computeJacobianMeanWrtCartesian();
584                     }
585                     cachedJacobian = jacobianMeanWrtCartesian;
586                     break;
587                 case ECCENTRIC :
588                     if (jacobianEccentricWrtCartesian == null) {
589                         // first call, we need to compute the Jacobian and cache it
590                         jacobianEccentricWrtCartesian = computeJacobianEccentricWrtCartesian();
591                     }
592                     cachedJacobian = jacobianEccentricWrtCartesian;
593                     break;
594                 case TRUE :
595                     if (jacobianTrueWrtCartesian == null) {
596                         // first call, we need to compute the Jacobian and cache it
597                         jacobianTrueWrtCartesian = computeJacobianTrueWrtCartesian();
598                     }
599                     cachedJacobian = jacobianTrueWrtCartesian;
600                     break;
601                 default :
602                     throw new OrekitInternalError(null);
603             }
604         }
605 
606         // fill the user provided array
607         for (int i = 0; i < cachedJacobian.length; ++i) {
608             System.arraycopy(cachedJacobian[i], 0, jacobian[i], 0, cachedJacobian[i].length);
609         }
610 
611     }
612 
613     /** Compute the Jacobian of the Cartesian parameters with respect to the orbital parameters.
614      * <p>
615      * Element {@code jacobian[i][j]} is the derivative of Cartesian coordinate i of the orbit with
616      * respect to orbital parameter j. This means each row corresponds to one Cartesian coordinate
617      * x, y, z, xdot, ydot, zdot whereas columns 0 to 5 correspond to the orbital parameters.
618      * </p>
619      * @param type type of the position angle to use
620      * @param jacobian placeholder 6x6 (or larger) matrix to be filled with the Jacobian, if matrix
621      * is larger than 6x6, only the 6x6 upper left corner will be modified
622      */
623     public void getJacobianWrtParameters(final PositionAngleType type, final double[][] jacobian) {
624 
625         final double[][] cachedJacobian;
626         synchronized (this) {
627             switch (type) {
628                 case MEAN :
629                     if (jacobianWrtParametersMean == null) {
630                         // first call, we need to compute the Jacobian and cache it
631                         jacobianWrtParametersMean = createInverseJacobian(type);
632                     }
633                     cachedJacobian = jacobianWrtParametersMean;
634                     break;
635                 case ECCENTRIC :
636                     if (jacobianWrtParametersEccentric == null) {
637                         // first call, we need to compute the Jacobian and cache it
638                         jacobianWrtParametersEccentric = createInverseJacobian(type);
639                     }
640                     cachedJacobian = jacobianWrtParametersEccentric;
641                     break;
642                 case TRUE :
643                     if (jacobianWrtParametersTrue == null) {
644                         // first call, we need to compute the Jacobian and cache it
645                         jacobianWrtParametersTrue = createInverseJacobian(type);
646                     }
647                     cachedJacobian = jacobianWrtParametersTrue;
648                     break;
649                 default :
650                     throw new OrekitInternalError(null);
651             }
652         }
653 
654         // fill the user-provided array
655         for (int i = 0; i < cachedJacobian.length; ++i) {
656             System.arraycopy(cachedJacobian[i], 0, jacobian[i], 0, cachedJacobian[i].length);
657         }
658 
659     }
660 
661     /** Create an inverse Jacobian.
662      * @param type type of the position angle to use
663      * @return inverse Jacobian
664      */
665     private double[][] createInverseJacobian(final PositionAngleType type) {
666 
667         // get the direct Jacobian
668         final double[][] directJacobian = new double[6][6];
669         getJacobianWrtCartesian(type, directJacobian);
670 
671         // invert the direct Jacobian
672         final RealMatrix matrix = MatrixUtils.createRealMatrix(directJacobian);
673         final DecompositionSolver solver = new QRDecomposition(matrix).getSolver();
674         return solver.getInverse().getData();
675 
676     }
677 
678     /** Compute the Jacobian of the orbital parameters with mean angle with respect to the Cartesian parameters.
679      * <p>
680      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
681      * respect to Cartesian coordinate j. This means each row correspond to one orbital parameter
682      * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.
683      * </p>
684      * <p>
685      * The array returned by this method will not be modified.
686      * </p>
687      * @return 6x6 Jacobian matrix
688      * @see #computeJacobianEccentricWrtCartesian()
689      * @see #computeJacobianTrueWrtCartesian()
690      */
691     protected abstract double[][] computeJacobianMeanWrtCartesian();
692 
693     /** Compute the Jacobian of the orbital parameters with eccentric angle with respect to the Cartesian parameters.
694      * <p>
695      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
696      * respect to Cartesian coordinate j. This means each row correspond to one orbital parameter
697      * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.
698      * </p>
699      * <p>
700      * The array returned by this method will not be modified.
701      * </p>
702      * @return 6x6 Jacobian matrix
703      * @see #computeJacobianMeanWrtCartesian()
704      * @see #computeJacobianTrueWrtCartesian()
705      */
706     protected abstract double[][] computeJacobianEccentricWrtCartesian();
707 
708     /** Compute the Jacobian of the orbital parameters with true angle with respect to the Cartesian parameters.
709      * <p>
710      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
711      * respect to Cartesian coordinate j. This means each row correspond to one orbital parameter
712      * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.
713      * </p>
714      * <p>
715      * The array returned by this method will not be modified.
716      * </p>
717      * @return 6x6 Jacobian matrix
718      * @see #computeJacobianMeanWrtCartesian()
719      * @see #computeJacobianEccentricWrtCartesian()
720      */
721     protected abstract double[][] computeJacobianTrueWrtCartesian();
722 
723     /** Add the contribution of the Keplerian motion to parameters derivatives
724      * <p>
725      * This method is used by integration-based propagators to evaluate the part of Keplerian
726      * motion to evolution of the orbital state.
727      * </p>
728      * @param type type of the position angle in the state
729      * @param gm attraction coefficient to use
730      * @param pDot array containing orbital state derivatives to update (the Keplerian
731      * part must be <em>added</em> to the array components, as the array may already
732      * contain some non-zero elements corresponding to non-Keplerian parts)
733      */
734     public abstract void addKeplerContribution(PositionAngleType type, double gm, double[] pDot);
735 
736         /** Fill a Jacobian half row with a single vector.
737      * @param a coefficient of the vector
738      * @param v vector
739      * @param row Jacobian matrix row
740      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
741      */
742     protected static void fillHalfRow(final double a, final Vector3D v, final double[] row, final int j) {
743         row[j]     = a * v.getX();
744         row[j + 1] = a * v.getY();
745         row[j + 2] = a * v.getZ();
746     }
747 
748     /** Fill a Jacobian half row with a linear combination of vectors.
749      * @param a1 coefficient of the first vector
750      * @param v1 first vector
751      * @param a2 coefficient of the second vector
752      * @param v2 second vector
753      * @param row Jacobian matrix row
754      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
755      */
756     protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2,
757                                       final double[] row, final int j) {
758         row[j]     = MathArrays.linearCombination(a1, v1.getX(), a2, v2.getX());
759         row[j + 1] = MathArrays.linearCombination(a1, v1.getY(), a2, v2.getY());
760         row[j + 2] = MathArrays.linearCombination(a1, v1.getZ(), a2, v2.getZ());
761     }
762 
763     /** Fill a Jacobian half row with a linear combination of vectors.
764      * @param a1 coefficient of the first vector
765      * @param v1 first vector
766      * @param a2 coefficient of the second vector
767      * @param v2 second vector
768      * @param a3 coefficient of the third vector
769      * @param v3 third vector
770      * @param row Jacobian matrix row
771      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
772      */
773     protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2,
774                                       final double a3, final Vector3D v3,
775                                       final double[] row, final int j) {
776         row[j]     = MathArrays.linearCombination(a1, v1.getX(), a2, v2.getX(), a3, v3.getX());
777         row[j + 1] = MathArrays.linearCombination(a1, v1.getY(), a2, v2.getY(), a3, v3.getY());
778         row[j + 2] = MathArrays.linearCombination(a1, v1.getZ(), a2, v2.getZ(), a3, v3.getZ());
779     }
780 
781     /** Fill a Jacobian half row with a linear combination of vectors.
782      * @param a1 coefficient of the first vector
783      * @param v1 first vector
784      * @param a2 coefficient of the second vector
785      * @param v2 second vector
786      * @param a3 coefficient of the third vector
787      * @param v3 third vector
788      * @param a4 coefficient of the fourth vector
789      * @param v4 fourth vector
790      * @param row Jacobian matrix row
791      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
792      */
793     protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2,
794                                       final double a3, final Vector3D v3, final double a4, final Vector3D v4,
795                                       final double[] row, final int j) {
796         row[j]     = MathArrays.linearCombination(a1, v1.getX(), a2, v2.getX(), a3, v3.getX(), a4, v4.getX());
797         row[j + 1] = MathArrays.linearCombination(a1, v1.getY(), a2, v2.getY(), a3, v3.getY(), a4, v4.getY());
798         row[j + 2] = MathArrays.linearCombination(a1, v1.getZ(), a2, v2.getZ(), a3, v3.getZ(), a4, v4.getZ());
799     }
800 
801     /** Fill a Jacobian half row with a linear combination of vectors.
802      * @param a1 coefficient of the first vector
803      * @param v1 first vector
804      * @param a2 coefficient of the second vector
805      * @param v2 second vector
806      * @param a3 coefficient of the third vector
807      * @param v3 third vector
808      * @param a4 coefficient of the fourth vector
809      * @param v4 fourth vector
810      * @param a5 coefficient of the fifth vector
811      * @param v5 fifth vector
812      * @param row Jacobian matrix row
813      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
814      */
815     protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2,
816                                       final double a3, final Vector3D v3, final double a4, final Vector3D v4,
817                                       final double a5, final Vector3D v5,
818                                       final double[] row, final int j) {
819         final double[] a = new double[] {
820             a1, a2, a3, a4, a5
821         };
822         row[j]     = MathArrays.linearCombination(a, new double[] {
823             v1.getX(), v2.getX(), v3.getX(), v4.getX(), v5.getX()
824         });
825         row[j + 1] = MathArrays.linearCombination(a, new double[] {
826             v1.getY(), v2.getY(), v3.getY(), v4.getY(), v5.getY()
827         });
828         row[j + 2] = MathArrays.linearCombination(a, new double[] {
829             v1.getZ(), v2.getZ(), v3.getZ(), v4.getZ(), v5.getZ()
830         });
831     }
832 
833     /** Fill a Jacobian half row with a linear combination of vectors.
834      * @param a1 coefficient of the first vector
835      * @param v1 first vector
836      * @param a2 coefficient of the second vector
837      * @param v2 second vector
838      * @param a3 coefficient of the third vector
839      * @param v3 third vector
840      * @param a4 coefficient of the fourth vector
841      * @param v4 fourth vector
842      * @param a5 coefficient of the fifth vector
843      * @param v5 fifth vector
844      * @param a6 coefficient of the sixth vector
845      * @param v6 sixth vector
846      * @param row Jacobian matrix row
847      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
848      */
849     protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2,
850                                       final double a3, final Vector3D v3, final double a4, final Vector3D v4,
851                                       final double a5, final Vector3D v5, final double a6, final Vector3D v6,
852                                       final double[] row, final int j) {
853         final double[] a = new double[] {
854             a1, a2, a3, a4, a5, a6
855         };
856         row[j]     = MathArrays.linearCombination(a, new double[] {
857             v1.getX(), v2.getX(), v3.getX(), v4.getX(), v5.getX(), v6.getX()
858         });
859         row[j + 1] = MathArrays.linearCombination(a, new double[] {
860             v1.getY(), v2.getY(), v3.getY(), v4.getY(), v5.getY(), v6.getY()
861         });
862         row[j + 2] = MathArrays.linearCombination(a, new double[] {
863             v1.getZ(), v2.getZ(), v3.getZ(), v4.getZ(), v5.getZ(), v6.getZ()
864         });
865     }
866 
867 }