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