File: mpicolltest.h

package info (click to toggle)
mpich 5.0.0-1
  • links: PTS, VCS
  • area: main
  • in suites: experimental
  • size: 251,828 kB
  • sloc: ansic: 1,323,147; cpp: 82,869; f90: 72,420; javascript: 40,763; perl: 28,296; sh: 19,399; python: 16,191; xml: 14,418; makefile: 9,468; fortran: 8,046; java: 4,635; pascal: 352; asm: 324; ruby: 176; awk: 27; lisp: 19; php: 8; sed: 4
file content (460 lines) | stat: -rw-r--r-- 18,427 bytes parent folder | download
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
/*
 * Copyright (C) by Argonne National Laboratory
 *     See COPYRIGHT in top-level directory
 */

#ifndef MPICOLLTEST_H_INCLUDED
#define MPICOLLTEST_H_INCLUDED

/* This file defines MTest wrapper functions for MPI collectives.
 * Test uses MTest_Collective calls that runs either the blocking or
 * non-blocking version of the collective depend on the is_blocking
 * parameter.
 */

#include "mpi.h"
#include "mpitestconf.h"

/* Wrap up MTest_Collective calls by non-blocking collectives.
 * It requires MPI_VERSION >= 3. */
#if MPI_VERSION >= 3
static inline int MTest_Barrier(int is_blocking, MPI_Comm comm)
{
    if (is_blocking) {
        return MPI_Barrier(comm);
    } else {
        int mpi_errno;
        MPI_Request req = MPI_REQUEST_NULL;

        mpi_errno = MPI_Ibarrier(comm, &req);
        if (mpi_errno != MPI_SUCCESS)
            return mpi_errno;
        mpi_errno = MPI_Wait(&req, MPI_STATUS_IGNORE);
        return mpi_errno;
    }
}

static inline int MTest_Bcast(int is_blocking, void *buffer, int count, MPI_Datatype datatype,
                              int root, MPI_Comm comm)
{
    if (is_blocking) {
        return MPI_Bcast(buffer, count, datatype, root, comm);
    } else {
        int mpi_errno;
        MPI_Request req = MPI_REQUEST_NULL;

        mpi_errno = MPI_Ibcast(buffer, count, datatype, root, comm, &req);
        if (mpi_errno != MPI_SUCCESS)
            return mpi_errno;
        mpi_errno = MPI_Wait(&req, MPI_STATUS_IGNORE);
        return mpi_errno;
    }
}

static inline int MTest_Gather(int is_blocking, const void *sendbuf, int sendcount,
                               MPI_Datatype sendtype, void *recvbuf, int recvcount,
                               MPI_Datatype recvtype, int root, MPI_Comm comm)
{
    if (is_blocking) {
        return MPI_Gather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm);
    } else {
        int mpi_errno;
        MPI_Request req = MPI_REQUEST_NULL;

        mpi_errno = MPI_Igather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype,
                                root, comm, &req);
        if (mpi_errno != MPI_SUCCESS)
            return mpi_errno;
        mpi_errno = MPI_Wait(&req, MPI_STATUS_IGNORE);
        return mpi_errno;
    }
}

static inline int MTest_Gatherv(int is_blocking, const void *sendbuf, int sendcount,
                                MPI_Datatype sendtype, void *recvbuf, const int *recvcounts,
                                const int *displs, MPI_Datatype recvtype, int root, MPI_Comm comm)
{
    if (is_blocking) {
        return MPI_Gatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs,
                           recvtype, root, comm);
    } else {
        int mpi_errno;
        MPI_Request req = MPI_REQUEST_NULL;

        mpi_errno = MPI_Igatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs,
                                 recvtype, root, comm, &req);
        if (mpi_errno != MPI_SUCCESS)
            return mpi_errno;
        mpi_errno = MPI_Wait(&req, MPI_STATUS_IGNORE);
        return mpi_errno;
    }
}

static inline int MTest_Scatter(int is_blocking, const void *sendbuf, int sendcount,
                                MPI_Datatype sendtype, void *recvbuf, int recvcount,
                                MPI_Datatype recvtype, int root, MPI_Comm comm)
{
    if (is_blocking) {
        return MPI_Scatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm);
    } else {
        int mpi_errno;
        MPI_Request req = MPI_REQUEST_NULL;

        mpi_errno = MPI_Iscatter(sendbuf, sendcount, sendtype, recvbuf, recvcount,
                                 recvtype, root, comm, &req);
        if (mpi_errno != MPI_SUCCESS)
            return mpi_errno;
        mpi_errno = MPI_Wait(&req, MPI_STATUS_IGNORE);
        return mpi_errno;
    }
}

