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