1   /* Copyright 2002-2019 CS Systèmes d'Information
2    * Licensed to CS Systèmes d'Information (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.DSFactory;
25  import org.hipparchus.analysis.differentiation.DerivativeStructure;
26  import org.hipparchus.analysis.interpolation.HermiteInterpolator;
27  import org.hipparchus.geometry.euclidean.threed.Vector3D;
28  import org.hipparchus.util.FastMath;
29  import org.hipparchus.util.MathUtils;
30  import org.orekit.errors.OrekitIllegalArgumentException;
31  import org.orekit.errors.OrekitInternalError;
32  import org.orekit.errors.OrekitMessages;
33  import org.orekit.frames.Frame;
34  import org.orekit.time.AbsoluteDate;
35  import org.orekit.utils.PVCoordinates;
36  import org.orekit.utils.TimeStampedPVCoordinates;
37  
38  
39  /**
40   * This class handles equinoctial orbital parameters, which can support both
41   * circular and equatorial orbits.
42   * <p>
43   * The parameters used internally are the equinoctial elements which can be
44   * related to Keplerian elements as follows:
45   *   <pre>
46   *     a
47   *     ex = e cos(ω + Ω)
48   *     ey = e sin(ω + Ω)
49   *     hx = tan(i/2) cos(Ω)
50   *     hy = tan(i/2) sin(Ω)
51   *     lv = v + ω + Ω
52   *   </pre>
53   * where ω stands for the Perigee Argument and Ω stands for the
54   * Right Ascension of the Ascending Node.
55   *
56   * <p>
57   * The conversion equations from and to Keplerian elements given above hold only
58   * when both sides are unambiguously defined, i.e. when orbit is neither equatorial
59   * nor circular. When orbit is either equatorial or circular, the equinoctial
60   * parameters are still unambiguously defined whereas some Keplerian elements
61   * (more precisely ω and Ω) become ambiguous. For this reason, equinoctial
62   * parameters are the recommended way to represent orbits.
63   * </p>
64   * <p>
65   * The instance <code>EquinoctialOrbit</code> is guaranteed to be immutable.
66   * </p>
67   * @see    Orbit
68   * @see    KeplerianOrbit
69   * @see    CircularOrbit
70   * @see    CartesianOrbit
71   * @author Mathieu Rom&eacute;ro
72   * @author Luc Maisonobe
73   * @author Guylaine Prat
74   * @author Fabien Maussion
75   * @author V&eacute;ronique Pommier-Maurussane
76   */
77  public class EquinoctialOrbit extends Orbit {
78  
79      /** Serializable UID. */
80      private static final long serialVersionUID = 20170414L;
81  
82      /** Factory for first time derivatives. */
83      private static final DSFactory FACTORY = new DSFactory(1, 1);
84  
85      /** Semi-major axis (m). */
86      private final double a;
87  
88      /** First component of the eccentricity vector. */
89      private final double ex;
90  
91      /** Second component of the eccentricity vector. */
92      private final double ey;
93  
94      /** First component of the inclination vector. */
95      private final double hx;
96  
97      /** Second component of the inclination vector. */
98      private final double hy;
99  
100     /** True longitude argument (rad). */
101     private final double lv;
102 
103     /** Semi-major axis derivative (m/s). */
104     private final double aDot;
105 
106     /** First component of the eccentricity vector derivative. */
107     private final double exDot;
108 
109     /** Second component of the eccentricity vector derivative. */
110     private final double eyDot;
111 
112     /** First component of the inclination vector derivative. */
113     private final double hxDot;
114 
115     /** Second component of the inclination vector derivative. */
116     private final double hyDot;
117 
118     /** True longitude argument derivative (rad/s). */
119     private final double lvDot;
120 
121     /** Partial Cartesian coordinates (position and velocity are valid, acceleration may be missing). */
122     private transient PVCoordinates partialPV;
123 
124     /** Creates a new instance.
125      * @param a  semi-major axis (m)
126      * @param ex e cos(ω + Ω), first component of eccentricity vector
127      * @param ey e sin(ω + Ω), second component of eccentricity vector
128      * @param hx tan(i/2) cos(Ω), first component of inclination vector
129      * @param hy tan(i/2) sin(Ω), second component of inclination vector
130      * @param l  (M or E or v) + ω + Ω, mean, eccentric or true longitude argument (rad)
131      * @param type type of longitude argument
132      * @param frame the frame in which the parameters are defined
133      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
134      * @param date date of the orbital parameters
135      * @param mu central attraction coefficient (m³/s²)
136      * @exception IllegalArgumentException if eccentricity is equal to 1 or larger or
137      * if frame is not a {@link Frame#isPseudoInertial pseudo-inertial frame}
138      */
139     public EquinoctialOrbit(final double a, final double ex, final double ey,
140                             final double hx, final double hy, final double l,
141                             final PositionAngle type,
142                             final Frame frame, final AbsoluteDate date, final double mu)
143         throws IllegalArgumentException {
144         this(a, ex, ey, hx, hy, l,
145              Double.NaN, Double.NaN, Double.NaN, Double.NaN, Double.NaN, Double.NaN,
146              type, frame, date, mu);
147     }
148 
149     /** Creates a new instance.
150      * @param a  semi-major axis (m)
151      * @param ex e cos(ω + Ω), first component of eccentricity vector
152      * @param ey e sin(ω + Ω), second component of eccentricity vector
153      * @param hx tan(i/2) cos(Ω), first component of inclination vector
154      * @param hy tan(i/2) sin(Ω), second component of inclination vector
155      * @param l  (M or E or v) + ω + Ω, mean, eccentric or true longitude argument (rad)
156      * @param aDot  semi-major axis derivative (m/s)
157      * @param exDot d(e cos(ω + Ω))/dt, first component of eccentricity vector derivative
158      * @param eyDot d(e sin(ω + Ω))/dt, second component of eccentricity vector derivative
159      * @param hxDot d(tan(i/2) cos(Ω))/dt, first component of inclination vector derivative
160      * @param hyDot d(tan(i/2) sin(Ω))/dt, second component of inclination vector derivative
161      * @param lDot  d(M or E or v) + ω + Ω)/dr, mean, eccentric or true longitude argument  derivative (rad/s)
162      * @param type type of longitude argument
163      * @param frame the frame in which the parameters are defined
164      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
165      * @param date date of the orbital parameters
166      * @param mu central attraction coefficient (m³/s²)
167      * @exception IllegalArgumentException if eccentricity is equal to 1 or larger or
168      * if frame is not a {@link Frame#isPseudoInertial pseudo-inertial frame}
169      */
170     public EquinoctialOrbit(final double a, final double ex, final double ey,
171                             final double hx, final double hy, final double l,
172                             final double aDot, final double exDot, final double eyDot,
173                             final double hxDot, final double hyDot, final double lDot,
174                             final PositionAngle type,
175                             final Frame frame, final AbsoluteDate date, final double mu)
176         throws IllegalArgumentException {
177         super(frame, date, mu);
178         if (ex * ex + ey * ey >= 1.0) {
179             throw new OrekitIllegalArgumentException(OrekitMessages.HYPERBOLIC_ORBIT_NOT_HANDLED_AS,
180                                                      getClass().getName());
181         }
182         this.a     = a;
183         this.aDot  = aDot;
184         this.ex    = ex;
185         this.exDot = exDot;
186         this.ey    = ey;
187         this.eyDot = eyDot;
188         this.hx    = hx;
189         this.hxDot = hxDot;
190         this.hy    = hy;
191         this.hyDot = hyDot;
192 
193         if (hasDerivatives()) {
194             final DerivativeStructure exDS = FACTORY.build(ex, exDot);
195             final DerivativeStructure eyDS = FACTORY.build(ey, eyDot);
196             final DerivativeStructure lDS  = FACTORY.build(l,  lDot);
197             final DerivativeStructure lvDS;
198             switch (type) {
199                 case MEAN :
200                     lvDS = FieldEquinoctialOrbit.eccentricToTrue(FieldEquinoctialOrbit.meanToEccentric(lDS, exDS, eyDS), exDS, eyDS);
201                     break;
202                 case ECCENTRIC :
203                     lvDS = FieldEquinoctialOrbit.eccentricToTrue(lDS, exDS, eyDS);
204                     break;
205                 case TRUE :
206                     lvDS = lDS;
207                     break;
208                 default : // this should never happen
209                     throw new OrekitInternalError(null);
210             }
211             this.lv    = lvDS.getValue();
212             this.lvDot = lvDS.getPartialDerivative(1);
213         } else {
214             switch (type) {
215                 case MEAN :
216                     this.lv = eccentricToTrue(meanToEccentric(l, ex, ey), ex, ey);
217                     break;
218                 case ECCENTRIC :
219                     this.lv = eccentricToTrue(l, ex, ey);
220                     break;
221                 case TRUE :
222                     this.lv = l;
223                     break;
224                 default : // this should never happen
225                     throw new OrekitInternalError(null);
226             }
227             this.lvDot = Double.NaN;
228         }
229 
230         this.partialPV = null;
231 
232     }
233 
234     /** Constructor from Cartesian parameters.
235      *
236      * <p> The acceleration provided in {@code pvCoordinates} is accessible using
237      * {@link #getPVCoordinates()} and {@link #getPVCoordinates(Frame)}. All other methods
238      * use {@code mu} and the position to compute the acceleration, including
239      * {@link #shiftedBy(double)} and {@link #getPVCoordinates(AbsoluteDate, Frame)}.
240      *
241      * @param pvCoordinates the position, velocity and acceleration
242      * @param frame the frame in which are defined the {@link PVCoordinates}
243      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
244      * @param mu central attraction coefficient (m³/s²)
245      * @exception IllegalArgumentException if eccentricity is equal to 1 or larger or
246      * if frame is not a {@link Frame#isPseudoInertial pseudo-inertial frame}
247      */
248     public EquinoctialOrbit(final TimeStampedPVCoordinates pvCoordinates,
249                             final Frame frame, final double mu)
250         throws IllegalArgumentException {
251         super(pvCoordinates, frame, mu);
252 
253         //  compute semi-major axis
254         final Vector3D pvP   = pvCoordinates.getPosition();
255         final Vector3D pvV   = pvCoordinates.getVelocity();
256         final Vector3D pvA   = pvCoordinates.getAcceleration();
257         final double r2      = pvP.getNormSq();
258         final double r       = FastMath.sqrt(r2);
259         final double V2      = pvV.getNormSq();
260         final double rV2OnMu = r * V2 / mu;
261 
262         if (rV2OnMu > 2) {
263             throw new OrekitIllegalArgumentException(OrekitMessages.HYPERBOLIC_ORBIT_NOT_HANDLED_AS,
264                                                      getClass().getName());
265         }
266 
267         // compute inclination vector
268         final Vector3D w = pvCoordinates.getMomentum().normalize();
269         final double d = 1.0 / (1 + w.getZ());
270         hx = -d * w.getY();
271         hy =  d * w.getX();
272 
273         // compute true longitude argument
274         final double cLv = (pvP.getX() - d * pvP.getZ() * w.getX()) / r;
275         final double sLv = (pvP.getY() - d * pvP.getZ() * w.getY()) / r;
276         lv = FastMath.atan2(sLv, cLv);
277 
278 
279         // compute semi-major axis
280         a = r / (2 - rV2OnMu);
281 
282         // compute eccentricity vector
283         final double eSE = Vector3D.dotProduct(pvP, pvV) / FastMath.sqrt(mu * a);
284         final double eCE = rV2OnMu - 1;
285         final double e2  = eCE * eCE + eSE * eSE;
286         final double f   = eCE - e2;
287         final double g   = FastMath.sqrt(1 - e2) * eSE;
288         ex = a * (f * cLv + g * sLv) / r;
289         ey = a * (f * sLv - g * cLv) / r;
290 
291         partialPV = pvCoordinates;
292 
293         if (hasNonKeplerianAcceleration(pvCoordinates, mu)) {
294             // we have a relevant acceleration, we can compute derivatives
295 
296             final double[][] jacobian = new double[6][6];
297             getJacobianWrtCartesian(PositionAngle.MEAN, jacobian);
298 
299             final Vector3D keplerianAcceleration    = new Vector3D(-mu / (r * r2), pvP);
300             final Vector3D nonKeplerianAcceleration = pvA.subtract(keplerianAcceleration);
301             final double   aX                       = nonKeplerianAcceleration.getX();
302             final double   aY                       = nonKeplerianAcceleration.getY();
303             final double   aZ                       = nonKeplerianAcceleration.getZ();
304             aDot  = jacobian[0][3] * aX + jacobian[0][4] * aY + jacobian[0][5] * aZ;
305             exDot = jacobian[1][3] * aX + jacobian[1][4] * aY + jacobian[1][5] * aZ;
306             eyDot = jacobian[2][3] * aX + jacobian[2][4] * aY + jacobian[2][5] * aZ;
307             hxDot = jacobian[3][3] * aX + jacobian[3][4] * aY + jacobian[3][5] * aZ;
308             hyDot = jacobian[4][3] * aX + jacobian[4][4] * aY + jacobian[4][5] * aZ;
309 
310             // in order to compute true anomaly derivative, we must compute
311             // mean anomaly derivative including Keplerian motion and convert to true anomaly
312             final double lMDot = getKeplerianMeanMotion() +
313                                  jacobian[5][3] * aX + jacobian[5][4] * aY + jacobian[5][5] * aZ;
314             final DerivativeStructure exDS = FACTORY.build(ex, exDot);
315             final DerivativeStructure eyDS = FACTORY.build(ey, eyDot);
316             final DerivativeStructure lMDS = FACTORY.build(getLM(), lMDot);
317             final DerivativeStructure lvDS = FieldEquinoctialOrbit.eccentricToTrue(FieldEquinoctialOrbit.meanToEccentric(lMDS, exDS, eyDS), exDS, eyDS);
318             lvDot = lvDS.getPartialDerivative(1);
319 
320         } else {
321             // acceleration is either almost zero or NaN,
322             // we assume acceleration was not known
323             // we don't set up derivatives
324             aDot  = Double.NaN;
325             exDot = Double.NaN;
326             eyDot = Double.NaN;
327             hxDot = Double.NaN;
328             hyDot = Double.NaN;
329             lvDot = Double.NaN;
330         }
331 
332     }
333 
334     /** Constructor from Cartesian parameters.
335      *
336      * <p> The acceleration provided in {@code pvCoordinates} is accessible using
337      * {@link #getPVCoordinates()} and {@link #getPVCoordinates(Frame)}. All other methods
338      * use {@code mu} and the position to compute the acceleration, including
339      * {@link #shiftedBy(double)} and {@link #getPVCoordinates(AbsoluteDate, Frame)}.
340      *
341      * @param pvCoordinates the position end velocity
342      * @param frame the frame in which are defined the {@link PVCoordinates}
343      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
344      * @param date date of the orbital parameters
345      * @param mu central attraction coefficient (m³/s²)
346      * @exception IllegalArgumentException if eccentricity is equal to 1 or larger or
347      * if frame is not a {@link Frame#isPseudoInertial pseudo-inertial frame}
348      */
349     public EquinoctialOrbit(final PVCoordinates pvCoordinates, final Frame frame,
350                             final AbsoluteDate date, final double mu)
351         throws IllegalArgumentException {
352         this(new TimeStampedPVCoordinates(date, pvCoordinates), frame, mu);
353     }
354 
355     /** Constructor from any kind of orbital parameters.
356      * @param op orbital parameters to copy
357      */
358     public EquinoctialOrbit(final Orbit op) {
359         super(op.getFrame(), op.getDate(), op.getMu());
360         a         = op.getA();
361         aDot      = op.getADot();
362         ex        = op.getEquinoctialEx();
363         exDot     = op.getEquinoctialExDot();
364         ey        = op.getEquinoctialEy();
365         eyDot     = op.getEquinoctialEyDot();
366         hx        = op.getHx();
367         hxDot     = op.getHxDot();
368         hy        = op.getHy();
369         hyDot     = op.getHyDot();
370         lv        = op.getLv();
371         lvDot     = op.getLvDot();
372         partialPV = null;
373     }
374 
375     /** {@inheritDoc} */
376     public OrbitType getType() {
377         return OrbitType.EQUINOCTIAL;
378     }
379 
380     /** {@inheritDoc} */
381     public double getA() {
382         return a;
383     }
384 
385     /** {@inheritDoc} */
386     public double getADot() {
387         return aDot;
388     }
389 
390     /** {@inheritDoc} */
391     public double getEquinoctialEx() {
392         return ex;
393     }
394 
395     /** {@inheritDoc} */
396     public double getEquinoctialExDot() {
397         return exDot;
398     }
399 
400     /** {@inheritDoc} */
401     public double getEquinoctialEy() {
402         return ey;
403     }
404 
405     /** {@inheritDoc} */
406     public double getEquinoctialEyDot() {
407         return eyDot;
408     }
409 
410     /** {@inheritDoc} */
411     public double getHx() {
412         return hx;
413     }
414 
415     /** {@inheritDoc} */
416     public double getHxDot() {
417         return hxDot;
418     }
419 
420     /** {@inheritDoc} */
421     public double getHy() {
422         return hy;
423     }
424 
425     /** {@inheritDoc} */
426     public double getHyDot() {
427         return hyDot;
428     }
429 
430     /** {@inheritDoc} */
431     public double getLv() {
432         return lv;
433     }
434 
435     /** {@inheritDoc} */
436     public double getLvDot() {
437         return lvDot;
438     }
439 
440     /** {@inheritDoc} */
441     public double getLE() {
442         return trueToEccentric(lv, ex, ey);
443     }
444 
445     /** {@inheritDoc} */
446     public double getLEDot() {
447         final DerivativeStructure lVDS = FACTORY.build(lv, lvDot);
448         final DerivativeStructure exDS = FACTORY.build(ex, exDot);
449         final DerivativeStructure eyDS = FACTORY.build(ey, eyDot);
450         final DerivativeStructure lEDS = FieldEquinoctialOrbit.trueToEccentric(lVDS, exDS, eyDS);
451         return lEDS.getPartialDerivative(1);
452     }
453 
454     /** {@inheritDoc} */
455     public double getLM() {
456         return eccentricToMean(trueToEccentric(lv, ex, ey), ex, ey);
457     }
458 
459     /** {@inheritDoc} */
460     public double getLMDot() {
461         final DerivativeStructure lVDS = FACTORY.build(lv, lvDot);
462         final DerivativeStructure exDS = FACTORY.build(ex, exDot);
463         final DerivativeStructure eyDS = FACTORY.build(ey, eyDot);
464         final DerivativeStructure lMDS = FieldEquinoctialOrbit.eccentricToMean(FieldEquinoctialOrbit.trueToEccentric(lVDS, exDS, eyDS), exDS, eyDS);
465         return lMDS.getPartialDerivative(1);
466     }
467 
468     /** Get the longitude argument.
469      * @param type type of the angle
470      * @return longitude argument (rad)
471      */
472     public double getL(final PositionAngle type) {
473         return (type == PositionAngle.MEAN) ? getLM() :
474                                               ((type == PositionAngle.ECCENTRIC) ? getLE() :
475                                                                                    getLv());
476     }
477 
478     /** Get the longitude argument derivative.
479      * @param type type of the angle
480      * @return longitude argument derivative (rad/s)
481      */
482     public double getLDot(final PositionAngle type) {
483         return (type == PositionAngle.MEAN) ? getLMDot() :
484                                               ((type == PositionAngle.ECCENTRIC) ? getLEDot() :
485                                                                                    getLvDot());
486     }
487 
488     /** Computes the true longitude argument from the eccentric longitude argument.
489      * @param lE = E + ω + Ω eccentric longitude argument (rad)
490      * @param ex first component of the eccentricity vector
491      * @param ey second component of the eccentricity vector
492      * @return the true longitude argument
493      */
494     public static double eccentricToTrue(final double lE, final double ex, final double ey) {
495         final double epsilon = FastMath.sqrt(1 - ex * ex - ey * ey);
496         final double cosLE   = FastMath.cos(lE);
497         final double sinLE   = FastMath.sin(lE);
498         final double num     = ex * sinLE - ey * cosLE;
499         final double den     = epsilon + 1 - ex * cosLE - ey * sinLE;
500         return lE + 2 * FastMath.atan(num / den);
501     }
502 
503     /** Computes the eccentric longitude argument from the true longitude argument.
504      * @param lv = v + ω + Ω true longitude argument (rad)
505      * @param ex first component of the eccentricity vector
506      * @param ey second component of the eccentricity vector
507      * @return the eccentric longitude argument
508      */
509     public static double trueToEccentric(final double lv, final double ex, final double ey) {
510         final double epsilon = FastMath.sqrt(1 - ex * ex - ey * ey);
511         final double cosLv   = FastMath.cos(lv);
512         final double sinLv   = FastMath.sin(lv);
513         final double num     = ey * cosLv - ex * sinLv;
514         final double den     = epsilon + 1 + ex * cosLv + ey * sinLv;
515         return lv + 2 * FastMath.atan(num / den);
516     }
517 
518     /** Computes the eccentric longitude argument from the mean longitude argument.
519      * @param lM = M + ω + Ω mean longitude argument (rad)
520      * @param ex first component of the eccentricity vector
521      * @param ey second component of the eccentricity vector
522      * @return the eccentric longitude argument
523      */
524     public static double meanToEccentric(final double lM, final double ex, final double ey) {
525         // Generalization of Kepler equation to equinoctial parameters
526         // with lE = PA + RAAN + E and
527         //      lM = PA + RAAN + M = lE - ex.sin(lE) + ey.cos(lE)
528         double lE = lM;
529         double shift = 0.0;
530         double lEmlM = 0.0;
531         double cosLE = FastMath.cos(lE);
532         double sinLE = FastMath.sin(lE);
533         int iter = 0;
534         do {
535             final double f2 = ex * sinLE - ey * cosLE;
536             final double f1 = 1.0 - ex * cosLE - ey * sinLE;
537             final double f0 = lEmlM - f2;
538 
539             final double f12 = 2.0 * f1;
540             shift = f0 * f12 / (f1 * f12 - f0 * f2);
541 
542             lEmlM -= shift;
543             lE     = lM + lEmlM;
544             cosLE  = FastMath.cos(lE);
545             sinLE  = FastMath.sin(lE);
546 
547         } while ((++iter < 50) && (FastMath.abs(shift) > 1.0e-12));
548 
549         return lE;
550 
551     }
552 
553     /** Computes the mean longitude argument from the eccentric longitude argument.
554      * @param lE = E + ω + Ω mean longitude argument (rad)
555      * @param ex first component of the eccentricity vector
556      * @param ey second component of the eccentricity vector
557      * @return the mean longitude argument
558      */
559     public static double eccentricToMean(final double lE, final double ex, final double ey) {
560         return lE - ex * FastMath.sin(lE) + ey * FastMath.cos(lE);
561     }
562 
563     /** {@inheritDoc} */
564     public double getE() {
565         return FastMath.sqrt(ex * ex + ey * ey);
566     }
567 
568     /** {@inheritDoc} */
569     public double getEDot() {
570         return (ex * exDot + ey * eyDot) / FastMath.sqrt(ex * ex + ey * ey);
571     }
572 
573     /** {@inheritDoc} */
574     public double getI() {
575         return 2 * FastMath.atan(FastMath.sqrt(hx * hx + hy * hy));
576     }
577 
578     /** {@inheritDoc} */
579     public double getIDot() {
580         final double h2 = hx * hx + hy * hy;
581         final double h  = FastMath.sqrt(h2);
582         return 2 * (hx * hxDot + hy * hyDot) / (h * (1 + h2));
583     }
584 
585     /** Compute position and velocity but not acceleration.
586      */
587     private void computePVWithoutA() {
588 
589         if (partialPV != null) {
590             // already computed
591             return;
592         }
593 
594         // get equinoctial parameters
595         final double lE = getLE();
596 
597         // inclination-related intermediate parameters
598         final double hx2   = hx * hx;
599         final double hy2   = hy * hy;
600         final double factH = 1. / (1 + hx2 + hy2);
601 
602         // reference axes defining the orbital plane
603         final double ux = (1 + hx2 - hy2) * factH;
604         final double uy =  2 * hx * hy * factH;
605         final double uz = -2 * hy * factH;
606 
607         final double vx = uy;
608         final double vy = (1 - hx2 + hy2) * factH;
609         final double vz =  2 * hx * factH;
610 
611         // eccentricity-related intermediate parameters
612         final double exey = ex * ey;
613         final double ex2  = ex * ex;
614         final double ey2  = ey * ey;
615         final double e2   = ex2 + ey2;
616         final double eta  = 1 + FastMath.sqrt(1 - e2);
617         final double beta = 1. / eta;
618 
619         // eccentric longitude argument
620         final double cLe    = FastMath.cos(lE);
621         final double sLe    = FastMath.sin(lE);
622         final double exCeyS = ex * cLe + ey * sLe;
623 
624         // coordinates of position and velocity in the orbital plane
625         final double x      = a * ((1 - beta * ey2) * cLe + beta * exey * sLe - ex);
626         final double y      = a * ((1 - beta * ex2) * sLe + beta * exey * cLe - ey);
627 
628         final double factor = FastMath.sqrt(getMu() / a) / (1 - exCeyS);
629         final double xdot   = factor * (-sLe + beta * ey * exCeyS);
630         final double ydot   = factor * ( cLe - beta * ex * exCeyS);
631 
632         final Vector3D position =
633                         new Vector3D(x * ux + y * vx, x * uy + y * vy, x * uz + y * vz);
634         final Vector3D velocity =
635                         new Vector3D(xdot * ux + ydot * vx, xdot * uy + ydot * vy, xdot * uz + ydot * vz);
636         partialPV = new PVCoordinates(position, velocity);
637 
638     }
639 
640     /** Compute non-Keplerian part of the acceleration from first time derivatives.
641      * <p>
642      * This method should be called only when {@link #hasDerivatives()} returns true.
643      * </p>
644      * @return non-Keplerian part of the acceleration
645      */
646     private Vector3D nonKeplerianAcceleration() {
647 
648         final double[][] dCdP = new double[6][6];
649         getJacobianWrtParameters(PositionAngle.MEAN, dCdP);
650 
651         final double nonKeplerianMeanMotion = getLMDot() - getKeplerianMeanMotion();
652         final double nonKeplerianAx = dCdP[3][0] * aDot  + dCdP[3][1] * exDot + dCdP[3][2] * eyDot +
653                                       dCdP[3][3] * hxDot + dCdP[3][4] * hyDot + dCdP[3][5] * nonKeplerianMeanMotion;
654         final double nonKeplerianAy = dCdP[4][0] * aDot  + dCdP[4][1] * exDot + dCdP[4][2] * eyDot +
655                                       dCdP[4][3] * hxDot + dCdP[4][4] * hyDot + dCdP[4][5] * nonKeplerianMeanMotion;
656         final double nonKeplerianAz = dCdP[5][0] * aDot  + dCdP[5][1] * exDot + dCdP[5][2] * eyDot +
657                                       dCdP[5][3] * hxDot + dCdP[5][4] * hyDot + dCdP[5][5] * nonKeplerianMeanMotion;
658 
659         return new Vector3D(nonKeplerianAx, nonKeplerianAy, nonKeplerianAz);
660 
661     }
662 
663     /** {@inheritDoc} */
664     protected TimeStampedPVCoordinates initPVCoordinates() {
665 
666         // position and velocity
667         computePVWithoutA();
668 
669         // acceleration
670         final double r2 = partialPV.getPosition().getNormSq();
671         final Vector3D keplerianAcceleration = new Vector3D(-getMu() / (r2 * FastMath.sqrt(r2)), partialPV.getPosition());
672         final Vector3D acceleration = hasDerivatives() ?
673                                       keplerianAcceleration.add(nonKeplerianAcceleration()) :
674                                       keplerianAcceleration;
675 
676         return new TimeStampedPVCoordinates(getDate(), partialPV.getPosition(), partialPV.getVelocity(), acceleration);
677 
678     }
679 
680     /** {@inheritDoc} */
681     public EquinoctialOrbit shiftedBy(final double dt) {
682 
683         // use Keplerian-only motion
684         final EquinoctialOrbitOrbit">EquinoctialOrbit keplerianShifted = new EquinoctialOrbit(a, ex, ey, hx, hy,
685                                                                        getLM() + getKeplerianMeanMotion() * dt,
686                                                                        PositionAngle.MEAN, getFrame(),
687                                                                        getDate().shiftedBy(dt), getMu());
688 
689         if (hasDerivatives()) {
690 
691             // extract non-Keplerian acceleration from first time derivatives
692             final Vector3D nonKeplerianAcceleration = nonKeplerianAcceleration();
693 
694             // add quadratic effect of non-Keplerian acceleration to Keplerian-only shift
695             keplerianShifted.computePVWithoutA();
696             final Vector3D fixedP   = new Vector3D(1, keplerianShifted.partialPV.getPosition(),
697                                                    0.5 * dt * dt, nonKeplerianAcceleration);
698             final double   fixedR2 = fixedP.getNormSq();
699             final double   fixedR  = FastMath.sqrt(fixedR2);
700             final Vector3D fixedV  = new Vector3D(1, keplerianShifted.partialPV.getVelocity(),
701                                                   dt, nonKeplerianAcceleration);
702             final Vector3D fixedA  = new Vector3D(-getMu() / (fixedR2 * fixedR), keplerianShifted.partialPV.getPosition(),
703                                                   1, nonKeplerianAcceleration);
704 
705             // build a new orbit, taking non-Keplerian acceleration into account
706             return new EquinoctialOrbit(new TimeStampedPVCoordinates(keplerianShifted.getDate(),
707                                                                      fixedP, fixedV, fixedA),
708                                         keplerianShifted.getFrame(), keplerianShifted.getMu());
709 
710         } else {
711             // Keplerian-only motion is all we can do
712             return keplerianShifted;
713         }
714 
715     }
716 
717     /** {@inheritDoc}
718      * <p>
719      * The interpolated instance is created by polynomial Hermite interpolation
720      * on equinoctial elements, without derivatives (which means the interpolation
721      * falls back to Lagrange interpolation only).
722      * </p>
723      * <p>
724      * As this implementation of interpolation is polynomial, it should be used only
725      * with small samples (about 10-20 points) in order to avoid <a
726      * href="http://en.wikipedia.org/wiki/Runge%27s_phenomenon">Runge's phenomenon</a>
727      * and numerical problems (including NaN appearing).
728      * </p>
729      * <p>
730      * If orbit interpolation on large samples is needed, using the {@link
731      * org.orekit.propagation.analytical.Ephemeris} class is a better way than using this
732      * low-level interpolation. The Ephemeris class automatically handles selection of
733      * a neighboring sub-sample with a predefined number of point from a large global sample
734      * in a thread-safe way.
735      * </p>
736      */
737     public EquinoctialOrbit interpolate(final AbsoluteDate date, final Stream<Orbit> sample) {
738 
739         // first pass to check if derivatives are available throughout the sample
740         final List<Orbit> list = sample.collect(Collectors.toList());
741         boolean useDerivatives = true;
742         for (final Orbit orbit : list) {
743             useDerivatives = useDerivatives && orbit.hasDerivatives();
744         }
745 
746         // set up an interpolator
747         final HermiteInterpolator interpolator = new HermiteInterpolator();
748 
749         // second pass to feed interpolator
750         AbsoluteDate previousDate = null;
751         double       previousLm   = Double.NaN;
752         for (final Orbit orbit : list) {
753             final EquinoctialOrbit/org/orekit/orbits/EquinoctialOrbit.html#EquinoctialOrbit">EquinoctialOrbit equi = (EquinoctialOrbit) OrbitType.EQUINOCTIAL.convertType(orbit);
754             final double continuousLm;
755             if (previousDate == null) {
756                 continuousLm = equi.getLM();
757             } else {
758                 final double dt       = equi.getDate().durationFrom(previousDate);
759                 final double keplerLm = previousLm + equi.getKeplerianMeanMotion() * dt;
760                 continuousLm = MathUtils.normalizeAngle(equi.getLM(), keplerLm);
761             }
762             previousDate = equi.getDate();
763             previousLm   = continuousLm;
764             if (useDerivatives) {
765                 interpolator.addSamplePoint(equi.getDate().durationFrom(date),
766                                             new double[] {
767                                                 equi.getA(),
768                                                 equi.getEquinoctialEx(),
769                                                 equi.getEquinoctialEy(),
770                                                 equi.getHx(),
771                                                 equi.getHy(),
772                                                 continuousLm
773                                             },
774                                             new double[] {
775                                                 equi.getADot(),
776                                                 equi.getEquinoctialExDot(),
777                                                 equi.getEquinoctialEyDot(),
778                                                 equi.getHxDot(),
779                                                 equi.getHyDot(),
780                                                 equi.getLMDot()
781                                             });
782             } else {
783                 interpolator.addSamplePoint(equi.getDate().durationFrom(date),
784                                             new double[] {
785                                                 equi.getA(),
786                                                 equi.getEquinoctialEx(),
787                                                 equi.getEquinoctialEy(),
788                                                 equi.getHx(),
789                                                 equi.getHy(),
790                                                 continuousLm
791                                             });
792             }
793         }
794 
795         // interpolate
796         final double[][] interpolated = interpolator.derivatives(0.0, 1);
797 
798         // build a new interpolated instance
799         return new EquinoctialOrbit(interpolated[0][0], interpolated[0][1], interpolated[0][2],
800                                     interpolated[0][3], interpolated[0][4], interpolated[0][5],
801                                     interpolated[1][0], interpolated[1][1], interpolated[1][2],
802                                     interpolated[1][3], interpolated[1][4], interpolated[1][5],
803                                     PositionAngle.MEAN, getFrame(), date, getMu());
804 
805     }
806 
807     /** {@inheritDoc} */
808     protected double[][] computeJacobianMeanWrtCartesian() {
809 
810         final double[][] jacobian = new double[6][6];
811 
812         // compute various intermediate parameters
813         computePVWithoutA();
814         final Vector3D position = partialPV.getPosition();
815         final Vector3D velocity = partialPV.getVelocity();
816         final double r2         = position.getNormSq();
817         final double r          = FastMath.sqrt(r2);
818         final double r3         = r * r2;
819 
820         final double mu         = getMu();
821         final double sqrtMuA    = FastMath.sqrt(a * mu);
822         final double a2         = a * a;
823 
824         final double e2         = ex * ex + ey * ey;
825         final double oMe2       = 1 - e2;
826         final double epsilon    = FastMath.sqrt(oMe2);
827         final double beta       = 1 / (1 + epsilon);
828         final double ratio      = epsilon * beta;
829 
830         final double hx2        = hx * hx;
831         final double hy2        = hy * hy;
832         final double hxhy       = hx * hy;
833 
834         // precomputing equinoctial frame unit vectors (f, g, w)
835         final Vector3D f  = new Vector3D(1 - hy2 + hx2, 2 * hxhy, -2 * hy).normalize();
836         final Vector3D g  = new Vector3D(2 * hxhy, 1 + hy2 - hx2, 2 * hx).normalize();
837         final Vector3D w  = Vector3D.crossProduct(position, velocity).normalize();
838 
839         // coordinates of the spacecraft in the equinoctial frame
840         final double x    = Vector3D.dotProduct(position, f);
841         final double y    = Vector3D.dotProduct(position, g);
842         final double xDot = Vector3D.dotProduct(velocity, f);
843         final double yDot = Vector3D.dotProduct(velocity, g);
844 
845         // drDot / dEx = dXDot / dEx * f + dYDot / dEx * g
846         final double c1 = a / (sqrtMuA * epsilon);
847         final double c2 = a * sqrtMuA * beta / r3;
848         final double c3 = sqrtMuA / (r3 * epsilon);
849         final Vector3D drDotSdEx = new Vector3D( c1 * xDot * yDot - c2 * ey * x - c3 * x * y, f,
850                                                 -c1 * xDot * xDot - c2 * ey * y + c3 * x * x, g);
851 
852         // drDot / dEy = dXDot / dEy * f + dYDot / dEy * g
853         final Vector3D drDotSdEy = new Vector3D( c1 * yDot * yDot + c2 * ex * x - c3 * y * y, f,
854                                                 -c1 * xDot * yDot + c2 * ex * y + c3 * x * y, g);
855 
856         // da
857         final Vector3D vectorAR = new Vector3D(2 * a2 / r3, position);
858         final Vector3D vectorARDot = new Vector3D(2 * a2 / mu, velocity);
859         fillHalfRow(1, vectorAR,    jacobian[0], 0);
860         fillHalfRow(1, vectorARDot, jacobian[0], 3);
861 
862         // dEx
863         final double d1 = -a * ratio / r3;
864         final double d2 = (hy * xDot - hx * yDot) / (sqrtMuA * epsilon);
865         final double d3 = (hx * y - hy * x) / sqrtMuA;
866         final Vector3D vectorExRDot =
867             new Vector3D((2 * x * yDot - xDot * y) / mu, g, -y * yDot / mu, f, -ey * d3 / epsilon, w);
868         fillHalfRow(ex * d1, position, -ey * d2, w, epsilon / sqrtMuA, drDotSdEy, jacobian[1], 0);
869         fillHalfRow(1, vectorExRDot, jacobian[1], 3);
870 
871         // dEy
872         final Vector3D vectorEyRDot =
873             new Vector3D((2 * xDot * y - x * yDot) / mu, f, -x * xDot / mu, g, ex * d3 / epsilon, w);
874         fillHalfRow(ey * d1, position, ex * d2, w, -epsilon / sqrtMuA, drDotSdEx, jacobian[2], 0);
875         fillHalfRow(1, vectorEyRDot, jacobian[2], 3);
876 
877         // dHx
878         final double h = (1 + hx2 + hy2) / (2 * sqrtMuA * epsilon);
879         fillHalfRow(-h * xDot, w, jacobian[3], 0);
880         fillHalfRow( h * x,    w, jacobian[3], 3);
881 
882         // dHy
883         fillHalfRow(-h * yDot, w, jacobian[4], 0);
884         fillHalfRow( h * y,    w, jacobian[4], 3);
885 
886         // dLambdaM
887         final double l = -ratio / sqrtMuA;
888         fillHalfRow(-1 / sqrtMuA, velocity, d2, w, l * ex, drDotSdEx, l * ey, drDotSdEy, jacobian[5], 0);
889         fillHalfRow(-2 / sqrtMuA, position, ex * beta, vectorEyRDot, -ey * beta, vectorExRDot, d3, w, jacobian[5], 3);
890 
891         return jacobian;
892 
893     }
894 
895     /** {@inheritDoc} */
896     protected double[][] computeJacobianEccentricWrtCartesian() {
897 
898         // start by computing the Jacobian with mean angle
899         final double[][] jacobian = computeJacobianMeanWrtCartesian();
900 
901         // Differentiating the Kepler equation lM = lE - ex sin lE + ey cos lE leads to:
902         // dlM = (1 - ex cos lE - ey sin lE) dE - sin lE dex + cos lE dey
903         // which is inverted and rewritten as:
904         // dlE = a/r dlM + sin lE a/r dex - cos lE a/r dey
905         final double le    = getLE();
906         final double cosLe = FastMath.cos(le);
907         final double sinLe = FastMath.sin(le);
908         final double aOr   = 1 / (1 - ex * cosLe - ey * sinLe);
909 
910         // update longitude row
911         final double[] rowEx = jacobian[1];
912         final double[] rowEy = jacobian[2];
913         final double[] rowL  = jacobian[5];
914         for (int j = 0; j < 6; ++j) {
915             rowL[j] = aOr * (rowL[j] + sinLe * rowEx[j] - cosLe * rowEy[j]);
916         }
917 
918         return jacobian;
919 
920     }
921 
922     /** {@inheritDoc} */
923     protected double[][] computeJacobianTrueWrtCartesian() {
924 
925         // start by computing the Jacobian with eccentric angle
926         final double[][] jacobian = computeJacobianEccentricWrtCartesian();
927 
928         // Differentiating the eccentric longitude equation
929         // tan((lV - lE)/2) = [ex sin lE - ey cos lE] / [sqrt(1-ex^2-ey^2) + 1 - ex cos lE - ey sin lE]
930         // leads to
931         // cT (dlV - dlE) = cE dlE + cX dex + cY dey
932         // with
933         // cT = [d^2 + (ex sin lE - ey cos lE)^2] / 2
934         // d  = 1 + sqrt(1-ex^2-ey^2) - ex cos lE - ey sin lE
935         // cE = (ex cos lE + ey sin lE) (sqrt(1-ex^2-ey^2) + 1) - ex^2 - ey^2
936         // cX =  sin lE (sqrt(1-ex^2-ey^2) + 1) - ey + ex (ex sin lE - ey cos lE) / sqrt(1-ex^2-ey^2)
937         // cY = -cos lE (sqrt(1-ex^2-ey^2) + 1) + ex + ey (ex sin lE - ey cos lE) / sqrt(1-ex^2-ey^2)
938         // which can be solved to find the differential of the true longitude
939         // dlV = (cT + cE) / cT dlE + cX / cT deX + cY / cT deX
940         final double le        = getLE();
941         final double cosLe     = FastMath.cos(le);
942         final double sinLe     = FastMath.sin(le);
943         final double eSinE     = ex * sinLe - ey * cosLe;
944         final double ecosE     = ex * cosLe + ey * sinLe;
945         final double e2        = ex * ex + ey * ey;
946         final double epsilon   = FastMath.sqrt(1 - e2);
947         final double onePeps   = 1 + epsilon;
948         final double d         = onePeps - ecosE;
949         final double cT        = (d * d + eSinE * eSinE) / 2;
950         final double cE        = ecosE * onePeps - e2;
951         final double cX        = ex * eSinE / epsilon - ey + sinLe * onePeps;
952         final double cY        = ey * eSinE / epsilon + ex - cosLe * onePeps;
953         final double factorLe  = (cT + cE) / cT;
954         final double factorEx  = cX / cT;
955         final double factorEy  = cY / cT;
956 
957         // update longitude row
958         final double[] rowEx = jacobian[1];
959         final double[] rowEy = jacobian[2];
960         final double[] rowL = jacobian[5];
961         for (int j = 0; j < 6; ++j) {
962             rowL[j] = factorLe * rowL[j] + factorEx * rowEx[j] + factorEy * rowEy[j];
963         }
964 
965         return jacobian;
966 
967     }
968 
969     /** {@inheritDoc} */
970     public void addKeplerContribution(final PositionAngle type, final double gm,
971                                       final double[] pDot) {
972         final double oMe2;
973         final double ksi;
974         final double n = FastMath.sqrt(gm / a) / a;
975         switch (type) {
976             case MEAN :
977                 pDot[5] += n;
978                 break;
979             case ECCENTRIC :
980                 oMe2 = 1 - ex * ex - ey * ey;
981                 ksi  = 1 + ex * FastMath.cos(lv) + ey * FastMath.sin(lv);
982                 pDot[5] += n * ksi / oMe2;
983                 break;
984             case TRUE :
985                 oMe2 = 1 - ex * ex - ey * ey;
986                 ksi  = 1 + ex * FastMath.cos(lv) + ey * FastMath.sin(lv);
987                 pDot[5] += n * ksi * ksi / (oMe2 * FastMath.sqrt(oMe2));
988                 break;
989             default :
990                 throw new OrekitInternalError(null);
991         }
992     }
993 
994     /**  Returns a string representation of this equinoctial parameters object.
995      * @return a string representation of this object
996      */
997     public String toString() {
998         return new StringBuffer().append("equinoctial parameters: ").append('{').
999                                   append("a: ").append(a).
1000                                   append("; ex: ").append(ex).append("; ey: ").append(ey).
1001                                   append("; hx: ").append(hx).append("; hy: ").append(hy).
1002                                   append("; lv: ").append(FastMath.toDegrees(lv)).
1003                                   append(";}").toString();
1004     }
1005 
1006     /** Replace the instance with a data transfer object for serialization.
1007      * @return data transfer object that will be serialized
1008      */
1009     private Object writeReplace() {
1010         return new DTO(this);
1011     }
1012 
1013     /** Internal class used only for serialization. */
1014     private static class DTO implements Serializable {
1015 
1016         /** Serializable UID. */
1017         private static final long serialVersionUID = 20170414L;
1018 
1019         /** Double values. */
1020         private double[] d;
1021 
1022         /** Frame in which are defined the orbital parameters. */
1023         private final Frame frame;
1024 
1025         /** Simple constructor.
1026          * @param orbit instance to serialize
1027          */
1028         private DTO(final EquinoctialOrbit orbit) {
1029 
1030             final TimeStampedPVCoordinates pv = orbit.getPVCoordinates();
1031 
1032             // decompose date
1033             final double epoch  = FastMath.floor(pv.getDate().durationFrom(AbsoluteDate.J2000_EPOCH));
1034             final double offset = pv.getDate().durationFrom(AbsoluteDate.J2000_EPOCH.shiftedBy(epoch));
1035 
1036             if (orbit.hasDerivatives()) {
1037                 // we have derivatives
1038                 this.d = new double[] {
1039                     epoch, offset, orbit.getMu(),
1040                     orbit.a, orbit.ex, orbit.ey,
1041                     orbit.hx, orbit.hy, orbit.lv,
1042                     orbit.aDot, orbit.exDot, orbit.eyDot,
1043                     orbit.hxDot, orbit.hyDot, orbit.lvDot
1044                 };
1045             } else {
1046                 // we don't have derivatives
1047                 this.d = new double[] {
1048                     epoch, offset, orbit.getMu(),
1049                     orbit.a, orbit.ex, orbit.ey,
1050                     orbit.hx, orbit.hy, orbit.lv
1051                 };
1052             }
1053 
1054             this.frame = orbit.getFrame();
1055 
1056         }
1057 
1058         /** Replace the deserialized data transfer object with a {@link EquinoctialOrbit}.
1059          * @return replacement {@link EquinoctialOrbit}
1060          */
1061         private Object readResolve() {
1062             if (d.length >= 15) {
1063                 // we have derivatives
1064                 return new EquinoctialOrbit(d[ 3], d[ 4], d[ 5], d[ 6], d[ 7], d[ 8],
1065                                             d[ 9], d[10], d[11], d[12], d[13], d[14],
1066                                             PositionAngle.TRUE,
1067                                             frame, AbsoluteDate.J2000_EPOCH.shiftedBy(d[0]).shiftedBy(d[1]),
1068                                             d[2]);
1069             } else {
1070                 // we don't have derivatives
1071                 return new EquinoctialOrbit(d[ 3], d[ 4], d[ 5], d[ 6], d[ 7], d[ 8], PositionAngle.TRUE,
1072                                             frame, AbsoluteDate.J2000_EPOCH.shiftedBy(d[0]).shiftedBy(d[1]),
1073                                             d[2]);
1074             }
1075         }
1076 
1077     }
1078 
1079 }