summaryrefslogtreecommitdiffstats
path: root/testpar
diff options
context:
space:
mode:
authorAlbert Cheng <acheng@hdfgroup.org>2004-01-03 02:54:50 (GMT)
committerAlbert Cheng <acheng@hdfgroup.org>2004-01-03 02:54:50 (GMT)
commitf65948d933ac5bc70ec5bd3ab43387b310722d97 (patch)
treee1ab74e1a87cd1b697d8cb2fbca01dfdf4f9fa2c /testpar
parenta59c045de10854aa53d69e2f27b52254138ee7c1 (diff)
downloadhdf5-f65948d933ac5bc70ec5bd3ab43387b310722d97.zip
hdf5-f65948d933ac5bc70ec5bd3ab43387b310722d97.tar.gz
hdf5-f65948d933ac5bc70ec5bd3ab43387b310722d97.tar.bz2
[svn-r8013] Description:
Added a test of fill value before any data is written to a dataset. Rename short_dataset() as dataset_fillvalue() as it reflects better the tests. Also removed the option of -S since the fill value test will be tested always. Platforms tested: "h5committested" Misc. update:
Diffstat (limited to 'testpar')
-rw-r--r--testpar/t_mdset.c101
-rw-r--r--testpar/testphdf5.c17
-rw-r--r--testpar/testphdf5.h2
3 files changed, 76 insertions, 44 deletions
diff --git a/testpar/t_mdset.c b/testpar/t_mdset.c
index 11f81af..9b4bfd1 100644
--- a/testpar/t_mdset.c
+++ b/testpar/t_mdset.c
@@ -306,11 +306,11 @@ void big_dataset(const char *filename)
VRFY((ret >= 0), "H5Pclose succeeded");
}
-/* Example of using PHDF5 to read "short" datasets. These datasets don't have
- * actual data written to the entire raw data area and rely on the "fill with
- * zeros" code in the VFL driver read routine to work correctly.
+/* Example of using PHDF5 to read a partial written dataset. The dataset does
+ * not have actual data written to the entire raw data area and relies on the
+ * default fill value of zeros to work correctly.
*/
-void short_dataset(const char *filename)
+void dataset_fillvalue(const char *filename)
{
int mpi_size, mpi_rank; /* MPI info */
hbool_t use_gpfs = FALSE; /* Don't use GPFS stuff for this test */
@@ -322,11 +322,12 @@ void short_dataset(const char *filename)
memspace, /* Memory dataspace ID */
filespace; /* Dataset's dataspace ID */
char dname[]="dataset"; /* Name of dataset */
- hsize_t dset_size[4] = {0, 6, 7, 8};
+ hsize_t dset_dims[4] = {0, 6, 7, 8};
hssize_t req_start[4] = {0, 0, 0, 0};
hsize_t req_count[4] = {1, 6, 7, 8};
+ hsize_t dset_size; /* Dataset size */
int *rdata, *wdata; /* Buffers for data to read and write */
- int *tdata, *tdata2; /* Temporary pointer into buffer */
+ int *twdata, *trdata; /* Temporary pointer into buffer */
int acc, i, j, k, l; /* Local index variables */
herr_t ret; /* Generic return value */
@@ -335,22 +336,17 @@ void short_dataset(const char *filename)
VRFY((mpi_size <= SIZE), "mpi_size <= SIZE");
+ /* Set the dataset dimension to be one row more than number of processes */
+ /* and calculate the actual dataset size. */
+ dset_dims[0]=mpi_size+1;
+ dset_size=dset_dims[0]*dset_dims[1]*dset_dims[2]*dset_dims[3];
+
/* Allocate space for the buffers */
- dset_size[0]=mpi_size+1;
- rdata=HDmalloc((size_t)(dset_size[0]*dset_size[1]*dset_size[2]*dset_size[3]*sizeof(int)));
+ rdata=HDmalloc((size_t)(dset_size*sizeof(int)));
VRFY((rdata != NULL), "HDcalloc succeeded for read buffer");
- wdata=HDmalloc((size_t)(dset_size[0]*dset_size[1]*dset_size[2]*dset_size[3]*sizeof(int)));
+ wdata=HDmalloc((size_t)(dset_size*sizeof(int)));
VRFY((wdata != NULL), "HDmalloc succeeded for write buffer");
- /* Initialize write buffer */
- HDmemset(rdata,2,(size_t)(dset_size[0]*dset_size[1]*dset_size[2]*dset_size[3]*sizeof(int)));
- tdata=wdata;
- for (i=0, acc=0; i<(int)dset_size[0]; i++)
- for (j=0; j<(int)dset_size[1]; j++)
- for (k=0; k<(int)dset_size[2]; k++)
- for (l=0; l<(int)dset_size[3]; l++)
- *tdata++ = acc++;
-
fapl = create_faccess_plist(MPI_COMM_WORLD, MPI_INFO_NULL, facc_type, use_gpfs);
VRFY((fapl >= 0), "create_faccess_plist succeeded");
@@ -360,15 +356,47 @@ void short_dataset(const char *filename)
iof = H5Fcreate (filename, H5F_ACC_TRUNC, H5P_DEFAULT, fapl);
VRFY((iof >= 0), "H5Fcreate succeeded");
- filespace = H5Screate_simple(4, dset_size, NULL);
+ filespace = H5Screate_simple(4, dset_dims, NULL);
VRFY((filespace >= 0), "File H5Screate_simple succeeded");
dataset = H5Dcreate(iof, dname, H5T_NATIVE_INT, filespace, H5P_DEFAULT);
VRFY((dataset >= 0), "H5Dcreate succeeded");
- memspace = H5Screate_simple(4, dset_size, NULL);
+ memspace = H5Screate_simple(4, dset_dims, NULL);
VRFY((memspace >= 0), "Memory H5Screate_simple succeeded");
+ /*
+ * Read dataset before any data is written.
+ */
+ /* set entire read buffer with the constant 2 */
+ HDmemset(rdata,2,(size_t)(dset_size*sizeof(int)));
+ /* Independently read the entire dataset back */
+ ret=H5Dread(dataset, H5T_NATIVE_INT, H5S_ALL, H5S_ALL, H5P_DEFAULT, rdata);
+ VRFY((ret >= 0), "H5Dread succeeded");
+
+ /* Verify all data read are the fill value 0 */
+ trdata=rdata;
+ err_num=0;
+ for (i=0; i<(int)dset_dims[0]; i++)
+ for (j=0; j<(int)dset_dims[1]; j++)
+ for (k=0; k<(int)dset_dims[2]; k++)
+ for (l=0; l<(int)dset_dims[3]; l++, twdata++, trdata++)
+ if( *trdata != 0)
+ if(err_num++ < MAX_ERR_REPORT || verbose)
+ printf("Dataset Verify failed at [%d][%d][%d][%d]: expect 0, got %d\n", i,j,k,l, *trdata);
+ if(err_num > MAX_ERR_REPORT && !verbose)
+ printf("[more errors ...]\n");
+ if(err_num){
+ printf("%d errors found in check_value\n", err_num);
+ nerrors++;
+ }
+
+ /* Barrier to ensure all processes have completed the above test. */
+ MPI_Barrier(MPI_COMM_WORLD);
+
+ /*
+ * Each process writes 1 row of data. Thus last row is not written.
+ */
/* Create hyperslabs in memory and file dataspaces */
req_start[0]=mpi_rank;
ret=H5Sselect_hyperslab(filespace, H5S_SELECT_SET, req_start, NULL, req_count, NULL);
@@ -383,6 +411,14 @@ void short_dataset(const char *filename)
ret=H5Pset_dxpl_mpio(dxpl, H5FD_MPIO_COLLECTIVE);
VRFY((ret >= 0), "H5Pset_dxpl_mpio succeeded");
+ /* Fill write buffer with some values */
+ twdata=wdata;
+ for (i=0, acc=0; i<(int)dset_dims[0]; i++)
+ for (j=0; j<(int)dset_dims[1]; j++)
+ for (k=0; k<(int)dset_dims[2]; k++)
+ for (l=0; l<(int)dset_dims[3]; l++)
+ *twdata++ = acc++;
+
/* Collectively write a hyperslab of data to the dataset */
ret=H5Dwrite(dataset, H5T_NATIVE_INT, memspace, filespace, dxpl, wdata);
VRFY((ret >= 0), "H5Dwrite succeeded");
@@ -390,27 +426,32 @@ void short_dataset(const char *filename)
/* Barrier here, to allow MPI-posix I/O to sync */
MPI_Barrier(MPI_COMM_WORLD);
+ /*
+ * Read dataset after partial write.
+ */
+ /* set entire read buffer with the constant 2 */
+ HDmemset(rdata,2,(size_t)(dset_size*sizeof(int)));
/* Independently read the entire dataset back */
ret=H5Dread(dataset, H5T_NATIVE_INT, H5S_ALL, H5S_ALL, H5P_DEFAULT, rdata);
VRFY((ret >= 0), "H5Dread succeeded");
/* Verify correct data read */
- tdata=wdata;
- tdata2=rdata;
+ twdata=wdata;
+ trdata=rdata;
err_num=0;
- for (i=0; i<(int)dset_size[0]; i++)
- for (j=0; j<(int)dset_size[1]; j++)
- for (k=0; k<(int)dset_size[2]; k++)
- for (l=0; l<(int)dset_size[3]; l++, tdata++, tdata2++)
+ for (i=0; i<(int)dset_dims[0]; i++)
+ for (j=0; j<(int)dset_dims[1]; j++)
+ for (k=0; k<(int)dset_dims[2]; k++)
+ for (l=0; l<(int)dset_dims[3]; l++, twdata++, trdata++)
if(i<mpi_size) {
- if( *tdata != *tdata2 )
+ if( *twdata != *trdata )
if(err_num++ < MAX_ERR_REPORT || verbose)
- printf("Dataset Verify failed at [%d][%d][%d][%d]: expect %d, got %d\n", i,j,k,l, *tdata, *tdata2);
+ printf("Dataset Verify failed at [%d][%d][%d][%d]: expect %d, got %d\n", i,j,k,l, *twdata, *trdata);
} /* end if */
else {
- if( *tdata2 != 0)
+ if( *trdata != 0)
if(err_num++ < MAX_ERR_REPORT || verbose)
- printf("Dataset Verify failed at [%d][%d][%d][%d]: expect 0, got %d\n", i,j,k,l, *tdata2);
+ printf("Dataset Verify failed at [%d][%d][%d][%d]: expect 0, got %d\n", i,j,k,l, *trdata);
} /* end else */
if(err_num > MAX_ERR_REPORT && !verbose)
printf("[more errors ...]\n");
diff --git a/testpar/testphdf5.c b/testpar/testphdf5.c
index 8394db4..4229876 100644
--- a/testpar/testphdf5.c
+++ b/testpar/testphdf5.c
@@ -45,7 +45,6 @@ int dowrite=1; /* write test */
int docompact=1; /* compact dataset test */
int doindependent=1; /* independent test */
unsigned dobig=0; /* "big" dataset tests */
-unsigned doshort=1; /* "short" dataset tests */
/* FILENAME and filenames must have the same number of names */
const char *FILENAME[10]={
@@ -57,7 +56,7 @@ const char *FILENAME[10]={
"ParaCompact",
"ParaIndividual",
"ParaBig",
- "ParaShort",
+ "ParaFill",
NULL};
char filenames[10][PATH_MAX];
hid_t fapl; /* file access property list */
@@ -130,7 +129,6 @@ usage(void)
printf("\t-o\t\tno compact dataset test\n");
printf("\t-i\t\tno independent read test\n");
printf("\t-b\t\trun big dataset test\n");
- printf("\t-S\t\tno short dataset test\n");
printf("\t-v\t\tverbose on\n");
printf("\t-f <prefix>\tfilename prefix\n");
printf("\t-s\t\tuse Split-file together with MPIO\n");
@@ -228,8 +226,6 @@ parse_options(int argc, char **argv)
break;
case 'h': /* print help message--return with nerrors set */
return(1);
- case 'S': doshort = 0;
- break;
default: nerrors++;
return(1);
}
@@ -484,15 +480,10 @@ int main(int argc, char **argv)
MPI_BANNER("big dataset test skipped");
}
- if (doshort) {
- MPI_BANNER("short dataset test...");
- short_dataset(filenames[8]);
- }
- else {
- MPI_BANNER("short dataset test skipped");
- }
+ MPI_BANNER("dataset fill value test...");
+ dataset_fillvalue(filenames[8]);
- if (!(dowrite || doread || ndatasets || ngroups || docompact || doindependent || dobig || doshort)){
+ if (!(dowrite || doread || ndatasets || ngroups || docompact || doindependent || dobig )){
usage();
nerrors++;
}
diff --git a/testpar/testphdf5.h b/testpar/testphdf5.h
index e82644a..c6069d5 100644
--- a/testpar/testphdf5.h
+++ b/testpar/testphdf5.h
@@ -145,7 +145,7 @@ void extend_readInd(char *filename);
void extend_readAll(char *filename);
void compact_dataset(char *filename);
void big_dataset(const char *filename);
-void short_dataset(const char *filename);
+void dataset_fillvalue(const char *filename);
int dataset_vrfy(hssize_t start[], hsize_t count[], hsize_t stride[],
hsize_t block[], DATATYPE *dataset, DATATYPE *original);