ESMF_MeshCreateFromUnstruct Function

private function ESMF_MeshCreateFromUnstruct(filename, convertToDual, fileformat, meshname, addUserArea, maskFlag, varname, rc)

Arguments

Type IntentOptional Attributes Name
character(len=*), intent(in) :: filename
logical, intent(in), optional :: convertToDual
type(ESMF_FileFormat_Flag), intent(in), optional :: fileformat
character(len=*), intent(in), optional :: meshname
logical, intent(in), optional :: addUserArea
type(ESMF_MeshLoc), intent(in), optional :: maskFlag
character(len=*), intent(in), optional :: varname
integer, intent(out), optional :: rc

Return Value type(ESMF_Mesh)


Source Code

    function ESMF_MeshCreateFromUnstruct(filename, convertToDual, fileformat, meshname, &
                        addUserArea, maskFlag, varname, rc)
!
!
! !RETURN VALUE:
    type(ESMF_Mesh)         :: ESMF_MeshCreateFromUnstruct
! !ARGUMENTS:
    character(len=*), intent(in)              :: filename
    logical, intent(in), optional               :: convertToDual
    type(ESMF_FileFormat_Flag), optional, intent(in) :: fileformat
    character(len=*), optional, intent(in)    :: meshname
    logical, intent(in), optional             :: addUserArea
    type(ESMF_MeshLoc), intent(in), optional  :: maskFlag
     character(len=*), optional, intent(in)    :: varname
   integer, intent(out), optional            :: rc
!
! !DESCRIPTION:
!   Create a mesh from a grid file defined in the ESMF Unstructured grid format.
!
!   \begin{description}
!   \item [filename]
!         The name of the grid file
!   \item[{[convertToDual]}]
!         if {\tt .true.}, the mesh will be converted to its dual. If not specified,
!         defaults to {\tt .false.}.
!   \item[{[addUserArea]}]
!         if {\tt .true.}, the cell area will be read in from the GRID file.  This feature is
!         only supported when the grid file is in the SCRIP or ESMF format.
!   \item [{[fileformat]}]
!         The type of grid file
!   \item[{[meshname]}]
!         The dummy variable for the mesh metadata in the UGRID file if the {\tt fileformat}
!         is {\tt ESMF\_FILEFORMAT\_UGRID}
!   \item[{[maskFlag]}]
!      If present, generate the mask using the missing\_value attribute defined in 'varname' on
!      the location defined by the flag, the accepted values are {\tt ESMF\_MESHLOC\_NODE} or
!      {\tt ESMF\_MESHLOC\_ELEMENT}
!   \item[{[varname]}]
!      If maskFlag is present, provide a variable name stored in the grid file and
!      the mask will be generated using the missing value of the data value of
!      this variable.  The first two dimensions of the variable has to be the
!      the longitude and the latitude dimension and the mask is derived from the
!      first 2D values of this variable even if this data is 3D, or 4D array.
!   \item [{[rc]}]
!         Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!   \end{description}
!
!EOPI
!------------------------------------------------------------------------------
    integer                             :: localrc      ! local return code
    integer                             :: PetNo, PetCnt
    real(ESMF_KIND_R8),pointer          :: nodeCoords(:,:), faceCoords(:,:)
    integer(ESMF_KIND_I4),pointer       :: elementConn(:)
    integer(ESMF_KIND_I4),pointer       :: elmtNum(:)
    integer                             :: startElmt
    integer                             :: NodeNo
    integer                             :: NodeCnt, total
    integer, allocatable                :: NodeId(:)
    integer, allocatable                :: NodeUsed(:)
    real(ESMF_KIND_R8), allocatable     :: NodeCoords1D(:)
    real(ESMF_KIND_R8), allocatable     :: NodeCoordsCart(:)
    real(ESMF_KIND_R8)                  :: coorX, coorY
    integer, allocatable                :: NodeOwners(:)
    integer, allocatable                :: NodeOwners1(:)
    integer, pointer                    :: glbNodeMask(:), NodeMask(:)

    integer                             :: ElemNo, TotalElements, startElemNo
    integer                             :: ElemCnt,i,j,k,dim, nedges
    integer                             :: localNodes, myStartElmt
    integer                             :: ConnNo, TotalConnects
    integer, allocatable                :: ElemId(:)
    integer, allocatable                :: ElemType(:)
    integer, allocatable                :: ElemConn(:)
    integer, pointer                    :: elementMask(:), ElemMask(:)
     real(ESMF_KIND_R8), pointer         :: elementArea(:), ElemArea(:)
    integer, allocatable                :: LocalElmTable(:)
    integer                             :: sndBuf(1)
    type(ESMF_VM)                       :: vm
    type(ESMF_Mesh)                     :: Mesh
    integer                             :: numPoly
