GD2GCE Subroutine

public subroutine GD2GCE(a, f, elong, phi, height, xyz, j)

Transform geodetic coordinates to geocentric for a reference ellipsoid of specified form.

Status: support routine.

Notes

  1. The equatorial radius, A, can be in any units, but meters is the conventional choice.

  2. The flattening, F, is (for the Earth) a value around 0.00335, i.e. around 1/298.

  3. The equatorial radius, A, and the height, HEIGHT, must be given in the same units, and determine the units of the returned geocentric vector, XYZ.

  4. No validation is performed on individual arguments. The error status J=-1 protects against (unrealistic) cases that would lead to arithmetic exceptions. If an error occurs, XYZ is unchanged.

  5. The inverse transformation is performed in the routine GC2GDE.

  6. The transformation for a standard ellipsoid (such as WGS84) can more conveniently be performed by calling GD2GC, which uses a numerical code (1 for WGS84) to identify the required A and F values.

References

  • Green, R.M., Spherical Astronomy, Cambridge University Press, (1985) Section 4.5, p96.

  • Explanatory Supplement to the Astronomical Almanac, P. Kenneth Seidelmann (ed), University Science Books (1992), Section 4.22, p202.

History

  • IAU SOFA revision: 2009 November 2

Arguments

TypeIntentOptionalAttributesName
real(kind=wp), intent(in) :: a

equatorial radius (Notes 1,4)

real(kind=wp), intent(in) :: f

flattening (Notes 2,4)

real(kind=wp), intent(in) :: elong

longitude (radians, east +ve)

real(kind=wp), intent(in) :: phi

latitude (geodetic, radians, Note 4)

real(kind=wp), intent(in) :: height

height above ellipsoid (geodetic, Notes 3,4)

real(kind=wp), intent(out), dimension(3):: xyz

geocentric vector (Note 3)

integer, intent(out) :: j
  • 0 = OK
  • -1 = illegal case (Note 4)

Called by

proc~~gd2gce~~CalledByGraph proc~gd2gce GD2GCE proc~gd2gc GD2GC proc~gd2gc->proc~gd2gce proc~pvtob PVTOB proc~pvtob->proc~gd2gc proc~apco APCO proc~apco->proc~pvtob proc~apio APIO proc~apio->proc~pvtob proc~apco13 APCO13 proc~apco13->proc~apco proc~apio13 APIO13 proc~apio13->proc~apio proc~atco13 ATCO13 proc~atco13->proc~apco13 proc~atoi13 ATOI13 proc~atoi13->proc~apio13 proc~atio13 ATIO13 proc~atio13->proc~apio13 proc~atoc13 ATOC13 proc~atoc13->proc~apco13

Contents

Source Code


Source Code

    subroutine GD2GCE ( a, f, elong, phi, height, xyz, j )

    implicit none

    real(wp),intent(in) :: a !! equatorial radius (Notes 1,4)
    real(wp),intent(in) :: f !! flattening (Notes 2,4)
    real(wp),intent(in) :: elong !! longitude (radians, east +ve)
    real(wp),intent(in) :: phi !! latitude (geodetic, radians, Note 4)
    real(wp),intent(in) :: height !! height above ellipsoid (geodetic, Notes 3,4)
    real(wp),dimension(3),intent(out) :: xyz !! geocentric vector (Note 3)
    integer,intent(out) :: j !! status:
                             !! * 0 = OK
                             !! * -1 = illegal case (Note 4)

    real(wp) :: sp, cp, w, d, ac, as, r

    !  Functions of geodetic latitude.
    sp = sin(phi)
    cp = cos(phi)
    w = 1.0_wp-f
    w = w*w
    d = cp*cp + w*sp*sp
    if ( d > 0.0_wp ) then
       ac = a / sqrt(d)
       as = w * ac

       !  Geocentric vector.
       r = ( ac + height ) * cp
       xyz(1) = r * cos(elong)
       xyz(2) = r * sin(elong)
       xyz(3) = ( as + height ) * sp

       !  Success.
       j = 0
    else

       !  Fail.
       j = -1
    end if

    end subroutine GD2GCE