1   /* Copyright 2002-2026 CS GROUP
2    * Licensed to CS GROUP (CS) under one or more
3    * contributor license agreements.  See the NOTICE file distributed with
4    * this work for additional information regarding copyright ownership.
5    * CS licenses this file to You under the Apache License, Version 2.0
6    * (the "License"); you may not use this file except in compliance with
7    * the License.  You may obtain a copy of the License at
8    *
9    *   http://www.apache.org/licenses/LICENSE-2.0
10   *
11   * Unless required by applicable law or agreed to in writing, software
12   * distributed under the License is distributed on an "AS IS" BASIS,
13   * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14   * See the License for the specific language governing permissions and
15   * limitations under the License.
16   */
17  package org.orekit.orbits;
18  
19  import org.hipparchus.analysis.differentiation.UnivariateDerivative1;
20  import org.hipparchus.geometry.euclidean.threed.Vector3D;
21  import org.hipparchus.util.FastMath;
22  import org.hipparchus.util.SinCos;
23  import org.orekit.errors.OrekitException;
24  import org.orekit.errors.OrekitIllegalArgumentException;
25  import org.orekit.errors.OrekitInternalError;
26  import org.orekit.errors.OrekitMessages;
27  import org.orekit.frames.Frame;
28  import org.orekit.frames.KinematicTransform;
29  import org.orekit.time.AbsoluteDate;
30  import org.orekit.time.TimeOffset;
31  import org.orekit.utils.PVCoordinates;
32  import org.orekit.utils.TimeStampedPVCoordinates;
33  
34  
35  /**
36   * This class handles traditional Keplerian orbital parameters.
37  
38   * <p>
39   * The parameters used internally are the classical Keplerian elements:
40   *   <pre>
41   *     a
42   *     e
43   *     i
44   *     ω
45   *     Ω
46   *     v
47   *   </pre>
48   * where ω stands for the Perigee Argument, Ω stands for the
49   * Right Ascension of the Ascending Node and v stands for the true anomaly.
50   *
51   * <p>
52   * This class supports hyperbolic orbits, using the convention that semi major
53   * axis is negative for such orbits (and of course eccentricity is greater than 1).
54   * </p>
55   * <p>
56   * When orbit is either equatorial or circular, some Keplerian elements
57   * (more precisely ω and Ω) become ambiguous so this class should not
58   * be used for such orbits. For this reason, {@link EquinoctialOrbit equinoctial
59   * orbits} is the recommended way to represent orbits.
60   * </p>
61  
62   * <p>
63   * The instance <code>KeplerianOrbit</code> is guaranteed to be immutable.
64   * </p>
65   * @see     Orbit
66   * @see    CircularOrbit
67   * @see    CartesianOrbit
68   * @see    EquinoctialOrbit
69   * @author Luc Maisonobe
70   * @author Guylaine Prat
71   * @author Fabien Maussion
72   * @author V&eacute;ronique Pommier-Maurussane
73   */
74  public class KeplerianOrbit extends Orbit implements PositionAngleBased<KeplerianOrbit> {
75  
76      /** Name of the eccentricity parameter. */
77      private static final String ECCENTRICITY = "eccentricity";
78  
79      /** Semi-major axis (m). */
80      private final double a;
81  
82      /** Eccentricity. */
83      private final double e;
84  
85      /** Inclination (rad). */
86      private final double i;
87  
88      /** Perigee Argument (rad). */
89      private final double pa;
90  
91      /** Right Ascension of Ascending Node (rad). */
92      private final double raan;
93  
94      /** Cached anomaly (rad). */
95      private final double cachedAnomaly;
96  
97      /** Semi-major axis derivative (m/s). */
98      private final double aDot;
99  
100     /** Eccentricity derivative. */
101     private final double eDot;
102 
103     /** Inclination derivative (rad/s). */
104     private final double iDot;
105 
106     /** Perigee Argument derivative (rad/s). */
107     private final double paDot;
108 
109     /** Right Ascension of Ascending Node derivative (rad/s). */
110     private final double raanDot;
111 
112     /** Derivative of cached anomaly (rad/s). */
113     private final double cachedAnomalyDot;
114 
115     /** Cached type of position angle. */
116     private final PositionAngleType cachedPositionAngleType;
117 
118     /** Partial Cartesian coordinates (position and velocity are valid, acceleration may be missing). */
119     private PVCoordinates partialPV;
120 
121     /** Creates a new instance.
122      * @param keplerianParameters Keplerian elements
123      * @param frame the frame in which the parameters are defined
124      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
125      * @param date date of the orbital parameters
126      * @param mu central attraction coefficient (m³/s²)
127      * @exception IllegalArgumentException if frame is not a {@link
128      * Frame#isPseudoInertial pseudo-inertial frame} or a and e don't match for hyperbolic orbits,
129      * or v is out of range for hyperbolic orbits
130      * @since 14.0
131      */
132     public KeplerianOrbit(final KeplerianParameters keplerianParameters, final Frame frame, final AbsoluteDate date,
133                           final double mu)
134             throws IllegalArgumentException {
135         this(keplerianParameters, 0., 0., 0., 0., 0.,
136                 computeKeplerianAnomalyDot(keplerianParameters.positionAngleType(), keplerianParameters.a(), keplerianParameters.e(), mu, keplerianParameters.anomaly(), keplerianParameters.positionAngleType()),
137                 keplerianParameters.positionAngleType(), frame, date, mu);
138     }
139 
140     /** Creates a new instance without derivatives and with cached position angle same as value inputted.
141      * @param a  semi-major axis (m), negative for hyperbolic orbits
142      * @param e eccentricity (positive or equal to 0)
143      * @param i inclination (rad)
144      * @param pa perigee argument (ω, rad)
145      * @param raan right ascension of ascending node (Ω, rad)
146      * @param anomaly mean, eccentric or true anomaly (rad)
147      * @param type type of anomaly
148      * @param frame the frame in which the parameters are defined
149      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
150      * @param date date of the orbital parameters
151      * @param mu central attraction coefficient (m³/s²)
152      * @exception IllegalArgumentException if frame is not a {@link
153      * Frame#isPseudoInertial pseudo-inertial frame} or a and e don't match for hyperbolic orbits,
154      * or v is out of range for hyperbolic orbits
155      */
156     public KeplerianOrbit(final double a, final double e, final double i, final double pa, final double raan,
157                           final double anomaly, final PositionAngleType type, final Frame frame,
158                           final AbsoluteDate date, final double mu)
159             throws IllegalArgumentException {
160         this(new KeplerianParameters(a, e, i, pa, raan, anomaly, type), frame, date, mu);
161     }
162 
163     /** Creates a new instance.
164      * @param keplerianParameters Keplerian elements
165      * @param aDot  semi-major axis derivative (m/s)
166      * @param eDot eccentricity derivative
167      * @param iDot inclination derivative (rad/s)
168      * @param paDot perigee argument derivative (rad/s)
169      * @param raanDot right ascension of ascending node derivative (rad/s)
170      * @param anomalyDot mean, eccentric or true anomaly derivative, with same type than elements (rad/s)
171      * @param cachedPositionAngleType type of cached anomaly
172      * @param frame the frame in which the parameters are defined
173      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
174      * @param date date of the orbital parameters
175      * @param mu central attraction coefficient (m³/s²)
176      * @exception IllegalArgumentException if frame is not a {@link
177      * Frame#isPseudoInertial pseudo-inertial frame} or a and e don't match for hyperbolic orbits,
178      * or v is out of range for hyperbolic orbits
179      * @since 14.0
180      */
181     public KeplerianOrbit(final KeplerianParameters keplerianParameters,
182                           final double aDot, final double eDot, final double iDot,
183                           final double paDot, final double raanDot, final double anomalyDot,
184                           final PositionAngleType cachedPositionAngleType,
185                           final Frame frame, final AbsoluteDate date, final double mu)
186             throws IllegalArgumentException {
187         super(frame, date, mu);
188         this.cachedPositionAngleType = cachedPositionAngleType;
189 
190         this.a = keplerianParameters.a();
191         this.e = keplerianParameters.e();
192         if (a * (1 - e) < 0) {
193             throw new OrekitIllegalArgumentException(OrekitMessages.ORBIT_A_E_MISMATCH_WITH_CONIC_TYPE, a, e);
194         }
195 
196         // Checking eccentricity range
197         checkParameterRangeInclusive(ECCENTRICITY, e, 0.0, Double.POSITIVE_INFINITY);
198 
199         this.aDot    = aDot;
200         this.eDot    = eDot;
201         this.i       = keplerianParameters.i();
202         this.iDot    = iDot;
203         this.pa      = keplerianParameters.pa();
204         this.paDot   = paDot;
205         this.raan    = keplerianParameters.raan();
206         this.raanDot = raanDot;
207 
208         final UnivariateDerivative1 cachedAnomalyUD = initializeCachedAnomaly(keplerianParameters.anomaly(), anomalyDot,
209                 keplerianParameters.positionAngleType());
210         this.cachedAnomaly = cachedAnomalyUD.getValue();
211         this.cachedAnomalyDot = cachedAnomalyUD.getFirstDerivative();
212 
213         // check true anomaly range
214         if (!isElliptical()) {
215             final double trueAnomaly = getTrueAnomaly();
216             if (1 + e * FastMath.cos(trueAnomaly) <= 0) {
217                 final double vMax = FastMath.acos(-1 / e);
218                 throw new OrekitIllegalArgumentException(OrekitMessages.ORBIT_ANOMALY_OUT_OF_HYPERBOLIC_RANGE,
219                         trueAnomaly, e, -vMax, vMax);
220             }
221         }
222 
223         this.partialPV = null;
224 
225     }
226 
227     /** Creates a new instance with cached position angle same as value inputted.
228      * @param a  semi-major axis (m), negative for hyperbolic orbits
229      * @param e eccentricity (positive or equal to 0)
230      * @param i inclination (rad)
231      * @param pa perigee argument (ω, rad)
232      * @param raan right ascension of ascending node (Ω, rad)
233      * @param anomaly mean, eccentric or true anomaly (rad)
234      * @param aDot  semi-major axis derivative (m/s)
235      * @param eDot eccentricity derivative
236      * @param iDot inclination derivative (rad/s)
237      * @param paDot perigee argument derivative (rad/s)
238      * @param raanDot right ascension of ascending node derivative (rad/s)
239      * @param anomalyDot mean, eccentric or true anomaly derivative (rad/s)
240      * @param type type of anomaly
241      * @param frame the frame in which the parameters are defined
242      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
243      * @param date date of the orbital parameters
244      * @param mu central attraction coefficient (m³/s²)
245      * @exception IllegalArgumentException if frame is not a {@link
246      * Frame#isPseudoInertial pseudo-inertial frame} or a and e don't match for hyperbolic orbits,
247      * or v is out of range for hyperbolic orbits
248      * @since 9.0
249      */
250     public KeplerianOrbit(final double a, final double e, final double i,
251                           final double pa, final double raan, final double anomaly,
252                           final double aDot, final double eDot, final double iDot,
253                           final double paDot, final double raanDot, final double anomalyDot,
254                           final PositionAngleType type,
255                           final Frame frame, final AbsoluteDate date, final double mu)
256             throws IllegalArgumentException {
257         this(new KeplerianParameters(a, e, i, pa, raan, anomaly, type), aDot, eDot, iDot, paDot, raanDot, anomalyDot, type,
258                 frame, date, mu);
259     }
260 
261     /** Constructor from Cartesian parameters.
262      *
263      * <p> The acceleration provided in {@code pvCoordinates} is accessible using
264      * {@link #getPVCoordinates()} and {@link #getPVCoordinates(Frame)}. All other methods
265      * use {@code mu} and the position to compute the acceleration, including
266      * {@link #shiftedBy(double)} and {@link #getPVCoordinates(AbsoluteDate, Frame)}.
267      *
268      * @param pvCoordinates the PVCoordinates of the satellite
269      * @param frame the frame in which are defined the {@link PVCoordinates}
270      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
271      * @param mu central attraction coefficient (m³/s²)
272      * @exception IllegalArgumentException if frame is not a {@link
273      * Frame#isPseudoInertial pseudo-inertial frame}
274      */
275     public KeplerianOrbit(final TimeStampedPVCoordinates pvCoordinates,
276                           final Frame frame, final double mu)
277             throws IllegalArgumentException {
278         this(pvCoordinates, frame, mu, hasNonKeplerianAcceleration(pvCoordinates, mu));
279     }
280 
281     /** Constructor from Cartesian parameters.
282      *
283      * <p> The acceleration provided in {@code pvCoordinates} is accessible using
284      * {@link #getPVCoordinates()} and {@link #getPVCoordinates(Frame)}. All other methods
285      * use {@code mu} and the position to compute the acceleration, including
286      * {@link #shiftedBy(double)} and {@link #getPVCoordinates(AbsoluteDate, Frame)}.
287      *
288      * @param pvCoordinates the PVCoordinates of the satellite
289      * @param frame the frame in which are defined the {@link PVCoordinates}
290      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
291      * @param mu central attraction coefficient (m³/s²)
292      * @param reliableAcceleration if true, the acceleration is considered to be reliable
293      * @exception IllegalArgumentException if frame is not a {@link
294      * Frame#isPseudoInertial pseudo-inertial frame}
295      */
296     private KeplerianOrbit(final TimeStampedPVCoordinates pvCoordinates,
297                            final Frame frame, final double mu,
298                            final boolean reliableAcceleration)
299             throws IllegalArgumentException {
300         super(pvCoordinates, frame, mu);
301 
302         // convert Cartesian to Keplerian elements
303         final KeplerianParametersConverter converter = new KeplerianParametersConverter(mu);
304         cachedPositionAngleType = PositionAngleType.ECCENTRIC;
305         final KeplerianParameters element = converter.toParameters(pvCoordinates, cachedPositionAngleType);
306         a = element.a();
307         e = element.e();
308         i = element.i();
309         raan = element.raan();
310         pa = element.pa();
311         cachedAnomaly = element.anomaly();
312 
313         // Checking eccentricity range
314         checkParameterRangeInclusive(ECCENTRICITY, e, 0.0, Double.POSITIVE_INFINITY);
315 
316         partialPV = pvCoordinates;
317 
318         if (reliableAcceleration) {
319             // we have a relevant acceleration, we can compute derivatives
320 
321             final double[][] jacobian = new double[6][6];
322             getJacobianWrtCartesian(PositionAngleType.MEAN, jacobian);
323 
324             final Vector3D pvP     = pvCoordinates.getPosition();
325             final double   r2      = pvP.getNorm2Sq();
326             final double   r       = FastMath.sqrt(r2);
327             final Vector3D keplerianAcceleration    = new Vector3D(-mu / (r * r2), pvP);
328             final Vector3D pvA     = pvCoordinates.getAcceleration();
329             final Vector3D nonKeplerianAcceleration = pvA.subtract(keplerianAcceleration);
330             final double   aX                       = nonKeplerianAcceleration.getX();
331             final double   aY                       = nonKeplerianAcceleration.getY();
332             final double   aZ                       = nonKeplerianAcceleration.getZ();
333             aDot    = jacobian[0][3] * aX + jacobian[0][4] * aY + jacobian[0][5] * aZ;
334             eDot    = jacobian[1][3] * aX + jacobian[1][4] * aY + jacobian[1][5] * aZ;
335             iDot    = jacobian[2][3] * aX + jacobian[2][4] * aY + jacobian[2][5] * aZ;
336             paDot   = jacobian[3][3] * aX + jacobian[3][4] * aY + jacobian[3][5] * aZ;
337             raanDot = jacobian[4][3] * aX + jacobian[4][4] * aY + jacobian[4][5] * aZ;
338 
339             // in order to compute cached anomaly derivative, we must compute
340             // mean anomaly derivative including Keplerian motion and convert to required anomaly
341             final double MDot = getKeplerianMeanMotion() +
342                     jacobian[5][3] * aX + jacobian[5][4] * aY + jacobian[5][5] * aZ;
343             final UnivariateDerivative1 eUD = new UnivariateDerivative1(e, eDot);
344             final UnivariateDerivative1 MUD = new UnivariateDerivative1(getMeanAnomaly(), MDot);
345             final UnivariateDerivative1 EUD = (a < 0) ?
346                     FieldKeplerianAnomalyUtility.hyperbolicMeanToEccentric(eUD, MUD) :
347                     FieldKeplerianAnomalyUtility.ellipticMeanToEccentric(eUD, MUD);
348             cachedAnomalyDot = EUD.getFirstDerivative();
349 
350         } else {
351             // acceleration is either almost zero or NaN,
352             // we assume acceleration was not known
353             aDot    = 0.;
354             eDot    = 0.;
355             iDot    = 0.;
356             paDot   = 0.;
357             raanDot = 0.;
358             cachedAnomalyDot = computeKeplerianAnomalyDot(cachedPositionAngleType, a, e, mu, cachedAnomaly, cachedPositionAngleType);
359         }
360 
361     }
362 
363     /** Constructor from Cartesian parameters.
364      *
365      * <p> The acceleration provided in {@code pvCoordinates} is accessible using
366      * {@link #getPVCoordinates()} and {@link #getPVCoordinates(Frame)}. All other methods
367      * use {@code mu} and the position to compute the acceleration, including
368      * {@link #shiftedBy(double)} and {@link #getPVCoordinates(AbsoluteDate, Frame)}.
369      *
370      * @param pvCoordinates the PVCoordinates of the satellite
371      * @param frame the frame in which are defined the {@link PVCoordinates}
372      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
373      * @param date date of the orbital parameters
374      * @param mu central attraction coefficient (m³/s²)
375      * @exception IllegalArgumentException if frame is not a {@link
376      * Frame#isPseudoInertial pseudo-inertial frame}
377      */
378     public KeplerianOrbit(final PVCoordinates pvCoordinates,
379                           final Frame frame, final AbsoluteDate date, final double mu)
380             throws IllegalArgumentException {
381         this(new TimeStampedPVCoordinates(date, pvCoordinates), frame, mu);
382     }
383 
384     /** Constructor from any kind of orbital parameters.
385      * @param op orbital parameters to copy
386      */
387     public KeplerianOrbit(final Orbit op) {
388         this(op.getPVCoordinates(), op.getFrame(), op.getMu(), op.hasNonKeplerianAcceleration());
389     }
390 
391     /** {@inheritDoc} */
392     @Override
393     public boolean hasNonKeplerianAcceleration() {
394         return aDot != 0. || eDot != 0. || paDot != 0. || iDot != 0. || raanDot != 0. ||
395                 FastMath.abs(cachedAnomalyDot - computeKeplerianAnomalyDot(cachedPositionAngleType, a, e, getMu(), cachedAnomaly, cachedPositionAngleType)) > TOLERANCE_POSITION_ANGLE_RATE;
396     }
397 
398     /** {@inheritDoc} */
399     @Override
400     public OrbitType getType() {
401         return OrbitType.KEPLERIAN;
402     }
403 
404     /** {@inheritDoc} */
405     @Override
406     public double getA() {
407         return a;
408     }
409 
410     /** {@inheritDoc} */
411     @Override
412     public double getADot() {
413         return aDot;
414     }
415 
416     /** {@inheritDoc} */
417     @Override
418     public double getE() {
419         return e;
420     }
421 
422     /** {@inheritDoc} */
423     @Override
424     public double getEDot() {
425         return eDot;
426     }
427 
428     /** {@inheritDoc} */
429     @Override
430     public double getI() {
431         return i;
432     }
433 
434     /** {@inheritDoc} */
435     @Override
436     public double getIDot() {
437         return iDot;
438     }
439 
440     /** Get the perigee argument.
441      * @return perigee argument (rad)
442      */
443     public double getPerigeeArgument() {
444         return pa;
445     }
446 
447     /** Get the perigee argument derivative.
448      * <p>
449      * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
450      * </p>
451      * @return perigee argument derivative (rad/s)
452      * @since 9.0
453      */
454     public double getPerigeeArgumentDot() {
455         return paDot;
456     }
457 
458     /** Get the right ascension of the ascending node.
459      * @return right ascension of the ascending node (rad)
460      */
461     public double getRightAscensionOfAscendingNode() {
462         return raan;
463     }
464 
465     /** Get the right ascension of the ascending node derivative.
466      * <p>
467      * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
468      * </p>
469      * @return right ascension of the ascending node derivative (rad/s)
470      * @since 9.0
471      */
472     public double getRightAscensionOfAscendingNodeDot() {
473         return raanDot;
474     }
475 
476     /** Get the true anomaly.
477      * @return true anomaly (rad)
478      */
479     public double getTrueAnomaly() {
480         return getAnomaly(PositionAngleType.TRUE);
481     }
482 
483     /** Get the true anomaly derivative.
484      * @return true anomaly derivative (rad/s)
485      */
486     public double getTrueAnomalyDot() {
487         return switch (cachedPositionAngleType) {
488             case MEAN -> {
489                 final UnivariateDerivative1 eUD = new UnivariateDerivative1(e, eDot);
490                 final UnivariateDerivative1 MUD = new UnivariateDerivative1(cachedAnomaly, cachedAnomalyDot);
491                 final UnivariateDerivative1 vUD = (a < 0) ?
492                         FieldKeplerianAnomalyUtility.hyperbolicMeanToTrue(eUD, MUD) :
493                         FieldKeplerianAnomalyUtility.ellipticMeanToTrue(eUD, MUD);
494                 yield vUD.getFirstDerivative();
495             }
496 
497             case TRUE -> cachedAnomalyDot;
498 
499             case ECCENTRIC -> {
500                 final UnivariateDerivative1 eUD2 = new UnivariateDerivative1(e, eDot);
501                 final UnivariateDerivative1 EUD = new UnivariateDerivative1(cachedAnomaly, cachedAnomalyDot);
502                 final UnivariateDerivative1 vUD2 = (a < 0) ?
503                         FieldKeplerianAnomalyUtility.hyperbolicEccentricToTrue(eUD2, EUD) :
504                         FieldKeplerianAnomalyUtility.ellipticEccentricToTrue(eUD2, EUD);
505                 yield vUD2.getFirstDerivative();
506             }
507         };
508     }
509 
510     /** Get the eccentric anomaly.
511      * @return eccentric anomaly (rad)
512      */
513     public double getEccentricAnomaly() {
514         return getAnomaly(PositionAngleType.ECCENTRIC);
515     }
516 
517     /** Get the eccentric anomaly derivative.
518      * @return eccentric anomaly derivative (rad/s)
519      * @since 9.0
520      */
521     public double getEccentricAnomalyDot() {
522         return switch (cachedPositionAngleType) {
523             case ECCENTRIC -> cachedAnomalyDot;
524 
525             case TRUE -> {
526                 final UnivariateDerivative1 eUD = new UnivariateDerivative1(e, eDot);
527                 final UnivariateDerivative1 vUD = new UnivariateDerivative1(cachedAnomaly, cachedAnomalyDot);
528                 final UnivariateDerivative1 EUD = (a < 0) ?
529                         FieldKeplerianAnomalyUtility.hyperbolicTrueToEccentric(eUD, vUD) :
530                         FieldKeplerianAnomalyUtility.ellipticTrueToEccentric(eUD, vUD);
531                 yield EUD.getFirstDerivative();
532             }
533 
534             case MEAN -> {
535                 final UnivariateDerivative1 eUD2 = new UnivariateDerivative1(e, eDot);
536                 final UnivariateDerivative1 MUD = new UnivariateDerivative1(cachedAnomaly, cachedAnomalyDot);
537                 final UnivariateDerivative1 EUD2 = (a < 0) ?
538                         FieldKeplerianAnomalyUtility.hyperbolicMeanToEccentric(eUD2, MUD) :
539                         FieldKeplerianAnomalyUtility.ellipticMeanToEccentric(eUD2, MUD);
540                 yield EUD2.getFirstDerivative();
541             }
542         };
543     }
544 
545     /** Get the mean anomaly.
546      * @return mean anomaly (rad)
547      */
548     public double getMeanAnomaly() {
549         return getAnomaly(PositionAngleType.MEAN);
550     }
551 
552     /** Get the mean anomaly derivative.
553      * @return mean anomaly derivative (rad/s)
554      * @since 9.0
555      */
556     public double getMeanAnomalyDot() {
557         return switch (cachedPositionAngleType) {
558             case MEAN -> cachedAnomalyDot;
559 
560             case ECCENTRIC -> {
561                 final UnivariateDerivative1 eUD = new UnivariateDerivative1(e, eDot);
562                 final UnivariateDerivative1 EUD = new UnivariateDerivative1(cachedAnomaly, cachedAnomalyDot);
563                 final UnivariateDerivative1 MUD = (a < 0) ?
564                         FieldKeplerianAnomalyUtility.hyperbolicEccentricToMean(eUD, EUD) :
565                         FieldKeplerianAnomalyUtility.ellipticEccentricToMean(eUD, EUD);
566                 yield MUD.getFirstDerivative();
567             }
568 
569             case TRUE -> {
570                 final UnivariateDerivative1 eUD2 = new UnivariateDerivative1(e, eDot);
571                 final UnivariateDerivative1 vUD = new UnivariateDerivative1(cachedAnomaly, cachedAnomalyDot);
572                 final UnivariateDerivative1 MUD2 = (a < 0) ?
573                         FieldKeplerianAnomalyUtility.hyperbolicTrueToMean(eUD2, vUD) :
574                         FieldKeplerianAnomalyUtility.ellipticTrueToMean(eUD2, vUD);
575                 yield MUD2.getFirstDerivative();
576             }
577         };
578     }
579 
580     /** Get the anomaly.
581      * @param type type of the angle
582      * @return anomaly (rad)
583      */
584     public double getAnomaly(final PositionAngleType type) {
585         return getKeplerianParameters().withPositionAngleType(type).anomaly();
586     }
587 
588     /** Get the anomaly derivative.
589      * @param type type of the angle
590      * @return anomaly derivative (rad/s)
591      * @since 9.0
592      */
593     public double getAnomalyDot(final PositionAngleType type) {
594         return switch (type) {
595             case MEAN -> getMeanAnomalyDot();
596             case ECCENTRIC -> getEccentricAnomalyDot();
597             case TRUE -> getTrueAnomalyDot();
598         };
599     }
600 
601     /**
602      * Method providing with the Keplerian elements, using the cached type for the anomaly.
603      * @return Keplerian elements
604      * @since 14.0
605      */
606     public KeplerianParameters getKeplerianParameters() {
607         return new KeplerianParameters(a, e, i, pa, raan, cachedAnomaly, cachedPositionAngleType);
608     }
609 
610     /** {@inheritDoc} */
611     @Override
612     public double getEquinoctialEx() {
613         return e * FastMath.cos(pa + raan);
614     }
615 
616     /** {@inheritDoc} */
617     @Override
618     public double getEquinoctialExDot() {
619         if (!hasNonKeplerianAcceleration()) {
620             return 0.;
621         }
622         final double paPraan = pa + raan;
623         final SinCos sc      = FastMath.sinCos(paPraan);
624         return eDot * sc.cos() - e * sc.sin() * (paDot + raanDot);
625     }
626 
627     /** {@inheritDoc} */
628     @Override
629     public double getEquinoctialEy() {
630         return e * FastMath.sin(pa + raan);
631     }
632 
633     /** {@inheritDoc} */
634     @Override
635     public double getEquinoctialEyDot() {
636         if (!hasNonKeplerianAcceleration()) {
637             return 0.;
638         }
639         final double paPraan = pa + raan;
640         final SinCos sc      = FastMath.sinCos(paPraan);
641         return eDot * sc.sin() + e * sc.cos() * (paDot + raanDot);
642     }
643 
644     /** {@inheritDoc} */
645     @Override
646     public double getHx() {
647         // Check for equatorial retrograde orbit
648         if (FastMath.abs(i - FastMath.PI) < 1.0e-10) {
649             return Double.NaN;
650         }
651         return FastMath.cos(raan) * FastMath.tan(0.5 * i);
652     }
653 
654     /** {@inheritDoc} */
655     @Override
656     public double getHxDot() {
657         // Check for equatorial retrograde orbit
658         if (FastMath.abs(i - FastMath.PI) < 1.0e-10) {
659             return Double.NaN;
660         }
661         if (!hasNonKeplerianAcceleration()) {
662             return 0.;
663         }
664         final SinCos sc      = FastMath.sinCos(raan);
665         final double tan     = FastMath.tan(0.5 * i);
666         return 0.5 * (1 + tan * tan) * sc.cos() * iDot - tan * sc.sin() * raanDot;
667     }
668 
669     /** {@inheritDoc} */
670     @Override
671     public double getHy() {
672         // Check for equatorial retrograde orbit
673         if (FastMath.abs(i - FastMath.PI) < 1.0e-10) {
674             return Double.NaN;
675         }
676         return FastMath.sin(raan) * FastMath.tan(0.5 * i);
677     }
678 
679     /** {@inheritDoc} */
680     @Override
681     public double getHyDot() {
682         // Check for equatorial retrograde orbit
683         if (FastMath.abs(i - FastMath.PI) < 1.0e-10) {
684             return Double.NaN;
685         }
686         if (!hasNonKeplerianAcceleration()) {
687             return 0.;
688         }
689         final SinCos sc      = FastMath.sinCos(raan);
690         final double tan     = FastMath.tan(0.5 * i);
691         return 0.5 * (1 + tan * tan) * sc.sin() * iDot + tan * sc.cos() * raanDot;
692     }
693 
694     /** {@inheritDoc} */
695     @Override
696     public double getLv() {
697         return pa + raan + getTrueAnomaly();
698     }
699 
700     /** {@inheritDoc} */
701     @Override
702     public double getLvDot() {
703         return paDot + raanDot + getTrueAnomalyDot();
704     }
705 
706     /** {@inheritDoc} */
707     @Override
708     public double getLE() {
709         return pa + raan + getEccentricAnomaly();
710     }
711 
712     /** {@inheritDoc} */
713     @Override
714     public double getLEDot() {
715         return paDot + raanDot + getEccentricAnomalyDot();
716     }
717 
718     /** {@inheritDoc} */
719     @Override
720     public double getLM() {
721         return pa + raan + getMeanAnomaly();
722     }
723 
724     /** {@inheritDoc} */
725     @Override
726     public double getLMDot() {
727         return paDot + raanDot + getMeanAnomalyDot();
728     }
729 
730     /** Initialize cached anomaly with rate.
731      * @param anomaly input anomaly
732      * @param anomalyDot rate of input anomaly
733      * @param inputType position angle type passed as input
734      * @return anomaly to cache with rate
735      * @since 12.1
736      */
737     private UnivariateDerivative1 initializeCachedAnomaly(final double anomaly, final double anomalyDot,
738                                                           final PositionAngleType inputType) {
739         if (cachedPositionAngleType == inputType) {
740             return new UnivariateDerivative1(anomaly, anomalyDot);
741 
742         } else {
743             final UnivariateDerivative1 eUD = new UnivariateDerivative1(e, eDot);
744             final UnivariateDerivative1 anomalyUD = new UnivariateDerivative1(anomaly, anomalyDot);
745 
746             if (a < 0) {
747                 switch (cachedPositionAngleType) {
748                     case MEAN:
749                         if (inputType == PositionAngleType.ECCENTRIC) {
750                             return FieldKeplerianAnomalyUtility.hyperbolicEccentricToMean(eUD, anomalyUD);
751                         } else {
752                             return FieldKeplerianAnomalyUtility.hyperbolicTrueToMean(eUD, anomalyUD);
753                         }
754 
755                     case ECCENTRIC:
756                         if (inputType == PositionAngleType.MEAN) {
757                             return FieldKeplerianAnomalyUtility.hyperbolicMeanToEccentric(eUD, anomalyUD);
758                         } else {
759                             return FieldKeplerianAnomalyUtility.hyperbolicTrueToEccentric(eUD, anomalyUD);
760                         }
761 
762                     case TRUE:
763                         if (inputType == PositionAngleType.MEAN) {
764                             return FieldKeplerianAnomalyUtility.hyperbolicMeanToTrue(eUD, anomalyUD);
765                         } else {
766                             return FieldKeplerianAnomalyUtility.hyperbolicEccentricToTrue(eUD, anomalyUD);
767                         }
768 
769                     default:
770                         break;
771                 }
772 
773             } else {
774                 switch (cachedPositionAngleType) {
775                     case MEAN:
776                         if (inputType == PositionAngleType.ECCENTRIC) {
777                             return FieldKeplerianAnomalyUtility.ellipticEccentricToMean(eUD, anomalyUD);
778                         } else {
779                             return FieldKeplerianAnomalyUtility.ellipticTrueToMean(eUD, anomalyUD);
780                         }
781 
782                     case ECCENTRIC:
783                         if (inputType == PositionAngleType.MEAN) {
784                             return FieldKeplerianAnomalyUtility.ellipticMeanToEccentric(eUD, anomalyUD);
785                         } else {
786                             return FieldKeplerianAnomalyUtility.ellipticTrueToEccentric(eUD, anomalyUD);
787                         }
788 
789                     case TRUE:
790                         if (inputType == PositionAngleType.MEAN) {
791                             return FieldKeplerianAnomalyUtility.ellipticMeanToTrue(eUD, anomalyUD);
792                         } else {
793                             return FieldKeplerianAnomalyUtility.ellipticEccentricToTrue(eUD, anomalyUD);
794                         }
795 
796                     default:
797                         break;
798                 }
799 
800             }
801             throw new OrekitInternalError(null);
802         }
803 
804     }
805 
806     /** Compute position and velocity but not acceleration.
807      */
808     private void computePVWithoutA() {
809 
810         if (partialPV != null) {
811             // already computed
812             return;
813         }
814 
815         final KeplerianParametersConverter converter = new KeplerianParametersConverter(getMu());
816         partialPV = converter.toCartesian(getKeplerianParameters());
817 
818     }
819 
820     /** {@inheritDoc} */
821     @Override
822     protected Vector3D initPosition() {
823 
824         final Vector3D[] axes = KeplerianParametersConverter.referenceAxes(i, pa, raan);
825 
826         if (isElliptical()) {
827 
828             // elliptical case
829 
830             // elliptic eccentric anomaly
831             final double uME2   = (1 - e) * (1 + e);
832             final double s1Me2  = FastMath.sqrt(uME2);
833             final SinCos scE    = FastMath.sinCos(getEccentricAnomaly());
834             final double cosE   = scE.cos();
835             final double sinE   = scE.sin();
836 
837             return new Vector3D(a * (cosE - e), axes[0], a * sinE * s1Me2, axes[1]);
838 
839         } else {
840 
841             // hyperbolic case
842 
843             // compute position and velocity factors
844             final SinCos scV       = FastMath.sinCos(getTrueAnomaly());
845             final double sinV      = scV.sin();
846             final double cosV      = scV.cos();
847             final double f         = a * (1 - e * e);
848             final double posFactor = f / (1 + e * cosV);
849 
850             return new Vector3D(posFactor * cosV, axes[0], posFactor * sinV, axes[1]);
851 
852         }
853 
854     }
855 
856     /** {@inheritDoc} */
857     @Override
858     protected TimeStampedPVCoordinates initPVCoordinates() {
859 
860         // position and velocity
861         computePVWithoutA();
862 
863         // acceleration
864         final double r2 = partialPV.getPosition().getNorm2Sq();
865         final Vector3D keplerianAcceleration = new Vector3D(-getMu() / (r2 * FastMath.sqrt(r2)), partialPV.getPosition());
866         final Vector3D acceleration = hasNonKeplerianAcceleration() ?
867                 keplerianAcceleration.add(nonKeplerianAcceleration()) :
868                 keplerianAcceleration;
869 
870         return new TimeStampedPVCoordinates(getDate(), partialPV.getPosition(), partialPV.getVelocity(), acceleration);
871 
872     }
873 
874     /** {@inheritDoc} */
875     @Override
876     public KeplerianOrbit inFrame(final Frame inertialFrame) {
877         final PVCoordinates pvCoordinates;
878         if (hasNonKeplerianAcceleration()) {
879             pvCoordinates = getPVCoordinates(inertialFrame);
880         } else {
881             final KinematicTransform transform = getFrame().getKinematicTransformTo(inertialFrame, getDate());
882             pvCoordinates = transform.transformOnlyPV(getPVCoordinates());
883         }
884         final KeplerianOrbit keplerianOrbit = new KeplerianOrbit(pvCoordinates, inertialFrame, getDate(), getMu());
885         if (keplerianOrbit.getCachedPositionAngleType() == getCachedPositionAngleType()) {
886             return keplerianOrbit;
887         } else {
888             return keplerianOrbit.withCachedPositionAngleType(getCachedPositionAngleType());
889         }
890     }
891 
892     /** {@inheritDoc} */
893     @Override
894     public KeplerianOrbit withCachedPositionAngleType(final PositionAngleType positionAngleType) {
895         return new KeplerianOrbit(a, e, i, pa, raan, getAnomaly(positionAngleType), aDot, eDot, iDot, paDot, raanDot,
896                 getAnomalyDot(positionAngleType), positionAngleType, getFrame(), getDate(), getMu());
897     }
898 
899     /** {@inheritDoc} */
900     @Override
901     public KeplerianOrbit shiftedBy(final double dt) {
902         return shiftedBy(new TimeOffset(dt));
903     }
904 
905     /** {@inheritDoc} */
906     @Override
907     public KeplerianOrbit shiftedBy(final TimeOffset dt) {
908 
909         final double dtS = dt.toDouble();
910 
911         // use Keplerian-only motion
912         final KeplerianParameters shiftedElements = new KeplerianParameters(a, e, i, pa, raan, getMeanAnomaly() + getKeplerianMeanMotion() * dtS,
913                 PositionAngleType.MEAN).withPositionAngleType(cachedPositionAngleType);
914         final KeplerianOrbit keplerianShifted = new KeplerianOrbit(shiftedElements, getFrame(),
915                 getDate().shiftedBy(dt), getMu());
916 
917         if (dtS != 0. && hasNonKeplerianAcceleration()) {
918 
919             // extract non-Keplerian acceleration from first time derivatives
920             final Vector3D nonKeplerianAcceleration = nonKeplerianAcceleration();
921 
922             // add quadratic effect of non-Keplerian acceleration to Keplerian-only shift
923             keplerianShifted.computePVWithoutA();
924             final Vector3D fixedP   = new Vector3D(1, keplerianShifted.partialPV.getPosition(),
925                     0.5 * dtS * dtS, nonKeplerianAcceleration);
926             final double   fixedR2 = fixedP.getNorm2Sq();
927             final double   fixedR  = FastMath.sqrt(fixedR2);
928             final Vector3D fixedV  = new Vector3D(1, keplerianShifted.partialPV.getVelocity(),
929                     dtS, nonKeplerianAcceleration);
930             final Vector3D fixedA  = new Vector3D(-getMu() / (fixedR2 * fixedR), keplerianShifted.partialPV.getPosition(),
931                     1, nonKeplerianAcceleration);
932 
933             // build a new orbit, taking non-Keplerian acceleration into account
934             return new KeplerianOrbit(new TimeStampedPVCoordinates(keplerianShifted.getDate(),
935                     fixedP, fixedV, fixedA),
936                     keplerianShifted.getFrame(), keplerianShifted.getMu());
937 
938         } else {
939             // Keplerian-only motion is all we can do
940             return keplerianShifted;
941         }
942 
943     }
944 
945     /** {@inheritDoc} */
946     @Override
947     protected double[][] computeJacobianMeanWrtCartesian() {
948         if (isElliptical()) {
949             return computeJacobianMeanWrtCartesianElliptical();
950         } else {
951             return computeJacobianMeanWrtCartesianHyperbolic();
952         }
953     }
954 
955     /** Compute the Jacobian of the orbital parameters with respect to the Cartesian parameters.
956      * <p>
957      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
958      * respect to Cartesian coordinate j (x for j=0, y for j=1, z for j=2, xDot for j=3,
959      * yDot for j=4, zDot for j=5).
960      * </p>
961      * @return 6x6 Jacobian matrix
962      */
963     private double[][] computeJacobianMeanWrtCartesianElliptical() {
964 
965         final double[][] jacobian = new double[6][6];
966 
967         // compute various intermediate parameters
968         computePVWithoutA();
969         final Vector3D position = partialPV.getPosition();
970         final Vector3D velocity = partialPV.getVelocity();
971         final Vector3D momentum = partialPV.getMomentum();
972         final double v2         = velocity.getNorm2Sq();
973         final double r2         = position.getNorm2Sq();
974         final double r          = FastMath.sqrt(r2);
975         final double r3         = r * r2;
976 
977         final double px         = position.getX();
978         final double py         = position.getY();
979         final double pz         = position.getZ();
980         final double vx         = velocity.getX();
981         final double vy         = velocity.getY();
982         final double vz         = velocity.getZ();
983         final double mx         = momentum.getX();
984         final double my         = momentum.getY();
985         final double mz         = momentum.getZ();
986 
987         final double mu         = getMu();
988         final double sqrtMuA    = FastMath.sqrt(a * mu);
989         final double sqrtAoMu   = FastMath.sqrt(a / mu);
990         final double a2         = a * a;
991         final double twoA       = 2 * a;
992         final double rOnA       = r / a;
993 
994         final double oMe2       = 1 - e * e;
995         final double epsilon    = FastMath.sqrt(oMe2);
996         final double sqrtRec    = 1 / epsilon;
997 
998         final SinCos scI        = FastMath.sinCos(i);
999         final SinCos scPA       = FastMath.sinCos(pa);
1000         final double cosI       = scI.cos();
1001         final double sinI       = scI.sin();
1002         final double cosPA      = scPA.cos();
1003         final double sinPA      = scPA.sin();
1004 
1005         final double pv         = Vector3D.dotProduct(position, velocity);
1006         final double cosE       = (a - r) / (a * e);
1007         final double sinE       = pv / (e * sqrtMuA);
1008 
1009         // da
1010         final Vector3D vectorAR = new Vector3D(2 * a2 / r3, position);
1011         final Vector3D vectorARDot = velocity.scalarMultiply(2 * a2 / mu);
1012         fillHalfRow(1, vectorAR,    jacobian[0], 0);
1013         fillHalfRow(1, vectorARDot, jacobian[0], 3);
1014 
1015         // de
1016         final double factorER3 = pv / twoA;
1017         final Vector3D vectorER   = new Vector3D(cosE * v2 / (r * mu), position,
1018                 sinE / sqrtMuA, velocity,
1019                 -factorER3 * sinE / sqrtMuA, vectorAR);
1020         final Vector3D vectorERDot = new Vector3D(sinE / sqrtMuA, position,
1021                 cosE * 2 * r / mu, velocity,
1022                 -factorER3 * sinE / sqrtMuA, vectorARDot);
1023         fillHalfRow(1, vectorER,    jacobian[1], 0);
1024         fillHalfRow(1, vectorERDot, jacobian[1], 3);
1025 
1026         // dE / dr (Eccentric anomaly)
1027         final double coefE = cosE / (e * sqrtMuA);
1028         final Vector3D  vectorEAnR =
1029                 new Vector3D(-sinE * v2 / (e * r * mu), position, coefE, velocity,
1030                         -factorER3 * coefE, vectorAR);
1031 
1032         // dE / drDot
1033         final Vector3D  vectorEAnRDot =
1034                 new Vector3D(-sinE * 2 * r / (e * mu), velocity, coefE, position,
1035                         -factorER3 * coefE, vectorARDot);
1036 
1037         // precomputing some more factors
1038         final double s1 = -sinE * pz / r - cosE * vz * sqrtAoMu;
1039         final double s2 = -cosE * pz / r3;
1040         final double s3 = -sinE * vz / (2 * sqrtMuA);
1041         final double t1 = sqrtRec * (cosE * pz / r - sinE * vz * sqrtAoMu);
1042         final double t2 = sqrtRec * (-sinE * pz / r3);
1043         final double t3 = sqrtRec * (cosE - e) * vz / (2 * sqrtMuA);
1044         final double t4 = sqrtRec * (e * sinI * cosPA * sqrtRec - vz * sqrtAoMu);
1045         final Vector3D s = new Vector3D(cosE / r, Vector3D.PLUS_K,
1046                 s1,       vectorEAnR,
1047                 s2,       position,
1048                 s3,       vectorAR);
1049         final Vector3D sDot = new Vector3D(-sinE * sqrtAoMu, Vector3D.PLUS_K,
1050                 s1,               vectorEAnRDot,
1051                 s3,               vectorARDot);
1052         final Vector3D t =
1053                 new Vector3D(sqrtRec * sinE / r, Vector3D.PLUS_K).add(new Vector3D(t1, vectorEAnR,
1054                         t2, position,
1055                         t3, vectorAR,
1056                         t4, vectorER));
1057         final Vector3D tDot = new Vector3D(sqrtRec * (cosE - e) * sqrtAoMu, Vector3D.PLUS_K,
1058                 t1,                              vectorEAnRDot,
1059                 t3,                              vectorARDot,
1060                 t4,                              vectorERDot);
1061 
1062         // di
1063         final double factorI1 = -sinI * sqrtRec / sqrtMuA;
1064         final double i1 =  factorI1;
1065         final double i2 = -factorI1 * mz / twoA;
1066         final double i3 =  factorI1 * mz * e / oMe2;
1067         final double i4 = cosI * sinPA;
1068         final double i5 = cosI * cosPA;
1069         fillHalfRow(i1, new Vector3D(vy, -vx, 0), i2, vectorAR, i3, vectorER, i4, s, i5, t,
1070                 jacobian[2], 0);
1071         fillHalfRow(i1, new Vector3D(-py, px, 0), i2, vectorARDot, i3, vectorERDot, i4, sDot, i5, tDot,
1072                 jacobian[2], 3);
1073 
1074         // dpa
1075         fillHalfRow(cosPA / sinI, s,    -sinPA / sinI, t,    jacobian[3], 0);
1076         fillHalfRow(cosPA / sinI, sDot, -sinPA / sinI, tDot, jacobian[3], 3);
1077 
1078         // dRaan
1079         final double factorRaanR = 1 / (mu * a * oMe2 * sinI * sinI);
1080         fillHalfRow(-factorRaanR * my, new Vector3D(  0, vz, -vy),
1081                 factorRaanR * mx, new Vector3D(-vz,  0,  vx),
1082                 jacobian[4], 0);
1083         fillHalfRow(-factorRaanR * my, new Vector3D( 0, -pz,  py),
1084                 factorRaanR * mx, new Vector3D(pz,   0, -px),
1085                 jacobian[4], 3);
1086 
1087         // dM
1088         fillHalfRow(rOnA, vectorEAnR,    -sinE, vectorER,    jacobian[5], 0);
1089         fillHalfRow(rOnA, vectorEAnRDot, -sinE, vectorERDot, jacobian[5], 3);
1090 
1091         return jacobian;
1092 
1093     }
1094 
1095     /** Compute the Jacobian of the orbital parameters with respect to the Cartesian parameters.
1096      * <p>
1097      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
1098      * respect to Cartesian coordinate j (x for j=0, y for j=1, z for j=2, xDot for j=3,
1099      * yDot for j=4, zDot for j=5).
1100      * </p>
1101      * @return 6x6 Jacobian matrix
1102      */
1103     private double[][] computeJacobianMeanWrtCartesianHyperbolic() {
1104 
1105         final double[][] jacobian = new double[6][6];
1106 
1107         // compute various intermediate parameters
1108         computePVWithoutA();
1109         final Vector3D position = partialPV.getPosition();
1110         final Vector3D velocity = partialPV.getVelocity();
1111         final Vector3D momentum = partialPV.getMomentum();
1112         final double r2         = position.getNorm2Sq();
1113         final double r          = FastMath.sqrt(r2);
1114         final double r3         = r * r2;
1115 
1116         final double x          = position.getX();
1117         final double y          = position.getY();
1118         final double z          = position.getZ();
1119         final double vx         = velocity.getX();
1120         final double vy         = velocity.getY();
1121         final double vz         = velocity.getZ();
1122         final double mx         = momentum.getX();
1123         final double my         = momentum.getY();
1124         final double mz         = momentum.getZ();
1125 
1126         final double mu         = getMu();
1127         final double absA       = -a;
1128         final double sqrtMuA    = FastMath.sqrt(absA * mu);
1129         final double a2         = a * a;
1130         final double rOa        = r / absA;
1131 
1132         final SinCos scI        = FastMath.sinCos(i);
1133         final double cosI       = scI.cos();
1134         final double sinI       = scI.sin();
1135 
1136         final double pv         = Vector3D.dotProduct(position, velocity);
1137 
1138         // da
1139         final Vector3D vectorAR = new Vector3D(-2 * a2 / r3, position);
1140         final Vector3D vectorARDot = velocity.scalarMultiply(-2 * a2 / mu);
1141         fillHalfRow(-1, vectorAR,    jacobian[0], 0);
1142         fillHalfRow(-1, vectorARDot, jacobian[0], 3);
1143 
1144         // differentials of the momentum
1145         final double m      = momentum.getNorm();
1146         final double oOm    = 1 / m;
1147         final Vector3D dcXP = new Vector3D(  0,  vz, -vy);
1148         final Vector3D dcYP = new Vector3D(-vz,   0,  vx);
1149         final Vector3D dcZP = new Vector3D( vy, -vx,   0);
1150         final Vector3D dcXV = new Vector3D(  0,  -z,   y);
1151         final Vector3D dcYV = new Vector3D(  z,   0,  -x);
1152         final Vector3D dcZV = new Vector3D( -y,   x,   0);
1153         final Vector3D dCP  = new Vector3D(mx * oOm, dcXP, my * oOm, dcYP, mz * oOm, dcZP);
1154         final Vector3D dCV  = new Vector3D(mx * oOm, dcXV, my * oOm, dcYV, mz * oOm, dcZV);
1155 
1156         // dp
1157         final double mOMu   = m / mu;
1158         final Vector3D dpP  = new Vector3D(2 * mOMu, dCP);
1159         final Vector3D dpV  = new Vector3D(2 * mOMu, dCV);
1160 
1161         // de
1162         final double p      = m * mOMu;
1163         final double moO2ae = 1 / (2 * absA * e);
1164         final double m2OaMu = -p / absA;
1165         fillHalfRow(moO2ae, dpP, m2OaMu * moO2ae, vectorAR,    jacobian[1], 0);
1166         fillHalfRow(moO2ae, dpV, m2OaMu * moO2ae, vectorARDot, jacobian[1], 3);
1167 
1168         // di
1169         final double cI1 = 1 / (m * sinI);
1170         final double cI2 = cosI * cI1;
1171         fillHalfRow(cI2, dCP, -cI1, dcZP, jacobian[2], 0);
1172         fillHalfRow(cI2, dCV, -cI1, dcZV, jacobian[2], 3);
1173 
1174         // dPA
1175         final double cP1     =  y * oOm;
1176         final double cP2     = -x * oOm;
1177         final double cP3     = -(mx * cP1 + my * cP2);
1178         final double cP4     = cP3 * oOm;
1179         final double cP5     = -1 / (r2 * sinI * sinI);
1180         final double cP6     = z  * cP5;
1181         final double cP7     = cP3 * cP5;
1182         final Vector3D dacP  = new Vector3D(cP1, dcXP, cP2, dcYP, cP4, dCP, oOm, new Vector3D(-my, mx, 0));
1183         final Vector3D dacV  = new Vector3D(cP1, dcXV, cP2, dcYV, cP4, dCV);
1184         final Vector3D dpoP  = new Vector3D(cP6, dacP, cP7, Vector3D.PLUS_K);
1185         final Vector3D dpoV  = new Vector3D(cP6, dacV);
1186 
1187         final double re2     = r2 * e * e;
1188         final double recOre2 = (p - r) / re2;
1189         final double resOre2 = (pv * mOMu) / re2;
1190         final Vector3D dreP  = new Vector3D(mOMu, velocity, pv / mu, dCP);
1191         final Vector3D dreV  = new Vector3D(mOMu, position, pv / mu, dCV);
1192         final Vector3D davP  = new Vector3D(-resOre2, dpP, recOre2, dreP, resOre2 / r, position);
1193         final Vector3D davV  = new Vector3D(-resOre2, dpV, recOre2, dreV);
1194         fillHalfRow(1, dpoP, -1, davP, jacobian[3], 0);
1195         fillHalfRow(1, dpoV, -1, davV, jacobian[3], 3);
1196 
1197         // dRAAN
1198         final double cO0 = cI1 * cI1;
1199         final double cO1 =  mx * cO0;
1200         final double cO2 = -my * cO0;
1201         fillHalfRow(cO1, dcYP, cO2, dcXP, jacobian[4], 0);
1202         fillHalfRow(cO1, dcYV, cO2, dcXV, jacobian[4], 3);
1203 
1204         // dM
1205         final double s2a    = pv / (2 * absA);
1206         final double oObux  = 1 / FastMath.sqrt(m * m + mu * absA);
1207         final double scasbu = pv * oObux;
1208         final Vector3D dauP = new Vector3D(1 / sqrtMuA, velocity, -s2a / sqrtMuA, vectorAR);
1209         final Vector3D dauV = new Vector3D(1 / sqrtMuA, position, -s2a / sqrtMuA, vectorARDot);
1210         final Vector3D dbuP = new Vector3D(oObux * mu / 2, vectorAR,    m * oObux, dCP);
1211         final Vector3D dbuV = new Vector3D(oObux * mu / 2, vectorARDot, m * oObux, dCV);
1212         final Vector3D dcuP = new Vector3D(oObux, velocity, -scasbu * oObux, dbuP);
1213         final Vector3D dcuV = new Vector3D(oObux, position, -scasbu * oObux, dbuV);
1214         fillHalfRow(1, dauP, -e / (1 + rOa), dcuP, jacobian[5], 0);
1215         fillHalfRow(1, dauV, -e / (1 + rOa), dcuV, jacobian[5], 3);
1216 
1217         return jacobian;
1218 
1219     }
1220 
1221     /** {@inheritDoc} */
1222     @Override
1223     protected double[][] computeJacobianEccentricWrtCartesian() {
1224         if (isElliptical()) {
1225             return computeJacobianEccentricWrtCartesianElliptical();
1226         } else {
1227             return computeJacobianEccentricWrtCartesianHyperbolic();
1228         }
1229     }
1230 
1231     /** Compute the Jacobian of the orbital parameters with respect to the Cartesian parameters.
1232      * <p>
1233      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
1234      * respect to Cartesian coordinate j (x for j=0, y for j=1, z for j=2, xDot for j=3,
1235      * yDot for j=4, zDot for j=5).
1236      * </p>
1237      * @return 6x6 Jacobian matrix
1238      */
1239     private double[][] computeJacobianEccentricWrtCartesianElliptical() {
1240 
1241         // start by computing the Jacobian with mean angle
1242         final double[][] jacobian = computeJacobianMeanWrtCartesianElliptical();
1243 
1244         // Differentiating the Kepler equation M = E - e sin E leads to:
1245         // dM = (1 - e cos E) dE - sin E de
1246         // which is inverted and rewritten as:
1247         // dE = a/r dM + sin E a/r de
1248         final SinCos scE              = FastMath.sinCos(getEccentricAnomaly());
1249         final double aOr              = 1 / (1 - e * scE.cos());
1250 
1251         // update anomaly row
1252         final double[] eRow           = jacobian[1];
1253         final double[] anomalyRow     = jacobian[5];
1254         for (int j = 0; j < anomalyRow.length; ++j) {
1255             anomalyRow[j] = aOr * (anomalyRow[j] + scE.sin() * eRow[j]);
1256         }
1257 
1258         return jacobian;
1259 
1260     }
1261 
1262     /** Compute the Jacobian of the orbital parameters with respect to the Cartesian parameters.
1263      * <p>
1264      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
1265      * respect to Cartesian coordinate j (x for j=0, y for j=1, z for j=2, xDot for j=3,
1266      * yDot for j=4, zDot for j=5).
1267      * </p>
1268      * @return 6x6 Jacobian matrix
1269      */
1270     private double[][] computeJacobianEccentricWrtCartesianHyperbolic() {
1271 
1272         // start by computing the Jacobian with mean angle
1273         final double[][] jacobian = computeJacobianMeanWrtCartesianHyperbolic();
1274 
1275         // Differentiating the Kepler equation M = e sinh H - H leads to:
1276         // dM = (e cosh H - 1) dH + sinh H de
1277         // which is inverted and rewritten as:
1278         // dH = 1 / (e cosh H - 1) dM - sinh H / (e cosh H - 1) de
1279         final double H      = getEccentricAnomaly();
1280         final double coshH  = FastMath.cosh(H);
1281         final double sinhH  = FastMath.sinh(H);
1282         final double absaOr = 1 / (e * coshH - 1);
1283 
1284         // update anomaly row
1285         final double[] eRow       = jacobian[1];
1286         final double[] anomalyRow = jacobian[5];
1287         for (int j = 0; j < anomalyRow.length; ++j) {
1288             anomalyRow[j] = absaOr * (anomalyRow[j] - sinhH * eRow[j]);
1289         }
1290 
1291         return jacobian;
1292 
1293     }
1294 
1295     /** {@inheritDoc} */
1296     @Override
1297     protected double[][] computeJacobianTrueWrtCartesian() {
1298         if (isElliptical()) {
1299             return computeJacobianTrueWrtCartesianElliptical();
1300         } else {
1301             return computeJacobianTrueWrtCartesianHyperbolic();
1302         }
1303     }
1304 
1305     /** Compute the Jacobian of the orbital parameters with respect to the Cartesian parameters.
1306      * <p>
1307      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
1308      * respect to Cartesian coordinate j (x for j=0, y for j=1, z for j=2, xDot for j=3,
1309      * yDot for j=4, zDot for j=5).
1310      * </p>
1311      * @return 6x6 Jacobian matrix
1312      */
1313     private double[][] computeJacobianTrueWrtCartesianElliptical() {
1314 
1315         // start by computing the Jacobian with eccentric angle
1316         final double[][] jacobian = computeJacobianEccentricWrtCartesianElliptical();
1317 
1318         // Differentiating the eccentric anomaly equation sin E = sqrt(1-e^2) sin v / (1 + e cos v)
1319         // and using cos E = (e + cos v) / (1 + e cos v) to get rid of cos E leads to:
1320         // dE = [sqrt (1 - e^2) / (1 + e cos v)] dv - [sin E / (1 - e^2)] de
1321         // which is inverted and rewritten as:
1322         // dv = sqrt (1 - e^2) a/r dE + [sin E / sqrt (1 - e^2)] a/r de
1323         final double e2           = e * e;
1324         final double oMe2         = 1 - e2;
1325         final double epsilon      = FastMath.sqrt(oMe2);
1326         final SinCos scE          = FastMath.sinCos(getEccentricAnomaly());
1327         final double aOr          = 1 / (1 - e * scE.cos());
1328         final double aFactor      = epsilon * aOr;
1329         final double eFactor      = scE.sin() * aOr / epsilon;
1330 
1331         // update anomaly row
1332         final double[] eRow       = jacobian[1];
1333         final double[] anomalyRow = jacobian[5];
1334         for (int j = 0; j < anomalyRow.length; ++j) {
1335             anomalyRow[j] = aFactor * anomalyRow[j] + eFactor * eRow[j];
1336         }
1337 
1338         return jacobian;
1339 
1340     }
1341 
1342     /** Compute the Jacobian of the orbital parameters with respect to the Cartesian parameters.
1343      * <p>
1344      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
1345      * respect to Cartesian coordinate j (x for j=0, y for j=1, z for j=2, xDot for j=3,
1346      * yDot for j=4, zDot for j=5).
1347      * </p>
1348      * @return 6x6 Jacobian matrix
1349      */
1350     private double[][] computeJacobianTrueWrtCartesianHyperbolic() {
1351 
1352         // start by computing the Jacobian with eccentric angle
1353         final double[][] jacobian = computeJacobianEccentricWrtCartesianHyperbolic();
1354 
1355         // Differentiating the eccentric anomaly equation sinh H = sqrt(e^2-1) sin v / (1 + e cos v)
1356         // and using cosh H = (e + cos v) / (1 + e cos v) to get rid of cosh H leads to:
1357         // dH = [sqrt (e^2 - 1) / (1 + e cos v)] dv + [sinh H / (e^2 - 1)] de
1358         // which is inverted and rewritten as:
1359         // dv = sqrt (1 - e^2) a/r dH - [sinh H / sqrt (e^2 - 1)] a/r de
1360         final double e2       = e * e;
1361         final double e2Mo     = e2 - 1;
1362         final double epsilon  = FastMath.sqrt(e2Mo);
1363         final double H        = getEccentricAnomaly();
1364         final double coshH    = FastMath.cosh(H);
1365         final double sinhH    = FastMath.sinh(H);
1366         final double aOr      = 1 / (e * coshH - 1);
1367         final double aFactor  = epsilon * aOr;
1368         final double eFactor  = sinhH * aOr / epsilon;
1369 
1370         // update anomaly row
1371         final double[] eRow           = jacobian[1];
1372         final double[] anomalyRow     = jacobian[5];
1373         for (int j = 0; j < anomalyRow.length; ++j) {
1374             anomalyRow[j] = aFactor * anomalyRow[j] - eFactor * eRow[j];
1375         }
1376 
1377         return jacobian;
1378 
1379     }
1380 
1381     /** {@inheritDoc} */
1382     @Override
1383     public void addKeplerContribution(final PositionAngleType type, final double gm,
1384                                       final double[] pDot) {
1385         pDot[5] += computeKeplerianAnomalyDot(type, a, e, gm, cachedAnomaly, cachedPositionAngleType);
1386     }
1387 
1388     /**
1389      * Compute rate of argument of latitude.
1390      * @param type position angle type of output rate
1391      * @param a semi major axis
1392      * @param e eccentricity
1393      * @param mu mu
1394      * @param anomaly anomaly
1395      * @param cachedType position angle type of input anomaly
1396      * @return first-order time derivative for anomaly
1397      * @since 12.2
1398      */
1399     private static double computeKeplerianAnomalyDot(final PositionAngleType type, final double a, final double e,
1400                                                      final double mu, final double anomaly, final PositionAngleType cachedType) {
1401         final double absA = FastMath.abs(a);
1402         final double n    = FastMath.sqrt(mu / absA) / absA;
1403         if (type == PositionAngleType.MEAN) {
1404             return n;
1405         }
1406         final double oMe2 = FastMath.abs(1 - e * e);
1407         final double ksi = 1 + e * FastMath.cos(KeplerianAnomalyUtility.convertAnomaly(cachedType, anomaly, e, PositionAngleType.TRUE));
1408         if (type == PositionAngleType.ECCENTRIC) {
1409             return n * ksi / oMe2;
1410         } else { // TRUE
1411             return n * ksi * ksi / (oMe2 * FastMath.sqrt(oMe2));
1412         }
1413     }
1414 
1415     /**  Returns a string representation of this Keplerian parameters object.
1416      * @return a string representation of this object
1417      */
1418     public String toString() {
1419         return "Keplerian parameters: " + '{' +
1420                 "a: " + a +
1421                 "; e: " + e +
1422                 "; i: " + FastMath.toDegrees(i) +
1423                 "; pa: " + FastMath.toDegrees(pa) +
1424                 "; raan: " + FastMath.toDegrees(raan) +
1425                 "; v: " + FastMath.toDegrees(getTrueAnomaly()) +
1426                 ";}";
1427     }
1428 
1429     /** {@inheritDoc} */
1430     @Override
1431     public PositionAngleType getCachedPositionAngleType() {
1432         return cachedPositionAngleType;
1433     }
1434 
1435     /** {@inheritDoc} */
1436     @Override
1437     public boolean hasNonKeplerianRates() {
1438         return hasNonKeplerianAcceleration();
1439     }
1440 
1441     /** {@inheritDoc} */
1442     @Override
1443     public KeplerianOrbit withKeplerianRates() {
1444         return new KeplerianOrbit(getKeplerianParameters(), getFrame(), getDate(), getMu());
1445     }
1446 
1447     /** Check if the given parameter is within an acceptable range.
1448      * The bounds are inclusive: an exception is raised when either of those conditions are met:
1449      * <ul>
1450      *     <li>The parameter is strictly greater than upperBound</li>
1451      *     <li>The parameter is strictly lower than lowerBound</li>
1452      * </ul>
1453      * <p>
1454      * In either of these cases, an OrekitException is raised.
1455      * </p>
1456      * @param parameterName name of the parameter
1457      * @param parameter value of the parameter
1458      * @param lowerBound lower bound of the acceptable range (inclusive)
1459      * @param upperBound upper bound of the acceptable range (inclusive)
1460      */
1461     private void checkParameterRangeInclusive(final String parameterName, final double parameter,
1462                                               final double lowerBound, final double upperBound) {
1463         if (parameter < lowerBound || parameter > upperBound) {
1464             throw new OrekitException(OrekitMessages.INVALID_PARAMETER_RANGE, parameterName,
1465                     parameter, lowerBound, upperBound);
1466         }
1467     }
1468 
1469 }