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