#if 0
    integer, parameter                  :: maxNumPoly=20
    real(ESMF_KIND_R8)                  :: polyCoords(3*maxNumPoly)
    real(ESMF_KIND_R8)                  :: polyDblBuf(3*maxNumPoly)
    real(ESMF_KIND_R8)                  :: area(maxNumPoly)
    integer                             :: polyIntBuf(maxNumPoly)
    integer                             :: triInd(3*(maxNumPoly-2))
#else
    integer                             :: maxNumPoly
    real(ESMF_KIND_R8),allocatable      :: polyCoords(:)
    real(ESMF_KIND_R8),allocatable      :: polyDblBuf(:)
    real(ESMF_KIND_R8),allocatable      :: area(:)
    integer,allocatable                 :: polyIntBuf(:)
    integer,allocatable                 :: triInd(:)
#endif
    real(ESMF_KIND_R8)                  :: totalarea
    integer                             :: spatialDim
    integer                             :: parametricDim
    integer                             :: lni,ti,tk
    type(ESMF_FileFormat_Flag)          :: fileformatlocal
    integer                             :: coordDim
    logical                             :: convertToDeg
    logical                             :: haveNodeMask, haveElmtMask
    logical                             :: haveMask
    logical                             :: localAddUserArea
    logical                             :: localConvertToDual
    type(ESMF_MeshLoc)                  :: localAddMask
    real(ESMF_KIND_R8), pointer         :: varbuffer(:)
    real(ESMF_KIND_R8)                  :: missingvalue
    type(ESMF_CoordSys_Flag)            :: coordSys
    integer                             :: maxEdges
    logical                             :: hasFaceCoords
    logical                             :: haveOrigGridDims
    integer                             :: origGridDims(2)
    integer                 :: poleVal, minPoleGid, maxPoleGid,poleObjType

    ! Initialize return code; assume failure until success is certain
    localrc = ESMF_RC_NOT_IMPL
    if (present(rc)) rc = ESMF_RC_NOT_IMPL

    ! set faceCoords to null
    faceCoords => NULL()
    hasFaceCoords = .false.

    if (present(addUserArea)) then
        localAddUserArea = addUserArea
    else
        localAddUserArea = .false.
    endif

    if (present(convertToDual)) then
        localConvertToDual = convertToDual
    else
        localConvertToDual = .false.
    endif

    if (present(maskFlag)) then
        localAddMask = maskFlag
     else
        localAddMask = ESMF_MESHLOC_NONE
    endif

    ! Read the mesh definition from the file
    if (present(fileformat)) then
        fileformatlocal = fileformat
    else
        fileformatlocal = ESMF_FILEFORMAT_ESMFMESH
    endif

#if 0
    if (fileformatlocal == ESMF_FILEFORMAT_UGRID) then
        if (.not. present(meshname)) then
           call ESMF_LogSetError(ESMF_RC_ARG_WRONG, &
                             msg="- meshname argument is missing", &
                             ESMF_CONTEXT, rcToReturn=rc)
        endif
    endif
