subroutine ESMF_GridspecGetVarByName(grid_filename, var_name, dimids, &
var_buffer, missing_value, start, count, rc)
character(len=*), intent(in) :: grid_filename
character(len=*), intent(in) :: var_name
integer, intent(in) :: dimids(:)
real(ESMF_KIND_R8), intent(out) :: var_buffer(:,:)
real(ESMF_KIND_R8), intent(out), optional:: missing_value
integer, intent(in), optional :: start(:), count(:)
integer, intent(out), optional:: rc
integer:: ncStatus
integer:: gridid, varid, ndims
integer:: vardimids(4)
integer:: len, i
integer, pointer:: lstart(:), lcount(:)
character(len=256) :: errmsg
integer, parameter :: nf90_noerror = 0
#ifdef ESMF_NETCDF
! Open the grid and mosaic 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_inq_varid( gridid, var_name, varid)
errmsg = "variable "//trim(var_name)// " in "//trim(grid_filename)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE,&
errmsg, &
rc)) return
! get missisng value
if (present(missing_value)) then
ncStatus = nf90_get_att(gridid, varid, "_FillValue", missing_value)
if (ncStatus /= nf90_noerror) then
ncStatus = nf90_get_att(gridid, varid, "missing_value", missing_value)
errmsg = "missing value attribute does not exist for "//trim(var_name)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE,&
errmsg, &
rc)) return
end if
end if
! inquire the dimension of the variable and make sure it matches with the dimids
! passed in
ncStatus = nf90_inquire_variable(gridid, varid, ndims=ndims, dimids=vardimids)
errmsg = "Variable "//trim(var_name)//" in "//trim(grid_filename)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE,&
errmsg, &
rc)) return
if (ndims < 2) then
call ESMF_LogSetError(rcToCheck=ESMF_FAILURE, &
msg="- variable dimension is less than 2", &
ESMF_CONTEXT, rcToReturn=rc)
return
else if (vardimids(1) /= dimids(1) .or. vardimids(2) /= dimids(2)) then
call ESMF_LogSetError(rcToCheck=ESMF_FAILURE, &
msg="- the first two dimensions does not match with the lat/lon dimension ids", &
ESMF_CONTEXT, rcToReturn=rc)
return
end if
! find the dimension length and check if it matches with the var array
! get variable
allocate(lstart(ndims))
allocate(lcount(ndims))
lstart(:)=1
lcount(:)=1
if (present(start) .and. .not. present(count)) then
call ESMF_LogSetError(rcToCheck=ESMF_FAILURE, &
msg="- Both start and count arguments have to be present at the same time", &
ESMF_CONTEXT, rcToReturn=rc)
return
endif
if (present(start)) then
lstart(1:2)=start
endif
if (present(count)) then
lcount(1:2)=count
else
ncStatus = nf90_inquire_dimension (gridid, dimids(1), len=lcount(1))
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE,&
errmsg,&
rc)) return
ncStatus = nf90_inquire_dimension (gridid, dimids(2), len=lcount(2))
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE,&
errmsg,&
rc)) return
endif
if (lcount(1) /= size(var_buffer, 1) .or. &
lcount(2) /= size(var_buffer, 2)) then
call ESMF_LogSetError(rcToCheck=ESMF_FAILURE, &
msg="- the variable array dimension does not match with dimension length", &
ESMF_CONTEXT, rcToReturn=rc)
return
endif
ncStatus = nf90_get_var(gridid, varid, var_buffer, start=lstart, &
count=lcount)
errmsg = "Variable "//trim(var_name)//" 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
deallocate(lstart, lcount)
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_GridspecGetVarByName