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