ESMF_UGridGetVar Subroutine

public subroutine ESMF_UGridGetVar(filename, meshId, nodeXcoords, nodeYcoords, faceXcoords, faceYcoords, faceNodeConnX, faceNodeConnY, rc)

Arguments

Type IntentOptional Attributes Name
character(len=*), intent(in) :: filename
integer :: meshId
real(kind=ESMF_KIND_R8), optional, pointer :: nodeXcoords(:)
real(kind=ESMF_KIND_R8), optional, pointer :: nodeYcoords(:)
real(kind=ESMF_KIND_R8), optional, pointer :: faceXcoords(:)
real(kind=ESMF_KIND_R8), optional, pointer :: faceYcoords(:)
real(kind=ESMF_KIND_R8), optional, pointer :: faceNodeConnX(:,:)
real(kind=ESMF_KIND_R8), optional, pointer :: faceNodeConnY(:,:)
integer, intent(out), optional :: rc

Source Code

subroutine ESMF_UGridGetVar (filename, meshId, &
                nodeXcoords, nodeYcoords, &
                faceXcoords, faceYcoords, &
                faceNodeConnX, faceNodeConnY, rc)
! !INTERFACE:
    character(len=*), intent(in)   :: filename
    integer                        :: meshId
    real(ESMF_KIND_R8), optional, pointer :: nodeXcoords(:), nodeYcoords(:)
    real(ESMF_KIND_R8), optional, pointer :: faceXcoords(:), faceYcoords(:)
    real(ESMF_KIND_R8), optional, pointer :: faceNodeConnX(:,:),faceNodeConnY(:,:)
    integer, intent(out), optional :: rc

    integer :: i,j,dim1,dim2
    integer, allocatable :: elemConn(:,:)
    integer:: localrc, ncStatus
    integer :: VarId, meshDim, pos1, pos2, len
    integer :: ncid, local_rank
    integer :: localFillValue, indexBase, offset
    real(ESMF_KIND_R8), pointer :: nodeXcoordsLocal(:), nodeYcoordsLocal(:)
    integer :: dimIds(1), nodeDim
    character(len=256):: errmsg, nodeCoordString, nodeCoordNames(2), elmtConnName
    integer, parameter :: nf90_noerror = 0
    character(len=80) :: varname

