Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
[TRACING] allGatherV -> allgatherv
[simgrid.git] / src / smpi / bindings / smpi_pmpi_coll.cpp
1 /* Copyright (c) 2007-2018. 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_datatype_derived.hpp"
10 #include "smpi_op.hpp"
11 #include "smpi_process.hpp"
12
13 XBT_LOG_EXTERNAL_DEFAULT_CATEGORY(smpi_pmpi);
14
15 /* PMPI User level calls */
16
17 int PMPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
18 {
19   int retval = 0;
20
21   smpi_bench_end();
22
23   if (comm == MPI_COMM_NULL) {
24     retval = MPI_ERR_COMM;
25   } else if (not datatype->is_valid()) {
26     retval = MPI_ERR_ARG;
27   } else {
28     int rank = simgrid::s4u::this_actor::get_pid();
29     TRACE_smpi_comm_in(rank, __func__,
30                        new simgrid::instr::CollTIData("bcast", root, -1.0,
31                                                       datatype->is_replayable() ? count : count * datatype->size(), -1,
32                                                       simgrid::smpi::Datatype::encode(datatype), ""));
33     if (comm->size() > 1)
34       simgrid::smpi::Colls::bcast(buf, count, datatype, root, comm);
35     retval = MPI_SUCCESS;
36
37     TRACE_smpi_comm_out(rank);
38   }
39   smpi_bench_begin();
40   return retval;
41 }
42
43 int PMPI_Barrier(MPI_Comm comm)
44 {
45   int retval = 0;
46
47   smpi_bench_end();
48
49   if (comm == MPI_COMM_NULL) {
50     retval = MPI_ERR_COMM;
51   } else {
52     int rank = simgrid::s4u::this_actor::get_pid();
53     TRACE_smpi_comm_in(rank, __func__, new simgrid::instr::NoOpTIData("barrier"));
54
55     simgrid::smpi::Colls::barrier(comm);
56
57     //Barrier can be used to synchronize RMA calls. Finish all requests from comm before.
58     comm->finish_rma_calls();
59
60     retval = MPI_SUCCESS;
61
62     TRACE_smpi_comm_out(rank);
63   }
64
65   smpi_bench_begin();
66   return retval;
67 }
68
69 int PMPI_Gather(void *sendbuf, int sendcount, MPI_Datatype sendtype,void *recvbuf, int recvcount, MPI_Datatype recvtype,
70                 int root, MPI_Comm comm)
71 {
72   int retval = 0;
73
74   smpi_bench_end();
75
76   if (comm == MPI_COMM_NULL) {
77     retval = MPI_ERR_COMM;
78   } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
79             ((comm->rank() == root) && (recvtype == MPI_DATATYPE_NULL))){
80     retval = MPI_ERR_TYPE;
81   } else if ((( sendbuf != MPI_IN_PLACE) && (sendcount <0)) || ((comm->rank() == root) && (recvcount <0))){
82     retval = MPI_ERR_COUNT;
83   } else {
84
85     char* sendtmpbuf = static_cast<char*>(sendbuf);
86     int sendtmpcount = sendcount;
87     MPI_Datatype sendtmptype = sendtype;
88     if( (comm->rank() == root) && (sendbuf == MPI_IN_PLACE )) {
89       sendtmpcount=0;
90       sendtmptype=recvtype;
91     }
92     int rank = simgrid::s4u::this_actor::get_pid();
93
94     TRACE_smpi_comm_in(
95         rank, __func__,
96         new simgrid::instr::CollTIData(
97             "gather", root, -1.0, sendtmptype->is_replayable() ? sendtmpcount : sendtmpcount * sendtmptype->size(),
98             (comm->rank() != root || recvtype->is_replayable()) ? recvcount : recvcount * recvtype->size(),
99             simgrid::smpi::Datatype::encode(sendtmptype), simgrid::smpi::Datatype::encode(recvtype)));
100
101     simgrid::smpi::Colls::gather(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcount, recvtype, root, comm);
102
103     retval = MPI_SUCCESS;
104     TRACE_smpi_comm_out(rank);
105   }
106
107   smpi_bench_begin();
108   return retval;
109 }
110
111 int PMPI_Gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, int *recvcounts, int *displs,
112                 MPI_Datatype recvtype, int root, MPI_Comm comm)
113 {
114   int retval = 0;
115
116   smpi_bench_end();
117
118   if (comm == MPI_COMM_NULL) {
119     retval = MPI_ERR_COMM;
120   } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
121             ((comm->rank() == root) && (recvtype == MPI_DATATYPE_NULL))){
122     retval = MPI_ERR_TYPE;
123   } else if (( sendbuf != MPI_IN_PLACE) && (sendcount <0)){
124     retval = MPI_ERR_COUNT;
125   } else if ((comm->rank() == root) && (recvcounts == nullptr || displs == nullptr)) {
126     retval = MPI_ERR_ARG;
127   } else {
128     char* sendtmpbuf = static_cast<char*>(sendbuf);
129     int sendtmpcount = sendcount;
130     MPI_Datatype sendtmptype = sendtype;
131     if( (comm->rank() == root) && (sendbuf == MPI_IN_PLACE )) {
132       sendtmpcount=0;
133       sendtmptype=recvtype;
134     }
135
136     int rank         = simgrid::s4u::this_actor::get_pid();
137     int dt_size_recv = recvtype->is_replayable() ? 1 : recvtype->size();
138
139     std::vector<int>* trace_recvcounts = new std::vector<int>;
140     if (comm->rank() == root) {
141       for (int i = 0; i < comm->size(); i++) // copy data to avoid bad free
142         trace_recvcounts->push_back(recvcounts[i] * dt_size_recv);
143     }
144
145     TRACE_smpi_comm_in(rank, __func__,
146                        new simgrid::instr::VarCollTIData(
147                            "gatherV", root,
148                            sendtmptype->is_replayable() ? sendtmpcount : sendtmpcount * sendtmptype->size(), nullptr,
149                            dt_size_recv, trace_recvcounts, simgrid::smpi::Datatype::encode(sendtmptype),
150                            simgrid::smpi::Datatype::encode(recvtype)));
151
152     retval = simgrid::smpi::Colls::gatherv(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcounts, displs, recvtype, root, comm);
153     TRACE_smpi_comm_out(rank);
154   }
155
156   smpi_bench_begin();
157   return retval;
158 }
159
160 int PMPI_Allgather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
161                    void *recvbuf, int recvcount, MPI_Datatype recvtype, MPI_Comm comm)
162 {
163   int retval = MPI_SUCCESS;
164
165   smpi_bench_end();
166
167   if (comm == MPI_COMM_NULL) {
168     retval = MPI_ERR_COMM;
169   } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
170             (recvtype == MPI_DATATYPE_NULL)){
171     retval = MPI_ERR_TYPE;
172   } else if ((( sendbuf != MPI_IN_PLACE) && (sendcount <0)) ||
173             (recvcount <0)){
174     retval = MPI_ERR_COUNT;
175   } else {
176     if(sendbuf == MPI_IN_PLACE) {
177       sendbuf=static_cast<char*>(recvbuf)+recvtype->get_extent()*recvcount*comm->rank();
178       sendcount=recvcount;
179       sendtype=recvtype;
180     }
181     int rank = simgrid::s4u::this_actor::get_pid();
182
183     TRACE_smpi_comm_in(rank, __func__,
184                        new simgrid::instr::CollTIData(
185                            "allGather", -1, -1.0, sendtype->is_replayable() ? sendcount : sendcount * sendtype->size(),
186                            recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(),
187                            simgrid::smpi::Datatype::encode(sendtype), simgrid::smpi::Datatype::encode(recvtype)));
188
189     simgrid::smpi::Colls::allgather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
190     TRACE_smpi_comm_out(rank);
191   }
192   smpi_bench_begin();
193   return retval;
194 }
195
196 int PMPI_Allgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
197                    void *recvbuf, int *recvcounts, int *displs, MPI_Datatype recvtype, MPI_Comm comm)
198 {
199   int retval = 0;
200
201   smpi_bench_end();
202
203   if (comm == MPI_COMM_NULL) {
204     retval = MPI_ERR_COMM;
205   } else if (((sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) || (recvtype == MPI_DATATYPE_NULL)) {
206     retval = MPI_ERR_TYPE;
207   } else if (( sendbuf != MPI_IN_PLACE) && (sendcount <0)){
208     retval = MPI_ERR_COUNT;
209   } else if (recvcounts == nullptr || displs == nullptr) {
210     retval = MPI_ERR_ARG;
211   } else {
212
213     if(sendbuf == MPI_IN_PLACE) {
214       sendbuf=static_cast<char*>(recvbuf)+recvtype->get_extent()*displs[comm->rank()];
215       sendcount=recvcounts[comm->rank()];
216       sendtype=recvtype;
217     }
218     int rank               = simgrid::s4u::this_actor::get_pid();
219     int dt_size_recv       = recvtype->is_replayable() ? 1 : recvtype->size();
220
221     std::vector<int>* trace_recvcounts = new std::vector<int>;
222     for (int i = 0; i < comm->size(); i++) // copy data to avoid bad free
223       trace_recvcounts->push_back(recvcounts[i] * dt_size_recv);
224
225     TRACE_smpi_comm_in(rank, __func__,
226                        new simgrid::instr::VarCollTIData(
227                            "allgatherv", -1, sendtype->is_replayable() ? sendcount : sendcount * sendtype->size(),
228                            nullptr, dt_size_recv, trace_recvcounts, simgrid::smpi::Datatype::encode(sendtype),
229                            simgrid::smpi::Datatype::encode(recvtype)));
230
231     simgrid::smpi::Colls::allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm);
232     retval = MPI_SUCCESS;
233     TRACE_smpi_comm_out(rank);
234   }
235
236   smpi_bench_begin();
237   return retval;
238 }
239
240 int PMPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
241                 void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm)
242 {
243   int retval = 0;
244
245   smpi_bench_end();
246
247   if (comm == MPI_COMM_NULL) {
248     retval = MPI_ERR_COMM;
249   } else if (((comm->rank() == root) && (not sendtype->is_valid())) ||
250              ((recvbuf != MPI_IN_PLACE) && (not recvtype->is_valid()))) {
251     retval = MPI_ERR_TYPE;
252   } else if ((sendbuf == recvbuf) ||
253       ((comm->rank()==root) && sendcount>0 && (sendbuf == nullptr))){
254     retval = MPI_ERR_BUFFER;
255   }else {
256
257     if (recvbuf == MPI_IN_PLACE) {
258       recvtype  = sendtype;
259       recvcount = sendcount;
260     }
261     int rank = simgrid::s4u::this_actor::get_pid();
262
263     TRACE_smpi_comm_in(
264         rank, __func__,
265         new simgrid::instr::CollTIData(
266             "scatter", root, -1.0,
267             (comm->rank() != root || sendtype->is_replayable()) ? sendcount : sendcount * sendtype->size(),
268             recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(),
269             simgrid::smpi::Datatype::encode(sendtype), simgrid::smpi::Datatype::encode(recvtype)));
270
271     simgrid::smpi::Colls::scatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm);
272     retval = MPI_SUCCESS;
273     TRACE_smpi_comm_out(rank);
274   }
275
276   smpi_bench_begin();
277   return retval;
278 }
279
280 int PMPI_Scatterv(void *sendbuf, int *sendcounts, int *displs,
281                  MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm)
282 {
283   int retval = 0;
284
285   smpi_bench_end();
286
287   if (comm == MPI_COMM_NULL) {
288     retval = MPI_ERR_COMM;
289   } else if (sendcounts == nullptr || displs == nullptr) {
290     retval = MPI_ERR_ARG;
291   } else if (((comm->rank() == root) && (sendtype == MPI_DATATYPE_NULL)) ||
292              ((recvbuf != MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL))) {
293     retval = MPI_ERR_TYPE;
294   } else {
295     if (recvbuf == MPI_IN_PLACE) {
296       recvtype  = sendtype;
297       recvcount = sendcounts[comm->rank()];
298     }
299     int rank               = simgrid::s4u::this_actor::get_pid();
300     int dt_size_send       = sendtype->is_replayable() ? 1 : sendtype->size();
301
302     std::vector<int>* trace_sendcounts = new std::vector<int>;
303     if (comm->rank() == root) {
304       for (int i = 0; i < comm->size(); i++) // copy data to avoid bad free
305         trace_sendcounts->push_back(sendcounts[i] * dt_size_send);
306     }
307
308     TRACE_smpi_comm_in(rank, __func__,
309                        new simgrid::instr::VarCollTIData(
310                            "scatterV", root, dt_size_send, trace_sendcounts,
311                            recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(), nullptr,
312                            simgrid::smpi::Datatype::encode(sendtype), simgrid::smpi::Datatype::encode(recvtype)));
313
314     retval = simgrid::smpi::Colls::scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm);
315
316     TRACE_smpi_comm_out(rank);
317   }
318
319   smpi_bench_begin();
320   return retval;
321 }
322
323 int PMPI_Reduce(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
324 {
325   int retval = 0;
326
327   smpi_bench_end();
328
329   if (comm == MPI_COMM_NULL) {
330     retval = MPI_ERR_COMM;
331   } else if (not datatype->is_valid() || op == MPI_OP_NULL) {
332     retval = MPI_ERR_ARG;
333   } else {
334     int rank = simgrid::s4u::this_actor::get_pid();
335
336     TRACE_smpi_comm_in(rank, __func__,
337                        new simgrid::instr::CollTIData("reduce", root, 0,
338                                                       datatype->is_replayable() ? count : count * datatype->size(), -1,
339                                                       simgrid::smpi::Datatype::encode(datatype), ""));
340
341     simgrid::smpi::Colls::reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
342
343     retval = MPI_SUCCESS;
344     TRACE_smpi_comm_out(rank);
345   }
346
347   smpi_bench_begin();
348   return retval;
349 }
350
351 int PMPI_Reduce_local(void *inbuf, void *inoutbuf, int count, MPI_Datatype datatype, MPI_Op op){
352   int retval = 0;
353
354   smpi_bench_end();
355   if (not datatype->is_valid() || op == MPI_OP_NULL) {
356     retval = MPI_ERR_ARG;
357   } else {
358     op->apply(inbuf, inoutbuf, &count, datatype);
359     retval = MPI_SUCCESS;
360   }
361   smpi_bench_begin();
362   return retval;
363 }
364
365 int PMPI_Allreduce(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
366 {
367   int retval = 0;
368
369   smpi_bench_end();
370
371   if (comm == MPI_COMM_NULL) {
372     retval = MPI_ERR_COMM;
373   } else if (not datatype->is_valid()) {
374     retval = MPI_ERR_TYPE;
375   } else if (op == MPI_OP_NULL) {
376     retval = MPI_ERR_OP;
377   } else {
378
379     char* sendtmpbuf = static_cast<char*>(sendbuf);
380     if( sendbuf == MPI_IN_PLACE ) {
381       sendtmpbuf = static_cast<char*>(xbt_malloc(count*datatype->get_extent()));
382       simgrid::smpi::Datatype::copy(recvbuf, count, datatype,sendtmpbuf, count, datatype);
383     }
384     int rank = simgrid::s4u::this_actor::get_pid();
385
386     TRACE_smpi_comm_in(rank, __func__,
387                        new simgrid::instr::CollTIData("allreduce", -1, 0,
388                                                       datatype->is_replayable() ? count : count * datatype->size(), -1,
389                                                       simgrid::smpi::Datatype::encode(datatype), ""));
390
391     simgrid::smpi::Colls::allreduce(sendtmpbuf, recvbuf, count, datatype, op, comm);
392
393     if( sendbuf == MPI_IN_PLACE )
394       xbt_free(sendtmpbuf);
395
396     retval = MPI_SUCCESS;
397     TRACE_smpi_comm_out(rank);
398   }
399
400   smpi_bench_begin();
401   return retval;
402 }
403
404 int PMPI_Scan(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
405 {
406   int retval = 0;
407
408   smpi_bench_end();
409
410   if (comm == MPI_COMM_NULL) {
411     retval = MPI_ERR_COMM;
412   } else if (not datatype->is_valid()) {
413     retval = MPI_ERR_TYPE;
414   } else if (op == MPI_OP_NULL) {
415     retval = MPI_ERR_OP;
416   } else {
417     int rank = simgrid::s4u::this_actor::get_pid();
418     void* sendtmpbuf = sendbuf;
419     if (sendbuf == MPI_IN_PLACE) {
420       sendtmpbuf = static_cast<void*>(xbt_malloc(count * datatype->size()));
421       memcpy(sendtmpbuf, recvbuf, count * datatype->size());
422     }
423     TRACE_smpi_comm_in(rank, __func__, new simgrid::instr::Pt2PtTIData(
424                                            "scan", -1, datatype->is_replayable() ? count : count * datatype->size(),
425                                            simgrid::smpi::Datatype::encode(datatype)));
426
427     retval = simgrid::smpi::Colls::scan(sendtmpbuf, recvbuf, count, datatype, op, comm);
428
429     TRACE_smpi_comm_out(rank);
430     if (sendbuf == MPI_IN_PLACE)
431       xbt_free(sendtmpbuf);
432   }
433
434   smpi_bench_begin();
435   return retval;
436 }
437
438 int PMPI_Exscan(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm){
439   int retval = 0;
440
441   smpi_bench_end();
442
443   if (comm == MPI_COMM_NULL) {
444     retval = MPI_ERR_COMM;
445   } else if (not datatype->is_valid()) {
446     retval = MPI_ERR_TYPE;
447   } else if (op == MPI_OP_NULL) {
448     retval = MPI_ERR_OP;
449   } else {
450     int rank         = simgrid::s4u::this_actor::get_pid();
451     void* sendtmpbuf = sendbuf;
452     if (sendbuf == MPI_IN_PLACE) {
453       sendtmpbuf = static_cast<void*>(xbt_malloc(count * datatype->size()));
454       memcpy(sendtmpbuf, recvbuf, count * datatype->size());
455     }
456
457     TRACE_smpi_comm_in(rank, __func__, new simgrid::instr::Pt2PtTIData(
458                                            "exscan", -1, datatype->is_replayable() ? count : count * datatype->size(),
459                                            simgrid::smpi::Datatype::encode(datatype)));
460
461     retval = simgrid::smpi::Colls::exscan(sendtmpbuf, recvbuf, count, datatype, op, comm);
462
463     TRACE_smpi_comm_out(rank);
464     if (sendbuf == MPI_IN_PLACE)
465       xbt_free(sendtmpbuf);
466   }
467
468   smpi_bench_begin();
469   return retval;
470 }
471
472 int PMPI_Reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
473 {
474   int retval = 0;
475   smpi_bench_end();
476
477   if (comm == MPI_COMM_NULL) {
478     retval = MPI_ERR_COMM;
479   } else if (not datatype->is_valid()) {
480     retval = MPI_ERR_TYPE;
481   } else if (op == MPI_OP_NULL) {
482     retval = MPI_ERR_OP;
483   } else if (recvcounts == nullptr) {
484     retval = MPI_ERR_ARG;
485   } else {
486     int rank                           = simgrid::s4u::this_actor::get_pid();
487     std::vector<int>* trace_recvcounts = new std::vector<int>;
488     int dt_send_size                   = datatype->is_replayable() ? 1 : datatype->size();
489     int totalcount    = 0;
490
491     for (int i = 0; i < comm->size(); i++) { // copy data to avoid bad free
492       trace_recvcounts->push_back(recvcounts[i] * dt_send_size);
493       totalcount += recvcounts[i];
494     }
495
496     void* sendtmpbuf = sendbuf;
497     if (sendbuf == MPI_IN_PLACE) {
498       sendtmpbuf = static_cast<void*>(xbt_malloc(totalcount * datatype->size()));
499       memcpy(sendtmpbuf, recvbuf, totalcount * datatype->size());
500     }
501
502     TRACE_smpi_comm_in(rank, __func__, new simgrid::instr::VarCollTIData(
503                                            "reduceScatter", -1, dt_send_size, nullptr, -1, trace_recvcounts,
504                                            simgrid::smpi::Datatype::encode(datatype), ""));
505
506     simgrid::smpi::Colls::reduce_scatter(sendtmpbuf, recvbuf, recvcounts, datatype, op, comm);
507     retval = MPI_SUCCESS;
508     TRACE_smpi_comm_out(rank);
509
510     if (sendbuf == MPI_IN_PLACE)
511       xbt_free(sendtmpbuf);
512   }
513
514   smpi_bench_begin();
515   return retval;
516 }
517
518 int PMPI_Reduce_scatter_block(void *sendbuf, void *recvbuf, int recvcount,
519                               MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
520 {
521   int retval;
522   smpi_bench_end();
523
524   if (comm == MPI_COMM_NULL) {
525     retval = MPI_ERR_COMM;
526   } else if (not datatype->is_valid()) {
527     retval = MPI_ERR_TYPE;
528   } else if (op == MPI_OP_NULL) {
529     retval = MPI_ERR_OP;
530   } else if (recvcount < 0) {
531     retval = MPI_ERR_ARG;
532   } else {
533     int count = comm->size();
534
535     int rank                           = simgrid::s4u::this_actor::get_pid();
536     int dt_send_size                   = datatype->is_replayable() ? 1 : datatype->size();
537     std::vector<int>* trace_recvcounts = new std::vector<int>(recvcount * dt_send_size); // copy data to avoid bad free
538
539     void* sendtmpbuf       = sendbuf;
540     if (sendbuf == MPI_IN_PLACE) {
541       sendtmpbuf = static_cast<void*>(xbt_malloc(recvcount * count * datatype->size()));
542       memcpy(sendtmpbuf, recvbuf, recvcount * count * datatype->size());
543     }
544
545     TRACE_smpi_comm_in(rank, __func__,
546                        new simgrid::instr::VarCollTIData("reduceScatter", -1, 0, nullptr, -1, trace_recvcounts,
547                                                          simgrid::smpi::Datatype::encode(datatype), ""));
548
549     int* recvcounts = new int[count];
550     for (int i      = 0; i < count; i++)
551       recvcounts[i] = recvcount;
552     simgrid::smpi::Colls::reduce_scatter(sendtmpbuf, recvbuf, recvcounts, datatype, op, comm);
553     delete[] recvcounts;
554     retval = MPI_SUCCESS;
555
556     TRACE_smpi_comm_out(rank);
557
558     if (sendbuf == MPI_IN_PLACE)
559       xbt_free(sendtmpbuf);
560   }
561
562   smpi_bench_begin();
563   return retval;
564 }
565
566 int PMPI_Alltoall(void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
567                   MPI_Datatype recvtype, MPI_Comm comm)
568 {
569   int retval = 0;
570   smpi_bench_end();
571
572   if (comm == MPI_COMM_NULL) {
573     retval = MPI_ERR_COMM;
574   } else if ((sendbuf != MPI_IN_PLACE && sendtype == MPI_DATATYPE_NULL) || recvtype == MPI_DATATYPE_NULL) {
575     retval = MPI_ERR_TYPE;
576   } else {
577     int rank                 = simgrid::s4u::this_actor::get_pid();
578     void* sendtmpbuf         = static_cast<char*>(sendbuf);
579     int sendtmpcount         = sendcount;
580     MPI_Datatype sendtmptype = sendtype;
581     if (sendbuf == MPI_IN_PLACE) {
582       sendtmpbuf = static_cast<void*>(xbt_malloc(recvcount * comm->size() * recvtype->size()));
583       memcpy(sendtmpbuf, recvbuf, recvcount * comm->size() * recvtype->size());
584       sendtmpcount = recvcount;
585       sendtmptype  = recvtype;
586     }
587
588     TRACE_smpi_comm_in(rank, __func__,
589                        new simgrid::instr::CollTIData(
590                            "allToAll", -1, -1.0,
591                            sendtmptype->is_replayable() ? sendtmpcount : sendtmpcount * sendtmptype->size(),
592                            recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(),
593                            simgrid::smpi::Datatype::encode(sendtmptype), simgrid::smpi::Datatype::encode(recvtype)));
594
595     retval = simgrid::smpi::Colls::alltoall(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcount, recvtype, comm);
596
597     TRACE_smpi_comm_out(rank);
598
599     if (sendbuf == MPI_IN_PLACE)
600       xbt_free(sendtmpbuf);
601   }
602
603   smpi_bench_begin();
604   return retval;
605 }
606
607 int PMPI_Alltoallv(void* sendbuf, int* sendcounts, int* senddisps, MPI_Datatype sendtype, void* recvbuf,
608                    int* recvcounts, int* recvdisps, MPI_Datatype recvtype, MPI_Comm comm)
609 {
610   int retval = 0;
611
612   smpi_bench_end();
613
614   if (comm == MPI_COMM_NULL) {
615     retval = MPI_ERR_COMM;
616   } else if (sendtype == MPI_DATATYPE_NULL || recvtype == MPI_DATATYPE_NULL) {
617     retval = MPI_ERR_TYPE;
618   } else if ((sendbuf != MPI_IN_PLACE && (sendcounts == nullptr || senddisps == nullptr)) || recvcounts == nullptr ||
619              recvdisps == nullptr) {
620     retval = MPI_ERR_ARG;
621   } else {
622     int rank                           = simgrid::s4u::this_actor::get_pid();
623     int size               = comm->size();
624     int send_size                      = 0;
625     int recv_size                      = 0;
626     std::vector<int>* trace_sendcounts = new std::vector<int>;
627     std::vector<int>* trace_recvcounts = new std::vector<int>;
628     int dt_size_recv       = recvtype->size();
629
630     void* sendtmpbuf         = static_cast<char*>(sendbuf);
631     int* sendtmpcounts       = sendcounts;
632     int* sendtmpdisps        = senddisps;
633     MPI_Datatype sendtmptype = sendtype;
634     int maxsize              = 0;
635     for (int i = 0; i < size; i++) { // copy data to avoid bad free
636       recv_size += recvcounts[i] * dt_size_recv;
637       trace_recvcounts->push_back(recvcounts[i] * dt_size_recv);
638       if (((recvdisps[i] + recvcounts[i]) * dt_size_recv) > maxsize)
639         maxsize = (recvdisps[i] + recvcounts[i]) * dt_size_recv;
640     }
641
642     if (sendbuf == MPI_IN_PLACE) {
643       sendtmpbuf = static_cast<void*>(xbt_malloc(maxsize));
644       memcpy(sendtmpbuf, recvbuf, maxsize);
645       sendtmpcounts = static_cast<int*>(xbt_malloc(size * sizeof(int)));
646       memcpy(sendtmpcounts, recvcounts, size * sizeof(int));
647       sendtmpdisps = static_cast<int*>(xbt_malloc(size * sizeof(int)));
648       memcpy(sendtmpdisps, recvdisps, size * sizeof(int));
649       sendtmptype = recvtype;
650     }
651
652     int dt_size_send = sendtmptype->size();
653
654     for (int i = 0; i < size; i++) { // copy data to avoid bad free
655       send_size += sendtmpcounts[i] * dt_size_send;
656       trace_sendcounts->push_back(sendtmpcounts[i] * dt_size_send);
657     }
658
659     TRACE_smpi_comm_in(rank, __func__,
660                        new simgrid::instr::VarCollTIData("allToAllV", -1, send_size, trace_sendcounts, recv_size,
661                                                          trace_recvcounts, simgrid::smpi::Datatype::encode(sendtype),
662                                                          simgrid::smpi::Datatype::encode(recvtype)));
663
664     retval = simgrid::smpi::Colls::alltoallv(sendtmpbuf, sendtmpcounts, sendtmpdisps, sendtmptype, recvbuf, recvcounts,
665                                     recvdisps, recvtype, comm);
666     TRACE_smpi_comm_out(rank);
667
668     if (sendbuf == MPI_IN_PLACE) {
669       xbt_free(sendtmpbuf);
670       xbt_free(sendtmpcounts);
671       xbt_free(sendtmpdisps);
672     }
673   }
674
675   smpi_bench_begin();
676   return retval;
677 }