IodLaplace.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.estimation.iod;

  18. import org.hipparchus.analysis.solvers.LaguerreSolver;
  19. import org.hipparchus.complex.Complex;
  20. import org.hipparchus.geometry.euclidean.threed.Vector3D;
  21. import org.hipparchus.linear.Array2DRowRealMatrix;
  22. import org.hipparchus.linear.LUDecomposition;
  23. import org.hipparchus.util.FastMath;
  24. import org.orekit.estimation.measurements.AngularRaDec;
  25. import org.orekit.frames.Frame;
  26. import org.orekit.orbits.CartesianOrbit;
  27. import org.orekit.time.AbsoluteDate;
  28. import org.orekit.utils.PVCoordinates;

  29. /**
  30.  * Laplace angles-only initial orbit determination, assuming Keplerian motion.
  31.  * An orbit is determined from three angular observations from the same site.
  32.  *
  33.  *
  34.  * Reference:
  35.  *    Bate, R., Mueller, D. D., & White, J. E. (1971). Fundamentals of astrodynamics.
  36.  *    New York: Dover Publications.
  37.  *
  38.  * @author Shiva Iyer
  39.  * @since 10.1
  40.  */
  41. public class IodLaplace {

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

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

  51.     /** Estimate the orbit from three angular observations at the same location.
  52.      *
  53.      * @param frame inertial frame for observer coordinates and orbit estimate
  54.      * @param obsPva Observer coordinates at time of raDec2
  55.      * @param raDec1 first angular observation
  56.      * @param raDec2 second angular observation
  57.      * @param raDec3 third angular observation
  58.      * @return estimate of the orbit at the central date obsDate2 or null if
  59.      *         no estimate is possible with the given data
  60.      * @since 11.0
  61.      */
  62.     public CartesianOrbit estimate(final Frame frame, final PVCoordinates obsPva,
  63.                                    final AngularRaDec raDec1, final AngularRaDec raDec2,
  64.                                    final AngularRaDec raDec3) {
  65.         return estimate(frame, obsPva,
  66.                         raDec1.getDate(), lineOfSight(raDec1),
  67.                         raDec2.getDate(), lineOfSight(raDec2),
  68.                         raDec3.getDate(), lineOfSight(raDec3));
  69.     }

  70.     /** Estimate orbit from three line of sight angles from the same location.
  71.      *
  72.      * @param frame inertial frame for observer coordinates and orbit estimate
  73.      * @param obsPva Observer coordinates at time obsDate2
  74.      * @param obsDate1 date of observation 1
  75.      * @param los1 line of sight unit vector 1
  76.      * @param obsDate2 date of observation 2
  77.      * @param los2 line of sight unit vector 2
  78.      * @param obsDate3 date of observation 3
  79.      * @param los3 line of sight unit vector 3
  80.      * @return estimate of the orbit at the central date obsDate2 or null if
  81.      *         no estimate is possible with the given data
  82.      */
  83.     public CartesianOrbit estimate(final Frame frame, final PVCoordinates obsPva,
  84.                                    final AbsoluteDate obsDate1, final Vector3D los1,
  85.                                    final AbsoluteDate obsDate2, final Vector3D los2,
  86.                                    final AbsoluteDate obsDate3, final Vector3D los3) {
  87.         // The first observation is taken as t1 = 0
  88.         final double t2 = obsDate2.durationFrom(obsDate1);
  89.         final double t3 = obsDate3.durationFrom(obsDate1);

  90.         // Calculate the first and second derivatives of the Line Of Sight vector at t2
  91.         final Vector3D Ldot = los1.scalarMultiply((t2 - t3) / (t2 * t3)).
  92.                         add(los2.scalarMultiply((2.0 * t2 - t3) / (t2 * (t2 - t3)))).
  93.                         add(los3.scalarMultiply(t2 / (t3 * (t3 - t2))));
  94.         final Vector3D Ldotdot = los1.scalarMultiply(2.0 / (t2 * t3)).
  95.                         add(los2.scalarMultiply(2.0 / (t2 * (t2 - t3)))).
  96.                         add(los3.scalarMultiply(2.0 / (t3 * (t3 - t2))));

  97.         // The determinant will vanish if the observer lies in the plane of the orbit at t2
  98.         final double D = 2.0 * getDeterminant(los2, Ldot, Ldotdot);
  99.         if (FastMath.abs(D) < 1.0E-14) {
  100.             return null;
  101.         }

  102.         final double Dsq = D * D;
  103.         final double R = obsPva.getPosition().getNorm();
  104.         final double RdotL = obsPva.getPosition().dotProduct(los2);

  105.         final double D1 = getDeterminant(los2, Ldot, obsPva.getAcceleration());
  106.         final double D2 = getDeterminant(los2, Ldot, obsPva.getPosition());

  107.         // Coefficients of the 8th order polynomial we need to solve to determine "r"
  108.         final double[] coeff = new double[] {-4.0 * mu * mu * D2 * D2 / Dsq,
  109.                                              0.0,
  110.                                              0.0,
  111.                                              4.0 * mu * D2 * (RdotL / D - 2.0 * D1 / Dsq),
  112.                                              0.0,
  113.                                              0.0,
  114.                                              4.0 * D1 * RdotL / D - 4.0 * D1 * D1 / Dsq - R * R, 0.0,
  115.                                              1.0};

  116.         // Use the Laguerre polynomial solver and take the initial guess to be
  117.         // 5 times the observer's position magnitude
  118.         final LaguerreSolver solver = new LaguerreSolver(1E-10, 1E-10, 1E-10);
  119.         final Complex[] roots = solver.solveAllComplex(coeff, 5.0 * R);

  120.         // We consider "r" to be the positive real root with the largest magnitude
  121.         double rMag = 0.0;
  122.         for (int i = 0; i < roots.length; i++) {
  123.             if (roots[i].getReal() > rMag &&
  124.                             FastMath.abs(roots[i].getImaginary()) < solver.getAbsoluteAccuracy()) {
  125.                 rMag = roots[i].getReal();
  126.             }
  127.         }
  128.         if (rMag == 0.0) {
  129.             return null;
  130.         }

  131.         // Calculate rho, the slant range from the observer to the satellite at t2.
  132.         // This yields the "r" vector, which is the satellite's position vector at t2.
  133.         final double rCubed = rMag * rMag * rMag;
  134.         final double rho = -2.0 * D1 / D - 2.0 * mu * D2 / (D * rCubed);
  135.         final Vector3D posVec = los2.scalarMultiply(rho).add(obsPva.getPosition());

  136.         // Calculate rho_dot at t2, which will yield the satellite's velocity vector at t2
  137.         final double D3 = getDeterminant(los2, obsPva.getAcceleration(), Ldotdot);
  138.         final double D4 = getDeterminant(los2, obsPva.getPosition(), Ldotdot);
  139.         final double rhoDot = -D3 / D - mu * D4 / (D * rCubed);
  140.         final Vector3D velVec = los2.scalarMultiply(rhoDot).
  141.                         add(Ldot.scalarMultiply(rho)).
  142.                         add(obsPva.getVelocity());

  143.         // Return the estimated orbit
  144.         return new CartesianOrbit(new PVCoordinates(posVec, velVec), frame, obsDate2, mu);
  145.     }

  146.     /**
  147.      * Calculates the line of sight vector.
  148.      * @param alpha right ascension angle, in radians
  149.      * @param delta declination angle, in radians
  150.      * @return the line of sight vector
  151.      * @since 11.0
  152.      */
  153.     public static Vector3D lineOfSight(final double alpha, final double delta) {
  154.         return new Vector3D(FastMath.cos(delta) * FastMath.cos(alpha),
  155.                             FastMath.cos(delta) * FastMath.sin(alpha),
  156.                             FastMath.sin(delta));
  157.     }

  158.     /**
  159.      * Calculate the line of sight vector from an AngularRaDec measurement.
  160.      * @param raDec measurement
  161.      * @return the line of sight vector
  162.      * @since 11.0
  163.      */
  164.     public static Vector3D lineOfSight(final AngularRaDec raDec) {

  165.         // Observed values
  166.         final double[] observed = raDec.getObservedValue();

  167.         // Return
  168.         return lineOfSight(observed[0], observed[1]);

  169.     }

  170.     /** Calculate the determinant of the matrix with given column vectors.
  171.      *
  172.      * @param col0 Matrix column 0
  173.      * @param col1 Matrix column 1
  174.      * @param col2 Matrix column 2
  175.      * @return matrix determinant
  176.      *
  177.      */
  178.     private double getDeterminant(final Vector3D col0, final Vector3D col1, final Vector3D col2) {
  179.         final Array2DRowRealMatrix mat = new Array2DRowRealMatrix(3, 3);
  180.         mat.setColumn(0, col0.toArray());
  181.         mat.setColumn(1, col1.toArray());
  182.         mat.setColumn(2, col2.toArray());
  183.         return new LUDecomposition(mat).getDeterminant();
  184.     }

  185. }