1   /* Copyright 2002-2023 CS GROUP
2    * Licensed to CS GROUP (CS) under one or more
3    * contributor license agreements.  See the NOTICE file distributed with
4    * this work for additional information regarding copyright ownership.
5    * CS licenses this file to You under the Apache License, Version 2.0
6    * (the "License"); you may not use this file except in compliance with
7    * the License.  You may obtain a copy of the License at
8    *
9    *   http://www.apache.org/licenses/LICENSE-2.0
10   *
11   * Unless required by applicable law or agreed to in writing, software
12   * distributed under the License is distributed on an "AS IS" BASIS,
13   * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14   * See the License for the specific language governing permissions and
15   * limitations under the License.
16   */
17  package org.orekit.orbits;
18  
19  import java.io.Serializable;
20  
21  import org.hipparchus.geometry.euclidean.threed.Vector3D;
22  import org.hipparchus.linear.DecompositionSolver;
23  import org.hipparchus.linear.MatrixUtils;
24  import org.hipparchus.linear.QRDecomposition;
25  import org.hipparchus.linear.RealMatrix;
26  import org.hipparchus.util.FastMath;
27  import org.hipparchus.util.MathArrays;
28  import org.orekit.errors.OrekitIllegalArgumentException;
29  import org.orekit.errors.OrekitInternalError;
30  import org.orekit.errors.OrekitMessages;
31  import org.orekit.frames.Frame;
32  import org.orekit.frames.StaticTransform;
33  import org.orekit.frames.Transform;
34  import org.orekit.time.AbsoluteDate;
35  import org.orekit.time.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>, Serializable, PVCoordinatesProvider {
66  
67      /** Serializable UID. */
68      private static final long serialVersionUID = 438733454597999578L;
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 transient Vector3D position;
83  
84      /** Computed PVCoordinates. */
85      private transient TimeStampedPVCoordinates pvCoordinates;
86  
87      /** Jacobian of the orbital parameters with mean angle with respect to the Cartesian coordinates. */
88      private transient double[][] jacobianMeanWrtCartesian;
89  
90      /** Jacobian of the Cartesian coordinates with respect to the orbital parameters with mean angle. */
91      private transient double[][] jacobianWrtParametersMean;
92  
93      /** Jacobian of the orbital parameters with eccentric angle with respect to the Cartesian coordinates. */
94      private transient double[][] jacobianEccentricWrtCartesian;
95  
96      /** Jacobian of the Cartesian coordinates with respect to the orbital parameters with eccentric angle. */
97      private transient double[][] jacobianWrtParametersEccentric;
98  
99      /** Jacobian of the orbital parameters with true angle with respect to the Cartesian coordinates. */
100     private transient double[][] jacobianTrueWrtCartesian;
101 
102     /** Jacobian of the Cartesian coordinates with respect to the orbital parameters with true angle. */
103     private transient 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     /** Check if Cartesian coordinates include non-Keplerian acceleration.
164      * @param pva Cartesian coordinates
165      * @param mu central attraction coefficient
166      * @return true if Cartesian coordinates include non-Keplerian acceleration
167      */
168     protected static boolean hasNonKeplerianAcceleration(final PVCoordinates pva, final double mu) {
169 
170         final Vector3D p = pva.getPosition();
171         final double r2 = p.getNormSq();
172         final double r  = FastMath.sqrt(r2);
173         final Vector3D keplerianAcceleration = new Vector3D(-mu / (r * r2), p);
174 
175         // Check if acceleration is null or relatively close to 0 compared to the keplerain acceleration
176         final Vector3D a = pva.getAcceleration();
177         if (a == null || a.getNorm() < 1e-9 * keplerianAcceleration.getNorm()) {
178             return false;
179         }
180 
181         final Vector3D nonKeplerianAcceleration = a.subtract(keplerianAcceleration);
182 
183         if ( nonKeplerianAcceleration.getNorm() > 1e-9 * keplerianAcceleration.getNorm()) {
184             // we have a relevant acceleration, we can compute derivatives
185             return true;
186         } else {
187             // the provided acceleration is either too small to be reliable (probably even 0), or NaN
188             return false;
189         }
190     }
191 
192     /** Returns true if and only if the orbit is elliptical i.e. has a non-negative semi-major axis.
193      * @return true if getA() is strictly greater than 0
194      * @since 12.0
195      */
196     public boolean isElliptical() {
197         return getA() > 0.;
198     }
199 
200     /** Get the orbit type.
201      * @return orbit type
202      */
203     public abstract OrbitType getType();
204 
205     /** Ensure the defining frame is a pseudo-inertial frame.
206      * @param frame frame to check
207      * @exception IllegalArgumentException if frame is not a {@link
208      * Frame#isPseudoInertial pseudo-inertial frame}
209      */
210     private static void ensurePseudoInertialFrame(final Frame frame)
211         throws IllegalArgumentException {
212         if (!frame.isPseudoInertial()) {
213             throw new OrekitIllegalArgumentException(OrekitMessages.NON_PSEUDO_INERTIAL_FRAME,
214                                                      frame.getName());
215         }
216     }
217 
218     /** Get the frame in which the orbital parameters are defined.
219      * @return frame in which the orbital parameters are defined
220      */
221     public Frame getFrame() {
222         return frame;
223     }
224 
225     /** Get the semi-major axis.
226      * <p>Note that the semi-major axis is considered negative for hyperbolic orbits.</p>
227      * @return semi-major axis (m)
228      */
229     public abstract double getA();
230 
231     /** Get the semi-major axis derivative.
232      * <p>Note that the semi-major axis is considered negative for hyperbolic orbits.</p>
233      * <p>
234      * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
235      * </p>
236      * @return semi-major axis  derivative (m/s)
237      * @see #hasDerivatives()
238      * @since 9.0
239      */
240     public abstract double getADot();
241 
242     /** Get the first component of the equinoctial eccentricity vector derivative.
243      * @return first component of the equinoctial eccentricity vector derivative
244      */
245     public abstract double getEquinoctialEx();
246 
247     /** Get the first component of the equinoctial eccentricity vector.
248      * <p>
249      * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
250      * </p>
251      * @return first component of the equinoctial eccentricity vector
252      * @see #hasDerivatives()
253      * @since 9.0
254      */
255     public abstract double getEquinoctialExDot();
256 
257     /** Get the second component of the equinoctial eccentricity vector derivative.
258      * @return second component of the equinoctial eccentricity vector derivative
259      */
260     public abstract double getEquinoctialEy();
261 
262     /** Get the second component of the equinoctial eccentricity vector.
263      * <p>
264      * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
265      * </p>
266      * @return second component of the equinoctial eccentricity vector
267      * @see #hasDerivatives()
268      * @since 9.0
269      */
270     public abstract double getEquinoctialEyDot();
271 
272     /** Get the first component of the inclination vector.
273      * @return first component of the inclination vector
274      */
275     public abstract double getHx();
276 
277     /** Get the first component of the inclination vector derivative.
278      * <p>
279      * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
280      * </p>
281      * @return first component of the inclination vector derivative
282      * @see #hasDerivatives()
283      * @since 9.0
284      */
285     public abstract double getHxDot();
286 
287     /** Get the second component of the inclination vector.
288      * @return second component of the inclination vector
289      */
290     public abstract double getHy();
291 
292     /** Get the second component of the inclination vector derivative.
293      * <p>
294      * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
295      * </p>
296      * @return second component of the inclination vector derivative
297      * @see #hasDerivatives()
298      * @since 9.0
299      */
300     public abstract double getHyDot();
301 
302     /** Get the eccentric longitude argument.
303      * @return E + ω + Ω eccentric longitude argument (rad)
304      */
305     public abstract double getLE();
306 
307     /** Get the eccentric longitude argument derivative.
308      * <p>
309      * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
310      * </p>
311      * @return d(E + ω + Ω)/dt eccentric longitude argument derivative (rad/s)
312      * @see #hasDerivatives()
313      * @since 9.0
314      */
315     public abstract double getLEDot();
316 
317     /** Get the true longitude argument.
318      * @return v + ω + Ω true longitude argument (rad)
319      */
320     public abstract double getLv();
321 
322     /** Get the true longitude argument derivative.
323      * <p>
324      * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
325      * </p>
326      * @return d(v + ω + Ω)/dt true longitude argument derivative (rad/s)
327      * @see #hasDerivatives()
328      * @since 9.0
329      */
330     public abstract double getLvDot();
331 
332     /** Get the mean longitude argument.
333      * @return M + ω + Ω mean longitude argument (rad)
334      */
335     public abstract double getLM();
336 
337     /** Get the mean longitude argument derivative.
338      * <p>
339      * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
340      * </p>
341      * @return d(M + ω + Ω)/dt mean longitude argument derivative (rad/s)
342      * @see #hasDerivatives()
343      * @since 9.0
344      */
345     public abstract double getLMDot();
346 
347     // Additional orbital elements
348 
349     /** Get the eccentricity.
350      * @return eccentricity
351      */
352     public abstract double getE();
353 
354     /** Get the eccentricity derivative.
355      * <p>
356      * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
357      * </p>
358      * @return eccentricity derivative
359      * @see #hasDerivatives()
360      * @since 9.0
361      */
362     public abstract double getEDot();
363 
364     /** Get the inclination.
365      * @return inclination (rad)
366      */
367     public abstract double getI();
368 
369     /** Get the inclination derivative.
370      * <p>
371      * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
372      * </p>
373      * @return inclination derivative (rad/s)
374      * @see #hasDerivatives()
375      * @since 9.0
376      */
377     public abstract double getIDot();
378 
379     /** Check if orbit includes derivatives.
380      * @return true if orbit includes derivatives
381      * @see #getADot()
382      * @see #getEquinoctialExDot()
383      * @see #getEquinoctialEyDot()
384      * @see #getHxDot()
385      * @see #getHyDot()
386      * @see #getLEDot()
387      * @see #getLvDot()
388      * @see #getLMDot()
389      * @see #getEDot()
390      * @see #getIDot()
391      * @since 9.0
392      */
393     public boolean hasDerivatives() {
394         return !Double.isNaN(getADot());
395     }
396 
397     /** Get the central acceleration constant.
398      * @return central acceleration constant
399      */
400     public double getMu() {
401         return mu;
402     }
403 
404     /** Get the Keplerian period.
405      * <p>The Keplerian period is computed directly from semi major axis
406      * and central acceleration constant.</p>
407      * @return Keplerian period in seconds, or positive infinity for hyperbolic orbits
408      */
409     public double getKeplerianPeriod() {
410         final double a = getA();
411         return isElliptical() ? 2.0 * FastMath.PI * a * FastMath.sqrt(a / mu) : Double.POSITIVE_INFINITY;
412     }
413 
414     /** Get the Keplerian mean motion.
415      * <p>The Keplerian mean motion is computed directly from semi major axis
416      * and central acceleration constant.</p>
417      * @return Keplerian mean motion in radians per second
418      */
419     public double getKeplerianMeanMotion() {
420         final double absA = FastMath.abs(getA());
421         return FastMath.sqrt(mu / absA) / absA;
422     }
423 
424     /** Get the derivative of the mean anomaly with respect to the semi major axis.
425      * @return derivative of the mean anomaly with respect to the semi major axis
426      */
427     public double getMeanAnomalyDotWrtA() {
428         return -1.5 * getKeplerianMeanMotion() / getA();
429     }
430 
431     /** Get the date of orbital parameters.
432      * @return date of the orbital parameters
433      */
434     public AbsoluteDate getDate() {
435         return date;
436     }
437 
438     /** Get the {@link TimeStampedPVCoordinates} in a specified frame.
439      * @param outputFrame frame in which the position/velocity coordinates shall be computed
440      * @return pvCoordinates in the specified output frame
441           * @see #getPVCoordinates()
442      */
443     public TimeStampedPVCoordinates getPVCoordinates(final Frame outputFrame) {
444         if (pvCoordinates == null) {
445             pvCoordinates = initPVCoordinates();
446         }
447 
448         // If output frame requested is the same as definition frame,
449         // PV coordinates are returned directly
450         if (outputFrame == frame) {
451             return pvCoordinates;
452         }
453 
454         // Else, PV coordinates are transformed to output frame
455         final Transform t = frame.getTransformTo(outputFrame, date);
456         return t.transformPVCoordinates(pvCoordinates);
457     }
458 
459     /** {@inheritDoc} */
460     public TimeStampedPVCoordinates getPVCoordinates(final AbsoluteDate otherDate, final Frame otherFrame) {
461         return shiftedBy(otherDate.durationFrom(getDate())).getPVCoordinates(otherFrame);
462     }
463 
464     /** Get the position in a specified frame.
465      * @param outputFrame frame in which the position coordinates shall be computed
466      * @return position in the specified output frame
467      * @see #getPosition()
468      * @since 12.0
469      */
470     public Vector3D getPosition(final Frame outputFrame) {
471         if (position == null) {
472             position = initPosition();
473         }
474 
475         // If output frame requested is the same as definition frame,
476         // Position vector is returned directly
477         if (outputFrame == frame) {
478             return position;
479         }
480 
481         // Else, position vector is transformed to output frame
482         final StaticTransform t = frame.getStaticTransformTo(outputFrame, date);
483         return t.transformPosition(position);
484 
485     }
486 
487     /** Get the position in definition frame.
488      * @return position in the definition frame
489      * @see #getPVCoordinates()
490      * @since 12.0
491      */
492     public Vector3D getPosition() {
493         if (position == null) {
494             position = initPosition();
495         }
496         return position;
497     }
498 
499     /** Get the {@link TimeStampedPVCoordinates} in definition frame.
500      * @return pvCoordinates in the definition frame
501      * @see #getPVCoordinates(Frame)
502      */
503     public TimeStampedPVCoordinates getPVCoordinates() {
504         if (pvCoordinates == null) {
505             pvCoordinates = initPVCoordinates();
506             position      = pvCoordinates.getPosition();
507         }
508         return pvCoordinates;
509     }
510 
511     /** Compute the position coordinates from the canonical parameters.
512      * @return computed position coordinates
513      * @since 12.0
514      */
515     protected abstract Vector3D initPosition();
516 
517     /** Compute the position/velocity coordinates from the canonical parameters.
518      * @return computed position/velocity coordinates
519      */
520     protected abstract TimeStampedPVCoordinates initPVCoordinates();
521 
522     /** Get a time-shifted orbit.
523      * <p>
524      * The orbit can be slightly shifted to close dates. The shifting model is a
525      * Keplerian one if no derivatives are available in the orbit, or Keplerian
526      * plus quadratic effect of the non-Keplerian acceleration if derivatives are
527      * available. Shifting is <em>not</em> intended as a replacement for proper
528      * orbit propagation but should be sufficient for small time shifts or coarse
529      * accuracy.
530      * </p>
531      * @param dt time shift in seconds
532      * @return a new orbit, shifted with respect to the instance (which is immutable)
533      */
534     public abstract Orbit shiftedBy(double dt);
535 
536     /** Compute the Jacobian of the orbital parameters with respect to the Cartesian parameters.
537      * <p>
538      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
539      * respect to Cartesian coordinate j. This means each row corresponds to one orbital parameter
540      * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.
541      * </p>
542      * @param type type of the position angle to use
543      * @param jacobian placeholder 6x6 (or larger) matrix to be filled with the Jacobian, if matrix
544      * is larger than 6x6, only the 6x6 upper left corner will be modified
545      */
546     public void getJacobianWrtCartesian(final PositionAngleType type, final double[][] jacobian) {
547 
548         final double[][] cachedJacobian;
549         synchronized (this) {
550             switch (type) {
551                 case MEAN :
552                     if (jacobianMeanWrtCartesian == null) {
553                         // first call, we need to compute the Jacobian and cache it
554                         jacobianMeanWrtCartesian = computeJacobianMeanWrtCartesian();
555                     }
556                     cachedJacobian = jacobianMeanWrtCartesian;
557                     break;
558                 case ECCENTRIC :
559                     if (jacobianEccentricWrtCartesian == null) {
560                         // first call, we need to compute the Jacobian and cache it
561                         jacobianEccentricWrtCartesian = computeJacobianEccentricWrtCartesian();
562                     }
563                     cachedJacobian = jacobianEccentricWrtCartesian;
564                     break;
565                 case TRUE :
566                     if (jacobianTrueWrtCartesian == null) {
567                         // first call, we need to compute the Jacobian and cache it
568                         jacobianTrueWrtCartesian = computeJacobianTrueWrtCartesian();
569                     }
570                     cachedJacobian = jacobianTrueWrtCartesian;
571                     break;
572                 default :
573                     throw new OrekitInternalError(null);
574             }
575         }
576 
577         // fill the user provided array
578         for (int i = 0; i < cachedJacobian.length; ++i) {
579             System.arraycopy(cachedJacobian[i], 0, jacobian[i], 0, cachedJacobian[i].length);
580         }
581 
582     }
583 
584     /** Compute the Jacobian of the Cartesian parameters with respect to the orbital parameters.
585      * <p>
586      * Element {@code jacobian[i][j]} is the derivative of Cartesian coordinate i of the orbit with
587      * respect to orbital parameter j. This means each row corresponds to one Cartesian coordinate
588      * x, y, z, xdot, ydot, zdot whereas columns 0 to 5 correspond to the orbital parameters.
589      * </p>
590      * @param type type of the position angle to use
591      * @param jacobian placeholder 6x6 (or larger) matrix to be filled with the Jacobian, if matrix
592      * is larger than 6x6, only the 6x6 upper left corner will be modified
593      */
594     public void getJacobianWrtParameters(final PositionAngleType type, final double[][] jacobian) {
595 
596         final double[][] cachedJacobian;
597         synchronized (this) {
598             switch (type) {
599                 case MEAN :
600                     if (jacobianWrtParametersMean == null) {
601                         // first call, we need to compute the Jacobian and cache it
602                         jacobianWrtParametersMean = createInverseJacobian(type);
603                     }
604                     cachedJacobian = jacobianWrtParametersMean;
605                     break;
606                 case ECCENTRIC :
607                     if (jacobianWrtParametersEccentric == null) {
608                         // first call, we need to compute the Jacobian and cache it
609                         jacobianWrtParametersEccentric = createInverseJacobian(type);
610                     }
611                     cachedJacobian = jacobianWrtParametersEccentric;
612                     break;
613                 case TRUE :
614                     if (jacobianWrtParametersTrue == null) {
615                         // first call, we need to compute the Jacobian and cache it
616                         jacobianWrtParametersTrue = createInverseJacobian(type);
617                     }
618                     cachedJacobian = jacobianWrtParametersTrue;
619                     break;
620                 default :
621                     throw new OrekitInternalError(null);
622             }
623         }
624 
625         // fill the user-provided array
626         for (int i = 0; i < cachedJacobian.length; ++i) {
627             System.arraycopy(cachedJacobian[i], 0, jacobian[i], 0, cachedJacobian[i].length);
628         }
629 
630     }
631 
632     /** Create an inverse Jacobian.
633      * @param type type of the position angle to use
634      * @return inverse Jacobian
635      */
636     private double[][] createInverseJacobian(final PositionAngleType type) {
637 
638         // get the direct Jacobian
639         final double[][] directJacobian = new double[6][6];
640         getJacobianWrtCartesian(type, directJacobian);
641 
642         // invert the direct Jacobian
643         final RealMatrix matrix = MatrixUtils.createRealMatrix(directJacobian);
644         final DecompositionSolver solver = new QRDecomposition(matrix).getSolver();
645         return solver.getInverse().getData();
646 
647     }
648 
649     /** Compute the Jacobian of the orbital parameters with mean angle with respect to the Cartesian parameters.
650      * <p>
651      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
652      * respect to Cartesian coordinate j. This means each row correspond to one orbital parameter
653      * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.
654      * </p>
655      * <p>
656      * The array returned by this method will not be modified.
657      * </p>
658      * @return 6x6 Jacobian matrix
659      * @see #computeJacobianEccentricWrtCartesian()
660      * @see #computeJacobianTrueWrtCartesian()
661      */
662     protected abstract double[][] computeJacobianMeanWrtCartesian();
663 
664     /** Compute the Jacobian of the orbital parameters with eccentric 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 #computeJacobianMeanWrtCartesian()
675      * @see #computeJacobianTrueWrtCartesian()
676      */
677     protected abstract double[][] computeJacobianEccentricWrtCartesian();
678 
679     /** Compute the Jacobian of the orbital parameters with true 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 #computeJacobianEccentricWrtCartesian()
691      */
692     protected abstract double[][] computeJacobianTrueWrtCartesian();
693 
694     /** Add the contribution of the Keplerian motion to parameters derivatives
695      * <p>
696      * This method is used by integration-based propagators to evaluate the part of Keplerian
697      * motion to evolution of the orbital state.
698      * </p>
699      * @param type type of the position angle in the state
700      * @param gm attraction coefficient to use
701      * @param pDot array containing orbital state derivatives to update (the Keplerian
702      * part must be <em>added</em> to the array components, as the array may already
703      * contain some non-zero elements corresponding to non-Keplerian parts)
704      */
705     public abstract void addKeplerContribution(PositionAngleType type, double gm, double[] pDot);
706 
707         /** Fill a Jacobian half row with a single vector.
708      * @param a coefficient of the vector
709      * @param v vector
710      * @param row Jacobian matrix row
711      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
712      */
713     protected static void fillHalfRow(final double a, final Vector3D v, final double[] row, final int j) {
714         row[j]     = a * v.getX();
715         row[j + 1] = a * v.getY();
716         row[j + 2] = a * v.getZ();
717     }
718 
719     /** Fill a Jacobian half row with a linear combination of vectors.
720      * @param a1 coefficient of the first vector
721      * @param v1 first vector
722      * @param a2 coefficient of the second vector
723      * @param v2 second vector
724      * @param row Jacobian matrix row
725      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
726      */
727     protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2,
728                                       final double[] row, final int j) {
729         row[j]     = MathArrays.linearCombination(a1, v1.getX(), a2, v2.getX());
730         row[j + 1] = MathArrays.linearCombination(a1, v1.getY(), a2, v2.getY());
731         row[j + 2] = MathArrays.linearCombination(a1, v1.getZ(), a2, v2.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 a3 coefficient of the third vector
740      * @param v3 third vector
741      * @param row Jacobian matrix row
742      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
743      */
744     protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2,
745                                       final double a3, final Vector3D v3,
746                                       final double[] row, final int j) {
747         row[j]     = MathArrays.linearCombination(a1, v1.getX(), a2, v2.getX(), a3, v3.getX());
748         row[j + 1] = MathArrays.linearCombination(a1, v1.getY(), a2, v2.getY(), a3, v3.getY());
749         row[j + 2] = MathArrays.linearCombination(a1, v1.getZ(), a2, v2.getZ(), a3, v3.getZ());
750     }
751 
752     /** Fill a Jacobian half row with a linear combination of vectors.
753      * @param a1 coefficient of the first vector
754      * @param v1 first vector
755      * @param a2 coefficient of the second vector
756      * @param v2 second vector
757      * @param a3 coefficient of the third vector
758      * @param v3 third vector
759      * @param a4 coefficient of the fourth vector
760      * @param v4 fourth 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, final double a4, final Vector3D v4,
766                                       final double[] row, final int j) {
767         row[j]     = MathArrays.linearCombination(a1, v1.getX(), a2, v2.getX(), a3, v3.getX(), a4, v4.getX());
768         row[j + 1] = MathArrays.linearCombination(a1, v1.getY(), a2, v2.getY(), a3, v3.getY(), a4, v4.getY());
769         row[j + 2] = MathArrays.linearCombination(a1, v1.getZ(), a2, v2.getZ(), a3, v3.getZ(), a4, v4.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 a5 coefficient of the fifth vector
782      * @param v5 fifth vector
783      * @param row Jacobian matrix row
784      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
785      */
786     protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2,
787                                       final double a3, final Vector3D v3, final double a4, final Vector3D v4,
788                                       final double a5, final Vector3D v5,
789                                       final double[] row, final int j) {
790         final double[] a = new double[] {
791             a1, a2, a3, a4, a5
792         };
793         row[j]     = MathArrays.linearCombination(a, new double[] {
794             v1.getX(), v2.getX(), v3.getX(), v4.getX(), v5.getX()
795         });
796         row[j + 1] = MathArrays.linearCombination(a, new double[] {
797             v1.getY(), v2.getY(), v3.getY(), v4.getY(), v5.getY()
798         });
799         row[j + 2] = MathArrays.linearCombination(a, new double[] {
800             v1.getZ(), v2.getZ(), v3.getZ(), v4.getZ(), v5.getZ()
801         });
802     }
803 
804     /** Fill a Jacobian half row with a linear combination of vectors.
805      * @param a1 coefficient of the first vector
806      * @param v1 first vector
807      * @param a2 coefficient of the second vector
808      * @param v2 second vector
809      * @param a3 coefficient of the third vector
810      * @param v3 third vector
811      * @param a4 coefficient of the fourth vector
812      * @param v4 fourth vector
813      * @param a5 coefficient of the fifth vector
814      * @param v5 fifth vector
815      * @param a6 coefficient of the sixth vector
816      * @param v6 sixth vector
817      * @param row Jacobian matrix row
818      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
819      */
820     protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2,
821                                       final double a3, final Vector3D v3, final double a4, final Vector3D v4,
822                                       final double a5, final Vector3D v5, final double a6, final Vector3D v6,
823                                       final double[] row, final int j) {
824         final double[] a = new double[] {
825             a1, a2, a3, a4, a5, a6
826         };
827         row[j]     = MathArrays.linearCombination(a, new double[] {
828             v1.getX(), v2.getX(), v3.getX(), v4.getX(), v5.getX(), v6.getX()
829         });
830         row[j + 1] = MathArrays.linearCombination(a, new double[] {
831             v1.getY(), v2.getY(), v3.getY(), v4.getY(), v5.getY(), v6.getY()
832         });
833         row[j + 2] = MathArrays.linearCombination(a, new double[] {
834             v1.getZ(), v2.getZ(), v3.getZ(), v4.getZ(), v5.getZ(), v6.getZ()
835         });
836     }
837 
838 }