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