Class Geoid

  • All Implemented Interfaces:
    Serializable, BodyShape, EarthShape

    public class Geoid
    extends Object
    implements EarthShape
    A geoid is a level surface of the gravity potential of a body. The gravity potential, W, is split so W = U + T, where U is the normal potential (defined by the ellipsoid) and T is the anomalous potential.[3](eq. 2-137)

    The getIntersectionPoint(Line, Vector3D, Frame, AbsoluteDate) method is tailored specifically for Earth's geoid. All of the other methods in this class are general and will work for an arbitrary body.

    There are several components that are needed to define a geoid[1]:

    • Geopotential field. These are the coefficients of the spherical harmonics: Sn,m and Cn,m
    • Reference Ellipsoid. The ellipsoid is used to define the undulation of the geoid (distance between ellipsoid and geoid) and U0 the value of the normal gravity potential at the surface of the ellipsoid.
    • W0, the potential at the geoid. The value of the potential on the level surface. This is taken to be U0, the normal gravity potential at the surface of the ReferenceEllipsoid.
    • Permanent Tide System. This implementation assumes that the geopotential field and the reference ellipsoid use the same permanent tide system. If the assumption is false it will produce errors of about 0.5 m. Conversion between tide systems is a possible improvement.[1,2]
    • Topographic Masses. That is mass outside of the geoid, e.g. mountains. This implementation ignores topographic masses, which causes up to 3m error in the Himalayas, and ~ 1.5m error in the Rockies. This could be improved through the use of DTED and calculating height anomalies or using the correction coefficients.[1]

    This implementation also assumes that the normal to the reference ellipsoid is the same as the normal to the geoid. This assumption enables the equation: (height above geoid) = (height above ellipsoid) - (undulation), which is used in transform(GeodeticPoint) and transform(Vector3D, Frame, AbsoluteDate).

    In testing, the error in the undulations calculated by this class were off by less than 3 meters, which matches the assumptions outlined above.

    References:

    1. Dru A. Smith. There is no such thing as "The" EGM96 geoid: Subtle points on the use of a global geopotential model. IGeS Bulletin No. 8:17-28, 1998. http ://www.ngs.noaa.gov/PUBS_LIB/EGM96_GEOID_PAPER/egm96_geoid_paper.html
    2. Martin Losch, Verena Seufer. How to Compute Geoid Undulations (Geoid Height Relative to a Given Reference Ellipsoid) from Spherical Harmonic Coefficients for Satellite Altimetry Applications. , 2003. mitgcm.org/~mlosch/geoidcookbook.pdf
    3. Weikko A. Heiskanen, Helmut Moritz. Physical Geodesy. W. H. Freeman and Company, 1967. (especially sections 2.13 and equation 2-144 Bruns Formula)
    4. S. A. Holmes, W. E. Featherstone. A unified approach to the Clenshaw summation and the recursive computation of very high degree and order normalised associated Legendre functions. Journal of Geodesy, 76(5):279, 2002.
    5. DMA TR 8350.2. 1984.
    6. Department of Defense World Geodetic System 1984. 2000. NIMA TR 8350.2 Third Edition, Amendment 1.
    Author:
    Evan Ward
    See Also:
    Serialized Form
    • Constructor Detail

      • Geoid

        public Geoid​(NormalizedSphericalHarmonicsProvider geopotential,
                     ReferenceEllipsoid referenceEllipsoid)
        Creates a geoid from the given geopotential, reference ellipsoid and the assumptions in the comment for Geoid.
        Parameters:
        geopotential - the gravity potential. Only the anomalous potential will be used. It is assumed that the geopotential and the referenceEllipsoid are defined in the same frame. Usually a constant geopotential is used to define a time-invariant Geoid.
        referenceEllipsoid - the normal gravity potential.
        Throws:
        NullPointerException - if geopotential == null || referenceEllipsoid == null
    • Method Detail

      • getBodyFrame

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

        public double getUndulation​(double geodeticLatitude,
                                    double longitude,
                                    AbsoluteDate date)
        Gets the Undulation of the Geoid, N at the given position. N is the distance between the reference ellipsoid and the geoid. The latitude and longitude parameters are both defined with respect to the reference ellipsoid. For EGM96 and the WGS84 ellipsoid the undulation is between -107m and +86m.

        NOTE: Restrictions are not put on the range of the arguments geodeticLatitude and longitude.

        Parameters:
        geodeticLatitude - geodetic latitude (angle between the local normal and the equatorial plane on the reference ellipsoid), in radians.
        longitude - on the reference ellipsoid, in radians.
        date - of evaluation. Used for time varying geopotential fields.
        Returns:
        the undulation in m, positive means the geoid is higher than the ellipsoid.
        See Also:
        Geoid, Geoid on Wikipedia
      • getEllipsoid

        public ReferenceEllipsoid getEllipsoid()
        Description copied from interface: EarthShape
        Get the underlying ellipsoid model that defines latitude and longitude. If the height component of a GeodeticPoint is not needed, then using the ellipsoid will provide the quickest transformation.
        Specified by:
        getEllipsoid in interface EarthShape
        Returns:
        the reference ellipsoid. May be this, but never null.
      • getIntersectionPoint

        public GeodeticPoint getIntersectionPoint​(org.hipparchus.geometry.euclidean.threed.Line lineInFrame,
                                                  org.hipparchus.geometry.euclidean.threed.Vector3D closeInFrame,
                                                  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.

        The intersection point is computed using a line search along the specified line. This is accurate when the geoid is slowly varying.

        Specified by:
        getIntersectionPoint in interface BodyShape
        Parameters:
        lineInFrame - test line (may intersect the body or not)
        closeInFrame - 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
      • projectToGround

        public org.hipparchus.geometry.euclidean.threed.Vector3D projectToGround​(org.hipparchus.geometry.euclidean.threed.Vector3D point,
                                                                                 AbsoluteDate date,
                                                                                 Frame frame)
        Description copied from interface: BodyShape
        Project a point to the ground.
        Specified by:
        projectToGround in interface BodyShape
        Parameters:
        point - point to project
        date - current date
        frame - frame in which moving point is expressed
        Returns:
        ground point exactly at the local vertical of specified point, in the same frame as specified point
        See Also:
        BodyShape.projectToGround(TimeStampedPVCoordinates, Frame)
      • getIntersectionPoint

        public <T extends org.hipparchus.RealFieldElement<T>> FieldGeodeticPoint<T> getIntersectionPoint​(org.hipparchus.geometry.euclidean.threed.FieldLine<T> lineInFrame,
                                                                                                         org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> closeInFrame,
                                                                                                         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.

        The intersection point is computed using a line search along the specified line. This is accurate when the geoid is slowly varying.

        Specified by:
        getIntersectionPoint in interface BodyShape
        Type Parameters:
        T - type of the field elements
        Parameters:
        lineInFrame - test line (may intersect the body or not)
        closeInFrame - 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 GeodeticPoint transform​(org.hipparchus.geometry.euclidean.threed.Vector3D point,
                                       Frame frame,
                                       AbsoluteDate date)
        Transform a Cartesian point to a surface-relative point.
        Specified by:
        transform in interface BodyShape
        Parameters:
        date - date of the conversion. Used for computing frame transformations and for time dependent geopotential.
        point - Cartesian point
        frame - frame in which Cartesian point is expressed
        Returns:
        The surface relative point at the same location. Altitude is orthometric height, that is height above the Geoid. Latitude and longitude are both geodetic and defined with respect to the reference ellipsoid.
        See Also:
        transform(GeodeticPoint), Orthometric_height
      • transform

        public <T extends org.hipparchus.RealFieldElement<T>> FieldGeodeticPoint<T> transform​(org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> point,
                                                                                              Frame frame,
                                                                                              FieldAbsoluteDate<T> date)
        Transform a Cartesian point to a surface-relative point.
        Specified by:
        transform in interface BodyShape
        Type Parameters:
        T - type fo the filed elements
        Parameters:
        date - date of the conversion. Used for computing frame transformations and for time dependent geopotential.
        point - Cartesian point
        frame - frame in which Cartesian point is expressed
        Returns:
        The surface relative point at the same location. Altitude is orthometric height, that is height above the Geoid. Latitude and longitude are both geodetic and defined with respect to the reference ellipsoid.
        See Also:
        transform(GeodeticPoint), Orthometric_height
      • transform

        public org.hipparchus.geometry.euclidean.threed.Vector3D transform​(GeodeticPoint point)
        Transform a surface-relative point to a Cartesian point.
        Specified by:
        transform in interface BodyShape
        Parameters:
        point - The surface relative point to transform. Altitude is orthometric height, that is height above the Geoid. Latitude and longitude are both geodetic and defined with respect to the reference ellipsoid.
        Returns:
        point at the same location but as a Cartesian point in the body frame.
        See Also:
        transform(Vector3D, Frame, AbsoluteDate)
      • transform

        public <T extends org.hipparchus.RealFieldElement<T>> org.hipparchus.geometry.euclidean.threed.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 field elements
        Parameters:
        point - The surface relative point to transform. Altitude is orthometric height, that is height above the Geoid. Latitude and longitude are both geodetic and defined with respect to the reference ellipsoid.
        Returns:
        point at the same location but as a Cartesian point in the body frame.
        Since:
        9.0
        See Also:
        transform(Vector3D, Frame, AbsoluteDate)