GC2GD Subroutine

public subroutine GC2GD(n, xyz, elong, phi, height, j)

Transform geocentric coordinates to geodetic using the specified reference ellipsoid.

Status: canonical transformation.

Notes

  1. The identifier N is a number that specifies the choice of reference ellipsoid. The following are supported:

    N   ellipsoid
    
    1    WGS84
    2    GRS80
    3    WGS72
    

    The number N has no significance outside the SOFA software.

  2. The geocentric vector (XYZ, given) and height (HEIGHT, returned) are in meters.

  3. An error status J=-1 means that the identifier N is illegal. An error status J=-2 is theoretically impossible. In all error cases, all three results are set to -1D9.

  4. The inverse transformation is performed in the routine GD2GC.

History

  • IAU SOFA revision: 2013 September 1

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: n

ellipsoid identifier (Note 1)

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

geocentric vector (Note 2)

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

longitude (radians, east +ve, Note 3)

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

latitude (geodetic, radians, Note 3)

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

height above ellipsoid (geodetic, Notes 2,3)

integer, intent(out) :: j
  • 0 = OK
  • -1 = illegal identifier (Note 3)
  • -2 = internal error (Note 3)

Calls

proc~~gc2gd~~CallsGraph proc~gc2gd GC2GD proc~eform EFORM proc~gc2gd->proc~eform proc~gc2gde GC2GDE proc~gc2gd->proc~gc2gde

Contents

Source Code


Source Code

    subroutine GC2GD ( n, xyz, elong, phi, height, j )

    implicit none

    integer,intent(in) :: n !! ellipsoid identifier (Note 1)
    real(wp),dimension(3),intent(in) :: xyz !! geocentric vector (Note 2)
    real(wp),intent(out) :: elong !! longitude (radians, east +ve, Note 3)
    real(wp),intent(out) :: phi !! latitude (geodetic, radians, Note 3)
    real(wp),intent(out) :: height !! height above ellipsoid (geodetic, Notes 2,3)
    integer,intent(out) :: j !! status:
                             !! * 0 = OK
                             !! * -1 = illegal identifier (Note 3)
                             !! * -2 = internal error (Note 3)

    real(wp) :: a, f

    !  Obtain reference ellipsoid parameters.
    call EFORM ( n, a, f, j )

    !  If OK, transform x,y,z to longitude, geodetic latitude, height.
    if ( j==0 ) then
       call GC2GDE ( a, f, xyz, elong, phi, height, j )
       if ( j<0 ) j=-2
    end if

    !  Deal with any errors.
    if ( j<0 ) then
       elong = -1.0e9_wp
       phi = -1.0e9_wp
       height = -1.0e9_wp
    end if

    end subroutine GC2GD