Logo AND Algorithmique Numérique Distribuée

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