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