ESMF_GridspecInq Subroutine

public subroutine ESMF_GridspecInq(grid_filename, ndims, grid_dims, coord_names, is3D, dimids, coordids, hasbound, units, rc)

Arguments

Type IntentOptional Attributes Name
character(len=*), intent(in) :: grid_filename
integer, intent(out) :: ndims
integer, intent(out) :: grid_dims(:)
character(len=*), intent(in), optional :: coord_names(:)
logical, intent(in), optional :: is3D
integer, intent(out), optional :: dimids(:)
integer, intent(out), optional :: coordids(:)
logical, intent(out), optional :: hasbound
character(len=*), intent(out), optional :: units
integer, intent(out), optional :: rc

Source Code

  subroutine ESMF_GridspecInq(grid_filename, ndims, grid_dims, coord_names, &
        is3D, dimids, coordids, hasbound, units, rc)

! !ARGUMENTS:

    character(len=*), intent(in)  :: grid_filename
    integer, intent(out)          :: ndims
    integer, intent(out)          :: grid_dims(:)
    logical, intent(in), optional :: is3D
    character(len=*), intent(in), optional :: coord_names(:)
    integer, intent(out), optional:: dimids(:)
    integer, intent(out), optional:: coordids(:)
    logical, intent(out), optional:: hasbound
    character(len=*), intent(out), optional ::units
    integer, intent(out), optional :: rc

    integer :: localrc, ncStatus
    integer :: varid, dimid
    integer :: varids(2)
    integer :: nvars, len
    integer :: gridid, i
    integer :: dimidslocal(2), localdimids(2)
    character (len=256) :: errmsg, boundvar
    character (len=80)  :: attstr, axisstr
    logical :: foundlon, foundlat
    logical :: useCoordName
    logical :: is3Dlocal;
    integer, parameter :: nf90_noerror = 0