#endif

    ! get global vm information
    !
    call ESMF_VMGetCurrent(vm, rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return

    ! set up local pet info
    call ESMF_VMGet(vm, localPet=PetNo, petCount=PetCnt, rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return

    ! define default coordinate system
    coordSys = ESMF_COORDSYS_SPH_DEG

    ! define default haveOrigGridDims
    haveOrigGridDims=.false.

    ! Get grid info
    if (fileformatlocal == ESMF_FILEFORMAT_ESMFMESH) then
       ! Get coordDim
       call ESMF_EsmfInq(filename,coordDim=coordDim, haveNodeMask=haveNodeMask, &
                    haveElmtMask=haveElmtMask, maxNodePElement=maxEdges, &
                    haveOrigGridDims=haveOrigGridDims, rc=localrc)
       if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return


       ! Don't convert if not 2D because that'll be cartesian right now
       if (coordDim .eq. 2) then
          convertToDeg = .true.
       else
          convertToDeg = .false.
       endif

       ! Get OrigGridDims
       if (haveOrigGridDims) then
          call ESMF_EsmfInq(filename, origGridDims=origGridDims, rc=localrc)
          if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
               ESMF_CONTEXT, rcToReturn=rc)) return
       endif

       ! Get information from file
       ! Need to return the coordinate system for the nodeCoords
       if (haveNodeMask) then
           call ESMF_EsmfGetNode(filename, nodeCoords, nodeMask=glbNodeMask,&
                                convertToDeg=convertToDeg, coordSys=coordSys, rc=localrc)
       else
           call ESMF_EsmfGetNode(filename, nodeCoords, &
                                convertToDeg=convertToDeg, coordSys=coordSys, rc=localrc)
       endif

       if (haveElmtMask .and. localAddUserArea) then
            call ESMF_EsmfGetElement(filename, elementConn, elmtNum, &
                                 startElmt, elementMask=elementMask, elementArea=elementArea, &
                                 centerCoords=faceCoords, &
                                 convertToDeg=convertToDeg, rc=localrc)
       elseif (haveElmtMask) then
            call ESMF_EsmfGetElement(filename, elementConn, elmtNum, &
                                 startElmt, elementMask=elementMask, &
                                 centerCoords=faceCoords, &
                                 convertToDeg=convertToDeg, rc=localrc)
       elseif (localAddUserArea) then
            call ESMF_EsmfGetElement(filename, elementConn, elmtNum, &
                                 startElmt, elementArea=elementArea, &
                                 centerCoords=faceCoords, &
                                 convertToDeg=convertToDeg, rc=localrc)
       else
            call ESMF_EsmfGetElement(filename, elementConn, elmtNum, startElmt, &
                                 centerCoords=faceCoords, &
                                 convertToDeg=convertToDeg, rc=localrc)
       endif

       ElemCnt = ubound (elmtNum, 1)
       totalConnects = ubound(elementConn, 1)
       if (associated(faceCoords)) then
            hasFaceCoords = .true.
       endif
       if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                              ESMF_CONTEXT, rcToReturn=rc)) return
    elseif (fileformatlocal == ESMF_FILEFORMAT_UGRID) then
       haveElmtMask = .false.
       haveNodeMask = .false.
       if (localAddMask == ESMF_MESHLOC_ELEMENT) then
          haveElmtMask = .true.
       elseif (localAddMask == ESMF_MESHLOC_NODE) then
          haveNodeMask = .true.
       endif
       ! Get information from file
       call ESMF_GetMeshFromUGridFile(filename, nodeCoords, elementConn, &
                           elmtNum, startElmt, convertToDeg=.true., &
                           faceCoords=faceCoords, rc=localrc)
       if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                   ESMF_CONTEXT, rcToReturn=rc)) return
       ! Chenk if the grid is 3D or 2D
       coordDim = ubound(nodeCoords,1)
       nodeCnt = ubound(nodeCoords,2)
       ElemCnt = ubound (elmtNum, 1)
       totalConnects = ubound(elementConn, 1)

       if ( associated(faceCoords)) then
                  hasFaceCoords = .true.
       endif

       if (coordDim == 2 .and. localAddMask == ESMF_MESHLOC_ELEMENT) then
          !Get the variable and the missing value attribute from file
          ! Total number of local elements
          allocate(varbuffer(ElemCnt))
          call ESMF_UGridGetVarByName(filename, varname, varbuffer, startind=startElmt, &
                count=ElemCnt, location="face", &
                missingvalue=missingvalue, rc=localrc)
          if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                   ESMF_CONTEXT, rcToReturn=rc)) return
          ! Create local mask
          allocate(elementMask(ElemCnt))
          elementMask(:)=1
          do i=1,ElemCnt
            if (varbuffer(i) == missingvalue) elementMask(i)=0
          enddo
          deallocate(varbuffer)
       elseif (coordDim == 2 .and. localAddMask == ESMF_MESHLOC_NODE) then
          !Get the variable and the missing value attribute from file
          ! Total number of total nodes
          allocate(varbuffer(nodeCnt))
          call ESMF_UGridGetVarByName(filename, varname, varbuffer, &
                location="node", &
                missingvalue=missingvalue, rc=localrc)
          if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                   ESMF_CONTEXT, rcToReturn=rc)) return
          ! Create local mask
          allocate(glbNodeMask(nodeCnt))
          glbNodeMask(:)=1
          do i=1,nodeCnt
            if (varbuffer(i) == missingvalue) glbNodeMask(i)=0
          enddo
          deallocate(varbuffer)
       endif
    else
       call ESMF_LogSetError(ESMF_RC_ARG_WRONG, &
                             msg="- unrecognized fileformat", &
                             ESMF_CONTEXT, rcToReturn=rc)
       return
    endif

    nodeCnt = ubound(nodeCoords,2)

   ! Figure out dimensions
    if (coordDim .eq. 2) then
       parametricDim = 2
       spatialDim = 2
    else if (coordDim .eq. 3) then
       parametricDim = 3
       spatialDim = 3
    else
       call ESMF_LogSetError(ESMF_RC_VAL_OUTOFRANGE, &
            msg="- only coordDim 2 or 3 is supported right now", &
            ESMF_CONTEXT, rcToReturn=rc)
       return
    endif

    ! create the mesh
    Mesh = ESMF_MeshCreate3part (parametricDim, spatialDim, &
               coordSys=coordSys, rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
         ESMF_CONTEXT, rcToReturn=rc)) return

    ! These two arrays are temp arrays
    ! NodeUsed() used for multiple purposes, first, find the owners of the node
    ! later, used to store the local Node ID to be used in the ElmtConn table
    ! NodeOwners1() is the receiving array for ESMF_VMAllReduce(), it will store
    ! the lowest PET number that stores the node.  That PET will become the owner
    ! of the node.
    allocate (NodeUsed(NodeCnt))
    allocate (NodeOwners1(NodeCnt))

    ! Set to a number > PetCnt because it will store the PetNo if this node is used by
    ! the local elements and we will do a global reduce to find the minimal values
    NodeUsed(:)=PetCnt+100

    ! Set the coorsponding NodeUsed(:) value to my PetNo if it is used by the local element
    ! Also calculate the total number of mesh elements based on elmtNum
    totalElements = ElemCnt
    maxNumPoly=0
    if (parametricDim .eq. 2) then
       j=1
       do ElemNo = 1, ElemCnt
          do i=1,elmtNum(ElemNo)
            if (elementConn(j) /= ESMF_MESH_POLYBREAK) then
                NodeUsed(elementConn(j))=PetNo
             endif
             j=j+1
          enddo
          if (elmtNum(ElemNo) > maxNumPoly) then
             maxNumPoly=elmtNum(ElemNo)
          endif
       end do
    else ! If not parametricDim==2, assuming parmetricDim==3
       j=1
       do ElemNo =1, ElemCnt
          do i=1,elmtNum(ElemNo)
             if (elementConn(j) /= ESMF_MESH_POLYBREAK) then
                NodeUsed(elementConn(j))=PetNo
             endif
             j=j+1
          enddo
       end do
    endif

    if (totalConnects /= (j-1)) then
         print *, 'Total number of connection mismatch:', j, ElemCnt, totalConnects
    endif
   ! write(*,*) "maxNumPoly=",maxNumPoly

    ! Do a global reduce to find out the lowest PET No that owns each node, the result is in
    ! NodeOwners1(:)
    call ESMF_VMAllReduce(vm, NodeUsed, NodeOwners1, NodeCnt, ESMF_REDUCE_MIN, rc=localrc)
   if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return


    ! count number of nodes used and convert NodeUsed values into local index
    localNodes = 0
    do NodeNo = 1, NodeCnt
       if (NodeUsed(NodeNo) == PetNo) then
         localNodes = localNodes+1
         NodeUsed(NodeNo) = localNodes
       else
         NodeUsed(NodeNo) = 0
       endif
    enddo

    ! allocate nodes arrays for ESMF_MeshAddNodes()
    allocate (NodeId(localNodes))
    allocate (NodeOwners(localNodes))

    if (parametricDim .eq. 2) then
       allocate (NodeCoords1D(localNodes*coordDim))
    else ! If not parametricDim==2, assuming parmetricDim==3
       allocate(NodeCoords1D(localNodes*3))
    endif

    ! copy vertex information into nodes, NodeUsed(:) now contains either 0 (not for me) or
    ! the local node index.  The owner of the node is stored in NodeOwners1(:)
    ! Also calculate how many nodes are "owned" by me -- total
    i = 1
    total = 0
    if (parametricDim .eq. 2) then
       do NodeNo = 1, NodeCnt
          if (NodeUsed(NodeNo) > 0) then
             NodeId(i) = NodeNo
             do dim = 1, coordDim
                NodeCoords1D ((i-1)*coordDim+dim) = nodeCoords (dim, NodeNo)
             end do
             NodeOwners (i) = NodeOwners1(NodeNo)
             if (NodeOwners1(NodeNo) == PetNo) total = total+1
             i = i+1
          endif
       end do
    else ! If not parametricDim==2, assuming parmetricDim==3
       do NodeNo = 1, NodeCnt
          if (NodeUsed(NodeNo) > 0) then
             NodeId(i) = NodeNo
             do dim = 1, 3
                NodeCoords1D ((i-1)*3+dim) = nodeCoords(dim, NodeNo)
             end do
             NodeOwners (i) = NodeOwners1(NodeNo)
             if (NodeOwners1(NodeNo) == PetNo) total = total+1
             i = i+1
          endif
       end do
    endif


    deallocate(nodeCoords)
    if (.not. haveNodeMask) then
       ! Add nodes
       call ESMF_MeshAddNodes (Mesh, NodeIds=NodeId, &
                            NodeCoords=NodeCoords1D, &
                            NodeOwners=NodeOwners, &
                            rc=localrc)
    else
       allocate(NodeMask(localNodes))
       do i=1,localNodes
         NodeMask(i)=glbNodeMask(NodeId(i))
       enddo
       call ESMF_MeshAddNodes (Mesh, NodeIds=NodeId, &
                            NodeCoords=NodeCoords1D, &
                            NodeOwners=NodeOwners, &
                            NodeMask = NodeMask, &
                            rc=localrc)
       deallocate(NodeMask)
       deallocate(glbNodeMask)
    endif

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

    ! Need to calculate the total number of ESMF_MESH objects and the start element ID
    ! Do a global gather to get all the local TotalElements

    allocate(localElmTable(PetCnt))
    sndBuf(1)=TotalElements
    call ESMF_VMAllGather(vm, sndBuf, localElmTable, 1, rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return

    ! Find out the start element ID
    myStartElmt=0
    do i=1,PetNo
      myStartElmt = myStartElmt+localElmTable(i)
    end do
    deallocate(localElmTable)

    ! allocate element arrays for the local elements
    allocate (ElemId(TotalElements))
    allocate (ElemType(TotalElements))
    allocate (ElemConn(TotalConnects))
    if (localAddUserArea) allocate(ElemArea(TotalElements))



    ! Allocate mask if the user wants one
    haveMask=.false.
    if (haveElmtMask) then
       allocate (ElemMask(TotalElements))
       haveMask=.true.
    endif


    ! The ElemId is the global ID.  The myStartElmt is the starting Element ID(-1), and the
    ! element IDs will be from startElmt to startElmt+ElemCnt-1
    ! The ElemConn() contains the four corner node IDs for each element and it is organized
    ! as a 1D array.  The node IDs are "local", which are stored in NodeUsed(:)
    ElemNo = 1
    ConnNo = 0
    if (parametricDim .eq. 2) then
       ! Loop through creating Mesh appropriate elements
       k=1
       do j = 1, ElemCnt
          if (elmtNum(j)==3) then
             ElemId(ElemNo) = myStartElmt+ElemNo
             ElemType (ElemNo) = ESMF_MESHELEMTYPE_TRI
             do i=1,3
                if (elementConn(k) /= ESMF_MESH_POLYBREAK) then
                   ElemConn (ConnNo+i) = NodeUsed(elementConn(k))
                else
                   ElemConn (ConnNo+i) = ESMF_MESH_POLYBREAK
                endif
                k=k+1
             end do
             if (haveElmtMask) ElemMask(ElemNo) = elementMask(j)
             if (localAddUserArea) ElemArea(ElemNo) = elementArea(j)
             ElemNo=ElemNo+1
             ConnNo=ConnNo+3
          elseif (elmtNum(j)==4) then
             ElemId(ElemNo) = myStartElmt+ElemNo
             ElemType (ElemNo) = ESMF_MESHELEMTYPE_QUAD
             do i=1,4
                if (elementConn(k) /= ESMF_MESH_POLYBREAK) then
                   ElemConn (ConnNo+i) = NodeUsed(elementConn(k))
                else
                   ElemConn (ConnNo+i) = ESMF_MESH_POLYBREAK
                endif
                k=k+1
             end do
             if (haveElmtMask) ElemMask(ElemNo) = elementMask(j)
             if (localAddUserArea) ElemArea(ElemNo) = elementArea(j)
             ElemNo=ElemNo+1
             ConnNo=ConnNo+4
          else
             ElemId(ElemNo) = myStartElmt+ElemNo
             ElemType (ElemNo) = elmtNum(j)
             do i=1,elmtNum(j)
                if (elementConn(k) /= ESMF_MESH_POLYBREAK) then
                   ElemConn (ConnNo+i) = NodeUsed(elementConn(k))
                else
                   ElemConn (ConnNo+i) = ESMF_MESH_POLYBREAK
                endif
                k=k+1
             end do
             if (haveElmtMask) ElemMask(ElemNo) = elementMask(j)
             if (localAddUserArea) ElemArea(ElemNo) = elementArea(j)
             ElemNo=ElemNo+1
             ConnNo=ConnNo+elmtNum(j)
          endif
       enddo
    else ! If not parametricDim==2, assuming parmetricDim==3
       k=1
       do j = 1, ElemCnt
          if (elmtNum(j)==4) then
             ElemType (ElemNo) = ESMF_MESHELEMTYPE_TETRA
          elseif (elmtNum(j)==8) then
             ElemType (ElemNo) = ESMF_MESHELEMTYPE_HEX
          else
             call ESMF_LogSetError(ESMF_RC_VAL_OUTOFRANGE, &
                  msg="- in 3D currently only support Tetra. (4 nodes) or Hexa. (8 nodes)", &
                  ESMF_CONTEXT, rcToReturn=rc)
             return
          endif

          do i=1,elmtNum(j)
             if (elementConn(k) /= ESMF_MESH_POLYBREAK) then
                ElemConn (ConnNo+i) = NodeUsed(elementConn(k))
             else
                ElemConn (ConnNo+i) = ESMF_MESH_POLYBREAK
             endif
             k=k+1
          end do
          ElemId(ElemNo) = myStartElmt+ElemNo
          if (haveElmtMask) ElemMask(ElemNo) = elementMask(j)
          if (localAddUserArea) ElemArea(ElemNo) = elementArea(j)
          ElemNo=ElemNo+1
          ConnNo=ConnNo+elmtNum(j)
       end do
    endif

    if ((ElemNo /= TotalElements+1) .or. (ConnNo /= TotalConnects)) then
        write (ESMF_UtilIOStdout,*)  &
            PetNo, ' TotalElements does not match ',ElemNo-1, TotalElements, ConnNo, TotalConnects
    end if
    ! Add elements

    if (hasFaceCoords) then
       if (haveMask .and. localAddUserArea) then
            call ESMF_MeshAddElements (Mesh, ElemId, ElemType, ElemConn, &
                        elementMask=ElemMask, elementArea=ElemArea, &
                        elementCoords=reshape(faceCoords,(/totalElements*coordDim/)), rc=localrc)
       elseif (haveMask) then
            call ESMF_MeshAddElements (Mesh, ElemId, ElemType, ElemConn, &
                        elementMask=ElemMask, &
                        elementCoords=reshape(faceCoords,(/totalElements*coordDim/)), rc=localrc)
       elseif (localAddUserArea) then
            call ESMF_MeshAddElements (Mesh, ElemId, ElemType, ElemConn, &
                        elementArea=ElemArea, &
                        elementCoords=reshape(faceCoords,(/totalElements*coordDim/)), rc=localrc)
       else
            call ESMF_MeshAddElements (Mesh, ElemId, ElemType, ElemConn, &
                        elementCoords=reshape(faceCoords,(/totalElements*coordDim/)), rc=localrc)
       end if
    else
       if (haveMask .and. localAddUserArea) then
            call ESMF_MeshAddElements (Mesh, ElemId, ElemType, ElemConn, &
                        elementMask=ElemMask, elementArea=ElemArea, rc=localrc)
       elseif (haveMask) then
            call ESMF_MeshAddElements (Mesh, ElemId, ElemType, ElemConn, &
                        elementMask=ElemMask, rc=localrc)
       elseif (localAddUserArea) then
            call ESMF_MeshAddElements (Mesh, ElemId, ElemType, ElemConn, &
                        elementArea=ElemArea, rc=localrc)
       else
            call ESMF_MeshAddElements (Mesh, ElemId, ElemType, ElemConn, rc=localrc)
       end if
    endif

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

    ! Add pole info if applicable
    if ((fileformatlocal == ESMF_FILEFORMAT_ESMFMESH) .and. &
         haveOrigGridDims) then
       poleVal=4
       poleObjType=1 ! Set elements
       minPoleGid=1
       maxPoleGid=origGridDims(1)
       call C_ESMC_MeshSetPoles(Mesh, poleObjType, &
            poleVal, minPoleGid, maxPoleGid, localrc)
       if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return

       poleVal=5
       poleObjType=1 ! Set elements
       minPoleGid=origGridDims(1)*origGridDims(2)-origGridDims(1)+1
       maxPoleGid=origGridDims(1)*origGridDims(2)
       call C_ESMC_MeshSetPoles(Mesh, poleObjType, &
            poleVal, minPoleGid, maxPoleGid, localrc)
       if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return
    endif


    deallocate(NodeUsed, NodeId, NodeCoords1D, NodeOwners, NodeOwners1)
    deallocate(ElemId, ElemType, ElemConn, elementConn, elmtNum)
    if (haveElmtMask) deallocate(elementMask)
    if (haveMask) deallocate(ElemMask)
    if (associated(faceCoords)) deallocate(faceCoords)
    if (localAddUserArea) deallocate(elementArea, ElemArea)

    if (localConvertToDual) then
       ESMF_MeshCreateFromUnstruct = ESMF_MeshCreateDual(Mesh, rc=localrc)
        if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
               ESMF_CONTEXT, rcToReturn=rc)) return
    else
        ESMF_MeshCreateFromUnstruct = Mesh
    endif
    if (present(rc)) rc=ESMF_SUCCESS
    return
end function ESMF_MeshCreateFromUnstruct