#ifdef ESMF_NETCDF
    if (present(rc)) rc=ESMF_SUCCESS
    ncStatus = nf90_open (path=trim(filename), mode=nf90_nowrite, ncid=ncid)
    if (CDFCheckError (ncStatus, &
      ESMF_METHOD, &
      ESMF_SRCLINE, trim(filename), rc)) return

    ! get dimension
    ncStatus = nf90_get_att (ncid, meshId, "topology_dimension", values=meshDim)
    if (ncStatus /= nf90_noerror) then
       ncStatus = nf90_get_att (ncid, meshId, "dimension", values=meshDim)
       errmsg = "Attribute topology_dimension or dimension in "//trim(filename)
       if (CDFCheckError (ncStatus, &
          ESMF_METHOD,  &
          ESMF_SRCLINE, errmsg, &
          rc)) return
    endif

    ! get number of nodes
    if (present(nodeXcoords) .or. present(nodeYcoords) .or. present(faceNodeConnX)) then
      ncStatus = nf90_inquire_attribute(ncid, meshId, "node_coordinates", len=len)
      errmsg = "Attribute node_coordinates in "//trim(filename)
      if (CDFCheckError (ncStatus, &
        ESMF_METHOD,  &
        ESMF_SRCLINE, errmsg, &
        rc)) return
      ncStatus = nf90_get_att (ncid, meshId, "node_coordinates", nodeCoordString)
      if (CDFCheckError (ncStatus, &
        ESMF_METHOD,  &
        ESMF_SRCLINE, errmsg, &
        rc)) return
      pos1 = index(nodeCoordString(1:)," ")
      nodeCoordNames(1) = nodeCoordString(1:pos1-1)
      pos2 = index(nodeCoordString(pos1+1:)," ")
      nodeCoordNames(2)=nodeCoordString(pos1+1:pos1+pos2-1)

      errmsg = "Variable "//trim(nodeCoordNames(1))//" in "//trim(filename)
      ncStatus = nf90_inq_varid (ncid, nodeCoordNames(1), VarId)
      if (CDFCheckError (ncStatus, &
        ESMF_METHOD,  &
        ESMF_SRCLINE, errmsg, &
        rc)) return

      ! find out the node dimension and allocate local nodeCoord arrays
      ncStatus = nf90_inquire_variable(ncid, VarId, dimids=dimIds)
      if (CDFCheckError (ncStatus, &
        ESMF_METHOD,  &
        ESMF_SRCLINE, errmsg, &
        rc)) return
      ncStatus = nf90_inquire_dimension (ncid, DimIds(1), len=nodeDim)
      errmsg = "Dimension in "//trim(filename)
      if (CDFCheckError (ncStatus, &
        ESMF_METHOD,  &
        ESMF_SRCLINE, errmsg, &
        rc)) return
      allocate(nodeXcoordsLocal(nodeDim), nodeYcoordsLocal(nodeDim))
      ncStatus = nf90_get_var (ncid, VarId, nodeXcoordsLocal)
      if (CDFCheckError (ncStatus, &
        ESMF_METHOD,  &
        ESMF_SRCLINE, errmsg, &
        rc)) return
      errmsg = "Variable "//trim(nodeCoordNames(2))//" in "//trim(filename)
      ncStatus = nf90_inq_varid (ncid, nodeCoordNames(2), VarId)
      if (CDFCheckError (ncStatus, &
        ESMF_METHOD,  &
        ESMF_SRCLINE, errmsg, &
        rc)) return
      ncStatus = nf90_get_var (ncid, VarId, nodeYcoordsLocal)
      if (CDFCheckError (ncStatus, &
        ESMF_METHOD,  &
        ESMF_SRCLINE, errmsg, &
        rc)) return
      if (present(nodeXcoords)) nodeXcoords = nodeXcoordsLocal
      if (present(nodeYcoords)) nodeYcoords = nodeYcoordsLocal

      ! faceNodeConnX and faceNodeConnY are the node coordinates for a given cell (element).
      ! In the UGRID standard, face_node_connectivity variable contains the indices
      ! to the nodes, rather than the actually coodinates.  This is used to write the weight
      ! file in the SCRIP format.
      if (meshDim == 2) then
        varname = "face_node_connectivity"
      else
        varname = "volume_node_connectivity"
      endif
      if (present(faceNodeConnX)) then
        errmsg = "Attribute "//trim(varname)//" in "//trim(filename)
        ncStatus = nf90_inquire_attribute(ncid, meshId, trim(varname), len=len)
        if (CDFCheckError (ncStatus, &
          ESMF_METHOD,  &
          ESMF_SRCLINE, errmsg, &
          rc)) return
        ncStatus = nf90_get_att (ncid, meshId, trim(varname), values=elmtConnName)
        if (CDFCheckError (ncStatus, &
          ESMF_METHOD,  &
          ESMF_SRCLINE, errmsg, &
          rc)) return
        ! Get element connectivity
        errmsg = "Variable "//elmtConnName(1:len)//" in "//trim(filename)
        ncStatus = nf90_inq_varid (ncid, elmtConnName(1:len), VarId)
        if (CDFCheckError (ncStatus, &
          ESMF_METHOD,  &
          ESMF_SRCLINE, errmsg, &
          rc)) return
        ! Get elmt conn fill value
        errmsg = "Attribute "//elmtConnName(1:len)//"_FillValue in "//trim(filename)
        ncStatus = nf90_get_att (ncid, VarId, "_FillValue", values=localFillValue)
        if (ncStatus /= nf90_noerror) localFillValue = -1
        ! Get start_index attribute to find out the index base (0 or 1)
        ncStatus = nf90_get_att (ncid, VarId, "start_index", values=indexBase)
        ! if not defined, default to 0-based
        if (ncStatus /= nf90_noerror) indexBase = 0

        dim1 = size(faceNodeConnX,1)
        dim2 = size(faceNodeConnX,2)
        allocate(elemConn(dim1,dim2) )
        ncStatus = nf90_get_var (ncid, VarId, elemConn)
        if (CDFCheckError (ncStatus, &
          ESMF_METHOD,  &
          ESMF_SRCLINE, errmsg, &
          rc)) return
        faceNodeConnX(:,:)=localFillValue
        faceNodeConnY(:,:)=localFillValue
        !adjust index to 1-based of the start_index is 0
        if (indexBase == 1) then
           offset = 0
        else
           offset = 1
        endif
        do i=1,dim1
          do j=1,dim2
                if (elemConn(i,j) /= localFillValue) then
                   faceNodeConnX(i,j) = nodeXcoordsLocal(elemConn(i,j)+offset)
                   faceNodeConnY(i,j) = nodeYcoordsLocal(elemConn(i,j)+offset)
                endif
          enddo
        enddo
        deallocate(elemConn)
      endif
      deallocate(nodeXcoordsLocal, nodeYcoordsLocal)
    endif

    ! get number of face
    if (present(faceXcoords) .or. present(faceYcoords)) then
      ncStatus = nf90_inquire_attribute(ncid, meshId, "face_coordinates", len=len)
      errmsg = "Attribute face_coordinates in "//trim(filename)
      if (CDFCheckError (ncStatus, &
        ESMF_METHOD,  &
        ESMF_SRCLINE, errmsg, &
        rc)) return
      ncStatus = nf90_get_att (ncid, meshId, "face_coordinates", nodeCoordString)
      if (CDFCheckError (ncStatus, &
        ESMF_METHOD,  &
        ESMF_SRCLINE, errmsg, &
        rc)) return
      pos1 = index(nodeCoordString(1:)," ")
      nodeCoordNames(1) = nodeCoordString(1:pos1-1)
      pos2 = index(nodeCoordString(pos1+1:)," ")
      nodeCoordNames(2)=nodeCoordString(pos1+1:pos1+pos2-1)

      errmsg = "Variable "//trim(nodeCoordNames(1))//" in "//trim(filename)
      ncStatus = nf90_inq_varid (ncid, nodeCoordNames(1), VarId)
      if (CDFCheckError (ncStatus, &
        ESMF_METHOD,  &
        ESMF_SRCLINE, errmsg, &
        rc)) return
      if (present(faceXcoords)) then
        ncStatus = nf90_get_var (ncid, VarId, faceXcoords)
        if (CDFCheckError (ncStatus, &
          ESMF_METHOD,  &
          ESMF_SRCLINE, errmsg, &
          rc)) return
      endif
      if (present(faceYcoords)) then
        errmsg = "Variable "//trim(nodeCoordNames(2))//" in "//trim(filename)
        ncStatus = nf90_inq_varid (ncid, nodeCoordNames(2), VarId)
        if (CDFCheckError (ncStatus, &
          ESMF_METHOD,  &
          ESMF_SRCLINE, errmsg, &
          rc)) return
        ncStatus = nf90_get_var (ncid, VarId, faceYcoords)
        if (CDFCheckError (ncStatus, &
          ESMF_METHOD,  &
          ESMF_SRCLINE, errmsg, &
          rc)) return
      endif
    endif

    ncStatus = nf90_close(ncid)

#else
    call ESMF_LogSetError(rcToCheck=ESMF_RC_LIB_NOT_PRESENT, &
                 msg="- ESMF_NETCDF not defined when lib was compiled", &
                 ESMF_CONTEXT, rcToReturn=rc)
#endif

    return
    end subroutine ESMF_UGridGetVar