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
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
1362
1363
1364
1365
1366
1367
1368
1369
1370
1371
1372
1373
1374
1375
1376
1377
1378
1379
1380
1381
1382
1383
1384
1385
1386
1387
1388
1389
1390
1391
1392
1393
1394
1395
1396
1397
1398
1399
1400
1401
1402
1403
1404
1405
1406
1407
1408
1409
1410
1411
1412
1413
1414
1415
1416
1417
1418
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
1429
1430
1431
1432
1433
1434
1435
1436
1437
1438
1439
1440
1441
1442
1443
1444
1445
1446
1447
1448
1449
1450
1451
1452
1453
1454
1455
1456
1457
1458
1459
1460
1461
1462
1463
1464
1465
1466
1467
1468
1469
1470
1471
1472
1473
1474
1475
1476
1477
1478
1479
1480
1481
1482
1483
1484
1485
1486
1487
1488
1489
1490
1491
1492
1493
1494
1495
1496
1497
1498
1499
1500
1501
|
/* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
* 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. *
* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
/* Pablo information */
/* (Put before include files to avoid problems with inline functions) */
#define PABLO_MASK H5FD_fphdf5_mask
#include "H5private.h" /* Generic Functions */
#include "H5ACprivate.h" /* Metadata cache */
#include "H5Dprivate.h" /* Datasets */
#include "H5Eprivate.h" /* Error handling */
#include "H5Fprivate.h" /* Files access */
#include "H5FDprivate.h" /* File drivers */
#include "H5FDmpi.h" /* MPI-based file drivers */
#include "H5Iprivate.h" /* IDs */
#include "H5MMprivate.h" /* Memory management */
#include "H5Pprivate.h" /* Property lists */
#ifdef H5_HAVE_FPHDF5
#include "H5FPprivate.h" /* Flexible PHDF5 */
/*
* The driver identification number, initialized at runtime if
* H5_HAVE_FPHDF5 is defined. This allows applications to still have
* the H5FD_FPHDF5 "constants" in their source code (it also makes this
* file strictly ANSI compliant when H5_HAVE_FPHDF5 isn't defined)
*/
static hid_t H5FD_FPHDF5_g = 0;
/*
* Private Prototypes
*/
/*
* Callbacks
*/
static void *H5FD_fphdf5_fapl_get(H5FD_t *_file);
static H5FD_t *H5FD_fphdf5_open(const char *name, unsigned flags,
hid_t fapl_id, haddr_t maxaddr);
static herr_t H5FD_fphdf5_close(H5FD_t *_file);
static herr_t H5FD_fphdf5_query(const H5FD_t *_f1, unsigned long *flags);
static haddr_t H5FD_fphdf5_get_eoa(H5FD_t *_file);
static herr_t H5FD_fphdf5_set_eoa(H5FD_t *_file, haddr_t addr);
static haddr_t H5FD_fphdf5_get_eof(H5FD_t *_file);
static herr_t H5FD_fphdf5_get_handle(H5FD_t *_file, hid_t fapl,
void **file_handle);
static herr_t H5FD_fphdf5_read(H5FD_t *_file, H5FD_mem_t mem_type, hid_t dxpl_id,
haddr_t addr, size_t size, void *buf);
static herr_t H5FD_fphdf5_write(H5FD_t *_file, H5FD_mem_t type, hid_t dxpl_id,
haddr_t addr, size_t size, const void *buf);
static herr_t H5FD_fphdf5_flush(H5FD_t *_file, hid_t dxpl_id, unsigned closing);
static int H5FD_fphdf5_mpi_rank(const H5FD_t *_file);
static int H5FD_fphdf5_mpi_size(const H5FD_t *_file);
static MPI_Comm H5FD_fphdf5_barrier_communicator(const H5FD_t *_file);
/*
* FPHDF5-specific file access properties
*/
typedef struct H5FD_fphdf5_fapl_t {
MPI_Comm comm; /*communicator */
MPI_Comm barrier_comm; /*barrier communicator */
MPI_Info info; /*file information */
unsigned sap_rank; /*SAP's rank */
unsigned capt_rank; /*captain rank */
} H5FD_fphdf5_fapl_t;
/*
* The FPHDF5 file driver information
*/
const H5FD_class_mpi_t H5FD_fphdf5_g = {
{ /* Start of superclass information */
"fphdf5", /*name */
HADDR_MAX, /*maxaddr */
H5F_CLOSE_SEMI, /*fc_degree */
NULL, /*sb_size */
NULL, /*sb_encode */
NULL, /*sb_decode */
sizeof(H5FD_fphdf5_fapl_t), /*fapl_size */
H5FD_fphdf5_fapl_get, /*fapl_get */
NULL, /*fapl_copy */
NULL, /*fapl_free */
0, /*dxpl_size */
NULL, /*dxpl_copy */
NULL, /*dxpl_free */
H5FD_fphdf5_open, /*open */
H5FD_fphdf5_close, /*close */
NULL, /*cmp */
H5FD_fphdf5_query, /*query */
NULL, /*alloc */
NULL, /*free */
H5FD_fphdf5_get_eoa, /*get_eoa */
H5FD_fphdf5_set_eoa, /*set_eoa */
H5FD_fphdf5_get_eof, /*get_eof */
H5FD_fphdf5_get_handle, /*get_handle */
H5FD_fphdf5_read, /*read */
H5FD_fphdf5_write, /*write */
H5FD_fphdf5_flush, /*flush */
NULL, /*lock */
NULL, /*unlock */
H5FD_FLMAP_SINGLE /*fl_map */
}, /* End of superclass information */
H5FD_fphdf5_mpi_rank, /*get_rank */
H5FD_fphdf5_mpi_size, /*get_size */
H5FD_fphdf5_barrier_communicator /*get_comm */
};
/* Interface initialization */
#define INTERFACE_INIT H5FD_fphdf5_init
static int interface_initialize_g = 0;
/*-------------------------------------------------------------------------
* Function: H5FD_fphdf5_init
* Purpose: Initialize this driver by registering the driver with the
* library.
* Return: Success: The driver ID for the FPHDF5 driver.
* Failure: Doesn't fail.
* Programmer: Bill Wendling
* 30. January 2003
* Modifications:
*-------------------------------------------------------------------------
*/
hid_t
H5FD_fphdf5_init(void)
{
hid_t ret_value;
FUNC_ENTER_NOAPI(H5FD_fphdf5_init, FAIL)
if (H5Iget_type(H5FD_FPHDF5_g) != H5I_VFL)
H5FD_FPHDF5_g = H5FD_register((const H5FD_class_t *)&H5FD_fphdf5_g,sizeof(H5FD_class_mpi_t));
/* Set return value */
ret_value = H5FD_FPHDF5_g;
done:
FUNC_LEAVE_NOAPI(ret_value)
}
/*-------------------------------------------------------------------------
* Function: H5Pset_fapl_fphdf5
* Purpose: Store the user supplied MPIO communicator COMM and INFO
* in the file access property list FAPL_ID which can then
* be used to create and/or open the file. This function is
* available only in the parallel HDF5 library and is not
* collective.
*
* COMM is the MPI communicator to be used for file open as
* defined in MPI_FILE_OPEN of MPI-2. This function does not
* make a duplicated communicator. Any modification to COMM
* after this function call returns may have an indeterminate
* effect on the access property list. Users should not
* modify the communicator while it is defined in a property
* list.
*
* INFO is the MPI info object to be used for file open as
* defined in MPI_FILE_OPEN of MPI-2. This function does not
* make a duplicated info. Any modification to info after
* this function call returns may have an indeterminate effect
* on the access property list. Users should not modify the
* info while it is defined in a property list.
* Return: Success: SUCCEED
* Failure: FAIL
* Programmer: Bill Wendling
* 30. January 2003
* Modifications:
*-------------------------------------------------------------------------
*/
herr_t
H5Pset_fapl_fphdf5(hid_t fapl_id, MPI_Comm comm, MPI_Comm barrier_comm,
MPI_Info info, unsigned sap_rank)
{
H5FD_fphdf5_fapl_t fa;
H5P_genplist_t *plist;
int mrc, comm_size;
herr_t ret_value;
FUNC_ENTER_API(H5Pset_fapl_fphdf5, FAIL)
H5TRACE5("e","iMcMcMiIu",fapl_id,comm,barrier_comm,info,sap_rank);
if (fapl_id == H5P_DEFAULT)
HGOTO_ERROR(H5E_PLIST, H5E_BADVALUE, FAIL, "can't set values in default property list")
/* Check arguments */
if ((plist = H5P_object_verify(fapl_id, H5P_FILE_ACCESS)) == NULL)
HGOTO_ERROR(H5E_PLIST, H5E_BADTYPE, FAIL, "not a file access list")
if ((mrc = MPI_Comm_size(comm, &comm_size)) != MPI_SUCCESS)
HMPI_GOTO_ERROR(FAIL, "MPI_Comm_size failed", mrc)
/* Initialize driver specific properties */
fa.comm = comm;
fa.barrier_comm = barrier_comm;
fa.info = info;
fa.sap_rank = sap_rank;
fa.capt_rank = (sap_rank + 1) % comm_size;
ret_value = H5P_set_driver(plist, H5FD_FPHDF5, &fa);
done:
FUNC_LEAVE_API(ret_value)
}
/*-------------------------------------------------------------------------
* Function: H5Pget_fapl_fphdf5
* Purpose: If the file access property list is set to the
* H5FD_FPHDF5 driver then this function returns the MPI
* communicator and information through the COMM and INFO
* pointers.
* Return: Success: SUCCEED with the communicator and information
* returned through the COMM and INFO arguments
* if non-null. Neither piece of information is
* copied and they are therefore valid only
* until the file access property list is
* modified or closed.
* Failure: FAIL
* Programmer: Bill Wendling
* 30. January 2003
* Modifications:
*-------------------------------------------------------------------------
*/
herr_t
H5Pget_fapl_fphdf5(hid_t fapl_id, MPI_Comm *comm, MPI_Comm *barrier_comm,
MPI_Info *info, unsigned *sap_rank, unsigned *capt_rank)
{
H5FD_fphdf5_fapl_t *fa;
H5P_genplist_t *plist;
herr_t ret_value = SUCCEED;
FUNC_ENTER_API(H5Pget_fapl_fphdf5, FAIL)
H5TRACE6("e","i*Mc*Mc*Mi*Iu*Iu",fapl_id,comm,barrier_comm,info,sap_rank,
capt_rank);
if ((plist = H5P_object_verify(fapl_id, H5P_FILE_ACCESS)) == NULL)
HGOTO_ERROR(H5E_PLIST, H5E_BADTYPE, FAIL, "not a file access list")
if (H5P_get_driver(plist) != H5FD_FPHDF5)
HGOTO_ERROR(H5E_PLIST, H5E_BADVALUE, FAIL, "incorrect VFL driver")
if ((fa = H5P_get_driver_info(plist)) == NULL)
HGOTO_ERROR(H5E_PLIST, H5E_BADVALUE, FAIL, "bad VFL driver info")
if (comm)
*comm = fa->comm;
if (barrier_comm)
*barrier_comm = fa->barrier_comm;
if (info)
*info = fa->info;
if (sap_rank)
*sap_rank = fa->sap_rank;
if (capt_rank)
*capt_rank = fa->capt_rank;
done:
FUNC_LEAVE_API(ret_value)
}
/*-------------------------------------------------------------------------
* Function: H5FD_fphdf5_file_id
* Purpose: Returns the file ID for the file.
* Return: Success: File ID
* Failure: Doesn't fail
* Programmer: Bill Wendling
* 19. February 2003
* Modifications:
*-------------------------------------------------------------------------
*/
unsigned
H5FD_fphdf5_file_id(H5FD_t *_file)
{
H5FD_fphdf5_t *file = (H5FD_fphdf5_t*)_file;
unsigned ret_value;
FUNC_ENTER_NOAPI(H5FD_fphdf5_file_id, 0)
/* check args */
assert(file);
assert(file->pub.driver_id == H5FD_FPHDF5);
/* Set return value */
ret_value = file->file_id;
done:
FUNC_LEAVE_NOAPI(ret_value)
}
/*-------------------------------------------------------------------------
* Function: H5FD_fphdf5_is_sap
* Purpose: Asks the question "Is this process the SAP?".
* Return: Success: Non-zero if it is, 0 otherwise
* Failure: Doesn't fails
* Programmer: Bill Wendling
* 19. February 2003
* Modifications:
*-------------------------------------------------------------------------
*/
hbool_t
H5FD_fphdf5_is_sap(H5FD_t *_file)
{
H5FD_fphdf5_t *file = (H5FD_fphdf5_t*)_file;
hbool_t ret_value = FALSE;
FUNC_ENTER_NOAPI(H5FD_fphdf5_is_sap, FALSE)
/* check args */
assert(file);
assert(file->pub.driver_id == H5FD_FPHDF5);
/* Set return value */
ret_value = ((unsigned)file->mpi_rank == H5FP_sap_rank);
done:
FUNC_LEAVE_NOAPI(ret_value)
}
/*-------------------------------------------------------------------------
* Function: H5FD_fphdf5_is_captain
* Purpose: Asks the question "Is this process the captain?".
* Return: Success: Non-zero if it is, 0 otherwise
* Failure: Doesn't fails
* Programmer: Bill Wendling
* 19. February 2003
* Modifications:
*-------------------------------------------------------------------------
*/
hbool_t
H5FD_fphdf5_is_captain(H5FD_t *_file)
{
H5FD_fphdf5_t *file = (H5FD_fphdf5_t*)_file;
hbool_t ret_value = FALSE;
FUNC_ENTER_NOAPI(H5FD_fphdf5_is_sap, FALSE)
/* check args */
assert(file);
assert(file->pub.driver_id == H5FD_FPHDF5);
/* Set return value */
ret_value = ((unsigned)file->mpi_rank == H5FP_capt_rank);
done:
FUNC_LEAVE_NOAPI(ret_value)
}
/*-------------------------------------------------------------------------
* Function: H5FD_is_fphdf5_driver
* Purpose: Asks the question "Is this an FPHDF5 driver?".
* Return: Success: Non-zero if it is, 0 otherwise
* Failure: Doesn't fails
* Programmer: Bill Wendling
* 26. February 2003
* Modifications:
*-------------------------------------------------------------------------
*/
hbool_t
H5FD_is_fphdf5_driver(H5FD_t *_file)
{
H5FD_fphdf5_t *file = (H5FD_fphdf5_t*)_file;
hbool_t ret_value = FALSE;
FUNC_ENTER_NOAPI(H5FD_is_fphdf5_driver, FALSE)
/* check args */
assert(file);
/* Set return value */
ret_value = file->pub.driver_id == H5FD_FPHDF5;
done:
FUNC_LEAVE_NOAPI(ret_value)
}
/*-------------------------------------------------------------------------
* Function: H5Pset_dxpl_fphdf5
* Purpose: Set the data transfer property list DXPL_ID to use
* transfer mode XFER_MODE. The property list can then be
* used to control the I/O transfer mode during data I/O
* operations. The valid transfer modes are:
*
* H5FD_MPIO_INDEPENDENT:
* Use independent I/O access (the default).
*
* H5FD_MPIO_COLLECTIVE:
* Use collective I/O access.
* Return: Success: SUCCEED
* Failure: FAIL
* Programmer: Bill Wendling
* 10. February 2003
* Modifications:
*-------------------------------------------------------------------------
*/
herr_t
H5Pset_dxpl_fphdf5(hid_t dxpl_id, H5FD_mpio_xfer_t xfer_mode)
{
H5P_genplist_t *plist;
herr_t ret_value = SUCCEED;
FUNC_ENTER_API(H5Pset_dxpl_fphdf5, FAIL)
H5TRACE2("e","iDt",dxpl_id,xfer_mode);
if (dxpl_id == H5P_DEFAULT)
HGOTO_ERROR(H5E_PLIST, H5E_BADVALUE, FAIL, "can't set values in default property list")
/* Check arguments */
if ((plist = H5P_object_verify(dxpl_id, H5P_DATASET_XFER)) == NULL)
HGOTO_ERROR(H5E_PLIST, H5E_BADTYPE, FAIL, "not a dxpl")
if (xfer_mode != H5FD_MPIO_INDEPENDENT && xfer_mode != H5FD_MPIO_COLLECTIVE)
HGOTO_ERROR(H5E_ARGS, H5E_BADVALUE, FAIL, "incorrect xfer_mode")
/* Set the transfer mode */
if (H5P_set(plist, H5D_XFER_IO_XFER_MODE_NAME, &xfer_mode) < 0)
HGOTO_ERROR(H5E_PLIST, H5E_CANTSET, FAIL, "unable to set value")
/* Initialize driver-specific properties */
ret_value = H5P_set_driver(plist, H5FD_FPHDF5, NULL);
done:
FUNC_LEAVE_API(ret_value)
}
/*-------------------------------------------------------------------------
* Function: H5Pget_dxpl_fphdf5
* Purpose: Queries the transfer mode current set in the data
* transfer property list DXPL_ID. This is not collective.
* Return: Success: SUCCEED - with the transfer mode returned
* through the XFER_MODE argument if
* it is non-null.
* Failure: FAIL
* Programmer: Bill Wendling
* 10. February 2003
* Modifications:
*-------------------------------------------------------------------------
*/
herr_t
H5Pget_dxpl_fphdf5(hid_t dxpl_id, H5FD_mpio_xfer_t *xfer_mode)
{
H5P_genplist_t *plist;
herr_t ret_value = SUCCEED;
FUNC_ENTER_API(H5Pget_dxpl_fphdf5, FAIL)
H5TRACE2("e","i*Dt",dxpl_id,xfer_mode);
if ((plist = H5P_object_verify(dxpl_id, H5P_DATASET_XFER)) == NULL)
HGOTO_ERROR(H5E_PLIST, H5E_BADTYPE, FAIL, "not a dxpl")
if (H5P_get_driver(plist) != H5FD_FPHDF5)
HGOTO_ERROR(H5E_PLIST, H5E_BADVALUE, FAIL, "incorrect VFL driver")
/* Get the transfer mode */
if (xfer_mode)
if (H5P_get(plist, H5D_XFER_IO_XFER_MODE_NAME, xfer_mode) < 0)
HGOTO_ERROR(H5E_PLIST, H5E_CANTSET, FAIL, "unable to get value")
done:
FUNC_LEAVE_API(ret_value)
}
/*-------------------------------------------------------------------------
* Function: H5FD_fphdf5_fapl_get
* Purpose: Returns a file access property list which could be used
* to create another file the same as this one.
* Return: Success: Ptr to new file access property list with all
* fields copied from the file pointer.
* Failure: NULL
* Programmer: Bill Wendling
* 07. February 2003
* Modifications:
*-------------------------------------------------------------------------
*/
static void *
H5FD_fphdf5_fapl_get(H5FD_t *_file)
{
H5FD_fphdf5_t *file = (H5FD_fphdf5_t*)_file;
H5FD_fphdf5_fapl_t *fa = NULL;
void *ret_value;
FUNC_ENTER_NOAPI(H5FD_fphdf5_fapl_get, NULL)
/* check args */
assert(file);
assert(file->pub.driver_id == H5FD_FPHDF5);
if ((fa = H5MM_calloc(sizeof(H5FD_fphdf5_fapl_t))) == NULL)
HGOTO_ERROR(H5E_RESOURCE, H5E_NOSPACE, NULL, "memory allocation failed")
/* These should both be copied. --rpm, 1999-08-13 */
fa->comm = file->comm;
fa->barrier_comm = file->barrier_comm;
fa->info = file->info;
/* Set return value */
ret_value = fa;
done:
FUNC_LEAVE_NOAPI(ret_value)
}
/*-------------------------------------------------------------------------
* Function: H5FD_fphdf5_open
* Purpose: Opens a file with name NAME. The FLAGS are a bit field with
* purpose similar to the second argument of open(2) and
* which are defined in H5Fpublic.h. The file access
* property list FAPL_ID contains the properties driver
* properties and MAXADDR is the largest address which this
* file will be expected to access. This is collective.
* Return: Success: A new file pointer.
* Failure: NULL
* Programmer: Bill Wendling
* 05. February 2003
* Modifications:
*-------------------------------------------------------------------------
*/
static H5FD_t *
H5FD_fphdf5_open(const char *name, unsigned flags, hid_t fapl_id, haddr_t maxaddr)
{
H5FD_fphdf5_t *file = NULL;
MPI_File fh;
int mpi_amode;
int mpi_rank;
int mpi_size;
int mrc;
MPI_Offset size;
H5FD_fphdf5_fapl_t _fa;
const H5FD_fphdf5_fapl_t *fa = NULL;
H5P_genplist_t *plist;
unsigned long feature_flags;
hsize_t meta_block_size = 0;
hsize_t sdata_block_size = 0;
hsize_t threshold;
hsize_t alignment;
unsigned file_id;
unsigned req_id = 0;
H5FD_t *ret_value = NULL;
/* Flag to indicate that the file was successfully opened */
unsigned file_opened = FALSE;
FUNC_ENTER_NOAPI(H5FD_fphdf5_open, NULL)
/* check args */
assert(name);
/* Obtain a pointer to fphdf5-specific file access properties */
if ((plist = H5P_object_verify(fapl_id, H5P_FILE_ACCESS)) == NULL)
HGOTO_ERROR(H5E_ARGS, H5E_BADTYPE, NULL, "not a file access property list")
if (fapl_id == H5P_FILE_ACCESS_DEFAULT || H5P_get_driver(plist) != H5FD_FPHDF5) {
_fa.comm = MPI_COMM_SELF; /*default*/
_fa.barrier_comm = MPI_COMM_SELF; /*default*/
_fa.info = MPI_INFO_NULL; /*default*/
fa = &_fa;
} else {
fa = H5P_get_driver_info(plist);
assert(fa);
}
/*
* Convert HDF5 flags to MPI-IO flags. Some combinations are illegal;
* let MPI-IO figure it out
*/
mpi_amode = (flags & H5F_ACC_RDWR) ? MPI_MODE_RDWR : MPI_MODE_RDONLY;
if (flags & H5F_ACC_CREAT) mpi_amode |= MPI_MODE_CREATE;
if (flags & H5F_ACC_EXCL) mpi_amode |= MPI_MODE_EXCL;
/* OKAY: CAST DISCARDS CONST */
if ((mrc = MPI_File_open(H5FP_SAP_BARRIER_COMM, (char *)name, mpi_amode,
fa->info, &fh)) != MPI_SUCCESS)
HMPI_GOTO_ERROR(NULL, "MPI_File_open failed", mrc)
file_opened = TRUE;
if (H5P_get(plist, H5F_ACS_META_BLOCK_SIZE_NAME, &meta_block_size) < 0)
HGOTO_ERROR(H5E_PLIST, H5E_CANTGET, NULL, "can't get meta data block size")
if (H5P_get(plist, H5F_ACS_SDATA_BLOCK_SIZE_NAME, &sdata_block_size) < 0)
HGOTO_ERROR(H5E_PLIST, H5E_CANTGET, NULL, "can't get 'small data' block size")
if (H5P_get(plist, H5F_ACS_ALIGN_THRHD_NAME, &threshold) < 0)
HGOTO_ERROR(H5E_PLIST, H5E_CANTGET, NULL, "can't get alignment threshold")
if (H5P_get(plist, H5F_ACS_ALIGN_NAME, &alignment) < 0)
HGOTO_ERROR(H5E_PLIST, H5E_CANTGET, NULL, "can't get alignment")
/* Retrieve the VFL driver feature flags */
H5FD_fphdf5_query(NULL, &feature_flags); /* doesn't fail */
/* Inform the SAP that the file was opened */
if (H5FP_request_open(H5FP_OBJ_FILE, maxaddr, feature_flags,
meta_block_size, sdata_block_size, threshold,
alignment, &file_id, &req_id) == FAIL)
HGOTO_ERROR(H5E_FPHDF5, H5E_CANTOPENFILE, NULL, "can't inform SAP of file open")
/* Broadcast the file ID */
if ((mrc = MPI_Bcast(&file_id, 1, MPI_UNSIGNED,
(int)H5FP_capt_barrier_rank,
H5FP_SAP_BARRIER_COMM)) != MPI_SUCCESS)
HMPI_GOTO_ERROR(NULL, "MPI_Bcast failed", mrc)
/* Grab the rank of this process */
if ((mrc = MPI_Comm_rank(H5FP_SAP_COMM, &mpi_rank)) != MPI_SUCCESS)
HMPI_GOTO_ERROR(NULL, "MPI_Comm_rank failed", mrc)
/* The captain rank will get the filesize and broadcast it. */
if ((unsigned)mpi_rank == H5FP_capt_rank)
/* Get current file size */
if ((mrc = MPI_File_get_size(fh, &size)) != MPI_SUCCESS)
HMPI_GOTO_ERROR(NULL, "MPI_File_get_size failed", mrc)
/* Broadcast file size */
if ((mrc = MPI_Bcast(&size, sizeof(MPI_Offset), MPI_BYTE,
(int)H5FP_capt_barrier_rank,
H5FP_SAP_BARRIER_COMM)) != MPI_SUCCESS)
HMPI_GOTO_ERROR(NULL, "MPI_Bcast failed", mrc)
/* Only if size > 0, truncate the file - if requested */
if (size && (flags & H5F_ACC_TRUNC)) {
if ((mrc = MPI_File_set_size(fh, (MPI_Offset)0)) != MPI_SUCCESS)
HMPI_GOTO_ERROR(NULL, "MPI_File_set_size (file truncation) failed", mrc)
/* Don't let any proc return until all have truncated the file. */
if ((mrc = MPI_Barrier(H5FP_SAP_BARRIER_COMM)) != MPI_SUCCESS)
HMPI_GOTO_ERROR(NULL, "MPI_Barrier failed", mrc)
size = 0;
}
/* Grab the size of the communicator */
if ((mrc = MPI_Comm_size(H5FP_SAP_COMM, &mpi_size)) != MPI_SUCCESS)
HMPI_GOTO_ERROR(NULL, "MPI_Comm_size failed", mrc)
/* Build the return value and initialize it */
if ((file = H5MM_calloc(sizeof(H5FD_fphdf5_t))) == NULL)
HGOTO_ERROR(H5E_RESOURCE, H5E_NOSPACE, NULL, "memory allocation failed")
file->file_id = file_id;
file->f = fh;
file->comm = fa->comm;
file->info = fa->info;
file->mpi_rank = mpi_rank;
file->mpi_size = mpi_size;
file->eof = H5FD_mpi_MPIOff_to_haddr(size);
/* Set return value */
ret_value = (H5FD_t *)file;
done:
if (!ret_value && file_opened)
MPI_File_close(&fh);
FUNC_LEAVE_NOAPI(ret_value)
}
/*-------------------------------------------------------------------------
* Function: H5FD_fphdf5_close
* Purpose: Closes a file. This is collective.
* Return: Success: SUCCEED
* Failure: FAIL
* Programmer: Bill Wendling
* 07. February 2003
* Modifications:
*-------------------------------------------------------------------------
*/
static herr_t
H5FD_fphdf5_close(H5FD_t *_file)
{
H5FD_fphdf5_t *file = (H5FD_fphdf5_t *)_file;
unsigned req_id = 0;
H5FP_status_t status = H5FP_STATUS_OK;
int mrc;
herr_t ret_value = SUCCEED;
FUNC_ENTER_NOAPI(H5FD_fphdf5_close, FAIL)
/* check args */
assert(file);
assert(file->pub.driver_id == H5FD_FPHDF5);
/* Tell the SAP that we're closing the file... */
if (H5FP_request_close(_file, file->file_id, &req_id, &status) == FAIL)
HGOTO_ERROR(H5E_IO, H5E_CANTCLOSEFILE, FAIL, "can't inform SAP of file close")
/* MPI_File_close sets argument to MPI_FILE_NULL */
if ((mrc = MPI_File_close(&file->f)) != MPI_SUCCESS)
HMPI_GOTO_ERROR(FAIL, "MPI_File_close failed", mrc)
/* Clean up other stuff */
H5MM_xfree(file);
done:
FUNC_LEAVE_NOAPI(ret_value)
}
/*-------------------------------------------------------------------------
* Function: H5FD_fphdf5_query
* Purpose: Set the flags that this VFL driver is capable of
* supporting. (listed in H5FDpublic.h)
* Return: Success: SUCCEED
* Failure: Doesn't fail.
* Programmer: Bill Wendling
* 07. February 2003
* Modifications:
*-------------------------------------------------------------------------
*/
static herr_t
H5FD_fphdf5_query(const H5FD_t UNUSED *_file, unsigned long *flags /* out */)
{
herr_t ret_value = SUCCEED;
FUNC_ENTER_NOAPI(H5FD_fphdf5_query, FAIL)
/* Set the VFL feature flags that this driver supports */
if (flags) {
*flags = 0;
/* OK to aggregate metadata allocations */
*flags |= H5FD_FEAT_AGGREGATE_METADATA;
/*
* Distinguish between updating the metadata accumulator on
* writes and reads. This is particularly (perhaps only, even)
* important for MPI-I/O where we guarantee that writes are
* collective, but reads may not be. If we were to allow the
* metadata accumulator to be written during a read operation,
* the application would hang.
*/
#if 0
/*
* FIXME: For now, having metadata accumulate causes problems for
* the SAP when it goes to allocate data (oddly enough, an
* allocation can result in a call to H5FD_free...which can
* result in a call to H5FD_write...which needs a data xfer
* property list...but only when metadata accumulation is turned
* on...go figure). Turn it off for now. -- BW 02/19/2003
*/
/* OK to accumulate metadata for faster writes */
*flags |= H5FD_FEAT_ACCUMULATE_METADATA_WRITE;
#endif /* 0 */
/* OK to aggregate "small" raw data allocations */
*flags |= H5FD_FEAT_AGGREGATE_SMALLDATA;
}
done:
FUNC_LEAVE_NOAPI(ret_value)
}
/*-------------------------------------------------------------------------
* Function: H5FD_fphdf5_get_eoa
* Purpose: Gets the end-of-address marker for the file. The EOA
* marker is the first address past the last byte allocated
* in the format address space.
* Return: Success: The end-of-address marker.
* Failure: HADDR_UNDEF
* Programmer: Bill Wendling
* 07. February 2003
* Modifications:
*-------------------------------------------------------------------------
*/
static haddr_t
H5FD_fphdf5_get_eoa(H5FD_t *_file)
{
H5FD_fphdf5_t *file = (H5FD_fphdf5_t *)_file;
unsigned req_id = 0;
H5FP_status_t status = H5FP_STATUS_OK;
haddr_t ret_value;
FUNC_ENTER_NOAPI(H5FD_fphdf5_get_eoa, HADDR_UNDEF)
/* check args */
assert(file);
assert(file->pub.driver_id == H5FD_FPHDF5);
/* The SAP's eoa field is correct */
if (!H5FD_fphdf5_is_sap(_file))
/* Retrieve the EOA information */
if (H5FP_request_get_eoa(_file, &file->eoa, &req_id, &status))
HGOTO_ERROR(H5E_FPHDF5, H5E_CANTRECV, HADDR_UNDEF, "can't retrieve eoa information")
/* Set return value */
ret_value = file->eoa;
done:
FUNC_LEAVE_NOAPI(ret_value)
}
/*-------------------------------------------------------------------------
* Function: H5FD_fphdf5_set_eoa
* Purpose: Set the end-of-address marker for the file. This function
* is called shortly after an existing HDF5 file is opened
* in order to tell the driver where the end of the HDF5
* data is located.
* Return: Success: SUCCEED
* Failure: FAIL
* Programmer: Bill Wendling
* 06. February 2003
* Modifications:
*-------------------------------------------------------------------------
*/
static herr_t
H5FD_fphdf5_set_eoa(H5FD_t *_file, haddr_t addr)
{
H5FD_fphdf5_t *file = (H5FD_fphdf5_t *)_file;
unsigned req_id = 0;
H5FP_status_t status = H5FP_STATUS_OK;
herr_t ret_value = SUCCEED;
FUNC_ENTER_NOAPI(H5FD_fphdf5_set_eoa, FAIL)
/* check args */
assert(file);
assert(file->pub.driver_id == H5FD_FPHDF5);
/* The SAP's eoa field is correct */
if (!H5FD_fphdf5_is_sap(_file))
/* Retrieve the EOA information */
if (H5FP_request_set_eoa(_file, addr, &req_id, &status))
HGOTO_ERROR(H5E_FPHDF5, H5E_CANTRECV, FAIL, "can't retrieve eoa information")
file->eoa = addr;
done:
FUNC_LEAVE_NOAPI(ret_value)
}
/*-------------------------------------------------------------------------
* Function: H5FD_fphdf5_get_eof
* Purpose: Gets the end-of-file marker for the file. The EOF marker
* is the real size of the file.
*
* The FPHDF5 driver doesn't bother keeping this field updated
* since that's a relatively expensive operation.
* Fortunately the library only needs the EOF just after the
* file is opened in order to determine whether the file is
* empty, truncated, or okay. Therefore, any MPIO I/O
* function will set its value to HADDR_UNDEF which is the
* error return value of this function.
* Return: Success: The end-of-address marker
* Failure: HADDR_UNDEF
* Programmer: Bill Wendling
* 06. February 2003
* Modifications:
*-------------------------------------------------------------------------
*/
static haddr_t
H5FD_fphdf5_get_eof(H5FD_t *_file)
{
H5FD_fphdf5_t *file = (H5FD_fphdf5_t*)_file;
haddr_t ret_value;
FUNC_ENTER_NOAPI(H5FD_fphdf5_get_eof, HADDR_UNDEF)
/* check args */
assert(file);
assert(file->pub.driver_id == H5FD_FPHDF5);
/* Set return value */
ret_value = file->eof;
done:
FUNC_LEAVE_NOAPI(ret_value)
}
/*-------------------------------------------------------------------------
* Function: H5FD_fphdf5_get_handle
* Purpose: Returns the file handle of MPIO file driver.
* Returns: Success: SUCCEED
* Failure: FAIL
* Programmer: Bill Wendling
* 06. February 2003
* Modifications:
*-------------------------------------------------------------------------
*/
static herr_t
H5FD_fphdf5_get_handle(H5FD_t *_file, hid_t UNUSED fapl, void** file_handle)
{
H5FD_fphdf5_t *file = (H5FD_fphdf5_t *)_file;
herr_t ret_value = SUCCEED;
FUNC_ENTER_NOAPI(H5FD_fphdf5_get_handle, FAIL)
/* check args */
assert(file);
if (!file_handle)
HGOTO_ERROR(H5E_ARGS, H5E_BADVALUE, FAIL, "file handle not valid")
*file_handle = &file->f;
done:
FUNC_LEAVE_NOAPI(ret_value)
}
/*-------------------------------------------------------------------------
* Function: H5FD_fphdf5_read
* Purpose: Reads SIZE bytes of data from FILE beginning at address
* ADDR into buffer BUF according to data transfer
* properties in DXPL_ID using potentially complex file and
* buffer types to effect the transfer.
*
* Reading past the end of the MPI file returns zeros
* instead of failing. MPI is able to coalesce requests
* from different processes (collective or independent).
* Return: Success: SUCCEED - Result is stored in caller-supplied
* buffer BUF
* Failure: FAIL - Contents of buffer BUF are undefined
* Programmer: Bill Wendling
* 10. February 2003
* Modifications:
*-------------------------------------------------------------------------
*/
static herr_t
H5FD_fphdf5_read(H5FD_t *_file, H5FD_mem_t mem_type, hid_t dxpl_id,
haddr_t addr, size_t size, void *buf)
{
H5FD_fphdf5_t *file = (H5FD_fphdf5_t*)_file;
MPI_Status status;
MPI_Offset mpi_off;
int size_i; /* Integer copy of 'size' to read */
int bytes_read; /* Number of bytes read in */
int mrc; /* MPI return code */
herr_t ret_value = SUCCEED;
FUNC_ENTER_NOAPI(H5FD_fphdf5_read, FAIL)
/* check args */
assert(file);
assert(file->pub.driver_id == H5FD_FPHDF5);
assert(buf);
/* Make certain we have the correct type of property list */
assert(H5I_get_type(dxpl_id) == H5I_GENPROP_LST);
assert(H5P_isa_class(dxpl_id, H5P_DATASET_XFER) == TRUE);
/* Portably initialize MPI status variable */
HDmemset(&status, 0, sizeof(MPI_Status));
/* Some numeric conversions */
if (H5FD_mpi_haddr_to_MPIOff(addr, &mpi_off) < 0)
HGOTO_ERROR(H5E_INTERNAL, H5E_BADRANGE, FAIL, "can't convert from haddr_t to MPI offset")
size_i = (int)size;
if ((hsize_t)size_i != size)
HGOTO_ERROR(H5E_INTERNAL, H5E_BADRANGE, FAIL, "can't convert from size_t to int")
/* Only do real MPI stuff for raw data transfers */
if(mem_type==H5FD_MEM_DRAW) {
H5P_genplist_t *plist;
H5FD_mpio_xfer_t xfer_mode; /* I/O tranfer mode */
MPI_Datatype buf_type; /* MPI datatype for I/O operation */
int n;
int type_size; /* MPI datatype used for I/O's size */
int io_size; /* Actual number of bytes requested */
unsigned use_view_this_time = 0; /* Flag to indicate we are using an MPI view */
/* Obtain the data transfer properties */
if ((plist = H5I_object(dxpl_id)) == NULL)
HGOTO_ERROR(H5E_ARGS, H5E_BADTYPE, FAIL, "not a file access property list")
xfer_mode = H5P_peek_unsigned(plist, H5D_XFER_IO_XFER_MODE_NAME);
/*
* Set up for a fancy xfer using complex types, or single byte block.
* We wouldn't need to rely on the use_view field if MPI semantics
* allowed us to test that btype == ftype == MPI_BYTE (or even
* MPI_TYPE_NULL, which could mean "use MPI_BYTE" by convention).
*/
if (xfer_mode == H5FD_MPIO_COLLECTIVE) {
MPI_Datatype file_type;
/* Remember that views are used */
use_view_this_time = TRUE;
/* Prepare for a full-blown xfer using btype, ftype, and disp */
if (H5P_get(plist, H5FD_MPI_XFER_MEM_MPI_TYPE_NAME, &buf_type) < 0)
HGOTO_ERROR(H5E_PLIST, H5E_CANTGET, FAIL, "can't get MPI-I/O type property")
if (H5P_get(plist, H5FD_MPI_XFER_FILE_MPI_TYPE_NAME, &file_type) < 0)
HGOTO_ERROR(H5E_PLIST, H5E_CANTGET, FAIL, "can't get MPI-I/O type property")
/* Set the file view when we are using MPI derived types */
if ((mrc = MPI_File_set_view(file->f, mpi_off, MPI_BYTE,
file_type, H5FD_mpi_native_g,
file->info)) != MPI_SUCCESS)
HMPI_GOTO_ERROR(FAIL, "MPI_File_set_view failed", mrc)
/*
* When using types, use the address as the displacement for
* MPI_File_set_view and reset the address for the read to zero
*/
mpi_off = 0;
/* Read the data. */
if ((mrc = MPI_File_read_at_all(file->f, mpi_off, buf, size_i,
buf_type, &status)) != MPI_SUCCESS)
HMPI_GOTO_ERROR(FAIL, "MPI_File_read_at_all failed", mrc)
/*
* Reset the file view when we used MPI derived types
*/
if ((mrc = MPI_File_set_view(file->f, (MPI_Offset)0, MPI_BYTE, MPI_BYTE,
H5FD_mpi_native_g, file->info)) != MPI_SUCCESS)
HMPI_GOTO_ERROR(FAIL, "MPI_File_set_view failed", mrc)
} else {
/*
* Prepare for a simple xfer of a contiguous block of bytes. The
* btype, ftype, and disp fields are not used.
*/
buf_type = MPI_BYTE;
/* Read the data. */
if ((mrc = MPI_File_read_at(file->f, mpi_off, buf, size_i,
buf_type, &status)) != MPI_SUCCESS)
HMPI_GOTO_ERROR(FAIL, "MPI_File_read_at failed", mrc)
}
/* How many bytes were actually read? */
/* [This works because the "basic elements" we use for all our MPI derived
* types are MPI_BYTE. We should be using the 'buf_type' for the MPI
* datatype in this call though... (We aren't because using it causes
* the LANL "qsc" machine to dump core - 12/19/03) - QAK]
*/
if (MPI_SUCCESS != (mrc=MPI_Get_elements(&status, MPI_BYTE, &bytes_read)))
HMPI_GOTO_ERROR(FAIL, "MPI_Get_elements failed", mrc)
/* Get the type's size */
if (MPI_SUCCESS != (mrc=MPI_Type_size(buf_type,&type_size)))
HMPI_GOTO_ERROR(FAIL, "MPI_Type_size failed", mrc)
/* Compute the actual number of bytes requested */
io_size=type_size*size_i;
/* Check for read failure */
if (bytes_read<0 || bytes_read>io_size)
HGOTO_ERROR(H5E_IO, H5E_READERROR, FAIL, "file read failed")
/*
* This gives us zeroes beyond end of physical MPI file.
*/
if ((n=(io_size-bytes_read)) > 0)
HDmemset((char*)buf+bytes_read, 0, (size_t)n);
} /* end if */
else {
/*
* This is metadata - we want to try to read it from the SAP
* first.
*/
unsigned req_id = 0;
H5FP_status_t sap_status = H5FP_STATUS_OK;
if (H5FP_request_read_metadata(_file, file->file_id, dxpl_id, mem_type,
addr, size, (uint8_t**)&buf,
&req_id, &sap_status) != SUCCEED) {
/* FIXME: The read failed, for some reason */
HDfprintf(stderr, "%s:%d: Metadata cache read failed!\n", FUNC, __LINE__);
}
if (sap_status == H5FP_STATUS_OK) {
/* WAH-HOO! We've found it! We can leave now */
} else if (sap_status == H5FP_STATUS_MDATA_NOT_CACHED) {
/* Metadata isn't cached, so grab it from the file */
if ((mrc = MPI_File_read_at(file->f, mpi_off, buf, size_i,
MPI_BYTE, &status)) != MPI_SUCCESS)
HMPI_GOTO_ERROR(FAIL, "MPI_File_read_at failed", mrc)
if ((mrc = MPI_Get_count(&status, MPI_BYTE, &bytes_read)) != MPI_SUCCESS)
HMPI_GOTO_ERROR(FAIL, "MPI_Get_count failed", mrc)
} else {
/* FIXME: something bad happened */
HDfprintf(stderr, "%s:%d: Metadata cache read failed!\n", FUNC, __LINE__);
}
} /* end else */
done:
FUNC_LEAVE_NOAPI(ret_value)
}
/*-------------------------------------------------------------------------
* Function: H5FD_fphdf5_write
* Purpose: Writes SIZE bytes of data to FILE beginning at address
* ADDR from buffer BUF according to data transfer
* properties in DXPL_ID using potentially complex file and
* buffer types to effect the transfer.
*
* MPI is able to coalesce requests from different processes
* (collective and independent).
* Return: Success: SUCCEED - USE_TYPES and OLD_USE_TYPES in the
* access params are altered.
* Failure: FAIL - USE_TYPES and OLD_USE_TYPES in the
* access params may be altered.
* Programmer: Bill Wendling
* 10. February 2003
* Modifications:
*-------------------------------------------------------------------------
*/
static herr_t
H5FD_fphdf5_write(H5FD_t *_file, H5FD_mem_t mem_type, hid_t dxpl_id,
haddr_t addr, size_t size, const void *buf)
{
H5FD_fphdf5_t *file = (H5FD_fphdf5_t*)_file;
int size_i;
unsigned dumping = 0;
H5P_genplist_t *plist;
herr_t ret_value;
FUNC_ENTER_NOAPI(H5FD_fphdf5_write, FAIL)
/* check args */
assert(file);
assert(file->pub.driver_id == H5FD_FPHDF5);
assert(buf);
/* Make certain we have the correct type of property list */
assert(H5I_get_type(dxpl_id) == H5I_GENPROP_LST);
assert(H5P_isa_class(dxpl_id, H5P_DATASET_XFER) == TRUE);
/* some numeric conversions */
size_i = (int)size;
if ((hsize_t)size_i != size)
HGOTO_ERROR(H5E_INTERNAL, H5E_BADRANGE, FAIL, "can't convert from size to size_i")
/* Obtain the data transfer properties */
if ((plist = H5I_object(dxpl_id)) == NULL)
HGOTO_ERROR(H5E_ARGS, H5E_BADTYPE, FAIL, "not a file access property list")
if (H5P_exist_plist(plist, H5FD_FPHDF5_XFER_DUMPING_METADATA) > 0)
if (H5P_get(plist, H5FD_FPHDF5_XFER_DUMPING_METADATA, &dumping) < 0)
HGOTO_ERROR(H5E_PLIST, H5E_CANTGET, FAIL, "can't get MPI-I/O type property")
/*
* If metadata, write to the metadata cache but only if we're not
* dumping the data from the SAP...
*/
if (mem_type != H5FD_MEM_DRAW && !dumping) {
unsigned req_id = 0;
H5FP_status_t sap_status = H5FP_STATUS_OK;
if (H5FP_request_write_metadata(_file, file->file_id, dxpl_id, mem_type,
addr, size_i, buf, &req_id,
&sap_status) != SUCCEED) {
/* FIXME: Couldn't write metadata. This is bad... */
HDfprintf(stderr, "%s:%d: Couldn't write metadata to SAP (%d)\n",
FUNC, __LINE__, sap_status);
}
switch (sap_status) {
case H5FP_STATUS_OK:
/* WAH-HOO! We've written it! We can leave now */
/* Forget the EOF value (see H5FD_fphdf5_get_eof()) */
file->eof = HADDR_UNDEF;
HGOTO_DONE(SUCCEED)
case H5FP_STATUS_FILE_CLOSING:
HGOTO_DONE(SUCCEED)
case H5FP_STATUS_DUMPING_FAILED:
case H5FP_STATUS_OOM:
case H5FP_STATUS_BAD_FILE_ID:
default:
/* FIXME: Something bad happened */
HDfprintf(stderr, "%s: Couldn't write metadata to SAP (%d)\n",
FUNC, sap_status);
break;
}
}
/* FIXME: Should I check this return value or just pass it on out? */
ret_value = H5FD_fphdf5_write_real(_file, mem_type, plist, addr, size_i, buf);
done:
FUNC_LEAVE_NOAPI(ret_value)
}
/*-------------------------------------------------------------------------
* Function: H5FD_fphdf5_write_real
* Purpose: Do the actual writing to a file. Split apart from the
* H5FD_fphdf5_write call since I need to write things
* directly if the SAP is dumping data to me.
* Return: Success: SUCCEED
* Failure: FAIL
* Programmer: Bill Wendling
* 12. February 2003
* Modifications:
*-------------------------------------------------------------------------
*/
herr_t
H5FD_fphdf5_write_real(H5FD_t *_file, H5FD_mem_t mem_type, H5P_genplist_t *plist,
haddr_t addr, int size, const void *buf)
{
H5FD_fphdf5_t *file = (H5FD_fphdf5_t*)_file;
MPI_Status status;
MPI_Datatype buf_type = MPI_BYTE;
MPI_Offset mpi_off;
int mrc;
int bytes_written;
int type_size; /* MPI datatype used for I/O's size */
int io_size; /* Actual number of bytes requested */
unsigned use_view_this_time = 0;
herr_t ret_value = SUCCEED;
FUNC_ENTER_NOAPI(H5FD_fphdf5_write_real, FAIL)
/* check args */
assert(file);
assert(file->pub.driver_id == H5FD_FPHDF5);
assert(plist);
assert(buf);
/* Portably initialize MPI status variable */
HDmemset(&status, 0, sizeof(MPI_Status));
/* some numeric conversions */
if (H5FD_mpi_haddr_to_MPIOff(addr, &mpi_off) < 0)
HGOTO_ERROR(H5E_INTERNAL, H5E_BADRANGE, FAIL, "can't convert from haddr to MPI off")
/* Only check for fancy transfers with raw data I/O */
if (mem_type == H5FD_MEM_DRAW) {
H5FD_mpio_xfer_t xfer_mode; /* I/O tranfer mode */
/* Obtain the data transfer properties */
xfer_mode = (H5FD_mpio_xfer_t)H5P_peek_unsigned(plist, H5D_XFER_IO_XFER_MODE_NAME);
/*
* Set up for a fancy xfer using complex types, or single byte block.
* We wouldn't need to rely on the use_view field if MPI semantics
* allowed us to test that btype == ftype == MPI_BYTE (or even
* MPI_TYPE_NULL, which could mean "use MPI_BYTE" by convention).
*/
if (xfer_mode == H5FD_MPIO_COLLECTIVE) {
MPI_Datatype file_type;
/* Remember that views are used */
use_view_this_time = TRUE;
/* Prepare for a full-blown xfer using btype, ftype, and disp */
if (H5P_get(plist, H5FD_MPI_XFER_MEM_MPI_TYPE_NAME, &buf_type) < 0)
HGOTO_ERROR(H5E_PLIST, H5E_CANTGET, FAIL, "can't get MPI-I/O type property")
if (H5P_get(plist, H5FD_MPI_XFER_FILE_MPI_TYPE_NAME, &file_type) < 0)
HGOTO_ERROR(H5E_PLIST, H5E_CANTGET, FAIL, "can't get MPI-I/O type property")
/* Set the file view when we are using MPI derived types */
if ((mrc = MPI_File_set_view(file->f, mpi_off, MPI_BYTE,
file_type, H5FD_mpi_native_g,
file->info)) != MPI_SUCCESS)
HMPI_GOTO_ERROR(FAIL, "MPI_File_set_view failed", mrc)
/*
* When using types, use the address as the displacement for
* MPI_File_set_view and reset the address for the read to zero
*/
mpi_off = 0;
}
} /* end if */
/* Write the data. */
if (use_view_this_time) {
/*OKAY: CAST DISCARDS CONST QUALIFIER*/
if ((mrc = MPI_File_write_at_all(file->f, mpi_off, (void*)buf,
size, buf_type, &status)) != MPI_SUCCESS)
HMPI_GOTO_ERROR(FAIL, "MPI_File_write_at_all failed", mrc)
/* Reset the file view when we used MPI derived types */
if ((mrc = MPI_File_set_view(file->f, (MPI_Offset)0, MPI_BYTE, MPI_BYTE,
H5FD_mpi_native_g, file->info)) != MPI_SUCCESS)
HMPI_GOTO_ERROR(FAIL, "MPI_File_set_view failed", mrc)
} /* end if */
else {
/*OKAY: CAST DISCARDS CONST QUALIFIER*/
if ((mrc = MPI_File_write_at(file->f, mpi_off, (void*)buf,
size, buf_type, &status)) != MPI_SUCCESS)
HMPI_GOTO_ERROR(FAIL, "MPI_File_write_at failed", mrc)
}
/* How many bytes were actually written? */
/* [This works because the "basic elements" we use for all our MPI derived
* types are MPI_BYTE. We should be using the 'buf_type' for the MPI
* datatype in this call though... (We aren't because using it causes
* the LANL "qsc" machine to dump core - 12/19/03) - QAK]
*/
if (MPI_SUCCESS != (mrc=MPI_Get_elements(&status, MPI_BYTE, &bytes_written)))
HMPI_GOTO_ERROR(FAIL, "MPI_Get_elements failed", mrc)
/* Get the type's size */
if (MPI_SUCCESS != (mrc=MPI_Type_size(buf_type,&type_size)))
HMPI_GOTO_ERROR(FAIL, "MPI_Type_size failed", mrc)
/* Compute the actual number of bytes requested */
io_size=type_size*size;
/* Check for write failure */
if (bytes_written!=io_size)
HGOTO_ERROR(H5E_IO, H5E_WRITEERROR, FAIL, "file write failed")
/* Forget the EOF value (see H5FD_fphdf5_get_eof()) */
file->eof = HADDR_UNDEF;
done:
FUNC_LEAVE_NOAPI(ret_value)
}
/*-------------------------------------------------------------------------
* Function: H5FD_fphdf5_flush
* Purpose: Makes sure that all data is on disk. This is collective.
* Return: Success: SUCCEED
* Failure: FAIL
* Programmer: Bill Wendling
* 12. February 2003
* Modifications:
*-------------------------------------------------------------------------
*/
static herr_t
H5FD_fphdf5_flush(H5FD_t *_file, hid_t dxpl_id, unsigned closing)
{
H5FD_fphdf5_t *file = (H5FD_fphdf5_t*)_file;
MPI_Offset mpi_off;
int mrc;
unsigned req_id = 0;
H5FP_status_t status = H5FP_STATUS_OK;
herr_t ret_value = SUCCEED;
FUNC_ENTER_NOAPI(H5FD_fphdf5_flush, FAIL)
/* check args */
assert(file);
assert(file->pub.driver_id == H5FD_FPHDF5);
/*
* Extend the file to make sure it's large enough, then sync.
* Unfortunately, keeping track of EOF is an expensive operation, so
* we can't just check whether EOF<EOA like with other drivers.
* Therefore we'll just read the byte at EOA-1 and then write it
* back.
*/
if (file->eoa > file->last_eoa) {
if (H5FD_mpi_haddr_to_MPIOff(file->eoa, &mpi_off) < 0)
HGOTO_ERROR(H5E_INTERNAL, H5E_BADRANGE, FAIL, "cannot convert from haddr_t to MPI_Offset")
/* Extend the file's size */
if ((mrc = MPI_File_set_size(file->f, mpi_off)) != MPI_SUCCESS)
HMPI_GOTO_ERROR(FAIL, "MPI_File_set_size failed", mrc)
/*
* Don't let any proc return until all have extended the file.
* (Prevents race condition where some processes go ahead and
* write more data to the file before all the processes have
* finished making it the shorter length, potentially truncating
* the file and dropping the new data written)
*/
if ((mrc = MPI_Barrier(H5FP_SAP_BARRIER_COMM)) != MPI_SUCCESS)
HMPI_GOTO_ERROR(FAIL, "MPI_Barrier failed", mrc)
/* Update the 'last' eoa value */
file->last_eoa = file->eoa;
}
/* Only the captain process needs to flush the metadata. */
if (H5FD_fphdf5_is_captain(_file)) {
if (H5FP_request_flush_metadata(_file, file->file_id, dxpl_id,
&req_id, &status) != SUCCEED)
HGOTO_ERROR(H5E_FPHDF5, H5E_CANTFLUSH, FAIL, "can't flush metadata")
/* Only sync the file if we are not going to immediately close it */
if (!closing)
if ((mrc = MPI_File_sync(file->f)) != MPI_SUCCESS)
HMPI_GOTO_ERROR(FAIL, "MPI_File_sync failed", mrc)
}
done:
FUNC_LEAVE_NOAPI(ret_value)
}
/*-------------------------------------------------------------------------
* Function: H5FD_fphdf5_mpi_rank
* Purpose: Returns the MPI rank for a process
* Return: Success: MPI rank
* Failure: Doesn't fail
* Programmer: Bill Wendling
* 30. January 2003
* Modifications:
*-------------------------------------------------------------------------
*/
static int
H5FD_fphdf5_mpi_rank(const H5FD_t *_file)
{
const H5FD_fphdf5_t *file = (const H5FD_fphdf5_t*)_file;
int ret_value;
FUNC_ENTER_NOAPI(H5FD_fphdf5_mpi_rank, FAIL)
/* check args */
assert(file);
assert(file->pub.driver_id == H5FD_FPHDF5);
/* Set return value */
ret_value = file->mpi_rank;
done:
FUNC_LEAVE_NOAPI(ret_value)
}
/*-------------------------------------------------------------------------
* Function: H5FD_fphdf5_mpi_size
* Purpose: Returns the number of MPI processes
* Return: Success: Number of MPI processes
* Failure: Doesn't fail
* Programmer: Bill Wendling
* 30. January 2003
* Modifications:
*-------------------------------------------------------------------------
*/
static int
H5FD_fphdf5_mpi_size(const H5FD_t *_file)
{
const H5FD_fphdf5_t *file = (const H5FD_fphdf5_t*)_file;
int ret_value;
FUNC_ENTER_NOAPI(H5FD_fphdf5_mpi_size, FAIL)
/* check args */
assert(file);
assert(file->pub.driver_id == H5FD_FPHDF5);
/* Set return value */
ret_value = file->mpi_size;
done:
FUNC_LEAVE_NOAPI(ret_value)
}
/*-------------------------------------------------------------------------
* Function: H5FD_fphdf5_barrier_communicator
* Purpose: Returns the MPI communicator for the file that can be
* used in an MPI_Barrier() statement for the client
* processes.
* Return: Success: The barrier communicator
* Failure: NULL
* Programmer: Bill Wendling
* 10. February 2003
* Modifications:
*-------------------------------------------------------------------------
*/
static MPI_Comm
H5FD_fphdf5_barrier_communicator(const H5FD_t *_file)
{
const H5FD_fphdf5_t *file = (const H5FD_fphdf5_t*)_file;
MPI_Comm ret_value;
FUNC_ENTER_NOAPI(H5FD_fphdf5_barrier_communicator, MPI_COMM_NULL)
/* check args */
assert(file);
assert(file->pub.driver_id == H5FD_FPHDF5);
/* Set return value */
ret_value = file->barrier_comm;
done:
FUNC_LEAVE_NOAPI(ret_value)
}
#endif /* H5_HAVE_FPHDF5 */
|