subroutine ESMF_GridspecGetVar1D(grid_filename, varids, loncoord, latcoord, &
cornerlon, cornerlat, rc)
!
! !ARGUMENTS:
character(len=*), intent(in) :: grid_filename
integer, intent(in) :: varids(:)
real(ESMF_KIND_R8), intent(out) :: loncoord(:)
real(ESMF_KIND_R8), intent(out) :: latcoord(:)
real(ESMF_KIND_R8), intent(out), optional :: cornerlon(:,:)
real(ESMF_KIND_R8), intent(out), optional :: cornerlat(:,:)
integer, intent(out), optional:: rc
integer:: ncStatus
integer:: gridid, boundId, VarId, totaldims
integer:: len, i
character(len=256) :: errmsg, boundvar
character (len=80) :: attstr
#ifdef ESMF_NETCDF
! Open the grid files
ncStatus = nf90_open (path=trim(grid_filename), mode=nf90_nowrite, ncid=gridid)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE,&
trim(grid_filename), &
rc)) return
ncStatus = nf90_get_var(gridid, varids(1), loncoord)
errmsg = "longitude variable in "//trim(grid_filename)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE,&
errmsg, &
rc)) return
ncStatus = nf90_get_var(gridid, varids(2), latcoord)
errmsg = "latitude variable in "//trim(grid_filename)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE,&
errmsg, &
rc)) return
if (present(cornerlon)) then
! find the bound variable for lon
ncStatus = nf90_inquire_attribute(gridid, varids(1), "bounds", len=len)
errmsg = "attribute bounds in "//trim(grid_filename)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE,&
errmsg, &
rc)) return
ncStatus = nf90_get_att(gridid, varids(1), "bounds", boundvar)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE,&
errmsg, &
rc)) return
ncStatus = nf90_inq_varid(gridid, boundvar(1:len), boundId)
errmsg = "longitude bound variable in "//trim(grid_filename)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE,&
errmsg, &
rc)) return
ncStatus = nf90_get_var(gridid, boundId, cornerlon)
errmsg = "longitude bound variable in "//trim(grid_filename)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE,&
errmsg, &
rc)) return
endif
if (present(cornerlat)) then
! find the bound variable for lat
ncStatus = nf90_inquire_attribute(gridid, varids(2), "bounds", len=len)
errmsg = "attribute bounds in "//trim(grid_filename)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE,&
errmsg, &
rc)) return
ncStatus = nf90_get_att(gridid, varids(2), "bounds", boundvar)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE,&
errmsg, &
rc)) return
ncStatus = nf90_inq_varid(gridid, boundvar(1:len), boundId)
errmsg = "latitude bound variable in "//trim(grid_filename)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE,&
errmsg, &
rc)) return
ncStatus = nf90_get_var(gridid, boundId, cornerlat)
errmsg = "longitude bound variable in "//trim(grid_filename)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE,&
errmsg, &
rc)) return
endif
ncStatus = nf90_close(gridid)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE,&
trim(grid_filename),&
rc)) return
if(present(rc)) rc = ESMF_SUCCESS
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_GridspecGetVar1D