static inline int MTest_Scatterv(int is_blocking, const void *sendbuf, const int *sendcounts,
                                 const int *displs, MPI_Datatype sendtype, void *recvbuf,
                                 int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm)
{
    if (is_blocking) {
        return MPI_Scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf,
                            recvcount, recvtype, root, comm);
    } else {
        int mpi_errno;
        MPI_Request req = MPI_REQUEST_NULL;

        mpi_errno = MPI_Iscatterv(sendbuf, sendcounts, displs, sendtype, recvbuf,
                                  recvcount, recvtype, root, comm, &req);
        if (mpi_errno != MPI_SUCCESS)
            return mpi_errno;
        mpi_errno = MPI_Wait(&req, MPI_STATUS_IGNORE);
        return mpi_errno;
    }
}

static inline int MTest_Allgather(int is_blocking, const void *sendbuf, int sendcount,
                                  MPI_Datatype sendtype, void *recvbuf, int recvcount,
                                  MPI_Datatype recvtype, MPI_Comm comm)
{
    if (is_blocking) {
        return MPI_Allgather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
    } else {
        int mpi_errno;
        MPI_Request req = MPI_REQUEST_NULL;

        mpi_errno = MPI_Iallgather(sendbuf, sendcount, sendtype, recvbuf, recvcount,
                                   recvtype, comm, &req);
        if (mpi_errno != MPI_SUCCESS)
            return mpi_errno;
        mpi_errno = MPI_Wait(&req, MPI_STATUS_IGNORE);
        return mpi_errno;
    }
}

static inline int MTest_Allgatherv(int is_blocking, const void *sendbuf, int sendcount,
                                   MPI_Datatype sendtype, void *recvbuf, const int *recvcounts,
                                   const int *displs, MPI_Datatype recvtype, MPI_Comm comm)
{
    if (is_blocking) {
        return MPI_Allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
                              displs, recvtype, comm);
    } else {
        int mpi_errno;
        MPI_Request req = MPI_REQUEST_NULL;

        mpi_errno = MPI_Iallgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
                                    displs, recvtype, comm, &req);
        if (mpi_errno != MPI_SUCCESS)
            return mpi_errno;
        mpi_errno = MPI_Wait(&req, MPI_STATUS_IGNORE);
        return mpi_errno;
    }
}

