IodLambert.java

  1. /* Copyright 2002-2020 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.estimation.iod;

  18. import org.hipparchus.geometry.euclidean.threed.Vector3D;
  19. import org.hipparchus.util.FastMath;
  20. import org.orekit.errors.OrekitException;
  21. import org.orekit.errors.OrekitMessages;
  22. import org.orekit.frames.Frame;
  23. import org.orekit.orbits.KeplerianOrbit;
  24. import org.orekit.time.AbsoluteDate;
  25. import org.orekit.utils.PVCoordinates;

  26. /**
  27.  * Lambert initial orbit determination, assuming Keplerian motion.
  28.  * An orbit is determined from two position vectors.
  29.  *
  30.  * References:
  31.  *  Battin, R.H., An Introduction to the Mathematics and Methods of Astrodynamics, AIAA Education, 1999.
  32.  *  Lancaster, E.R. and Blanchard, R.C., A Unified Form of Lambert’s Theorem, Goddard Space Flight Center, 1968.
  33.  *
  34.  * @author Joris Olympio
  35.  * @since 8.0
  36.  */
  37. public class IodLambert {

  38.     /** gravitational constant. */
  39.     private final double mu;

  40.     /** Creator.
  41.      *
  42.      * @param mu gravitational constant
  43.      */
  44.     public IodLambert(final double mu) {
  45.         this.mu = mu;
  46.     }

  47.     /** Estimate a Keplerian orbit given two position vectors and a duration.
  48.      * <p>
  49.      * The logic for setting {@code posigrade} and {@code nRev} is that the
  50.      * sweep angle Δυ travelled by the object between {@code t1} and {@code t2} is
  51.      * 2π {@code nRev +1} - α if {@code posigrade} is false and 2π {@code nRev} + α
  52.      * if {@code posigrade} is true, where α is the separation angle between
  53.      * {@code p1} and {@code p2}, which is always computed between 0 and π
  54.      * (because in 3D without a normal reference, vector angles cannot go past π).
  55.      * </p>
  56.      * <p>
  57.      * This implies that {@code posigrade} should be set to true if {@code p2} is
  58.      * located in the half orbit starting at {@code p1} and it should be set to
  59.      * false if {@code p2} is located in the half orbit ending at {@code p1},
  60.      * regardless of the number of periods between {@code t1} and {@code t2},
  61.      * and {@code nRev} should be set accordingly.
  62.      * </p>
  63.      * <p>
  64.      * As an example, if {@code t2} is less than half a period after {@code t1},
  65.      * then {@code posigrade} should be {@code true} and {@code nRev} should be 0.
  66.      * If {@code t2} is more than half a period after {@code t1} but less than
  67.      * one period after {@code t1}, {@code posigrade} should be {@code false} and
  68.      * {@code nRev} should be 0.
  69.      * </p>
  70.      * @param frame     frame
  71.      * @param posigrade flag indicating the direction of motion
  72.      * @param nRev      number of revolutions
  73.      * @param p1        position vector 1
  74.      * @param t1        date of observation 1
  75.      * @param p2        position vector 2
  76.      * @param t2        date of observation 2
  77.      * @return  an initial Keplerian orbit estimate
  78.      */
  79.     public KeplerianOrbit estimate(final Frame frame, final boolean posigrade,
  80.                                    final int nRev,
  81.                                    final Vector3D p1, final AbsoluteDate t1,
  82.                                    final Vector3D p2, final AbsoluteDate t2) {

  83.         final double r1 = p1.getNorm();
  84.         final double r2 = p2.getNorm();
  85.         final double tau = t2.durationFrom(t1); // in seconds

  86.         // Exception if t2 < t1
  87.         if (tau < 0.0) {
  88.             throw new OrekitException(OrekitMessages.NON_CHRONOLOGICAL_DATES_FOR_OBSERVATIONS, t1, t2);
  89.         }

  90.         // normalizing constants
  91.         final double R = FastMath.max(r1, r2); // in m
  92.         final double V = FastMath.sqrt(mu / R);  // in m/s
  93.         final double T = R / V; // in seconds

  94.         // sweep angle
  95.         double dth = Vector3D.angle(p1, p2);
  96.         // compute the number of revolutions
  97.         if (!posigrade) {
  98.             dth = 2 * FastMath.PI - dth;
  99.         }

  100.         // velocity vectors in the orbital plane, in the R-T frame
  101.         final double[] Vdep = new double[2];

  102.         // call Lambert's problem solver
  103.         final boolean exitflag = solveLambertPb(r1 / R, r2 / R, dth, tau / T, nRev, Vdep);

  104.         if (exitflag) {
  105.             // basis vectors
  106.             // normal to the orbital arc plane
  107.             final Vector3D Pn = p1.crossProduct(p2);
  108.             // perpendicular to the radius vector, in the orbital arc plane
  109.             final Vector3D Pt = Pn.crossProduct(p1);

  110.             // tangential velocity norm
  111.             double RT = Pt.getNorm();
  112.             if (!posigrade) {
  113.                 RT = -RT;
  114.             }

  115.             // velocity vector at P1
  116.             final Vector3D Vel1 = new Vector3D(V * Vdep[0] / r1, p1,
  117.                                                V * Vdep[1] / RT, Pt);

  118.             // compute the equivalent Keplerian orbit
  119.             return new KeplerianOrbit(new PVCoordinates(p1, Vel1), frame, t1, mu);
  120.         }

  121.         return null;
  122.     }

  123.     /** Lambert's solver.
  124.      * Assume mu=1.
  125.      *
  126.      * @param r1 radius 1
  127.      * @param r2  radius 2
  128.      * @param dth sweep angle
  129.      * @param tau time of flight
  130.      * @param mRev number of revs
  131.      * @param V1 velocity at departure in (T, N) basis
  132.      * @return something
  133.      */
  134.     boolean solveLambertPb(final double r1, final double r2, final double dth, final double tau,
  135.                            final int mRev, final double[] V1) {
  136.         // decide whether to use the left or right branch (for multi-revolution
  137.         // problems), and the long- or short way.
  138.         final boolean leftbranch = dth < FastMath.PI;
  139.         int longway = 1;
  140.         if (dth > FastMath.PI) {
  141.             longway = -1;
  142.         }

  143.         final int m = FastMath.abs(mRev);
  144.         final double rtof = FastMath.abs(tau);
  145.         final double theta = dth;

  146.         // non-dimensional chord ||r2-r1||
  147.         final double chord = FastMath.sqrt(r1 * r1 + r2 * r2 - 2 * r1 * r2 * FastMath.cos(theta));

  148.         // non-dimensional semi-perimeter of the triangle
  149.         final double speri = 0.5 * (r1 + r2 + chord);

  150.         // minimum energy ellipse semi-major axis
  151.         final double minSma = speri / 2.;

  152.         // lambda parameter (Eq 7.6)
  153.         final double lambda = longway * FastMath.sqrt(1 - chord / speri);

  154.         // reference tof value for the Newton solver
  155.         final double logt = FastMath.log(rtof);

  156.         // initialisation of the iterative root finding process (secant method)
  157.         // initial values
  158.         //  -1 < x < 1  =>  Elliptical orbits
  159.         //  x = 1           Parabolic orbit
  160.         //  x > 1           Hyperbolic orbits
  161.         final double in1;
  162.         final double in2;
  163.         double x1;
  164.         double x2;
  165.         if (m == 0) {
  166.             // one revolution, one solution. Define the left and right asymptotes.
  167.             in1 = -0.6523333;
  168.             in2 = 0.6523333;
  169.             x1   = FastMath.log(1 + in1);
  170.             x2   = FastMath.log(1 + in2);
  171.         } else {
  172.             // select initial values, depending on the branch
  173.             if (!leftbranch) {
  174.                 in1 = -0.523334;
  175.                 in2 = -0.223334;
  176.             } else {
  177.                 in1 = 0.723334;
  178.                 in2 = 0.523334;
  179.             }
  180.             x1 = FastMath.tan(in1 * 0.5 * FastMath.PI);
  181.             x2 = FastMath.tan(in2 * 0.5 * FastMath.PI);
  182.         }

  183.         // initial estimates for the tof
  184.         final double tof1 = timeOfFlight(in1, longway, m, minSma, speri, chord);
  185.         final double tof2 = timeOfFlight(in2, longway, m, minSma, speri, chord);

  186.         // initial bounds for y
  187.         double y1;
  188.         double y2;
  189.         if (m == 0) {
  190.             y1 = FastMath.log(tof1) - logt;
  191.             y2 = FastMath.log(tof2) - logt;
  192.         } else {
  193.             y1 = tof1 - rtof;
  194.             y2 = tof2 - rtof;
  195.         }

  196.         // Solve for x with the secant method
  197.         double err = 1e20;
  198.         int iterations = 0;
  199.         final double tol = 1e-13;
  200.         final int maxiter = 50;
  201.         double xnew = 0;
  202.         while ((err > tol) && (iterations < maxiter)) {
  203.             // new x
  204.             xnew = (x1 * y2 - y1 * x2) / (y2 - y1);

  205.             // evaluate new time of flight
  206.             final double x;
  207.             if (m == 0) {
  208.                 x = FastMath.exp(xnew) - 1;
  209.             } else {
  210.                 x = FastMath.atan(xnew) * 2 / FastMath.PI;
  211.             }

  212.             final double tof = timeOfFlight(x, longway, m, minSma, speri, chord);

  213.             // new value of y
  214.             final double ynew;
  215.             if (m == 0) {
  216.                 ynew = FastMath.log(tof) - logt;
  217.             } else {
  218.                 ynew = tof - rtof;
  219.             }

  220.             // save previous and current values for the next iteration
  221.             x1 = x2;
  222.             x2 = xnew;
  223.             y1 = y2;
  224.             y2 = ynew;

  225.             // update error
  226.             err = FastMath.abs(x1 - xnew);

  227.             // increment number of iterations
  228.             ++iterations;
  229.         }

  230.         // failure to converge
  231.         if (err > tol) {
  232.             return false;
  233.         }

  234.         // convert converged value of x
  235.         final double x;
  236.         if (m == 0) {
  237.             x = FastMath.exp(xnew) - 1;
  238.         } else {
  239.             x = FastMath.atan(xnew) * 2 / FastMath.PI;
  240.         }

  241.         // Solution for the semi-major axis (Eq. 7.20)
  242.         final double sma = minSma / (1 - x * x);

  243.         // compute velocities
  244.         final double eta;
  245.         if (x < 1) {
  246.             // ellipse, Eqs. 7.7, 7.17
  247.             final double alfa = 2 * FastMath.acos(x);
  248.             final double beta = longway * 2 * FastMath.asin(FastMath.sqrt((speri - chord) / (2. * sma)));
  249.             final double psi  = (alfa - beta) / 2;
  250.             // Eq. 7.21
  251.             final double sinPsi = FastMath.sin(psi);
  252.             final double etaSq = 2 * sma * sinPsi * sinPsi / speri;
  253.             eta  = FastMath.sqrt(etaSq);
  254.         } else {
  255.             // hyperbola
  256.             final double gamma = 2 * FastMath.acosh(x);
  257.             final double delta = longway * 2 * FastMath.asinh(FastMath.sqrt((chord - speri) / (2 * sma)));
  258.             //
  259.             final double psi  = (gamma - delta) / 2.;
  260.             final double sinhPsi = FastMath.sinh(psi);
  261.             final double etaSq = -2 * sma * sinhPsi * sinhPsi / speri;
  262.             eta  = FastMath.sqrt(etaSq);
  263.         }

  264.         // radial and tangential directions for departure velocity (Eq. 7.36)
  265.         final double VR1 = (1. / eta) * FastMath.sqrt(1. / minSma) * (2 * lambda * minSma / r1 - (lambda + x * eta));
  266.         final double VT1 = (1. / eta) * FastMath.sqrt(1. / minSma) * FastMath.sqrt(r2 / r1) * FastMath.sin(dth / 2);
  267.         V1[0] = VR1;
  268.         V1[1] = VT1;

  269.         return true;
  270.     }

  271.     /** Compute the time of flight of a given arc of orbit.
  272.      * The time of flight is evaluated via the Lagrange expression.
  273.      *
  274.      * @param x          x
  275.      * @param longway    solution number; the long way or the short war
  276.      * @param mrev       number of revolutions of the arc of orbit
  277.      * @param minSma     minimum possible semi-major axis
  278.      * @param speri      semi-parameter of the arc of orbit
  279.      * @param chord      chord of the arc of orbit
  280.      * @return the time of flight for the given arc of orbit
  281.      */
  282.     private double timeOfFlight(final double x, final int longway, final int mrev, final double minSma,
  283.                                 final double speri, final double chord) {

  284.         final double a = minSma / (1 - x * x);

  285.         final double tof;
  286.         if (FastMath.abs(x) < 1) {
  287.             // Lagrange form of the time of flight equation Eq. (7.9)
  288.             // elliptical orbit (note: mu = 1)
  289.             final double beta = longway * 2 * FastMath.asin(FastMath.sqrt((speri - chord) / (2. * a)));
  290.             final double alfa = 2 * FastMath.acos(x);
  291.             tof = a * FastMath.sqrt(a) * ((alfa - FastMath.sin(alfa)) - (beta - FastMath.sin(beta)) + 2 * FastMath.PI * mrev);
  292.         } else {
  293.             // hyperbolic orbit
  294.             final double alfa = 2 * FastMath.acosh(x);
  295.             final double beta = longway * 2 * FastMath.asinh(FastMath.sqrt((speri - chord) / (-2. * a)));
  296.             tof = -a * FastMath.sqrt(-a) * ((FastMath.sinh(alfa) - alfa) - (FastMath.sinh(beta) - beta));
  297.         }

  298.         return tof;
  299.     }
  300. }