summaryrefslogtreecommitdiffstats
path: root/testpar
diff options
context:
space:
mode:
authorQuincey Koziol <koziol@hdfgroup.org>2003-09-16 17:33:00 (GMT)
committerQuincey Koziol <koziol@hdfgroup.org>2003-09-16 17:33:00 (GMT)
commitd7bde16f45fac765f45172d88a1a9cd44a1f95fa (patch)
treebe73c954ea0220b752cfbb49597e35954aa47cd9 /testpar
parentbf1c2f0e8bf9788c2e47a1b1ac963cc321afab0e (diff)
downloadhdf5-d7bde16f45fac765f45172d88a1a9cd44a1f95fa.zip
hdf5-d7bde16f45fac765f45172d88a1a9cd44a1f95fa.tar.gz
hdf5-d7bde16f45fac765f45172d88a1a9cd44a1f95fa.tar.bz2
[svn-r7480] Purpose:
Bug fix Description: The MPI_File_set_size() routine on ASCI Red is not able to extend files so that they are larger than 2GB. Solution: Add an extra macro which controls whether MPI_File_set_size() can handle >2GB offsets or if our "older" way of reading a byte, then writing a byte at the appropriate offset should be used. Platforms tested: FreeBSD 4.9 (sleipnir) h5committest
Diffstat (limited to 'testpar')
-rw-r--r--testpar/t_mdset.c118
-rw-r--r--testpar/testphdf5.c39
-rw-r--r--testpar/testphdf5.h2
3 files changed, 156 insertions, 3 deletions
diff --git a/testpar/t_mdset.c b/testpar/t_mdset.c
index b8bf4b6..3d40307 100644
--- a/testpar/t_mdset.c
+++ b/testpar/t_mdset.c
@@ -211,6 +211,124 @@ void compact_dataset(char *filename)
H5Fclose(iof);
}
+/* Example of using PHDF5 to create "large" datasets. (>2GB, >4GB, >8GB)
+ * Actual data is _not_ written to these datasets. Dataspaces are exact
+ * sizes (2GB, 4GB, etc.), but the metadata for the file pushes the file over
+ * the boundary of interest.
+ */
+void big_dataset(const char *filename)
+{
+ int mpi_size, mpi_rank; /* MPI info */
+ hbool_t use_gpfs = FALSE; /* Don't use GPFS stuff for this test */
+ hid_t iof, /* File ID */
+ fapl, /* File access property list ID */
+ dataset, /* Dataset ID */
+ filespace; /* Dataset's dataspace ID */
+ hsize_t file_dims [4]; /* Dimensions of dataspace */
+ char dname[]="dataset"; /* Name of dataset */
+ MPI_Offset file_size; /* Size of file on disk */
+ herr_t ret; /* Generic return value */
+
+ MPI_Comm_rank (MPI_COMM_WORLD, &mpi_rank);
+ MPI_Comm_size (MPI_COMM_WORLD, &mpi_size);
+
+ VRFY((mpi_size <= SIZE), "mpi_size <= SIZE");
+
+ fapl = create_faccess_plist(MPI_COMM_WORLD, MPI_INFO_NULL, facc_type, use_gpfs);
+ VRFY((fapl >= 0), "create_faccess_plist succeeded");
+
+ /*
+ * Create >2GB HDF5 file
+ */
+ iof = H5Fcreate (filename, H5F_ACC_TRUNC, H5P_DEFAULT, fapl);
+ VRFY((iof >= 0), "H5Fcreate succeeded");
+
+ /* Define dataspace for 2GB dataspace */
+ file_dims[0]= 2;
+ file_dims[1]= 1024;
+ file_dims[2]= 1024;
+ file_dims[3]= 1024;
+ filespace = H5Screate_simple (4, file_dims, NULL);
+ VRFY((filespace >= 0), "H5Screate_simple succeeded");
+
+ dataset = H5Dcreate (iof, dname, H5T_NATIVE_UCHAR, filespace, H5P_DEFAULT);
+ VRFY((dataset >= 0), "H5Dcreate succeeded");
+
+ /* Close all file objects */
+ ret=H5Dclose (dataset);
+ VRFY((ret >= 0), "H5Dclose succeeded");
+ ret=H5Sclose (filespace);
+ VRFY((ret >= 0), "H5Sclose succeeded");
+ ret=H5Fclose (iof);
+ VRFY((ret >= 0), "H5Fclose succeeded");
+
+ /* Check that file of the correct size was created */
+ file_size=h5_mpi_get_file_size(filename, MPI_COMM_WORLD, MPI_INFO_NULL);
+ VRFY((file_size == 2147485696ULL), "File is correct size");
+
+ /*
+ * Create >4GB HDF5 file
+ */
+ iof = H5Fcreate (filename, H5F_ACC_TRUNC, H5P_DEFAULT, fapl);
+ VRFY((iof >= 0), "H5Fcreate succeeded");
+
+ /* Define dataspace for 4GB dataspace */
+ file_dims[0]= 4;
+ file_dims[1]= 1024;
+ file_dims[2]= 1024;
+ file_dims[3]= 1024;
+ filespace = H5Screate_simple (4, file_dims, NULL);
+ VRFY((filespace >= 0), "H5Screate_simple succeeded");
+
+ dataset = H5Dcreate (iof, dname, H5T_NATIVE_UCHAR, filespace, H5P_DEFAULT);
+ VRFY((dataset >= 0), "H5Dcreate succeeded");
+
+ /* Close all file objects */
+ ret=H5Dclose (dataset);
+ VRFY((ret >= 0), "H5Dclose succeeded");
+ ret=H5Sclose (filespace);
+ VRFY((ret >= 0), "H5Sclose succeeded");
+ ret=H5Fclose (iof);
+ VRFY((ret >= 0), "H5Fclose succeeded");
+
+ /* Check that file of the correct size was created */
+ file_size=h5_mpi_get_file_size(filename, MPI_COMM_WORLD, MPI_INFO_NULL);
+ VRFY((file_size == 4294969344ULL), "File is correct size");
+
+ /*
+ * Create >8GB HDF5 file
+ */
+ iof = H5Fcreate (filename, H5F_ACC_TRUNC, H5P_DEFAULT, fapl);
+ VRFY((iof >= 0), "H5Fcreate succeeded");
+
+ /* Define dataspace for 8GB dataspace */
+ file_dims[0]= 8;
+ file_dims[1]= 1024;
+ file_dims[2]= 1024;
+ file_dims[3]= 1024;
+ filespace = H5Screate_simple (4, file_dims, NULL);
+ VRFY((filespace >= 0), "H5Screate_simple succeeded");
+
+ dataset = H5Dcreate (iof, dname, H5T_NATIVE_UCHAR, filespace, H5P_DEFAULT);
+ VRFY((dataset >= 0), "H5Dcreate succeeded");
+
+ /* Close all file objects */
+ ret=H5Dclose (dataset);
+ VRFY((ret >= 0), "H5Dclose succeeded");
+ ret=H5Sclose (filespace);
+ VRFY((ret >= 0), "H5Sclose succeeded");
+ ret=H5Fclose (iof);
+ VRFY((ret >= 0), "H5Fclose succeeded");
+
+ /* Check that file of the correct size was created */
+ file_size=h5_mpi_get_file_size(filename, MPI_COMM_WORLD, MPI_INFO_NULL);
+ VRFY((file_size == 8589936640ULL), "File is correct size");
+
+ /* Close fapl */
+ ret=H5Pclose (fapl);
+ VRFY((ret >= 0), "H5Pclose succeeded");
+}
+
/* Write multiple groups with a chunked dataset in each group collectively.
* These groups and datasets are for testing independent read later.
*/
diff --git a/testpar/testphdf5.c b/testpar/testphdf5.c
index 6a4a0b2..ff44424 100644
--- a/testpar/testphdf5.c
+++ b/testpar/testphdf5.c
@@ -44,9 +44,10 @@ int doread=1; /* read test */
int dowrite=1; /* write test */
int docompact=1; /* compact dataset test */
int doindependent=1; /* independent test */
+unsigned dobig=1; /* "big" dataset tests */
/* FILENAME and filenames must have the same number of names */
-const char *FILENAME[8]={
+const char *FILENAME[9]={
"ParaEg1",
"ParaEg2",
"ParaEg3",
@@ -54,8 +55,9 @@ const char *FILENAME[8]={
"ParaMgroup",
"ParaCompact",
"ParaIndividual",
+ "ParaBig",
NULL};
-char filenames[8][PATH_MAX];
+char filenames[9][PATH_MAX];
hid_t fapl; /* file access property list */
#ifdef USE_PAUSE
@@ -177,6 +179,8 @@ parse_options(int argc, char **argv)
break;
case 'i': doindependent = 0;
break;
+ case 'b': dobig = 0;
+ break;
case 'v': verbose = 1;
break;
case 'f': if (--argc < 1) {
@@ -326,6 +330,27 @@ create_faccess_plist(MPI_Comm comm, MPI_Info info, int l_facc_type,
return (ret_pl);
}
+/*
+ * Check the size of a file using MPI routines
+ */
+MPI_Offset
+h5_mpi_get_file_size(const char *filename, MPI_Comm comm, MPI_Info info)
+{
+ MPI_File fh; /* MPI file handle */
+ MPI_Offset size=0; /* File size to return */
+
+ if (MPI_SUCCESS != MPI_File_open(comm, (char*)filename, MPI_MODE_RDONLY, info, &fh))
+ goto done;
+
+ if (MPI_SUCCESS != (MPI_File_get_size(fh, &size)))
+ goto done;
+
+ if (MPI_SUCCESS != MPI_File_close(&fh))
+ size=0;
+
+done:
+ return(size);
+}
int main(int argc, char **argv)
{
@@ -445,7 +470,15 @@ int main(int argc, char **argv)
MPI_BANNER("Independent test skipped");
}
- if (!(dowrite || doread || ndatasets || ngroups || docompact)){
+ if (dobig && sizeof(MPI_Offset)>4){
+ MPI_BANNER("big dataset test...");
+ big_dataset(filenames[7]);
+ }
+ else {
+ MPI_BANNER("big dataset test skipped");
+ }
+
+ if (!(dowrite || doread || ndatasets || ngroups || docompact || doindependent || dobig)){
usage();
nerrors++;
}
diff --git a/testpar/testphdf5.h b/testpar/testphdf5.h
index afbce3e..b6f5dc2 100644
--- a/testpar/testphdf5.h
+++ b/testpar/testphdf5.h
@@ -123,6 +123,7 @@ extern int facc_type; /*Test file access type */
/* prototypes */
hid_t create_faccess_plist(MPI_Comm comm, MPI_Info info, int l_facc_type, hbool_t use_gpfs);
+MPI_Offset h5_mpi_get_file_size(const char *filename, MPI_Comm comm, MPI_Info info);
void multiple_dset_write(char *filename, int ndatasets);
void multiple_group_write(char *filename, int ngroups);
void multiple_group_read(char *filename, int ngroups);
@@ -140,6 +141,7 @@ void dataset_readAll(char *filename);
void extend_readInd(char *filename);
void extend_readAll(char *filename);
void compact_dataset(char *filename);
+void big_dataset(const char *filename);
int dataset_vrfy(hssize_t start[], hsize_t count[], hsize_t stride[], hsize_t block[], DATATYPE *dataset, DATATYPE *original);
#endif /* PHDF5TEST_H */