#ifdef ESMF_NETCDF
    foundlat = .false.
    foundlon = .false.
    if (present(is3D)) then
      is3Dlocal = is3D
    else
      is3Dlocal = .false.
    endif
    ncStatus = nf90_open (path=trim(grid_filename), mode=nf90_nowrite, ncid=gridid)
    if (CDFCheckError (ncStatus, &
        ESMF_METHOD,  &
        ESMF_SRCLINE, &
        trim(grid_filename), &
        rc)) return

    if (present(coord_names)) then
        if (size(coord_names) /= 2) then
            call ESMF_LogSetError(ESMF_FAILURE, &
                 msg="- coord_names has to be of size 2", &
                 ESMF_CONTEXT, rcToReturn=rc)
            return
        endif
      do i=1,2
        ncStatus = nf90_inq_varid(gridid, coord_names(i), varid)
        errmsg ="variable "//trim(coord_names(1))//" in "//trim(grid_filename)
        if (CDFCheckError (ncStatus, &
                ESMF_METHOD,  &
                ESMF_SRCLINE, &
                errmsg, &
                rc)) return
        ! The coordinate variables are not in particular order.  It could be lon, lat
        ! or lat, lon.  We need to use its units value to decide whether it is longitude or
        ! latitude.
        ! First check units attribute for longitude, it has to be one of degrees_east,
        ! degree_east, degree_E, degreeE or degreesE
        !
        ! Added support for Cartesian coordinates in v8.0.0.  The Catesian coordinates
        ! should have attribute "axis" set to "X" or "Y" and units are "m" or
        ! "meters" for meters and "km" for kilometers.

        ncStatus = nf90_inquire_attribute(gridid, varid, "units", len=len)
        errmsg ="attribute units for "//trim(coord_names(1))//" in "//trim(grid_filename)
        if (CDFCheckError (ncStatus, &
                ESMF_METHOD,  &
                ESMF_SRCLINE, &
                errmsg, &
                rc)) return
        ncStatus = nf90_get_att(gridid, varid, 'units',attstr)
        if (CDFCheckError (ncStatus, &
                ESMF_METHOD,  &
                ESMF_SRCLINE, &
                errmsg, &
                rc)) return
        if (attstr(len:len) .eq. achar(0)) len = len-1
        if (attstr(1:len) .eq. 'degrees_east' .or. &
                attstr(1:len) .eq. 'degree_east' .or. &
                attstr(1:len) .eq. 'degree_E' .or. &
                attstr(1:len) .eq. 'degrees_E' .or. &
                attstr(1:len) .eq. 'degreeE' .or. &
                attstr(1:len) .eq. 'degreesE')  then
                foundlon = .true.
                varids(1)=varid
                if (present(units)) units = "degrees"
        endif
        ! It is not longitude, check for Latitude units
        if (attstr(1:len) .eq. 'degrees_north' .or. &
                attstr(1:len) .eq. 'degree_north' .or. &
                attstr(1:len) .eq. 'degree_N' .or. &
                attstr(1:len) .eq. 'degrees_N' .or. &
                attstr(1:len) .eq. 'degreeN' .or. &
                attstr(1:len) .eq. 'degreesN')  then
                varids(2)=varid
                foundlat = .true.
                if (present(units)) units = "degrees"
        endif
        if (attstr(1:len) .eq. 'm' .or. attstr(1:len) .eq. 'meters' .or. &
            attstr(1:len) .eq. 'km' .or. attstr(1:len) .eq. 'kilometers') then
            if (present(units)) units = attstr(1:len)
            ! check if the axis attribute exists
            ncStatus = nf90_get_att(gridid, varid, 'axis',axisstr)
            if (trim(axisstr) .eq. 'X') then
                foundlon = .true.
                varids(1)=varid
            elseif (trim(axisstr) .eq. 'Y') then
                varids(2)=varids(1)
                foundlat = .true.
                varids(2)=varid
            else
                call ESMF_LogSetError(ESMF_FAILURE, &
                   msg="- Not a valid coordinate variable.", &
                   ESMF_CONTEXT, rcToReturn=rc)
                return
            endif
        endif
      enddo
      if (.not. (foundlon .and. foundlat)) then
        call ESMF_LogSetError(ESMF_FAILURE, &
           msg="- Not a valid coordinate variable.", &
           ESMF_CONTEXT, rcToReturn=rc)
        return
      endif
    else
        ! get a list of variables in the file, inquire its standard_name and/or long_name to
        ! find out which one is longitude and which one is latitude variable
        ncStatus = nf90_inquire(gridid, nVariables = nvars)
        if (CDFCheckError (ncStatus, &
            ESMF_METHOD,  &
            ESMF_SRCLINE, &
            trim(grid_filename), &
            rc)) return

        foundlon = .false.
        foundlat  = .false.
        do i=1,nvars
           ! check if units attribute exists. If not, skip.  If yes, check if
           ! the value is degrees_east or degrees_north
           !
           ncStatus = nf90_inquire_attribute(gridid, i, "units", len=len)
           if (ncStatus /= nf90_noerror) CYCLE
           ncStatus = nf90_get_att(gridid, i, 'units',attstr)
           if (CDFCheckError (ncStatus, &
                ESMF_METHOD,  &
                ESMF_SRCLINE, &
                errmsg, &
                rc)) return

           if (attstr(len:len) .eq. achar(0)) len = len-1
           if (len >= 6 .and. (attstr(1:6) .eq. "degree")) then
              if (attstr(1:len) .eq. "degrees_east" .or. &
                  attstr(1:len) .eq. "degree_east" .or. &
                  attstr(1:len) .eq. "degree_E" .or. &
                  attstr(1:len) .eq. "degrees_E" .or. &
                  attstr(1:len) .eq. "degreeE" .or. &
                  attstr(1:len) .eq. "degreesE")  then
                  if (foundlon) then
                   call ESMF_LogSetError(ESMF_FAILURE, &
                      msg="- Duplicate longitude variables defined", &
                      ESMF_CONTEXT, rcToReturn=rc)
                  return
                else
                   varids(1)=i
                   foundlon = .true.
                   if (present(units)) units = "degrees"
                endif
              else if (attstr(1:len) .eq. "degrees_north" .or. &
                   attstr(1:len) .eq. "degree_north" .or. &
                   attstr(1:len) .eq. "degree_N" .or. &
                   attstr(1:len) .eq. "degrees_N" .or. &
                   attstr(1:len) .eq. "degreeN" .or. &
                   attstr(1:len) .eq. "degreesN")  then
                if (foundlat) then
                  call ESMF_LogSetError(ESMF_FAILURE, &
                     msg="- Duplicate latitude variables defined", &
                     ESMF_CONTEXT, rcToReturn=rc)
                  return
                else
                  varids(2)=i
                  foundlat = .true.
                  if (present(units)) units = "degrees"
                endif
              !else
          !        print *, "not the right units :", attstr(1:len), len
              endif
           endif
         enddo
    endif

    ! Check if "bounds" attribute is defined for the coordinate variable
    if (present(hasbound)) then
      ncStatus = nf90_get_att(gridid, varids(1), "bounds", boundvar)
      if (ncStatus /= nf90_noerror) then
         hasbound = .false.
      else
         hasbound = .true.
      endif
    endif

    ! find the dimension of the coordinate variables
    ncStatus = nf90_inquire_variable(gridid, varids(1), ndims=ndims, dimids=dimidslocal)
    errmsg ="longitude variable in "//trim(grid_filename)
    if (CDFCheckError (ncStatus, &
            ESMF_METHOD, &
            ESMF_SRCLINE,&
            errmsg,&
            rc)) return

    ! if the longitude is a 1D array, need to get the latitude dimension from the
    ! latitude variable
    if (ndims == 1) then
         if (present(dimids)) dimids(1)=dimidslocal(1)
         localdimids(1)=dimidslocal(1)
         ncStatus = nf90_inquire_variable(gridid, varids(2), ndims=ndims, dimids=dimidslocal)
         errmsg ="latitude variable in "//trim(grid_filename)
         if (CDFCheckError (ncStatus, &
            ESMF_METHOD, &
            ESMF_SRCLINE,&
            errmsg,&
            rc)) return
         if (present(dimids)) dimids(2)=dimidslocal(1)
         localdimids(2)=dimidslocal(1)
    else        
         if (present(dimids)) dimids = dimidslocal
         localdimids = dimidslocal
    endif
    ! find the dimension values
    ncStatus = nf90_inquire_dimension (gridid, localdimids(1), len=grid_dims(1))
    errmsg ="grid 1st dimension in "//trim(grid_filename)
    if (CDFCheckError (ncStatus, &
            ESMF_METHOD, &
            ESMF_SRCLINE,&
            errmsg,&
            rc)) return
    ncStatus = nf90_inquire_dimension (gridid, localdimids(2), len=grid_dims(2))
    errmsg ="grid 2nd dimension in "//trim(grid_filename)
    if (CDFCheckError (ncStatus, &
            ESMF_METHOD, &
            ESMF_SRCLINE,&
            errmsg,&
            rc)) return

    ncStatus = nf90_close(gridid)
    if (CDFCheckError (ncStatus, &
      ESMF_METHOD,  &
      ESMF_SRCLINE, &
      trim(grid_filename), &
      rc)) return

    if (present(coordids)) coordids = varids
    if (present(rc)) rc=ESMF_SUCCESS
    return
#else
    ndims = 0
    grid_dims = 0
    call ESMF_LogSetError(ESMF_RC_LIB_NOT_PRESENT, &
                 msg="- ESMF_NETCDF not defined when lib was compiled", &
                 ESMF_CONTEXT, rcToReturn=rc)
    return
#endif

end subroutine ESMF_GridspecInq