1   /* Copyright 2002-2016 CS Systèmes d'Information
2    * Licensed to CS Systèmes d'Information (CS) under one or more
3    * contributor license agreements.  See the NOTICE file distributed with
4    * this work for additional information regarding copyright ownership.
5    * CS licenses this file to You under the Apache License, Version 2.0
6    * (the "License"); you may not use this file except in compliance with
7    * the License.  You may obtain a copy of the License at
8    *
9    *   http://www.apache.org/licenses/LICENSE-2.0
10   *
11   * Unless required by applicable law or agreed to in writing, software
12   * distributed under the License is distributed on an "AS IS" BASIS,
13   * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14   * See the License for the specific language governing permissions and
15   * limitations under the License.
16   */
17  package org.orekit.orbits;
18  
19  import java.io.Serializable;
20  import java.util.Collection;
21  
22  import org.apache.commons.math3.analysis.interpolation.HermiteInterpolator;
23  import org.apache.commons.math3.geometry.euclidean.threed.Vector3D;
24  import org.apache.commons.math3.util.FastMath;
25  import org.apache.commons.math3.util.MathUtils;
26  import org.orekit.errors.OrekitIllegalArgumentException;
27  import org.orekit.errors.OrekitInternalError;
28  import org.orekit.errors.OrekitMessages;
29  import org.orekit.frames.Frame;
30  import org.orekit.time.AbsoluteDate;
31  import org.orekit.utils.PVCoordinates;
32  import org.orekit.utils.TimeStampedPVCoordinates;
33  
34  
35  /**
36   * This class handles equinoctial orbital parameters, which can support both
37   * circular and equatorial orbits.
38   * <p>
39   * The parameters used internally are the equinoctial elements which can be
40   * related to keplerian elements as follows:
41   *   <pre>
42   *     a
43   *     ex = e cos(ω + Ω)
44   *     ey = e sin(ω + Ω)
45   *     hx = tan(i/2) cos(Ω)
46   *     hy = tan(i/2) sin(Ω)
47   *     lv = v + ω + Ω
48   *   </pre>
49   * where ω stands for the Perigee Argument and Ω stands for the
50   * Right Ascension of the Ascending Node.
51   * </p>
52   * <p>
53   * The conversion equations from and to keplerian elements given above hold only
54   * when both sides are unambiguously defined, i.e. when orbit is neither equatorial
55   * nor circular. When orbit is either equatorial or circular, the equinoctial
56   * parameters are still unambiguously defined whereas some keplerian elements
57   * (more precisely ω and Ω) become ambiguous. For this reason, equinoctial
58   * parameters are the recommended way to represent orbits.
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 {
74  
75      /** Serializable UID. */
76      private static final long serialVersionUID = 20141228L;
77  
78      /** Semi-major axis (m). */
79      private final double a;
80  
81      /** First component of the eccentricity vector. */
82      private final double ex;
83  
84      /** Second component of the eccentricity vector. */
85      private final double ey;
86  
87      /** First component of the inclination vector. */
88      private final double hx;
89  
90      /** Second component of the inclination vector. */
91      private final double hy;
92  
93      /** True longitude argument (rad). */
94      private final double lv;
95  
96      /** Creates a new instance.
97       * @param a  semi-major axis (m)
98       * @param ex e cos(ω + Ω), first component of eccentricity vector
99       * @param ey e sin(ω + Ω), second component of eccentricity vector
100      * @param hx tan(i/2) cos(Ω), first component of inclination vector
101      * @param hy tan(i/2) sin(Ω), second component of inclination vector
102      * @param l  (M or E or v) + ω + Ω, mean, eccentric or true longitude argument (rad)
103      * @param type type of longitude argument
104      * @param frame the frame in which the parameters are defined
105      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
106      * @param date date of the orbital parameters
107      * @param mu central attraction coefficient (m³/s²)
108      * @exception IllegalArgumentException if eccentricity is equal to 1 or larger or
109      * if frame is not a {@link Frame#isPseudoInertial pseudo-inertial frame}
110      */
111     public EquinoctialOrbit(final double a, final double ex, final double ey,
112                             final double hx, final double hy,
113                             final double l, final PositionAngle type,
114                             final Frame frame, final AbsoluteDate date, final double mu)
115         throws IllegalArgumentException {
116         super(frame, date, mu);
117         if (ex * ex + ey * ey >= 1.0) {
118             throw new OrekitIllegalArgumentException(OrekitMessages.HYPERBOLIC_ORBIT_NOT_HANDLED_AS,
119                                                      getClass().getName());
120         }
121         this.a  =  a;
122         this.ex = ex;
123         this.ey = ey;
124         this.hx = hx;
125         this.hy = hy;
126 
127         switch (type) {
128             case MEAN :
129                 this.lv = eccentricToTrue(meanToEccentric(l));
130                 break;
131             case ECCENTRIC :
132                 this.lv = eccentricToTrue(l);
133                 break;
134             case TRUE :
135                 this.lv = l;
136                 break;
137             default :
138                 throw new OrekitInternalError(null);
139         }
140 
141     }
142 
143     /** Constructor from cartesian parameters.
144      *
145      * <p> The acceleration provided in {@code pvCoordinates} is accessible using
146      * {@link #getPVCoordinates()} and {@link #getPVCoordinates(Frame)}. All other methods
147      * use {@code mu} and the position to compute the acceleration, including
148      * {@link #shiftedBy(double)} and {@link #getPVCoordinates(AbsoluteDate, Frame)}.
149      *
150      * @param pvCoordinates the position, velocity and acceleration
151      * @param frame the frame in which are defined the {@link PVCoordinates}
152      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
153      * @param mu central attraction coefficient (m³/s²)
154      * @exception IllegalArgumentException if eccentricity is equal to 1 or larger or
155      * if frame is not a {@link Frame#isPseudoInertial pseudo-inertial frame}
156      */
157     public EquinoctialOrbit(final TimeStampedPVCoordinates pvCoordinates,
158                             final Frame frame, final double mu)
159         throws IllegalArgumentException {
160         super(pvCoordinates, frame, mu);
161 
162         //  compute semi-major axis
163         final Vector3D pvP = pvCoordinates.getPosition();
164         final Vector3D pvV = pvCoordinates.getVelocity();
165         final double r = pvP.getNorm();
166         final double V2 = pvV.getNormSq();
167         final double rV2OnMu = r * V2 / mu;
168 
169         if (rV2OnMu > 2) {
170             throw new OrekitIllegalArgumentException(OrekitMessages.HYPERBOLIC_ORBIT_NOT_HANDLED_AS,
171                                                      getClass().getName());
172         }
173 
174         // compute inclination vector
175         final Vector3D w = pvCoordinates.getMomentum().normalize();
176         final double d = 1.0 / (1 + w.getZ());
177         hx = -d * w.getY();
178         hy =  d * w.getX();
179 
180         // compute true longitude argument
181         final double cLv = (pvP.getX() - d * pvP.getZ() * w.getX()) / r;
182         final double sLv = (pvP.getY() - d * pvP.getZ() * w.getY()) / r;
183         lv = FastMath.atan2(sLv, cLv);
184 
185 
186         // compute semi-major axis
187         a = r / (2 - rV2OnMu);
188 
189         // compute eccentricity vector
190         final double eSE = Vector3D.dotProduct(pvP, pvV) / FastMath.sqrt(mu * a);
191         final double eCE = rV2OnMu - 1;
192         final double e2  = eCE * eCE + eSE * eSE;
193         final double f   = eCE - e2;
194         final double g   = FastMath.sqrt(1 - e2) * eSE;
195         ex = a * (f * cLv + g * sLv) / r;
196         ey = a * (f * sLv - g * cLv) / r;
197 
198     }
199 
200     /** Constructor from cartesian parameters.
201      *
202      * <p> The acceleration provided in {@code pvCoordinates} is accessible using
203      * {@link #getPVCoordinates()} and {@link #getPVCoordinates(Frame)}. All other methods
204      * use {@code mu} and the position to compute the acceleration, including
205      * {@link #shiftedBy(double)} and {@link #getPVCoordinates(AbsoluteDate, Frame)}.
206      *
207      * @param pvCoordinates the position end velocity
208      * @param frame the frame in which are defined the {@link PVCoordinates}
209      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
210      * @param date date of the orbital parameters
211      * @param mu central attraction coefficient (m³/s²)
212      * @exception IllegalArgumentException if eccentricity is equal to 1 or larger or
213      * if frame is not a {@link Frame#isPseudoInertial pseudo-inertial frame}
214      */
215     public EquinoctialOrbit(final PVCoordinates pvCoordinates, final Frame frame,
216                             final AbsoluteDate date, final double mu)
217         throws IllegalArgumentException {
218         this(new TimeStampedPVCoordinates(date, pvCoordinates), frame, mu);
219     }
220 
221     /** Constructor from any kind of orbital parameters.
222      * @param op orbital parameters to copy
223      */
224     public EquinoctialOrbit(final Orbit op) {
225         super(op.getFrame(), op.getDate(), op.getMu());
226         a  = op.getA();
227         ex = op.getEquinoctialEx();
228         ey = op.getEquinoctialEy();
229         hx = op.getHx();
230         hy = op.getHy();
231         lv = op.getLv();
232     }
233 
234     /** {@inheritDoc} */
235     public OrbitType getType() {
236         return OrbitType.EQUINOCTIAL;
237     }
238 
239     /** {@inheritDoc} */
240     public double getA() {
241         return a;
242     }
243 
244     /** {@inheritDoc} */
245     public double getEquinoctialEx() {
246         return ex;
247     }
248 
249     /** {@inheritDoc} */
250     public double getEquinoctialEy() {
251         return ey;
252     }
253 
254     /** {@inheritDoc} */
255     public double getHx() {
256         return hx;
257     }
258 
259     /** {@inheritDoc} */
260     public double getHy() {
261         return hy;
262     }
263 
264     /** Get the longitude argument.
265      * @param type type of the angle
266      * @return longitude argument (rad)
267      */
268     public double getL(final PositionAngle type) {
269         return (type == PositionAngle.MEAN) ? getLM() :
270                                               ((type == PositionAngle.ECCENTRIC) ? getLE() :
271                                                                                    getLv());
272     }
273 
274     /** {@inheritDoc} */
275     public double getLv() {
276         return lv;
277     }
278 
279     /** {@inheritDoc} */
280     public double getLE() {
281         final double epsilon = FastMath.sqrt(1 - ex * ex - ey * ey);
282         final double cosLv   = FastMath.cos(lv);
283         final double sinLv   = FastMath.sin(lv);
284         final double num     = ey * cosLv - ex * sinLv;
285         final double den     = epsilon + 1 + ex * cosLv + ey * sinLv;
286         return lv + 2 * FastMath.atan(num / den);
287     }
288 
289     /** Computes the true longitude argument from the eccentric longitude argument.
290      * @param lE = E + ω + Ω eccentric longitude argument (rad)
291      * @return the true longitude argument
292      */
293     private double eccentricToTrue(final double lE) {
294         final double epsilon = FastMath.sqrt(1 - ex * ex - ey * ey);
295         final double cosLE   = FastMath.cos(lE);
296         final double sinLE   = FastMath.sin(lE);
297         final double num     = ex * sinLE - ey * cosLE;
298         final double den     = epsilon + 1 - ex * cosLE - ey * sinLE;
299         return lE + 2 * FastMath.atan(num / den);
300     }
301 
302     /** {@inheritDoc} */
303     public double getLM() {
304         final double lE = getLE();
305         return lE - ex * FastMath.sin(lE) + ey * FastMath.cos(lE);
306     }
307 
308     /** Computes the eccentric longitude argument from the mean longitude argument.
309      * @param lM = M + ω + Ω mean longitude argument (rad)
310      * @return the eccentric longitude argument
311      */
312     private double meanToEccentric(final double lM) {
313         // Generalization of Kepler equation to equinoctial parameters
314         // with lE = PA + RAAN + E and
315         //      lM = PA + RAAN + M = lE - ex.sin(lE) + ey.cos(lE)
316         double lE = lM;
317         double shift = 0.0;
318         double lEmlM = 0.0;
319         double cosLE = FastMath.cos(lE);
320         double sinLE = FastMath.sin(lE);
321         int iter = 0;
322         do {
323             final double f2 = ex * sinLE - ey * cosLE;
324             final double f1 = 1.0 - ex * cosLE - ey * sinLE;
325             final double f0 = lEmlM - f2;
326 
327             final double f12 = 2.0 * f1;
328             shift = f0 * f12 / (f1 * f12 - f0 * f2);
329 
330             lEmlM -= shift;
331             lE     = lM + lEmlM;
332             cosLE  = FastMath.cos(lE);
333             sinLE  = FastMath.sin(lE);
334 
335         } while ((++iter < 50) && (FastMath.abs(shift) > 1.0e-12));
336 
337         return lE;
338 
339     }
340 
341     /** {@inheritDoc} */
342     public double getE() {
343         return FastMath.sqrt(ex * ex + ey * ey);
344     }
345 
346     /** {@inheritDoc} */
347     public double getI() {
348         return 2 * FastMath.atan(FastMath.sqrt(hx * hx + hy * hy));
349     }
350 
351     /** {@inheritDoc} */
352     protected TimeStampedPVCoordinates initPVCoordinates() {
353 
354         // get equinoctial parameters
355         final double lE = getLE();
356 
357         // inclination-related intermediate parameters
358         final double hx2   = hx * hx;
359         final double hy2   = hy * hy;
360         final double factH = 1. / (1 + hx2 + hy2);
361 
362         // reference axes defining the orbital plane
363         final double ux = (1 + hx2 - hy2) * factH;
364         final double uy =  2 * hx * hy * factH;
365         final double uz = -2 * hy * factH;
366 
367         final double vx = uy;
368         final double vy = (1 - hx2 + hy2) * factH;
369         final double vz =  2 * hx * factH;
370 
371         // eccentricity-related intermediate parameters
372         final double exey = ex * ey;
373         final double ex2  = ex * ex;
374         final double ey2  = ey * ey;
375         final double e2   = ex2 + ey2;
376         final double eta  = 1 + FastMath.sqrt(1 - e2);
377         final double beta = 1. / eta;
378 
379         // eccentric longitude argument
380         final double cLe    = FastMath.cos(lE);
381         final double sLe    = FastMath.sin(lE);
382         final double exCeyS = ex * cLe + ey * sLe;
383 
384         // coordinates of position and velocity in the orbital plane
385         final double x      = a * ((1 - beta * ey2) * cLe + beta * exey * sLe - ex);
386         final double y      = a * ((1 - beta * ex2) * sLe + beta * exey * cLe - ey);
387 
388         final double factor = FastMath.sqrt(getMu() / a) / (1 - exCeyS);
389         final double xdot   = factor * (-sLe + beta * ey * exCeyS);
390         final double ydot   = factor * ( cLe - beta * ex * exCeyS);
391 
392         final Vector3D position =
393             new Vector3D(x * ux + y * vx, x * uy + y * vy, x * uz + y * vz);
394         final double r2 = position.getNormSq();
395         final Vector3D velocity =
396             new Vector3D(xdot * ux + ydot * vx, xdot * uy + ydot * vy, xdot * uz + ydot * vz);
397 
398         final Vector3D acceleration = new Vector3D(-getMu() / (r2 * FastMath.sqrt(r2)), position);
399 
400         return new TimeStampedPVCoordinates(getDate(), position, velocity, acceleration);
401 
402     }
403 
404     /** {@inheritDoc} */
405     public EquinoctialOrbit shiftedBy(final double dt) {
406         return new EquinoctialOrbit(a, ex, ey, hx, hy,
407                                     getLM() + getKeplerianMeanMotion() * dt,
408                                     PositionAngle.MEAN, getFrame(),
409                                     getDate().shiftedBy(dt), getMu());
410     }
411 
412     /** {@inheritDoc}
413      * <p>
414      * The interpolated instance is created by polynomial Hermite interpolation
415      * on equinoctial elements, without derivatives (which means the interpolation
416      * falls back to Lagrange interpolation only).
417      * </p>
418      * <p>
419      * As this implementation of interpolation is polynomial, it should be used only
420      * with small samples (about 10-20 points) in order to avoid <a
421      * href="http://en.wikipedia.org/wiki/Runge%27s_phenomenon">Runge's phenomenon</a>
422      * and numerical problems (including NaN appearing).
423      * </p>
424      * <p>
425      * If orbit interpolation on large samples is needed, using the {@link
426      * org.orekit.propagation.analytical.Ephemeris} class is a better way than using this
427      * low-level interpolation. The Ephemeris class automatically handles selection of
428      * a neighboring sub-sample with a predefined number of point from a large global sample
429      * in a thread-safe way.
430      * </p>
431      */
432     public EquinoctialOrbit interpolate(final AbsoluteDate date, final Collection<Orbit> sample) {
433 
434         // set up an interpolator
435         final HermiteInterpolator interpolator = new HermiteInterpolator();
436 
437         // add sample points
438         AbsoluteDate previousDate = null;
439         double previousLm = Double.NaN;
440         for (final Orbit orbit : sample) {
441             final EquinoctialOrbit equi = (EquinoctialOrbit) OrbitType.EQUINOCTIAL.convertType(orbit);
442             final double continuousLm;
443             if (previousDate == null) {
444                 continuousLm = equi.getLM();
445             } else {
446                 final double dt       = equi.getDate().durationFrom(previousDate);
447                 final double keplerLm = previousLm + equi.getKeplerianMeanMotion() * dt;
448                 continuousLm = MathUtils.normalizeAngle(equi.getLM(), keplerLm);
449             }
450             previousDate = equi.getDate();
451             previousLm   = continuousLm;
452             interpolator.addSamplePoint(equi.getDate().durationFrom(date),
453                                         new double[] {
454                                             equi.getA(),
455                                             equi.getEquinoctialEx(),
456                                             equi.getEquinoctialEy(),
457                                             equi.getHx(),
458                                             equi.getHy(),
459                                             continuousLm
460                                         });
461         }
462 
463         // interpolate
464         final double[] interpolated = interpolator.value(0);
465 
466         // build a new interpolated instance
467         return new EquinoctialOrbit(interpolated[0], interpolated[1], interpolated[2],
468                                     interpolated[3], interpolated[4], interpolated[5],
469                                     PositionAngle.MEAN, getFrame(), date, getMu());
470 
471     }
472 
473     /** {@inheritDoc} */
474     protected double[][] computeJacobianMeanWrtCartesian() {
475 
476         final double[][] jacobian = new double[6][6];
477 
478         // compute various intermediate parameters
479         final Vector3D position = getPVCoordinates().getPosition();
480         final Vector3D velocity = getPVCoordinates().getVelocity();
481         final double r2         = position.getNormSq();
482         final double r          = FastMath.sqrt(r2);
483         final double r3         = r * r2;
484 
485         final double mu         = getMu();
486         final double sqrtMuA    = FastMath.sqrt(a * mu);
487         final double a2         = a * a;
488 
489         final double e2         = ex * ex + ey * ey;
490         final double oMe2       = 1 - e2;
491         final double epsilon    = FastMath.sqrt(oMe2);
492         final double beta       = 1 / (1 + epsilon);
493         final double ratio      = epsilon * beta;
494 
495         final double hx2        = hx * hx;
496         final double hy2        = hy * hy;
497         final double hxhy       = hx * hy;
498 
499         // precomputing equinoctial frame unit vectors (f,g,w)
500         final Vector3D f  = new Vector3D(1 - hy2 + hx2, 2 * hxhy, -2 * hy).normalize();
501         final Vector3D g  = new Vector3D(2 * hxhy, 1 + hy2 - hx2, 2 * hx).normalize();
502         final Vector3D w  = Vector3D.crossProduct(position, velocity).normalize();
503 
504         // coordinates of the spacecraft in the equinoctial frame
505         final double x    = Vector3D.dotProduct(position, f);
506         final double y    = Vector3D.dotProduct(position, g);
507         final double xDot = Vector3D.dotProduct(velocity, f);
508         final double yDot = Vector3D.dotProduct(velocity, g);
509 
510         // drDot / dEx = dXDot / dEx * f + dYDot / dEx * g
511         final double c1 = a / (sqrtMuA * epsilon);
512         final double c2 = a * sqrtMuA * beta / r3;
513         final double c3 = sqrtMuA / (r3 * epsilon);
514         final Vector3D drDotSdEx = new Vector3D( c1 * xDot * yDot - c2 * ey * x - c3 * x * y, f,
515                                                 -c1 * xDot * xDot - c2 * ey * y + c3 * x * x, g);
516 
517         // drDot / dEy = dXDot / dEy * f + dYDot / dEy * g
518         final Vector3D drDotSdEy = new Vector3D( c1 * yDot * yDot + c2 * ex * x - c3 * y * y, f,
519                                                 -c1 * xDot * yDot + c2 * ex * y + c3 * x * y, g);
520 
521         // da
522         final Vector3D vectorAR = new Vector3D(2 * a2 / r3, position);
523         final Vector3D vectorARDot = new Vector3D(2 * a2 / mu, velocity);
524         fillHalfRow(1, vectorAR,    jacobian[0], 0);
525         fillHalfRow(1, vectorARDot, jacobian[0], 3);
526 
527         // dEx
528         final double d1 = -a * ratio / r3;
529         final double d2 = (hy * xDot - hx * yDot) / (sqrtMuA * epsilon);
530         final double d3 = (hx * y - hy * x) / sqrtMuA;
531         final Vector3D vectorExRDot =
532             new Vector3D((2 * x * yDot - xDot * y) / mu, g, -y * yDot / mu, f, -ey * d3 / epsilon, w);
533         fillHalfRow(ex * d1, position, -ey * d2, w, epsilon / sqrtMuA, drDotSdEy, jacobian[1], 0);
534         fillHalfRow(1, vectorExRDot, jacobian[1], 3);
535 
536         // dEy
537         final Vector3D vectorEyRDot =
538             new Vector3D((2 * xDot * y - x * yDot) / mu, f, -x * xDot / mu, g, ex * d3 / epsilon, w);
539         fillHalfRow(ey * d1, position, ex * d2, w, -epsilon / sqrtMuA, drDotSdEx, jacobian[2], 0);
540         fillHalfRow(1, vectorEyRDot, jacobian[2], 3);
541 
542         // dHx
543         final double h = (1 + hx2 + hy2) / (2 * sqrtMuA * epsilon);
544         fillHalfRow(-h * xDot, w, jacobian[3], 0);
545         fillHalfRow( h * x,    w, jacobian[3], 3);
546 
547        // dHy
548         fillHalfRow(-h * yDot, w, jacobian[4], 0);
549         fillHalfRow( h * y,    w, jacobian[4], 3);
550 
551         // dLambdaM
552         final double l = -ratio / sqrtMuA;
553         fillHalfRow(-1 / sqrtMuA, velocity, d2, w, l * ex, drDotSdEx, l * ey, drDotSdEy, jacobian[5], 0);
554         fillHalfRow(-2 / sqrtMuA, position, ex * beta, vectorEyRDot, -ey * beta, vectorExRDot, d3, w, jacobian[5], 3);
555 
556         return jacobian;
557 
558     }
559 
560     /** {@inheritDoc} */
561     protected double[][] computeJacobianEccentricWrtCartesian() {
562 
563         // start by computing the Jacobian with mean angle
564         final double[][] jacobian = computeJacobianMeanWrtCartesian();
565 
566         // Differentiating the Kepler equation lM = lE - ex sin lE + ey cos lE leads to:
567         // dlM = (1 - ex cos lE - ey sin lE) dE - sin lE dex + cos lE dey
568         // which is inverted and rewritten as:
569         // dlE = a/r dlM + sin lE a/r dex - cos lE a/r dey
570         final double le    = getLE();
571         final double cosLe = FastMath.cos(le);
572         final double sinLe = FastMath.sin(le);
573         final double aOr   = 1 / (1 - ex * cosLe - ey * sinLe);
574 
575         // update longitude row
576         final double[] rowEx = jacobian[1];
577         final double[] rowEy = jacobian[2];
578         final double[] rowL  = jacobian[5];
579         for (int j = 0; j < 6; ++j) {
580             rowL[j] = aOr * (rowL[j] + sinLe * rowEx[j] - cosLe * rowEy[j]);
581         }
582 
583         return jacobian;
584 
585     }
586 
587     /** {@inheritDoc} */
588     protected double[][] computeJacobianTrueWrtCartesian() {
589 
590         // start by computing the Jacobian with eccentric angle
591         final double[][] jacobian = computeJacobianEccentricWrtCartesian();
592 
593         // Differentiating the eccentric longitude equation
594         // tan((lV - lE)/2) = [ex sin lE - ey cos lE] / [sqrt(1-ex^2-ey^2) + 1 - ex cos lE - ey sin lE]
595         // leads to
596         // cT (dlV - dlE) = cE dlE + cX dex + cY dey
597         // with
598         // cT = [d^2 + (ex sin lE - ey cos lE)^2] / 2
599         // d  = 1 + sqrt(1-ex^2-ey^2) - ex cos lE - ey sin lE
600         // cE = (ex cos lE + ey sin lE) (sqrt(1-ex^2-ey^2) + 1) - ex^2 - ey^2
601         // cX =  sin lE (sqrt(1-ex^2-ey^2) + 1) - ey + ex (ex sin lE - ey cos lE) / sqrt(1-ex^2-ey^2)
602         // cY = -cos lE (sqrt(1-ex^2-ey^2) + 1) + ex + ey (ex sin lE - ey cos lE) / sqrt(1-ex^2-ey^2)
603         // which can be solved to find the differential of the true longitude
604         // dlV = (cT + cE) / cT dlE + cX / cT deX + cY / cT deX
605         final double le        = getLE();
606         final double cosLe     = FastMath.cos(le);
607         final double sinLe     = FastMath.sin(le);
608         final double eSinE     = ex * sinLe - ey * cosLe;
609         final double ecosE     = ex * cosLe + ey * sinLe;
610         final double e2        = ex * ex + ey * ey;
611         final double epsilon   = FastMath.sqrt(1 - e2);
612         final double onePeps   = 1 + epsilon;
613         final double d         = onePeps - ecosE;
614         final double cT        = (d * d + eSinE * eSinE) / 2;
615         final double cE        = ecosE * onePeps - e2;
616         final double cX        = ex * eSinE / epsilon - ey + sinLe * onePeps;
617         final double cY        = ey * eSinE / epsilon + ex - cosLe * onePeps;
618         final double factorLe  = (cT + cE) / cT;
619         final double factorEx  = cX / cT;
620         final double factorEy  = cY / cT;
621 
622         // update longitude row
623         final double[] rowEx = jacobian[1];
624         final double[] rowEy = jacobian[2];
625         final double[] rowL = jacobian[5];
626         for (int j = 0; j < 6; ++j) {
627             rowL[j] = factorLe * rowL[j] + factorEx * rowEx[j] + factorEy * rowEy[j];
628         }
629 
630         return jacobian;
631 
632     }
633 
634     /** {@inheritDoc} */
635     public void addKeplerContribution(final PositionAngle type, final double gm,
636                                       final double[] pDot) {
637         final double oMe2;
638         final double ksi;
639         final double n = FastMath.sqrt(gm / a) / a;
640         switch (type) {
641             case MEAN :
642                 pDot[5] += n;
643                 break;
644             case ECCENTRIC :
645                 oMe2 = 1 - ex * ex - ey * ey;
646                 ksi  = 1 + ex * FastMath.cos(lv) + ey * FastMath.sin(lv);
647                 pDot[5] += n * ksi / oMe2;
648                 break;
649             case TRUE :
650                 oMe2 = 1 - ex * ex - ey * ey;
651                 ksi  = 1 + ex * FastMath.cos(lv) + ey * FastMath.sin(lv);
652                 pDot[5] += n * ksi * ksi / (oMe2 * FastMath.sqrt(oMe2));
653                 break;
654             default :
655                 throw new OrekitInternalError(null);
656         }
657     }
658 
659     /**  Returns a string representation of this equinoctial parameters object.
660      * @return a string representation of this object
661      */
662     public String toString() {
663         return new StringBuffer().append("equinoctial parameters: ").append('{').
664                                   append("a: ").append(a).
665                                   append("; ex: ").append(ex).append("; ey: ").append(ey).
666                                   append("; hx: ").append(hx).append("; hy: ").append(hy).
667                                   append("; lv: ").append(FastMath.toDegrees(lv)).
668                                   append(";}").toString();
669     }
670 
671     /** Replace the instance with a data transfer object for serialization.
672      * @return data transfer object that will be serialized
673      */
674     private Object writeReplace() {
675         return new DTO(this);
676     }
677 
678     /** Internal class used only for serialization. */
679     private static class DTO implements Serializable {
680 
681         /** Serializable UID. */
682         private static final long serialVersionUID = 20140617L;
683 
684         /** Double values. */
685         private double[] d;
686 
687         /** Frame in which are defined the orbital parameters. */
688         private final Frame frame;
689 
690         /** Simple constructor.
691          * @param orbit instance to serialize
692          */
693         private DTO(final EquinoctialOrbit orbit) {
694 
695             final TimeStampedPVCoordinates pv = orbit.getPVCoordinates();
696 
697             // decompose date
698             final double epoch  = FastMath.floor(pv.getDate().durationFrom(AbsoluteDate.J2000_EPOCH));
699             final double offset = pv.getDate().durationFrom(AbsoluteDate.J2000_EPOCH.shiftedBy(epoch));
700 
701             this.d = new double[] {
702                 epoch, offset, orbit.getMu(),
703                 orbit.a, orbit.ex, orbit.ey,
704                 orbit.hx, orbit.hy, orbit.lv
705             };
706 
707             this.frame = orbit.getFrame();
708 
709         }
710 
711         /** Replace the deserialized data transfer object with a {@link EquinoctialOrbit}.
712          * @return replacement {@link EquinoctialOrbit}
713          */
714         private Object readResolve() {
715             return new EquinoctialOrbit(d[3], d[4], d[5], d[6], d[7], d[8], PositionAngle.TRUE,
716                                         frame, AbsoluteDate.J2000_EPOCH.shiftedBy(d[0]).shiftedBy(d[1]),
717                                         d[2]);
718         }
719 
720     }
721 
722 }