ESMF_UGridInq Subroutine

public subroutine ESMF_UGridInq(filename, meshname, nodeCount, elementCount, maxNodePElement, units, fillvalue, nodeCoordDim, faceCoordFlag, meshid, rc)

Arguments

Type IntentOptional Attributes Name
character(len=*), intent(in) :: filename
character(len=*), intent(in), optional :: meshname
integer, intent(out), optional :: nodeCount
integer, intent(out), optional :: elementCount
integer, intent(out), optional :: maxNodePElement
character(len=*), intent(out), optional :: units
integer, intent(out), optional :: fillvalue
integer, intent(out), optional :: nodeCoordDim
logical, intent(out), optional :: faceCoordFlag
integer, intent(out), optional :: meshid
integer, intent(out), optional :: rc

Source Code

subroutine ESMF_UGridInq(filename, meshname, nodeCount, elementCount, &
                        maxNodePElement, units, fillvalue, nodeCoordDim, &
                        faceCoordFlag, meshid, rc)

! !ARGUMENTS:

    character(len=*), intent(in)   :: filename
    character(len=*), intent(in), optional :: meshname
    integer, intent(out), optional :: nodeCount
    integer, intent(out), optional :: elementCount
    integer, intent(out), optional :: maxNodePElement
    character(len=*), intent(out), optional :: units
    integer, intent(out), optional :: fillvalue
    integer, intent(out), optional :: nodeCoordDim
    logical, intent(out), optional :: faceCoordFlag
    integer, intent(out), optional :: meshid
    integer, intent(out), optional :: rc

    integer:: localrc, ncStatus
    integer :: DimIds(2), VarId, dummyId
    integer :: ncid, local_rank, len
    character(len=256):: errmsg, nodeCoordString, nodeCoordNames(2), elmtConnName
    integer :: pos, meshDim
    character(len=80):: varname, attvalue
    integer :: i, nvars

    integer, parameter :: nf90_noerror = 0

