1   /* Copyright 2002-2013 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.ArrayList;
21  import java.util.Collection;
22  import java.util.List;
23  
24  import org.apache.commons.math3.exception.ConvergenceException;
25  import org.apache.commons.math3.geometry.euclidean.threed.Rotation;
26  import org.apache.commons.math3.geometry.euclidean.threed.Vector3D;
27  import org.apache.commons.math3.util.FastMath;
28  import org.apache.commons.math3.util.Pair;
29  import org.orekit.errors.OrekitMessages;
30  import org.orekit.frames.Frame;
31  import org.orekit.time.AbsoluteDate;
32  import org.orekit.utils.PVCoordinates;
33  
34  
35  /** This class holds cartesian orbital parameters.
36  
37   * <p>
38   * The parameters used internally are the cartesian coordinates:
39   *   <ul>
40   *     <li>x</li>
41   *     <li>y</li>
42   *     <li>z</li>
43   *     <li>xDot</li>
44   *     <li>yDot</li>
45   *     <li>zDot</li>
46   *   </ul>
47   * contained in {@link PVCoordinates}.
48   * </p>
49  
50   * <p>
51   * Note that the implementation of this class delegates all non-cartesian related
52   * computations ({@link #getA()}, {@link #getEquinoctialEx()}, ...) to an underlying
53   * instance of the {@link EquinoctialOrbit} class. This implies that using this class
54   * only for analytical computations which are always based on non-cartesian
55   * parameters is perfectly possible but somewhat sub-optimal.
56   * </p>
57   * <p>
58   * The instance <code>CartesianOrbit</code> is guaranteed to be immutable.
59   * </p>
60   * @see    Orbit
61   * @see    KeplerianOrbit
62   * @see    CircularOrbit
63   * @see    EquinoctialOrbit
64   * @author Luc Maisonobe
65   * @author Guylaine Prat
66   * @author Fabien Maussion
67   * @author V&eacute;ronique Pommier-Maurussane
68   */
69  public class CartesianOrbit extends Orbit {
70  
71      /** Serializable UID. */
72      private static final long serialVersionUID = -5411308212620896302L;
73  
74      /** Underlying equinoctial orbit to which high-level methods are delegated. */
75      private transient EquinoctialOrbit equinoctial;
76  
77      /** Constructor from cartesian parameters.
78       * @param pvCoordinates the position and velocity of the satellite.
79       * @param frame the frame in which the {@link PVCoordinates} are defined
80       * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
81       * @param date date of the orbital parameters
82       * @param mu central attraction coefficient (m<sup>3</sup>/s<sup>2</sup>)
83       * @exception IllegalArgumentException if frame is not a {@link
84       * Frame#isPseudoInertial pseudo-inertial frame}
85       */
86      public CartesianOrbit(final PVCoordinates pvCoordinates, final Frame frame,
87                            final AbsoluteDate date, final double mu)
88          throws IllegalArgumentException {
89          super(pvCoordinates, frame, date, mu);
90          equinoctial = null;
91      }
92  
93      /** Constructor from any kind of orbital parameters.
94       * @param op orbital parameters to copy
95       */
96      public CartesianOrbit(final Orbit op) {
97          super(op.getPVCoordinates(), op.getFrame(), op.getDate(), op.getMu());
98          if (op instanceof EquinoctialOrbit) {
99              equinoctial = (EquinoctialOrbit) op;
100         } else if (op instanceof CartesianOrbit) {
101             equinoctial = ((CartesianOrbit) op).equinoctial;
102         } else {
103             equinoctial = null;
104         }
105     }
106 
107     /** {@inheritDoc} */
108     public OrbitType getType() {
109         return OrbitType.CARTESIAN;
110     }
111 
112     /** Lazy evaluation of equinoctial parameters. */
113     private void initEquinoctial() {
114         if (equinoctial == null) {
115             equinoctial = new EquinoctialOrbit(getPVCoordinates(), getFrame(), getDate(), getMu());
116         }
117     }
118 
119     /** {@inheritDoc} */
120     public double getA() {
121         // lazy evaluation of semi-major axis
122         final double r  = getPVCoordinates().getPosition().getNorm();
123         final double V2 = getPVCoordinates().getVelocity().getNormSq();
124         return r / (2 - r * V2 / getMu());
125     }
126 
127     /** {@inheritDoc} */
128     public double getE() {
129         final Vector3D pvP   = getPVCoordinates().getPosition();
130         final Vector3D pvV   = getPVCoordinates().getVelocity();
131         final double rV2OnMu = pvP.getNorm() * pvV.getNormSq() / getMu();
132         final double eSE     = Vector3D.dotProduct(pvP, pvV) / FastMath.sqrt(getMu() * getA());
133         final double eCE     = rV2OnMu - 1;
134         return FastMath.sqrt(eCE * eCE + eSE * eSE);
135     }
136 
137     /** {@inheritDoc} */
138     public double getI() {
139         return Vector3D.angle(Vector3D.PLUS_K, getPVCoordinates().getMomentum());
140     }
141 
142     /** {@inheritDoc} */
143     public double getEquinoctialEx() {
144         initEquinoctial();
145         return equinoctial.getEquinoctialEx();
146     }
147 
148     /** {@inheritDoc} */
149     public double getEquinoctialEy() {
150         initEquinoctial();
151         return equinoctial.getEquinoctialEy();
152     }
153 
154     /** {@inheritDoc} */
155     public double getHx() {
156         final Vector3D w = getPVCoordinates().getMomentum().normalize();
157         // Check for equatorial retrograde orbit
158         if (((w.getX() * w.getX() + w.getY() * w.getY()) == 0) && w.getZ() < 0) {
159             return Double.NaN;
160         }
161         return -w.getY() / (1 + w.getZ());
162     }
163 
164     /** {@inheritDoc} */
165     public double getHy() {
166         final Vector3D w = getPVCoordinates().getMomentum().normalize();
167         // Check for equatorial retrograde orbit
168         if (((w.getX() * w.getX() + w.getY() * w.getY()) == 0) && w.getZ() < 0) {
169             return Double.NaN;
170         }
171         return  w.getX() / (1 + w.getZ());
172     }
173 
174     /** {@inheritDoc} */
175     public double getLv() {
176         initEquinoctial();
177         return equinoctial.getLv();
178     }
179 
180     /** {@inheritDoc} */
181     public double getLE() {
182         initEquinoctial();
183         return equinoctial.getLE();
184     }
185 
186     /** {@inheritDoc} */
187     public double getLM() {
188         initEquinoctial();
189         return equinoctial.getLM();
190     }
191 
192     /** {@inheritDoc} */
193     protected PVCoordinates initPVCoordinates() {
194         // nothing to do here, as the canonical elements are already the cartesian ones
195         return getPVCoordinates();
196     }
197 
198     /** {@inheritDoc} */
199     public CartesianOrbit shiftedBy(final double dt) {
200         final PVCoordinates shiftedPV = (getA() < 0) ? shiftPVHyperbolic(dt) : shiftPVElliptic(dt);
201         return new CartesianOrbit(shiftedPV, getFrame(), getDate().shiftedBy(dt), getMu());
202     }
203 
204     /** {@inheritDoc}
205      * <p>
206      * The interpolated instance is created by polynomial Hermite interpolation
207      * ensuring velocity remains the exact derivative of position.
208      * </p>
209      * <p>
210      * As this implementation of interpolation is polynomial, it should be used only
211      * with small samples (about 10-20 points) in order to avoid <a
212      * href="http://en.wikipedia.org/wiki/Runge%27s_phenomenon">Runge's phenomenon</a>
213      * and numerical problems (including NaN appearing).
214      * </p>
215      * <p>
216      * If orbit interpolation on large samples is needed, using the {@link
217      * org.orekit.propagation.analytical.Ephemeris} class is a better way than using this
218      * low-level interpolation. The Ephemeris class automatically handles selection of
219      * a neighboring sub-sample with a predefined number of point from a large global sample
220      * in a thread-safe way.
221      * </p>
222      */
223     public CartesianOrbit interpolate(final AbsoluteDate date, final Collection<Orbit> sample) {
224         final List<Pair<AbsoluteDate, PVCoordinates>> datedPV =
225                 new ArrayList<Pair<AbsoluteDate, PVCoordinates>>(sample.size());
226         for (final Orbit orbit : sample) {
227             datedPV.add(new Pair<AbsoluteDate, PVCoordinates>(orbit.getDate(), orbit.getPVCoordinates()));
228         }
229         final PVCoordinates interpolated = PVCoordinates.interpolate(date, true, datedPV);
230         return new CartesianOrbit(interpolated, getFrame(), date, getMu());
231     }
232 
233     /** Compute shifted position and velocity in elliptic case.
234      * @param dt time shift
235      * @return shifted position and velocity
236      */
237     private PVCoordinates shiftPVElliptic(final double dt) {
238 
239         // preliminary computation
240         final Vector3D pvP   = getPVCoordinates().getPosition();
241         final Vector3D pvV   = getPVCoordinates().getVelocity();
242         final double r       = pvP.getNorm();
243         final double rV2OnMu = r * pvV.getNormSq() / getMu();
244         final double a       = getA();
245         final double eSE     = Vector3D.dotProduct(pvP, pvV) / FastMath.sqrt(getMu() * a);
246         final double eCE     = rV2OnMu - 1;
247         final double e2      = eCE * eCE + eSE * eSE;
248 
249         // we can use any arbitrary reference 2D frame in the orbital plane
250         // in order to simplify some equations below, we use the current position as the u axis
251         final Vector3D u     = pvP.normalize();
252         final Vector3D v     = Vector3D.crossProduct(getPVCoordinates().getMomentum(), u).normalize();
253 
254         // the following equations rely on the specific choice of u explained above,
255         // some coefficients that vanish to 0 in this case have already been removed here
256         final double ex      = (eCE - e2) * a / r;
257         final double ey      = -FastMath.sqrt(1 - e2) * eSE * a / r;
258         final double beta    = 1 / (1 + FastMath.sqrt(1 - e2));
259         final double thetaE0 = FastMath.atan2(ey + eSE * beta * ex, r / a + ex - eSE * beta * ey);
260         final double thetaM0 = thetaE0 - ex * FastMath.sin(thetaE0) + ey * FastMath.cos(thetaE0);
261 
262         // compute in-plane shifted eccentric argument
263         final double thetaM1 = thetaM0 + getKeplerianMeanMotion() * dt;
264         final double thetaE1 = meanToEccentric(thetaM1, ex, ey);
265         final double cTE     = FastMath.cos(thetaE1);
266         final double sTE     = FastMath.sin(thetaE1);
267 
268         // compute shifted in-plane cartesian coordinates
269         final double exey   = ex * ey;
270         final double exCeyS = ex * cTE + ey * sTE;
271         final double x      = a * ((1 - beta * ey * ey) * cTE + beta * exey * sTE - ex);
272         final double y      = a * ((1 - beta * ex * ex) * sTE + beta * exey * cTE - ey);
273         final double factor = FastMath.sqrt(getMu() / a) / (1 - exCeyS);
274         final double xDot   = factor * (-sTE + beta * ey * exCeyS);
275         final double yDot   = factor * ( cTE - beta * ex * exCeyS);
276 
277         return new PVCoordinates(new Vector3D(x, u, y, v), new Vector3D(xDot, u, yDot, v));
278 
279     }
280 
281     /** Compute shifted position and velocity in hyperbolic case.
282      * @param dt time shift
283      * @return shifted position and velocity
284      */
285     private PVCoordinates shiftPVHyperbolic(final double dt) {
286 
287         final PVCoordinates pv = getPVCoordinates();
288         final Vector3D pvP   = pv.getPosition();
289         final Vector3D pvV   = pv.getVelocity();
290         final Vector3D pvM   = pv.getMomentum();
291         final double r       = pvP.getNorm();
292         final double rV2OnMu = r * pvV.getNormSq() / getMu();
293         final double a       = getA();
294         final double muA     = getMu() * a;
295         final double e       = FastMath.sqrt(1 - Vector3D.dotProduct(pvM, pvM) / muA);
296         final double sqrt    = FastMath.sqrt((e + 1) / (e - 1));
297 
298         // compute mean anomaly
299         final double eSH     = Vector3D.dotProduct(pvP, pvV) / FastMath.sqrt(-muA);
300         final double eCH     = rV2OnMu - 1;
301         final double H0      = FastMath.log((eCH + eSH) / (eCH - eSH)) / 2;
302         final double M0      = e * FastMath.sinh(H0) - H0;
303 
304         // find canonical 2D frame with p pointing to perigee
305         final double v0      = 2 * FastMath.atan(sqrt * FastMath.tanh(H0 / 2));
306         final Vector3D p     = new Rotation(pvM, -v0).applyTo(pvP).normalize();
307         final Vector3D q     = Vector3D.crossProduct(pvM, p).normalize();
308 
309         // compute shifted eccentric anomaly
310         final double M1      = M0 + getKeplerianMeanMotion() * dt;
311         final double H1      = meanToHyperbolicEccentric(M1, e);
312 
313         // compute shifted in-plane cartesian coordinates
314         final double cH     = FastMath.cosh(H1);
315         final double sH     = FastMath.sinh(H1);
316         final double sE2m1  = FastMath.sqrt((e - 1) * (e + 1));
317 
318         // coordinates of position and velocity in the orbital plane
319         final double x      = a * (cH - e);
320         final double y      = -a * sE2m1 * sH;
321         final double factor = FastMath.sqrt(getMu() / -a) / (e * cH - 1);
322         final double xDot   = -factor * sH;
323         final double yDot   =  factor * sE2m1 * cH;
324 
325         return new PVCoordinates(new Vector3D(x, p, y, q), new Vector3D(xDot, p, yDot, q));
326 
327     }
328 
329     /** Computes the eccentric in-plane argument from the mean in-plane argument.
330      * @param thetaM = mean in-plane argument (rad)
331      * @param ex first component of eccentricity vector
332      * @param ey second component of eccentricity vector
333      * @return the eccentric in-plane argument.
334      */
335     private double meanToEccentric(final double thetaM, final double ex, final double ey) {
336         // Generalization of Kepler equation to in-plane parameters
337         // with thetaE = eta + E and
338         //      thetaM = eta + M = thetaE - ex.sin(thetaE) + ey.cos(thetaE)
339         // and eta being counted from an arbitrary reference in the orbital plane
340         double thetaE        = thetaM;
341         double thetaEMthetaM = 0.0;
342         int    iter          = 0;
343         do {
344             final double cosThetaE = FastMath.cos(thetaE);
345             final double sinThetaE = FastMath.sin(thetaE);
346 
347             final double f2 = ex * sinThetaE - ey * cosThetaE;
348             final double f1 = 1.0 - ex * cosThetaE - ey * sinThetaE;
349             final double f0 = thetaEMthetaM - f2;
350 
351             final double f12 = 2.0 * f1;
352             final double shift = f0 * f12 / (f1 * f12 - f0 * f2);
353 
354             thetaEMthetaM -= shift;
355             thetaE         = thetaM + thetaEMthetaM;
356 
357             if (FastMath.abs(shift) <= 1.0e-12) {
358                 return thetaE;
359             }
360 
361         } while (++iter < 50);
362 
363         throw new ConvergenceException();
364 
365     }
366 
367     /** Computes the hyperbolic eccentric anomaly from the mean anomaly.
368      * <p>
369      * The algorithm used here for solving hyperbolic Kepler equation is
370      * Danby's iterative method (3rd order) with Vallado's initial guess.
371      * </p>
372      * @param M mean anomaly (rad)
373      * @param ecc eccentricity
374      * @return the hyperbolic eccentric anomaly
375      */
376     private double meanToHyperbolicEccentric(final double M, final double ecc) {
377 
378         // Resolution of hyperbolic Kepler equation for keplerian parameters
379 
380         // Initial guess
381         double H;
382         if (ecc < 1.6) {
383             if ((-FastMath.PI < M && M < 0.) || M > FastMath.PI) {
384                 H = M - ecc;
385             } else {
386                 H = M + ecc;
387             }
388         } else {
389             if (ecc < 3.6 && FastMath.abs(M) > FastMath.PI) {
390                 H = M - FastMath.copySign(ecc, M);
391             } else {
392                 H = M / (ecc - 1.);
393             }
394         }
395 
396         // Iterative computation
397         int iter = 0;
398         do {
399             final double f3  = ecc * FastMath.cosh(H);
400             final double f2  = ecc * FastMath.sinh(H);
401             final double f1  = f3 - 1.;
402             final double f0  = f2 - H - M;
403             final double f12 = 2. * f1;
404             final double d   = f0 / f12;
405             final double fdf = f1 - d * f2;
406             final double ds  = f0 / fdf;
407 
408             final double shift = f0 / (fdf + ds * ds * f3 / 6.);
409 
410             H -= shift;
411 
412             if (FastMath.abs(shift) <= 1.0e-12) {
413                 return H;
414             }
415 
416         } while (++iter < 50);
417 
418         throw new ConvergenceException(OrekitMessages.UNABLE_TO_COMPUTE_HYPERBOLIC_ECCENTRIC_ANOMALY,
419                                        iter);
420     }
421 
422     @Override
423     public void getJacobianWrtCartesian(final PositionAngle type, final double[][] jacobian) {
424         // this is the fastest way to set the 6x6 identity matrix
425         for (int i = 0; i < 6; i++) {
426             for (int j = 0; j < 6; j++) {
427                 jacobian[i][j] = 0;
428             }
429             jacobian[i][i] = 1;
430         }
431     }
432 
433     @Override
434     protected double[][] computeJacobianMeanWrtCartesian() {
435         // not used
436         return null;
437     }
438     @Override
439     protected double[][] computeJacobianEccentricWrtCartesian() {
440         // not used
441         return null;
442     }
443 
444     @Override
445     protected double[][] computeJacobianTrueWrtCartesian() {
446         // not used
447         return null;
448     }
449 
450     /** {@inheritDoc} */
451     public void addKeplerContribution(final PositionAngle type, final double gm,
452                                       final double[] pDot) {
453 
454         final PVCoordinates pv = getPVCoordinates();
455 
456         // position derivative is velocity
457         final Vector3D velocity = pv.getVelocity();
458         pDot[0] += velocity.getX();
459         pDot[1] += velocity.getY();
460         pDot[2] += velocity.getZ();
461 
462         // velocity derivative is Newtonian acceleration
463         final Vector3D position = pv.getPosition();
464         final double r2         = position.getNormSq();
465         final double coeff      = -gm / (r2 * FastMath.sqrt(r2));
466         pDot[3] += coeff * position.getX();
467         pDot[4] += coeff * position.getY();
468         pDot[5] += coeff * position.getZ();
469 
470     }
471 
472     /**  Returns a string representation of this Orbit object.
473      * @return a string representation of this object
474      */
475     public String toString() {
476         return "cartesian parameters: " + getPVCoordinates().toString();
477     }
478 
479     /** Replace the instance with a data transfer object for serialization.
480      * <p>
481      * This intermediate class serializes all needed parameters,
482      * including position-velocity which are <em>not</em> serialized by parent
483      * {@link Orbit} class.
484      * </p>
485      * @return data transfer object that will be serialized
486      */
487     private Object writeReplace() {
488         return new DataTransferObject(getPVCoordinates(), getFrame(), getDate(), getMu());
489     }
490 
491     /** Internal class used only for serialization. */
492     private static class DataTransferObject implements Serializable {
493 
494         /** Serializable UID. */
495         private static final long serialVersionUID = 4184412866917874790L;
496 
497         /** Computed PVCoordinates. */
498         private PVCoordinates pvCoordinates;
499 
500         /** Frame in which are defined the orbital parameters. */
501         private final Frame frame;
502 
503         /** Date of the orbital parameters. */
504         private final AbsoluteDate date;
505 
506         /** Value of mu used to compute position and velocity (m<sup>3</sup>/s<sup>2</sup>). */
507         private final double mu;
508 
509         /** Simple constructor.
510          * @param pvCoordinates the position and velocity of the satellite.
511          * @param frame the frame in which the {@link PVCoordinates} are defined
512          * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
513          * @param date date of the orbital parameters
514          * @param mu central attraction coefficient (m<sup>3</sup>/s<sup>2</sup>)
515          */
516         private DataTransferObject(final PVCoordinates pvCoordinates, final Frame frame,
517                                    final AbsoluteDate date, final double mu) {
518             this.pvCoordinates = pvCoordinates;
519             this.frame         = frame;
520             this.date          = date;
521             this.mu            = mu;
522         }
523 
524         /** Replace the deserialized data transfer object with a {@link CartesianOrbit}.
525          * @return replacement {@link CartesianOrbit}
526          */
527         private Object readResolve() {
528             // build a new provider, with an empty cache
529             return new CartesianOrbit(pvCoordinates, frame, date, mu);
530         }
531 
532     }
533 
534 }