summaryrefslogtreecommitdiffstats
path: root/testpar/t_filter_read.c
blob: f001cc9e1445d9984ff6f42af4025af2270c7395 (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
/* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
 * Copyright by The HDF Group.                                               *
 * All rights reserved.                                                      *
 *                                                                           *
 * This file is part of HDF5.  The full HDF5 copyright notice, including     *
 * terms governing use, modification, and redistribution, is contained in    *
 * the COPYING file, which can be found at the root of the source code       *
 * distribution tree, or in https://www.hdfgroup.org/licenses.               *
 * If you do not have access to either file, you may request a copy from     *
 * help@hdfgroup.org.                                                        *
 * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */

/*
 * This verifies the correctness of parallel reading of a dataset that has been
 * written serially using filters.
 *
 */

#include "testphdf5.h"

#ifdef H5_HAVE_SZLIB_H
#include "szlib.h"
#endif

static int mpi_size, mpi_rank;

/* Chunk sizes */
#define CHUNK_DIM1 7
#define CHUNK_DIM2 27

/* Sizes of the vertical hyperslabs. Total dataset size is
   {HS_DIM1, HS_DIM2 * mpi_size } */
#define HS_DIM1 200
#define HS_DIM2 100

/*-------------------------------------------------------------------------
 * Function:    filter_read_internal
 *
 * Purpose:     Tests parallel reading of a 2D dataset written serially using
 *              filters. During the parallel reading phase, the dataset is
 *              divided evenly among the processors in vertical hyperslabs.
 *-------------------------------------------------------------------------
 */
static void
filter_read_internal(const char *filename, hid_t dcpl, hsize_t *dset_size)
{
    hid_t   file, dataset; /* HDF5 IDs */
    hid_t   access_plist;  /* Access property list ID */
    hid_t   sid, memspace; /* Dataspace IDs */
    hsize_t size[2];       /* Dataspace dimensions */
    hsize_t hs_offset[2];  /* Hyperslab offset */
    hsize_t hs_size[2];    /* Hyperslab size */
    size_t  i, j;          /* Local index variables */
    char    name[32] = "dataset";
    herr_t  hrc; /* Error status */
    bool    vol_is_native;
    int    *points = NULL; /* Writing buffer for entire dataset */
    int    *check  = NULL; /* Reading buffer for selected hyperslab */

    /* set up MPI parameters */
    MPI_Comm_size(MPI_COMM_WORLD, &mpi_size);
    MPI_Comm_rank(MPI_COMM_WORLD, &mpi_rank);

    /* set sizes for dataset and hyperslabs */
    hs_size[0] = size[0] = HS_DIM1;
    hs_size[1]           = HS_DIM2;

    size[1] = hs_size[1] * (hsize_t)mpi_size;

    hs_offset[0] = 0;
    hs_offset[1] = hs_size[1] * (hsize_t)mpi_rank;

    /* Create the data space */
    sid = H5Screate_simple(2, size, NULL);
    VRFY(sid >= 0, "H5Screate_simple");

    /* Create buffers */
    points = (int *)malloc(size[0] * size[1] * sizeof(int));
    VRFY(points != NULL, "malloc");

    check = (int *)malloc(hs_size[0] * hs_size[1] * sizeof(int));
    VRFY(check != NULL, "malloc");

    /* Initialize writing buffer with random data */
    for (i = 0; i < size[0]; i++)
        for (j = 0; j < size[1]; j++)
            points[i * size[1] + j] = (int)(i + j + 7);

    VRFY(H5Pall_filters_avail(dcpl), "Incorrect filter availability");

    /* Serial write phase */
    if (MAINPROCESS) {

        file = H5Fcreate(h5_rmprefix(filename), H5F_ACC_TRUNC, H5P_DEFAULT, H5P_DEFAULT);
        VRFY(file >= 0, "H5Fcreate");

        /* Check if native VOL is being used */
        VRFY((h5_using_native_vol(H5P_DEFAULT, file, &vol_is_native) >= 0), "h5_using_native_vol");

        /* Create the dataset */
        dataset = H5Dcreate2(file, name, H5T_NATIVE_INT, sid, H5P_DEFAULT, dcpl, H5P_DEFAULT);
        VRFY(dataset >= 0, "H5Dcreate2");

        hrc = H5Dwrite(dataset, H5T_NATIVE_INT, H5S_ALL, H5S_ALL, H5P_DEFAULT, points);
        VRFY(hrc >= 0, "H5Dwrite");

        if (vol_is_native) {
            *dset_size = H5Dget_storage_size(dataset);
            VRFY(*dset_size > 0, "H5Dget_storage_size");
        }

        hrc = H5Dclose(dataset);
        VRFY(hrc >= 0, "H5Dclose");

        hrc = H5Fclose(file);
        VRFY(hrc >= 0, "H5Fclose");
    }

    MPI_Barrier(MPI_COMM_WORLD);

    /* Parallel read phase */
    /* Set up MPIO file access property lists */
    access_plist = H5Pcreate(H5P_FILE_ACCESS);
    VRFY((access_plist >= 0), "H5Pcreate");

    hrc = H5Pset_fapl_mpio(access_plist, MPI_COMM_WORLD, MPI_INFO_NULL);
    VRFY((hrc >= 0), "H5Pset_fapl_mpio");

    /* Open the file */
    file = H5Fopen(filename, H5F_ACC_RDWR, access_plist);
    VRFY((file >= 0), "H5Fopen");

    /* Check if native VOL is being used */
    VRFY((h5_using_native_vol(H5P_DEFAULT, file, &vol_is_native) >= 0), "h5_using_native_vol");

    dataset = H5Dopen2(file, name, H5P_DEFAULT);
    VRFY((dataset >= 0), "H5Dopen2");

    hrc = H5Sselect_hyperslab(sid, H5S_SELECT_SET, hs_offset, NULL, hs_size, NULL);
    VRFY(hrc >= 0, "H5Sselect_hyperslab");

    memspace = H5Screate_simple(2, hs_size, NULL);
    VRFY(memspace >= 0, "H5Screate_simple");

    hrc = H5Dread(dataset, H5T_NATIVE_INT, memspace, sid, H5P_DEFAULT, check);
    VRFY(hrc >= 0, "H5Dread");

    /* Check that the values read are the same as the values written */
    for (i = 0; i < hs_size[0]; i++) {
        for (j = 0; j < hs_size[1]; j++) {
            if (points[i * size[1] + (size_t)hs_offset[1] + j] != check[i * hs_size[1] + j]) {
                fprintf(stderr, "    Read different values than written.\n");
                fprintf(stderr, "    At index %lu,%lu\n", (unsigned long)(i),
                        (unsigned long)(hs_offset[1] + j));
                fprintf(stderr, "    At original: %d\n", (int)points[i * size[1] + (size_t)hs_offset[1] + j]);
                fprintf(stderr, "    At returned: %d\n", (int)check[i * hs_size[1] + j]);
                VRFY(false, "");
            }
        }
    }

    if (vol_is_native) {
        /* Get the storage size of the dataset */
        *dset_size = H5Dget_storage_size(dataset);
        VRFY(*dset_size != 0, "H5Dget_storage_size");
    }

    /* Clean up objects used for this test */
    hrc = H5Dclose(dataset);
    VRFY(hrc >= 0, "H5Dclose");

    hrc = H5Sclose(sid);
    VRFY(hrc >= 0, "H5Sclose");

    hrc = H5Sclose(memspace);
    VRFY(hrc >= 0, "H5Sclose");

    hrc = H5Pclose(access_plist);
    VRFY(hrc >= 0, "H5Pclose");

    hrc = H5Fclose(file);
    VRFY(hrc >= 0, "H5Fclose");

    free(points);
    free(check);

    MPI_Barrier(MPI_COMM_WORLD);
}

/*-------------------------------------------------------------------------
 * Function:    test_filter_read
 *
 * Purpose:    Tests parallel reading of datasets written serially using
 *              several (combinations of) filters.
 *-------------------------------------------------------------------------
 */

void
test_filter_read(void)
{
    hid_t         dc;                                       /* HDF5 IDs */
    const hsize_t chunk_size[2] = {CHUNK_DIM1, CHUNK_DIM2}; /* Chunk dimensions */
    hsize_t       null_size;                                /* Size of dataset without filters */
    unsigned      chunk_opts;                               /* Chunk options */
    unsigned      disable_partial_chunk_filters; /* Whether filters are disabled on partial chunks */
    herr_t        hrc;
    const char   *filename;
    bool          vol_is_native;
    hsize_t       fletcher32_size; /* Size of dataset with Fletcher32 checksum */

#ifdef H5_HAVE_FILTER_DEFLATE
    hsize_t deflate_size; /* Size of dataset with deflate filter */
#endif                    /* H5_HAVE_FILTER_DEFLATE */

#ifdef H5_HAVE_FILTER_SZIP
    hsize_t  szip_size; /* Size of dataset with szip filter */
    unsigned szip_options_mask     = H5_SZIP_NN_OPTION_MASK;
    unsigned szip_pixels_per_block = 4;
#endif /* H5_HAVE_FILTER_SZIP */

    hsize_t shuffle_size = 0; /* Size of dataset with shuffle filter */

#if (defined H5_HAVE_FILTER_DEFLATE || defined H5_HAVE_FILTER_SZIP)
    hsize_t combo_size; /* Size of dataset with multiple filters */
#endif                  /* H5_HAVE_FILTER_DEFLATE || H5_HAVE_FILTER_SZIP */

    filename = GetTestParameters();

    if (VERBOSE_MED)
        printf("Parallel reading of dataset written with filters %s\n", filename);

    MPI_Comm_rank(MPI_COMM_WORLD, &mpi_rank);

    /* Make sure the connector supports the API functions being tested */
    if (!(vol_cap_flags_g & H5VL_CAP_FLAG_FILE_BASIC) || !(vol_cap_flags_g & H5VL_CAP_FLAG_DATASET_BASIC) ||
        !(vol_cap_flags_g & H5VL_CAP_FLAG_FILTERS)) {
        if (MAINPROCESS) {
            puts("SKIPPED");
            printf(
                "    API functions for basic file, dataset or filter aren't supported with this connector\n");
            fflush(stdout);
        }

        return;
    }

    /* Check if native VOL is being used */
    VRFY(h5_using_native_vol(H5P_DEFAULT, H5I_INVALID_HID, &vol_is_native) >= 0, "h5_using_native_vol");

    /*----------------------------------------------------------
     * STEP 0: Test without filters.
     *----------------------------------------------------------
     */
    dc = H5Pcreate(H5P_DATASET_CREATE);
    VRFY(dc >= 0, "H5Pcreate");

    hrc = H5Pset_chunk(dc, 2, chunk_size);
    VRFY(hrc >= 0, "H5Pset_chunk");

    filter_read_internal(filename, dc, &null_size);

    /* Clean up objects used for this test */
    hrc = H5Pclose(dc);
    VRFY(hrc >= 0, "H5Pclose");

    /* Run steps 1-3 both with and without filters disabled on partial chunks */
    for (disable_partial_chunk_filters = 0; disable_partial_chunk_filters <= 1;
         disable_partial_chunk_filters++) {
        /* Set chunk options appropriately */
        dc = H5Pcreate(H5P_DATASET_CREATE);
        VRFY(dc >= 0, "H5Pcreate");

        hrc = H5Pset_chunk(dc, 2, chunk_size);
        VRFY(hrc >= 0, "H5Pset_filter");

        hrc = H5Pget_chunk_opts(dc, &chunk_opts);
        VRFY(hrc >= 0, "H5Pget_chunk_opts");

        if (disable_partial_chunk_filters)
            chunk_opts |= H5D_CHUNK_DONT_FILTER_PARTIAL_CHUNKS;

        hrc = H5Pclose(dc);
        VRFY(hrc >= 0, "H5Pclose");

        /*----------------------------------------------------------
         * STEP 1: Test Fletcher32 Checksum by itself.
         *----------------------------------------------------------
         */

        dc = H5Pcreate(H5P_DATASET_CREATE);
        VRFY(dc >= 0, "H5Pset_filter");

        hrc = H5Pset_chunk(dc, 2, chunk_size);
        VRFY(hrc >= 0, "H5Pset_filter");

        hrc = H5Pset_chunk_opts(dc, chunk_opts);
        VRFY(hrc >= 0, "H5Pset_chunk_opts");

        hrc = H5Pset_filter(dc, H5Z_FILTER_FLETCHER32, 0, 0, NULL);
        VRFY(hrc >= 0, "H5Pset_filter");

        filter_read_internal(filename, dc, &fletcher32_size);

        if (vol_is_native)
            VRFY(fletcher32_size > null_size, "Size after checksumming is incorrect.");

        /* Clean up objects used for this test */
        hrc = H5Pclose(dc);
        VRFY(hrc >= 0, "H5Pclose");

        /*----------------------------------------------------------
         * STEP 2: Test deflation by itself.
         *----------------------------------------------------------
         */
#ifdef H5_HAVE_FILTER_DEFLATE

        dc = H5Pcreate(H5P_DATASET_CREATE);
        VRFY(dc >= 0, "H5Pcreate");

        hrc = H5Pset_chunk(dc, 2, chunk_size);
        VRFY(hrc >= 0, "H5Pset_chunk");

        hrc = H5Pset_chunk_opts(dc, chunk_opts);
        VRFY(hrc >= 0, "H5Pset_chunk_opts");

        hrc = H5Pset_deflate(dc, 6);
        VRFY(hrc >= 0, "H5Pset_deflate");

        filter_read_internal(filename, dc, &deflate_size);

        /* Clean up objects used for this test */
        hrc = H5Pclose(dc);
        VRFY(hrc >= 0, "H5Pclose");

#endif /* H5_HAVE_FILTER_DEFLATE */

        /*----------------------------------------------------------
         * STEP 3: Test szip compression by itself.
         *----------------------------------------------------------
         */
#ifdef H5_HAVE_FILTER_SZIP
        if (h5_szip_can_encode() == 1) {
            dc = H5Pcreate(H5P_DATASET_CREATE);
            VRFY(dc >= 0, "H5Pcreate");

            hrc = H5Pset_chunk(dc, 2, chunk_size);
            VRFY(hrc >= 0, "H5Pset_chunk");

            hrc = H5Pset_chunk_opts(dc, chunk_opts);
            VRFY(hrc >= 0, "H5Pset_chunk_opts");

            hrc = H5Pset_szip(dc, szip_options_mask, szip_pixels_per_block);
            VRFY(hrc >= 0, "H5Pset_szip");

            filter_read_internal(filename, dc, &szip_size);

            /* Clean up objects used for this test */
            hrc = H5Pclose(dc);
            VRFY(hrc >= 0, "H5Pclose");
        }
#endif /* H5_HAVE_FILTER_SZIP */
    }  /* end for */

    /*----------------------------------------------------------
     * STEP 4: Test shuffling by itself.
     *----------------------------------------------------------
     */

    dc = H5Pcreate(H5P_DATASET_CREATE);
    VRFY(dc >= 0, "H5Pcreate");

    hrc = H5Pset_chunk(dc, 2, chunk_size);
    VRFY(hrc >= 0, "H5Pset_chunk");

    hrc = H5Pset_shuffle(dc);
    VRFY(hrc >= 0, "H5Pset_shuffle");

    filter_read_internal(filename, dc, &shuffle_size);

    if (vol_is_native)
        VRFY(shuffle_size == null_size, "Shuffled size not the same as uncompressed size.");

    /* Clean up objects used for this test */
    hrc = H5Pclose(dc);
    VRFY(hrc >= 0, "H5Pclose");

    /*----------------------------------------------------------
     * STEP 5: Test shuffle + deflate + checksum in any order.
     *----------------------------------------------------------
     */
#ifdef H5_HAVE_FILTER_DEFLATE
    /* Testing shuffle+deflate+checksum filters (checksum first) */
    dc = H5Pcreate(H5P_DATASET_CREATE);
    VRFY(dc >= 0, "H5Pcreate");

    hrc = H5Pset_chunk(dc, 2, chunk_size);
    VRFY(hrc >= 0, "H5Pset_chunk");

    hrc = H5Pset_fletcher32(dc);
    VRFY(hrc >= 0, "H5Pset_fletcher32");

    hrc = H5Pset_shuffle(dc);
    VRFY(hrc >= 0, "H5Pset_shuffle");

    hrc = H5Pset_deflate(dc, 6);
    VRFY(hrc >= 0, "H5Pset_deflate");

    filter_read_internal(filename, dc, &combo_size);

    /* Clean up objects used for this test */
    hrc = H5Pclose(dc);
    VRFY(hrc >= 0, "H5Pclose");

    /* Testing shuffle+deflate+checksum filters (checksum last) */
    dc = H5Pcreate(H5P_DATASET_CREATE);
    VRFY(dc >= 0, "H5Pcreate");

    hrc = H5Pset_chunk(dc, 2, chunk_size);
    VRFY(hrc >= 0, "H5Pset_chunk");

    hrc = H5Pset_shuffle(dc);
    VRFY(hrc >= 0, "H5Pset_shuffle");

    hrc = H5Pset_deflate(dc, 6);
    VRFY(hrc >= 0, "H5Pset_deflate");

    hrc = H5Pset_fletcher32(dc);
    VRFY(hrc >= 0, "H5Pset_fletcher32");

    filter_read_internal(filename, dc, &combo_size);

    /* Clean up objects used for this test */
    hrc = H5Pclose(dc);
    VRFY(hrc >= 0, "H5Pclose");

#endif /* H5_HAVE_FILTER_DEFLATE */

    /*----------------------------------------------------------
     * STEP 6: Test shuffle + szip + checksum in any order.
     *----------------------------------------------------------
     */
#ifdef H5_HAVE_FILTER_SZIP

    /* Testing shuffle+szip(with encoder)+checksum filters(checksum first) */
    dc = H5Pcreate(H5P_DATASET_CREATE);
    VRFY(dc >= 0, "H5Pcreate");

    hrc = H5Pset_chunk(dc, 2, chunk_size);
    VRFY(hrc >= 0, "H5Pset_chunk");

    hrc = H5Pset_fletcher32(dc);
    VRFY(hrc >= 0, "H5Pset_fletcher32");

    hrc = H5Pset_shuffle(dc);
    VRFY(hrc >= 0, "H5Pset_shuffle");

    /* Make sure encoding is enabled */
    if (h5_szip_can_encode() == 1) {
        hrc = H5Pset_szip(dc, szip_options_mask, szip_pixels_per_block);
        VRFY(hrc >= 0, "H5Pset_szip");

        filter_read_internal(filename, dc, &combo_size);
    }

    /* Clean up objects used for this test */
    hrc = H5Pclose(dc);
    VRFY(hrc >= 0, "H5Pclose");

    /* Testing shuffle+szip(with encoder)+checksum filters(checksum last) */
    /* Make sure encoding is enabled */
    if (h5_szip_can_encode() == 1) {
        dc = H5Pcreate(H5P_DATASET_CREATE);
        VRFY(dc >= 0, "H5Pcreate");

        hrc = H5Pset_chunk(dc, 2, chunk_size);
        VRFY(hrc >= 0, "H5Pset_chunk");

        hrc = H5Pset_shuffle(dc);
        VRFY(hrc >= 0, "H5Pset_shuffle");

        hrc = H5Pset_szip(dc, szip_options_mask, szip_pixels_per_block);
        VRFY(hrc >= 0, "H5Pset_szip");

        hrc = H5Pset_fletcher32(dc);
        VRFY(hrc >= 0, "H5Pset_fletcher32");

        filter_read_internal(filename, dc, &combo_size);

        /* Clean up objects used for this test */
        hrc = H5Pclose(dc);
        VRFY(hrc >= 0, "H5Pclose");
    }

#endif /* H5_HAVE_FILTER_SZIP */
}