ESMF_EsmfInq Subroutine

public subroutine ESMF_EsmfInq(filename, nodeCount, elementCount, maxNodePElement, coordDim, haveNodeMask, haveElmtMask, haveArea, haveOrigGridDims, origGridDims, rc)

Arguments

Type IntentOptional Attributes Name
character(len=*), intent(in) :: filename
integer, intent(out), optional :: nodeCount
integer, intent(out), optional :: elementCount
integer, intent(out), optional :: maxNodePElement
integer, intent(out), optional :: coordDim
logical, intent(out), optional :: haveNodeMask
logical, intent(out), optional :: haveElmtMask
logical, intent(out), optional :: haveArea
logical, intent(out), optional :: haveOrigGridDims
integer, intent(out), optional :: origGridDims(2)
integer, intent(out), optional :: rc

Source Code

subroutine ESMF_EsmfInq(filename, nodeCount, elementCount, &
                        maxNodePElement, coordDim,  &
                        haveNodeMask, haveElmtMask, haveArea, &
                        haveOrigGridDims, origGridDims, rc)

! !ARGUMENTS:

    character(len=*), intent(in)   :: filename
    integer, intent(out), optional :: nodeCount
    integer, intent(out), optional :: elementCount
    integer, intent(out), optional :: maxNodePElement
    integer, intent(out), optional :: coordDim
    logical, intent(out), optional :: haveNodeMask
    logical, intent(out), optional :: haveElmtMask
    logical, intent(out), optional :: haveArea
    logical, intent(out), optional :: haveOrigGridDims
    integer, intent(out), optional :: origGridDims(2)
    integer, intent(out), optional :: rc

    integer:: localrc, ncStatus
    integer :: DimId, VarId
    integer :: ncid, local_rank
    character(len=256):: errmsg
    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 number of elements
    if (present(nodeCount)) then
      ncStatus = nf90_inq_dimid (ncid, "nodeCount", DimId)
      errmsg = "Dimension nodeCount in "//trim(filename)
      if (CDFCheckError (ncStatus, &
        ESMF_METHOD, &
        ESMF_SRCLINE,errmsg,&
        rc)) return

      ncStatus = nf90_inquire_dimension (ncid, DimId, len=nodeCount)
      if (CDFCheckError (ncStatus, &
        ESMF_METHOD, &
        ESMF_SRCLINE,errmsg,&
        rc)) return
    end if

    ! get number of elements
    if (present(elementCount)) then
      ncStatus = nf90_inq_dimid (ncid, "elementCount", DimId)
      errmsg = "Dimension elementCount in "//trim(filename)
      if (CDFCheckError (ncStatus, &
        ESMF_METHOD, &
        ESMF_SRCLINE,errmsg,&
        rc)) return

      ncStatus = nf90_inquire_dimension (ncid, DimId, len=elementCount)
      if (CDFCheckError (ncStatus, &
        ESMF_METHOD, &
        ESMF_SRCLINE,errmsg,&
        rc)) return
    end if

    ! get number of elements
    if (present(maxNodePElement)) then
      ncStatus = nf90_inq_dimid (ncid, "maxNodePElement", DimId)
      errmsg = "Dimension maxNodePElement in "//trim(filename)
      if (ncStatus /= nf90_noerror) then
         maxNodePElement = -1
      else
         ncStatus = nf90_inquire_dimension (ncid, DimId, len=maxNodePElement)
         if (CDFCheckError (ncStatus, &
           ESMF_METHOD, &
           ESMF_SRCLINE,errmsg,&
           rc)) return
      end if
    end if

    ! get number of elements
    if (present(coordDim)) then
      ncStatus = nf90_inq_dimid (ncid, "coordDim", DimId)
      errmsg = "Dimension coordDim in "//trim(filename)
      if (CDFCheckError (ncStatus, &
        ESMF_METHOD, &
        ESMF_SRCLINE,errmsg,&
        rc)) return

      ncStatus = nf90_inquire_dimension (ncid, DimId, len=coordDim)
      if (CDFCheckError (ncStatus, &
        ESMF_METHOD, &
        ESMF_SRCLINE,errmsg,&
        rc)) return
    end if

    ! check if elementMask exit
    if (present(haveElmtMask)) then
      ncStatus = nf90_inq_varid (ncid, "elementMask", VarId)
      if (ncStatus == nf90_noerror) then
           haveElmtMask = .true.
      else
           haveElmtMask = .false.
      end if
    end if

    ! check if nodeMask exit
    if (present(haveNodeMask)) then
      ncStatus = nf90_inq_varid (ncid, "nodeMask", VarId)
      if (ncStatus == nf90_noerror) then
           haveNodeMask = .true.
      else
           haveNodeMask = .false.
      end if
    end if

    ! check if elementArea exit
    if (present(haveArea)) then
      ncStatus = nf90_inq_varid (ncid, "elementArea", VarId)
      if (ncStatus == nf90_noerror) then
           haveArea = .true.
      else
           haveArea = .false.
      end if
    end if

    ! check if haveOrigGridDims
    if (present(haveOrigGridDims)) then
      ncStatus = nf90_inq_varid (ncid, "origGridDims", VarId)
      if (ncStatus == nf90_noerror) then
           haveOrigGridDims = .true.
      else
           haveOrigGridDims = .false.
      end if
    end if

! XMRKX !
    ! Get origGridDims
    if (present(origGridDims)) then
      ncStatus = nf90_inq_varid (ncid, "origGridDims", VarId)
      errmsg = "Variable origGridDims in "//trim(filename)
      if (CDFCheckError (ncStatus, &
        ESMF_METHOD, &
        ESMF_SRCLINE,errmsg,&
        rc)) return

      ncStatus=nf90_get_var(ncid, VarId, origGridDims)
      if (CDFCheckError (ncStatus, &
           ESMF_METHOD, &
           ESMF_SRCLINE,errmsg,&
           rc)) return
    end if

    if (present(rc)) rc=ESMF_SUCCESS
    return
#else
    if (ESMF_LogFoundError(ESMF_RC_LIB_NOT_PRESENT, &
                msg="- ESMF_NETCDF not defined when lib was compiled", &
                ESMF_CONTEXT, rcToReturn=rc)) return
#endif

end subroutine ESMF_EsmfInq