Logo AND Algorithmique Numérique Distribuée

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