Logo AND Algorithmique Numérique Distribuée

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