Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Merge branch 'trace_smpi_execute_flops' into 'master'
[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 && sendcount > 0) || ((comm->rank() == root) && recvbuf == nullptr && recvcount > 0))
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   CHECK_ARGS(comm == MPI_COMM_NULL, MPI_ERR_COMM, "Iallgather: the communicator cannot be MPI_COMM_NULL");
217   CHECK_ARGS(recvbuf == nullptr && recvcount > 0, MPI_ERR_BUFFER, "Iallgather: param 4 recvbuf cannot be NULL");
218   CHECK_ARGS(sendbuf == nullptr && sendcount > 0, MPI_ERR_BUFFER,
219              "Iallgather: param 1 sendbuf cannot be NULL when sendcound > 0");
220   CHECK_ARGS((sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL), MPI_ERR_TYPE,
221              "Iallgather: param 3 sendtype cannot be MPI_DATATYPE_NULL when sendbuff is not MPI_IN_PLACE");
222   CHECK_ARGS(recvtype == MPI_DATATYPE_NULL, MPI_ERR_TYPE, "Iallgather: param 6 recvtype cannot be MPI_DATATYPE_NULL");
223   CHECK_ARGS(recvcount < 0, MPI_ERR_COUNT, "Iallgather: param 5 recvcount cannot be negative");
224   CHECK_ARGS((sendbuf != MPI_IN_PLACE) && (sendcount < 0), MPI_ERR_COUNT,
225              "Iallgather: param 2 sendcount cannot be negative when sendbuf is not MPI_IN_PLACE");
226   CHECK_ARGS(request == nullptr, MPI_ERR_ARG, "Iallgather: param 8 request cannot be NULL");
227
228   smpi_bench_end();
229   if (sendbuf == MPI_IN_PLACE) {
230     sendbuf   = static_cast<char*>(recvbuf) + recvtype->get_extent() * recvcount * comm->rank();
231     sendcount = recvcount;
232     sendtype  = recvtype;
233   }
234   int rank = simgrid::s4u::this_actor::get_pid();
235
236   TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Allgather" : "PMPI_Iallggather",
237                      new simgrid::instr::CollTIData(
238                          request == MPI_REQUEST_IGNORED ? "allgather" : "iallgather", -1, -1.0,
239                          sendtype->is_replayable() ? sendcount : sendcount * sendtype->size(),
240                          recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(),
241                          simgrid::smpi::Datatype::encode(sendtype), simgrid::smpi::Datatype::encode(recvtype)));
242   if (request == MPI_REQUEST_IGNORED)
243     simgrid::smpi::Colls::allgather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
244   else
245     simgrid::smpi::Colls::iallgather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm, request);
246
247   TRACE_smpi_comm_out(rank);
248   smpi_bench_begin();
249   return MPI_SUCCESS;
250 }
251
252 int PMPI_Allgatherv(const void *sendbuf, int sendcount, MPI_Datatype sendtype,
253                    void *recvbuf, const int *recvcounts, const int *displs, MPI_Datatype recvtype, MPI_Comm comm){
254   return PMPI_Iallgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm, MPI_REQUEST_IGNORED);
255 }
256
257 int PMPI_Iallgatherv(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, const int* recvcounts, const int* displs,
258                      MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request)
259 {
260   if (comm == MPI_COMM_NULL)
261     return MPI_ERR_COMM;
262   if (sendbuf == nullptr && sendcount > 0)
263     return MPI_ERR_BUFFER;
264   if (((sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) || (recvtype == MPI_DATATYPE_NULL))
265     return MPI_ERR_TYPE;
266   if ((sendbuf != MPI_IN_PLACE) && (sendcount < 0))
267     return MPI_ERR_COUNT;
268   if (recvcounts == nullptr || displs == nullptr)
269     return MPI_ERR_ARG;
270   if (request == nullptr)
271     return MPI_ERR_ARG;
272
273   for (int i = 0; i < comm->size(); i++) {
274     if (recvcounts[i] < 0)
275       return MPI_ERR_COUNT;
276     else if (recvcounts[i] > 0 && recvbuf == nullptr)
277       return MPI_ERR_BUFFER;
278   }
279
280   smpi_bench_end();
281   if (sendbuf == MPI_IN_PLACE) {
282     sendbuf   = static_cast<char*>(recvbuf) + recvtype->get_extent() * displs[comm->rank()];
283     sendcount = recvcounts[comm->rank()];
284     sendtype  = recvtype;
285   }
286   int rank         = simgrid::s4u::this_actor::get_pid();
287   int dt_size_recv = recvtype->is_replayable() ? 1 : recvtype->size();
288
289   std::vector<int>* trace_recvcounts = new std::vector<int>;
290   for (int i = 0; i < comm->size(); i++) { // copy data to avoid bad free
291     trace_recvcounts->push_back(recvcounts[i] * dt_size_recv);
292   }
293
294   TRACE_smpi_comm_in(
295       rank, request == MPI_REQUEST_IGNORED ? "PMPI_Allgatherv" : "PMPI_Iallgatherv",
296       new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "allgatherv" : "iallgatherv", -1,
297                                         sendtype->is_replayable() ? sendcount : sendcount * sendtype->size(), nullptr,
298                                         dt_size_recv, trace_recvcounts, simgrid::smpi::Datatype::encode(sendtype),
299                                         simgrid::smpi::Datatype::encode(recvtype)));
300   if (request == MPI_REQUEST_IGNORED)
301     simgrid::smpi::Colls::allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm);
302   else
303     simgrid::smpi::Colls::iallgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm,
304                                       request);
305
306   TRACE_smpi_comm_out(rank);
307   smpi_bench_begin();
308   return MPI_SUCCESS;
309 }
310
311 int PMPI_Scatter(const void *sendbuf, int sendcount, MPI_Datatype sendtype,
312                 void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm){
313   return PMPI_Iscatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm, MPI_REQUEST_IGNORED);
314 }
315
316 int PMPI_Iscatter(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
317                   MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request* request)
318 {
319   if (comm == MPI_COMM_NULL)
320     return MPI_ERR_COMM;
321   if (((comm->rank() == root) && (sendtype == MPI_DATATYPE_NULL || not sendtype->is_valid())) ||
322       ((recvbuf != MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL || not recvtype->is_valid())))
323     return MPI_ERR_TYPE;
324   if (((comm->rank() == root) && (sendcount < 0)) || ((recvbuf != MPI_IN_PLACE) && (recvcount < 0)))
325     return MPI_ERR_COUNT;
326   if ((sendbuf == recvbuf) || ((comm->rank() == root) && sendcount > 0 && (sendbuf == nullptr)) ||
327       (recvcount > 0 && recvbuf == nullptr))
328     return MPI_ERR_BUFFER;
329   if (root < 0 || root >= comm->size())
330     return MPI_ERR_ROOT;
331   if (request == nullptr)
332     return MPI_ERR_ARG;
333
334   smpi_bench_end();
335   if (recvbuf == MPI_IN_PLACE) {
336     recvtype  = sendtype;
337     recvcount = sendcount;
338   }
339   int rank = simgrid::s4u::this_actor::get_pid();
340
341   TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Scatter" : "PMPI_Iscatter",
342                      new simgrid::instr::CollTIData(
343                          request == MPI_REQUEST_IGNORED ? "scatter" : "iscatter", root, -1.0,
344                          (comm->rank() != root || sendtype->is_replayable()) ? sendcount : sendcount * sendtype->size(),
345                          recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(),
346                          simgrid::smpi::Datatype::encode(sendtype), simgrid::smpi::Datatype::encode(recvtype)));
347   if (request == MPI_REQUEST_IGNORED)
348     simgrid::smpi::Colls::scatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm);
349   else
350     simgrid::smpi::Colls::iscatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm, request);
351
352   TRACE_smpi_comm_out(rank);
353   smpi_bench_begin();
354   return MPI_SUCCESS;
355 }
356
357 int PMPI_Scatterv(const void *sendbuf, const int *sendcounts, const int *displs,
358                  MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm){
359   return PMPI_Iscatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm, MPI_REQUEST_IGNORED);
360 }
361
362 int PMPI_Iscatterv(const void* sendbuf, const int* sendcounts, const int* displs, MPI_Datatype sendtype, void* recvbuf, int recvcount,
363                    MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request* request)
364 {
365   CHECK_ARGS(comm == MPI_COMM_NULL, MPI_ERR_COMM, "Iscatterv: the communicator cannot be MPI_COMM_NULL");
366   CHECK_ARGS((comm->rank() == root) && (sendcounts == nullptr), MPI_ERR_ARG,
367              "Iscatterv: param 2 sendcounts cannot be NULL on the root rank");
368   CHECK_ARGS((comm->rank() == root) && (displs == nullptr), MPI_ERR_ARG,
369              "Iscatterv: param 3 displs cannot be NULL on the root rank");
370   CHECK_ARGS((comm->rank() == root) && (sendtype == MPI_DATATYPE_NULL), MPI_ERR_TYPE,
371              "Iscatterv: The sendtype cannot be NULL on the root rank");
372   CHECK_ARGS((recvbuf != MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL), MPI_ERR_TYPE,
373              "Iscatterv: the recvtype cannot be NULL when not receiving in place");
374   CHECK_ARGS(request == nullptr, MPI_ERR_ARG, "Iscatterv: param 10 request cannot be NULL");
375   CHECK_ARGS(recvbuf != MPI_IN_PLACE && recvcount < 0, MPI_ERR_COUNT,
376              "Iscatterv: When not receiving in place, the recvcound cannot be negative");
377   CHECK_ARGS(root < 0, MPI_ERR_ROOT, "Iscatterv: root cannot be negative");
378   CHECK_ARGS(root >= comm->size(), MPI_ERR_ROOT, "Iscatterv: root (=%d) is larger than communicator size (=%d)", root,
379              comm->size());
380
381   if (comm->rank() == root) {
382     if (recvbuf == MPI_IN_PLACE) {
383       recvtype  = sendtype;
384       recvcount = sendcounts[comm->rank()];
385     }
386     for (int i = 0; i < comm->size(); i++)
387       CHECK_ARGS(sendcounts[i] < 0, MPI_ERR_COUNT, "Iscatterv: sendcounts[%d]=%d but this cannot be negative", i,
388                  sendcounts[i]);
389   }
390
391   smpi_bench_end();
392
393   int rank         = simgrid::s4u::this_actor::get_pid();
394   int dt_size_send = sendtype->is_replayable() ? 1 : sendtype->size();
395
396   std::vector<int>* trace_sendcounts = new std::vector<int>;
397   if (comm->rank() == root) {
398     for (int i = 0; i < comm->size(); i++) { // copy data to avoid bad free
399       trace_sendcounts->push_back(sendcounts[i] * dt_size_send);
400     }
401   }
402
403   TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Scatterv" : "PMPI_Iscatterv",
404                      new simgrid::instr::VarCollTIData(
405                          request == MPI_REQUEST_IGNORED ? "scatterv" : "iscatterv", root, dt_size_send,
406                          trace_sendcounts, recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(),
407                          nullptr, simgrid::smpi::Datatype::encode(sendtype),
408                          simgrid::smpi::Datatype::encode(recvtype)));
409   if (request == MPI_REQUEST_IGNORED)
410     simgrid::smpi::Colls::scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm);
411   else
412     simgrid::smpi::Colls::iscatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm,
413                                     request);
414
415   TRACE_smpi_comm_out(rank);
416   smpi_bench_begin();
417   return MPI_SUCCESS;
418 }
419
420 int PMPI_Reduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
421 {
422   return PMPI_Ireduce(sendbuf, recvbuf, count, datatype, op, root, comm, MPI_REQUEST_IGNORED);
423 }
424
425 int PMPI_Ireduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm, MPI_Request* request)
426 {
427   if (comm == MPI_COMM_NULL)
428     return MPI_ERR_COMM;
429   if ((sendbuf == nullptr && count > 0) || ((comm->rank() == root) && recvbuf == nullptr))
430     return MPI_ERR_BUFFER;
431   if (datatype == MPI_DATATYPE_NULL || not datatype->is_valid())
432     return MPI_ERR_TYPE;
433   if (op == MPI_OP_NULL)
434     return MPI_ERR_OP;
435   if (request == nullptr)
436     return MPI_ERR_ARG;
437   if (root < 0 || root >= comm->size())
438     return MPI_ERR_ROOT;
439   if (count < 0)
440     return MPI_ERR_COUNT;
441
442   smpi_bench_end();
443   int rank = simgrid::s4u::this_actor::get_pid();
444
445   TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Reduce" : "PMPI_Ireduce",
446                      new simgrid::instr::CollTIData(request == MPI_REQUEST_IGNORED ? "reduce" : "ireduce", root, 0,
447                                                     datatype->is_replayable() ? count : count * datatype->size(), -1,
448                                                     simgrid::smpi::Datatype::encode(datatype), ""));
449   if (request == MPI_REQUEST_IGNORED)
450     simgrid::smpi::Colls::reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
451   else
452     simgrid::smpi::Colls::ireduce(sendbuf, recvbuf, count, datatype, op, root, comm, request);
453
454   TRACE_smpi_comm_out(rank);
455   smpi_bench_begin();
456   return MPI_SUCCESS;
457 }
458
459 int PMPI_Reduce_local(const void* inbuf, void* inoutbuf, int count, MPI_Datatype datatype, MPI_Op op)
460 {
461   if (datatype == MPI_DATATYPE_NULL || not datatype->is_valid())
462     return MPI_ERR_TYPE;
463   if (op == MPI_OP_NULL)
464     return MPI_ERR_OP;
465   if (count < 0)
466     return MPI_ERR_COUNT;
467
468   smpi_bench_end();
469   op->apply(inbuf, inoutbuf, &count, datatype);
470   smpi_bench_begin();
471   return MPI_SUCCESS;
472 }
473
474 int PMPI_Allreduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
475 {
476   return PMPI_Iallreduce(sendbuf, recvbuf, count, datatype, op, comm, MPI_REQUEST_IGNORED);
477 }
478
479 int PMPI_Iallreduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request *request)
480 {
481   if (comm == MPI_COMM_NULL)
482     return MPI_ERR_COMM;
483   if ((sendbuf == nullptr && count > 0) || (recvbuf == nullptr))
484     return MPI_ERR_BUFFER;
485   if (datatype == MPI_DATATYPE_NULL || not datatype->is_valid())
486     return MPI_ERR_TYPE;
487   if (count < 0)
488     return MPI_ERR_COUNT;
489   if (op == MPI_OP_NULL)
490     return MPI_ERR_OP;
491   if (request == nullptr)
492     return MPI_ERR_ARG;
493
494   smpi_bench_end();
495   const void* real_sendbuf = sendbuf;
496   std::unique_ptr<unsigned char[]> tmp_sendbuf;
497   if (sendbuf == MPI_IN_PLACE) {
498     tmp_sendbuf.reset(new unsigned char[count * datatype->get_extent()]);
499     simgrid::smpi::Datatype::copy(recvbuf, count, datatype, tmp_sendbuf.get(), count, datatype);
500     real_sendbuf = tmp_sendbuf.get();
501   }
502   int rank = simgrid::s4u::this_actor::get_pid();
503
504   TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Allreduce" : "PMPI_Iallreduce",
505                      new simgrid::instr::CollTIData(request == MPI_REQUEST_IGNORED ? "allreduce" : "iallreduce", -1, 0,
506                                                     datatype->is_replayable() ? count : count * datatype->size(), -1,
507                                                     simgrid::smpi::Datatype::encode(datatype), ""));
508
509   if (request == MPI_REQUEST_IGNORED)
510     simgrid::smpi::Colls::allreduce(real_sendbuf, recvbuf, count, datatype, op, comm);
511   else
512     simgrid::smpi::Colls::iallreduce(real_sendbuf, recvbuf, count, datatype, op, comm, request);
513
514   TRACE_smpi_comm_out(rank);
515   smpi_bench_begin();
516   return MPI_SUCCESS;
517 }
518
519 int PMPI_Scan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
520 {
521   return PMPI_Iscan(sendbuf, recvbuf, count, datatype, op, comm, MPI_REQUEST_IGNORED);
522 }
523
524 int PMPI_Iscan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request* request)
525 {
526   if (comm == MPI_COMM_NULL)
527     return MPI_ERR_COMM;
528   if (datatype == MPI_DATATYPE_NULL || not datatype->is_valid())
529     return MPI_ERR_TYPE;
530   if (op == MPI_OP_NULL)
531     return MPI_ERR_OP;
532   if (request == nullptr)
533     return MPI_ERR_ARG;
534   if (count < 0)
535     return MPI_ERR_COUNT;
536   if (sendbuf == nullptr || recvbuf == nullptr)
537     return MPI_ERR_BUFFER;
538
539   smpi_bench_end();
540   int rank         = simgrid::s4u::this_actor::get_pid();
541   const void* real_sendbuf = sendbuf;
542   std::unique_ptr<unsigned char[]> tmp_sendbuf;
543   if (sendbuf == MPI_IN_PLACE) {
544     tmp_sendbuf.reset(new unsigned char[count * datatype->size()]);
545     real_sendbuf = memcpy(tmp_sendbuf.get(), recvbuf, count * datatype->size());
546   }
547   TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Scan" : "PMPI_Iscan",
548                      new simgrid::instr::Pt2PtTIData(request == MPI_REQUEST_IGNORED ? "scan" : "iscan", -1,
549                                                      datatype->is_replayable() ? count : count * datatype->size(),
550                                                      simgrid::smpi::Datatype::encode(datatype)));
551
552   int retval;
553   if (request == MPI_REQUEST_IGNORED)
554     retval = simgrid::smpi::Colls::scan(real_sendbuf, recvbuf, count, datatype, op, comm);
555   else
556     retval = simgrid::smpi::Colls::iscan(real_sendbuf, recvbuf, count, datatype, op, comm, request);
557
558   TRACE_smpi_comm_out(rank);
559   smpi_bench_begin();
560   return retval;
561 }
562
563 int PMPI_Exscan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
564 {
565   return PMPI_Iexscan(sendbuf, recvbuf, count, datatype, op, comm, MPI_REQUEST_IGNORED);
566 }
567
568 int PMPI_Iexscan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request* request){
569   if (comm == MPI_COMM_NULL)
570     return MPI_ERR_COMM;
571   if (not datatype->is_valid())
572     return MPI_ERR_TYPE;
573   if (op == MPI_OP_NULL)
574     return MPI_ERR_OP;
575   if (request == nullptr)
576     return MPI_ERR_ARG;
577   if (count < 0)
578     return MPI_ERR_COUNT;
579   if (sendbuf == nullptr || recvbuf == nullptr)
580     return MPI_ERR_BUFFER;
581
582   smpi_bench_end();
583   int rank         = simgrid::s4u::this_actor::get_pid();
584   const void* real_sendbuf = sendbuf;
585   std::unique_ptr<unsigned char[]> tmp_sendbuf;
586   if (sendbuf == MPI_IN_PLACE) {
587     tmp_sendbuf.reset(new unsigned char[count * datatype->size()]);
588     real_sendbuf = memcpy(tmp_sendbuf.get(), recvbuf, count * datatype->size());
589   }
590
591   TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Exscan" : "PMPI_Iexscan",
592                      new simgrid::instr::Pt2PtTIData(request == MPI_REQUEST_IGNORED ? "exscan" : "iexscan", -1,
593                                                      datatype->is_replayable() ? count : count * datatype->size(),
594                                                      simgrid::smpi::Datatype::encode(datatype)));
595
596   int retval;
597   if (request == MPI_REQUEST_IGNORED)
598     retval = simgrid::smpi::Colls::exscan(real_sendbuf, recvbuf, count, datatype, op, comm);
599   else
600     retval = simgrid::smpi::Colls::iexscan(real_sendbuf, recvbuf, count, datatype, op, comm, request);
601
602   TRACE_smpi_comm_out(rank);
603   smpi_bench_begin();
604   return retval;
605 }
606
607 int PMPI_Reduce_scatter(const void *sendbuf, void *recvbuf, const int *recvcounts, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
608 {
609   return PMPI_Ireduce_scatter(sendbuf, recvbuf, recvcounts, datatype, op, comm, MPI_REQUEST_IGNORED);
610 }
611
612 int PMPI_Ireduce_scatter(const void *sendbuf, void *recvbuf, const int *recvcounts, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request *request)
613 {
614   if (comm == MPI_COMM_NULL)
615     return MPI_ERR_COMM;
616   if ((sendbuf == nullptr) || (recvbuf == nullptr))
617     return MPI_ERR_BUFFER;
618   if (datatype == MPI_DATATYPE_NULL || not datatype->is_valid())
619     return MPI_ERR_TYPE;
620   if (op == MPI_OP_NULL)
621     return MPI_ERR_OP;
622   if (recvcounts == nullptr)
623     return MPI_ERR_ARG;
624   if (request == nullptr)
625     return MPI_ERR_ARG;
626
627   for (int i = 0; i < comm->size(); i++) { // copy data to avoid bad free
628     if (recvcounts[i] < 0)
629       return MPI_ERR_COUNT;
630   }
631
632   smpi_bench_end();
633   int rank                           = simgrid::s4u::this_actor::get_pid();
634   std::vector<int>* trace_recvcounts = new std::vector<int>;
635   int dt_send_size                   = datatype->is_replayable() ? 1 : datatype->size();
636   int totalcount                     = 0;
637
638   for (int i = 0; i < comm->size(); i++) { // copy data to avoid bad free
639     trace_recvcounts->push_back(recvcounts[i] * dt_send_size);
640     totalcount += recvcounts[i];
641   }
642
643   const void* real_sendbuf = sendbuf;
644   std::unique_ptr<unsigned char[]> tmp_sendbuf;
645   if (sendbuf == MPI_IN_PLACE) {
646     tmp_sendbuf.reset(new unsigned char[totalcount * datatype->size()]);
647     real_sendbuf = memcpy(tmp_sendbuf.get(), recvbuf, totalcount * datatype->size());
648   }
649
650   TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Reduce_scatter" : "PMPI_Ireduce_scatter",
651                      new simgrid::instr::VarCollTIData(
652                          request == MPI_REQUEST_IGNORED ? "reducescatter" : "ireducescatter", -1, dt_send_size, nullptr,
653                          -1, trace_recvcounts, simgrid::smpi::Datatype::encode(datatype), ""));
654
655   if (request == MPI_REQUEST_IGNORED)
656     simgrid::smpi::Colls::reduce_scatter(real_sendbuf, recvbuf, recvcounts, datatype, op, comm);
657   else
658     simgrid::smpi::Colls::ireduce_scatter(real_sendbuf, recvbuf, recvcounts, datatype, op, comm, request);
659
660   TRACE_smpi_comm_out(rank);
661   smpi_bench_begin();
662   return MPI_SUCCESS;
663 }
664
665 int PMPI_Reduce_scatter_block(const void *sendbuf, void *recvbuf, int recvcount,
666                               MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
667 {
668   return PMPI_Ireduce_scatter_block(sendbuf, recvbuf, recvcount, datatype, op, comm, MPI_REQUEST_IGNORED);
669 }
670
671 int PMPI_Ireduce_scatter_block(const void* sendbuf, void* recvbuf, int recvcount, MPI_Datatype datatype, MPI_Op op,
672                                MPI_Comm comm, MPI_Request* request)
673 {
674   if (comm == MPI_COMM_NULL)
675     return MPI_ERR_COMM;
676   if (not datatype->is_valid())
677     return MPI_ERR_TYPE;
678   if (op == MPI_OP_NULL)
679     return MPI_ERR_OP;
680   if (recvcount < 0)
681     return MPI_ERR_ARG;
682   if (request == nullptr)
683     return MPI_ERR_ARG;
684
685   smpi_bench_end();
686   int count = comm->size();
687
688   int rank                           = simgrid::s4u::this_actor::get_pid();
689   int dt_send_size                   = datatype->is_replayable() ? 1 : datatype->size();
690   std::vector<int>* trace_recvcounts = new std::vector<int>(recvcount * dt_send_size); // copy data to avoid bad free
691
692   const void* real_sendbuf = sendbuf;
693   std::unique_ptr<unsigned char[]> tmp_sendbuf;
694   if (sendbuf == MPI_IN_PLACE) {
695     tmp_sendbuf.reset(new unsigned char[recvcount * count * datatype->size()]);
696     real_sendbuf = memcpy(tmp_sendbuf.get(), recvbuf, recvcount * count * datatype->size());
697   }
698
699   TRACE_smpi_comm_in(
700       rank, request == MPI_REQUEST_IGNORED ? "PMPI_Reduce_scatter_block" : "PMPI_Ireduce_scatter_block",
701       new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "reducescatter" : "ireducescatter", -1, 0,
702                                         nullptr, -1, trace_recvcounts, simgrid::smpi::Datatype::encode(datatype), ""));
703
704   int* recvcounts = new int[count];
705   for (int i      = 0; i < count; i++)
706     recvcounts[i] = recvcount;
707   if (request == MPI_REQUEST_IGNORED)
708     simgrid::smpi::Colls::reduce_scatter(real_sendbuf, recvbuf, recvcounts, datatype, op, comm);
709   else
710     simgrid::smpi::Colls::ireduce_scatter(real_sendbuf, recvbuf, recvcounts, datatype, op, comm, request);
711   delete[] recvcounts;
712
713   TRACE_smpi_comm_out(rank);
714   smpi_bench_begin();
715   return MPI_SUCCESS;
716 }
717
718 int PMPI_Alltoall(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
719                   MPI_Datatype recvtype, MPI_Comm comm){
720   return PMPI_Ialltoall(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm, MPI_REQUEST_IGNORED);
721 }
722
723 int PMPI_Ialltoall(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
724                    MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request)
725 {
726   if (comm == MPI_COMM_NULL)
727     return MPI_ERR_COMM;
728   if ((sendbuf == nullptr && sendcount > 0) || (recvbuf == nullptr && recvcount > 0))
729     return MPI_ERR_BUFFER;
730   if ((sendbuf != MPI_IN_PLACE && sendtype == MPI_DATATYPE_NULL) || recvtype == MPI_DATATYPE_NULL)
731     return MPI_ERR_TYPE;
732   if ((sendbuf != MPI_IN_PLACE && sendcount < 0) || recvcount < 0)
733     return MPI_ERR_COUNT;
734   if (request == nullptr)
735     return MPI_ERR_ARG;
736
737   smpi_bench_end();
738   int rank                 = simgrid::s4u::this_actor::get_pid();
739   const void* real_sendbuf = sendbuf;
740   int real_sendcount         = sendcount;
741   MPI_Datatype real_sendtype = sendtype;
742   std::unique_ptr<unsigned char[]> tmp_sendbuf;
743   if (sendbuf == MPI_IN_PLACE) {
744     tmp_sendbuf.reset(new unsigned char[recvcount * comm->size() * recvtype->size()]);
745     // memcpy(??,nullptr,0) is actually undefined behavor, even if harmless.
746     if (recvbuf != nullptr)
747       memcpy(tmp_sendbuf.get(), recvbuf, recvcount * comm->size() * recvtype->size());
748     real_sendbuf = tmp_sendbuf.get();
749     real_sendcount = recvcount;
750     real_sendtype  = recvtype;
751   }
752
753   TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Alltoall" : "PMPI_Ialltoall",
754                      new simgrid::instr::CollTIData(
755                          request == MPI_REQUEST_IGNORED ? "alltoall" : "ialltoall", -1, -1.0,
756                          real_sendtype->is_replayable() ? real_sendcount : real_sendcount * real_sendtype->size(),
757                          recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(),
758                          simgrid::smpi::Datatype::encode(real_sendtype), simgrid::smpi::Datatype::encode(recvtype)));
759   int retval;
760   if (request == MPI_REQUEST_IGNORED)
761     retval =
762         simgrid::smpi::Colls::alltoall(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcount, recvtype, comm);
763   else
764     retval = simgrid::smpi::Colls::ialltoall(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcount, recvtype,
765                                              comm, request);
766
767   TRACE_smpi_comm_out(rank);
768   smpi_bench_begin();
769   return retval;
770 }
771
772 int PMPI_Alltoallv(const void* sendbuf, const int* sendcounts, const int* senddisps, MPI_Datatype sendtype, void* recvbuf,
773                    const int* recvcounts, const int* recvdisps, MPI_Datatype recvtype, MPI_Comm comm)
774 {
775   return PMPI_Ialltoallv(sendbuf, sendcounts, senddisps, sendtype, recvbuf, recvcounts, recvdisps, recvtype, comm, MPI_REQUEST_IGNORED);
776 }
777
778 int PMPI_Ialltoallv(const void* sendbuf, const int* sendcounts, const int* senddisps, MPI_Datatype sendtype, void* recvbuf,
779                     const int* recvcounts, const int* recvdisps, MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request)
780 {
781   if (comm == MPI_COMM_NULL)
782     return MPI_ERR_COMM;
783   if (sendbuf == nullptr || recvbuf == nullptr)
784     return MPI_ERR_BUFFER;
785   if ((sendbuf != MPI_IN_PLACE && sendtype == MPI_DATATYPE_NULL) || recvtype == MPI_DATATYPE_NULL)
786     return MPI_ERR_TYPE;
787   if ((sendbuf != MPI_IN_PLACE && (sendcounts == nullptr || senddisps == nullptr)) || recvcounts == nullptr ||
788       recvdisps == nullptr)
789     return MPI_ERR_ARG;
790   if (request == nullptr)
791     return MPI_ERR_ARG;
792
793   int rank = simgrid::s4u::this_actor::get_pid();
794   int size = comm->size();
795   for (int i = 0; i < size; i++) {
796     if (recvcounts[i] < 0 || (sendbuf != MPI_IN_PLACE && sendcounts[i] < 0))
797       return MPI_ERR_COUNT;
798   }
799
800   smpi_bench_end();
801   int send_size                      = 0;
802   int recv_size                      = 0;
803   std::vector<int>* trace_sendcounts = new std::vector<int>;
804   std::vector<int>* trace_recvcounts = new std::vector<int>;
805   int dt_size_recv                   = recvtype->size();
806
807   const void* real_sendbuf   = sendbuf;
808   const int* real_sendcounts = sendcounts;
809   const int* real_senddisps  = senddisps;
810   MPI_Datatype real_sendtype = sendtype;
811   int maxsize              = 0;
812   for (int i = 0; i < size; i++) { // copy data to avoid bad free
813     recv_size += recvcounts[i] * dt_size_recv;
814     trace_recvcounts->push_back(recvcounts[i] * dt_size_recv);
815     if (((recvdisps[i] + recvcounts[i]) * dt_size_recv) > maxsize)
816       maxsize = (recvdisps[i] + recvcounts[i]) * dt_size_recv;
817   }
818
819   std::unique_ptr<unsigned char[]> tmp_sendbuf;
820   std::unique_ptr<int[]> tmp_sendcounts;
821   std::unique_ptr<int[]> tmp_senddisps;
822   if (sendbuf == MPI_IN_PLACE) {
823     tmp_sendbuf.reset(new unsigned char[maxsize]);
824     real_sendbuf = memcpy(tmp_sendbuf.get(), recvbuf, maxsize);
825     tmp_sendcounts.reset(new int[size]);
826     std::copy(recvcounts, recvcounts + size, tmp_sendcounts.get());
827     real_sendcounts = tmp_sendcounts.get();
828     tmp_senddisps.reset(new int[size]);
829     std::copy(recvdisps, recvdisps + size, tmp_senddisps.get());
830     real_senddisps = tmp_senddisps.get();
831     real_sendtype  = recvtype;
832   }
833
834   int dt_size_send = real_sendtype->size();
835
836   for (int i = 0; i < size; i++) { // copy data to avoid bad free
837     send_size += real_sendcounts[i] * dt_size_send;
838     trace_sendcounts->push_back(real_sendcounts[i] * dt_size_send);
839   }
840
841   TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Alltoallv" : "PMPI_Ialltoallv",
842                      new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "alltoallv" : "ialltoallv", -1,
843                                                        send_size, trace_sendcounts, recv_size, trace_recvcounts,
844                                                        simgrid::smpi::Datatype::encode(real_sendtype),
845                                                        simgrid::smpi::Datatype::encode(recvtype)));
846
847   int retval;
848   if (request == MPI_REQUEST_IGNORED)
849     retval = simgrid::smpi::Colls::alltoallv(real_sendbuf, real_sendcounts, real_senddisps, real_sendtype, recvbuf,
850                                              recvcounts, recvdisps, recvtype, comm);
851   else
852     retval = simgrid::smpi::Colls::ialltoallv(real_sendbuf, real_sendcounts, real_senddisps, real_sendtype, recvbuf,
853                                               recvcounts, recvdisps, recvtype, comm, request);
854
855   TRACE_smpi_comm_out(rank);
856   smpi_bench_begin();
857   return retval;
858 }
859
860 int PMPI_Alltoallw(const void* sendbuf, const int* sendcounts, const int* senddisps, const MPI_Datatype* sendtypes, void* recvbuf,
861                    const int* recvcounts, const int* recvdisps, const MPI_Datatype* recvtypes, MPI_Comm comm)
862 {
863   return PMPI_Ialltoallw(sendbuf, sendcounts, senddisps, sendtypes, recvbuf, recvcounts, recvdisps, recvtypes, comm, MPI_REQUEST_IGNORED);
864 }
865
866 int PMPI_Ialltoallw(const void* sendbuf, const int* sendcounts, const int* senddisps, const MPI_Datatype* sendtypes, void* recvbuf,
867                     const int* recvcounts, const int* recvdisps, const MPI_Datatype* recvtypes, MPI_Comm comm, MPI_Request* request)
868 {
869   if (comm == MPI_COMM_NULL)
870     return MPI_ERR_COMM;
871   if (sendbuf == nullptr || recvbuf == nullptr)
872     return MPI_ERR_BUFFER;
873   if ((sendbuf != MPI_IN_PLACE && sendtypes == nullptr) || recvtypes == nullptr)
874     return MPI_ERR_TYPE;
875   if ((sendbuf != MPI_IN_PLACE && (sendcounts == nullptr || senddisps == nullptr)) || recvcounts == nullptr ||
876       recvdisps == nullptr)
877     return MPI_ERR_ARG;
878   if (request == nullptr)
879     return MPI_ERR_ARG;
880
881   smpi_bench_end();
882   int rank = simgrid::s4u::this_actor::get_pid();
883   int size = comm->size();
884   for (int i = 0; i < size; i++) {
885     if (recvcounts[i] < 0 || (sendbuf != MPI_IN_PLACE && sendcounts[i] < 0))
886       return MPI_ERR_COUNT;
887   }
888   int send_size                      = 0;
889   int recv_size                      = 0;
890   std::vector<int>* trace_sendcounts = new std::vector<int>;
891   std::vector<int>* trace_recvcounts = new std::vector<int>;
892
893   const void* real_sendbuf           = sendbuf;
894   const int* real_sendcounts         = sendcounts;
895   const int* real_senddisps          = senddisps;
896   const MPI_Datatype* real_sendtypes = sendtypes;
897   unsigned long maxsize      = 0;
898   for (int i = 0; i < size; i++) { // copy data to avoid bad free
899     if (recvtypes[i] == MPI_DATATYPE_NULL) {
900       delete trace_recvcounts;
901       delete trace_sendcounts;
902       return MPI_ERR_TYPE;
903     }
904     recv_size += recvcounts[i] * recvtypes[i]->size();
905     trace_recvcounts->push_back(recvcounts[i] * recvtypes[i]->size());
906     if ((recvdisps[i] + (recvcounts[i] * recvtypes[i]->size())) > maxsize)
907       maxsize = recvdisps[i] + (recvcounts[i] * recvtypes[i]->size());
908   }
909
910   std::unique_ptr<unsigned char[]> tmp_sendbuf;
911   std::unique_ptr<int[]> tmp_sendcounts;
912   std::unique_ptr<int[]> tmp_senddisps;
913   std::unique_ptr<MPI_Datatype[]> tmp_sendtypes;
914   if (sendbuf == MPI_IN_PLACE) {
915     tmp_sendbuf.reset(new unsigned char[maxsize]);
916     real_sendbuf = memcpy(tmp_sendbuf.get(), recvbuf, maxsize);
917     tmp_sendcounts.reset(new int[size]);
918     std::copy(recvcounts, recvcounts + size, tmp_sendcounts.get());
919     real_sendcounts = tmp_sendcounts.get();
920     tmp_senddisps.reset(new int[size]);
921     std::copy(recvdisps, recvdisps + size, tmp_senddisps.get());
922     real_senddisps = tmp_senddisps.get();
923     tmp_sendtypes.reset(new MPI_Datatype[size]);
924     std::copy(recvtypes, recvtypes + size, tmp_sendtypes.get());
925     real_sendtypes = tmp_sendtypes.get();
926   }
927
928   for (int i = 0; i < size; i++) { // copy data to avoid bad free
929     send_size += real_sendcounts[i] * real_sendtypes[i]->size();
930     trace_sendcounts->push_back(real_sendcounts[i] * real_sendtypes[i]->size());
931   }
932
933   TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Alltoallw" : "PMPI_Ialltoallw",
934                      new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "alltoallv" : "ialltoallv", -1,
935                                                        send_size, trace_sendcounts, recv_size, trace_recvcounts,
936                                                        simgrid::smpi::Datatype::encode(real_sendtypes[0]),
937                                                        simgrid::smpi::Datatype::encode(recvtypes[0])));
938
939   int retval;
940   if (request == MPI_REQUEST_IGNORED)
941     retval = simgrid::smpi::Colls::alltoallw(real_sendbuf, real_sendcounts, real_senddisps, real_sendtypes, recvbuf,
942                                              recvcounts, recvdisps, recvtypes, comm);
943   else
944     retval = simgrid::smpi::Colls::ialltoallw(real_sendbuf, real_sendcounts, real_senddisps, real_sendtypes, recvbuf,
945                                               recvcounts, recvdisps, recvtypes, comm, request);
946
947   TRACE_smpi_comm_out(rank);
948   smpi_bench_begin();
949   return retval;
950 }