1   /* Copyright 2002-2018 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.stream.Stream;
21  
22  import org.hipparchus.analysis.differentiation.DSFactory;
23  import org.hipparchus.analysis.differentiation.DerivativeStructure;
24  import org.hipparchus.exception.LocalizedCoreFormats;
25  import org.hipparchus.exception.MathIllegalStateException;
26  import org.hipparchus.geometry.euclidean.threed.FieldVector3D;
27  import org.hipparchus.geometry.euclidean.threed.Rotation;
28  import org.hipparchus.geometry.euclidean.threed.RotationConvention;
29  import org.hipparchus.geometry.euclidean.threed.Vector3D;
30  import org.hipparchus.util.FastMath;
31  import org.orekit.errors.OrekitMessages;
32  import org.orekit.frames.Frame;
33  import org.orekit.time.AbsoluteDate;
34  import org.orekit.utils.CartesianDerivativesFilter;
35  import org.orekit.utils.PVCoordinates;
36  import org.orekit.utils.TimeStampedPVCoordinates;
37  
38  
39  /** This class holds Cartesian orbital parameters.
40  
41   * <p>
42   * The parameters used internally are the Cartesian coordinates:
43   *   <ul>
44   *     <li>x</li>
45   *     <li>y</li>
46   *     <li>z</li>
47   *     <li>xDot</li>
48   *     <li>yDot</li>
49   *     <li>zDot</li>
50   *   </ul>
51   * contained in {@link PVCoordinates}.
52   *
53  
54   * <p>
55   * Note that the implementation of this class delegates all non-Cartesian related
56   * computations ({@link #getA()}, {@link #getEquinoctialEx()}, ...) to an underlying
57   * instance of the {@link EquinoctialOrbit} class. This implies that using this class
58   * only for analytical computations which are always based on non-Cartesian
59   * parameters is perfectly possible but somewhat sub-optimal.
60   * </p>
61   * <p>
62   * The instance <code>CartesianOrbit</code> is guaranteed to be immutable.
63   * </p>
64   * @see    Orbit
65   * @see    KeplerianOrbit
66   * @see    CircularOrbit
67   * @see    EquinoctialOrbit
68   * @author Luc Maisonobe
69   * @author Guylaine Prat
70   * @author Fabien Maussion
71   * @author V&eacute;ronique Pommier-Maurussane
72   */
73  public class CartesianOrbit extends Orbit {
74  
75      /** Serializable UID. */
76      private static final long serialVersionUID = 20170414L;
77  
78      /** Factory for first time derivatives. */
79      private static final DSFactory FACTORY = new DSFactory(1, 1);
80  
81      /** Indicator for non-Keplerian derivatives. */
82      private final transient boolean hasNonKeplerianAcceleration;
83  
84      /** Underlying equinoctial orbit to which high-level methods are delegated. */
85      private transient EquinoctialOrbit equinoctial;
86  
87      /** Constructor from Cartesian parameters.
88       *
89       * <p> The acceleration provided in {@code pvCoordinates} is accessible using
90       * {@link #getPVCoordinates()} and {@link #getPVCoordinates(Frame)}. All other methods
91       * use {@code mu} and the position to compute the acceleration, including
92       * {@link #shiftedBy(double)} and {@link #getPVCoordinates(AbsoluteDate, Frame)}.
93       *
94       * @param pvaCoordinates the position, velocity and acceleration of the satellite.
95       * @param frame the frame in which the {@link PVCoordinates} are defined
96       * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
97       * @param mu central attraction coefficient (m³/s²)
98       * @exception IllegalArgumentException if frame is not a {@link
99       * Frame#isPseudoInertial pseudo-inertial frame}
100      */
101     public CartesianOrbit(final TimeStampedPVCoordinates pvaCoordinates,
102                           final Frame frame, final double mu)
103         throws IllegalArgumentException {
104         super(pvaCoordinates, frame, mu);
105         hasNonKeplerianAcceleration = hasNonKeplerianAcceleration(pvaCoordinates, mu);
106         equinoctial = null;
107     }
108 
109     /** Constructor from Cartesian parameters.
110      *
111      * <p> The acceleration provided in {@code pvCoordinates} is accessible using
112      * {@link #getPVCoordinates()} and {@link #getPVCoordinates(Frame)}. All other methods
113      * use {@code mu} and the position to compute the acceleration, including
114      * {@link #shiftedBy(double)} and {@link #getPVCoordinates(AbsoluteDate, Frame)}.
115      *
116      * @param pvaCoordinates the position and velocity of the satellite.
117      * @param frame the frame in which the {@link PVCoordinates} are defined
118      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
119      * @param date date of the orbital parameters
120      * @param mu central attraction coefficient (m³/s²)
121      * @exception IllegalArgumentException if frame is not a {@link
122      * Frame#isPseudoInertial pseudo-inertial frame}
123      */
124     public CartesianOrbit(final PVCoordinates pvaCoordinates, final Frame frame,
125                           final AbsoluteDate date, final double mu)
126         throws IllegalArgumentException {
127         this(new TimeStampedPVCoordinates(date, pvaCoordinates), frame, mu);
128     }
129 
130     /** Constructor from any kind of orbital parameters.
131      * @param op orbital parameters to copy
132      */
133     public CartesianOrbit(final Orbit op) {
134         super(op.getPVCoordinates(), op.getFrame(), op.getMu());
135         hasNonKeplerianAcceleration = op.hasDerivatives();
136         if (op instanceof EquinoctialOrbit) {
137             equinoctial = (EquinoctialOrbit) op;
138         } else if (op instanceof CartesianOrbit) {
139             equinoctial = ((CartesianOrbit) op).equinoctial;
140         } else {
141             equinoctial = null;
142         }
143     }
144 
145     /** {@inheritDoc} */
146     public OrbitType getType() {
147         return OrbitType.CARTESIAN;
148     }
149 
150     /** Lazy evaluation of equinoctial parameters. */
151     private void initEquinoctial() {
152         if (equinoctial == null) {
153             if (hasDerivatives()) {
154                 // getPVCoordinates includes accelerations that will be interpreted as derivatives
155                 equinoctial = new EquinoctialOrbit(getPVCoordinates(), getFrame(), getDate(), getMu());
156             } else {
157                 // get rid of Keplerian acceleration so we don't assume
158                 // we have derivatives when in fact we don't have them
159                 equinoctial = new EquinoctialOrbit(new PVCoordinates(getPVCoordinates().getPosition(),
160                                                                      getPVCoordinates().getVelocity()),
161                                                    getFrame(), getDate(), getMu());
162             }
163         }
164     }
165 
166     /** Get position with derivatives.
167      * @return position with derivatives
168      */
169     private FieldVector3D<DerivativeStructure> getPositionDS() {
170         final Vector3D p = getPVCoordinates().getPosition();
171         final Vector3D v = getPVCoordinates().getVelocity();
172         return new FieldVector3D<>(FACTORY.build(p.getX(), v.getX()),
173                                    FACTORY.build(p.getY(), v.getY()),
174                                    FACTORY.build(p.getZ(), v.getZ()));
175     }
176 
177     /** Get velocity with derivatives.
178      * @return velocity with derivatives
179      */
180     private FieldVector3D<DerivativeStructure> getVelocityDS() {
181         final Vector3D v = getPVCoordinates().getVelocity();
182         final Vector3D a = getPVCoordinates().getAcceleration();
183         return new FieldVector3D<>(FACTORY.build(v.getX(), a.getX()),
184                                    FACTORY.build(v.getY(), a.getY()),
185                                    FACTORY.build(v.getZ(), a.getZ()));
186     }
187 
188     /** {@inheritDoc} */
189     public double getA() {
190         final double r  = getPVCoordinates().getPosition().getNorm();
191         final double V2 = getPVCoordinates().getVelocity().getNormSq();
192         return r / (2 - r * V2 / getMu());
193     }
194 
195     /** {@inheritDoc} */
196     public double getADot() {
197         if (hasDerivatives()) {
198             final DerivativeStructure r  = getPositionDS().getNorm();
199             final DerivativeStructure V2 = getVelocityDS().getNormSq();
200             final DerivativeStructure a  = r.divide(r.multiply(V2).divide(getMu()).subtract(2).negate());
201             return a.getPartialDerivative(1);
202         } else {
203             return Double.NaN;
204         }
205     }
206 
207     /** {@inheritDoc} */
208     public double getE() {
209         final double muA = getMu() * getA();
210         if (muA > 0) {
211             // elliptic or circular orbit
212             final Vector3D pvP   = getPVCoordinates().getPosition();
213             final Vector3D pvV   = getPVCoordinates().getVelocity();
214             final double rV2OnMu = pvP.getNorm() * pvV.getNormSq() / getMu();
215             final double eSE     = Vector3D.dotProduct(pvP, pvV) / FastMath.sqrt(muA);
216             final double eCE     = rV2OnMu - 1;
217             return FastMath.sqrt(eCE * eCE + eSE * eSE);
218         } else {
219             // hyperbolic orbit
220             final Vector3D pvM = getPVCoordinates().getMomentum();
221             return FastMath.sqrt(1 - pvM.getNormSq() / muA);
222         }
223     }
224 
225     /** {@inheritDoc} */
226     public double getEDot() {
227         if (hasDerivatives()) {
228             final FieldVector3D<DerivativeStructure> pvP   = getPositionDS();
229             final FieldVector3D<DerivativeStructure> pvV   = getVelocityDS();
230             final DerivativeStructure r       = getPositionDS().getNorm();
231             final DerivativeStructure V2      = getVelocityDS().getNormSq();
232             final DerivativeStructure rV2OnMu = r.multiply(V2).divide(getMu());
233             final DerivativeStructure a       = r.divide(rV2OnMu.negate().add(2));
234             final DerivativeStructure eSE     = FieldVector3D.dotProduct(pvP, pvV).divide(a.multiply(getMu()).sqrt());
235             final DerivativeStructure eCE     = rV2OnMu.subtract(1);
236             final DerivativeStructure e       = eCE.multiply(eCE).add(eSE.multiply(eSE)).sqrt();
237             return e.getPartialDerivative(1);
238         } else {
239             return Double.NaN;
240         }
241     }
242 
243     /** {@inheritDoc} */
244     public double getI() {
245         return Vector3D.angle(Vector3D.PLUS_K, getPVCoordinates().getMomentum());
246     }
247 
248     /** {@inheritDoc} */
249     public double getIDot() {
250         if (hasDerivatives()) {
251             final FieldVector3D<DerivativeStructure> momentum =
252                             FieldVector3D.crossProduct(getPositionDS(), getVelocityDS());
253             final DerivativeStructure i = FieldVector3D.angle(Vector3D.PLUS_K, momentum);
254             return i.getPartialDerivative(1);
255         } else {
256             return Double.NaN;
257         }
258     }
259 
260     /** {@inheritDoc} */
261     public double getEquinoctialEx() {
262         initEquinoctial();
263         return equinoctial.getEquinoctialEx();
264     }
265 
266     /** {@inheritDoc} */
267     public double getEquinoctialExDot() {
268         initEquinoctial();
269         return equinoctial.getEquinoctialExDot();
270     }
271 
272     /** {@inheritDoc} */
273     public double getEquinoctialEy() {
274         initEquinoctial();
275         return equinoctial.getEquinoctialEy();
276     }
277 
278     /** {@inheritDoc} */
279     public double getEquinoctialEyDot() {
280         initEquinoctial();
281         return equinoctial.getEquinoctialEyDot();
282     }
283 
284     /** {@inheritDoc} */
285     public double getHx() {
286         final Vector3D w = getPVCoordinates().getMomentum().normalize();
287         // Check for equatorial retrograde orbit
288         if (((w.getX() * w.getX() + w.getY() * w.getY()) == 0) && w.getZ() < 0) {
289             return Double.NaN;
290         }
291         return -w.getY() / (1 + w.getZ());
292     }
293 
294     /** {@inheritDoc} */
295     public double getHxDot() {
296         if (hasDerivatives()) {
297             final FieldVector3D<DerivativeStructure> w =
298                             FieldVector3D.crossProduct(getPositionDS(), getVelocityDS()).normalize();
299             // Check for equatorial retrograde orbit
300             final double x = w.getX().getValue();
301             final double y = w.getY().getValue();
302             final double z = w.getZ().getValue();
303             if (((x * x + y * y) == 0) && z < 0) {
304                 return Double.NaN;
305             }
306             final DerivativeStructure hx = w.getY().negate().divide(w.getZ().add(1));
307             return hx.getPartialDerivative(1);
308         } else {
309             return Double.NaN;
310         }
311     }
312 
313     /** {@inheritDoc} */
314     public double getHy() {
315         final Vector3D w = getPVCoordinates().getMomentum().normalize();
316         // Check for equatorial retrograde orbit
317         if (((w.getX() * w.getX() + w.getY() * w.getY()) == 0) && w.getZ() < 0) {
318             return Double.NaN;
319         }
320         return  w.getX() / (1 + w.getZ());
321     }
322 
323     /** {@inheritDoc} */
324     public double getHyDot() {
325         if (hasDerivatives()) {
326             final FieldVector3D<DerivativeStructure> w =
327                             FieldVector3D.crossProduct(getPositionDS(), getVelocityDS()).normalize();
328             // Check for equatorial retrograde orbit
329             final double x = w.getX().getValue();
330             final double y = w.getY().getValue();
331             final double z = w.getZ().getValue();
332             if (((x * x + y * y) == 0) && z < 0) {
333                 return Double.NaN;
334             }
335             final DerivativeStructure hy = w.getX().divide(w.getZ().add(1));
336             return hy.getPartialDerivative(1);
337         } else {
338             return Double.NaN;
339         }
340     }
341 
342     /** {@inheritDoc} */
343     public double getLv() {
344         initEquinoctial();
345         return equinoctial.getLv();
346     }
347 
348     /** {@inheritDoc} */
349     public double getLvDot() {
350         initEquinoctial();
351         return equinoctial.getLvDot();
352     }
353 
354     /** {@inheritDoc} */
355     public double getLE() {
356         initEquinoctial();
357         return equinoctial.getLE();
358     }
359 
360     /** {@inheritDoc} */
361     public double getLEDot() {
362         initEquinoctial();
363         return equinoctial.getLEDot();
364     }
365 
366     /** {@inheritDoc} */
367     public double getLM() {
368         initEquinoctial();
369         return equinoctial.getLM();
370     }
371 
372     /** {@inheritDoc} */
373     public double getLMDot() {
374         initEquinoctial();
375         return equinoctial.getLMDot();
376     }
377 
378     /** {@inheritDoc} */
379     public boolean hasDerivatives() {
380         return hasNonKeplerianAcceleration;
381     }
382 
383     /** {@inheritDoc} */
384     protected TimeStampedPVCoordinates initPVCoordinates() {
385         // nothing to do here, as the canonical elements are already the Cartesian ones
386         return getPVCoordinates();
387     }
388 
389     /** {@inheritDoc} */
390     public CartesianOrbit shiftedBy(final double dt) {
391         final PVCoordinates shiftedPV = (getA() < 0) ? shiftPVHyperbolic(dt) : shiftPVElliptic(dt);
392         return new CartesianOrbit(shiftedPV, getFrame(), getDate().shiftedBy(dt), getMu());
393     }
394 
395     /** {@inheritDoc}
396      * <p>
397      * The interpolated instance is created by polynomial Hermite interpolation
398      * ensuring velocity remains the exact derivative of position.
399      * </p>
400      * <p>
401      * As this implementation of interpolation is polynomial, it should be used only
402      * with small samples (about 10-20 points) in order to avoid <a
403      * href="http://en.wikipedia.org/wiki/Runge%27s_phenomenon">Runge's phenomenon</a>
404      * and numerical problems (including NaN appearing).
405      * </p>
406      * <p>
407      * If orbit interpolation on large samples is needed, using the {@link
408      * org.orekit.propagation.analytical.Ephemeris} class is a better way than using this
409      * low-level interpolation. The Ephemeris class automatically handles selection of
410      * a neighboring sub-sample with a predefined number of point from a large global sample
411      * in a thread-safe way.
412      * </p>
413      */
414     public CartesianOrbit interpolate(final AbsoluteDate date, final Stream<Orbit> sample) {
415         final TimeStampedPVCoordinates interpolated =
416                 TimeStampedPVCoordinates.interpolate(date, CartesianDerivativesFilter.USE_PVA,
417                                                      sample.map(orbit -> orbit.getPVCoordinates()));
418         return new CartesianOrbit(interpolated, getFrame(), date, getMu());
419     }
420 
421     /** Compute shifted position and velocity in elliptic case.
422      * @param dt time shift
423      * @return shifted position and velocity
424      */
425     private PVCoordinates shiftPVElliptic(final double dt) {
426 
427         // preliminary computation
428         final Vector3D pvP   = getPVCoordinates().getPosition();
429         final Vector3D pvV   = getPVCoordinates().getVelocity();
430         final double r2      = pvP.getNormSq();
431         final double r       = FastMath.sqrt(r2);
432         final double rV2OnMu = r * pvV.getNormSq() / getMu();
433         final double a       = getA();
434         final double eSE     = Vector3D.dotProduct(pvP, pvV) / FastMath.sqrt(getMu() * a);
435         final double eCE     = rV2OnMu - 1;
436         final double e2      = eCE * eCE + eSE * eSE;
437 
438         // we can use any arbitrary reference 2D frame in the orbital plane
439         // in order to simplify some equations below, we use the current position as the u axis
440         final Vector3D u     = pvP.normalize();
441         final Vector3D v     = Vector3D.crossProduct(getPVCoordinates().getMomentum(), u).normalize();
442 
443         // the following equations rely on the specific choice of u explained above,
444         // some coefficients that vanish to 0 in this case have already been removed here
445         final double ex      = (eCE - e2) * a / r;
446         final double ey      = -FastMath.sqrt(1 - e2) * eSE * a / r;
447         final double beta    = 1 / (1 + FastMath.sqrt(1 - e2));
448         final double thetaE0 = FastMath.atan2(ey + eSE * beta * ex, r / a + ex - eSE * beta * ey);
449         final double thetaM0 = thetaE0 - ex * FastMath.sin(thetaE0) + ey * FastMath.cos(thetaE0);
450 
451         // compute in-plane shifted eccentric argument
452         final double thetaM1 = thetaM0 + getKeplerianMeanMotion() * dt;
453         final double thetaE1 = meanToEccentric(thetaM1, ex, ey);
454         final double cTE     = FastMath.cos(thetaE1);
455         final double sTE     = FastMath.sin(thetaE1);
456 
457         // compute shifted in-plane Cartesian coordinates
458         final double exey   = ex * ey;
459         final double exCeyS = ex * cTE + ey * sTE;
460         final double x      = a * ((1 - beta * ey * ey) * cTE + beta * exey * sTE - ex);
461         final double y      = a * ((1 - beta * ex * ex) * sTE + beta * exey * cTE - ey);
462         final double factor = FastMath.sqrt(getMu() / a) / (1 - exCeyS);
463         final double xDot   = factor * (-sTE + beta * ey * exCeyS);
464         final double yDot   = factor * ( cTE - beta * ex * exCeyS);
465 
466         final Vector3D shiftedP = new Vector3D(x, u, y, v);
467         final Vector3D shiftedV = new Vector3D(xDot, u, yDot, v);
468         if (hasNonKeplerianAcceleration) {
469 
470             // extract non-Keplerian part of the initial acceleration
471             final Vector3D nonKeplerianAcceleration = new Vector3D(1, getPVCoordinates().getAcceleration(),
472                                                                    getMu() / (r2 * r), pvP);
473 
474             // add the quadratic motion due to the non-Keplerian acceleration to the Keplerian motion
475             final Vector3D fixedP   = new Vector3D(1, shiftedP,
476                                                    0.5 * dt * dt, nonKeplerianAcceleration);
477             final double   fixedR2 = fixedP.getNormSq();
478             final double   fixedR  = FastMath.sqrt(fixedR2);
479             final Vector3D fixedV  = new Vector3D(1, shiftedV,
480                                                   dt, nonKeplerianAcceleration);
481             final Vector3D fixedA  = new Vector3D(-getMu() / (fixedR2 * fixedR), shiftedP,
482                                                   1, nonKeplerianAcceleration);
483 
484             return new PVCoordinates(fixedP, fixedV, fixedA);
485 
486         } else {
487             // don't include acceleration,
488             // so the shifted orbit is not considered to have derivatives
489             return new PVCoordinates(shiftedP, shiftedV);
490         }
491 
492     }
493 
494     /** Compute shifted position and velocity in hyperbolic case.
495      * @param dt time shift
496      * @return shifted position and velocity
497      */
498     private PVCoordinates shiftPVHyperbolic(final double dt) {
499 
500         final PVCoordinates pv = getPVCoordinates();
501         final Vector3D pvP   = pv.getPosition();
502         final Vector3D pvV   = pv.getVelocity();
503         final Vector3D pvM   = pv.getMomentum();
504         final double r2      = pvP.getNormSq();
505         final double r       = FastMath.sqrt(r2);
506         final double rV2OnMu = r * pvV.getNormSq() / getMu();
507         final double a       = getA();
508         final double muA     = getMu() * a;
509         final double e       = FastMath.sqrt(1 - Vector3D.dotProduct(pvM, pvM) / muA);
510         final double sqrt    = FastMath.sqrt((e + 1) / (e - 1));
511 
512         // compute mean anomaly
513         final double eSH     = Vector3D.dotProduct(pvP, pvV) / FastMath.sqrt(-muA);
514         final double eCH     = rV2OnMu - 1;
515         final double H0      = FastMath.log((eCH + eSH) / (eCH - eSH)) / 2;
516         final double M0      = e * FastMath.sinh(H0) - H0;
517 
518         // find canonical 2D frame with p pointing to perigee
519         final double v0      = 2 * FastMath.atan(sqrt * FastMath.tanh(H0 / 2));
520         final Vector3D p     = new Rotation(pvM, v0, RotationConvention.FRAME_TRANSFORM).applyTo(pvP).normalize();
521         final Vector3D q     = Vector3D.crossProduct(pvM, p).normalize();
522 
523         // compute shifted eccentric anomaly
524         final double M1      = M0 + getKeplerianMeanMotion() * dt;
525         final double H1      = meanToHyperbolicEccentric(M1, e);
526 
527         // compute shifted in-plane Cartesian coordinates
528         final double cH     = FastMath.cosh(H1);
529         final double sH     = FastMath.sinh(H1);
530         final double sE2m1  = FastMath.sqrt((e - 1) * (e + 1));
531 
532         // coordinates of position and velocity in the orbital plane
533         final double x      = a * (cH - e);
534         final double y      = -a * sE2m1 * sH;
535         final double factor = FastMath.sqrt(getMu() / -a) / (e * cH - 1);
536         final double xDot   = -factor * sH;
537         final double yDot   =  factor * sE2m1 * cH;
538 
539         final Vector3D shiftedP = new Vector3D(x, p, y, q);
540         final Vector3D shiftedV = new Vector3D(xDot, p, yDot, q);
541         if (hasNonKeplerianAcceleration) {
542 
543             // extract non-Keplerian part of the initial acceleration
544             final Vector3D nonKeplerianAcceleration = new Vector3D(1, getPVCoordinates().getAcceleration(),
545                                                                    getMu() / (r2 * r), pvP);
546 
547             // add the quadratic motion due to the non-Keplerian acceleration to the Keplerian motion
548             final Vector3D fixedP   = new Vector3D(1, shiftedP,
549                                                    0.5 * dt * dt, nonKeplerianAcceleration);
550             final double   fixedR2 = fixedP.getNormSq();
551             final double   fixedR  = FastMath.sqrt(fixedR2);
552             final Vector3D fixedV  = new Vector3D(1, shiftedV,
553                                                   dt, nonKeplerianAcceleration);
554             final Vector3D fixedA  = new Vector3D(-getMu() / (fixedR2 * fixedR), shiftedP,
555                                                   1, nonKeplerianAcceleration);
556 
557             return new PVCoordinates(fixedP, fixedV, fixedA);
558 
559         } else {
560             // don't include acceleration,
561             // so the shifted orbit is not considered to have derivatives
562             return new PVCoordinates(shiftedP, shiftedV);
563         }
564 
565     }
566 
567     /** Computes the eccentric in-plane argument from the mean in-plane argument.
568      * @param thetaM = mean in-plane argument (rad)
569      * @param ex first component of eccentricity vector
570      * @param ey second component of eccentricity vector
571      * @return the eccentric in-plane argument.
572      */
573     private double meanToEccentric(final double thetaM, final double ex, final double ey) {
574         // Generalization of Kepler equation to in-plane parameters
575         // with thetaE = eta + E and
576         //      thetaM = eta + M = thetaE - ex.sin(thetaE) + ey.cos(thetaE)
577         // and eta being counted from an arbitrary reference in the orbital plane
578         double thetaE        = thetaM;
579         double thetaEMthetaM = 0.0;
580         int    iter          = 0;
581         do {
582             final double cosThetaE = FastMath.cos(thetaE);
583             final double sinThetaE = FastMath.sin(thetaE);
584 
585             final double f2 = ex * sinThetaE - ey * cosThetaE;
586             final double f1 = 1.0 - ex * cosThetaE - ey * sinThetaE;
587             final double f0 = thetaEMthetaM - f2;
588 
589             final double f12 = 2.0 * f1;
590             final double shift = f0 * f12 / (f1 * f12 - f0 * f2);
591 
592             thetaEMthetaM -= shift;
593             thetaE         = thetaM + thetaEMthetaM;
594 
595             if (FastMath.abs(shift) <= 1.0e-12) {
596                 return thetaE;
597             }
598 
599         } while (++iter < 50);
600 
601         throw new MathIllegalStateException(LocalizedCoreFormats.CONVERGENCE_FAILED);
602 
603     }
604 
605     /** Computes the hyperbolic eccentric anomaly from the mean anomaly.
606      * <p>
607      * The algorithm used here for solving hyperbolic Kepler equation is
608      * Danby's iterative method (3rd order) with Vallado's initial guess.
609      * </p>
610      * @param M mean anomaly (rad)
611      * @param ecc eccentricity
612      * @return the hyperbolic eccentric anomaly
613      */
614     private double meanToHyperbolicEccentric(final double M, final double ecc) {
615 
616         // Resolution of hyperbolic Kepler equation for Keplerian parameters
617 
618         // Initial guess
619         double H;
620         if (ecc < 1.6) {
621             if ((-FastMath.PI < M && M < 0.) || M > FastMath.PI) {
622                 H = M - ecc;
623             } else {
624                 H = M + ecc;
625             }
626         } else {
627             if (ecc < 3.6 && FastMath.abs(M) > FastMath.PI) {
628                 H = M - FastMath.copySign(ecc, M);
629             } else {
630                 H = M / (ecc - 1.);
631             }
632         }
633 
634         // Iterative computation
635         int iter = 0;
636         do {
637             final double f3  = ecc * FastMath.cosh(H);
638             final double f2  = ecc * FastMath.sinh(H);
639             final double f1  = f3 - 1.;
640             final double f0  = f2 - H - M;
641             final double f12 = 2. * f1;
642             final double d   = f0 / f12;
643             final double fdf = f1 - d * f2;
644             final double ds  = f0 / fdf;
645 
646             final double shift = f0 / (fdf + ds * ds * f3 / 6.);
647 
648             H -= shift;
649 
650             if (FastMath.abs(shift) <= 1.0e-12) {
651                 return H;
652             }
653 
654         } while (++iter < 50);
655 
656         throw new MathIllegalStateException(OrekitMessages.UNABLE_TO_COMPUTE_HYPERBOLIC_ECCENTRIC_ANOMALY,
657                                             iter);
658     }
659 
660     /** Create a 6x6 identity matrix.
661      * @return 6x6 identity matrix
662      */
663     private double[][] create6x6Identity() {
664         // this is the fastest way to set the 6x6 identity matrix
665         final double[][] identity = new double[6][6];
666         for (int i = 0; i < 6; i++) {
667             identity[i][i] = 1.0;
668         }
669         return identity;
670     }
671 
672     @Override
673     protected double[][] computeJacobianMeanWrtCartesian() {
674         return create6x6Identity();
675     }
676 
677     @Override
678     protected double[][] computeJacobianEccentricWrtCartesian() {
679         return create6x6Identity();
680     }
681 
682     @Override
683     protected double[][] computeJacobianTrueWrtCartesian() {
684         return create6x6Identity();
685     }
686 
687     /** {@inheritDoc} */
688     public void addKeplerContribution(final PositionAngle type, final double gm,
689                                       final double[] pDot) {
690 
691         final PVCoordinates pv = getPVCoordinates();
692 
693         // position derivative is velocity
694         final Vector3D velocity = pv.getVelocity();
695         pDot[0] += velocity.getX();
696         pDot[1] += velocity.getY();
697         pDot[2] += velocity.getZ();
698 
699         // velocity derivative is Newtonian acceleration
700         final Vector3D position = pv.getPosition();
701         final double r2         = position.getNormSq();
702         final double coeff      = -gm / (r2 * FastMath.sqrt(r2));
703         pDot[3] += coeff * position.getX();
704         pDot[4] += coeff * position.getY();
705         pDot[5] += coeff * position.getZ();
706 
707     }
708 
709     /**  Returns a string representation of this Orbit object.
710      * @return a string representation of this object
711      */
712     public String toString() {
713         return "Cartesian parameters: " + getPVCoordinates().toString();
714     }
715 
716     /** Replace the instance with a data transfer object for serialization.
717      * <p>
718      * This intermediate class serializes all needed parameters,
719      * including position-velocity which are <em>not</em> serialized by parent
720      * {@link Orbit} class.
721      * </p>
722      * @return data transfer object that will be serialized
723      */
724     private Object writeReplace() {
725         return new DTO(this);
726     }
727 
728     /** Internal class used only for serialization. */
729     private static class DTO implements Serializable {
730 
731         /** Serializable UID. */
732         private static final long serialVersionUID = 20170414L;
733 
734         /** Double values. */
735         private double[] d;
736 
737         /** Frame in which are defined the orbital parameters. */
738         private final Frame frame;
739 
740         /** Simple constructor.
741          * @param orbit instance to serialize
742          */
743         private DTO(final CartesianOrbit orbit) {
744 
745             final TimeStampedPVCoordinates pv = orbit.getPVCoordinates();
746 
747             // decompose date
748             final double epoch  = FastMath.floor(pv.getDate().durationFrom(AbsoluteDate.J2000_EPOCH));
749             final double offset = pv.getDate().durationFrom(AbsoluteDate.J2000_EPOCH.shiftedBy(epoch));
750 
751             if (orbit.hasDerivatives()) {
752                 this.d = new double[] {
753                     epoch, offset, orbit.getMu(),
754                     pv.getPosition().getX(),     pv.getPosition().getY(),     pv.getPosition().getZ(),
755                     pv.getVelocity().getX(),     pv.getVelocity().getY(),     pv.getVelocity().getZ(),
756                     pv.getAcceleration().getX(), pv.getAcceleration().getY(), pv.getAcceleration().getZ()
757                 };
758             } else {
759                 this.d = new double[] {
760                     epoch, offset, orbit.getMu(),
761                     pv.getPosition().getX(),     pv.getPosition().getY(),     pv.getPosition().getZ(),
762                     pv.getVelocity().getX(),     pv.getVelocity().getY(),     pv.getVelocity().getZ()
763                 };
764             }
765 
766             this.frame = orbit.getFrame();
767 
768         }
769 
770         /** Replace the deserialized data transfer object with a {@link CartesianOrbit}.
771          * @return replacement {@link CartesianOrbit}
772          */
773         private Object readResolve() {
774             if (d.length >= 12) {
775                 // we have derivatives
776                 return new CartesianOrbit(new TimeStampedPVCoordinates(AbsoluteDate.J2000_EPOCH.shiftedBy(d[0]).shiftedBy(d[1]),
777                                                                        new Vector3D(d[3], d[ 4], d[ 5]),
778                                                                        new Vector3D(d[6], d[ 7], d[ 8]),
779                                                                        new Vector3D(d[9], d[10], d[11])),
780                                           frame, d[2]);
781             } else {
782                 // we don't have derivatives
783                 return new CartesianOrbit(new TimeStampedPVCoordinates(AbsoluteDate.J2000_EPOCH.shiftedBy(d[0]).shiftedBy(d[1]),
784                                                                        new Vector3D(d[3], d[ 4], d[ 5]),
785                                                                        new Vector3D(d[6], d[ 7], d[ 8])),
786                                           frame, d[2]);
787             }
788         }
789 
790     }
791 
792 }