static inline int MTest_Alltoall(int is_blocking, const void *sendbuf, int sendcount,
                                 MPI_Datatype sendtype, void *recvbuf, int recvcount,
                                 MPI_Datatype recvtype, MPI_Comm comm)
{
    if (is_blocking) {
        return MPI_Alltoall(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
    } else {
        int mpi_errno;
        MPI_Request req = MPI_REQUEST_NULL;

        mpi_errno = MPI_Ialltoall(sendbuf, sendcount, sendtype, recvbuf, recvcount,
                                  recvtype, comm, &req);
        if (mpi_errno != MPI_SUCCESS)
            return mpi_errno;
        mpi_errno = MPI_Wait(&req, MPI_STATUS_IGNORE);
        return mpi_errno;
    }
}

static inline int MTest_Alltoallv(int is_blocking, const void *sendbuf, const int *sendcounts,
                                  const int *sdispls, MPI_Datatype sendtype, void *recvbuf,
                                  const int *recvcounts, const int *rdispls, MPI_Datatype recvtype,
                                  MPI_Comm comm)
{
    if (is_blocking) {
        return MPI_Alltoallv(sendbuf, sendcounts, sdispls, sendtype, recvbuf,
                             recvcounts, rdispls, recvtype, comm);
    } else {
        int mpi_errno;
        MPI_Request req = MPI_REQUEST_NULL;

        mpi_errno = MPI_Ialltoallv(sendbuf, sendcounts, sdispls, sendtype, recvbuf,
                                   recvcounts, rdispls, recvtype, comm, &req);
        if (mpi_errno != MPI_SUCCESS)
            return mpi_errno;
        mpi_errno = MPI_Wait(&req, MPI_STATUS_IGNORE);
        return mpi_errno;
    }
}

static inline int MTest_Alltoallw(int is_blocking, const void *sendbuf, const int *sendcounts,
                                  const int *sdispls, const MPI_Datatype * sendtypes, void *recvbuf,
                                  const int *recvcounts, const int *rdispls,
                                  const MPI_Datatype * recvtypes, MPI_Comm comm)
{
    if (is_blocking) {
        return MPI_Alltoallw(sendbuf, sendcounts, sdispls, sendtypes, recvbuf,
                             recvcounts, rdispls, recvtypes, comm);
    } else {
        int mpi_errno;
        MPI_Request req = MPI_REQUEST_NULL;

        mpi_errno = MPI_Ialltoallw(sendbuf, sendcounts, sdispls, sendtypes, recvbuf,
                                   recvcounts, rdispls, recvtypes, comm, &req);
        if (mpi_errno != MPI_SUCCESS)
            return mpi_errno;
        mpi_errno = MPI_Wait(&req, MPI_STATUS_IGNORE);
        return mpi_errno;
    }
}

static inline int MTest_Reduce(int is_blocking, const void *sendbuf, void *recvbuf, int count,
                               MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
{
    if (is_blocking) {
        return MPI_Reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
    } else {
        int mpi_errno;
        MPI_Request req = MPI_REQUEST_NULL;

        mpi_errno = MPI_Ireduce(sendbuf, recvbuf, count, datatype, op, root, comm, &req);
        if (mpi_errno != MPI_SUCCESS)
            return mpi_errno;
        mpi_errno = MPI_Wait(&req, MPI_STATUS_IGNORE);
        return mpi_errno;
    }
}

static inline int MTest_Allreduce(int is_blocking, const void *sendbuf, void *recvbuf, int count,
                                  MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
{
    if (is_blocking) {
        return MPI_Allreduce(sendbuf, recvbuf, count, datatype, op, comm);
    } else {
        int mpi_errno;
        MPI_Request req = MPI_REQUEST_NULL;

        mpi_errno = MPI_Iallreduce(sendbuf, recvbuf, count, datatype, op, comm, &req);
        if (mpi_errno != MPI_SUCCESS)
            return mpi_errno;
        mpi_errno = MPI_Wait(&req, MPI_STATUS_IGNORE);
        return mpi_errno;
    }
}

static inline int MTest_Reduce_scatter(int is_blocking, const void *sendbuf, void *recvbuf,
                                       const int *recvcounts, MPI_Datatype datatype, MPI_Op op,
                                       MPI_Comm comm)
{
    if (is_blocking) {
        return MPI_Reduce_scatter(sendbuf, recvbuf, recvcounts, datatype, op, comm);
    } else {
        int mpi_errno;
        MPI_Request req = MPI_REQUEST_NULL;

        mpi_errno = MPI_Ireduce_scatter(sendbuf, recvbuf, recvcounts, datatype, op, comm, &req);
        if (mpi_errno != MPI_SUCCESS)
            return mpi_errno;
        mpi_errno = MPI_Wait(&req, MPI_STATUS_IGNORE);
        return mpi_errno;
    }
}

static inline int MTest_Reduce_scatter_block(int is_blocking, const void *sendbuf, void *recvbuf,
                                             int recvcount, MPI_Datatype datatype, MPI_Op op,
                                             MPI_Comm comm)
{
    if (is_blocking) {
        return MPI_Reduce_scatter_block(sendbuf, recvbuf, recvcount, datatype, op, comm);
    } else {
        int mpi_errno;
        MPI_Request req = MPI_REQUEST_NULL;

        mpi_errno =
            MPI_Ireduce_scatter_block(sendbuf, recvbuf, recvcount, datatype, op, comm, &req);
        if (mpi_errno != MPI_SUCCESS)
            return mpi_errno;
        mpi_errno = MPI_Wait(&req, MPI_STATUS_IGNORE);
        return mpi_errno;
    }
}

static inline int MTest_Scan(int is_blocking, const void *sendbuf, void *recvbuf, int count,
                             MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
{
    if (is_blocking) {
        return MPI_Scan(sendbuf, recvbuf, count, datatype, op, comm);
    } else {
        int mpi_errno;
        MPI_Request req = MPI_REQUEST_NULL;

        mpi_errno = MPI_Iscan(sendbuf, recvbuf, count, datatype, op, comm, &req);
        if (mpi_errno != MPI_SUCCESS)
            return mpi_errno;
        mpi_errno = MPI_Wait(&req, MPI_STATUS_IGNORE);
        return mpi_errno;
    }
}

static inline int MTest_Exscan(int is_blocking, const void *sendbuf, void *recvbuf, int count,
                               MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
{
    if (is_blocking) {
        return MPI_Exscan(sendbuf, recvbuf, count, datatype, op, comm);
    } else {
        int mpi_errno;
        MPI_Request req = MPI_REQUEST_NULL;

        mpi_errno = MPI_Iexscan(sendbuf, recvbuf, count, datatype, op, comm, &req);
        if (mpi_errno != MPI_SUCCESS)
            return mpi_errno;
        mpi_errno = MPI_Wait(&req, MPI_STATUS_IGNORE);
        return mpi_errno;
    }
}
#else
/* If MPI_VERSION is less than 3,
 * wrap up MTest_Collective calls by traditional blocking collectives.*/

static inline int MTest_Barrier(int is_blocking, MPI_Comm comm)
{
    return MPI_Barrier(comm);
}

static inline int MTest_Bcast(int is_blocking, void *buffer, int count, MPI_Datatype datatype,
                              int root, MPI_Comm comm)
{
    return MPI_Bcast(buffer, count, datatype, root, comm);
}

static inline int MTest_Gather(int is_blocking, const void *sendbuf, int sendcount,
                               MPI_Datatype sendtype, void *recvbuf, int recvcount,
                               MPI_Datatype recvtype, int root, MPI_Comm comm)
{
    return MPI_Gather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm);
}

static inline int MTest_Gatherv(int is_blocking, const void *sendbuf, int sendcount,
                                MPI_Datatype sendtype, void *recvbuf, const int *recvcounts,
                                const int *displs, MPI_Datatype recvtype, int root, MPI_Comm comm)
{
    return MPI_Gatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs,
                       recvtype, root, comm);
}

static inline int MTest_Scatter(int is_blocking, const void *sendbuf, int sendcount,
                                MPI_Datatype sendtype, void *recvbuf, int recvcount,
                                MPI_Datatype recvtype, int root, MPI_Comm comm)
{
    return MPI_Scatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm);
}

static inline int MTest_Scatterv(int is_blocking, const void *sendbuf, const int *sendcounts,
                                 const int *displs, MPI_Datatype sendtype, void *recvbuf,
                                 int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm)
{
    return MPI_Scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf,
                        recvcount, recvtype, root, comm);
}

static inline int MTest_Allgather(int is_blocking, const void *sendbuf, int sendcount,
                                  MPI_Datatype sendtype, void *recvbuf, int recvcount,
                                  MPI_Datatype recvtype, MPI_Comm comm)
{
    return MPI_Allgather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
}

static inline int MTest_Allgatherv(int is_blocking, const void *sendbuf, int sendcount,
                                   MPI_Datatype sendtype, void *recvbuf, const int *recvcounts,
                                   const int *displs, MPI_Datatype recvtype, MPI_Comm comm)
{
    return MPI_Allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
                          displs, recvtype, comm);
}

