|
| 1 | +/* This is part of the netCDF package. Copyright 2020 University |
| 2 | + Corporation for Atmospheric Research/Unidata See COPYRIGHT file for |
| 3 | + conditions of use. |
| 4 | +
|
| 5 | + Test HDF5 file code. These are not intended to be exhaustive tests, |
| 6 | + but they use HDF5 the same way that netCDF-4 does, so if these |
| 7 | + tests don't work, than netCDF-4 won't work either. |
| 8 | +
|
| 9 | + This files tests parallel I/O using compression filters. This |
| 10 | + functionality is only available in HDF5-1.10.3 and later versions. |
| 11 | +
|
| 12 | + Ed Hartnett |
| 13 | +*/ |
| 14 | +#include <nc_tests.h> |
| 15 | +#include "err_macros.h" |
| 16 | +#include <hdf5.h> |
| 17 | + |
| 18 | +#define FILE_NAME "tst_h_par_compress.h5" |
| 19 | +#define VAR_NAME "HALs_memory" |
| 20 | +#define NDIMS 1 |
| 21 | +#define MILLION 1000000 |
| 22 | +#define DIM2_LEN 16000000 |
| 23 | +#define SC1 100000 /* slice count. */ |
| 24 | + |
| 25 | +/* The following code, when uncommented, adds szip testing for |
| 26 | + * parallel I/O. However, this currently fails. I have a support |
| 27 | + * request in to HDF5 about this. Ed 7/8/20 */ |
| 28 | +/* #ifdef USE_SZIP */ |
| 29 | +/* #define NUM_COMPRESS_FILTERS 2 */ |
| 30 | +/* #else */ |
| 31 | +/* #define NUM_COMPRESS_FILTERS 1 */ |
| 32 | +/* #endif /\* USE_SZIP *\/ */ |
| 33 | +#define NUM_COMPRESS_FILTERS 1 |
| 34 | + |
| 35 | +int |
| 36 | +main(int argc, char **argv) |
| 37 | +{ |
| 38 | + int p, my_rank; |
| 39 | + |
| 40 | + MPI_Init(&argc, &argv); |
| 41 | + MPI_Comm_rank(MPI_COMM_WORLD, &my_rank); |
| 42 | + MPI_Comm_size(MPI_COMM_WORLD, &p); |
| 43 | + |
| 44 | + /* For builds with HDF5 prior to 1.10.3, just return success. */ |
| 45 | +#ifdef HDF5_SUPPORTS_PAR_FILTERS |
| 46 | + for (int cf = 0; cf < NUM_COMPRESS_FILTERS; cf++) |
| 47 | + { |
| 48 | + if (!my_rank) |
| 49 | + printf("*** Testing parallel I/O with %s compression...", cf ? "szip" : "zlib"); |
| 50 | + { |
| 51 | + hid_t fapl_id, fileid, whole_spaceid, dsid, slice_spaceid, whole_spaceid1, xferid; |
| 52 | + hid_t plistid; |
| 53 | + hsize_t start[NDIMS], count[NDIMS]; |
| 54 | + hsize_t dims[1], chunksize = SC1; |
| 55 | + int data[SC1], data_in[SC1]; |
| 56 | + int num_steps; |
| 57 | + int deflate_level = 4; |
| 58 | + int i, s; |
| 59 | + |
| 60 | + /* We will write the same slice of random data over and over to |
| 61 | + * fill the file. */ |
| 62 | + for (i = 0; i < SC1; i++) |
| 63 | + data[i] = rand(); |
| 64 | + |
| 65 | + /* Create file. */ |
| 66 | + if ((fapl_id = H5Pcreate(H5P_FILE_ACCESS)) < 0) ERR; |
| 67 | + if (H5Pset_fapl_mpio(fapl_id, MPI_COMM_WORLD, MPI_INFO_NULL) < 0) ERR; |
| 68 | + if ((fileid = H5Fcreate(FILE_NAME, H5F_ACC_TRUNC, H5P_DEFAULT, |
| 69 | + fapl_id)) < 0) ERR; |
| 70 | + |
| 71 | + /* Create a space to deal with one slice in memory. */ |
| 72 | + dims[0] = SC1; |
| 73 | + if ((slice_spaceid = H5Screate_simple(NDIMS, dims, NULL)) < 0) ERR; |
| 74 | + |
| 75 | + /* Create a space to write all slices. */ |
| 76 | + dims[0] = DIM2_LEN; |
| 77 | + if ((whole_spaceid = H5Screate_simple(NDIMS, dims, NULL)) < 0) ERR; |
| 78 | + |
| 79 | + /* Create property list for dataset. */ |
| 80 | + if ((plistid = H5Pcreate(H5P_DATASET_CREATE)) < 0) ERR; |
| 81 | + |
| 82 | + /* Turn off object tracking times in HDF5 (as is done in nc4hdf.c). */ |
| 83 | + if (H5Pset_obj_track_times(plistid, 0) < 0) ERR; |
| 84 | + |
| 85 | + /* Required to truly turn HDF5 fill values off */ |
| 86 | + if (H5Pset_fill_time(plistid, H5D_FILL_TIME_NEVER) < 0) ERR; |
| 87 | + |
| 88 | + /* Set compression, either deflate or szip. */ |
| 89 | + if (cf == 0) |
| 90 | + { |
| 91 | + if (H5Pset_deflate(plistid, deflate_level) < 0) ERR; |
| 92 | + } |
| 93 | + else |
| 94 | + { |
| 95 | + int options_mask = 32; |
| 96 | + int bits_per_pixel = 32; |
| 97 | + if (H5Pset_szip(plistid, options_mask, bits_per_pixel)) ERR; |
| 98 | + } |
| 99 | + |
| 100 | + /* Set chunking. */ |
| 101 | + if (H5Pset_chunk(plistid, NDIMS, &chunksize) < 0) ERR; |
| 102 | + |
| 103 | + /* Turn on creation order tracking. */ |
| 104 | + if (H5Pset_attr_creation_order(plistid, H5P_CRT_ORDER_TRACKED| |
| 105 | + H5P_CRT_ORDER_INDEXED) < 0) ERR; |
| 106 | + |
| 107 | + /* Create dataset. */ |
| 108 | + if ((dsid = H5Dcreate2(fileid, VAR_NAME, H5T_NATIVE_INT, |
| 109 | + whole_spaceid, H5P_DEFAULT, plistid, H5P_DEFAULT)) < 0) ERR; |
| 110 | + |
| 111 | + /* Use collective write operations. */ |
| 112 | + if ((xferid = H5Pcreate(H5P_DATASET_XFER)) < 0) ERR; |
| 113 | + if (H5Pset_dxpl_mpio(xferid, H5FD_MPIO_COLLECTIVE) < 0) ERR; |
| 114 | + |
| 115 | + /* Write the data in num_step steps. */ |
| 116 | + num_steps = (DIM2_LEN/SC1) / p; |
| 117 | + for (s = 0; s < num_steps; s++) |
| 118 | + { |
| 119 | + /* Select hyperslab for write of one slice. */ |
| 120 | + start[0] = s * SC1 * p + my_rank * SC1; |
| 121 | + count[0] = SC1; |
| 122 | + if (H5Sselect_hyperslab(whole_spaceid, H5S_SELECT_SET, |
| 123 | + start, NULL, count, NULL) < 0) ERR; |
| 124 | + |
| 125 | + if (H5Dwrite(dsid, H5T_NATIVE_INT, slice_spaceid, whole_spaceid, |
| 126 | + xferid, data) < 0) ERR; |
| 127 | + |
| 128 | + } |
| 129 | + |
| 130 | + /* Close. These collective operations will allow every process |
| 131 | + * to catch up. */ |
| 132 | + if (H5Dclose(dsid) < 0 || |
| 133 | + H5Sclose(whole_spaceid) < 0 || |
| 134 | + H5Sclose(slice_spaceid) < 0 || |
| 135 | + H5Pclose(fapl_id) < 0 || |
| 136 | + H5Pclose(plistid) < 0 || |
| 137 | + H5Fclose(fileid) < 0) |
| 138 | + ERR; |
| 139 | + |
| 140 | + /* Open the file. */ |
| 141 | + if ((fapl_id = H5Pcreate(H5P_FILE_ACCESS)) < 0) ERR; |
| 142 | + if (H5Pset_fapl_mpio(fapl_id, MPI_COMM_WORLD, MPI_INFO_NULL) < 0) ERR; |
| 143 | + |
| 144 | + |
| 145 | + if (H5Pset_libver_bounds(fapl_id, H5F_LIBVER_LATEST, H5F_LIBVER_LATEST) < 0) ERR; |
| 146 | + if ((fileid = H5Fopen(FILE_NAME, H5F_ACC_RDONLY, fapl_id)) < 0) ERR; |
| 147 | + |
| 148 | + /* Create a space to deal with one slice in memory. */ |
| 149 | + dims[0] = SC1; |
| 150 | + if ((slice_spaceid = H5Screate_simple(NDIMS, dims, NULL)) < 0) ERR; |
| 151 | + |
| 152 | + /* Open the dataset. */ |
| 153 | + if ((dsid = H5Dopen(fileid, VAR_NAME)) < 0) ERR; |
| 154 | + if ((whole_spaceid1 = H5Dget_space(dsid)) < 0) ERR; |
| 155 | + |
| 156 | + /* Read the data, a slice at a time. */ |
| 157 | + for (s = 0; s < num_steps; s++) |
| 158 | + { |
| 159 | + /* Select hyperslab for read of one slice. */ |
| 160 | + start[0] = s * SC1 * p + my_rank * SC1; |
| 161 | + count[0] = SC1; |
| 162 | + if (H5Sselect_hyperslab(whole_spaceid1, H5S_SELECT_SET, |
| 163 | + start, NULL, count, NULL) < 0) |
| 164 | + { |
| 165 | + ERR; |
| 166 | + return 2; |
| 167 | + } |
| 168 | + |
| 169 | + if (H5Dread(dsid, H5T_NATIVE_INT, slice_spaceid, whole_spaceid1, |
| 170 | + H5P_DEFAULT, data_in) < 0) |
| 171 | + { |
| 172 | + ERR; |
| 173 | + return 2; |
| 174 | + } |
| 175 | + |
| 176 | + /* Check the slice of data. */ |
| 177 | + for (i = 0; i < SC1; i++) |
| 178 | + if (data[i] != data_in[i]) |
| 179 | + { |
| 180 | + ERR; |
| 181 | + return 2; |
| 182 | + } |
| 183 | + } |
| 184 | + |
| 185 | + /* Close down. */ |
| 186 | + if (H5Dclose(dsid) < 0 || |
| 187 | + H5Sclose(slice_spaceid) < 0 || |
| 188 | + H5Sclose(whole_spaceid1) < 0 || |
| 189 | + H5Pclose(fapl_id) < 0 || |
| 190 | + H5Fclose(fileid) < 0) |
| 191 | + ERR; |
| 192 | + } |
| 193 | + if (!my_rank) |
| 194 | + SUMMARIZE_ERR; |
| 195 | + |
| 196 | + } /* next cf */ |
| 197 | +#else |
| 198 | + { |
| 199 | + if (!my_rank) |
| 200 | + printf("*** HDF5 1.10.3 or greater required for this test.\n"); |
| 201 | + } |
| 202 | + |
| 203 | +#endif /* HDF5_SUPPORTS_PAR_FILTERS */ |
| 204 | + |
| 205 | + MPI_Finalize(); |
| 206 | + |
| 207 | + if (!my_rank) |
| 208 | + FINAL_RESULTS; |
| 209 | + return 0; |
| 210 | +} |
0 commit comments