Class FieldKeplerianOrbit<T extends org.hipparchus.RealFieldElement<T>>

  • All Implemented Interfaces:
    FieldTimeInterpolable<FieldOrbit<T>,​T>, FieldTimeShiftable<FieldOrbit<T>,​T>, FieldTimeStamped<T>, FieldPVCoordinatesProvider<T>

    public class FieldKeplerianOrbit<T extends org.hipparchus.RealFieldElement<T>>
    extends FieldOrbit<T>
    This class handles traditional Keplerian orbital parameters.

    The parameters used internally are the classical Keplerian elements:

         a
         e
         i
         ω
         Ω
         v
       
    where ω stands for the Perigee Argument, Ω stands for the Right Ascension of the Ascending Node and v stands for the true anomaly.

    This class supports hyperbolic orbits, using the convention that semi major axis is negative for such orbits (and of course eccentricity is greater than 1).

    When orbit is either equatorial or circular, some Keplerian elements (more precisely ω and Ω) become ambiguous so this class should not be used for such orbits. For this reason, equinoctial orbits is the recommended way to represent orbits.

    The instance KeplerianOrbit is guaranteed to be immutable.

    Since:
    9.0
    Author:
    Luc Maisonobe, Guylaine Prat, Fabien Maussion, Véronique Pommier-Maurussane, Andrea Antolino
    See Also:
    Orbit, CircularOrbit, CartesianOrbit, EquinoctialOrbit
    • Method Detail

      • getType

        public OrbitType getType()
        Get the orbit type.
        Specified by:
        getType in class FieldOrbit<T extends org.hipparchus.RealFieldElement<T>>
        Returns:
        orbit type
      • getA

        public T getA()
        Get the semi-major axis.

        Note that the semi-major axis is considered negative for hyperbolic orbits.

        Specified by:
        getA in class FieldOrbit<T extends org.hipparchus.RealFieldElement<T>>
        Returns:
        semi-major axis (m)
      • getADot

        public T getADot()
        Get the semi-major axis derivative.

        Note that the semi-major axis is considered negative for hyperbolic orbits.

        If the orbit was created without derivatives, the value returned is null.

        Specified by:
        getADot in class FieldOrbit<T extends org.hipparchus.RealFieldElement<T>>
        Returns:
        semi-major axis derivative (m/s)
      • getE

        public T getE()
        Get the eccentricity.
        Specified by:
        getE in class FieldOrbit<T extends org.hipparchus.RealFieldElement<T>>
        Returns:
        eccentricity
      • getEDot

        public T getEDot()
        Get the eccentricity derivative.

        If the orbit was created without derivatives, the value returned is null.

        Specified by:
        getEDot in class FieldOrbit<T extends org.hipparchus.RealFieldElement<T>>
        Returns:
        eccentricity derivative
      • getI

        public T getI()
        Get the inclination.

        If the orbit was created without derivatives, the value returned is null.

        Specified by:
        getI in class FieldOrbit<T extends org.hipparchus.RealFieldElement<T>>
        Returns:
        inclination (rad)
      • getIDot

        public T getIDot()
        Get the inclination derivative.
        Specified by:
        getIDot in class FieldOrbit<T extends org.hipparchus.RealFieldElement<T>>
        Returns:
        inclination derivative (rad/s)
      • getPerigeeArgument

        public T getPerigeeArgument()
        Get the perigee argument.
        Returns:
        perigee argument (rad)
      • getPerigeeArgumentDot

        public T getPerigeeArgumentDot()
        Get the perigee argument derivative.

        If the orbit was created without derivatives, the value returned is null.

        Returns:
        perigee argument derivative (rad/s)
      • getRightAscensionOfAscendingNode

        public T getRightAscensionOfAscendingNode()
        Get the right ascension of the ascending node.
        Returns:
        right ascension of the ascending node (rad)
      • getRightAscensionOfAscendingNodeDot

        public T getRightAscensionOfAscendingNodeDot()
        Get the right ascension of the ascending node derivative.

        If the orbit was created without derivatives, the value returned is null.

        Returns:
        right ascension of the ascending node derivative (rad/s)
      • getTrueAnomaly

        public T getTrueAnomaly()
        Get the true anomaly.
        Returns:
        true anomaly (rad)
      • getTrueAnomalyDot

        public T getTrueAnomalyDot()
        Get the true anomaly derivative.

        If the orbit was created without derivatives, the value returned is null.

        Returns:
        true anomaly derivative (rad/s)
      • getEccentricAnomaly

        public T getEccentricAnomaly()
        Get the eccentric anomaly.
        Returns:
        eccentric anomaly (rad)
      • getEccentricAnomalyDot

        public T getEccentricAnomalyDot()
        Get the eccentric anomaly derivative.

        If the orbit was created without derivatives, the value returned is null.

        Returns:
        eccentric anomaly derivative (rad/s)
      • getMeanAnomaly

        public T getMeanAnomaly()
        Get the mean anomaly.
        Returns:
        mean anomaly (rad)
      • getMeanAnomalyDot

        public T getMeanAnomalyDot()
        Get the mean anomaly derivative.

        If the orbit was created without derivatives, the value returned is null.

        Returns:
        mean anomaly derivative (rad/s)
      • getAnomaly

        public T getAnomaly​(PositionAngle type)
        Get the anomaly.
        Parameters:
        type - type of the angle
        Returns:
        anomaly (rad)
      • getAnomalyDot

        public T getAnomalyDot​(PositionAngle type)
        Get the anomaly derivative.

        If the orbit was created without derivatives, the value returned is null.

        Parameters:
        type - type of the angle
        Returns:
        anomaly derivative (rad/s)
      • ellipticEccentricToTrue

        public static <T extends org.hipparchus.RealFieldElement<T>> T ellipticEccentricToTrue​(T E,
                                                                                               T e)
        Computes the true anomaly from the elliptic eccentric anomaly.
        Type Parameters:
        T - type of the field elements
        Parameters:
        E - eccentric anomaly (rad)
        e - eccentricity
        Returns:
        v the true anomaly
      • trueToEllipticEccentric

        public static <T extends org.hipparchus.RealFieldElement<T>> T trueToEllipticEccentric​(T v,
                                                                                               T e)
        Computes the elliptic eccentric anomaly from the true anomaly.
        Type Parameters:
        T - type of the field elements
        Parameters:
        v - true anomaly (rad)
        e - eccentricity
        Returns:
        E the elliptic eccentric anomaly
      • hyperbolicEccentricToTrue

        public static <T extends org.hipparchus.RealFieldElement<T>> T hyperbolicEccentricToTrue​(T H,
                                                                                                 T e)
        Computes the true anomaly from the hyperbolic eccentric anomaly.
        Type Parameters:
        T - type of the field elements
        Parameters:
        H - hyperbolic eccentric anomaly (rad)
        e - eccentricity
        Returns:
        v the true anomaly
      • trueToHyperbolicEccentric

        public static <T extends org.hipparchus.RealFieldElement<T>> T trueToHyperbolicEccentric​(T v,
                                                                                                 T e)
        Computes the hyperbolic eccentric anomaly from the true anomaly.
        Type Parameters:
        T - type of the field elements
        Parameters:
        v - true anomaly (rad)
        e - eccentricity
        Returns:
        H the hyperbolic eccentric anomaly
      • hyperbolicEccentricToMean

        public static <T extends org.hipparchus.RealFieldElement<T>> T hyperbolicEccentricToMean​(T H,
                                                                                                 T e)
        Computes the mean anomaly from the hyperbolic eccentric anomaly.
        Type Parameters:
        T - type of the field elements
        Parameters:
        H - hyperbolic eccentric anomaly (rad)
        e - eccentricity
        Returns:
        M the mean anomaly
      • meanToEllipticEccentric

        public static <T extends org.hipparchus.RealFieldElement<T>> T meanToEllipticEccentric​(T M,
                                                                                               T e)
        Computes the elliptic eccentric anomaly from the mean anomaly.

        The algorithm used here for solving Kepler equation has been published in: "Procedures for solving Kepler's Equation", A. W. Odell and R. H. Gooding, Celestial Mechanics 38 (1986) 307-334

        Type Parameters:
        T - type of the field elements
        Parameters:
        M - mean anomaly (rad)
        e - eccentricity
        Returns:
        E the eccentric anomaly
      • meanToHyperbolicEccentric

        public static <T extends org.hipparchus.RealFieldElement<T>> T meanToHyperbolicEccentric​(T M,
                                                                                                 T e)
        Computes the hyperbolic eccentric anomaly from the mean anomaly.

        The algorithm used here for solving hyperbolic Kepler equation is Danby's iterative method (3rd order) with Vallado's initial guess.

        Type Parameters:
        T - Type of the field elements
        Parameters:
        M - mean anomaly (rad)
        e - eccentricity
        Returns:
        H the hyperbolic eccentric anomaly
      • ellipticEccentricToMean

        public static <T extends org.hipparchus.RealFieldElement<T>> T ellipticEccentricToMean​(T E,
                                                                                               T e)
        Computes the mean anomaly from the elliptic eccentric anomaly.
        Type Parameters:
        T - type of the field elements
        Parameters:
        E - eccentric anomaly (rad)
        e - eccentricity
        Returns:
        M the mean anomaly
      • getEquinoctialEx

        public T getEquinoctialEx()
        Get the first component of the equinoctial eccentricity vector.
        Specified by:
        getEquinoctialEx in class FieldOrbit<T extends org.hipparchus.RealFieldElement<T>>
        Returns:
        first component of the equinoctial eccentricity vector
      • getEquinoctialExDot

        public T getEquinoctialExDot()
        Get the first component of the equinoctial eccentricity vector.

        If the orbit was created without derivatives, the value returned is null.

        Specified by:
        getEquinoctialExDot in class FieldOrbit<T extends org.hipparchus.RealFieldElement<T>>
        Returns:
        first component of the equinoctial eccentricity vector
      • getEquinoctialEy

        public T getEquinoctialEy()
        Get the second component of the equinoctial eccentricity vector.
        Specified by:
        getEquinoctialEy in class FieldOrbit<T extends org.hipparchus.RealFieldElement<T>>
        Returns:
        second component of the equinoctial eccentricity vector
      • getEquinoctialEyDot

        public T getEquinoctialEyDot()
        Get the second component of the equinoctial eccentricity vector.

        If the orbit was created without derivatives, the value returned is null.

        Specified by:
        getEquinoctialEyDot in class FieldOrbit<T extends org.hipparchus.RealFieldElement<T>>
        Returns:
        second component of the equinoctial eccentricity vector
      • getHx

        public T getHx()
        Get the first component of the inclination vector.
        Specified by:
        getHx in class FieldOrbit<T extends org.hipparchus.RealFieldElement<T>>
        Returns:
        first component of the inclination vector
      • getHxDot

        public T getHxDot()
        Get the first component of the inclination vector derivative.

        If the orbit was created without derivatives, the value returned is null.

        Specified by:
        getHxDot in class FieldOrbit<T extends org.hipparchus.RealFieldElement<T>>
        Returns:
        first component of the inclination vector derivative
      • getHy

        public T getHy()
        Get the second component of the inclination vector.
        Specified by:
        getHy in class FieldOrbit<T extends org.hipparchus.RealFieldElement<T>>
        Returns:
        second component of the inclination vector
      • getHyDot

        public T getHyDot()
        Get the second component of the inclination vector derivative.
        Specified by:
        getHyDot in class FieldOrbit<T extends org.hipparchus.RealFieldElement<T>>
        Returns:
        second component of the inclination vector derivative
      • getLv

        public T getLv()
        Get the true longitude argument.
        Specified by:
        getLv in class FieldOrbit<T extends org.hipparchus.RealFieldElement<T>>
        Returns:
        v + ω + Ω true longitude argument (rad)
      • getLvDot

        public T getLvDot()
        Get the true longitude argument derivative.

        If the orbit was created without derivatives, the value returned is null.

        Specified by:
        getLvDot in class FieldOrbit<T extends org.hipparchus.RealFieldElement<T>>
        Returns:
        d(v + ω + Ω)/dt true longitude argument derivative (rad/s)
      • getLE

        public T getLE()
        Get the eccentric longitude argument.
        Specified by:
        getLE in class FieldOrbit<T extends org.hipparchus.RealFieldElement<T>>
        Returns:
        E + ω + Ω eccentric longitude argument (rad)
      • getLEDot

        public T getLEDot()
        Get the eccentric longitude argument derivative.

        If the orbit was created without derivatives, the value returned is null.

        Specified by:
        getLEDot in class FieldOrbit<T extends org.hipparchus.RealFieldElement<T>>
        Returns:
        d(E + ω + Ω)/dt eccentric longitude argument derivative (rad/s)
      • getLM

        public T getLM()
        Get the mean longitude argument.
        Specified by:
        getLM in class FieldOrbit<T extends org.hipparchus.RealFieldElement<T>>
        Returns:
        M + ω + Ω mean longitude argument (rad)
      • getLMDot

        public T getLMDot()
        Get the mean longitude argument derivative.

        If the orbit was created without derivatives, the value returned is null.

        Specified by:
        getLMDot in class FieldOrbit<T extends org.hipparchus.RealFieldElement<T>>
        Returns:
        d(M + ω + Ω)/dt mean longitude argument derivative (rad/s)
      • initPVCoordinates

        protected TimeStampedFieldPVCoordinates<T> initPVCoordinates()
        Compute the position/velocity coordinates from the canonical parameters.
        Specified by:
        initPVCoordinates in class FieldOrbit<T extends org.hipparchus.RealFieldElement<T>>
        Returns:
        computed position/velocity coordinates
      • shiftedBy

        public FieldKeplerianOrbit<T> shiftedBy​(double dt)
        Get a time-shifted instance.
        Parameters:
        dt - time shift in seconds
        Returns:
        a new instance, shifted with respect to instance (which is not changed)
      • shiftedBy

        public FieldKeplerianOrbit<T> shiftedBy​(T dt)
        Get a time-shifted orbit.

        The orbit can be slightly shifted to close dates. This shift is based on a simple Keplerian model. It is not intended as a replacement for proper orbit and attitude propagation but should be sufficient for small time shifts or coarse accuracy.

        Specified by:
        shiftedBy in interface FieldTimeShiftable<FieldOrbit<T extends org.hipparchus.RealFieldElement<T>>,​T extends org.hipparchus.RealFieldElement<T>>
        Specified by:
        shiftedBy in class FieldOrbit<T extends org.hipparchus.RealFieldElement<T>>
        Parameters:
        dt - time shift in seconds
        Returns:
        a new orbit, shifted with respect to the instance (which is immutable)
      • interpolate

        public FieldKeplerianOrbit<T> interpolate​(FieldAbsoluteDate<T> date,
                                                  Stream<FieldOrbit<T>> sample)
        Get an interpolated instance.

        Note that the state of the current instance may not be used in the interpolation process, only its type and non interpolable fields are used (for example central attraction coefficient or frame when interpolating orbits). The interpolable fields taken into account are taken only from the states of the sample points. So if the state of the instance must be used, the instance should be included in the sample points.

        Note that this method is designed for small samples only (say up to about 10-20 points) so it can be implemented using polynomial interpolation (typically Hermite interpolation). Using too much points may induce Runge's phenomenon and numerical problems (including NaN appearing).

        The interpolated instance is created by polynomial Hermite interpolation on Keplerian elements, without derivatives (which means the interpolation falls back to Lagrange interpolation only).

        As this implementation of interpolation is polynomial, it should be used only with small samples (about 10-20 points) in order to avoid Runge's phenomenon and numerical problems (including NaN appearing).

        If orbit interpolation on large samples is needed, using the Ephemeris class is a better way than using this low-level interpolation. The Ephemeris class automatically handles selection of a neighboring sub-sample with a predefined number of point from a large global sample in a thread-safe way.

        Parameters:
        date - interpolation date
        sample - sample points on which interpolation should be done
        Returns:
        a new instance, interpolated at specified date
      • computeJacobianMeanWrtCartesian

        protected T[][] computeJacobianMeanWrtCartesian()
        Compute the Jacobian of the orbital parameters with mean angle with respect to the Cartesian parameters.

        Element jacobian[i][j] is the derivative of parameter i of the orbit with respect to Cartesian coordinate j. This means each row correspond to one orbital parameter whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.

        Specified by:
        computeJacobianMeanWrtCartesian in class FieldOrbit<T extends org.hipparchus.RealFieldElement<T>>
        Returns:
        6x6 Jacobian matrix
        See Also:
        FieldOrbit.computeJacobianEccentricWrtCartesian(), FieldOrbit.computeJacobianTrueWrtCartesian()
      • computeJacobianEccentricWrtCartesian

        protected T[][] computeJacobianEccentricWrtCartesian()
        Compute the Jacobian of the orbital parameters with eccentric angle with respect to the Cartesian parameters.

        Element jacobian[i][j] is the derivative of parameter i of the orbit with respect to Cartesian coordinate j. This means each row correspond to one orbital parameter whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.

        Specified by:
        computeJacobianEccentricWrtCartesian in class FieldOrbit<T extends org.hipparchus.RealFieldElement<T>>
        Returns:
        6x6 Jacobian matrix
        See Also:
        FieldOrbit.computeJacobianMeanWrtCartesian(), FieldOrbit.computeJacobianTrueWrtCartesian()
      • computeJacobianTrueWrtCartesian

        protected T[][] computeJacobianTrueWrtCartesian()
        Compute the Jacobian of the orbital parameters with true angle with respect to the Cartesian parameters.

        Element jacobian[i][j] is the derivative of parameter i of the orbit with respect to Cartesian coordinate j. This means each row correspond to one orbital parameter whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.

        Specified by:
        computeJacobianTrueWrtCartesian in class FieldOrbit<T extends org.hipparchus.RealFieldElement<T>>
        Returns:
        6x6 Jacobian matrix
        See Also:
        FieldOrbit.computeJacobianMeanWrtCartesian(), FieldOrbit.computeJacobianEccentricWrtCartesian()
      • addKeplerContribution

        public void addKeplerContribution​(PositionAngle type,
                                          double gm,
                                          T[] pDot)
        Add the contribution of the Keplerian motion to parameters derivatives

        This method is used by integration-based propagators to evaluate the part of Keplerian motion to evolution of the orbital state.

        Specified by:
        addKeplerContribution in class FieldOrbit<T extends org.hipparchus.RealFieldElement<T>>
        Parameters:
        type - type of the position angle in the state
        gm - attraction coefficient to use
        pDot - array containing orbital state derivatives to update (the Keplerian part must be added to the array components, as the array may already contain some non-zero elements corresponding to non-Keplerian parts)
      • toString

        public String toString()
        Returns a string representation of this Keplerian parameters object.
        Overrides:
        toString in class Object
        Returns:
        a string representation of this object
      • normalizeAngle

        public static <T extends org.hipparchus.RealFieldElement<T>> T normalizeAngle​(T a,
                                                                                      T center)
        Normalize an angle in a 2π wide interval around a center value.

        This method has three main uses:

        • normalize an angle between 0 and 2π:
          a = MathUtils.normalizeAngle(a, FastMath.PI);
        • normalize an angle between -π and +π
          a = MathUtils.normalizeAngle(a, 0.0);
        • compute the angle between two defining angular positions:
          angle = MathUtils.normalizeAngle(end, start) - start;

        Note that due to numerical accuracy and since π cannot be represented exactly, the result interval is closed, it cannot be half-closed as would be more satisfactory in a purely mathematical view.

        Type Parameters:
        T - the type of the field elements
        Parameters:
        a - angle to normalize
        center - center of the desired 2π interval for the result
        Returns:
        a-2kπ with integer k and center-π <= a-2kπ <= center+π
      • toOrbit

        public KeplerianOrbit toOrbit()
        Description copied from class: FieldOrbit
        Transforms the FieldOrbit instance into an Orbit instance.
        Specified by:
        toOrbit in class FieldOrbit<T extends org.hipparchus.RealFieldElement<T>>
        Returns:
        Orbit instance with same properties