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