IodLambert.java

  1. /* Copyright 2002-2025 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.estimation.measurements.PV;
  23. import org.orekit.estimation.measurements.Position;
  24. import org.orekit.frames.Frame;
  25. import org.orekit.orbits.CartesianOrbit;
  26. import org.orekit.orbits.Orbit;
  27. import org.orekit.time.AbsoluteDate;
  28. import org.orekit.utils.PVCoordinates;

  29. /**
  30.  * Lambert position-based Initial Orbit Determination (IOD) algorithm, assuming Keplerian motion.
  31.  * <p>
  32.  * An orbit is determined from two position vectors.
  33.  *
  34.  * References:
  35.  *  Battin, R.H., An Introduction to the Mathematics and Methods of Astrodynamics, AIAA Education, 1999.
  36.  *  Lancaster, E.R. and Blanchard, R.C., A Unified Form of Lambert’s Theorem, Goddard Space Flight Center, 1968.
  37.  * </p>
  38.  * @author Joris Olympio
  39.  * @since 8.0
  40.  */
  41. public class IodLambert {

  42.     /** gravitational constant. */
  43.     private final double mu;

  44.     /** Creator.
  45.      *
  46.      * @param mu gravitational constant
  47.      */
  48.     public IodLambert(final double mu) {
  49.         this.mu = mu;
  50.     }

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

  87.     /** Estimate an initial orbit from two PV measurements.
  88.      * <p>
  89.      * The logic for setting {@code posigrade} and {@code nRev} is that the
  90.      * sweep angle Δυ travelled by the object between {@code t1} and {@code t2} is
  91.      * 2π {@code nRev +1} - α if {@code posigrade} is false and 2π {@code nRev} + α
  92.      * if {@code posigrade} is true, where α is the separation angle between
  93.      * {@code p1} and {@code p2}, which is always computed between 0 and π
  94.      * (because in 3D without a normal reference, vector angles cannot go past π).
  95.      * </p>
  96.      * <p>
  97.      * This implies that {@code posigrade} should be set to true if {@code p2} is
  98.      * located in the half orbit starting at {@code p1} and it should be set to
  99.      * false if {@code p2} is located in the half orbit ending at {@code p1},
  100.      * regardless of the number of periods between {@code t1} and {@code t2},
  101.      * and {@code nRev} should be set accordingly.
  102.      * </p>
  103.      * <p>
  104.      * As an example, if {@code t2} is less than half a period after {@code t1},
  105.      * then {@code posigrade} should be {@code true} and {@code nRev} should be 0.
  106.      * If {@code t2} is more than half a period after {@code t1} but less than
  107.      * one period after {@code t1}, {@code posigrade} should be {@code false} and
  108.      * {@code nRev} should be 0.
  109.      * </p>
  110.      * @param frame     measurements frame
  111.      * @param posigrade flag indicating the direction of motion
  112.      * @param nRev      number of revolutions
  113.      * @param pv1       first PV measurement
  114.      * @param pv2       second PV measurement
  115.      * @return an initial Keplerian orbit estimation at the first observation date t1
  116.      * @since 12.0
  117.      */
  118.     public Orbit estimate(final Frame frame, final boolean posigrade,
  119.                           final int nRev, final PV pv1,  final PV pv2) {
  120.         return estimate(frame, posigrade, nRev,
  121.                         pv1.getPosition(), pv1.getDate(), pv2.getPosition(), pv2.getDate());
  122.     }

  123.     /** Estimate a Keplerian orbit given two position vectors and a duration.
  124.      * <p>
  125.      * The logic for setting {@code posigrade} and {@code nRev} is that the
  126.      * sweep angle Δυ travelled by the object between {@code t1} and {@code t2} is
  127.      * 2π {@code nRev +1} - α if {@code posigrade} is false and 2π {@code nRev} + α
  128.      * if {@code posigrade} is true, where α is the separation angle between
  129.      * {@code p1} and {@code p2}, which is always computed between 0 and π
  130.      * (because in 3D without a normal reference, vector angles cannot go past π).
  131.      * </p>
  132.      * <p>
  133.      * This implies that {@code posigrade} should be set to true if {@code p2} is
  134.      * located in the half orbit starting at {@code p1} and it should be set to
  135.      * false if {@code p2} is located in the half orbit ending at {@code p1},
  136.      * regardless of the number of periods between {@code t1} and {@code t2},
  137.      * and {@code nRev} should be set accordingly.
  138.      * </p>
  139.      * <p>
  140.      * As an example, if {@code t2} is less than half a period after {@code t1},
  141.      * then {@code posigrade} should be {@code true} and {@code nRev} should be 0.
  142.      * If {@code t2} is more than half a period after {@code t1} but less than
  143.      * one period after {@code t1}, {@code posigrade} should be {@code false} and
  144.      * {@code nRev} should be 0.
  145.      * </p>
  146.      * @param frame     frame
  147.      * @param posigrade flag indicating the direction of motion
  148.      * @param nRev      number of revolutions
  149.      * @param p1        position vector 1
  150.      * @param t1        date of observation 1
  151.      * @param p2        position vector 2
  152.      * @param t2        date of observation 2
  153.      * @return  an initial Keplerian orbit estimate at the first observation date t1
  154.      */
  155.     public Orbit estimate(final Frame frame, final boolean posigrade,
  156.                           final int nRev,
  157.                           final Vector3D p1, final AbsoluteDate t1,
  158.                           final Vector3D p2, final AbsoluteDate t2) {

  159.         final double r1 = p1.getNorm();
  160.         final double r2 = p2.getNorm();
  161.         final double tau = t2.durationFrom(t1); // in seconds

  162.         // Exception if t2 < t1
  163.         if (tau < 0.0) {
  164.             throw new OrekitException(OrekitMessages.NON_CHRONOLOGICAL_DATES_FOR_OBSERVATIONS, t1, t2, -tau);
  165.         }

  166.         // normalizing constants
  167.         final double R = FastMath.max(r1, r2); // in m
  168.         final double V = FastMath.sqrt(mu / R);  // in m/s
  169.         final double T = R / V; // in seconds

  170.         // sweep angle
  171.         double dth = Vector3D.angle(p1, p2);
  172.         // compute the number of revolutions
  173.         if (!posigrade) {
  174.             dth = 2 * FastMath.PI - dth;
  175.         }

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

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

  180.         if (exitflag) {
  181.             // basis vectors
  182.             // normal to the orbital arc plane
  183.             final Vector3D Pn = p1.crossProduct(p2);
  184.             // perpendicular to the radius vector, in the orbital arc plane
  185.             final Vector3D Pt = Pn.crossProduct(p1);

  186.             // tangential velocity norm
  187.             double RT = Pt.getNorm();
  188.             if (!posigrade) {
  189.                 RT = -RT;
  190.             }

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

  194.             // compute the equivalent Cartesian orbit
  195.             return new CartesianOrbit(new PVCoordinates(p1, Vel1), frame, t1, mu);
  196.         }

  197.         return null;
  198.     }

  199.     /** Lambert's solver.
  200.      * Assume mu=1.
  201.      *
  202.      * @param r1 radius 1
  203.      * @param r2  radius 2
  204.      * @param dth sweep angle
  205.      * @param tau time of flight
  206.      * @param mRev number of revs
  207.      * @param V1 velocity at departure in (T, N) basis
  208.      * @return something
  209.      */
  210.     boolean solveLambertPb(final double r1, final double r2, final double dth, final double tau,
  211.                            final int mRev, final double[] V1) {
  212.         // decide whether to use the left or right branch (for multi-revolution
  213.         // problems), and the long- or short way.
  214.         final boolean leftbranch = dth < FastMath.PI;
  215.         int longway = 1;
  216.         if (dth > FastMath.PI) {
  217.             longway = -1;
  218.         }

  219.         final int m = FastMath.abs(mRev);
  220.         final double rtof = FastMath.abs(tau);

  221.         // non-dimensional chord ||r2-r1||
  222.         final double chord = FastMath.sqrt(r1 * r1 + r2 * r2 - 2 * r1 * r2 * FastMath.cos(dth));

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

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

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

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

  231.         // initialisation of the iterative root finding process (secant method)
  232.         // initial values
  233.         //  -1 < x < 1  =>  Elliptical orbits
  234.         //  x = 1           Parabolic orbit
  235.         //  x > 1           Hyperbolic orbits
  236.         final double in1;
  237.         final double in2;
  238.         double x1;
  239.         double x2;
  240.         if (m == 0) {
  241.             // one revolution, one solution. Define the left and right asymptotes.
  242.             in1 = -0.6523333;
  243.             in2 = 0.6523333;
  244.             x1   = FastMath.log(1 + in1);
  245.             x2   = FastMath.log(1 + in2);
  246.         } else {
  247.             // select initial values, depending on the branch
  248.             if (!leftbranch) {
  249.                 in1 = -0.523334;
  250.                 in2 = -0.223334;
  251.             } else {
  252.                 in1 = 0.723334;
  253.                 in2 = 0.523334;
  254.             }
  255.             x1 = FastMath.tan(in1 * 0.5 * FastMath.PI);
  256.             x2 = FastMath.tan(in2 * 0.5 * FastMath.PI);
  257.         }

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

  261.         // initial bounds for y
  262.         double y1;
  263.         double y2;
  264.         if (m == 0) {
  265.             y1 = FastMath.log(tof1) - logt;
  266.             y2 = FastMath.log(tof2) - logt;
  267.         } else {
  268.             y1 = tof1 - rtof;
  269.             y2 = tof2 - rtof;
  270.         }

  271.         // Solve for x with the secant method
  272.         double err = 1e20;
  273.         int iterations = 0;
  274.         final double tol = 1e-13;
  275.         final int maxiter = 50;
  276.         double xnew = 0;
  277.         while (err > tol && iterations < maxiter) {
  278.             // new x
  279.             xnew = (x1 * y2 - y1 * x2) / (y2 - y1);

  280.             // evaluate new time of flight
  281.             final double x;
  282.             if (m == 0) {
  283.                 x = FastMath.exp(xnew) - 1;
  284.             } else {
  285.                 x = FastMath.atan(xnew) * 2 / FastMath.PI;
  286.             }

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

  288.             // new value of y
  289.             final double ynew;
  290.             if (m == 0) {
  291.                 ynew = FastMath.log(tof) - logt;
  292.             } else {
  293.                 ynew = tof - rtof;
  294.             }

  295.             // save previous and current values for the next iteration
  296.             x1 = x2;
  297.             x2 = xnew;
  298.             y1 = y2;
  299.             y2 = ynew;

  300.             // update error
  301.             err = FastMath.abs(x1 - xnew);

  302.             // increment number of iterations
  303.             ++iterations;
  304.         }

  305.         // failure to converge
  306.         if (err > tol) {
  307.             return false;
  308.         }

  309.         // convert converged value of x
  310.         final double x;
  311.         if (m == 0) {
  312.             x = FastMath.exp(xnew) - 1;
  313.         } else {
  314.             x = FastMath.atan(xnew) * 2 / FastMath.PI;
  315.         }

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

  318.         // compute velocities
  319.         final double eta;
  320.         if (x < 1) {
  321.             // ellipse, Eqs. 7.7, 7.17
  322.             final double alfa = 2 * FastMath.acos(x);
  323.             final double beta = longway * 2 * FastMath.asin(FastMath.sqrt((speri - chord) / (2. * sma)));
  324.             final double psi  = (alfa - beta) / 2;
  325.             // Eq. 7.21
  326.             final double sinPsi = FastMath.sin(psi);
  327.             final double etaSq = 2 * sma * sinPsi * sinPsi / speri;
  328.             eta  = FastMath.sqrt(etaSq);
  329.         } else {
  330.             // hyperbola
  331.             final double gamma = 2 * FastMath.acosh(x);
  332.             final double delta = longway * 2 * FastMath.asinh(FastMath.sqrt((chord - speri) / (2 * sma)));
  333.             //
  334.             final double psi  = (gamma - delta) / 2.;
  335.             final double sinhPsi = FastMath.sinh(psi);
  336.             final double etaSq = -2 * sma * sinhPsi * sinhPsi / speri;
  337.             eta  = FastMath.sqrt(etaSq);
  338.         }

  339.         // radial and tangential directions for departure velocity (Eq. 7.36)
  340.         final double VR1 = (1. / eta) * FastMath.sqrt(1. / minSma) * (2 * lambda * minSma / r1 - (lambda + x * eta));
  341.         final double VT1 = (1. / eta) * FastMath.sqrt(1. / minSma) * FastMath.sqrt(r2 / r1) * FastMath.sin(dth / 2);
  342.         V1[0] = VR1;
  343.         V1[1] = VT1;

  344.         return true;
  345.     }

  346.     /** Compute the time of flight of a given arc of orbit.
  347.      * The time of flight is evaluated via the Lagrange expression.
  348.      *
  349.      * @param x          x
  350.      * @param longway    solution number; the long way or the short war
  351.      * @param mrev       number of revolutions of the arc of orbit
  352.      * @param minSma     minimum possible semi-major axis
  353.      * @param speri      semi-parameter of the arc of orbit
  354.      * @param chord      chord of the arc of orbit
  355.      * @return the time of flight for the given arc of orbit
  356.      */
  357.     private double timeOfFlight(final double x, final int longway, final int mrev, final double minSma,
  358.                                 final double speri, final double chord) {

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

  360.         final double tof;
  361.         if (FastMath.abs(x) < 1) {
  362.             // Lagrange form of the time of flight equation Eq. (7.9)
  363.             // elliptical orbit (note: mu = 1)
  364.             final double beta = longway * 2 * FastMath.asin(FastMath.sqrt((speri - chord) / (2. * a)));
  365.             final double alfa = 2 * FastMath.acos(x);
  366.             tof = a * FastMath.sqrt(a) * ((alfa - FastMath.sin(alfa)) - (beta - FastMath.sin(beta)) + 2 * FastMath.PI * mrev);
  367.         } else {
  368.             // hyperbolic orbit
  369.             final double alfa = 2 * FastMath.acosh(x);
  370.             final double beta = longway * 2 * FastMath.asinh(FastMath.sqrt((speri - chord) / (-2. * a)));
  371.             tof = -a * FastMath.sqrt(-a) * ((FastMath.sinh(alfa) - alfa) - (FastMath.sinh(beta) - beta));
  372.         }

  373.         return tof;
  374.     }
  375. }