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