Function nf90_get_att gets the value(s) of a netCDF attribute, given its variable ID and name.
function nf90_get_att(ncid, varid, name, values)
integer, intent( in) :: ncid, varid
character(len = *), intent( in) :: name
any valid type, scalar or array of rank 1, &
intent(out) :: values
integer :: nf90_get_att
ncidvaridnamevaluesNF90_GET_ATT_ type returns the value NF90_NOERR if no errors occurred. Otherwise, the returned status indicates an error. Possible causes of errors include:
Here is an example using NF90_GET_ATT to determine the values of an attribute named valid_range for a netCDF variable named rh and a global attribute named title in an existing netCDF dataset named foo.nc. In this example, it is assumed that we don't know how many values will be returned, so we first inquire about the length of the attributes to make sure we have enough space to store them:
use netcdf
implicit none
integer :: ncid, status
integer :: RHVarID ! Variable ID
integer :: validRangeLength, titleLength ! Attribute lengths
real, dimension(:), allocatable, &
:: validRange
character (len = 80) :: title
...
status = nf90_open("foo.nc", nf90_nowrite, ncid)
if (status /= nf90_noerr) call handle_err(status)
...
! Find the lengths of the attributes
status = nf90_inq_varid(ncid, "rh", RHVarID)
if (status /= nf90_noerr) call handle_err(status)
status = nf90_inquire_attribute(ncid, RHVarID, "valid_range", &
len = validRangeLength)
if (status /= nf90_noerr) call handle_err(status)
status = nf90_inquire_attribute(ncid, nf90_global, "title", len = titleLength)
if (status /= nf90_noerr) call handle_err(status)
...
!Allocate space to hold attribute values, check string lengths
allocate(validRange(validRangeLength), stat = status)
if(status /= 0 .or. len(title) < titleLength)
print *, "Not enough space to put attribute values."
exit
end if
! Read the attributes.
status = nf90_get_att(ncid, RHVarID, "valid_range", validRange)
if (status /= nf90_noerr) call handle_err(status)
status = nf90_get_att(ncid, nf90_global, "title", title)
if (status /= nf90_noerr) call handle_err(status)