GD2GC Subroutine

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

Transform geodetic coordinates to geocentric 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 height (HEIGHT, given) and the geocentric vector (XYZ, returned) are in meters.

  3. No validation is performed on the arguments ELONG, PHI and HEIGHT. An error status J=-1 means that the identifier N is illegal. An error status J=-2 protects against cases that would lead to arithmetic exceptions. In all error cases, XYZ is set to zeros.

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

History

  • IAU SOFA revision: 2010 January 18

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: n

ellipsoid identifier (Note 1)

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

longitude (radians, east +ve)

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

latitude (geodetic, radians, Note 3)

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

height above ellipsoid (geodetic, Notes 2,3)

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

geocentric vector (Note 2)

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

Calls

proc~~gd2gc~~CallsGraph proc~gd2gc GD2GC proc~eform EFORM proc~gd2gc->proc~eform proc~gd2gce GD2GCE proc~gd2gc->proc~gd2gce proc~zp ZP proc~gd2gc->proc~zp

Called by

proc~~gd2gc~~CalledByGraph proc~gd2gc GD2GC 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 GD2GC ( n, elong, phi, height, xyz, j )

    implicit none

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

    real(wp) :: a, f

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

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

    !  Deal with any errors.
    if ( j/=0 ) call ZP ( xyz )

    end subroutine GD2GC