Logo AND Algorithmique Numérique Distribuée

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