mirror of
https://github.com/Unidata/netcdf-c.git
synced 2024-12-21 08:39:46 +08:00
f3e711e2b8
re: https://github.com/Unidata/netcdf-c/issues/2177 re: https://github.com/Unidata/netcdf-c/pull/2178 Provide get/set functions to store global data alignment information and apply it when a file is created. The api is as follows: ```` int nc_set_alignment(int threshold, int alignment); int nc_get_alignment(int* thresholdp, int* alignmentp); ```` If defined, then for every file created opened after the call to nc_set_alignment, for every new variable added to the file, the most recently set threshold and alignment values will be applied to that variable. The nc_get_alignment function return the last values set by nc_set_alignment. If nc_set_alignment has not been called, then it returns the value 0 for both threshold and alignment. The alignment parameters are stored in the NCglobalstate object (see below) for use as needed. Repeated calls to nc_set_alignment will overwrite any existing values in NCglobalstate. The alignment parameters are applied in libhdf5/hdf5create.c and libhdf5/hdf5open.c The set/get alignment functions are defined in libsrc4/nc4internal.c. A test program was added as nc_test4/tst_alignment.c. ## Misc. Changes Unrelated to Alignment * The NCRCglobalstate type was renamed to NCglobalstate to indicate that it represented more general global state than just .rc data. It was also moved to nc4internal.h. This led to a large number of small changes: mostly renaming. The global state management functions were moved to nc4internal.c. * The global chunk cache variables have been moved into NCglobalstate. As warranted, other global state will be moved as well. * Some misc. problems with the nczarr performance tests were corrected.
97 lines
2.8 KiB
C
97 lines
2.8 KiB
C
/* This is part of the netCDF package.
|
|
Copyright 2018 University Corporation for Atmospheric Research/Unidata
|
|
See COPYRIGHT file for conditions of use.
|
|
|
|
Test HDF5 alignment
|
|
Dennis Heimbigner
|
|
*/
|
|
|
|
#include <nc_tests.h>
|
|
#include "err_macros.h"
|
|
|
|
#include <hdf5.h>
|
|
#include <H5DSpublic.h>
|
|
|
|
#undef DEBUG
|
|
|
|
#define THRESHOLD 512
|
|
#define ALIGNMENT 4096
|
|
#define CHUNKSIZE 513
|
|
|
|
int
|
|
main(int argc, char **argv)
|
|
{
|
|
int i,ncid, varid, dimids[1];
|
|
size_t chunks[1];
|
|
unsigned char data[CHUNKSIZE];
|
|
hid_t fileid, grpid, datasetid;
|
|
hid_t dxpl_id = H5P_DEFAULT; /*data transfer property list */
|
|
unsigned int filter_mask = 0;
|
|
hsize_t hoffset[1];
|
|
haddr_t addr;
|
|
hsize_t size ;
|
|
hid_t fspace;
|
|
|
|
H5Eset_auto2(H5E_DEFAULT,(H5E_auto2_t)H5Eprint1,stderr);
|
|
|
|
printf("\n*** Testing HDF5 alignment.\n");
|
|
|
|
printf("chunksize=%d threshold=%d alignment=%d\n",CHUNKSIZE,THRESHOLD,ALIGNMENT);
|
|
|
|
if(nc_set_alignment(THRESHOLD,ALIGNMENT)) ERR;
|
|
if (nc_create("tst_alignment.nc", NC_NETCDF4, &ncid)) ERR;
|
|
if (nc_def_dim(ncid, "d0", CHUNKSIZE, &dimids[0])) ERR;
|
|
if (nc_def_var(ncid, "var", NC_UBYTE, 1, dimids, &varid)) ERR;
|
|
chunks[0] = CHUNKSIZE;
|
|
if (nc_def_var_chunking(ncid, varid, NC_CHUNKED, chunks)) ERR;
|
|
if (nc_enddef(ncid)) ERR;
|
|
|
|
for(i=0;i<CHUNKSIZE;i++) data[i] = (i)%2;
|
|
if(nc_put_var(ncid,varid,data)) ERR;
|
|
|
|
if (nc_close(ncid)) ERR;
|
|
|
|
/* Use HDF5 API to verify */
|
|
|
|
if ((fileid = H5Fopen("tst_alignment.nc", H5F_ACC_RDONLY, H5P_DEFAULT)) < 0) ERR;
|
|
if ((grpid = H5Gopen1(fileid, "/")) < 0) ERR;
|
|
if ((datasetid = H5Dopen1(grpid, "var")) < 0) ERR;
|
|
|
|
/* Test the offset */
|
|
if((fspace = H5Dget_space(datasetid)) < 0) ERR;
|
|
if(H5Dget_chunk_info(datasetid, fspace, 0, hoffset, &filter_mask, &addr, &size) < 0) ERR;
|
|
#ifdef DEBUG
|
|
fprintf(stderr,"H5Dget_chunk_info: offset=%lu addr=%lu size=%lu\n",(unsigned long)hoffset[0], (unsigned long)addr, (unsigned long)size);
|
|
fprintf(stderr,"\t%% alignment: offset=%lu, addr=%lu\n",
|
|
(((unsigned long)hoffset[0]) % ALIGNMENT),
|
|
(((unsigned long)addr) % ALIGNMENT)
|
|
);
|
|
#endif
|
|
|
|
printf("H5Dget_chunk_info: addr=%lu (addr %% alignment)=%lu\n", (unsigned long)addr, (((unsigned long)hoffset[0]) % ALIGNMENT));
|
|
if((addr % ALIGNMENT) != 0) ERR;
|
|
|
|
/* Test chunk content */
|
|
memset(data,0,sizeof(data));
|
|
hoffset[0] = 0;
|
|
|
|
if(H5Dread_chunk(datasetid, dxpl_id, hoffset, &filter_mask, data) < 0) ERR;
|
|
|
|
#ifdef DEBUG
|
|
fprintf(stderr,"H5Dread_chunk: offset=%lu offset %% alignment=%lu\n",(unsigned long)hoffset[0], (((unsigned long)hoffset) % ALIGNMENT));
|
|
#endif
|
|
for(i=0;i<CHUNKSIZE;i++) {
|
|
if(data[i] != (i)%2) {
|
|
fprintf(stderr,"data[%d] mismatch\n",i);
|
|
ERR;
|
|
}
|
|
}
|
|
|
|
if (H5Dclose(datasetid) < 0) ERR;
|
|
if (H5Gclose(grpid) < 0) ERR;
|
|
if (H5Fclose(fileid) < 0) ERR;
|
|
|
|
SUMMARIZE_ERR;
|
|
FINAL_RESULTS;
|
|
}
|