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