Logo AND Algorithmique Numérique Distribuée

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