static inline int MTest_Alltoall(int is_blocking, const void *sendbuf, int sendcount,
                                 MPI_Datatype sendtype, void *recvbuf, int recvcount,
                                 MPI_Datatype recvtype, MPI_Comm comm)
{
    return MPI_Alltoall(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
}

static inline int MTest_Alltoallv(int is_blocking, const void *sendbuf, const int *sendcounts,
                                  const int *sdispls, MPI_Datatype sendtype, void *recvbuf,
                                  const int *recvcounts, const int *rdispls, MPI_Datatype recvtype,
                                  MPI_Comm comm)
{
    return MPI_Alltoallv(sendbuf, sendcounts, sdispls, sendtype, recvbuf,
                         recvcounts, rdispls, recvtype, comm);
}

static inline int MTest_Alltoallw(int is_blocking, const void *sendbuf, const int *sendcounts,
                                  const int *sdispls, const MPI_Datatype * sendtypes, void *recvbuf,
                                  const int *recvcounts, const int *rdispls,
                                  const MPI_Datatype * recvtypes, MPI_Comm comm)
{
    return MPI_Alltoallw(sendbuf, sendcounts, sdispls, sendtypes, recvbuf,
                         recvcounts, rdispls, recvtypes, comm);
}

static inline int MTest_Reduce(int is_blocking, const void *sendbuf, void *recvbuf, int count,
                               MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
{
    return MPI_Reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
}

static inline int MTest_Allreduce(int is_blocking, const void *sendbuf, void *recvbuf, int count,
                                  MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
{
    return MPI_Allreduce(sendbuf, recvbuf, count, datatype, op, comm);
}

static inline int MTest_Reduce_scatter(int is_blocking, const void *sendbuf, void *recvbuf,
                                       const int *recvcounts, MPI_Datatype datatype, MPI_Op op,
                                       MPI_Comm comm)
{
    return MPI_Reduce_scatter(sendbuf, recvbuf, recvcounts, datatype, op, comm);
}

static inline int MTest_Reduce_scatter_block(int is_blocking, const void *sendbuf, void *recvbuf,
                                             int recvcount, MPI_Datatype datatype, MPI_Op op,
                                             MPI_Comm comm)
{
    return MPI_Reduce_scatter_block(sendbuf, recvbuf, recvcount, datatype, op, comm);
}

static inline int MTest_Scan(int is_blocking, const void *sendbuf, void *recvbuf, int count,
                             MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
{
    return MPI_Scan(sendbuf, recvbuf, count, datatype, op, comm);
}

static inline int MTest_Exscan(int is_blocking, const void *sendbuf, void *recvbuf, int count,
                               MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
{
    return MPI_Exscan(sendbuf, recvbuf, count, datatype, op, comm);
}

#endif /* MPI_VERSION >= 3 */
#endif /* MPICOLLTEST_H_INCLUDED */