Next: NF_INQ_VAR_CHUNKING, Previous: NF_DEF_VAR, Up: Variables
NF_DEF_VAR_CHUNKINGThe function NF_DEF_VAR_CHUNKING sets the storage parameters for a variable in a netCDF-4 file. It can set the chunk sizes to get chunked storage, or it can set the contiguous flag to get contiguous storage.
The total size of a chunk must be less than 4 GiB. That is, the product of all chunksizes and the size of the data (or the size of nc_vlen_t for VLEN types) must be less than 4 GiB.
This function may only be called after the variable is defined, but before nf_enddef is called. Once the chunking parameters are set for a variable, they cannot be changed.
NF_DEF_VAR_CHUNKING(INTEGER NCID, INTEGER VARID, INTEGER STORAGE, INTEGER CHUNKSIZES)
ncidvaridstorageIf NF_CHUNKED, then chunked storage is used for this variable.
Chunk sizes may be specified with the chunksizes parameter.
Default sizes will be used if chunking is required and this function
is not called.
chunksizesNF_DEF_VAR_CHUNKING returns the value NF_NOERR if no errors occurred. Otherwise, the returned status indicates an error.
Possible return codes include:
NF_NOERRNF_BADIDNF_ENOTNC4NF_ENOTVARNF_ELATEDEFNF_ENOTINDEFINENF_ESTRICTNC3In this example from nf_test/ftst_vars.F, a file is created, two dimensions and a variable are defined, and the chunksizes of the data are set to the size of the data (that is, data will be written in one chunk).
C Create the netCDF file.
retval = nf_create(FILE_NAME, NF_NETCDF4, ncid)
if (retval .ne. nf_noerr) call handle_err(retval)
C Define the dimensions.
retval = nf_def_dim(ncid, "x", NX, x_dimid)
if (retval .ne. nf_noerr) call handle_err(retval)
retval = nf_def_dim(ncid, "y", NY, y_dimid)
if (retval .ne. nf_noerr) call handle_err(retval)
C Define the variable.
dimids(1) = y_dimid
dimids(2) = x_dimid
retval = NF_DEF_VAR(ncid, "data", NF_INT, NDIMS, dimids, varid)
if (retval .ne. nf_noerr) call handle_err(retval)
C Turn on chunking.
chunks(1) = NY
chunks(2) = NX
retval = NF_DEF_VAR_chunking(ncid, varid, NF_CHUNKED, chunks)
if (retval .ne. nf_noerr) call handle_err(retval)