#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 the dummy variable meshname which contains the topology data as attributes
    ! If the meshname is not given, find it in the file using its attribute cf_role or
    ! standard_name
    dummyId = 0
    if (present(meshname)) then
      ncStatus = nf90_inq_varid (ncid, trim(meshname), dummyId)
      errmsg = "Dummy Variable "//trim(meshname)//" does not exist in "//trim(filename)
      if (CDFCheckError (ncStatus, &
        ESMF_METHOD,  &
        ESMF_SRCLINE, errmsg, &
        rc)) return
    else
       ! find the mesh name using attribute "cf_role"
       ncStatus = nf90_inquire(ncid, nVariables=nvars)
       errmsg = "inquiry error with "//trim(filename)
       if (CDFCheckError (ncStatus, &
           ESMF_METHOD,  &
           ESMF_SRCLINE, errmsg, &
           rc)) return
       do i=1,nvars
          ncStatus=nf90_get_att(ncid, i, 'cf_role', attvalue)
          if (ncStatus == nf90_noerror) then
               ncStatus = nf90_inquire_attribute(ncid, i, 'cf_role', len=len)
          else
               ncStatus = nf90_get_att (ncid, i, "standard_name", attvalue)
               if (ncStatus == nf90_noerror) then
                  ncStatus = nf90_inquire_attribute(ncid, i, 'standard_name', len=len)
               endif
          endif
          if (ncStatus == nf90_noerror) then
               if (attvalue(len:len) .eq. achar(0)) len = len-1
               if (attvalue(1:len) .eq. 'mesh_topology') then
                 dummyId=i
               endif
          endif
       enddo
    endif
    if (dummyId == 0) then
       ! Mesh variable not found, return error
       errmsg = "- dummy mesh variable not found in "//trim(filename)
       call ESMF_LogSetError(rcToCheck=ESMF_FAILURE, &
             msg=errmsg, ESMF_CONTEXT, rcToReturn=rc)
       return
    endif

    if (present(meshId)) meshId=dummyId
    ! get dimension
    ! Change to topology_dimension based on the update on 2/28/2013 at
    ! http://publicwiki.deltares.nl/display/NETCDF/Deltares+CF+proposal+for+Unstructured+Grid+data+model
    ! for backward compatibility, use dimension if topology_dimension does not exist

    ncStatus = nf90_get_att (ncid, dummyId, "topology_dimension", values=meshDim)
    if (ncStatus /= nf90_noerror) then
       ncStatus = nf90_get_att (ncid, dummyId, "dimension", values=meshDim)
       errmsg = "Attribute topology_dimension or dimension in "//trim(filename)
       if (CDFCheckError (ncStatus, &
          ESMF_METHOD,  &
          ESMF_SRCLINE, errmsg, &
          rc)) return
    endif
    if (present(NodeCoordDim)) NodeCoordDim=meshDim
    ! get number of nodes
    if (present(nodeCount) .or. present(units)) then
      ncStatus = nf90_inquire_attribute(ncid, dummyId, "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, dummyId, "node_coordinates", nodeCoordString)
      if (CDFCheckError (ncStatus, &
        ESMF_METHOD,  &
        ESMF_SRCLINE, errmsg, &
        rc)) return
      pos = index(nodeCoordString(1:)," ")
      nodeCoordNames(1) = nodeCoordString(1:pos-1)

      ! Get dimension (# nodes) used to define the node coordinates
      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(nodeCount)) then
        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=nodeCount)
        if (CDFCheckError (ncStatus, &
          ESMF_METHOD,  &
          ESMF_SRCLINE, errmsg, &
          rc)) return
      endif
      if (present(units)) then
        ! get the attribute 'units'
        ncStatus = nf90_inquire_attribute(ncid, VarId, "units", len=len)
        errmsg = "Attribute units for "//nodeCoordNames(1)//" in "//trim(filename)
        if (CDFCheckError (ncStatus, &
            ESMF_METHOD, &
            ESMF_SRCLINE,&
            errmsg,&
            rc)) return
        ncStatus = nf90_get_att(ncid, VarId, "units", units)
        if (CDFCheckError (ncStatus, &
          ESMF_METHOD, &
          ESMF_SRCLINE,&
          errmsg,&
          rc)) return
        units = ESMF_UtilStringLowerCase(units(1:len))
        if (units(len:len) .eq. achar(0)) len = len-1
        units = units(1:7)
      endif
    end if

    if (present(elementCount) .or. present(maxNodePElement) .or. present(fillvalue)) then
      if (meshDim == 1) then
          errmsg = "- 1D network topology does not have elements" 
          call ESMF_LogSetError(rcToCheck=ESMF_FAILURE, &
                 msg=errmsg, ESMF_CONTEXT, rcToReturn=rc) 
          return 
      elseif (meshDim == 2) then
        varname = "face_node_connectivity"
      else
        varname = "volume_node_connectivity"
      endif
      ncStatus = nf90_inquire_attribute(ncid, dummyId, varname, len=len)
      errmsg = "Attribute face_node_connectivity in "//trim(filename)
      if (CDFCheckError (ncStatus, &
        ESMF_METHOD,  &
        ESMF_SRCLINE, errmsg, &
        rc)) return
      ncStatus = nf90_get_att (ncid, dummyId, 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 dimensions of element connectivity
      ncStatus = nf90_inquire_variable (ncid, VarId, dimids=DimIds)
      if (CDFCheckError (ncStatus, &
        ESMF_METHOD,  &
        ESMF_SRCLINE, errmsg, &
        rc)) return
      if (present(elementCount)) then   
        ncStatus = nf90_inquire_dimension (ncid, DimIds(2), len=elementCount)
        if (CDFCheckError (ncStatus, &
          ESMF_METHOD,  &
          ESMF_SRCLINE, errmsg, &
          rc)) return
      endif
      if (present(maxNodePElement)) then
        ncStatus = nf90_inquire_dimension (ncid, DimIds(1), len=maxNodePElement)
        if (CDFCheckError (ncStatus, &
          ESMF_METHOD,  &
          ESMF_SRCLINE, errmsg, &
          rc)) return
      endif
      if (present(fillvalue)) then
        ! Get elmt conn fill value (if it is mixed topology with different number of node
        ! per element.)
        errmsg = "Attribute "//trim(elmtConnName)//"_FillValue in "//trim(filename)
        ncStatus = nf90_get_att (ncid, VarId, "_FillValue", values=fillvalue)
        if (CDFCheckError (ncStatus, &
          ESMF_METHOD,  &
          ESMF_SRCLINE, errmsg, &
          rc)) return
      endif
   endif

   if (present(faceCoordFlag)) then
      ncStatus = nf90_inquire_attribute(ncid, dummyId, "face_coordinates", len=len)
      if (ncStatus /= nf90_noerror) then
          faceCoordFlag = .FALSE.
      else
          faceCoordFlag = .TRUE.
      endif
   endif

   ncStatus = nf90_close(ncid)
   if (CDFCheckError (ncStatus, &
      ESMF_METHOD,  &
      ESMF_SRCLINE, errmsg, &
      rc)) return
#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_UGridInq