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