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