Logo AND Algorithmique Numérique Distribuée

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