ESMF_EsmfGetNode Subroutine

public subroutine ESMF_EsmfGetNode(filename, nodeCoords, nodeMask, convertToDeg, coordSys, rc)

Arguments

Type IntentOptional Attributes Name
character(len=*), intent(in) :: filename
real(kind=ESMF_KIND_R8), pointer :: nodeCoords(:,:)
integer(kind=ESMF_KIND_I4), optional, pointer :: nodeMask(:)
logical, intent(in), optional :: convertToDeg
type(ESMF_CoordSys_Flag), optional :: coordSys
integer, intent(out), optional :: rc

Source Code

subroutine ESMF_EsmfGetNode (filename, nodeCoords, nodeMask, &
                            convertToDeg, coordSys, rc)

    character(len=*), intent(in)   :: filename
    real(ESMF_KIND_R8), pointer    :: nodeCoords (:,:)
    integer(ESMF_KIND_I4), pointer, optional :: nodeMask (:)
    logical, intent(in), optional  :: convertToDeg
    type(ESMF_CoordSys_Flag), optional :: coordSys
    integer, intent(out), optional :: rc

    type(ESMF_VM) :: vm
    integer :: PetNo, PetCnt

    integer :: ncid
    integer :: ncStatus
    integer :: RecCnt (2)

    integer :: DimId
    integer :: nodeCnt, ElmtCount, MaxNodePerElmt, NodeDim
    integer :: localCount, remain

    integer :: VarNo
    character(len=256)::errmsg
    character(len=80) :: units
    integer :: len
    logical :: convertToDegLocal
    type(ESMF_CoordSys_Flag) :: coordSysLocal
    integer :: memstat

#ifdef ESMF_NETCDF
    coordSysLocal = ESMF_COORDSYS_SPH_DEG
    convertToDegLocal = .false.
    if (present(convertToDeg)) convertToDegLocal = convertToDeg

    call ESMF_VMGetCurrent(vm, rc=rc)
    if (rc /= ESMF_SUCCESS) return
    ! set up local pet info
    call ESMF_VMGet(vm, localPet=PetNo, petCount=PetCnt, rc=rc)
    if (rc /= ESMF_SUCCESS) return

    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 vertices
    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=nodeCnt)
    if (CDFCheckError (ncStatus, &
      ESMF_METHOD,  &
      ESMF_SRCLINE, errmsg, &
      rc)) return

    ! Get vertex dimension
    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=NodeDim)
    if (CDFCheckError (ncStatus, &
      ESMF_METHOD,  &
      ESMF_SRCLINE, errmsg, &
      rc)) return

    ! allocate memory for verticies
    allocate (nodeCoords (NodeDim, nodeCnt), stat=memstat)
    if (ESMF_LogFoundAllocError(memstat,  &
        ESMF_CONTEXT, rcToReturn=rc)) return


    RecCnt(:) = ubound(nodeCoords)
    !print *, "nodeCoords:",nodeCnt, NodeDim

    ! read vertex data
    ncStatus = nf90_inq_varid (ncid, "nodeCoords", VarNo)
    errmsg = "Variable nodeCoords in "//trim(filename)
    if (CDFCheckError (ncStatus, &
      ESMF_METHOD,  &
      ESMF_SRCLINE, errmsg, &
      rc)) return

    ncStatus = nf90_get_var (ncid, VarNo, nodeCoords, start=(/1,1/), count=(/NodeDim, nodeCnt/))
    errmsg = "Variable nodeCoords in "//trim(filename)
    if (CDFCheckError (ncStatus, &
      ESMF_METHOD,  &
      ESMF_SRCLINE, errmsg, &
      rc)) return

    ! Check units, and set the coordSys output, if nodeDim==3, normalize the height
    ! field
    !if (nodeDim==2) then
       ncStatus = nf90_inquire_attribute(ncid, VarNo, "units", len=len)
       if (CDFCheckError (ncStatus, &
            ESMF_METHOD, &
            ESMF_SRCLINE,errmsg,&
            rc)) return

       ncStatus = nf90_get_att(ncid, VarNo, "units", units)
       if (CDFCheckError (ncStatus, &
            ESMF_METHOD, &
            ESMF_SRCLINE,errmsg,&
            rc)) return
       ! if len != 7, something is wrong, check the value.  If it starts
       ! with Degres/degrees/Radians/radians, ignore the garbage after the
       ! word.  Otherwise, return the whole thing
       if (units(len:len) .eq. achar(0)) len = len-1
       units = ESMF_UtilStringLowerCase(units(1:len))
       ! if the units is meters, kilometers, or km, make it Cartisian 2D
       if (units(1:len) .eq. "meters" .or. units(1:len) .eq. "m" .or. &
           units(1:len) .eq. "km" .or. units(1:len) .eq. "kilometers") then     
           coordSysLocal = ESMF_COORDSYS_CART
       elseif (units(1:len) .eq. 'radians' .and. .not. convertToDegLocal) then
            coordSysLocal = ESMF_COORDSYS_SPH_RAD
       elseif (units(1:len) .ne. 'degrees' .and. &
            units(1:len) .ne. 'radians') then
          call ESMF_LogSetError(rcToCheck=ESMF_FAILURE, &
                 msg="- units attribute is not degrees or radians", &
                 ESMF_CONTEXT, rcToReturn=rc)
          return
       endif

       ! if units is "radians", convert it to degree
       ! if nodeDim=3, only convert the first two dimensions, leave the 3rd dim unchanged
       if (convertToDegLocal) then
          if (units(1:len) .eq. "radians") then
             nodeCoords(1:2,:) = &
                 nodeCoords(1:2,:)*ESMF_COORDSYS_RAD2DEG
          endif
       endif

    !endif

    ! get nodeMask
    if (present(nodeMask)) then
       allocate(nodeMask(nodeCnt), stat=memstat)
       if (ESMF_LogFoundAllocError(memstat,  &
           ESMF_CONTEXT, rcToReturn=rc)) return

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

       ncStatus = nf90_get_var (ncid, VarNo, nodeMask, start=(/1/), count=(/nodeCnt/))
       if (CDFCheckError (ncStatus, &
          ESMF_METHOD,  &
          ESMF_SRCLINE, errmsg, &
          rc)) return
    endif

    ncStatus = nf90_close (ncid=ncid)
    if (CDFCheckError (ncStatus, &
      ESMF_METHOD,  &
      ESMF_SRCLINE, trim(filename), &
      rc)) return

    if (present(coordSys)) then
       coordSys = coordSysLocal
    endif
#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_EsmfGetNode