ESMF_PointListCreateFrmMesh Function

private function ESMF_PointListCreateFrmMesh(mesh, meshLoc, maskValues, addOrigCoords, rc)

Arguments

Type IntentOptional Attributes Name
type(ESMF_Mesh), intent(in) :: mesh
type(ESMF_MeshLoc), intent(in) :: meshLoc
integer(kind=ESMF_KIND_I4), optional :: maskValues(:)
logical, optional :: addOrigCoords
integer, intent(out), optional :: rc

Return Value type(ESMF_PointList)


Source Code

  function ESMF_PointListCreateFrmMesh(mesh, meshLoc, maskValues, addOrigCoords, rc)
!
! !RETURN VALUE:
    type(ESMF_PointList) :: ESMF_PointListCreateFrmMesh
!
! !ARGUMENTS:
    type(ESMF_Mesh), intent(in)     :: mesh
    type(ESMF_MeshLoc), intent(in)  :: meshLoc
    integer(ESMF_KIND_I4), optional :: maskValues(:)
    logical, optional               :: addOrigCoords
    integer, intent(out), optional  :: rc
!
! !DESCRIPTION:
!   Allocates memory for a new {\tt ESMF\_PointList} object and
!   constructs its internals from input Mesh.
!
!   The arguments are:
!   \begin{description}
!   \item[{mesh}]
!     The mesh to get the information from to create the PointList.
!   \item[{meshLoc}]
!     mesh location
!   \item[{maskValues}]
!     Values to set as masked
!   \item[{[addOrigCoords]}]
!     Also put the original (spherical) coordinate values in the PointList
!     Defaults to false.    
!   \item[{[rc]}]
!     Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!   \end{description}
!
!EOPI
!------------------------------------------------------------------------------
    integer             :: localrc      ! local return code
    type(ESMF_PointList)  :: pointlist
    type(ESMF_InterArray) :: maskValuesArg
    type(ESMF_Logical) :: opt_addOrigCoords
    
    ! initialize return code; assume routine not implemented
    localrc = ESMF_RC_NOT_IMPL
    if (present(rc)) rc = ESMF_RC_NOT_IMPL
    pointlist%this = ESMF_NULL_POINTER

    ! convert mask values
    maskValuesArg = ESMF_InterArrayCreate(maskValues, rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
        ESMF_CONTEXT, rcToReturn=rc)) return

    ! Convert addOrigCoords
    opt_addOrigCoords=ESMF_FALSE
    if (present(addOrigCoords)) opt_addOrigCoords=addOrigCoords
    
    ! Call C++ create code
    call c_ESMC_PointListCreateFrmMesh(mesh, meshLoc, maskValuesArg, opt_addOrigCoords, &         
         pointlist, localrc)
    if (ESMF_LogFoundError(localrc, &
      ESMF_ERR_PASSTHRU, &
      ESMF_CONTEXT, rcToReturn=rc)) return

    call ESMF_InterArrayDestroy(maskValuesArg, rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return

    ! Set return values
    ESMF_PointListCreateFrmMesh = pointlist

    ESMF_INIT_SET_CREATED(ESMF_PointListCreateFrmMesh)

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

  end function ESMF_PointListCreateFrmMesh