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