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