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