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