Logo AND Algorithmique Numérique Distribuée

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