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