subroutine ESMF_GridspecQueryTileSize(filename, nx, ny, units, rc)
character(len=*), intent(in) :: filename
integer, intent(out) :: nx, ny
character(len=*), intent(out), optional :: units
integer, intent(out), optional :: rc
#ifdef ESMF_NETCDF
integer :: ncid, nvars, attlen, i
integer :: ncStatus
integer :: ndims, dimids(2)
character(len=128) :: attstr
character(len=1024) :: errmsg
if (present(rc)) rc=ESMF_SUCCESS
ncStatus = nf90_open(path=filename, mode=nf90_nowrite, ncid=ncid)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE, &
trim(filename), &
rc)) return
ncStatus = nf90_inquire(ncid, nVariables=nvars)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE, &
trim(filename), &
rc)) return
do i=1,nvars
! Check its standard_name attribute
ncStatus = nf90_inquire_attribute(ncid, i, 'standard_name', len=attlen)
if (ncStatus == nf90_noerr) then
ncStatus = nf90_get_att(ncid, i, 'standard_name', values=attstr)
if (attstr(1:attlen) .eq. 'geographic_longitude' .or. &
attstr(1:attlen) .eq. 'geographic_latitude') then
! read the longitude or latitude variable
! First find the dimension of this variable
ncStatus = nf90_inquire_variable(ncid, i, ndims=ndims, dimids=dimids)
if (ndims /= 2) then
call ESMF_LogSetError(rcToCheck=ESMF_FAILURE, &
msg="- The longitude variable should have dimension 2", &
ESMF_CONTEXT, rcToReturn=rc)
return
endif
! find out the dimenison size
ncStatus = nf90_inquire_dimension(ncid, dimids(1), len=nx)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE, &
"contact dimension inquire", &
rc)) return
ncStatus = nf90_inquire_dimension(ncid, dimids(2), len=ny)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE, &
"contact dimension inquire", &
rc)) return
! return the dimension of the center grid
nx = (nx-1)/2
ny = (ny-1)/2
! find out units attribute
if (present(units)) then
ncStatus = nf90_inquire_attribute(ncid, i, "units", len=attlen)
errmsg ="attribute units for coordinate variable" //" in "//trim(filename)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE, &
errmsg, &
rc)) return
ncStatus = nf90_get_att(ncid, i, 'units',attstr)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE, &
errmsg, &
rc)) return
if (attstr(1:6) .eq. 'degree') then
units = 'degrees'
elseif (attstr(1:6) .eq. 'radian') then
units = 'radians'
endif
endif
goto 20
endif
endif
enddo
20 continue
ncStatus = nf90_close(ncid)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE, &
"close tile file", &
rc)) return
return
#else
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_GridspecQueryTileSize