Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Merge branch 'master' of https://framagit.org/simgrid/simgrid
[simgrid.git] / src / smpi / bindings / smpi_pmpi_coll.cpp
1 /* Copyright (c) 2007-2019. The SimGrid Team. All rights reserved.          */
2
3 /* This program is free software; you can redistribute it and/or modify it
4  * under the terms of the license (GNU LGPL) which comes with this package. */
5
6 #include "private.hpp"
7 #include "smpi_coll.hpp"
8 #include "smpi_comm.hpp"
9 #include "smpi_request.hpp"
10 #include "smpi_datatype_derived.hpp"
11 #include "smpi_op.hpp"
12 #include "src/smpi/include/smpi_actor.hpp"
13
14 XBT_LOG_EXTERNAL_DEFAULT_CATEGORY(smpi_pmpi);
15
16 #define CHECK_ARGS(test, errcode, ...)                                                                                 \
17   if (test) {                                                                                                          \
18     XBT_WARN(__VA_ARGS__);                                                                                             \
19     return (errcode);                                                                                                  \
20   }
21
22   static const void* smpi_get_in_place_buf(const void* inplacebuf, const void* otherbuf,std::unique_ptr<unsigned char[]>& tmp_sendbuf, int count, MPI_Datatype datatype){
23   if (inplacebuf == MPI_IN_PLACE) {
24       tmp_sendbuf.reset(new unsigned char[count * datatype->get_extent()]);
25       simgrid::smpi::Datatype::copy(otherbuf, count, datatype, tmp_sendbuf.get(), count, datatype);
26     return tmp_sendbuf.get();
27   }else{
28     return inplacebuf;
29   }
30 }
31 /* PMPI User level calls */
32
33 int PMPI_Barrier(MPI_Comm comm)
34 {
35   return PMPI_Ibarrier(comm, MPI_REQUEST_IGNORED);
36 }
37
38 int PMPI_Ibarrier(MPI_Comm comm, MPI_Request *request)
39 {
40   if (comm == MPI_COMM_NULL)
41     return MPI_ERR_COMM;
42   if (request == nullptr)
43     return MPI_ERR_ARG;
44
45   smpi_bench_end();
46   int rank = simgrid::s4u::this_actor::get_pid();
47   TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Barrier" : "PMPI_Ibarrier",
48                      new simgrid::instr::NoOpTIData(request == MPI_REQUEST_IGNORED ? "barrier" : "ibarrier"));
49   if (request == MPI_REQUEST_IGNORED) {
50     simgrid::smpi::colls::barrier(comm);
51     // Barrier can be used to synchronize RMA calls. Finish all requests from comm before.
52     comm->finish_rma_calls();
53   } else
54     simgrid::smpi::colls::ibarrier(comm, request);
55
56   TRACE_smpi_comm_out(rank);
57   smpi_bench_begin();
58   return MPI_SUCCESS;
59 }
60
61 int PMPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
62 {
63   return PMPI_Ibcast(buf, count, datatype, root, comm, MPI_REQUEST_IGNORED);
64 }
65
66 int PMPI_Ibcast(void *buf, int count, MPI_Datatype datatype, 
67                    int root, MPI_Comm comm, MPI_Request* request)
68 {
69   if (comm == MPI_COMM_NULL)
70     return MPI_ERR_COMM;
71   if (buf == nullptr && count > 0)
72     return MPI_ERR_BUFFER;
73   if (datatype == MPI_DATATYPE_NULL || not datatype->is_valid())
74     return MPI_ERR_TYPE;
75   if (count < 0)
76     return MPI_ERR_COUNT;
77   if (root < 0 || root >= comm->size())
78     return MPI_ERR_ROOT;
79   if (request == nullptr)
80     return MPI_ERR_ARG;
81
82   smpi_bench_end();
83   int rank = simgrid::s4u::this_actor::get_pid();
84   TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Bcast" : "PMPI_Ibcast",
85                      new simgrid::instr::CollTIData(request == MPI_REQUEST_IGNORED ? "bcast" : "ibcast", root, -1.0,
86                                                     datatype->is_replayable() ? count : count * datatype->size(), -1,
87                                                     simgrid::smpi::Datatype::encode(datatype), ""));
88   if (comm->size() > 1) {
89     if (request == MPI_REQUEST_IGNORED)
90       simgrid::smpi::colls::bcast(buf, count, datatype, root, comm);
91     else
92       simgrid::smpi::colls::ibcast(buf, count, datatype, root, comm, request);
93   } else {
94     if (request != MPI_REQUEST_IGNORED)
95       *request = MPI_REQUEST_NULL;
96   }
97
98   TRACE_smpi_comm_out(rank);
99   smpi_bench_begin();
100   return MPI_SUCCESS;
101 }
102
103 int PMPI_Gather(const void *sendbuf, int sendcount, MPI_Datatype sendtype,void *recvbuf, int recvcount, MPI_Datatype recvtype,
104                 int root, MPI_Comm comm){
105   return PMPI_Igather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm, MPI_REQUEST_IGNORED);
106 }
107
108 int PMPI_Igather(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
109                  MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request* request)
110 {
111   if (comm == MPI_COMM_NULL)
112     return MPI_ERR_COMM;
113   if ((sendbuf == nullptr && sendcount > 0) || ((comm->rank() == root) && recvbuf == nullptr && recvcount > 0))
114     return MPI_ERR_BUFFER;
115   if (((sendbuf != MPI_IN_PLACE && sendcount > 0) && (sendtype == MPI_DATATYPE_NULL)) ||
116       ((comm->rank() == root) && (recvtype == MPI_DATATYPE_NULL)))
117     return MPI_ERR_TYPE;
118   if (((sendbuf != MPI_IN_PLACE) && (sendcount < 0)) || ((comm->rank() == root) && (recvcount < 0)))
119     return MPI_ERR_COUNT;
120   if (root < 0 || root >= comm->size())
121     return MPI_ERR_ROOT;
122   if (request == nullptr)
123     return MPI_ERR_ARG;
124
125   smpi_bench_end();
126   const void* real_sendbuf   = sendbuf;
127   int real_sendcount         = sendcount;
128   MPI_Datatype real_sendtype = sendtype;
129   if ((comm->rank() == root) && (sendbuf == MPI_IN_PLACE)) {
130     real_sendcount = 0;
131     real_sendtype  = recvtype;
132   }
133   int rank = simgrid::s4u::this_actor::get_pid();
134
135   TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Gather" : "PMPI_Igather",
136                      new simgrid::instr::CollTIData(
137                          request == MPI_REQUEST_IGNORED ? "gather" : "igather", root, -1.0,
138                          real_sendtype->is_replayable() ? real_sendcount : real_sendcount * real_sendtype->size(),
139                          (comm->rank() != root || recvtype->is_replayable()) ? recvcount : recvcount * recvtype->size(),
140                          simgrid::smpi::Datatype::encode(real_sendtype), simgrid::smpi::Datatype::encode(recvtype)));
141   if (request == MPI_REQUEST_IGNORED)
142     simgrid::smpi::colls::gather(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcount, recvtype, root, comm);
143   else
144     simgrid::smpi::colls::igather(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcount, recvtype, root, comm,
145                                   request);
146
147   TRACE_smpi_comm_out(rank);
148   smpi_bench_begin();
149   return MPI_SUCCESS;
150 }
151
152 int PMPI_Gatherv(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, const int *recvcounts, const int *displs,
153                 MPI_Datatype recvtype, int root, MPI_Comm comm){
154   return PMPI_Igatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, root, comm, MPI_REQUEST_IGNORED);
155 }
156
157 int PMPI_Igatherv(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, const int* recvcounts, const int* displs,
158                   MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request* request)
159 {
160   if (comm == MPI_COMM_NULL)
161     return MPI_ERR_COMM;
162   if ((sendbuf == nullptr && sendcount > 0) || ((comm->rank() == root) && recvbuf == nullptr))
163     return MPI_ERR_BUFFER;
164   if (((sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
165       ((comm->rank() == root) && (recvtype == MPI_DATATYPE_NULL)))
166     return MPI_ERR_TYPE;
167   if ((sendbuf != MPI_IN_PLACE) && (sendcount < 0))
168     return MPI_ERR_COUNT;
169   if ((comm->rank() == root) && (recvcounts == nullptr || displs == nullptr))
170     return MPI_ERR_ARG;
171   if (root < 0 || root >= comm->size())
172     return MPI_ERR_ROOT;
173   if (request == nullptr)
174     return MPI_ERR_ARG;
175
176   for (int i = 0; i < comm->size(); i++) {
177     if ((comm->rank() == root) && (recvcounts[i] < 0))
178       return MPI_ERR_COUNT;
179   }
180
181   smpi_bench_end();
182   const void* real_sendbuf   = sendbuf;
183   int real_sendcount         = sendcount;
184   MPI_Datatype real_sendtype = sendtype;
185   if ((comm->rank() == root) && (sendbuf == MPI_IN_PLACE)) {
186     real_sendcount = 0;
187     real_sendtype  = recvtype;
188   }
189
190   int rank         = simgrid::s4u::this_actor::get_pid();
191   int dt_size_recv = recvtype->is_replayable() ? 1 : recvtype->size();
192
193   std::vector<int>* trace_recvcounts = new std::vector<int>;
194   if (comm->rank() == root) {
195     for (int i = 0; i < comm->size(); i++) // copy data to avoid bad free
196       trace_recvcounts->push_back(recvcounts[i] * dt_size_recv);
197   }
198
199   TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Gatherv" : "PMPI_Igatherv",
200                      new simgrid::instr::VarCollTIData(
201                          request == MPI_REQUEST_IGNORED ? "gatherv" : "igatherv", root,
202                          real_sendtype->is_replayable() ? real_sendcount : real_sendcount * real_sendtype->size(),
203                          nullptr, dt_size_recv, trace_recvcounts, simgrid::smpi::Datatype::encode(real_sendtype),
204                          simgrid::smpi::Datatype::encode(recvtype)));
205   if (request == MPI_REQUEST_IGNORED)
206     simgrid::smpi::colls::gatherv(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcounts, displs, recvtype,
207                                   root, comm);
208   else
209     simgrid::smpi::colls::igatherv(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcounts, displs, recvtype,
210                                    root, comm, request);
211
212   TRACE_smpi_comm_out(rank);
213   smpi_bench_begin();
214   return MPI_SUCCESS;
215 }
216
217 int PMPI_Allgather(const void *sendbuf, int sendcount, MPI_Datatype sendtype,
218                    void *recvbuf, int recvcount, MPI_Datatype recvtype, MPI_Comm comm){
219   return PMPI_Iallgather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm, MPI_REQUEST_IGNORED);
220 }
221
222 int PMPI_Iallgather(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
223                     MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request)
224 {
225   CHECK_ARGS(comm == MPI_COMM_NULL, MPI_ERR_COMM, "Iallgather: the communicator cannot be MPI_COMM_NULL");
226   CHECK_ARGS(recvbuf == nullptr && recvcount > 0, MPI_ERR_BUFFER, "Iallgather: param 4 recvbuf cannot be NULL");
227   CHECK_ARGS(sendbuf == nullptr && sendcount > 0, MPI_ERR_BUFFER,
228              "Iallgather: param 1 sendbuf cannot be NULL when sendcound > 0");
229   CHECK_ARGS((sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL), MPI_ERR_TYPE,
230              "Iallgather: param 3 sendtype cannot be MPI_DATATYPE_NULL when sendbuff is not MPI_IN_PLACE");
231   CHECK_ARGS(recvtype == MPI_DATATYPE_NULL, MPI_ERR_TYPE, "Iallgather: param 6 recvtype cannot be MPI_DATATYPE_NULL");
232   CHECK_ARGS(recvcount < 0, MPI_ERR_COUNT, "Iallgather: param 5 recvcount cannot be negative");
233   CHECK_ARGS((sendbuf != MPI_IN_PLACE) && (sendcount < 0), MPI_ERR_COUNT,
234              "Iallgather: param 2 sendcount cannot be negative when sendbuf is not MPI_IN_PLACE");
235   CHECK_ARGS(request == nullptr, MPI_ERR_ARG, "Iallgather: param 8 request cannot be NULL");
236
237   smpi_bench_end();
238   if (sendbuf == MPI_IN_PLACE) {
239     sendbuf   = static_cast<char*>(recvbuf) + recvtype->get_extent() * recvcount * comm->rank();
240     sendcount = recvcount;
241     sendtype  = recvtype;
242   }
243   int rank = simgrid::s4u::this_actor::get_pid();
244
245   TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Allgather" : "PMPI_Iallggather",
246                      new simgrid::instr::CollTIData(
247                          request == MPI_REQUEST_IGNORED ? "allgather" : "iallgather", -1, -1.0,
248                          sendtype->is_replayable() ? sendcount : sendcount * sendtype->size(),
249                          recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(),
250                          simgrid::smpi::Datatype::encode(sendtype), simgrid::smpi::Datatype::encode(recvtype)));
251   if (request == MPI_REQUEST_IGNORED)
252     simgrid::smpi::colls::allgather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
253   else
254     simgrid::smpi::colls::iallgather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm, request);
255
256   TRACE_smpi_comm_out(rank);
257   smpi_bench_begin();
258   return MPI_SUCCESS;
259 }
260
261 int PMPI_Allgatherv(const void *sendbuf, int sendcount, MPI_Datatype sendtype,
262                    void *recvbuf, const int *recvcounts, const int *displs, MPI_Datatype recvtype, MPI_Comm comm){
263   return PMPI_Iallgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm, MPI_REQUEST_IGNORED);
264 }
265
266 int PMPI_Iallgatherv(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, const int* recvcounts, const int* displs,
267                      MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request)
268 {
269   if (comm == MPI_COMM_NULL)
270     return MPI_ERR_COMM;
271   if (sendbuf == nullptr && sendcount > 0)
272     return MPI_ERR_BUFFER;
273   if (((sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) || (recvtype == MPI_DATATYPE_NULL))
274     return MPI_ERR_TYPE;
275   if ((sendbuf != MPI_IN_PLACE) && (sendcount < 0))
276     return MPI_ERR_COUNT;
277   if (recvcounts == nullptr || displs == nullptr)
278     return MPI_ERR_ARG;
279   if (request == nullptr)
280     return MPI_ERR_ARG;
281
282   for (int i = 0; i < comm->size(); i++) {
283     if (recvcounts[i] < 0)
284       return MPI_ERR_COUNT;
285     else if (recvcounts[i] > 0 && recvbuf == nullptr)
286       return MPI_ERR_BUFFER;
287   }
288
289   smpi_bench_end();
290   if (sendbuf == MPI_IN_PLACE) {
291     sendbuf   = static_cast<char*>(recvbuf) + recvtype->get_extent() * displs[comm->rank()];
292     sendcount = recvcounts[comm->rank()];
293     sendtype  = recvtype;
294   }
295   int rank         = simgrid::s4u::this_actor::get_pid();
296   int dt_size_recv = recvtype->is_replayable() ? 1 : recvtype->size();
297
298   std::vector<int>* trace_recvcounts = new std::vector<int>;
299   for (int i = 0; i < comm->size(); i++) { // copy data to avoid bad free
300     trace_recvcounts->push_back(recvcounts[i] * dt_size_recv);
301   }
302
303   TRACE_smpi_comm_in(
304       rank, request == MPI_REQUEST_IGNORED ? "PMPI_Allgatherv" : "PMPI_Iallgatherv",
305       new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "allgatherv" : "iallgatherv", -1,
306                                         sendtype->is_replayable() ? sendcount : sendcount * sendtype->size(), nullptr,
307                                         dt_size_recv, trace_recvcounts, simgrid::smpi::Datatype::encode(sendtype),
308                                         simgrid::smpi::Datatype::encode(recvtype)));
309   if (request == MPI_REQUEST_IGNORED)
310     simgrid::smpi::colls::allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm);
311   else
312     simgrid::smpi::colls::iallgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm,
313                                       request);
314
315   TRACE_smpi_comm_out(rank);
316   smpi_bench_begin();
317   return MPI_SUCCESS;
318 }
319
320 int PMPI_Scatter(const void *sendbuf, int sendcount, MPI_Datatype sendtype,
321                 void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm){
322   return PMPI_Iscatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm, MPI_REQUEST_IGNORED);
323 }
324
325 int PMPI_Iscatter(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
326                   MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request* request)
327 {
328   if (comm == MPI_COMM_NULL)
329     return MPI_ERR_COMM;
330   if (((comm->rank() == root) && (sendtype == MPI_DATATYPE_NULL || not sendtype->is_valid())) ||
331       ((recvbuf != MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL || not recvtype->is_valid())))
332     return MPI_ERR_TYPE;
333   if (((comm->rank() == root) && (sendcount < 0)) || ((recvbuf != MPI_IN_PLACE) && (recvcount < 0)))
334     return MPI_ERR_COUNT;
335   if ((sendbuf == recvbuf) || ((comm->rank() == root) && sendcount > 0 && (sendbuf == nullptr)) ||
336       (recvcount > 0 && recvbuf == nullptr))
337     return MPI_ERR_BUFFER;
338   if (root < 0 || root >= comm->size())
339     return MPI_ERR_ROOT;
340   if (request == nullptr)
341     return MPI_ERR_ARG;
342
343   smpi_bench_end();
344   if (recvbuf == MPI_IN_PLACE) {
345     recvtype  = sendtype;
346     recvcount = sendcount;
347   }
348   int rank = simgrid::s4u::this_actor::get_pid();
349
350   TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Scatter" : "PMPI_Iscatter",
351                      new simgrid::instr::CollTIData(
352                          request == MPI_REQUEST_IGNORED ? "scatter" : "iscatter", root, -1.0,
353                          (comm->rank() != root || sendtype->is_replayable()) ? sendcount : sendcount * sendtype->size(),
354                          recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(),
355                          simgrid::smpi::Datatype::encode(sendtype), simgrid::smpi::Datatype::encode(recvtype)));
356   if (request == MPI_REQUEST_IGNORED)
357     simgrid::smpi::colls::scatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm);
358   else
359     simgrid::smpi::colls::iscatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm, request);
360
361   TRACE_smpi_comm_out(rank);
362   smpi_bench_begin();
363   return MPI_SUCCESS;
364 }
365
366 int PMPI_Scatterv(const void *sendbuf, const int *sendcounts, const int *displs,
367                  MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm){
368   return PMPI_Iscatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm, MPI_REQUEST_IGNORED);
369 }
370
371 int PMPI_Iscatterv(const void* sendbuf, const int* sendcounts, const int* displs, MPI_Datatype sendtype, void* recvbuf, int recvcount,
372                    MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request* request)
373 {
374   CHECK_ARGS(comm == MPI_COMM_NULL, MPI_ERR_COMM, "Iscatterv: the communicator cannot be MPI_COMM_NULL");
375   CHECK_ARGS((comm->rank() == root) && (sendcounts == nullptr), MPI_ERR_ARG,
376              "Iscatterv: param 2 sendcounts cannot be NULL on the root rank");
377   CHECK_ARGS((comm->rank() == root) && (displs == nullptr), MPI_ERR_ARG,
378              "Iscatterv: param 3 displs cannot be NULL on the root rank");
379   CHECK_ARGS((comm->rank() == root) && (sendtype == MPI_DATATYPE_NULL), MPI_ERR_TYPE,
380              "Iscatterv: The sendtype cannot be NULL on the root rank");
381   CHECK_ARGS((recvbuf != MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL), MPI_ERR_TYPE,
382              "Iscatterv: the recvtype cannot be NULL when not receiving in place");
383   CHECK_ARGS(request == nullptr, MPI_ERR_ARG, "Iscatterv: param 10 request cannot be NULL");
384   CHECK_ARGS(recvbuf != MPI_IN_PLACE && recvcount < 0, MPI_ERR_COUNT,
385              "Iscatterv: When not receiving in place, the recvcount cannot be negative");
386   CHECK_ARGS(root < 0, MPI_ERR_ROOT, "Iscatterv: root cannot be negative");
387   CHECK_ARGS(root >= comm->size(), MPI_ERR_ROOT, "Iscatterv: root (=%d) is larger than communicator size (=%d)", root,
388              comm->size());
389
390   if (comm->rank() == root) {
391     if (recvbuf == MPI_IN_PLACE) {
392       recvtype  = sendtype;
393       recvcount = sendcounts[comm->rank()];
394     }
395     for (int i = 0; i < comm->size(); i++)
396       CHECK_ARGS(sendcounts[i] < 0, MPI_ERR_COUNT, "Iscatterv: sendcounts[%d]=%d but this cannot be negative", i,
397                  sendcounts[i]);
398   }
399
400   smpi_bench_end();
401
402   int rank         = simgrid::s4u::this_actor::get_pid();
403   int dt_size_send = sendtype->is_replayable() ? 1 : sendtype->size();
404
405   std::vector<int>* trace_sendcounts = new std::vector<int>;
406   if (comm->rank() == root) {
407     for (int i = 0; i < comm->size(); i++) { // copy data to avoid bad free
408       trace_sendcounts->push_back(sendcounts[i] * dt_size_send);
409     }
410   }
411
412   TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Scatterv" : "PMPI_Iscatterv",
413                      new simgrid::instr::VarCollTIData(
414                          request == MPI_REQUEST_IGNORED ? "scatterv" : "iscatterv", root, dt_size_send,
415                          trace_sendcounts, recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(),
416                          nullptr, simgrid::smpi::Datatype::encode(sendtype),
417                          simgrid::smpi::Datatype::encode(recvtype)));
418   if (request == MPI_REQUEST_IGNORED)
419     simgrid::smpi::colls::scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm);
420   else
421     simgrid::smpi::colls::iscatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm,
422                                     request);
423
424   TRACE_smpi_comm_out(rank);
425   smpi_bench_begin();
426   return MPI_SUCCESS;
427 }
428
429 int PMPI_Reduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
430 {
431   return PMPI_Ireduce(sendbuf, recvbuf, count, datatype, op, root, comm, MPI_REQUEST_IGNORED);
432 }
433
434 int PMPI_Ireduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm, MPI_Request* request)
435 {
436   if (comm == MPI_COMM_NULL)
437     return MPI_ERR_COMM;
438   if ((sendbuf == nullptr && count > 0) || ((comm->rank() == root) && recvbuf == nullptr))
439     return MPI_ERR_BUFFER;
440   if (datatype == MPI_DATATYPE_NULL || not datatype->is_valid())
441     return MPI_ERR_TYPE;
442   if (op == MPI_OP_NULL)
443     return MPI_ERR_OP;
444   if (request == nullptr)
445     return MPI_ERR_ARG;
446   if (root < 0 || root >= comm->size())
447     return MPI_ERR_ROOT;
448   if (count < 0)
449     return MPI_ERR_COUNT;
450
451   smpi_bench_end();
452   int rank = simgrid::s4u::this_actor::get_pid();
453
454   TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Reduce" : "PMPI_Ireduce",
455                      new simgrid::instr::CollTIData(request == MPI_REQUEST_IGNORED ? "reduce" : "ireduce", root, 0,
456                                                     datatype->is_replayable() ? count : count * datatype->size(), -1,
457                                                     simgrid::smpi::Datatype::encode(datatype), ""));
458   if (request == MPI_REQUEST_IGNORED)
459     simgrid::smpi::colls::reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
460   else
461     simgrid::smpi::colls::ireduce(sendbuf, recvbuf, count, datatype, op, root, comm, request);
462
463   TRACE_smpi_comm_out(rank);
464   smpi_bench_begin();
465   return MPI_SUCCESS;
466 }
467
468 int PMPI_Reduce_local(const void* inbuf, void* inoutbuf, int count, MPI_Datatype datatype, MPI_Op op)
469 {
470   if (datatype == MPI_DATATYPE_NULL || not datatype->is_valid())
471     return MPI_ERR_TYPE;
472   if (op == MPI_OP_NULL)
473     return MPI_ERR_OP;
474   if (count < 0)
475     return MPI_ERR_COUNT;
476
477   smpi_bench_end();
478   op->apply(inbuf, inoutbuf, &count, datatype);
479   smpi_bench_begin();
480   return MPI_SUCCESS;
481 }
482
483 int PMPI_Allreduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
484 {
485   return PMPI_Iallreduce(sendbuf, recvbuf, count, datatype, op, comm, MPI_REQUEST_IGNORED);
486 }
487
488 int PMPI_Iallreduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request *request)
489 {
490   if (comm == MPI_COMM_NULL)
491     return MPI_ERR_COMM;
492   if ((sendbuf == nullptr && count > 0) || (recvbuf == nullptr))
493     return MPI_ERR_BUFFER;
494   if (datatype == MPI_DATATYPE_NULL || not datatype->is_valid())
495     return MPI_ERR_TYPE;
496   if (count < 0)
497     return MPI_ERR_COUNT;
498   if (op == MPI_OP_NULL)
499     return MPI_ERR_OP;
500   if (request == nullptr)
501     return MPI_ERR_ARG;
502
503   smpi_bench_end();
504   std::unique_ptr<unsigned char[]> tmp_sendbuf;
505   const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, count, datatype);
506
507   int rank = simgrid::s4u::this_actor::get_pid();
508
509   TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Allreduce" : "PMPI_Iallreduce",
510                      new simgrid::instr::CollTIData(request == MPI_REQUEST_IGNORED ? "allreduce" : "iallreduce", -1, 0,
511                                                     datatype->is_replayable() ? count : count * datatype->size(), -1,
512                                                     simgrid::smpi::Datatype::encode(datatype), ""));
513
514   if (request == MPI_REQUEST_IGNORED)
515     simgrid::smpi::colls::allreduce(real_sendbuf, recvbuf, count, datatype, op, comm);
516   else
517     simgrid::smpi::colls::iallreduce(real_sendbuf, recvbuf, count, datatype, op, comm, request);
518
519   TRACE_smpi_comm_out(rank);
520   smpi_bench_begin();
521   return MPI_SUCCESS;
522 }
523
524 int PMPI_Scan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
525 {
526   return PMPI_Iscan(sendbuf, recvbuf, count, datatype, op, comm, MPI_REQUEST_IGNORED);
527 }
528
529 int PMPI_Iscan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request* request)
530 {
531   if (comm == MPI_COMM_NULL)
532     return MPI_ERR_COMM;
533   if (datatype == MPI_DATATYPE_NULL || not datatype->is_valid())
534     return MPI_ERR_TYPE;
535   if (op == MPI_OP_NULL)
536     return MPI_ERR_OP;
537   if (request == nullptr)
538     return MPI_ERR_ARG;
539   if (count < 0)
540     return MPI_ERR_COUNT;
541   if (sendbuf == nullptr || recvbuf == nullptr)
542     return MPI_ERR_BUFFER;
543
544   smpi_bench_end();
545   int rank         = simgrid::s4u::this_actor::get_pid();
546   std::unique_ptr<unsigned char[]> tmp_sendbuf;
547   const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, count, datatype);
548
549   TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Scan" : "PMPI_Iscan",
550                      new simgrid::instr::Pt2PtTIData(request == MPI_REQUEST_IGNORED ? "scan" : "iscan", -1,
551                                                      datatype->is_replayable() ? count : count * datatype->size(),
552                                                      simgrid::smpi::Datatype::encode(datatype)));
553
554   int retval;
555   if (request == MPI_REQUEST_IGNORED)
556     retval = simgrid::smpi::colls::scan(real_sendbuf, recvbuf, count, datatype, op, comm);
557   else
558     retval = simgrid::smpi::colls::iscan(real_sendbuf, recvbuf, count, datatype, op, comm, request);
559
560   TRACE_smpi_comm_out(rank);
561   smpi_bench_begin();
562   return retval;
563 }
564
565 int PMPI_Exscan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
566 {
567   return PMPI_Iexscan(sendbuf, recvbuf, count, datatype, op, comm, MPI_REQUEST_IGNORED);
568 }
569
570 int PMPI_Iexscan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request* request){
571   if (comm == MPI_COMM_NULL)
572     return MPI_ERR_COMM;
573   if (not datatype->is_valid())
574     return MPI_ERR_TYPE;
575   if (op == MPI_OP_NULL)
576     return MPI_ERR_OP;
577   if (request == nullptr)
578     return MPI_ERR_ARG;
579   if (count < 0)
580     return MPI_ERR_COUNT;
581   if (sendbuf == nullptr || recvbuf == nullptr)
582     return MPI_ERR_BUFFER;
583
584   smpi_bench_end();
585   int rank         = simgrid::s4u::this_actor::get_pid();
586   std::unique_ptr<unsigned char[]> tmp_sendbuf;
587   const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, count, datatype);
588
589   TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Exscan" : "PMPI_Iexscan",
590                      new simgrid::instr::Pt2PtTIData(request == MPI_REQUEST_IGNORED ? "exscan" : "iexscan", -1,
591                                                      datatype->is_replayable() ? count : count * datatype->size(),
592                                                      simgrid::smpi::Datatype::encode(datatype)));
593
594   int retval;
595   if (request == MPI_REQUEST_IGNORED)
596     retval = simgrid::smpi::colls::exscan(real_sendbuf, recvbuf, count, datatype, op, comm);
597   else
598     retval = simgrid::smpi::colls::iexscan(real_sendbuf, recvbuf, count, datatype, op, comm, request);
599
600   TRACE_smpi_comm_out(rank);
601   smpi_bench_begin();
602   return retval;
603 }
604
605 int PMPI_Reduce_scatter(const void *sendbuf, void *recvbuf, const int *recvcounts, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
606 {
607   return PMPI_Ireduce_scatter(sendbuf, recvbuf, recvcounts, datatype, op, comm, MPI_REQUEST_IGNORED);
608 }
609
610 int PMPI_Ireduce_scatter(const void *sendbuf, void *recvbuf, const int *recvcounts, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request *request)
611 {
612   if (comm == MPI_COMM_NULL)
613     return MPI_ERR_COMM;
614   if ((sendbuf == nullptr) || (recvbuf == nullptr))
615     return MPI_ERR_BUFFER;
616   if (datatype == MPI_DATATYPE_NULL || not datatype->is_valid())
617     return MPI_ERR_TYPE;
618   if (op == MPI_OP_NULL)
619     return MPI_ERR_OP;
620   if (recvcounts == nullptr)
621     return MPI_ERR_ARG;
622   if (request == nullptr)
623     return MPI_ERR_ARG;
624
625   for (int i = 0; i < comm->size(); i++) { // copy data to avoid bad free
626     if (recvcounts[i] < 0)
627       return MPI_ERR_COUNT;
628   }
629
630   smpi_bench_end();
631   int rank                           = simgrid::s4u::this_actor::get_pid();
632   std::vector<int>* trace_recvcounts = new std::vector<int>;
633   int dt_send_size                   = datatype->is_replayable() ? 1 : datatype->size();
634   int totalcount                     = 0;
635
636   for (int i = 0; i < comm->size(); i++) { // copy data to avoid bad free
637     trace_recvcounts->push_back(recvcounts[i] * dt_send_size);
638     totalcount += recvcounts[i];
639   }
640   std::unique_ptr<unsigned char[]> tmp_sendbuf;
641   const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, totalcount, datatype);
642
643   TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Reduce_scatter" : "PMPI_Ireduce_scatter",
644                      new simgrid::instr::VarCollTIData(
645                          request == MPI_REQUEST_IGNORED ? "reducescatter" : "ireducescatter", -1, dt_send_size, nullptr,
646                          -1, trace_recvcounts, simgrid::smpi::Datatype::encode(datatype), ""));
647
648   if (request == MPI_REQUEST_IGNORED)
649     simgrid::smpi::colls::reduce_scatter(real_sendbuf, recvbuf, recvcounts, datatype, op, comm);
650   else
651     simgrid::smpi::colls::ireduce_scatter(real_sendbuf, recvbuf, recvcounts, datatype, op, comm, request);
652
653   TRACE_smpi_comm_out(rank);
654   smpi_bench_begin();
655   return MPI_SUCCESS;
656 }
657
658 int PMPI_Reduce_scatter_block(const void *sendbuf, void *recvbuf, int recvcount,
659                               MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
660 {
661   return PMPI_Ireduce_scatter_block(sendbuf, recvbuf, recvcount, datatype, op, comm, MPI_REQUEST_IGNORED);
662 }
663
664 int PMPI_Ireduce_scatter_block(const void* sendbuf, void* recvbuf, int recvcount, MPI_Datatype datatype, MPI_Op op,
665                                MPI_Comm comm, MPI_Request* request)
666 {
667   if (comm == MPI_COMM_NULL)
668     return MPI_ERR_COMM;
669   if (not datatype->is_valid())
670     return MPI_ERR_TYPE;
671   if (op == MPI_OP_NULL)
672     return MPI_ERR_OP;
673   if (recvcount < 0)
674     return MPI_ERR_ARG;
675   if (request == nullptr)
676     return MPI_ERR_ARG;
677
678   smpi_bench_end();
679   int count = comm->size();
680
681   int rank                           = simgrid::s4u::this_actor::get_pid();
682   int dt_send_size                   = datatype->is_replayable() ? 1 : datatype->size();
683   std::vector<int>* trace_recvcounts = new std::vector<int>(recvcount * dt_send_size); // copy data to avoid bad free
684   std::unique_ptr<unsigned char[]> tmp_sendbuf;
685   const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, recvcount * count, datatype);
686
687   TRACE_smpi_comm_in(
688       rank, request == MPI_REQUEST_IGNORED ? "PMPI_Reduce_scatter_block" : "PMPI_Ireduce_scatter_block",
689       new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "reducescatter" : "ireducescatter", -1, 0,
690                                         nullptr, -1, trace_recvcounts, simgrid::smpi::Datatype::encode(datatype), ""));
691
692   int* recvcounts = new int[count];
693   for (int i      = 0; i < count; i++)
694     recvcounts[i] = recvcount;
695   if (request == MPI_REQUEST_IGNORED)
696     simgrid::smpi::colls::reduce_scatter(real_sendbuf, recvbuf, recvcounts, datatype, op, comm);
697   else
698     simgrid::smpi::colls::ireduce_scatter(real_sendbuf, recvbuf, recvcounts, datatype, op, comm, request);
699   delete[] recvcounts;
700
701   TRACE_smpi_comm_out(rank);
702   smpi_bench_begin();
703   return MPI_SUCCESS;
704 }
705
706 int PMPI_Alltoall(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
707                   MPI_Datatype recvtype, MPI_Comm comm){
708   return PMPI_Ialltoall(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm, MPI_REQUEST_IGNORED);
709 }
710
711 int PMPI_Ialltoall(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
712                    MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request)
713 {
714   if (comm == MPI_COMM_NULL)
715     return MPI_ERR_COMM;
716   if ((sendbuf == nullptr && sendcount > 0) || (recvbuf == nullptr && recvcount > 0))
717     return MPI_ERR_BUFFER;
718   if ((sendbuf != MPI_IN_PLACE && sendtype == MPI_DATATYPE_NULL) || recvtype == MPI_DATATYPE_NULL)
719     return MPI_ERR_TYPE;
720   if ((sendbuf != MPI_IN_PLACE && sendcount < 0) || recvcount < 0)
721     return MPI_ERR_COUNT;
722   if (request == nullptr)
723     return MPI_ERR_ARG;
724
725   smpi_bench_end();
726   int rank                 = simgrid::s4u::this_actor::get_pid();
727   int real_sendcount         = sendcount;
728   MPI_Datatype real_sendtype = sendtype;
729   
730   std::unique_ptr<unsigned char[]> tmp_sendbuf;
731   const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, recvcount * comm->size(), recvtype);
732
733   if (sendbuf == MPI_IN_PLACE) {
734     real_sendcount = recvcount;
735     real_sendtype  = recvtype;
736   }
737
738   TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Alltoall" : "PMPI_Ialltoall",
739                      new simgrid::instr::CollTIData(
740                          request == MPI_REQUEST_IGNORED ? "alltoall" : "ialltoall", -1, -1.0,
741                          real_sendtype->is_replayable() ? real_sendcount : real_sendcount * real_sendtype->size(),
742                          recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(),
743                          simgrid::smpi::Datatype::encode(real_sendtype), simgrid::smpi::Datatype::encode(recvtype)));
744   int retval;
745   if (request == MPI_REQUEST_IGNORED)
746     retval =
747         simgrid::smpi::colls::alltoall(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcount, recvtype, comm);
748   else
749     retval = simgrid::smpi::colls::ialltoall(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcount, recvtype,
750                                              comm, request);
751
752   TRACE_smpi_comm_out(rank);
753   smpi_bench_begin();
754   return retval;
755 }
756
757 int PMPI_Alltoallv(const void* sendbuf, const int* sendcounts, const int* senddisps, MPI_Datatype sendtype, void* recvbuf,
758                    const int* recvcounts, const int* recvdisps, MPI_Datatype recvtype, MPI_Comm comm)
759 {
760   return PMPI_Ialltoallv(sendbuf, sendcounts, senddisps, sendtype, recvbuf, recvcounts, recvdisps, recvtype, comm, MPI_REQUEST_IGNORED);
761 }
762
763 int PMPI_Ialltoallv(const void* sendbuf, const int* sendcounts, const int* senddisps, MPI_Datatype sendtype, void* recvbuf,
764                     const int* recvcounts, const int* recvdisps, MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request)
765 {
766   if (comm == MPI_COMM_NULL)
767     return MPI_ERR_COMM;
768   if (sendbuf == nullptr || recvbuf == nullptr)
769     return MPI_ERR_BUFFER;
770   if ((sendbuf != MPI_IN_PLACE && sendtype == MPI_DATATYPE_NULL) || recvtype == MPI_DATATYPE_NULL)
771     return MPI_ERR_TYPE;
772   if ((sendbuf != MPI_IN_PLACE && (sendcounts == nullptr || senddisps == nullptr)) || recvcounts == nullptr ||
773       recvdisps == nullptr)
774     return MPI_ERR_ARG;
775   if (request == nullptr)
776     return MPI_ERR_ARG;
777
778   int rank = simgrid::s4u::this_actor::get_pid();
779   int size = comm->size();
780   for (int i = 0; i < size; i++) {
781     if (recvcounts[i] < 0 || (sendbuf != MPI_IN_PLACE && sendcounts[i] < 0))
782       return MPI_ERR_COUNT;
783   }
784
785   smpi_bench_end();
786   int send_size                      = 0;
787   int recv_size                      = 0;
788   std::vector<int>* trace_sendcounts = new std::vector<int>;
789   std::vector<int>* trace_recvcounts = new std::vector<int>;
790   int dt_size_recv                   = recvtype->size();
791
792   const int* real_sendcounts = sendcounts;
793   const int* real_senddisps  = senddisps;
794   MPI_Datatype real_sendtype = sendtype;
795   int maxsize              = 0;
796   for (int i = 0; i < size; i++) { // copy data to avoid bad free
797     recv_size += recvcounts[i] * dt_size_recv;
798     trace_recvcounts->push_back(recvcounts[i] * dt_size_recv);
799     if (((recvdisps[i] + recvcounts[i]) * dt_size_recv) > maxsize)
800       maxsize = (recvdisps[i] + recvcounts[i]) * dt_size_recv;
801   }
802
803   std::unique_ptr<unsigned char[]> tmp_sendbuf;
804   std::unique_ptr<int[]> tmp_sendcounts;
805   std::unique_ptr<int[]> tmp_senddisps;
806   const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, maxsize, MPI_CHAR);
807   if (sendbuf == MPI_IN_PLACE) {
808     tmp_sendcounts.reset(new int[size]);
809     std::copy(recvcounts, recvcounts + size, tmp_sendcounts.get());
810     real_sendcounts = tmp_sendcounts.get();
811     tmp_senddisps.reset(new int[size]);
812     std::copy(recvdisps, recvdisps + size, tmp_senddisps.get());
813     real_senddisps = tmp_senddisps.get();
814     real_sendtype  = recvtype;
815   }
816
817   int dt_size_send = real_sendtype->size();
818
819   for (int i = 0; i < size; i++) { // copy data to avoid bad free
820     send_size += real_sendcounts[i] * dt_size_send;
821     trace_sendcounts->push_back(real_sendcounts[i] * dt_size_send);
822   }
823
824   TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Alltoallv" : "PMPI_Ialltoallv",
825                      new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "alltoallv" : "ialltoallv", -1,
826                                                        send_size, trace_sendcounts, recv_size, trace_recvcounts,
827                                                        simgrid::smpi::Datatype::encode(real_sendtype),
828                                                        simgrid::smpi::Datatype::encode(recvtype)));
829
830   int retval;
831   if (request == MPI_REQUEST_IGNORED)
832     retval = simgrid::smpi::colls::alltoallv(real_sendbuf, real_sendcounts, real_senddisps, real_sendtype, recvbuf,
833                                              recvcounts, recvdisps, recvtype, comm);
834   else
835     retval = simgrid::smpi::colls::ialltoallv(real_sendbuf, real_sendcounts, real_senddisps, real_sendtype, recvbuf,
836                                               recvcounts, recvdisps, recvtype, comm, request);
837
838   TRACE_smpi_comm_out(rank);
839   smpi_bench_begin();
840   return retval;
841 }
842
843 int PMPI_Alltoallw(const void* sendbuf, const int* sendcounts, const int* senddisps, const MPI_Datatype* sendtypes, void* recvbuf,
844                    const int* recvcounts, const int* recvdisps, const MPI_Datatype* recvtypes, MPI_Comm comm)
845 {
846   return PMPI_Ialltoallw(sendbuf, sendcounts, senddisps, sendtypes, recvbuf, recvcounts, recvdisps, recvtypes, comm, MPI_REQUEST_IGNORED);
847 }
848
849 int PMPI_Ialltoallw(const void* sendbuf, const int* sendcounts, const int* senddisps, const MPI_Datatype* sendtypes, void* recvbuf,
850                     const int* recvcounts, const int* recvdisps, const MPI_Datatype* recvtypes, MPI_Comm comm, MPI_Request* request)
851 {
852   if (comm == MPI_COMM_NULL)
853     return MPI_ERR_COMM;
854   if (sendbuf == nullptr || recvbuf == nullptr)
855     return MPI_ERR_BUFFER;
856   if ((sendbuf != MPI_IN_PLACE && sendtypes == nullptr) || recvtypes == nullptr)
857     return MPI_ERR_TYPE;
858   if ((sendbuf != MPI_IN_PLACE && (sendcounts == nullptr || senddisps == nullptr)) || recvcounts == nullptr ||
859       recvdisps == nullptr)
860     return MPI_ERR_ARG;
861   if (request == nullptr)
862     return MPI_ERR_ARG;
863
864   smpi_bench_end();
865   int rank = simgrid::s4u::this_actor::get_pid();
866   int size = comm->size();
867   for (int i = 0; i < size; i++) {
868     if (recvcounts[i] < 0 || (sendbuf != MPI_IN_PLACE && sendcounts[i] < 0))
869       return MPI_ERR_COUNT;
870   }
871   int send_size                      = 0;
872   int recv_size                      = 0;
873   std::vector<int>* trace_sendcounts = new std::vector<int>;
874   std::vector<int>* trace_recvcounts = new std::vector<int>;
875
876   const int* real_sendcounts         = sendcounts;
877   const int* real_senddisps          = senddisps;
878   const MPI_Datatype* real_sendtypes = sendtypes;
879   unsigned long maxsize      = 0;
880   for (int i = 0; i < size; i++) { // copy data to avoid bad free
881     if (recvtypes[i] == MPI_DATATYPE_NULL) {
882       delete trace_recvcounts;
883       delete trace_sendcounts;
884       return MPI_ERR_TYPE;
885     }
886     recv_size += recvcounts[i] * recvtypes[i]->size();
887     trace_recvcounts->push_back(recvcounts[i] * recvtypes[i]->size());
888     if ((recvdisps[i] + (recvcounts[i] * recvtypes[i]->size())) > maxsize)
889       maxsize = recvdisps[i] + (recvcounts[i] * recvtypes[i]->size());
890   }
891
892   std::unique_ptr<unsigned char[]> tmp_sendbuf;
893   std::unique_ptr<int[]> tmp_sendcounts;
894   std::unique_ptr<int[]> tmp_senddisps;
895   std::unique_ptr<MPI_Datatype[]> tmp_sendtypes;
896   const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, maxsize, MPI_CHAR);
897   if (sendbuf == MPI_IN_PLACE) {
898     tmp_sendcounts.reset(new int[size]);
899     std::copy(recvcounts, recvcounts + size, tmp_sendcounts.get());
900     real_sendcounts = tmp_sendcounts.get();
901     tmp_senddisps.reset(new int[size]);
902     std::copy(recvdisps, recvdisps + size, tmp_senddisps.get());
903     real_senddisps = tmp_senddisps.get();
904     tmp_sendtypes.reset(new MPI_Datatype[size]);
905     std::copy(recvtypes, recvtypes + size, tmp_sendtypes.get());
906     real_sendtypes = tmp_sendtypes.get();
907   }
908
909   for (int i = 0; i < size; i++) { // copy data to avoid bad free
910     send_size += real_sendcounts[i] * real_sendtypes[i]->size();
911     trace_sendcounts->push_back(real_sendcounts[i] * real_sendtypes[i]->size());
912   }
913
914   TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Alltoallw" : "PMPI_Ialltoallw",
915                      new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "alltoallv" : "ialltoallv", -1,
916                                                        send_size, trace_sendcounts, recv_size, trace_recvcounts,
917                                                        simgrid::smpi::Datatype::encode(real_sendtypes[0]),
918                                                        simgrid::smpi::Datatype::encode(recvtypes[0])));
919
920   int retval;
921   if (request == MPI_REQUEST_IGNORED)
922     retval = simgrid::smpi::colls::alltoallw(real_sendbuf, real_sendcounts, real_senddisps, real_sendtypes, recvbuf,
923                                              recvcounts, recvdisps, recvtypes, comm);
924   else
925     retval = simgrid::smpi::colls::ialltoallw(real_sendbuf, real_sendcounts, real_senddisps, real_sendtypes, recvbuf,
926                                               recvcounts, recvdisps, recvtypes, comm, request);
927
928   TRACE_smpi_comm_out(rank);
929   smpi_bench_begin();
930   return retval;
931 }