Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
oops, maxsize was needed here.
[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   int rank = comm->rank();
103   if(sendbuf != MPI_IN_PLACE){
104     CHECK_COUNT(2, sendcount)
105     CHECK_TYPE(3, sendtype)
106     CHECK_BUFFER(1,sendbuf, sendcount, sendtype)
107   }
108   if(rank == root){
109     SET_BUF2(recvbuf)
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   int rank = comm->rank();
163   if(sendbuf != MPI_IN_PLACE){
164     CHECK_TYPE(3, sendtype)
165     CHECK_COUNT(2, sendcount)
166   }
167   CHECK_BUFFER(1, sendbuf, sendcount, sendtype)
168   if(rank == root){
169     SET_BUF2(recvbuf)
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     trace_recvcounts->insert(trace_recvcounts->end(), &recvcounts[0], &recvcounts[comm->size()]);
201   else //this is not significant outside of root, put 0 as we don't know if recvcounts is initialized
202     trace_recvcounts->insert(trace_recvcounts->end(), comm->size(), 0);
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                          sendcount,
208                          nullptr, -1, 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   trace_recvcounts->insert(trace_recvcounts->end(), &recvcounts[0], &recvcounts[comm->size()]);
311
312   TRACE_smpi_comm_in(
313       pid, request == MPI_REQUEST_IGNORED ? "PMPI_Allgatherv" : "PMPI_Iallgatherv",
314       new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "allgatherv" : "iallgatherv", -1,
315                                         sendcount, nullptr,
316                                         -1, trace_recvcounts, simgrid::smpi::Datatype::encode(sendtype),
317                                         simgrid::smpi::Datatype::encode(recvtype)));
318   if (request == MPI_REQUEST_IGNORED)
319     simgrid::smpi::colls::allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm);
320   else
321     simgrid::smpi::colls::iallgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm,
322                                       request);
323
324   TRACE_smpi_comm_out(pid);
325   return MPI_SUCCESS;
326 }
327
328 int PMPI_Scatter(const void *sendbuf, int sendcount, MPI_Datatype sendtype,
329                 void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm){
330   return PMPI_Iscatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm, MPI_REQUEST_IGNORED);
331 }
332
333 int PMPI_Iscatter(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
334                   MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request* request)
335 {
336   CHECK_COMM(8)
337   SET_BUF2(recvbuf)
338   int rank = comm->rank();
339   if(rank == root){
340     SET_BUF1(sendbuf)
341     CHECK_NOT_IN_PLACE_ROOT(1, sendbuf)
342     CHECK_COUNT(2, sendcount)
343     CHECK_TYPE(3, sendtype)
344     CHECK_BUFFER(1, sendbuf, sendcount, sendtype)
345   } else {
346     CHECK_NOT_IN_PLACE_ROOT(4, recvbuf)
347   }
348   if(recvbuf != MPI_IN_PLACE){
349     CHECK_COUNT(5, recvcount)
350     CHECK_TYPE(6, recvtype)
351     CHECK_BUFFER(4, recvbuf, recvcount, recvtype)
352   }
353   CHECK_ROOT(8)
354   CHECK_REQUEST(9)
355
356   if (recvbuf == MPI_IN_PLACE) {
357     recvtype  = sendtype;
358     recvcount = sendcount;
359   }
360
361   if((rank == root) && (recvtype->size() * recvcount !=  sendtype->size() * sendcount)){
362     XBT_WARN("MPI_(I)Scatter : sent size to each process differs from receive size");
363     return MPI_ERR_TRUNCATE;
364   }
365
366   const SmpiBenchGuard suspend_bench;
367
368   aid_t pid = simgrid::s4u::this_actor::get_pid();
369
370   TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Scatter" : "PMPI_Iscatter",
371                      new simgrid::instr::CollTIData(
372                          request == MPI_REQUEST_IGNORED ? "scatter" : "iscatter", root, -1.0,
373                          sendcount, recvcount,
374                          simgrid::smpi::Datatype::encode(sendtype), simgrid::smpi::Datatype::encode(recvtype)));
375   if (request == MPI_REQUEST_IGNORED)
376     simgrid::smpi::colls::scatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm);
377   else
378     simgrid::smpi::colls::iscatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm, request);
379
380   TRACE_smpi_comm_out(pid);
381   return MPI_SUCCESS;
382 }
383
384 int PMPI_Scatterv(const void *sendbuf, const int *sendcounts, const int *displs,
385                  MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm){
386   return PMPI_Iscatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm, MPI_REQUEST_IGNORED);
387 }
388
389 int PMPI_Iscatterv(const void* sendbuf, const int* sendcounts, const int* displs, MPI_Datatype sendtype, void* recvbuf, int recvcount,
390                    MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request* request)
391 {
392   SET_BUF2(recvbuf)
393   CHECK_COMM(9)
394   int rank = comm->rank();
395   if(recvbuf != MPI_IN_PLACE){
396     CHECK_COUNT(5, recvcount)
397     CHECK_TYPE(7, recvtype)
398     CHECK_BUFFER(4, recvbuf, recvcount, recvtype)
399   }
400   CHECK_ROOT(9)
401   CHECK_REQUEST(10)
402   if (rank == root) {
403     SET_BUF1(sendbuf)
404     CHECK_NOT_IN_PLACE_ROOT(1, sendbuf)
405     CHECK_NULL(2, MPI_ERR_COUNT, sendcounts)
406     CHECK_NULL(3, MPI_ERR_ARG, displs)
407     CHECK_TYPE(4, sendtype)
408     for (int i = 0; i < comm->size(); i++){
409       CHECK_COUNT(2, sendcounts[i])
410       CHECK_BUFFER(1, sendbuf, sendcounts[i], sendtype)
411     }
412     if (recvbuf == MPI_IN_PLACE) {
413       recvtype  = sendtype;
414       recvcount = sendcounts[rank];
415     }
416   } else {
417     CHECK_NOT_IN_PLACE_ROOT(4, recvbuf)
418   }
419
420   const SmpiBenchGuard suspend_bench;
421
422   aid_t pid        = simgrid::s4u::this_actor::get_pid();
423
424   auto trace_sendcounts = std::make_shared<std::vector<int>>();
425   if (rank == root)
426     trace_sendcounts->insert(trace_sendcounts->end(), &sendcounts[0], &sendcounts[comm->size()]);
427   else //this is not significant outside of root, put 0 as we don't know if sendcounts is initialized
428     trace_sendcounts->insert(trace_sendcounts->end(), comm->size(), 0);
429
430
431   TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Scatterv" : "PMPI_Iscatterv",
432                      new simgrid::instr::VarCollTIData(
433                          request == MPI_REQUEST_IGNORED ? "scatterv" : "iscatterv", root, -1,
434                          trace_sendcounts, recvcount,
435                          nullptr, simgrid::smpi::Datatype::encode(sendtype),
436                          simgrid::smpi::Datatype::encode(recvtype)));
437   if (request == MPI_REQUEST_IGNORED)
438     simgrid::smpi::colls::scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm);
439   else
440     simgrid::smpi::colls::iscatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm,
441                                     request);
442
443   TRACE_smpi_comm_out(pid);
444   return MPI_SUCCESS;
445 }
446
447 int PMPI_Reduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
448 {
449   return PMPI_Ireduce(sendbuf, recvbuf, count, datatype, op, root, comm, MPI_REQUEST_IGNORED);
450 }
451
452 int PMPI_Ireduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm, MPI_Request* request)
453 {
454   CHECK_COMM(7)
455   SET_BUF1(sendbuf)
456   int rank = comm->rank();
457   CHECK_TYPE(4, datatype)
458   CHECK_COUNT(3, count)
459   CHECK_BUFFER(1, sendbuf, count, datatype)
460   if(rank == root){
461     SET_BUF2(recvbuf)
462     CHECK_NOT_IN_PLACE(2, recvbuf)
463     CHECK_BUFFER(5, recvbuf, count, datatype)
464   }
465   CHECK_OP(5, op, datatype)
466   CHECK_ROOT(7)
467   CHECK_REQUEST(8)
468
469   const SmpiBenchGuard suspend_bench;
470   aid_t pid = simgrid::s4u::this_actor::get_pid();
471
472   TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Reduce" : "PMPI_Ireduce",
473                      new simgrid::instr::CollTIData(request == MPI_REQUEST_IGNORED ? "reduce" : "ireduce", root, 0,
474                                                     count, 0,
475                                                     simgrid::smpi::Datatype::encode(datatype), ""));
476   if (request == MPI_REQUEST_IGNORED)
477     simgrid::smpi::colls::reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
478   else
479     simgrid::smpi::colls::ireduce(sendbuf, recvbuf, count, datatype, op, root, comm, request);
480
481   TRACE_smpi_comm_out(pid);
482   return MPI_SUCCESS;
483 }
484
485 int PMPI_Reduce_local(const void* inbuf, void* inoutbuf, int count, MPI_Datatype datatype, MPI_Op op)
486 {
487   SET_BUF1(inbuf)
488   SET_BUF2(inoutbuf)
489   CHECK_TYPE(4, datatype)
490   CHECK_COUNT(3, count)
491   CHECK_BUFFER(1, inbuf, count, datatype)
492   CHECK_BUFFER(2, inoutbuf, count, datatype)
493   CHECK_OP(5, op, datatype)
494
495   const SmpiBenchGuard suspend_bench;
496   op->apply(inbuf, inoutbuf, &count, datatype);
497   return MPI_SUCCESS;
498 }
499
500 int PMPI_Allreduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
501 {
502   return PMPI_Iallreduce(sendbuf, recvbuf, count, datatype, op, comm, MPI_REQUEST_IGNORED);
503 }
504
505 int PMPI_Iallreduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request *request)
506 {
507   CHECK_COMM(6)
508   SET_BUF1(sendbuf)
509   SET_BUF2(recvbuf)
510   int rank = comm->rank();
511   CHECK_NOT_IN_PLACE(2, recvbuf)
512   CHECK_TYPE(4, datatype)
513   CHECK_OP(5, op, datatype)
514   CHECK_COUNT(3, count)
515   CHECK_BUFFER(1, sendbuf, count, datatype)
516   CHECK_BUFFER(2, recvbuf, count, datatype)
517   CHECK_REQUEST(7)
518
519   const SmpiBenchGuard suspend_bench;
520   std::vector<unsigned char> tmp_sendbuf;
521   const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, count, datatype);
522
523   aid_t pid = simgrid::s4u::this_actor::get_pid();
524
525   TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Allreduce" : "PMPI_Iallreduce",
526                      new simgrid::instr::CollTIData(request == MPI_REQUEST_IGNORED ? "allreduce" : "iallreduce", -1, 0,
527                                                     count, 0,
528                                                     simgrid::smpi::Datatype::encode(datatype), ""));
529
530   if (request == MPI_REQUEST_IGNORED)
531     simgrid::smpi::colls::allreduce(real_sendbuf, recvbuf, count, datatype, op, comm);
532   else
533     simgrid::smpi::colls::iallreduce(real_sendbuf, recvbuf, count, datatype, op, comm, request);
534
535   TRACE_smpi_comm_out(pid);
536   return MPI_SUCCESS;
537 }
538
539 int PMPI_Scan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
540 {
541   return PMPI_Iscan(sendbuf, recvbuf, count, datatype, op, comm, MPI_REQUEST_IGNORED);
542 }
543
544 int PMPI_Iscan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request* request)
545 {
546   CHECK_COMM(6)
547   SET_BUF1(sendbuf)
548   SET_BUF2(recvbuf)
549   CHECK_TYPE(4, datatype)
550   CHECK_COUNT(3, count)
551   CHECK_BUFFER(1,sendbuf,count, datatype)
552   CHECK_BUFFER(2,recvbuf,count, datatype)
553   CHECK_REQUEST(7)
554   CHECK_OP(5, op, datatype)
555
556   const SmpiBenchGuard suspend_bench;
557   aid_t pid = simgrid::s4u::this_actor::get_pid();
558   std::vector<unsigned char> tmp_sendbuf;
559   const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, count, datatype);
560
561   TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Scan" : "PMPI_Iscan",
562                      new simgrid::instr::CollTIData(request == MPI_REQUEST_IGNORED ? "scan" : "iscan", -1, 0.0,
563                                                     count, 0, simgrid::smpi::Datatype::encode(datatype), ""));
564
565   int retval;
566   if (request == MPI_REQUEST_IGNORED)
567     retval = simgrid::smpi::colls::scan(real_sendbuf, recvbuf, count, datatype, op, comm);
568   else
569     retval = simgrid::smpi::colls::iscan(real_sendbuf, recvbuf, count, datatype, op, comm, request);
570
571   TRACE_smpi_comm_out(pid);
572   return retval;
573 }
574
575 int PMPI_Exscan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
576 {
577   return PMPI_Iexscan(sendbuf, recvbuf, count, datatype, op, comm, MPI_REQUEST_IGNORED);
578 }
579
580 int PMPI_Iexscan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request* request){
581   CHECK_COMM(6)
582   SET_BUF1(sendbuf)
583   SET_BUF2(recvbuf)
584   CHECK_TYPE(4, datatype)
585   CHECK_COUNT(3, count)
586   CHECK_BUFFER(1, sendbuf, count, datatype)
587   CHECK_BUFFER(2, recvbuf, count, datatype)
588   CHECK_REQUEST(7)
589   CHECK_OP(5, op, datatype)
590
591   const SmpiBenchGuard suspend_bench;
592   aid_t pid = simgrid::s4u::this_actor::get_pid();
593   std::vector<unsigned char> tmp_sendbuf;
594   const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, count, datatype);
595
596   TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Exscan" : "PMPI_Iexscan",
597                      new simgrid::instr::CollTIData(request == MPI_REQUEST_IGNORED ? "exscan" : "iexscan", -1, 0.0,
598                                                     count, 0, simgrid::smpi::Datatype::encode(datatype), ""));
599
600   int retval;
601   if (request == MPI_REQUEST_IGNORED)
602     retval = simgrid::smpi::colls::exscan(real_sendbuf, recvbuf, count, datatype, op, comm);
603   else
604     retval = simgrid::smpi::colls::iexscan(real_sendbuf, recvbuf, count, datatype, op, comm, request);
605
606   TRACE_smpi_comm_out(pid);
607   return retval;
608 }
609
610 int PMPI_Reduce_scatter(const void *sendbuf, void *recvbuf, const int *recvcounts, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
611 {
612   return PMPI_Ireduce_scatter(sendbuf, recvbuf, recvcounts, datatype, op, comm, MPI_REQUEST_IGNORED);
613 }
614
615 int PMPI_Ireduce_scatter(const void *sendbuf, void *recvbuf, const int *recvcounts, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request *request)
616 {
617   CHECK_COMM(6)
618   SET_BUF1(sendbuf)
619   SET_BUF2(recvbuf)
620   int rank = comm->rank();
621   CHECK_NOT_IN_PLACE(2, recvbuf)
622   CHECK_TYPE(4, datatype)
623   CHECK_NULL(3, MPI_ERR_COUNT, recvcounts)
624   CHECK_REQUEST(7)
625   CHECK_OP(5, op, datatype)
626   for (int i = 0; i < comm->size(); i++) {
627     CHECK_COUNT(3, recvcounts[i])
628     CHECK_BUFFER(1, sendbuf, recvcounts[i], datatype)
629     CHECK_BUFFER(2, recvbuf, recvcounts[i], datatype)
630   }
631
632   const SmpiBenchGuard suspend_bench;
633   aid_t pid                          = simgrid::s4u::this_actor::get_pid();
634   auto trace_recvcounts              = std::make_shared<std::vector<int>>();
635   trace_recvcounts->insert(trace_recvcounts->end(), &recvcounts[0], &recvcounts[comm->size()]);
636
637   int totalcount                     = 0;
638
639   for (int i = 0; i < comm->size(); i++) { // copy data to avoid bad free
640     totalcount += recvcounts[i];
641   }
642   std::vector<unsigned char> tmp_sendbuf;
643   const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, totalcount, datatype);
644
645   TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Reduce_scatter" : "PMPI_Ireduce_scatter",
646                      new simgrid::instr::VarCollTIData(
647                          request == MPI_REQUEST_IGNORED ? "reducescatter" : "ireducescatter", -1, -1, nullptr,
648                          0, trace_recvcounts, simgrid::smpi::Datatype::encode(datatype), ""));
649
650   if (request == MPI_REQUEST_IGNORED)
651     simgrid::smpi::colls::reduce_scatter(real_sendbuf, recvbuf, recvcounts, datatype, op, comm);
652   else
653     simgrid::smpi::colls::ireduce_scatter(real_sendbuf, recvbuf, recvcounts, datatype, op, comm, request);
654
655   TRACE_smpi_comm_out(pid);
656   return MPI_SUCCESS;
657 }
658
659 int PMPI_Reduce_scatter_block(const void *sendbuf, void *recvbuf, int recvcount,
660                               MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
661 {
662   return PMPI_Ireduce_scatter_block(sendbuf, recvbuf, recvcount, datatype, op, comm, MPI_REQUEST_IGNORED);
663 }
664
665 int PMPI_Ireduce_scatter_block(const void* sendbuf, void* recvbuf, int recvcount, MPI_Datatype datatype, MPI_Op op,
666                                MPI_Comm comm, MPI_Request* request)
667 {
668   CHECK_COMM(6)
669   SET_BUF1(sendbuf)
670   SET_BUF2(recvbuf)
671   CHECK_TYPE(4, datatype)
672   CHECK_COUNT(3, recvcount)
673   CHECK_BUFFER(1, sendbuf, recvcount, datatype)
674   CHECK_BUFFER(2, recvbuf, recvcount, datatype)
675   CHECK_REQUEST(7)
676   CHECK_OP(5, op, datatype)
677
678   const SmpiBenchGuard suspend_bench;
679   int count = comm->size();
680
681   aid_t pid                          = simgrid::s4u::this_actor::get_pid();
682   auto trace_recvcounts = std::make_shared<std::vector<int>>(recvcount);
683
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, -1,
690                                         nullptr, -1, 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   trace_recvcounts->insert(trace_recvcounts->end(), &recvcounts[0], &recvcounts[size]);
801
802   int dt_size_recv                   = recvtype->size();
803
804   const int* real_sendcounts = sendcounts;
805   const int* real_senddispls  = senddispls;
806   MPI_Datatype real_sendtype = sendtype;
807   int maxsize              = 0;
808   std::vector<unsigned char> tmp_sendbuf;
809   std::vector<int> tmp_sendcounts;
810   std::vector<int> tmp_senddispls;
811   const void* real_sendbuf;
812
813   if (sendbuf == MPI_IN_PLACE) {
814     tmp_sendcounts.assign(recvcounts, recvcounts + size);
815     real_sendcounts = tmp_sendcounts.data();
816     tmp_senddispls.assign(recvdispls, recvdispls + size);
817     real_senddispls = tmp_senddispls.data();
818     real_sendtype  = recvtype;
819   }
820
821   for (int i = 0; i < size; i++) { // copy data to avoid bad free
822     send_size += real_sendcounts[i] ;
823     recv_size += recvcounts[i];
824     if (((recvdispls[i] + recvcounts[i]) * dt_size_recv) > maxsize)
825       maxsize = (recvdispls[i] + recvcounts[i]) * dt_size_recv;
826   }
827   real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, maxsize, MPI_CHAR);
828
829   if(recvtype->size() * recvcounts[comm->rank()] !=  real_sendtype->size() * real_sendcounts[comm->rank()]){
830     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()]);
831     return MPI_ERR_TRUNCATE;
832   }
833
834   trace_sendcounts->insert(trace_sendcounts->end(), &real_sendcounts[0], &real_sendcounts[size]);
835
836   TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Alltoallv" : "PMPI_Ialltoallv",
837                      new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "alltoallv" : "ialltoallv", -1,
838                                                        send_size, trace_sendcounts, recv_size, trace_recvcounts,
839                                                        simgrid::smpi::Datatype::encode(real_sendtype),
840                                                        simgrid::smpi::Datatype::encode(recvtype)));
841
842   int retval;
843   if (request == MPI_REQUEST_IGNORED)
844     retval = simgrid::smpi::colls::alltoallv(real_sendbuf, real_sendcounts, real_senddispls, real_sendtype, recvbuf,
845                                              recvcounts, recvdispls, recvtype, comm);
846   else
847     retval = simgrid::smpi::colls::ialltoallv(real_sendbuf, real_sendcounts, real_senddispls, real_sendtype, recvbuf,
848                                               recvcounts, recvdispls, recvtype, comm, request);
849
850   TRACE_smpi_comm_out(pid);
851   return retval;
852 }
853
854 int PMPI_Alltoallw(const void* sendbuf, const int* sendcounts, const int* senddispls, const MPI_Datatype* sendtypes, void* recvbuf,
855                    const int* recvcounts, const int* recvdispls, const MPI_Datatype* recvtypes, MPI_Comm comm)
856 {
857   return PMPI_Ialltoallw(sendbuf, sendcounts, senddispls, sendtypes, recvbuf, recvcounts, recvdispls, recvtypes, comm, MPI_REQUEST_IGNORED);
858 }
859
860 int PMPI_Ialltoallw(const void* sendbuf, const int* sendcounts, const int* senddispls, const MPI_Datatype* sendtypes, void* recvbuf,
861                     const int* recvcounts, const int* recvdispls, const MPI_Datatype* recvtypes, MPI_Comm comm, MPI_Request* request)
862 {
863   CHECK_COMM(9)
864   SET_BUF1(sendbuf)
865   SET_BUF2(recvbuf)
866   if(sendbuf != MPI_IN_PLACE){
867     CHECK_NULL(2, MPI_ERR_COUNT, sendcounts)
868     CHECK_NULL(3, MPI_ERR_ARG, senddispls)
869     CHECK_NULL(4, MPI_ERR_TYPE, sendtypes)
870   }
871   CHECK_NULL(6, MPI_ERR_COUNT, recvcounts)
872   CHECK_NULL(7, MPI_ERR_ARG, recvdispls)
873   CHECK_NULL(8, MPI_ERR_TYPE, recvtypes)
874   CHECK_REQUEST(10)
875   aid_t pid = simgrid::s4u::this_actor::get_pid();
876   int size = comm->size();
877   for (int i = 0; i < size; i++) {
878     if(sendbuf != MPI_IN_PLACE){
879       CHECK_COUNT(2, sendcounts[i])
880       CHECK_TYPE(4, sendtypes[i])
881       CHECK_BUFFER(1, sendbuf, sendcounts[i], sendtypes[i])
882     }
883     CHECK_COUNT(6, recvcounts[i])
884     CHECK_TYPE(8, recvtypes[i])
885     CHECK_BUFFER(5, recvbuf, recvcounts[i], recvtypes[i])
886   }
887
888   const SmpiBenchGuard suspend_bench;
889
890   int send_size                      = 0;
891   int recv_size                      = 0;
892   auto trace_sendcounts              = std::make_shared<std::vector<int>>();
893   auto trace_recvcounts              = std::make_shared<std::vector<int>>();
894   trace_recvcounts->insert(trace_recvcounts->end(), &recvcounts[0], &recvcounts[size]);
895
896   const int* real_sendcounts         = sendcounts;
897   const int* real_senddispls          = senddispls;
898   const MPI_Datatype* real_sendtypes = sendtypes;
899
900   unsigned long maxsize      = 0;
901   for (int i = 0; i < size; i++) { // copy data to avoid bad free
902     if (recvtypes[i] == MPI_DATATYPE_NULL)
903       return MPI_ERR_TYPE;
904     recv_size += recvcounts[i] * recvtypes[i]->size();
905     if ((recvdispls[i] + (recvcounts[i] * recvtypes[i]->size())) > maxsize)
906       maxsize = recvdispls[i] + (recvcounts[i] * recvtypes[i]->size());
907   }
908
909   std::vector<unsigned char> tmp_sendbuf;
910   std::vector<int> tmp_sendcounts;
911   std::vector<int> tmp_senddispls;
912   std::vector<MPI_Datatype> tmp_sendtypes;
913   const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, maxsize, MPI_CHAR);
914   if (sendbuf == MPI_IN_PLACE) {
915     tmp_sendcounts.assign(recvcounts, recvcounts + size);
916     real_sendcounts = tmp_sendcounts.data();
917     tmp_senddispls.assign(recvdispls, recvdispls + size);
918     real_senddispls = tmp_senddispls.data();
919     tmp_sendtypes.assign(recvtypes, recvtypes + size);
920     real_sendtypes = tmp_sendtypes.data();
921   }
922
923
924   if(recvtypes[comm->rank()]->size() * recvcounts[comm->rank()] !=  real_sendtypes[comm->rank()]->size() * real_sendcounts[comm->rank()]){
925     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()]);
926     return MPI_ERR_TRUNCATE;
927   }
928
929   trace_sendcounts->insert(trace_sendcounts->end(), &real_sendcounts[0], &real_sendcounts[size]);
930   for (int i = 0; i < size; i++) { // copy data to avoid bad free
931     send_size += real_sendcounts[i] * real_sendtypes[i]->size();
932   }
933
934   TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Alltoallw" : "PMPI_Ialltoallw",
935                      new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "alltoallv" : "ialltoallv", -1,
936                                                        send_size, trace_sendcounts, recv_size, trace_recvcounts,
937                                                        simgrid::smpi::Datatype::encode(real_sendtypes[0]),
938                                                        simgrid::smpi::Datatype::encode(recvtypes[0])));
939
940   int retval;
941   if (request == MPI_REQUEST_IGNORED)
942     retval = simgrid::smpi::colls::alltoallw(real_sendbuf, real_sendcounts, real_senddispls, real_sendtypes, recvbuf,
943                                              recvcounts, recvdispls, recvtypes, comm);
944   else
945     retval = simgrid::smpi::colls::ialltoallw(real_sendbuf, real_sendcounts, real_senddispls, real_sendtypes, recvbuf,
946                                               recvcounts, recvdispls, recvtypes, comm, request);
947
948   TRACE_smpi_comm_out(pid);
949   return retval;
950 }