Logo AND Algorithmique Numérique Distribuée

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