ESMF_PointListAdd Subroutine

public subroutine ESMF_PointListAdd(pointlist, id, loc_coords, loc_orig_coords, rc)

Arguments

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

Source Code

  subroutine ESMF_PointListAdd(pointlist, id, loc_coords, loc_orig_coords, rc)
!
! !ARGUMENTS:
    type(ESMF_PointList), intent(in)                      :: pointlist
    integer,            intent(in)                        :: id
    real(ESMF_KIND_R8), intent(in), dimension(:)          :: loc_coords
    real(ESMF_KIND_R8), intent(in), dimension(:),optional :: loc_orig_coords
    integer,            intent(out), optional             :: rc
!
! !DESCRIPTION:
!   Add a point to an {\tt ESMF\_PointList} with the given values.
!
!     The arguments are:
!     \begin{description}
!     \item[pointlist]
!          {\tt ESMF\_PointList} to be queried.
!     \item[{[id]}]
!          The id associated with point to add.
!     \item[{[loc_coords]}]
!          The array of coordinates associated with point to add.
!     \item[{[loc_orig_coords]}]
!          The array of orig. coordinates associated with point to add.
!          (Only valid if the pointlist has been created to hold original
!           coords (e.g. by creating with origCoordDim != 0.0))
!     \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)

    ! Add based on presence of orig coords.
    if (present(loc_orig_coords)) then
       call c_ESMC_PointListAddWOrig(pointlist, id, loc_coords, loc_orig_coords, localrc)
       if (ESMF_LogFoundError(localrc, &
            ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return
    else
       call c_ESMC_PointListAdd(pointlist, id, 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_PointListAdd