Class OneAxisEllipsoid

  • All Implemented Interfaces:
    Serializable, BodyShape
    Direct Known Subclasses:
    ReferenceEllipsoid

    public class OneAxisEllipsoid
    extends Ellipsoid
    implements BodyShape
    Modeling of a one-axis ellipsoid.

    One-axis ellipsoids is a good approximate model for most planet-size and larger natural bodies. It is the equilibrium shape reached by a fluid body under its own gravity field when it rotates. The symmetry axis is the rotation or polar axis.

    Author:
    Luc Maisonobe, Guylaine Prat
    See Also:
    Serialized Form
    • Method Detail

      • setAngularThreshold

        public void setAngularThreshold​(double angularThreshold)
        Set the angular convergence threshold.

        The angular threshold is used both to identify points close to the ellipse axes and as the convergence threshold used to stop the iterations in the transform(Vector3D, Frame, AbsoluteDate) method.

        If this method is not called, the default value is set to 10-12.

        Parameters:
        angularThreshold - angular convergence threshold (rad)
      • getEquatorialRadius

        public double getEquatorialRadius()
        Get the equatorial radius of the body.
        Returns:
        equatorial radius of the body (m)
      • getFlattening

        public double getFlattening()
        Get the flattening of the body: f = (a-b)/a.
        Returns:
        the flattening
      • getEccentricitySquared

        public double getEccentricitySquared()
        Get the first eccentricity squared of the ellipsoid: e^2 = f * (2.0 - f).
        Returns:
        the eccentricity squared
      • getEccentricity

        public double getEccentricity()
        Get the first eccentricity of the ellipsoid: e = sqrt(f * (2.0 - f)).
        Returns:
        the eccentricity
      • getBodyFrame

        public Frame getBodyFrame()
        Get body frame related to body shape.
        Specified by:
        getBodyFrame in interface BodyShape
        Returns:
        body frame related to body shape
      • getCartesianIntersectionPoint

        public Vector3D getCartesianIntersectionPoint​(Line line,
                                                      Vector3D close,
                                                      Frame frame,
                                                      AbsoluteDate date)
        Get the intersection point of a line with the surface of the body.

        A line may have several intersection points with a closed surface (we consider the one point case as a degenerated two points case). The close parameter is used to select which of these points should be returned. The selected point is the one that is closest to the close point.

        Parameters:
        line - test line (may intersect the body or not)
        close - point used for intersections selection
        frame - frame in which line is expressed
        date - date of the line in given frame
        Returns:
        intersection point at altitude zero or null if the line does not intersect the surface
        Since:
        9.3
      • getIntersectionPoint

        public GeodeticPoint getIntersectionPoint​(Line line,
                                                  Vector3D close,
                                                  Frame frame,
                                                  AbsoluteDate date)
        Get the intersection point of a line with the surface of the body.

        A line may have several intersection points with a closed surface (we consider the one point case as a degenerated two points case). The close parameter is used to select which of these points should be returned. The selected point is the one that is closest to the close point.

        Specified by:
        getIntersectionPoint in interface BodyShape
        Parameters:
        line - test line (may intersect the body or not)
        close - point used for intersections selection
        frame - frame in which line is expressed
        date - date of the line in given frame
        Returns:
        intersection point at altitude zero or null if the line does not intersect the surface
      • getCartesianIntersectionPoint

        public <T extends CalculusFieldElement<T>> FieldVector3D<T> getCartesianIntersectionPoint​(FieldLine<T> line,
                                                                                                  FieldVector3D<T> close,
                                                                                                  Frame frame,
                                                                                                  FieldAbsoluteDate<T> date)
        Get the intersection point of a line with the surface of the body.

        A line may have several intersection points with a closed surface (we consider the one point case as a degenerated two points case). The close parameter is used to select which of these points should be returned. The selected point is the one that is closest to the close point.

        Type Parameters:
        T - type of the field elements
        Parameters:
        line - test line (may intersect the body or not)
        close - point used for intersections selection
        frame - frame in which line is expressed
        date - date of the line in given frame
        Returns:
        intersection point at altitude zero or null if the line does not intersect the surface
        Since:
        9.3
      • getIntersectionPoint

        public <T extends CalculusFieldElement<T>> FieldGeodeticPoint<T> getIntersectionPoint​(FieldLine<T> line,
                                                                                              FieldVector3D<T> close,
                                                                                              Frame frame,
                                                                                              FieldAbsoluteDate<T> date)
        Get the intersection point of a line with the surface of the body.

        A line may have several intersection points with a closed surface (we consider the one point case as a degenerated two points case). The close parameter is used to select which of these points should be returned. The selected point is the one that is closest to the close point.

        Specified by:
        getIntersectionPoint in interface BodyShape
        Type Parameters:
        T - type of the field elements
        Parameters:
        line - test line (may intersect the body or not)
        close - point used for intersections selection
        frame - frame in which line is expressed
        date - date of the line in given frame
        Returns:
        intersection point at altitude zero or null if the line does not intersect the surface
      • transform

        public Vector3D transform​(GeodeticPoint point)
        Transform a surface-relative point to a Cartesian point.
        Specified by:
        transform in interface BodyShape
        Parameters:
        point - surface-relative point
        Returns:
        point at the same location but as a Cartesian point
      • transform

        public <T extends CalculusFieldElement<T>> FieldVector3D<T> transform​(FieldGeodeticPoint<T> point)
        Transform a surface-relative point to a Cartesian point.
        Specified by:
        transform in interface BodyShape
        Type Parameters:
        T - type of the filed elements
        Parameters:
        point - surface-relative point
        Returns:
        point at the same location but as a Cartesian point
      • transform

        public GeodeticPoint transform​(Vector3D point,
                                       Frame frame,
                                       AbsoluteDate date)
        Transform a Cartesian point to a surface-relative point.

        This method is based on Toshio Fukushima's algorithm which uses Halley's method. transformation from Cartesian to Geodetic Coordinates Accelerated by Halley's Method, Toshio Fukushima, Journal of Geodesy 9(12):689-693, February 2006

        Some changes have been added to the original method:

        • in order to handle more accurately corner cases near the pole
        • in order to handle properly corner cases near the equatorial plane, even far inside the ellipsoid
        • in order to handle very flat ellipsoids

        In some rare cases (for example very flat ellipsoid, or points close to ellipsoid center), the loop may fail to converge. As this seems to happen only in degenerate cases, a design choice was to return an approximate point corresponding to last iteration. This point may be incorrect and fail to give the initial point back if doing roundtrip by calling transform(GeodeticPoint). This design choice was made to avoid NaNs appearing for example in inter-satellites visibility checks when two satellites are almost on opposite sides of Earth. The intermediate points far within the Earth should not prevent the detection algorithm to find visibility start/end.

        Specified by:
        transform in interface BodyShape
        Parameters:
        point - Cartesian point
        frame - frame in which Cartesian point is expressed
        date - date of the computation (used for frames conversions)
        Returns:
        point at the same location but as a surface-relative point
      • transform

        public <T extends CalculusFieldElement<T>> FieldGeodeticPoint<T> transform​(FieldVector3D<T> point,
                                                                                   Frame frame,
                                                                                   FieldAbsoluteDate<T> date)
        Transform a Cartesian point to a surface-relative point.

        This method is based on Toshio Fukushima's algorithm which uses Halley's method. transformation from Cartesian to Geodetic Coordinates Accelerated by Halley's Method, Toshio Fukushima, Journal of Geodesy 9(12):689-693, February 2006

        Some changes have been added to the original method:

        • in order to handle more accurately corner cases near the pole
        • in order to handle properly corner cases near the equatorial plane, even far inside the ellipsoid
        • in order to handle very flat ellipsoids

        In some rare cases (for example very flat ellipsoid, or points close to ellipsoid center), the loop may fail to converge. As this seems to happen only in degenerate cases, a design choice was to return an approximate point corresponding to last iteration. This point may be incorrect and fail to give the initial point back if doing roundtrip by calling transform(GeodeticPoint). This design choice was made to avoid NaNs appearing for example in inter-satellites visibility checks when two satellites are almost on opposite sides of Earth. The intermediate points far within the Earth should not prevent the detection algorithm to find visibility start/end.

        Specified by:
        transform in interface BodyShape
        Type Parameters:
        T - type of the filed elements
        Parameters:
        point - Cartesian point
        frame - frame in which Cartesian point is expressed
        date - date of the computation (used for frames conversions)
        Returns:
        point at the same location but as a surface-relative point
      • transform

        public FieldGeodeticPoint<DerivativeStructure> transform​(PVCoordinates point,
                                                                 Frame frame,
                                                                 AbsoluteDate date)
        Transform a Cartesian point to a surface-relative point.
        Parameters:
        point - Cartesian point
        frame - frame in which Cartesian point is expressed
        date - date of the computation (used for frames conversions)
        Returns:
        point at the same location but as a surface-relative point, using time as the single derivation parameter
      • azimuthBetweenPoints

        public double azimuthBetweenPoints​(GeodeticPoint origin,
                                           GeodeticPoint destination)
        Compute the azimuth angle from local north between the two points. The angle is calculated clockwise from local north at the origin point and follows the rhumb line to the destination point.
        Parameters:
        origin - the origin point, at which the azimuth angle will be computed (non-null)
        destination - the destination point, to which the angle is defined (non-null)
        Returns:
        the resulting azimuth angle (radians, [0-2pi))
        Since:
        11.3
      • azimuthBetweenPoints

        public <T extends CalculusFieldElement<T>> T azimuthBetweenPoints​(FieldGeodeticPoint<T> origin,
                                                                          FieldGeodeticPoint<T> destination)
        Compute the azimuth angle from local north between the two points. The angle is calculated clockwise from local north at the origin point and follows the rhumb line to the destination point.
        Type Parameters:
        T - the type of field elements
        Parameters:
        origin - the origin point, at which the azimuth angle will be computed (non-null)
        destination - the destination point, to which the angle is defined (non-null)
        Returns:
        the resulting azimuth angle (radians, [0-2pi))
        Since:
        11.3
      • geodeticToIsometricLatitude

        public double geodeticToIsometricLatitude​(double geodeticLatitude)
        Compute the isometric latitude corresponding to the provided latitude.
        Parameters:
        geodeticLatitude - the latitude (radians, within interval [-pi/2, +pi/2])
        Returns:
        the isometric latitude (radians)
        Since:
        11.3
      • geodeticToIsometricLatitude

        public <T extends CalculusFieldElement<T>> T geodeticToIsometricLatitude​(T geodeticLatitude)
        Compute the isometric latitude corresponding to the provided latitude.
        Type Parameters:
        T - the type of field elements
        Parameters:
        geodeticLatitude - the latitude (radians, within interval [-pi/2, +pi/2])
        Returns:
        the isometric latitude (radians)
        Since:
        11.3
      • lowestAltitudeIntermediate

        public GeodeticPoint lowestAltitudeIntermediate​(Vector3D endpoint1,
                                                        Vector3D endpoint2)
        Find intermediate point of lowest altitude along a line between two endpoints.
        Parameters:
        endpoint1 - first endpoint, in body frame
        endpoint2 - second endpoint, in body frame
        Returns:
        point with lowest altitude between endpoint1 and endpoint2.
        Since:
        12.0
      • lowestAltitudeIntermediate

        public <T extends CalculusFieldElement<T>> FieldGeodeticPoint<T> lowestAltitudeIntermediate​(FieldVector3D<T> endpoint1,
                                                                                                    FieldVector3D<T> endpoint2)
        Find intermediate point of lowest altitude along a line between two endpoints.
        Type Parameters:
        T - type of the field elements
        Parameters:
        endpoint1 - first endpoint, in body frame
        endpoint2 - second endpoint, in body frame
        Returns:
        point with lowest altitude between endpoint1 and endpoint2.
        Since:
        12.0