Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Merge branch 'depencencies' of https://framagit.org/simgrid/simgrid into depencencies
[simgrid.git] / src / smpi / include / smpi_coll.hpp
1 /* High level handling of collective algorithms */
2 /* Copyright (c) 2009-2020. The SimGrid Team. All rights reserved.          */
3
4 /* This program is free software; you can redistribute it and/or modify it
5  * under the terms of the license (GNU LGPL) which comes with this package. */
6
7 #ifndef SMPI_COLL_HPP
8 #define SMPI_COLL_HPP
9
10 #include "private.hpp"
11 #include "xbt/base.h"
12
13 /** @brief MPI collective description */
14
15 #define COLL_DEFS(cat, ret, args, args2)                                                                               \
16   extern int(*cat) args;
17
18 #define COLL_UNPAREN(...)  __VA_ARGS__
19 #define COLL_APPLY(action, sig, name) action(sig, name)
20
21 #define COLL_GATHER_SIG gather, int,                                                                                   \
22     (const void *send_buff, int send_count, MPI_Datatype send_type,                                                    \
23            void *recv_buff, int recv_count, MPI_Datatype recv_type,  int root, MPI_Comm comm)
24 #define COLL_ALLGATHER_SIG allgather, int,                                                                             \
25     (const void *send_buff, int send_count, MPI_Datatype send_type,                                                    \
26             void *recv_buff, int recv_count, MPI_Datatype recv_type,  MPI_Comm comm)
27 #define COLL_ALLGATHERV_SIG allgatherv, int,                                                                           \
28     (const void *send_buff, int send_count, MPI_Datatype send_type,                                                    \
29            void *recv_buff, const int *recv_count, const int *recv_disps,  MPI_Datatype recv_type, MPI_Comm comm)
30 #define COLL_ALLTOALL_SIG alltoall, int,                                                                               \
31     (const void *send_buff, int send_count, MPI_Datatype send_type,                                                    \
32            void *recv_buff, int recv_count, MPI_Datatype recv_type,  MPI_Comm comm)
33 #define COLL_ALLTOALLV_SIG alltoallv, int,                                                                             \
34     (const void *send_buff, const int *send_counts, const int *send_disps, MPI_Datatype send_type,                     \
35            void *recv_buff, const int *recv_counts, const int *recv_disps, MPI_Datatype recv_type, MPI_Comm comm)
36 #define COLL_BCAST_SIG bcast, int,                                                                                     \
37     (void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
38 #define COLL_REDUCE_SIG reduce, int,                                                                                   \
39     (const void *buf, void *rbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
40 #define COLL_ALLREDUCE_SIG allreduce, int,                                                                             \
41     (const void *sbuf, void *rbuf, int rcount,  MPI_Datatype dtype, MPI_Op op, MPI_Comm comm)
42 #define COLL_REDUCE_SCATTER_SIG reduce_scatter, int,                                                                   \
43     (const void *sbuf, void *rbuf, const int *rcounts, MPI_Datatype dtype,MPI_Op  op,MPI_Comm  comm)
44 #define COLL_SCATTER_SIG scatter, int,                                                                                 \
45     (const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm)
46 #define COLL_BARRIER_SIG barrier, int,                                                                                 \
47     (MPI_Comm comm)
48
49 namespace simgrid {
50 namespace smpi {
51
52 struct s_mpi_coll_description_t {
53   std::string name;
54   std::string description;
55   void *coll;
56 };
57
58 namespace colls {
59 void set_collectives();
60 XBT_PUBLIC std::vector<s_mpi_coll_description_t>* get_smpi_coll_descriptions(const std::string& name);
61
62 void set_gather(const std::string& name);
63 void set_allgather(const std::string& name);
64 void set_allgatherv(const std::string& name);
65 void set_reduce(const std::string& name);
66 void set_allreduce(const std::string& name);
67 void set_reduce_scatter(const std::string& name);
68 void set_scatter(const std::string& name);
69 void set_barrier(const std::string& name);
70 void set_bcast(const std::string& name);
71 void set_alltoall(const std::string& name);
72 void set_alltoallv(const std::string& name);
73
74 // for each collective type, create the function pointer
75 //  extern int(*gather)(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count,
76 //                      MPI_Datatype recv_type, int root, MPI_Comm comm);
77 COLL_APPLY(COLL_DEFS, COLL_GATHER_SIG, "")
78 COLL_APPLY(COLL_DEFS, COLL_ALLGATHER_SIG, "")
79 COLL_APPLY(COLL_DEFS, COLL_ALLGATHERV_SIG, "")
80 COLL_APPLY(COLL_DEFS, COLL_REDUCE_SIG, "")
81 COLL_APPLY(COLL_DEFS, COLL_ALLREDUCE_SIG, "")
82 COLL_APPLY(COLL_DEFS, COLL_REDUCE_SCATTER_SIG, "")
83 COLL_APPLY(COLL_DEFS, COLL_SCATTER_SIG, "")
84 COLL_APPLY(COLL_DEFS, COLL_BARRIER_SIG, "")
85 COLL_APPLY(COLL_DEFS, COLL_BCAST_SIG, "")
86 COLL_APPLY(COLL_DEFS, COLL_ALLTOALL_SIG, "")
87 COLL_APPLY(COLL_DEFS, COLL_ALLTOALLV_SIG, "")
88
89 // These fairly unused collectives only have one implementation in SMPI
90 int gatherv(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, const int* recvcounts,
91             const int* displs, MPI_Datatype recvtype, int root, MPI_Comm comm);
92 int scatterv(const void* sendbuf, const int* sendcounts, const int* displs, MPI_Datatype sendtype, void* recvbuf,
93              int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm);
94 int scan(const void* sendbuf, void* recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm);
95 int exscan(const void* sendbuf, void* recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm);
96 int alltoallw(const void* sendbuf, const int* sendcounts, const int* senddisps, const MPI_Datatype* sendtypes,
97               void* recvbuf, const int* recvcounts, const int* recvdisps, const MPI_Datatype* recvtypes, MPI_Comm comm);
98
99 // async collectives
100 int ibarrier(MPI_Comm comm, MPI_Request* request, int external = 1);
101 int ibcast(void* buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm, MPI_Request* request,
102            int external = 1);
103 int igather(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
104             MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request* request, int external = 1);
105 int igatherv(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, const int* recvcounts,
106              const int* displs, MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request* request, int external = 1);
107 int iallgather(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
108                MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request, int external = 1);
109 int iallgatherv(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, const int* recvcounts,
110                 const int* displs, MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request, int external = 1);
111 int iscatter(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
112              MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request* request, int external = 1);
113 int iscatterv(const void* sendbuf, const int* sendcounts, const int* displs, MPI_Datatype sendtype, void* recvbuf,
114               int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request* request, int external = 1);
115 int ireduce(const void* sendbuf, void* recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm,
116             MPI_Request* request, int external = 1);
117 int iallreduce(const void* sendbuf, void* recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm,
118                MPI_Request* request, int external = 1);
119 int iscan(const void* sendbuf, void* recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm,
120           MPI_Request* request, int external = 1);
121 int iexscan(const void* sendbuf, void* recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm,
122             MPI_Request* request, int external = 1);
123 int ireduce_scatter(const void* sendbuf, void* recvbuf, const int* recvcounts, MPI_Datatype datatype, MPI_Op op,
124                     MPI_Comm comm, MPI_Request* request, int external = 1);
125 int ireduce_scatter_block(const void* sendbuf, void* recvbuf, int recvcount, MPI_Datatype datatype, MPI_Op op,
126                           MPI_Comm comm, MPI_Request* request, int external = 1);
127 int ialltoall(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
128               MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request, int external = 1);
129 int ialltoallv(const void* sendbuf, const int* sendcounts, const int* senddisps, MPI_Datatype sendtype, void* recvbuf,
130                const int* recvcounts, const int* recvdisps, MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request,
131                int external = 1);
132 int ialltoallw(const void* sendbuf, const int* sendcounts, const int* senddisps, const MPI_Datatype* sendtypes,
133                void* recvbuf, const int* recvcounts, const int* recvdisps, const MPI_Datatype* recvtypes, MPI_Comm comm,
134                MPI_Request* request, int external = 1);
135
136 extern void (*smpi_coll_cleanup_callback)();
137 }
138
139 /***********************************************
140  * Prototypes of each and every implementation *
141  ***********************************************/
142
143 int gather__default(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, int root, MPI_Comm comm);
144 int gather__ompi(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, int root, MPI_Comm comm);
145 int gather__ompi_basic_linear(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, int root, MPI_Comm comm);
146 int gather__ompi_binomial(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, int root, MPI_Comm comm);
147 int gather__ompi_linear_sync(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, int root, MPI_Comm comm);
148 int gather__mpich(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, int root, MPI_Comm comm);
149 int gather__mvapich2(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, int root, MPI_Comm comm);
150 int gather__mvapich2_two_level(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, int root, MPI_Comm comm);
151 int gather__impi(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, int root, MPI_Comm comm);
152 int gather__automatic(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, int root, MPI_Comm comm);
153
154 int allgather__default(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
155 int allgather__2dmesh(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
156 int allgather__3dmesh(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
157 int allgather__bruck(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
158 int allgather__GB(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
159 int allgather__loosely_lr(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
160 int allgather__NTSLR(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
161 int allgather__NTSLR_NB(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
162 int allgather__pair(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
163 int allgather__rdb(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
164 int allgather__rhv(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
165 int allgather__ring(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
166 int allgather__SMP_NTS(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
167 int allgather__smp_simple(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
168 int allgather__spreading_simple(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
169 int allgather__ompi(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
170 int allgather__ompi_neighborexchange (const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
171 int allgather__mvapich2(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
172 int allgather__mvapich2_smp(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
173 int allgather__mpich(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
174 int allgather__impi(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
175 int allgather__automatic(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
176
177 int allgatherv__default(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, const int *recv_count, const int *recv_disps, MPI_Datatype recv_type, MPI_Comm comm);
178 int allgatherv__GB(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, const int *recv_count, const int *recv_disps, MPI_Datatype recv_type, MPI_Comm comm);
179 int allgatherv__pair(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, const int *recv_count, const int *recv_disps, MPI_Datatype recv_type, MPI_Comm comm);
180 int allgatherv__ring(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, const int *recv_count, const int *recv_disps, MPI_Datatype recv_type, MPI_Comm comm);
181 int allgatherv__ompi(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, const int *recv_count, const int *recv_disps, MPI_Datatype recv_type, MPI_Comm comm);
182 int allgatherv__ompi_neighborexchange(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, const int *recv_count, const int *recv_disps, MPI_Datatype recv_type, MPI_Comm comm);
183 int allgatherv__ompi_bruck (const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, const int *recv_count, const int *recv_disps, MPI_Datatype recv_type, MPI_Comm comm);
184 int allgatherv__mpich(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, const int *recv_count, const int *recv_disps, MPI_Datatype recv_type, MPI_Comm comm);
185 int allgatherv__mpich_rdb(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, const int *recv_count, const int *recv_disps, MPI_Datatype recv_type, MPI_Comm comm);
186 int allgatherv__mpich_ring(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, const int *recv_count, const int *recv_disps, MPI_Datatype recv_type, MPI_Comm comm);
187 int allgatherv__mvapich2(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, const int *recv_count, const int *recv_disps, MPI_Datatype recv_type, MPI_Comm comm);
188 int allgatherv__impi(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, const int *recv_count, const int *recv_disps, MPI_Datatype recv_type, MPI_Comm comm);
189 int allgatherv__automatic(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, const int *recv_count, const int *recv_disps, MPI_Datatype recv_type, MPI_Comm comm);
190
191 int allreduce__default(const void *sbuf, void *rbuf, int rcount, MPI_Datatype dtype, MPI_Op op, MPI_Comm comm);
192 int allreduce__lr(const void *sbuf, void *rbuf, int rcount, MPI_Datatype dtype, MPI_Op op, MPI_Comm comm);
193 int allreduce__rab1(const void *sbuf, void *rbuf, int rcount, MPI_Datatype dtype, MPI_Op op, MPI_Comm comm);
194 int allreduce__rab2(const void *sbuf, void *rbuf, int rcount, MPI_Datatype dtype, MPI_Op op, MPI_Comm comm);
195 int allreduce__rab_rdb(const void *sbuf, void *rbuf, int rcount, MPI_Datatype dtype, MPI_Op op, MPI_Comm comm);
196 int allreduce__rdb(const void *sbuf, void *rbuf, int rcount, MPI_Datatype dtype, MPI_Op op, MPI_Comm comm);
197 int allreduce__smp_binomial(const void *sbuf, void *rbuf, int rcount, MPI_Datatype dtype, MPI_Op op, MPI_Comm comm);
198 int allreduce__smp_binomial_pipeline(const void *sbuf, void *rbuf, int rcount, MPI_Datatype dtype, MPI_Op op, MPI_Comm comm);
199 int allreduce__smp_rdb(const void *sbuf, void *rbuf, int rcount, MPI_Datatype dtype, MPI_Op op, MPI_Comm comm);
200 int allreduce__smp_rsag(const void *sbuf, void *rbuf, int rcount, MPI_Datatype dtype, MPI_Op op, MPI_Comm comm);
201 int allreduce__smp_rsag_lr(const void *sbuf, void *rbuf, int rcount, MPI_Datatype dtype, MPI_Op op, MPI_Comm comm);
202 int allreduce__smp_rsag_rab(const void *sbuf, void *rbuf, int rcount, MPI_Datatype dtype, MPI_Op op, MPI_Comm comm);
203 int allreduce__redbcast(const void *sbuf, void *rbuf, int rcount, MPI_Datatype dtype, MPI_Op op, MPI_Comm comm);
204 int allreduce__ompi(const void *sbuf, void *rbuf, int rcount, MPI_Datatype dtype, MPI_Op op, MPI_Comm comm);
205 int allreduce__ompi_ring_segmented(const void *sbuf, void *rbuf, int rcount, MPI_Datatype dtype, MPI_Op op, MPI_Comm comm);
206 int allreduce__mpich(const void *sbuf, void *rbuf, int rcount, MPI_Datatype dtype, MPI_Op op, MPI_Comm comm);
207 int allreduce__mvapich2(const void *sbuf, void *rbuf, int rcount, MPI_Datatype dtype, MPI_Op op, MPI_Comm comm);
208 int allreduce__mvapich2_rs(const void *sbuf, void *rbuf, int rcount, MPI_Datatype dtype, MPI_Op op, MPI_Comm comm);
209 int allreduce__mvapich2_two_level(const void *sbuf, void *rbuf, int rcount, MPI_Datatype dtype, MPI_Op op, MPI_Comm comm);
210 int allreduce__impi(const void *sbuf, void *rbuf, int rcount, MPI_Datatype dtype, MPI_Op op, MPI_Comm comm);
211 int allreduce__rab(const void *sbuf, void *rbuf, int rcount, MPI_Datatype dtype, MPI_Op op, MPI_Comm comm);
212 int allreduce__automatic(const void *sbuf, void *rbuf, int rcount, MPI_Datatype dtype, MPI_Op op, MPI_Comm comm);
213
214 int alltoall__default(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
215 int alltoall__2dmesh(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
216 int alltoall__3dmesh(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
217 int alltoall__basic_linear(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
218 int alltoall__bruck(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
219 int alltoall__pair(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
220 int alltoall__pair_rma(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
221 int alltoall__pair_light_barrier(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
222 int alltoall__pair_mpi_barrier(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
223 int alltoall__pair_one_barrier(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
224 int alltoall__rdb(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
225 int alltoall__ring(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
226 int alltoall__ring_light_barrier(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
227 int alltoall__ring_mpi_barrier(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
228 int alltoall__ring_one_barrier(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
229 int alltoall__mvapich2(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
230 int alltoall__mvapich2_scatter_dest(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
231 int alltoall__ompi(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
232 int alltoall__mpich(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
233 int alltoall__impi(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
234 int alltoall__automatic(const void *send_buff, int send_count, MPI_Datatype send_type, void *recv_buff, int recv_count, MPI_Datatype recv_type, MPI_Comm comm);
235
236 int alltoallv__default(const void *send_buff, const int *send_counts, const int *send_disps, MPI_Datatype send_type, void *recv_buff, const int *recv_counts, const int *recv_disps, MPI_Datatype recv_type, MPI_Comm comm);
237 int alltoallv__bruck(const void *send_buff, const int *send_counts, const int *send_disps, MPI_Datatype send_type, void *recv_buff, const int *recv_counts, const int *recv_disps, MPI_Datatype recv_type, MPI_Comm comm);
238 int alltoallv__pair(const void *send_buff, const int *send_counts, const int *send_disps, MPI_Datatype send_type, void *recv_buff, const int *recv_counts, const int *recv_disps, MPI_Datatype recv_type, MPI_Comm comm);
239 int alltoallv__pair_light_barrier(const void *send_buff, const int *send_counts, const int *send_disps, MPI_Datatype send_type, void *recv_buff, const int *recv_counts, const int *recv_disps, MPI_Datatype recv_type, MPI_Comm comm);
240 int alltoallv__pair_mpi_barrier(const void *send_buff, const int *send_counts, const int *send_disps, MPI_Datatype send_type, void *recv_buff, const int *recv_counts, const int *recv_disps, MPI_Datatype recv_type, MPI_Comm comm);
241 int alltoallv__pair_one_barrier(const void *send_buff, const int *send_counts, const int *send_disps, MPI_Datatype send_type, void *recv_buff, const int *recv_counts, const int *recv_disps, MPI_Datatype recv_type, MPI_Comm comm);
242 int alltoallv__ring(const void *send_buff, const int *send_counts, const int *send_disps, MPI_Datatype send_type, void *recv_buff, const int *recv_counts, const int *recv_disps, MPI_Datatype recv_type, MPI_Comm comm);
243 int alltoallv__ring_light_barrier(const void *send_buff, const int *send_counts, const int *send_disps, MPI_Datatype send_type, void *recv_buff, const int *recv_counts, const int *recv_disps, MPI_Datatype recv_type, MPI_Comm comm);
244 int alltoallv__ring_mpi_barrier(const void *send_buff, const int *send_counts, const int *send_disps, MPI_Datatype send_type, void *recv_buff, const int *recv_counts, const int *recv_disps, MPI_Datatype recv_type, MPI_Comm comm);
245 int alltoallv__ring_one_barrier(const void *send_buff, const int *send_counts, const int *send_disps, MPI_Datatype send_type, void *recv_buff, const int *recv_counts, const int *recv_disps, MPI_Datatype recv_type, MPI_Comm comm);
246 int alltoallv__ompi(const void *send_buff, const int *send_counts, const int *send_disps, MPI_Datatype send_type, void *recv_buff, const int *recv_counts, const int *recv_disps, MPI_Datatype recv_type, MPI_Comm comm);
247 int alltoallv__mpich(const void *send_buff, const int *send_counts, const int *send_disps, MPI_Datatype send_type, void *recv_buff, const int *recv_counts, const int *recv_disps, MPI_Datatype recv_type, MPI_Comm comm);
248 int alltoallv__ompi_basic_linear(const void *send_buff, const int *send_counts, const int *send_disps, MPI_Datatype send_type, void *recv_buff, const int *recv_counts, const int *recv_disps, MPI_Datatype recv_type, MPI_Comm comm);
249 int alltoallv__mvapich2(const void *send_buff, const int *send_counts, const int *send_disps, MPI_Datatype send_type, void *recv_buff, const int *recv_counts, const int *recv_disps, MPI_Datatype recv_type, MPI_Comm comm);
250 int alltoallv__impi(const void *send_buff, const int *send_counts, const int *send_disps, MPI_Datatype send_type, void *recv_buff, const int *recv_counts, const int *recv_disps, MPI_Datatype recv_type, MPI_Comm comm);
251 int alltoallv__automatic(const void *send_buff, const int *send_counts, const int *send_disps, MPI_Datatype send_type, void *recv_buff, const int *recv_counts, const int *recv_disps, MPI_Datatype recv_type, MPI_Comm comm);
252
253 int bcast__default(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm);
254 int bcast__arrival_pattern_aware(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm);
255 int bcast__arrival_pattern_aware_wait(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm);
256 int bcast__arrival_scatter(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm);
257 int bcast__binomial_tree(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm);
258 int bcast__flattree(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm);
259 int bcast__flattree_pipeline(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm);
260 int bcast__NTSB(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm);
261 int bcast__NTSL(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm);
262 int bcast__NTSL_Isend(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm);
263 int bcast__scatter_LR_allgather(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm);
264 int bcast__scatter_rdb_allgather(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm);
265 int bcast__SMP_binary(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm);
266 int bcast__SMP_binomial(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm);
267 int bcast__SMP_linear(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm);
268 int bcast__ompi(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm);
269 int bcast__ompi_split_bintree(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm);
270 int bcast__ompi_pipeline(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm);
271 int bcast__mpich(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm);
272 int bcast__mvapich2(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm);
273 int bcast__mvapich2_inter_node(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm);
274 int bcast__mvapich2_intra_node(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm);
275 int bcast__mvapich2_knomial_intra_node(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm);
276 int bcast__impi(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm);
277 int bcast__automatic(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm);
278
279 int reduce__default(const void *buf, void *rbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm);
280 int reduce__arrival_pattern_aware(const void *buf, void *rbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm);
281 int reduce__binomial(const void *buf, void *rbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm);
282 int reduce__flat_tree(const void *buf, void *rbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm);
283 int reduce__NTSL(const void *buf, void *rbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm);
284 int reduce__scatter_gather(const void *buf, void *rbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm);
285 int reduce__ompi(const void *buf, void *rbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm);
286 int reduce__ompi_chain(const void *buf, void *rbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm);
287 int reduce__ompi_pipeline(const void *buf, void *rbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm);
288 int reduce__ompi_basic_linear(const void *buf, void *rbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm);
289 int reduce__ompi_in_order_binary(const void *buf, void *rbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm);
290 int reduce__ompi_binary(const void *buf, void *rbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm);
291 int reduce__ompi_binomial(const void *buf, void *rbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm);
292 int reduce__mpich(const void *buf, void *rbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm);
293 int reduce__mvapich2(const void *buf, void *rbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm);
294 int reduce__mvapich2_knomial(const void *buf, void *rbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm);
295 int reduce__mvapich2_two_level(const void *buf, void *rbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm);
296 int reduce__impi(const void *buf, void *rbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm);
297 int reduce__rab(const void *buf, void *rbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm);
298 int reduce__automatic(const void *buf, void *rbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm);
299
300 int reduce_scatter__default(const void *sbuf, void *rbuf, const int *rcounts, MPI_Datatype dtype,MPI_Op  op,MPI_Comm  comm);
301 int reduce_scatter__ompi(const void *sbuf, void *rbuf, const int *rcounts, MPI_Datatype dtype,MPI_Op  op,MPI_Comm  comm);
302 int reduce_scatter__ompi_basic_recursivehalving(const void *sbuf, void *rbuf, const int *rcounts, MPI_Datatype dtype,MPI_Op  op,MPI_Comm  comm);
303 int reduce_scatter__ompi_ring(const void *sbuf, void *rbuf, const int *rcounts, MPI_Datatype dtype,MPI_Op  op,MPI_Comm  comm);
304 int reduce_scatter__mpich(const void *sbuf, void *rbuf, const int *rcounts, MPI_Datatype dtype,MPI_Op  op,MPI_Comm  comm);
305 int reduce_scatter__mpich_pair(const void *sbuf, void *rbuf, const int *rcounts, MPI_Datatype dtype,MPI_Op  op,MPI_Comm  comm);
306 int reduce_scatter__mpich_rdb(const void *sbuf, void *rbuf, const int *rcounts, MPI_Datatype dtype,MPI_Op  op,MPI_Comm  comm);
307 int reduce_scatter__mpich_noncomm(const void *sbuf, void *rbuf, const int *rcounts, MPI_Datatype dtype,MPI_Op  op,MPI_Comm  comm);
308 int reduce_scatter__mvapich2(const void *sbuf, void *rbuf, const int *rcounts, MPI_Datatype dtype,MPI_Op  op,MPI_Comm  comm);
309 int reduce_scatter__impi(const void *sbuf, void *rbuf, const int *rcounts, MPI_Datatype dtype,MPI_Op  op,MPI_Comm  comm);
310 int reduce_scatter__automatic(const void *sbuf, void *rbuf, const int *rcounts, MPI_Datatype dtype,MPI_Op  op,MPI_Comm  comm);
311
312 int scatter__default(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm);
313 int scatter__ompi(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm);
314 int scatter__ompi_basic_linear(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm);
315 int scatter__ompi_binomial(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm);
316 int scatter__mpich(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm);
317 int scatter__mvapich2(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm);
318 int scatter__mvapich2_two_level_binomial(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm);
319 int scatter__mvapich2_two_level_direct(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm);
320 int scatter__impi (const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm);
321 int scatter__automatic (const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm);
322
323 int barrier__default(MPI_Comm comm);
324 int barrier__ompi(MPI_Comm comm);
325 int barrier__ompi_basic_linear(MPI_Comm comm);
326 int barrier__ompi_two_procs(MPI_Comm comm);
327 int barrier__ompi_tree(MPI_Comm comm);
328 int barrier__ompi_bruck(MPI_Comm comm);
329 int barrier__ompi_recursivedoubling(MPI_Comm comm);
330 int barrier__ompi_doublering(MPI_Comm comm);
331 int barrier__mpich_smp(MPI_Comm comm);
332 int barrier__mpich(MPI_Comm comm);
333 int barrier__mvapich2_pair(MPI_Comm comm);
334 int barrier__mvapich2 (MPI_Comm comm);
335 int barrier__impi(MPI_Comm comm);
336 int barrier__automatic(MPI_Comm comm);
337
338 }
339 }
340 #endif