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
|
/* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
* 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. *
* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
/*
* Programmer: Quincey Koziol
* Friday, January 30, 2004
*
* Purpose: Common routines for all MPI-based VFL drivers.
*
*/
#include "H5private.h" /* Generic Functions */
#include "H5CXprivate.h" /* API Contexts */
#include "H5Eprivate.h" /* Error handling */
#include "H5Fprivate.h" /* File access */
#include "H5FDprivate.h" /* File drivers */
#include "H5FDmpi.h" /* Common MPI file driver */
#include "H5Pprivate.h" /* Property lists */
#ifdef H5_HAVE_PARALLEL
/*-------------------------------------------------------------------------
* Function: H5FD_mpi_get_rank
*
* Purpose: Retrieves the rank of an MPI process.
*
* Return: Success: The rank (non-negative)
*
* Failure: Negative
*
* Programmer: Quincey Koziol
* Friday, January 30, 2004
*
*-------------------------------------------------------------------------
*/
int
H5FD_mpi_get_rank(const H5FD_t *file)
{
const H5FD_class_mpi_t *cls;
int ret_value;
FUNC_ENTER_NOAPI(FAIL)
HDassert(file);
cls = (const H5FD_class_mpi_t *)(file->cls);
HDassert(cls);
HDassert(cls->get_rank); /* All MPI drivers must implement this */
/* Dispatch to driver */
if ((ret_value = (cls->get_rank)(file)) < 0)
HGOTO_ERROR(H5E_VFL, H5E_CANTGET, FAIL, "driver get_rank request failed")
done:
FUNC_LEAVE_NOAPI(ret_value)
} /* end H5FD_mpi_get_rank() */
/*-------------------------------------------------------------------------
* Function: H5FD_mpi_get_size
*
* Purpose: Retrieves the size of the communicator used for the file
*
* Return: Success: The communicator size (non-negative)
*
* Failure: Negative
*
* Programmer: Quincey Koziol
* Friday, January 30, 2004
*
*-------------------------------------------------------------------------
*/
int
H5FD_mpi_get_size(const H5FD_t *file)
{
const H5FD_class_mpi_t *cls;
int ret_value;
FUNC_ENTER_NOAPI(FAIL)
HDassert(file);
cls = (const H5FD_class_mpi_t *)(file->cls);
HDassert(cls);
HDassert(cls->get_size); /* All MPI drivers must implement this */
/* Dispatch to driver */
if ((ret_value = (cls->get_size)(file)) < 0)
HGOTO_ERROR(H5E_VFL, H5E_CANTGET, FAIL, "driver get_size request failed")
done:
FUNC_LEAVE_NOAPI(ret_value)
} /* end H5FD_mpi_get_size() */
/*-------------------------------------------------------------------------
* Function: H5FD_mpi_get_comm
*
* Purpose: Retrieves the file's communicator
*
* Return: Success: The communicator (non-negative)
*
* Failure: Negative
*
* Programmer: Quincey Koziol
* Friday, January 30, 2004
*
*-------------------------------------------------------------------------
*/
MPI_Comm
H5FD_mpi_get_comm(const H5FD_t *file)
{
const H5FD_class_mpi_t *cls;
MPI_Comm ret_value;
FUNC_ENTER_NOAPI(MPI_COMM_NULL)
HDassert(file);
cls = (const H5FD_class_mpi_t *)(file->cls);
HDassert(cls);
HDassert(cls->get_comm); /* All MPI drivers must implement this */
/* Dispatch to driver */
if ((ret_value = (cls->get_comm)(file)) == MPI_COMM_NULL)
HGOTO_ERROR(H5E_VFL, H5E_CANTGET, MPI_COMM_NULL, "driver get_comm request failed")
done:
FUNC_LEAVE_NOAPI(ret_value)
} /* end H5FD_mpi_get_comm() */
/*-------------------------------------------------------------------------
* Function: H5FD_mpi_MPIOff_to_haddr
*
* Purpose: Convert an MPI_Offset value to haddr_t.
*
* Return: Success: The haddr_t equivalent of the MPI_OFF
* argument.
*
* Failure: HADDR_UNDEF
*
* Programmer: Unknown
* January 30, 1998
*
*-------------------------------------------------------------------------
*/
haddr_t
H5FD_mpi_MPIOff_to_haddr(MPI_Offset mpi_off)
{
haddr_t ret_value = HADDR_UNDEF;
FUNC_ENTER_NOAPI_NOINIT_NOERR
if (mpi_off != (MPI_Offset)(haddr_t)mpi_off)
ret_value = HADDR_UNDEF;
else
ret_value = (haddr_t)mpi_off;
FUNC_LEAVE_NOAPI(ret_value)
}
/*-------------------------------------------------------------------------
* Function: H5FD_mpi_haddr_to_MPIOff
*
* Purpose: Convert an haddr_t value to MPI_Offset.
*
* Return: Success: Non-negative, the MPI_OFF argument contains
* the converted value.
*
* Failure: Negative, MPI_OFF is undefined.
*
* Programmer: Unknown
* January 30, 1998
*
*-------------------------------------------------------------------------
*/
herr_t
H5FD_mpi_haddr_to_MPIOff(haddr_t addr, MPI_Offset *mpi_off /*out*/)
{
herr_t ret_value = FAIL;
FUNC_ENTER_NOAPI_NOINIT_NOERR
HDassert(mpi_off);
/* Convert the HDF5 address into an MPI offset */
*mpi_off = (MPI_Offset)addr;
if (addr != (haddr_t)((MPI_Offset)addr))
ret_value = FAIL;
else
ret_value = SUCCEED;
FUNC_LEAVE_NOAPI(ret_value)
}
/*-------------------------------------------------------------------------
* Function: H5FD_mpi_comm_info_dup
*
* Purpose: Make duplicates of communicator and Info object.
* If the Info object is in fact MPI_INFO_NULL, no duplicate
* is made but the same value assigned to the new Info object
* handle.
*
* Return: Success: Non-negative. The new communicator and Info
* object handles are returned via comm_new and
* info_new pointers.
*
* Failure: Negative.
*
* Programmer: Albert Cheng
* Jan 8, 2003
*-------------------------------------------------------------------------
*/
herr_t
H5FD_mpi_comm_info_dup(MPI_Comm comm, MPI_Info info, MPI_Comm *comm_new, MPI_Info *info_new)
{
herr_t ret_value = SUCCEED;
MPI_Comm comm_dup = MPI_COMM_NULL;
MPI_Info info_dup = MPI_INFO_NULL;
int mpi_code;
FUNC_ENTER_NOAPI(FAIL)
/* Check arguments */
if (MPI_COMM_NULL == comm)
HGOTO_ERROR(H5E_INTERNAL, H5E_BADVALUE, FAIL, "not a valid argument")
if (!comm_new || !info_new)
HGOTO_ERROR(H5E_INTERNAL, H5E_BADVALUE, FAIL, "bad pointers")
/* Dup them. Using temporary variables for error recovery cleanup. */
if (MPI_SUCCESS != (mpi_code = MPI_Comm_dup(comm, &comm_dup)))
HMPI_GOTO_ERROR(FAIL, "MPI_Comm_dup failed", mpi_code)
if (MPI_INFO_NULL != info) {
if (MPI_SUCCESS != (mpi_code = MPI_Info_dup(info, &info_dup)))
HMPI_GOTO_ERROR(FAIL, "MPI_Info_dup failed", mpi_code)
}
else {
/* No dup, just copy it. */
info_dup = info;
}
/* Set MPI_ERRORS_RETURN on comm_dup so that MPI failures are not fatal,
and return codes can be checked and handled. May 23, 2017 FTW */
if (MPI_SUCCESS != (mpi_code = MPI_Comm_set_errhandler(comm_dup, MPI_ERRORS_RETURN)))
HMPI_GOTO_ERROR(FAIL, "MPI_Errhandler_set failed", mpi_code)
/* copy them to the return arguments */
*comm_new = comm_dup;
*info_new = info_dup;
done:
if (FAIL == ret_value) {
/* need to free anything created here */
if (MPI_COMM_NULL != comm_dup)
MPI_Comm_free(&comm_dup);
if (MPI_INFO_NULL != info_dup)
MPI_Info_free(&info_dup);
}
FUNC_LEAVE_NOAPI(ret_value)
}
/*-------------------------------------------------------------------------
* Function: H5FD_mpi_comm_info_free
*
* Purpose: Free the communicator and Info object.
* If comm or info is in fact MPI_COMM_NULL or MPI_INFO_NULL
* respectively, no action occurs to it.
*
* Return: Success: Non-negative. The values the pointers refer
* to will be set to the corresponding NULL
* handles.
*
* Failure: Negative.
*
* Programmer: Albert Cheng
* Jan 8, 2003
*-------------------------------------------------------------------------
*/
herr_t
H5FD_mpi_comm_info_free(MPI_Comm *comm, MPI_Info *info)
{
herr_t ret_value = SUCCEED;
FUNC_ENTER_NOAPI(FAIL)
/* Check arguments */
if (!comm || !info)
HGOTO_ERROR(H5E_INTERNAL, H5E_BADVALUE, FAIL, "not a valid argument")
if (MPI_COMM_NULL != *comm)
MPI_Comm_free(comm);
if (MPI_INFO_NULL != *info)
MPI_Info_free(info);
done:
FUNC_LEAVE_NOAPI(ret_value)
}
#ifdef NOT_YET
/*-------------------------------------------------------------------------
* Function: H5FD_mpio_wait_for_left_neighbor
*
* Purpose: Blocks until (empty) msg is received from immediately
* lower-rank neighbor. In conjunction with
* H5FD_mpio_signal_right_neighbor, useful for enforcing
* 1-process-at-at-time access to critical regions to avoid race
* conditions (though it is overkill to require that the
* processes be allowed to proceed strictly in order of their
* rank).
*
* Note: This routine doesn't read or write any file, just performs
* interprocess coordination. It really should reside in a
* separate package of such routines.
*
* Return: Success: 0
* Failure: -1
*
* Programmer: rky
* 19981207
*
*-------------------------------------------------------------------------
*/
herr_t
H5FD_mpio_wait_for_left_neighbor(H5FD_t *_file)
{
H5FD_mpio_t *file = (H5FD_mpio_t *)_file;
char msgbuf[1];
MPI_Status rcvstat;
int mpi_code; /* mpi return code */
herr_t ret_value = SUCCEED; /* Return value */
FUNC_ENTER_NOAPI(FAIL)
HDassert(file);
HDassert(H5FD_MPIO == file->pub.driver_id);
/* Portably initialize MPI status variable */
HDmemset(&rcvstat, 0, sizeof(MPI_Status));
/* p0 has no left neighbor; all other procs wait for msg */
if (file->mpi_rank != 0) {
if (MPI_SUCCESS != (mpi_code = MPI_Recv(&msgbuf, 1, MPI_CHAR, file->mpi_rank - 1, MPI_ANY_TAG,
file->comm, &rcvstat)))
HMPI_GOTO_ERROR(FAIL, "MPI_Recv failed", mpi_code)
}
done:
FUNC_LEAVE_NOAPI(ret_value)
}
/*-------------------------------------------------------------------------
* Function: H5FD_mpio_signal_right_neighbor
*
* Purpose: Blocks until (empty) msg is received from immediately
* lower-rank neighbor. In conjunction with
* H5FD_mpio_wait_for_left_neighbor, useful for enforcing
* 1-process-at-at-time access to critical regions to avoid race
* conditions (though it is overkill to require that the
* processes be allowed to proceed strictly in order of their
* rank).
*
* Note: This routine doesn't read or write any file, just performs
* interprocess coordination. It really should reside in a
* separate package of such routines.
*
* Return: Success: 0
* Failure: -1
*
* Programmer: rky
* 19981207
*
*-------------------------------------------------------------------------
*/
herr_t
H5FD_mpio_signal_right_neighbor(H5FD_t *_file)
{
H5FD_mpio_t *file = (H5FD_mpio_t *)_file;
char msgbuf[1];
int mpi_code; /* mpi return code */
herr_t ret_value = SUCCEED; /* Return value */
FUNC_ENTER_NOAPI(FAIL)
HDassert(file);
HDassert(H5FD_MPIO == file->pub.driver_id);
if (file->mpi_rank != (file->mpi_size - 1))
if (MPI_SUCCESS !=
(mpi_code = MPI_Send(&msgbuf, 0 /*empty msg*/, MPI_CHAR, file->mpi_rank + 1, 0, file->comm)))
HMPI_GOTO_ERROR(FAIL, "MPI_Send failed", mpi_code)
done:
FUNC_LEAVE_NOAPI(ret_value)
}
#endif /* NOT_YET */
#endif /* H5_HAVE_PARALLEL */
|