Logo AND Algorithmique Numérique Distribuée

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