Next: , Previous: NF90_INQ_VAR_FILL, Up: Variables


6.6 Get Information about a Variable from Its ID: NF90_INQUIRE_VARIABLE

NF90_INQUIRE_VARIABLE returns information about a netCDF variable given its ID. Information about a variable includes its name, type, number of dimensions, a list of dimension IDs describing the shape of the variable, and the number of variable attributes that have been assigned to the variable.

All parameters after nAtts are optional, and only supported if netCDF was built with netCDF-4 features enabled, and if the variable is in a netCDF-4/HDF5 file.

Usage

       function nf90_inquire_variable(ncid, varid, name, xtype, ndims, dimids, nAtts, &
            contiguous, chunksizes, deflate_level, shuffle, fletcher32, endianness)
         integer, intent(in) :: ncid, varid
         character (len = *), optional, intent(out) :: name
         integer, optional, intent(out) :: xtype, ndims
         integer, dimension(:), optional, intent(out) :: dimids
         integer, optional, intent(out) :: nAtts
         logical, optional, intent(out) :: contiguous
         integer, optional, dimension(:), intent(out) :: chunksizes
         integer, optional, intent(out) :: deflate_level
         logical, optional, intent(out) :: shuffle, fletcher32
         integer, optional, intent(out) :: endianness
         integer :: nf90_inquire_variable
ncid
NetCDF ID, from a previous call to NF90_OPEN or NF90_CREATE.
varid
Variable ID.
name
Returned variable name. The caller must allocate space for the returned name. The maximum possible length, in characters, of a variable name is given by the predefined constant NF90_MAX_NAME.
xtype
Returned variable type, one of the set of predefined netCDF external data types. The type of this parameter, NF90_TYPE, is defined in the netCDF header file. The valid netCDF external data types are NF90_BYTE, NF90_CHAR, NF90_SHORT, NF90_INT, NF90_FLOAT, AND NF90_DOUBLE.
ndims
Returned number of dimensions the variable was defined as using. For example, 2 indicates a matrix, 1 indicates a vector, and 0 means the variable is a scalar with no dimensions.
dimids
Returned vector of *ndimsp dimension IDs corresponding to the variable dimensions. The caller must allocate enough space for a vector of at least *ndimsp integers to be returned. The maximum possible number of dimensions for a variable is given by the predefined constant NF90_MAX_VAR_DIMS.
natts
Returned number of variable attributes assigned to this variable.
contiguous
On return, set to NF90_CONTIGUOUS if this variable uses contiguous storage, NF90_CHUNKED if it uses chunked storage.
chunksizes
An array of chunk sizes. The array must have the one element for each dimension in the variable.
shuffle
True if the shuffle filter is turned on for this variable.
deflate_level
The deflate_level from 0 to 9. A value of zero indicates no deflation is in use.
fletcher32
Set to true if the fletcher32 checksum filter is turned on for this variable.
endianness
Will be set to NF90_ENDIAN_LITTLE if this variable is stored in little-endian format, NF90_ENDIAN_BIG if it is stored in big-endian format, and NF90_ENDIAN_NATIVE if the endianness is not set, and the variable is not created yet.

These functions return the value NF90_NOERR if no errors occurred. Otherwise, the returned status indicates an error. Possible causes of errors include:

Example

Here is an example using NF90_INQ_VAR to find out about a variable named rh in an existing netCDF dataset named foo.nc:

         use netcdf
         implicit none
         integer                            :: status, ncid, &
                                               RhVarId       &
                                               numDims, numAtts
      integer, dimension(nf90_max_var_dims) :: rhDimIds
      ...
      status = nf90_open("foo.nc", nf90_NoWrite, ncid)
      if(status /= nf90_NoErr) call handle_error(status)
      ...
      status = nf90_inq_varid(ncid, "rh", RhVarId)
      if(status /= nf90_NoErr) call handle_err(status)
      status = nf90_inquire_variable(ncid, RhVarId, ndims = numDims, natts = numAtts)
      if(status /= nf90_NoErr) call handle_err(status)
      status = nf90_inquire_variable(ncid, RhVarId, dimids = rhDimIds(:numDims))
      if(status /= nf90_NoErr) call handle_err(status)