Logo AND Algorithmique Numérique Distribuée

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