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
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
|
/* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
* Copyright by the Board of Trustees of the University of Illinois. *
* 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 files COPYING and Copyright.html. COPYING can be found at the root *
* of the source code distribution tree; Copyright.html can be found at the *
* root level of an installed copy of the electronic HDF5 document set and *
* is linked from the top-level documents page. It can also be found at *
* http://hdf.ncsa.uiuc.edu/HDF5/doc/Copyright.html. If you do not have *
* access to either file, you may request a copy from hdfhelp@ncsa.uiuc.edu. *
* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
/*
This program will test irregular hyperslab selections with collective write and read.
The way to test whether collective write and read works is to use independent IO
output to verify the collective output.
1) We will write two datasets with the same hyperslab selection settings;
one in independent mode,
one in collective mode,
2) We will read two datasets with the same hyperslab selection settings,
1. independent read to read independent output,
independent read to read collecive output,
Compare the result,
If the result is the same, then collective write succeeds.
2. collective read to read independent output,
independent read to read independent output,
Compare the result,
If the result is the same, then collective read succeeds.
*/
#include "hdf5.h"
#include "H5private.h"
#include "testphdf5.h"
void coll_irregular_cont_write()
{
char *filename;
hid_t acc_plist,xfer_plist;
hid_t file, datasetc,dataseti; /* File and dataset identifiers */
hid_t mspaceid1, mspaceid, fspaceid,fspaceid1; /* Dataspace identifiers */
hid_t plist; /* Dataset property list identifier */
hsize_t mdim1[] = {MSPACE1_DIM}; /* Dimension size of the first dataset
(in memory) */
hsize_t fsdim[] = {FSPACE_DIM1, FSPACE_DIM2};
/* Dimension sizes of the dataset (on disk) */
hsize_t mdim[] = {MSPACE_DIM1, MSPACE_DIM2}; /* Dimension sizes of the
dataset in memory when we
read selection from the
dataset on the disk */
hssize_t start[2]; /* Start of hyperslab */
hsize_t stride[2]; /* Stride of hyperslab */
hsize_t count[2]; /* Block count */
hsize_t block[2]; /* Block sizes */
herr_t ret;
unsigned i,j;
int fillvalue = 0; /* Fill value for the dataset */
int matrix_out[MSPACE_DIM1][MSPACE_DIM2];
int matrix_out1[MSPACE_DIM1][MSPACE_DIM2]; /* Buffer to read from the
dataset */
int vector[MSPACE1_DIM];
int mpi_size,mpi_rank;
MPI_Comm comm = MPI_COMM_WORLD;
MPI_Info info = MPI_INFO_NULL;
/*set up MPI parameters */
MPI_Comm_size(comm,&mpi_size);
MPI_Comm_rank(comm,&mpi_rank);
/* Obtain file name */
filename = (char *) GetTestParameters();
/*
* Buffers' initialization.
*/
vector[0] = vector[MSPACE1_DIM - 1] = -1;
for (i = 1; i < MSPACE1_DIM - 1; i++) vector[i] = i;
acc_plist = H5Pcreate(H5P_FILE_ACCESS);
VRFY((acc_plist >= 0),"");
ret = H5Pset_fapl_mpio(acc_plist,comm,info);
VRFY((ret >= 0),"MPIO creation property list succeeded");
/*
* Create a file.
*/
file = H5Fcreate(filename, H5F_ACC_TRUNC, H5P_DEFAULT, acc_plist);
VRFY((file >= 0),"H5Fcreate succeeded");
/*
* Create property list for a dataset and set up fill values.
*/
plist = H5Pcreate(H5P_DATASET_CREATE);
VRFY((acc_plist >= 0),"");
ret = H5Pset_fill_value(plist, H5T_NATIVE_INT, &fillvalue);
VRFY((ret >= 0),"Fill value creation property list succeeded");
/*
* Create dataspace for the dataset in the file.
*/
fspaceid = H5Screate_simple(FSPACE_RANK, fsdim, NULL);
VRFY((fspaceid >= 0),"file dataspace created succeeded");
/*
* Create dataset in the file. Notice that creation
* property list plist is used.
*/
datasetc = H5Dcreate(file, "collect_write", H5T_NATIVE_INT, fspaceid, plist);
VRFY((datasetc >= 0),"dataset created succeeded");
dataseti = H5Dcreate(file, "independ_write", H5T_NATIVE_INT, fspaceid, plist);
VRFY((dataseti >= 0),"dataset created succeeded");
/*
* Select hyperslab for the dataset in the file, using 3x2 blocks,
* (4,3) stride and (1,4) count starting at the position (0,1)
for the first selection
*/
start[0] = FHSTART0;
start[1] = FHSTART1+mpi_rank*FHSTRIDE1*FHCOUNT1/mpi_size;
stride[0] = FHSTRIDE0;
stride[1] = FHSTRIDE1;
count[0] = FHCOUNT0;
count[1] = FHCOUNT1/mpi_size;
block[0] = FHBLOCK0;
block[1] = FHBLOCK1;
ret = H5Sselect_hyperslab(fspaceid, H5S_SELECT_SET, start, stride, count, block);
VRFY((ret >= 0),"hyperslab selection succeeded");
/*
* Select hyperslab for the dataset in the file, using 3x2*4 blocks,
* stride 1 and (1,1) count starting at the position (4,0).
*/
start[0] = SHSTART0;
start[1] = SHSTART1+SHCOUNT1*SHBLOCK1*mpi_rank/mpi_size;
stride[0] = SHSTRIDE0;
stride[1] = SHSTRIDE1;
count[0] = SHCOUNT0;
count[1] = SHCOUNT1;
block[0] = SHBLOCK0;
block[1] = SHBLOCK1/mpi_size;
ret = H5Sselect_hyperslab(fspaceid, H5S_SELECT_OR, start, stride, count, block);
VRFY((ret >= 0),"hyperslab selection succeeded");
/*
* Create dataspace for the first dataset.
*/
mspaceid1 = H5Screate_simple(MSPACE1_RANK, mdim1, NULL);
VRFY((mspaceid1 >= 0),"memory dataspace created succeeded");
/*
* Select hyperslab.
* We will use 48 elements of the vector buffer starting at the second element.
* Selected elements are 1 2 3 . . . 48
*/
start[0] = MHSTART0;
stride[0] = MHSTRIDE0;
count[0] = MHCOUNT0/mpi_size;
block[0] = MHBLOCK0;
ret = H5Sselect_hyperslab(mspaceid1, H5S_SELECT_SET, start, stride, count, block);
VRFY((ret >= 0),"hyperslab selection succeeded");
ret = H5Dwrite(dataseti, H5T_NATIVE_INT, mspaceid1, fspaceid, H5P_DEFAULT, vector);
VRFY((ret >= 0),"dataset independent write succeed");
xfer_plist = H5Pcreate(H5P_DATASET_XFER);
VRFY((xfer_plist >= 0),"");
ret = H5Pset_dxpl_mpio(xfer_plist, H5FD_MPIO_COLLECTIVE);
VRFY((ret >= 0),"MPIO data transfer property list succeed");
ret = H5Dwrite(datasetc, H5T_NATIVE_INT, mspaceid1, fspaceid, xfer_plist, vector);
VRFY((ret >= 0),"dataset collective write succeed");
ret = H5Sclose(mspaceid1);
VRFY((ret >= 0),"");
ret = H5Sclose(fspaceid);
VRFY((ret >= 0),"");
/*
* Close dataset.
*/
ret = H5Dclose(datasetc);
VRFY((ret >= 0),"");
ret = H5Dclose(dataseti);
VRFY((ret >= 0),"");
/*
* Close the file.
*/
ret = H5Fclose(file);
VRFY((ret >= 0),"");
/*
* Close property list
*/
ret = H5Pclose(acc_plist);
VRFY((ret >= 0),"");
ret = H5Pclose(xfer_plist);
VRFY((ret >= 0),"");
ret = H5Pclose(plist);
VRFY((ret >= 0),"");
/*
* Open the file.
*/
/*** For testing collective hyperslab selection write ***/
acc_plist = H5Pcreate(H5P_FILE_ACCESS);
VRFY((acc_plist >= 0),"");
ret = H5Pset_fapl_mpio(acc_plist,comm,info);
VRFY((ret >= 0),"MPIO creation property list succeeded");
file = H5Fopen(filename, H5F_ACC_RDONLY, acc_plist);
VRFY((file >= 0),"H5Fopen succeeded");
/*
* Open the dataset.
*/
datasetc = H5Dopen(file,"collect_write");
VRFY((datasetc >= 0),"H5Dopen succeeded");
dataseti = H5Dopen(file,"independ_write");
VRFY((dataseti >= 0),"H5Dopen succeeded");
/*
* Get dataspace of the open dataset.
*/
fspaceid = H5Dget_space(datasetc);
VRFY((fspaceid >= 0),"file dataspace obtained succeeded");
fspaceid1 = H5Dget_space(dataseti);
VRFY((fspaceid1 >= 0),"file dataspace obtained succeeded");
start[0] = RFFHSTART0;
start[1] = RFFHSTART1+mpi_rank*RFFHCOUNT1/mpi_size;
block[0] = RFFHBLOCK0;
block[1] = RFFHBLOCK1;
stride[0] = RFFHSTRIDE0;
stride[1] = RFFHSTRIDE1;
count[0] = RFFHCOUNT0;
count[1] = RFFHCOUNT1/mpi_size;
ret = H5Sselect_hyperslab(fspaceid, H5S_SELECT_SET, start, stride, count, block);
VRFY((ret >= 0),"hyperslab selection succeeded");
ret = H5Sselect_hyperslab(fspaceid1, H5S_SELECT_SET, start, stride, count, block);
VRFY((ret >= 0),"hyperslab selection succeeded");
/*start[0] = RFSHSTART0+mpi_rank*RFSHCOUNT1/mpi_size; */
start[0] = RFSHSTART0;
start[1] = RFSHSTART1+RFSHCOUNT1*mpi_rank/mpi_size;
block[0] = RFSHBLOCK0;
block[1] = RFSHBLOCK1;
stride[0] = RFSHSTRIDE0;
stride[1] = RFSHSTRIDE0;
count[0] = RFSHCOUNT0;
count[1] = RFSHCOUNT1/mpi_size;
ret = H5Sselect_hyperslab(fspaceid, H5S_SELECT_OR, start, stride, count, block);
VRFY((ret >= 0),"hyperslab selection succeeded");
ret = H5Sselect_hyperslab(fspaceid1, H5S_SELECT_OR, start, stride, count, block);
VRFY((ret >= 0),"hyperslab selection succeeded");
/*
* Create memory dataspace.
*/
mspaceid = H5Screate_simple(MSPACE_RANK, mdim, NULL);
/*
* Select two hyperslabs in memory. Hyperslabs has the same
* size and shape as the selected hyperslabs for the file dataspace.
*/
start[0] = RMFHSTART0;
start[1] = RMFHSTART1+mpi_rank*RMFHCOUNT1/mpi_size;
block[0] = RMFHBLOCK0;
block[1] = RMFHBLOCK1;
stride[0] = RMFHSTRIDE0;
stride[1] = RMFHSTRIDE1;
count[0] = RMFHCOUNT0;
count[1] = RMFHCOUNT1/mpi_size;
ret = H5Sselect_hyperslab(mspaceid, H5S_SELECT_SET, start, stride, count, block);
VRFY((ret >= 0),"hyperslab selection succeeded");
start[0] = RMSHSTART0;
start[1] = RMSHSTART1+mpi_rank*RMSHCOUNT1/mpi_size;
block[0] = RMSHBLOCK0;
block[1] = RMSHBLOCK1;
stride[0] = RMSHSTRIDE0;
stride[1] = RMSHSTRIDE1;
count[0] = RMSHCOUNT0;
count[1] = RMSHCOUNT1/mpi_size;
ret = H5Sselect_hyperslab(mspaceid, H5S_SELECT_OR, start, stride, count, block);
VRFY((ret >= 0),"hyperslab selection succeeded");
/*
* Initialize data buffer.
*/
for (i = 0; i < MSPACE_DIM1; i++) {
for (j = 0; j < MSPACE_DIM2; j++)
matrix_out[i][j] = 0;
}
/*
* Read data back to the buffer matrix_out.
*/
ret = H5Dread(datasetc, H5T_NATIVE_INT, mspaceid, fspaceid,
H5P_DEFAULT, matrix_out);
VRFY((ret >= 0),"H5D independent read succeed");
for (i = 0; i < MSPACE_DIM1; i++) {
for (j = 0; j < MSPACE_DIM2; j++)
matrix_out1[i][j] = 0;
}
ret = H5Dread(dataseti, H5T_NATIVE_INT, mspaceid, fspaceid,
H5P_DEFAULT, matrix_out1);
VRFY((ret >= 0),"H5D independent read succeed");
ret = 0;
for (i = 0; i < MSPACE_DIM1; i++){
for (j = 0; j < MSPACE_DIM2; j++){
if(matrix_out[i][j]!=matrix_out1[i][j]) ret = -1;
if(ret < 0) break;
}
}
VRFY((ret >= 0),"H5D contiguous irregular collective write succeed");
/*
* Close memory file and memory dataspaces.
*/
ret = H5Sclose(mspaceid);
VRFY((ret >= 0),"");
ret = H5Sclose(fspaceid);
VRFY((ret >= 0),"");
/*
* Close dataset.
*/
ret = H5Dclose(dataseti);
VRFY((ret >= 0),"");
ret = H5Dclose(datasetc);
VRFY((ret >= 0),"");
/*
* Close property list
*/
ret = H5Pclose(acc_plist);
VRFY((ret >= 0),"");
/*
* Close the file.
*/
ret = H5Fclose(file);
VRFY((ret >= 0),"");
return ;
}
void coll_irregular_cont_read()
{
char *filename;
hid_t acc_plist,xfer_plist;
hid_t file, dataseti; /* File and dataset identifiers */
hid_t mspaceid, fspaceid1; /* Dataspace identifiers */
/* Dimension sizes of the dataset (on disk) */
hsize_t mdim[] = {MSPACE_DIM1, MSPACE_DIM2}; /* Dimension sizes of the
dataset in memory when we
read selection from the
dataset on the disk */
hssize_t start[2]; /* Start of hyperslab */
hsize_t stride[2]; /* Stride of hyperslab */
hsize_t count[2]; /* Block count */
hsize_t block[2]; /* Block sizes */
herr_t ret;
unsigned i,j;
int matrix_out[MSPACE_DIM1][MSPACE_DIM2];
int matrix_out1[MSPACE_DIM1][MSPACE_DIM2]; /* Buffer to read from the
dataset */
int mpi_size,mpi_rank;
MPI_Comm comm = MPI_COMM_WORLD;
MPI_Info info = MPI_INFO_NULL;
/*set up MPI parameters */
MPI_Comm_size(comm,&mpi_size);
MPI_Comm_rank(comm,&mpi_rank);
/* Obtain file name */
filename = (char *) GetTestParameters();
/*
* Buffers' initialization.
*/
/*
* Open the file.
*/
/*** For testing collective hyperslab selection read ***/
acc_plist = H5Pcreate(H5P_FILE_ACCESS);
VRFY((acc_plist >= 0),"");
ret = H5Pset_fapl_mpio(acc_plist,comm,info);
VRFY((ret >= 0),"MPIO creation property list succeeded");
file = H5Fopen(filename, H5F_ACC_RDONLY, acc_plist);
VRFY((file >= 0),"H5Fopen succeeded");
/*
* Open the dataset.
*/
dataseti = H5Dopen(file,"independ_write");
VRFY((dataseti >= 0),"H5Dopen succeeded");
/*
* Get dataspace of the open dataset.
*/
fspaceid1 = H5Dget_space(dataseti);
VRFY((fspaceid1 >= 0),"file dataspace obtained succeeded");
start[0] = RFFHSTART0;
start[1] = RFFHSTART1+mpi_rank*RFFHCOUNT1/mpi_size;
block[0] = RFFHBLOCK0;
block[1] = RFFHBLOCK1;
stride[0] = RFFHSTRIDE0;
stride[1] = RFFHSTRIDE1;
count[0] = RFFHCOUNT0;
count[1] = RFFHCOUNT1/mpi_size;
ret = H5Sselect_hyperslab(fspaceid1, H5S_SELECT_SET, start, stride, count, block);
VRFY((ret >= 0),"hyperslab selection succeeded");
start[0] = RFSHSTART0;
start[1] = RFSHSTART1+RFSHCOUNT1*mpi_rank/mpi_size;
block[0] = RFSHBLOCK0;
block[1] = RFSHBLOCK1;
stride[0] = RFSHSTRIDE0;
stride[1] = RFSHSTRIDE0;
count[0] = RFSHCOUNT0;
count[1] = RFSHCOUNT1/mpi_size;
ret = H5Sselect_hyperslab(fspaceid1, H5S_SELECT_OR, start, stride, count, block);
VRFY((ret >= 0),"hyperslab selection succeeded");
/*
* Create memory dataspace.
*/
mspaceid = H5Screate_simple(MSPACE_RANK, mdim, NULL);
/*
* Select two hyperslabs in memory. Hyperslabs has the same
* size and shape as the selected hyperslabs for the file dataspace.
*/
start[0] = RMFHSTART0;
start[1] = RMFHSTART1+mpi_rank*RMFHCOUNT1/mpi_size;
block[0] = RMFHBLOCK0;
block[1] = RMFHBLOCK1;
stride[0] = RMFHSTRIDE0;
stride[1] = RMFHSTRIDE1;
count[0] = RMFHCOUNT0;
count[1] = RMFHCOUNT1/mpi_size;
ret = H5Sselect_hyperslab(mspaceid, H5S_SELECT_SET, start, stride, count, block);
VRFY((ret >= 0),"hyperslab selection succeeded");
start[0] = RMSHSTART0;
start[1] = RMSHSTART1+mpi_rank*RMSHCOUNT1/mpi_size;
block[0] = RMSHBLOCK0;
block[1] = RMSHBLOCK1;
stride[0] = RMSHSTRIDE0;
stride[1] = RMSHSTRIDE1;
count[0] = RMSHCOUNT0;
count[1] = RMSHCOUNT1/mpi_size;
ret = H5Sselect_hyperslab(mspaceid, H5S_SELECT_OR, start, stride, count, block);
VRFY((ret >= 0),"hyperslab selection succeeded");
/*
* Initialize data buffer.
*/
for (i = 0; i < MSPACE_DIM1; i++) {
for (j = 0; j < MSPACE_DIM2; j++)
matrix_out[i][j] = 0;
}
/*
* Read data back to the buffer matrix_out.
*/
xfer_plist = H5Pcreate(H5P_DATASET_XFER);
VRFY((xfer_plist >= 0),"");
ret = H5Pset_dxpl_mpio(xfer_plist, H5FD_MPIO_COLLECTIVE);
VRFY((ret >= 0),"MPIO data transfer property list succeed");
ret = H5Dread(dataseti, H5T_NATIVE_INT, mspaceid, fspaceid1,
xfer_plist, matrix_out);
VRFY((ret >= 0),"H5D collecive read succeed");
ret = H5Pclose(xfer_plist);
VRFY((ret >= 0),"");
for (i = 0; i < MSPACE_DIM1; i++) {
for (j = 0; j < MSPACE_DIM2; j++)
matrix_out1[i][j] = 0;
}
ret = H5Dread(dataseti, H5T_NATIVE_INT, mspaceid, fspaceid1,
H5P_DEFAULT, matrix_out1);
VRFY((ret >= 0),"H5D independent read succeed");
ret = 0;
for (i = 0; i < MSPACE_DIM1; i++){
for (j = 0; j < MSPACE_DIM2; j++){
if(matrix_out[i][j]!=matrix_out1[i][j])ret = -1;
if(ret < 0) break;
}
}
VRFY((ret >= 0),"H5D contiguous irregular collective read succeed");
/*
* Close memory file and memory dataspaces.
*/
ret = H5Sclose(mspaceid);
VRFY((ret >= 0),"");
ret = H5Sclose(fspaceid1);
VRFY((ret >= 0),"");
/*
* Close dataset.
*/
ret = H5Dclose(dataseti);
VRFY((ret >= 0),"");
/*
* Close property list
*/
ret = H5Pclose(acc_plist);
VRFY((ret >= 0),"");
/*
* Close the file.
*/
ret = H5Fclose(file);
VRFY((ret >= 0),"");
return ;
}
|