Logo AND Algorithmique Numérique Distribuée

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