FieldKeplerianAnomalyUtility.java

  1. /* Copyright 2002-2022 CS GROUP
  2.  * Licensed to CS GROUP (CS) under one or more
  3.  * contributor license agreements.  See the NOTICE file distributed with
  4.  * this work for additional information regarding copyright ownership.
  5.  * CS licenses this file to You under the Apache License, Version 2.0
  6.  * (the "License"); you may not use this file except in compliance with
  7.  * the License.  You may obtain a copy of the License at
  8.  *
  9.  *   http://www.apache.org/licenses/LICENSE-2.0
  10.  *
  11.  * Unless required by applicable law or agreed to in writing, software
  12.  * distributed under the License is distributed on an "AS IS" BASIS,
  13.  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
  14.  * See the License for the specific language governing permissions and
  15.  * limitations under the License.
  16.  */
  17. package org.orekit.orbits;

  18. import org.hipparchus.CalculusFieldElement;
  19. import org.hipparchus.Field;
  20. import org.hipparchus.exception.MathIllegalStateException;
  21. import org.hipparchus.util.FastMath;
  22. import org.hipparchus.util.FieldSinCos;
  23. import org.hipparchus.util.MathUtils;
  24. import org.hipparchus.util.Precision;
  25. import org.orekit.errors.OrekitMessages;

  26. /**
  27.  * Utility methods for converting between different Keplerian anomalies.
  28.  * @author Luc Maisonobe
  29.  * @author Andrea Antolino
  30.  * @author Andrew Goetz
  31.  */
  32. public class FieldKeplerianAnomalyUtility {

  33.     /** First coefficient to compute elliptic Kepler equation solver starter. */
  34.     private static final double A;

  35.     /** Second coefficient to compute elliptic Kepler equation solver starter. */
  36.     private static final double B;

  37.     static {
  38.         final double k1 = 3 * FastMath.PI + 2;
  39.         final double k2 = FastMath.PI - 1;
  40.         final double k3 = 6 * FastMath.PI - 1;
  41.         A = 3 * k2 * k2 / k1;
  42.         B = k3 * k3 / (6 * k1);
  43.     }

  44.     private FieldKeplerianAnomalyUtility() {
  45.     }

  46.     /**
  47.      * Computes the elliptic true anomaly from the elliptic mean anomaly.
  48.      * @param <T> field type
  49.      * @param e eccentricity such that 0 &le; e &lt; 1
  50.      * @param M elliptic mean anomaly (rad)
  51.      * @return elliptic true anomaly (rad)
  52.      */
  53.     public static <T extends CalculusFieldElement<T>> T ellipticMeanToTrue(final T e, final T M) {
  54.         final T E = ellipticMeanToEccentric(e, M);
  55.         final T v = ellipticEccentricToTrue(e, E);
  56.         return v;
  57.     }

  58.     /**
  59.      * Computes the elliptic mean anomaly from the elliptic true anomaly.
  60.      * @param <T> field type
  61.      * @param e eccentricity such that 0 &le; e &lt; 1
  62.      * @param v elliptic true anomaly (rad)
  63.      * @return elliptic mean anomaly (rad)
  64.      */
  65.     public static <T extends CalculusFieldElement<T>> T ellipticTrueToMean(final T e, final T v) {
  66.         final T E = ellipticTrueToEccentric(e, v);
  67.         final T M = ellipticEccentricToMean(e, E);
  68.         return M;
  69.     }

  70.     /**
  71.      * Computes the elliptic true anomaly from the elliptic eccentric anomaly.
  72.      * @param <T> field type
  73.      * @param e eccentricity such that 0 &le; e &lt; 1
  74.      * @param E elliptic eccentric anomaly (rad)
  75.      * @return elliptic true anomaly (rad)
  76.      */
  77.     public static <T extends CalculusFieldElement<T>> T ellipticEccentricToTrue(final T e, final T E) {
  78.         final T beta = e.divide(e.multiply(e).negate().add(1).sqrt().add(1));
  79.         final FieldSinCos<T> scE = FastMath.sinCos(E);
  80.         return E.add(beta.multiply(scE.sin()).divide(beta.multiply(scE.cos()).subtract(1).negate()).atan().multiply(2));
  81.     }

  82.     /**
  83.      * Computes the elliptic eccentric anomaly from the elliptic true anomaly.
  84.      * @param <T> field type
  85.      * @param e eccentricity such that 0 &le; e &lt; 1
  86.      * @param v elliptic true anomaly (rad)
  87.      * @return elliptic eccentric anomaly (rad)
  88.      */
  89.     public static <T extends CalculusFieldElement<T>> T ellipticTrueToEccentric(final T e, final T v) {
  90.         final T beta = e.divide(e.multiply(e).negate().add(1).sqrt().add(1));
  91.         final FieldSinCos<T> scv = FastMath.sinCos(v);
  92.         return v.subtract((beta.multiply(scv.sin()).divide(beta.multiply(scv.cos()).add(1))).atan().multiply(2));
  93.     }

  94.     /**
  95.      * Computes the elliptic eccentric anomaly from the elliptic mean anomaly.
  96.      * <p>
  97.      * The algorithm used here for solving hyperbolic Kepler equation is from Odell,
  98.      * A.W., Gooding, R.H. "Procedures for solving Kepler's equation." Celestial
  99.      * Mechanics 38, 307–334 (1986). https://doi.org/10.1007/BF01238923
  100.      * </p>
  101.      * @param <T> field type
  102.      * @param e eccentricity such that 0 &le; e &lt; 1
  103.      * @param M elliptic mean anomaly (rad)
  104.      * @return elliptic eccentric anomaly (rad)
  105.      */
  106.     public static <T extends CalculusFieldElement<T>> T ellipticMeanToEccentric(final T e, final T M) {
  107.         // reduce M to [-PI PI) interval
  108.         final T reducedM = MathUtils.normalizeAngle(M, M.getField().getZero());

  109.         // compute start value according to A. W. Odell and R. H. Gooding S12 starter
  110.         T E;
  111.         if (reducedM.abs().getReal() < 1.0 / 6.0) {
  112.             if (FastMath.abs(reducedM.getReal()) < Precision.SAFE_MIN) {
  113.                 // this is an Orekit change to the S12 starter, mainly used when T is some kind
  114.                 // of derivative structure. If reducedM is 0.0, the derivative of cbrt is
  115.                 // infinite which induces NaN appearing later in the computation. As in this
  116.                 // case E and M are almost equal, we initialize E with reducedM
  117.                 E = reducedM;
  118.             } else {
  119.                 // this is the standard S12 starter
  120.                 E = reducedM.add(e.multiply((reducedM.multiply(6).cbrt()).subtract(reducedM)));
  121.             }
  122.         } else {
  123.             final T pi = e.getPi();
  124.             if (reducedM.getReal() < 0) {
  125.                 final T w = reducedM.add(pi);
  126.                 E = reducedM.add(e.multiply(w.multiply(A).divide(w.negate().add(B)).subtract(pi).subtract(reducedM)));
  127.             } else {
  128.                 final T w = reducedM.negate().add(pi);
  129.                 E = reducedM
  130.                         .add(e.multiply(w.multiply(A).divide(w.negate().add(B)).negate().subtract(reducedM).add(pi)));
  131.             }
  132.         }

  133.         final T e1 = e.negate().add(1);
  134.         final boolean noCancellationRisk = (e1.getReal() + E.getReal() * E.getReal() / 6) >= 0.1;

  135.         // perform two iterations, each consisting of one Halley step and one
  136.         // Newton-Raphson step
  137.         for (int j = 0; j < 2; ++j) {

  138.             final T f;
  139.             T fd;
  140.             final FieldSinCos<T> scE = FastMath.sinCos(E);
  141.             final T fdd = e.multiply(scE.sin());
  142.             final T fddd = e.multiply(scE.cos());

  143.             if (noCancellationRisk) {

  144.                 f = (E.subtract(fdd)).subtract(reducedM);
  145.                 fd = fddd.negate().add(1);
  146.             } else {

  147.                 f = eMeSinE(e, E).subtract(reducedM);
  148.                 final T s = E.multiply(0.5).sin();
  149.                 fd = e1.add(e.multiply(s).multiply(s).multiply(2));
  150.             }
  151.             final T dee = f.multiply(fd).divide(f.multiply(fdd).multiply(0.5).subtract(fd.multiply(fd)));

  152.             // update eccentric anomaly, using expressions that limit underflow problems
  153.             final T w = fd.add(dee.multiply(fdd.add(dee.multiply(fddd.divide(3)))).multiply(0.5));
  154.             fd = fd.add(dee.multiply(fdd.add(dee.multiply(fddd).multiply(0.5))));
  155.             E = E.subtract(f.subtract(dee.multiply(fd.subtract(w))).divide(fd));

  156.         }

  157.         // expand the result back to original range
  158.         E = E.add(M).subtract(reducedM);
  159.         return E;
  160.     }

  161.     /**
  162.      * Accurate computation of E - e sin(E).
  163.      * <p>
  164.      * This method is used when E is close to 0 and e close to 1, i.e. near the
  165.      * perigee of almost parabolic orbits
  166.      * </p>
  167.      * @param <T> field type
  168.      * @param e eccentricity
  169.      * @param E eccentric anomaly (rad)
  170.      * @return E - e sin(E)
  171.      */
  172.     private static <T extends CalculusFieldElement<T>> T eMeSinE(final T e, final T E) {
  173.         T x = (e.negate().add(1)).multiply(E.sin());
  174.         final T mE2 = E.negate().multiply(E);
  175.         T term = E;
  176.         double d = 0;
  177.         // the inequality test below IS intentional and should NOT be replaced by a
  178.         // check with a small tolerance
  179.         for (T x0 = E.getField().getZero().add(Double.NaN); !Double.valueOf(x.getReal())
  180.                 .equals(Double.valueOf(x0.getReal()));) {
  181.             d += 2;
  182.             term = term.multiply(mE2.divide(d * (d + 1)));
  183.             x0 = x;
  184.             x = x.subtract(term);
  185.         }
  186.         return x;
  187.     }

  188.     /**
  189.      * Computes the elliptic mean anomaly from the elliptic eccentric anomaly.
  190.      * @param <T> field type
  191.      * @param e eccentricity such that 0 &le; e &lt; 1
  192.      * @param E elliptic eccentric anomaly (rad)
  193.      * @return elliptic mean anomaly (rad)
  194.      */
  195.     public static <T extends CalculusFieldElement<T>> T ellipticEccentricToMean(final T e, final T E) {
  196.         return E.subtract(e.multiply(E.sin()));
  197.     }

  198.     /**
  199.      * Computes the hyperbolic true anomaly from the hyperbolic mean anomaly.
  200.      * @param <T> field type
  201.      * @param e eccentricity &gt; 1
  202.      * @param M hyperbolic mean anomaly
  203.      * @return hyperbolic true anomaly (rad)
  204.      */
  205.     public static <T extends CalculusFieldElement<T>> T hyperbolicMeanToTrue(final T e, final T M) {
  206.         final T H = hyperbolicMeanToEccentric(e, M);
  207.         final T v = hyperbolicEccentricToTrue(e, H);
  208.         return v;
  209.     }

  210.     /**
  211.      * Computes the hyperbolic mean anomaly from the hyperbolic true anomaly.
  212.      * @param <T> field type
  213.      * @param e eccentricity &gt; 1
  214.      * @param v hyperbolic true anomaly (rad)
  215.      * @return hyperbolic mean anomaly
  216.      */
  217.     public static <T extends CalculusFieldElement<T>> T hyperbolicTrueToMean(final T e, final T v) {
  218.         final T H = hyperbolicTrueToEccentric(e, v);
  219.         final T M = hyperbolicEccentricToMean(e, H);
  220.         return M;
  221.     }

  222.     /**
  223.      * Computes the hyperbolic true anomaly from the hyperbolic eccentric anomaly.
  224.      * @param <T> field type
  225.      * @param e eccentricity &gt; 1
  226.      * @param H hyperbolic eccentric anomaly
  227.      * @return hyperbolic true anomaly (rad)
  228.      */
  229.     public static <T extends CalculusFieldElement<T>> T hyperbolicEccentricToTrue(final T e, final T H) {
  230.         final T s = e.add(1).divide(e.subtract(1)).sqrt();
  231.         final T tanH = H.multiply(0.5).tanh();
  232.         return s.multiply(tanH).atan().multiply(2);
  233.     }

  234.     /**
  235.      * Computes the hyperbolic eccentric anomaly from the hyperbolic true anomaly.
  236.      * @param <T> field type
  237.      * @param e eccentricity &gt; 1
  238.      * @param v hyperbolic true anomaly (rad)
  239.      * @return hyperbolic eccentric anomaly
  240.      */
  241.     public static <T extends CalculusFieldElement<T>> T hyperbolicTrueToEccentric(final T e, final T v) {
  242.         final FieldSinCos<T> scv = FastMath.sinCos(v);
  243.         final T sinhH = e.multiply(e).subtract(1).sqrt().multiply(scv.sin()).divide(e.multiply(scv.cos()).add(1));
  244.         return sinhH.asinh();
  245.     }

  246.     /**
  247.      * Computes the hyperbolic eccentric anomaly from the hyperbolic mean anomaly.
  248.      * <p>
  249.      * The algorithm used here for solving hyperbolic Kepler equation is from
  250.      * Gooding, R.H., Odell, A.W. "The hyperbolic Kepler equation (and the elliptic
  251.      * equation revisited)." Celestial Mechanics 44, 267–282 (1988).
  252.      * https://doi.org/10.1007/BF01235540
  253.      * </p>
  254.      * @param <T> field type
  255.      * @param e eccentricity &gt; 1
  256.      * @param M hyperbolic mean anomaly
  257.      * @return hyperbolic eccentric anomaly
  258.      */
  259.     public static <T extends CalculusFieldElement<T>> T hyperbolicMeanToEccentric(final T e, final T M) {
  260.         final Field<T> field = e.getField();
  261.         final T zero = field.getZero();
  262.         final T one = field.getOne();
  263.         final T two = zero.add(2.0);
  264.         final T three = zero.add(3.0);
  265.         final T half = zero.add(0.5);
  266.         final T onePointFive = zero.add(1.5);
  267.         final T fourThirds = zero.add(4.0).divide(zero.add(3.0));

  268.         // Solve L = S - g * asinh(S) for S = sinh(H).
  269.         final T L = M.divide(e);
  270.         final T g = e.reciprocal();
  271.         final T g1 = one.subtract(g);

  272.         // Starter based on Lagrange's theorem.
  273.         T S = L;
  274.         if (L.isZero()) {
  275.             return M.getField().getZero();
  276.         }
  277.         final T cl = L.multiply(L).add(one).sqrt();
  278.         final T al = L.asinh();
  279.         final T w = g.multiply(g).multiply(al).divide(cl.multiply(cl).multiply(cl));
  280.         S = one.subtract(g.divide(cl));
  281.         S = L.add(g.multiply(al).divide(S.multiply(S).multiply(S)
  282.                 .add(w.multiply(L).multiply(onePointFive.subtract(fourThirds.multiply(g)))).cbrt()));

  283.         // Two iterations (at most) of Halley-then-Newton process.
  284.         for (int i = 0; i < 2; ++i) {
  285.             final T s0 = S.multiply(S);
  286.             final T s1 = s0.add(one);
  287.             final T s2 = s1.sqrt();
  288.             final T s3 = s1.multiply(s2);
  289.             final T fdd = g.multiply(S).divide(s3);
  290.             final T fddd = g.multiply(one.subtract(two.multiply(s0))).divide(s1.multiply(s3));

  291.             T f;
  292.             T fd;
  293.             if (s0.divide(zero.add(6.0)).add(g1).getReal() >= 0.5) {
  294.                 f = S.subtract(g.multiply(S.asinh())).subtract(L);
  295.                 fd = one.subtract(g.divide(s2));
  296.             } else {
  297.                 // Accurate computation of S - (1 - g1) * asinh(S)
  298.                 // when (g1, S) is close to (0, 0).
  299.                 final T t = S.divide(one.add(one.add(S.multiply(S)).sqrt()));
  300.                 final T tsq = t.multiply(t);
  301.                 T x = S.multiply(g1.add(g.multiply(tsq)));
  302.                 T term = two.multiply(g).multiply(t);
  303.                 T twoI1 = one;
  304.                 T x0;
  305.                 int j = 0;
  306.                 do {
  307.                     if (++j == 1000000) {
  308.                         // This isn't expected to happen, but it protects against an infinite loop.
  309.                         throw new MathIllegalStateException(
  310.                                 OrekitMessages.UNABLE_TO_COMPUTE_HYPERBOLIC_ECCENTRIC_ANOMALY, j);
  311.                     }
  312.                     twoI1 = twoI1.add(2.0);
  313.                     term = term.multiply(tsq);
  314.                     x0 = x;
  315.                     x = x.subtract(term.divide(twoI1));
  316.                 } while (x.getReal() != x0.getReal());
  317.                 f = x.subtract(L);
  318.                 fd = s0.divide(s2.add(one)).add(g1).divide(s2);
  319.             }
  320.             final T ds = f.multiply(fd).divide(half.multiply(f).multiply(fdd).subtract(fd.multiply(fd)));
  321.             final T stemp = S.add(ds);
  322.             if (S.getReal() == stemp.getReal()) {
  323.                 break;
  324.             }
  325.             f = f.add(ds.multiply(fd.add(half.multiply(ds.multiply(fdd.add(ds.divide(three).multiply(fddd)))))));
  326.             fd = fd.add(ds.multiply(fdd.add(half.multiply(ds).multiply(fddd))));
  327.             S = stemp.subtract(f.divide(fd));
  328.         }

  329.         final T H = S.asinh();
  330.         return H;
  331.     }

  332.     /**
  333.      * Computes the hyperbolic mean anomaly from the hyperbolic eccentric anomaly.
  334.      * @param <T> field type
  335.      * @param e eccentricity &gt; 1
  336.      * @param H hyperbolic eccentric anomaly
  337.      * @return hyperbolic mean anomaly
  338.      */
  339.     public static <T extends CalculusFieldElement<T>> T hyperbolicEccentricToMean(final T e, final T H) {
  340.         return e.multiply(H.sinh()).subtract(H);
  341.     }

  342. }