ESMF_PointListGetForLoc Subroutine

public subroutine ESMF_PointListGetForLoc(pointlist, loc, id, loc_coords, rc)

Arguments

Type IntentOptional Attributes Name
type(ESMF_PointList), intent(in) :: pointlist
integer, intent(in) :: loc
integer, intent(out), optional :: id
real(kind=ESMF_KIND_R8), intent(out), optional, dimension(:) :: loc_coords
integer, intent(out), optional :: rc

Source Code

  subroutine ESMF_PointListGetForLoc(pointlist, loc, id, loc_coords, rc)
!
! !ARGUMENTS:
    type(ESMF_PointList), intent(in)                          :: pointlist
    integer,            intent(in)                          :: loc
    integer,            intent(out), optional               :: id
    real(ESMF_KIND_R8), intent(out), optional, dimension(:) :: loc_coords
    integer,            intent(out), optional               :: rc
!
! !DESCRIPTION:
!     Returns information about an {\tt ESMF\_PointList}.
!
!     The arguments are:
!     \begin{description}
!     \item[pointlist]
!          {\tt ESMF\_PointList} to be queried.
!     \item[loc]
!          Location within Pointlist to be queried. Locations values begin with zero.
!     \item[{[id]}]
!          Returns the id associated with location.
!     \item[{[loc_coords]}]
!          Returns array of coordinates associated with location.
!     \item[{[rc]}]
!          Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!     \end{description}
!
!EOPI
!------------------------------------------------------------------------------
    integer                 :: localrc      ! local return code

    ! initialize return code; assume routine not implemented
    localrc = ESMF_RC_NOT_IMPL
    if (present(rc)) rc = ESMF_RC_NOT_IMPL

    ESMF_INIT_CHECK_DEEP(ESMF_PointListGetInit,pointlist,rc)

    if (present(id)) then
      call c_ESMC_PointListGetId(pointlist, loc, id, localrc)
      if (ESMF_LogFoundError(localrc, &
        ESMF_ERR_PASSTHRU, &
        ESMF_CONTEXT, rcToReturn=rc)) return
    endif

    if (present(loc_coords)) then
      call c_ESMC_PointListGetCoords(pointlist, loc, loc_coords, localrc)

      if (ESMF_LogFoundError(localrc, &
        ESMF_ERR_PASSTHRU, &
        ESMF_CONTEXT, rcToReturn=rc)) return
    endif

    ! Return successfully
    if (present(rc)) rc = ESMF_SUCCESS

  end subroutine ESMF_PointListGetForLoc