subroutine ESMF_UGridGetVarByName (filename, varname, varbuffer, &
startind, count, location, missingvalue, rc)
character(len=*), intent(in) :: filename
character(len=*), intent(in) :: varname
real(ESMF_KIND_R8), pointer :: varbuffer (:)
integer, intent(in), optional :: startind
integer, intent(in), optional :: count
character(len=*), intent(in), optional:: location
real(ESMF_KIND_R8), intent(out), optional:: missingvalue
integer :: rc
integer :: ncid, VarId
integer :: ncStatus, localrc
character(len=256) :: errmsg
character(len=80) :: attstr
integer :: ndims, dimids(3), dimsize, length
integer, pointer :: starts(:), counts(:)
integer :: locflag
integer, parameter :: nf90_noerror=0
#ifdef ESMF_NETCDF
if (present(startind)) then
if (.not. present(count)) then
call ESMF_LogSetError(rcToCheck=ESMF_FAILURE, &
msg="- optional argument count is missing", &
ESMF_CONTEXT, rcToReturn=rc)
return
endif
endif
ncStatus = nf90_open (path=trim(filename), mode=nf90_nowrite, ncid=ncid)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE, trim(filename), &
rc)) return
errmsg = "Variable "//varname//" in "//trim(filename)
ncStatus = nf90_inq_varid (ncid, varname, VarId)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE, errmsg, &
rc)) return
! if location is given, check if the location attribute of the variable match with the
! input value
if (present(location)) then
call ESMF_UGridInqVarLoc(ncid, VarId, varname, locflag, localrc)
if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
ESMF_CONTEXT, rcToReturn=rc)) return
if (((location .eq. 'node') .and. (locflag /= 1)) .or. &
((location .eq. 'face') .and. (locflag /= 2))) then
errmsg = "- variable "//varname//" is not defined on location "//location
call ESMF_LogSetError(rcToCheck=ESMF_FAILURE, &
msg=errmsg, ESMF_CONTEXT, rcToReturn=rc)
return
endif
endif
! check if the variable dimension matches with the allocated array
errmsg = "Variable "//varname//" in "//trim(filename)
ncStatus = nf90_inquire_variable(ncid, VarId, ndims=ndims, dimids=dimids)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE, errmsg, &
rc)) return
! get the first dimension length of the variable
ncStatus = nf90_inquire_dimension (ncid, dimids(1), len=dimsize)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE,&
errmsg,&
rc)) return
if (present(startind)) then
if (size(varbuffer) < count) then
call ESMF_LogSetError(rcToCheck=ESMF_FAILURE, &
msg="- the varbuffer is too small", &
ESMF_CONTEXT, rcToReturn=rc)
return
endif
else if (dimsize /= size(varbuffer)) then
call ESMF_LogSetError(rcToCheck=ESMF_FAILURE, &
msg="- the variable array dimension does not match with dimension length", &
ESMF_CONTEXT, rcToReturn=rc)
return
endif
allocate(starts(ndims), counts(ndims))
starts(:)=1
counts(:)=1
counts(1)=dimsize
if (present(startind)) then
starts(1)=startind
counts(1)=count
endif
ncStatus = nf90_get_var (ncid, VarId, varbuffer, start=starts, count=counts)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE, errmsg, &
rc)) return
deallocate(starts, counts)
! get missisng value
if (present(missingvalue)) then
ncStatus = nf90_get_att(ncid, VarId, "_FillValue", missingvalue)
if (ncStatus /= nf90_noerror) then
ncStatus = nf90_get_att(ncid, varid, "missing_value", missingvalue)
errmsg = "missing value attribute does not exist for "//trim(varname)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE,&
errmsg, &
rc)) return
end if
end if
ncStatus = nf90_close (ncid=ncid)
if (CDFCheckError (ncStatus, &
ESMF_METHOD, &
ESMF_SRCLINE, trim(filename), &
rc)) return
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_UGridGetVarByName