Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
delete unused parameters
[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.h"
7 #include "smpi_comm.hpp"
8 #include "smpi_coll.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
16 /* PMPI User level calls */
17 extern "C" { // Obviously, the C MPI interface should use the C linkage
18
19 int PMPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
20 {
21   int retval = 0;
22
23   smpi_bench_end();
24
25   if (comm == MPI_COMM_NULL) {
26     retval = MPI_ERR_COMM;
27   } else if (not datatype->is_valid()) {
28     retval = MPI_ERR_ARG;
29   } else {
30     int rank        = comm != MPI_COMM_NULL ? smpi_process()->index() : -1;
31     int root_traced = comm->group()->index(root);
32
33     instr_extra_data extra = xbt_new0(s_instr_extra_data_t, 1);
34     extra->type            = TRACING_BCAST;
35     extra->root            = root_traced;
36     int known              = 0;
37     extra->datatype1       = encode_datatype(datatype, &known);
38     int dt_size_send       = 1;
39     if (known == 0)
40       dt_size_send   = datatype->size();
41     extra->send_size = count * dt_size_send;
42     TRACE_smpi_collective_in(rank, __FUNCTION__, extra);
43     if (comm->size() > 1)
44       simgrid::smpi::Colls::bcast(buf, count, datatype, root, comm);
45     retval = MPI_SUCCESS;
46
47     TRACE_smpi_collective_out(rank, __FUNCTION__);
48   }
49   smpi_bench_begin();
50   return retval;
51 }
52
53 int PMPI_Barrier(MPI_Comm comm)
54 {
55   int retval = 0;
56
57   smpi_bench_end();
58
59   if (comm == MPI_COMM_NULL) {
60     retval = MPI_ERR_COMM;
61   } else {
62     int rank               = comm != MPI_COMM_NULL ? smpi_process()->index() : -1;
63     instr_extra_data extra = xbt_new0(s_instr_extra_data_t, 1);
64     extra->type            = TRACING_BARRIER;
65     TRACE_smpi_collective_in(rank, __FUNCTION__, extra);
66
67     simgrid::smpi::Colls::barrier(comm);
68
69     //Barrier can be used to synchronize RMA calls. Finish all requests from comm before.
70     comm->finish_rma_calls();
71
72     retval = MPI_SUCCESS;
73
74     TRACE_smpi_collective_out(rank, __FUNCTION__);
75   }
76
77   smpi_bench_begin();
78   return retval;
79 }
80
81 int PMPI_Gather(void *sendbuf, int sendcount, MPI_Datatype sendtype,void *recvbuf, int recvcount, MPI_Datatype recvtype,
82                 int root, MPI_Comm comm)
83 {
84   int retval = 0;
85
86   smpi_bench_end();
87
88   if (comm == MPI_COMM_NULL) {
89     retval = MPI_ERR_COMM;
90   } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
91             ((comm->rank() == root) && (recvtype == MPI_DATATYPE_NULL))){
92     retval = MPI_ERR_TYPE;
93   } else if ((( sendbuf != MPI_IN_PLACE) && (sendcount <0)) || ((comm->rank() == root) && (recvcount <0))){
94     retval = MPI_ERR_COUNT;
95   } else {
96
97     char* sendtmpbuf = static_cast<char*>(sendbuf);
98     int sendtmpcount = sendcount;
99     MPI_Datatype sendtmptype = sendtype;
100     if( (comm->rank() == root) && (sendbuf == MPI_IN_PLACE )) {
101       sendtmpcount=0;
102       sendtmptype=recvtype;
103     }
104     int rank               = comm != MPI_COMM_NULL ? smpi_process()->index() : -1;
105     int root_traced        = comm->group()->index(root);
106     instr_extra_data extra = xbt_new0(s_instr_extra_data_t, 1);
107     extra->type            = TRACING_GATHER;
108     extra->root            = root_traced;
109     int known              = 0;
110     extra->datatype1       = encode_datatype(sendtmptype, &known);
111     int dt_size_send       = 1;
112     if (known == 0)
113       dt_size_send   = sendtmptype->size();
114     extra->send_size = sendtmpcount * dt_size_send;
115     extra->datatype2 = encode_datatype(recvtype, &known);
116     int dt_size_recv = 1;
117     if ((comm->rank() == root) && known == 0)
118       dt_size_recv   = recvtype->size();
119     extra->recv_size = recvcount * dt_size_recv;
120
121     TRACE_smpi_collective_in(rank, __FUNCTION__, extra);
122
123     simgrid::smpi::Colls::gather(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcount, recvtype, root, comm);
124
125     retval = MPI_SUCCESS;
126     TRACE_smpi_collective_out(rank, __FUNCTION__);
127   }
128
129   smpi_bench_begin();
130   return retval;
131 }
132
133 int PMPI_Gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, int *recvcounts, int *displs,
134                 MPI_Datatype recvtype, int root, MPI_Comm comm)
135 {
136   int retval = 0;
137
138   smpi_bench_end();
139
140   if (comm == MPI_COMM_NULL) {
141     retval = MPI_ERR_COMM;
142   } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
143             ((comm->rank() == root) && (recvtype == MPI_DATATYPE_NULL))){
144     retval = MPI_ERR_TYPE;
145   } else if (( sendbuf != MPI_IN_PLACE) && (sendcount <0)){
146     retval = MPI_ERR_COUNT;
147   } else if (recvcounts == nullptr || displs == nullptr) {
148     retval = MPI_ERR_ARG;
149   } else {
150     char* sendtmpbuf = static_cast<char*>(sendbuf);
151     int sendtmpcount = sendcount;
152     MPI_Datatype sendtmptype = sendtype;
153     if( (comm->rank() == root) && (sendbuf == MPI_IN_PLACE )) {
154       sendtmpcount=0;
155       sendtmptype=recvtype;
156     }
157
158     int rank               = comm != MPI_COMM_NULL ? smpi_process()->index() : -1;
159     int root_traced        = comm->group()->index(root);
160     int size               = comm->size();
161     instr_extra_data extra = xbt_new0(s_instr_extra_data_t, 1);
162     extra->type            = TRACING_GATHERV;
163     extra->num_processes   = size;
164     extra->root            = root_traced;
165     int known              = 0;
166     extra->datatype1       = encode_datatype(sendtmptype, &known);
167     int dt_size_send       = 1;
168     if (known == 0)
169       dt_size_send   = sendtype->size();
170     extra->send_size = sendtmpcount * dt_size_send;
171     extra->datatype2 = encode_datatype(recvtype, &known);
172     int dt_size_recv = 1;
173     if (known == 0)
174       dt_size_recv = recvtype->size();
175     if (comm->rank() == root) {
176       extra->recvcounts = xbt_new(int, size);
177       for (int i = 0; i < size; i++) // copy data to avoid bad free
178         extra->recvcounts[i] = recvcounts[i] * dt_size_recv;
179     }
180     TRACE_smpi_collective_in(rank, __FUNCTION__, extra);
181
182     retval = simgrid::smpi::Colls::gatherv(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcounts, displs, recvtype, root, comm);
183     TRACE_smpi_collective_out(rank, __FUNCTION__);
184   }
185
186   smpi_bench_begin();
187   return retval;
188 }
189
190 int PMPI_Allgather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
191                    void *recvbuf, int recvcount, MPI_Datatype recvtype, MPI_Comm comm)
192 {
193   int retval = 0;
194
195   smpi_bench_end();
196
197   if (comm == MPI_COMM_NULL) {
198     retval = MPI_ERR_COMM;
199   } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
200             (recvtype == MPI_DATATYPE_NULL)){
201     retval = MPI_ERR_TYPE;
202   } else if ((( sendbuf != MPI_IN_PLACE) && (sendcount <0)) ||
203             (recvcount <0)){
204     retval = MPI_ERR_COUNT;
205   } else {
206     if(sendbuf == MPI_IN_PLACE) {
207       sendbuf=static_cast<char*>(recvbuf)+recvtype->get_extent()*recvcount*comm->rank();
208       sendcount=recvcount;
209       sendtype=recvtype;
210     }
211     int rank               = comm != MPI_COMM_NULL ? smpi_process()->index() : -1;
212     instr_extra_data extra = xbt_new0(s_instr_extra_data_t, 1);
213     extra->type            = TRACING_ALLGATHER;
214     int known              = 0;
215     extra->datatype1       = encode_datatype(sendtype, &known);
216     int dt_size_send       = 1;
217     if (known == 0)
218       dt_size_send   = sendtype->size();
219     extra->send_size = sendcount * dt_size_send;
220     extra->datatype2 = encode_datatype(recvtype, &known);
221     int dt_size_recv = 1;
222     if (known == 0)
223       dt_size_recv   = recvtype->size();
224     extra->recv_size = recvcount * dt_size_recv;
225
226     TRACE_smpi_collective_in(rank, __FUNCTION__, extra);
227
228     simgrid::smpi::Colls::allgather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
229     retval = MPI_SUCCESS;
230     TRACE_smpi_collective_out(rank, __FUNCTION__);
231   }
232   smpi_bench_begin();
233   return retval;
234 }
235
236 int PMPI_Allgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
237                    void *recvbuf, int *recvcounts, int *displs, MPI_Datatype recvtype, MPI_Comm comm)
238 {
239   int retval = 0;
240
241   smpi_bench_end();
242
243   if (comm == MPI_COMM_NULL) {
244     retval = MPI_ERR_COMM;
245   } else if (((sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) || (recvtype == MPI_DATATYPE_NULL)) {
246     retval = MPI_ERR_TYPE;
247   } else if (( sendbuf != MPI_IN_PLACE) && (sendcount <0)){
248     retval = MPI_ERR_COUNT;
249   } else if (recvcounts == nullptr || displs == nullptr) {
250     retval = MPI_ERR_ARG;
251   } else {
252
253     if(sendbuf == MPI_IN_PLACE) {
254       sendbuf=static_cast<char*>(recvbuf)+recvtype->get_extent()*displs[comm->rank()];
255       sendcount=recvcounts[comm->rank()];
256       sendtype=recvtype;
257     }
258     int rank               = comm != MPI_COMM_NULL ? smpi_process()->index() : -1;
259     int i                  = 0;
260     int size               = comm->size();
261     instr_extra_data extra = xbt_new0(s_instr_extra_data_t, 1);
262     extra->type            = TRACING_ALLGATHERV;
263     extra->num_processes   = size;
264     int known              = 0;
265     extra->datatype1       = encode_datatype(sendtype, &known);
266     int dt_size_send       = 1;
267     if (known == 0)
268       dt_size_send   = sendtype->size();
269     extra->send_size = sendcount * dt_size_send;
270     extra->datatype2 = encode_datatype(recvtype, &known);
271     int dt_size_recv = 1;
272     if (known == 0)
273       dt_size_recv    = recvtype->size();
274     extra->recvcounts = xbt_new(int, size);
275     for (i                 = 0; i < size; i++) // copy data to avoid bad free
276       extra->recvcounts[i] = recvcounts[i] * dt_size_recv;
277
278     TRACE_smpi_collective_in(rank, __FUNCTION__, extra);
279
280     simgrid::smpi::Colls::allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm);
281     retval = MPI_SUCCESS;
282     TRACE_smpi_collective_out(rank, __FUNCTION__);
283   }
284
285   smpi_bench_begin();
286   return retval;
287 }
288
289 int PMPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
290                 void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm)
291 {
292   int retval = 0;
293
294   smpi_bench_end();
295
296   if (comm == MPI_COMM_NULL) {
297     retval = MPI_ERR_COMM;
298   } else if (((comm->rank() == root) && (not sendtype->is_valid())) ||
299              ((recvbuf != MPI_IN_PLACE) && (not recvtype->is_valid()))) {
300     retval = MPI_ERR_TYPE;
301   } else if ((sendbuf == recvbuf) ||
302       ((comm->rank()==root) && sendcount>0 && (sendbuf == nullptr))){
303     retval = MPI_ERR_BUFFER;
304   }else {
305
306     if (recvbuf == MPI_IN_PLACE) {
307       recvtype  = sendtype;
308       recvcount = sendcount;
309     }
310     int rank               = comm != MPI_COMM_NULL ? smpi_process()->index() : -1;
311     int root_traced        = comm->group()->index(root);
312     instr_extra_data extra = xbt_new0(s_instr_extra_data_t, 1);
313     extra->type            = TRACING_SCATTER;
314     extra->root            = root_traced;
315     int known              = 0;
316     extra->datatype1       = encode_datatype(sendtype, &known);
317     int dt_size_send       = 1;
318     if ((comm->rank() == root) && known == 0)
319       dt_size_send   = sendtype->size();
320     extra->send_size = sendcount * dt_size_send;
321     extra->datatype2 = encode_datatype(recvtype, &known);
322     int dt_size_recv = 1;
323     if (known == 0)
324       dt_size_recv   = recvtype->size();
325     extra->recv_size = recvcount * dt_size_recv;
326     TRACE_smpi_collective_in(rank, __FUNCTION__, extra);
327
328     simgrid::smpi::Colls::scatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm);
329     retval = MPI_SUCCESS;
330     TRACE_smpi_collective_out(rank, __FUNCTION__);
331   }
332
333   smpi_bench_begin();
334   return retval;
335 }
336
337 int PMPI_Scatterv(void *sendbuf, int *sendcounts, int *displs,
338                  MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm)
339 {
340   int retval = 0;
341
342   smpi_bench_end();
343
344   if (comm == MPI_COMM_NULL) {
345     retval = MPI_ERR_COMM;
346   } else if (sendcounts == nullptr || displs == nullptr) {
347     retval = MPI_ERR_ARG;
348   } else if (((comm->rank() == root) && (sendtype == MPI_DATATYPE_NULL)) ||
349              ((recvbuf != MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL))) {
350     retval = MPI_ERR_TYPE;
351   } else {
352     if (recvbuf == MPI_IN_PLACE) {
353       recvtype  = sendtype;
354       recvcount = sendcounts[comm->rank()];
355     }
356     int rank               = comm != MPI_COMM_NULL ? smpi_process()->index() : -1;
357     int root_traced        = comm->group()->index(root);
358     int size               = comm->size();
359     instr_extra_data extra = xbt_new0(s_instr_extra_data_t, 1);
360     extra->type            = TRACING_SCATTERV;
361     extra->num_processes   = size;
362     extra->root            = root_traced;
363     int known              = 0;
364     extra->datatype1       = encode_datatype(sendtype, &known);
365     int dt_size_send       = 1;
366     if (known == 0)
367       dt_size_send = sendtype->size();
368     if (comm->rank() == root) {
369       extra->sendcounts = xbt_new(int, size);
370       for (int i = 0; i < size; i++) // copy data to avoid bad free
371         extra->sendcounts[i] = sendcounts[i] * dt_size_send;
372     }
373     extra->datatype2 = encode_datatype(recvtype, &known);
374     int dt_size_recv = 1;
375     if (known == 0)
376       dt_size_recv   = recvtype->size();
377     extra->recv_size = recvcount * dt_size_recv;
378     TRACE_smpi_collective_in(rank, __FUNCTION__, extra);
379
380     retval = simgrid::smpi::Colls::scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm);
381
382     TRACE_smpi_collective_out(rank, __FUNCTION__);
383   }
384
385   smpi_bench_begin();
386   return retval;
387 }
388
389 int PMPI_Reduce(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
390 {
391   int retval = 0;
392
393   smpi_bench_end();
394
395   if (comm == MPI_COMM_NULL) {
396     retval = MPI_ERR_COMM;
397   } else if (not datatype->is_valid() || op == MPI_OP_NULL) {
398     retval = MPI_ERR_ARG;
399   } else {
400     int rank               = comm != MPI_COMM_NULL ? smpi_process()->index() : -1;
401     int root_traced        = comm->group()->index(root);
402     instr_extra_data extra = xbt_new0(s_instr_extra_data_t, 1);
403     extra->type            = TRACING_REDUCE;
404     int known              = 0;
405     extra->datatype1       = encode_datatype(datatype, &known);
406     int dt_size_send       = 1;
407     if (known == 0)
408       dt_size_send   = datatype->size();
409     extra->send_size = count * dt_size_send;
410     extra->root      = root_traced;
411
412     TRACE_smpi_collective_in(rank, __FUNCTION__, extra);
413
414     simgrid::smpi::Colls::reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
415
416     retval = MPI_SUCCESS;
417     TRACE_smpi_collective_out(rank, __FUNCTION__);
418   }
419
420   smpi_bench_begin();
421   return retval;
422 }
423
424 int PMPI_Reduce_local(void *inbuf, void *inoutbuf, int count, MPI_Datatype datatype, MPI_Op op){
425   int retval = 0;
426
427   smpi_bench_end();
428   if (not datatype->is_valid() || op == MPI_OP_NULL) {
429     retval = MPI_ERR_ARG;
430   } else {
431     op->apply(inbuf, inoutbuf, &count, datatype);
432     retval = MPI_SUCCESS;
433   }
434   smpi_bench_begin();
435   return retval;
436 }
437
438 int PMPI_Allreduce(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
439 {
440   int retval = 0;
441
442   smpi_bench_end();
443
444   if (comm == MPI_COMM_NULL) {
445     retval = MPI_ERR_COMM;
446   } else if (not datatype->is_valid()) {
447     retval = MPI_ERR_TYPE;
448   } else if (op == MPI_OP_NULL) {
449     retval = MPI_ERR_OP;
450   } else {
451
452     char* sendtmpbuf = static_cast<char*>(sendbuf);
453     if( sendbuf == MPI_IN_PLACE ) {
454       sendtmpbuf = static_cast<char*>(xbt_malloc(count*datatype->get_extent()));
455       simgrid::smpi::Datatype::copy(recvbuf, count, datatype,sendtmpbuf, count, datatype);
456     }
457     int rank               = comm != MPI_COMM_NULL ? smpi_process()->index() : -1;
458     instr_extra_data extra = xbt_new0(s_instr_extra_data_t, 1);
459     extra->type            = TRACING_ALLREDUCE;
460     int known              = 0;
461     extra->datatype1       = encode_datatype(datatype, &known);
462     int dt_size_send       = 1;
463     if (known == 0)
464       dt_size_send   = datatype->size();
465     extra->send_size = count * dt_size_send;
466
467     TRACE_smpi_collective_in(rank, __FUNCTION__, extra);
468
469     simgrid::smpi::Colls::allreduce(sendtmpbuf, recvbuf, count, datatype, op, comm);
470
471     if( sendbuf == MPI_IN_PLACE )
472       xbt_free(sendtmpbuf);
473
474     retval = MPI_SUCCESS;
475     TRACE_smpi_collective_out(rank, __FUNCTION__);
476   }
477
478   smpi_bench_begin();
479   return retval;
480 }
481
482 int PMPI_Scan(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
483 {
484   int retval = 0;
485
486   smpi_bench_end();
487
488   if (comm == MPI_COMM_NULL) {
489     retval = MPI_ERR_COMM;
490   } else if (not datatype->is_valid()) {
491     retval = MPI_ERR_TYPE;
492   } else if (op == MPI_OP_NULL) {
493     retval = MPI_ERR_OP;
494   } else {
495     int rank               = comm != MPI_COMM_NULL ? smpi_process()->index() : -1;
496     instr_extra_data extra = xbt_new0(s_instr_extra_data_t, 1);
497     extra->type            = TRACING_SCAN;
498     int known              = 0;
499     extra->datatype1       = encode_datatype(datatype, &known);
500     int dt_size_send       = 1;
501     if (known == 0)
502       dt_size_send   = datatype->size();
503     extra->send_size = count * dt_size_send;
504
505     TRACE_smpi_collective_in(rank, __FUNCTION__, extra);
506
507     retval = simgrid::smpi::Colls::scan(sendbuf, recvbuf, count, datatype, op, comm);
508
509     TRACE_smpi_collective_out(rank, __FUNCTION__);
510   }
511
512   smpi_bench_begin();
513   return retval;
514 }
515
516 int PMPI_Exscan(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm){
517   int retval = 0;
518
519   smpi_bench_end();
520
521   if (comm == MPI_COMM_NULL) {
522     retval = MPI_ERR_COMM;
523   } else if (not datatype->is_valid()) {
524     retval = MPI_ERR_TYPE;
525   } else if (op == MPI_OP_NULL) {
526     retval = MPI_ERR_OP;
527   } else {
528     int rank               = comm != MPI_COMM_NULL ? smpi_process()->index() : -1;
529     instr_extra_data extra = xbt_new0(s_instr_extra_data_t, 1);
530     extra->type            = TRACING_EXSCAN;
531     int known              = 0;
532     extra->datatype1       = encode_datatype(datatype, &known);
533     int dt_size_send       = 1;
534     if (known == 0)
535       dt_size_send   = datatype->size();
536     extra->send_size = count * dt_size_send;
537     void* sendtmpbuf = sendbuf;
538     if (sendbuf == MPI_IN_PLACE) {
539       sendtmpbuf = static_cast<void*>(xbt_malloc(count * datatype->size()));
540       memcpy(sendtmpbuf, recvbuf, count * datatype->size());
541     }
542     TRACE_smpi_collective_in(rank, __FUNCTION__, extra);
543
544     retval = simgrid::smpi::Colls::exscan(sendtmpbuf, recvbuf, count, datatype, op, comm);
545
546     TRACE_smpi_collective_out(rank, __FUNCTION__);
547     if (sendbuf == MPI_IN_PLACE)
548       xbt_free(sendtmpbuf);
549   }
550
551   smpi_bench_begin();
552   return retval;
553 }
554
555 int PMPI_Reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
556 {
557   int retval = 0;
558   smpi_bench_end();
559
560   if (comm == MPI_COMM_NULL) {
561     retval = MPI_ERR_COMM;
562   } else if (not datatype->is_valid()) {
563     retval = MPI_ERR_TYPE;
564   } else if (op == MPI_OP_NULL) {
565     retval = MPI_ERR_OP;
566   } else if (recvcounts == nullptr) {
567     retval = MPI_ERR_ARG;
568   } else {
569     int rank               = comm != MPI_COMM_NULL ? smpi_process()->index() : -1;
570     int i                  = 0;
571     int size               = comm->size();
572     instr_extra_data extra = xbt_new0(s_instr_extra_data_t, 1);
573     extra->type            = TRACING_REDUCE_SCATTER;
574     extra->num_processes   = size;
575     int known              = 0;
576     extra->datatype1       = encode_datatype(datatype, &known);
577     int dt_size_send       = 1;
578     if (known == 0)
579       dt_size_send    = datatype->size();
580     extra->send_size  = 0;
581     extra->recvcounts = xbt_new(int, size);
582     int totalcount    = 0;
583     for (i = 0; i < size; i++) { // copy data to avoid bad free
584       extra->recvcounts[i] = recvcounts[i] * dt_size_send;
585       totalcount += recvcounts[i];
586     }
587     void* sendtmpbuf = sendbuf;
588     if (sendbuf == MPI_IN_PLACE) {
589       sendtmpbuf = static_cast<void*>(xbt_malloc(totalcount * datatype->size()));
590       memcpy(sendtmpbuf, recvbuf, totalcount * datatype->size());
591     }
592
593     TRACE_smpi_collective_in(rank, __FUNCTION__, extra);
594
595     simgrid::smpi::Colls::reduce_scatter(sendtmpbuf, recvbuf, recvcounts, datatype, op, comm);
596     retval = MPI_SUCCESS;
597     TRACE_smpi_collective_out(rank, __FUNCTION__);
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_Reduce_scatter_block(void *sendbuf, void *recvbuf, int recvcount,
608                               MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
609 {
610   int retval;
611   smpi_bench_end();
612
613   if (comm == MPI_COMM_NULL) {
614     retval = MPI_ERR_COMM;
615   } else if (not datatype->is_valid()) {
616     retval = MPI_ERR_TYPE;
617   } else if (op == MPI_OP_NULL) {
618     retval = MPI_ERR_OP;
619   } else if (recvcount < 0) {
620     retval = MPI_ERR_ARG;
621   } else {
622     int count = comm->size();
623
624     int rank               = comm != MPI_COMM_NULL ? smpi_process()->index() : -1;
625     instr_extra_data extra = xbt_new0(s_instr_extra_data_t, 1);
626     extra->type            = TRACING_REDUCE_SCATTER;
627     extra->num_processes   = count;
628     int known              = 0;
629     extra->datatype1       = encode_datatype(datatype, &known);
630     int dt_size_send       = 1;
631     if (known == 0)
632       dt_size_send    = datatype->size();
633     extra->send_size  = 0;
634     extra->recvcounts = xbt_new(int, count);
635     for (int i             = 0; i < count; i++) // copy data to avoid bad free
636       extra->recvcounts[i] = recvcount * dt_size_send;
637     void* sendtmpbuf       = sendbuf;
638     if (sendbuf == MPI_IN_PLACE) {
639       sendtmpbuf = static_cast<void*>(xbt_malloc(recvcount * count * datatype->size()));
640       memcpy(sendtmpbuf, recvbuf, recvcount * count * datatype->size());
641     }
642
643     TRACE_smpi_collective_in(rank, __FUNCTION__, extra);
644
645     int* recvcounts = static_cast<int*>(xbt_malloc(count * sizeof(int)));
646     for (int i      = 0; i < count; i++)
647       recvcounts[i] = recvcount;
648     simgrid::smpi::Colls::reduce_scatter(sendtmpbuf, recvbuf, recvcounts, datatype, op, comm);
649     xbt_free(recvcounts);
650     retval = MPI_SUCCESS;
651
652     TRACE_smpi_collective_out(rank, __FUNCTION__);
653
654     if (sendbuf == MPI_IN_PLACE)
655       xbt_free(sendtmpbuf);
656   }
657
658   smpi_bench_begin();
659   return retval;
660 }
661
662 int PMPI_Alltoall(void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
663                   MPI_Datatype recvtype, MPI_Comm comm)
664 {
665   int retval = 0;
666   smpi_bench_end();
667
668   if (comm == MPI_COMM_NULL) {
669     retval = MPI_ERR_COMM;
670   } else if ((sendbuf != MPI_IN_PLACE && sendtype == MPI_DATATYPE_NULL) || recvtype == MPI_DATATYPE_NULL) {
671     retval = MPI_ERR_TYPE;
672   } else {
673     int rank               = comm != MPI_COMM_NULL ? smpi_process()->index() : -1;
674     instr_extra_data extra = xbt_new0(s_instr_extra_data_t, 1);
675     extra->type            = TRACING_ALLTOALL;
676
677     void* sendtmpbuf         = static_cast<char*>(sendbuf);
678     int sendtmpcount         = sendcount;
679     MPI_Datatype sendtmptype = sendtype;
680     if (sendbuf == MPI_IN_PLACE) {
681       sendtmpbuf = static_cast<void*>(xbt_malloc(recvcount * comm->size() * recvtype->size()));
682       memcpy(sendtmpbuf, recvbuf, recvcount * comm->size() * recvtype->size());
683       sendtmpcount = recvcount;
684       sendtmptype  = recvtype;
685     }
686
687     int known        = 0;
688     extra->datatype1 = encode_datatype(sendtmptype, &known);
689     if (known == 0)
690       extra->send_size = sendtmpcount * sendtmptype->size();
691     else
692       extra->send_size = sendtmpcount;
693     extra->datatype2   = encode_datatype(recvtype, &known);
694     if (known == 0)
695       extra->recv_size = recvcount * recvtype->size();
696     else
697       extra->recv_size = recvcount;
698
699     TRACE_smpi_collective_in(rank, __FUNCTION__, extra);
700
701     retval = simgrid::smpi::Colls::alltoall(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcount, recvtype, comm);
702
703     TRACE_smpi_collective_out(rank, __FUNCTION__);
704
705     if (sendbuf == MPI_IN_PLACE)
706       xbt_free(sendtmpbuf);
707   }
708
709   smpi_bench_begin();
710   return retval;
711 }
712
713 int PMPI_Alltoallv(void* sendbuf, int* sendcounts, int* senddisps, MPI_Datatype sendtype, void* recvbuf,
714                    int* recvcounts, int* recvdisps, MPI_Datatype recvtype, MPI_Comm comm)
715 {
716   int retval = 0;
717
718   smpi_bench_end();
719
720   if (comm == MPI_COMM_NULL) {
721     retval = MPI_ERR_COMM;
722   } else if (sendtype == MPI_DATATYPE_NULL || recvtype == MPI_DATATYPE_NULL) {
723     retval = MPI_ERR_TYPE;
724   } else if ((sendbuf != MPI_IN_PLACE && (sendcounts == nullptr || senddisps == nullptr)) || recvcounts == nullptr ||
725              recvdisps == nullptr) {
726     retval = MPI_ERR_ARG;
727   } else {
728     int rank               = comm != MPI_COMM_NULL ? smpi_process()->index() : -1;
729     int i                  = 0;
730     int size               = comm->size();
731     instr_extra_data extra = xbt_new0(s_instr_extra_data_t, 1);
732     extra->type            = TRACING_ALLTOALLV;
733     extra->send_size       = 0;
734     extra->recv_size       = 0;
735     extra->recvcounts      = xbt_new(int, size);
736     extra->sendcounts      = xbt_new(int, size);
737     int known              = 0;
738     extra->datatype2       = encode_datatype(recvtype, &known);
739     int dt_size_recv       = recvtype->size();
740
741     void* sendtmpbuf         = static_cast<char*>(sendbuf);
742     int* sendtmpcounts       = sendcounts;
743     int* sendtmpdisps        = senddisps;
744     MPI_Datatype sendtmptype = sendtype;
745     int maxsize              = 0;
746     for (i = 0; i < size; i++) { // copy data to avoid bad free
747       extra->recv_size += recvcounts[i] * dt_size_recv;
748       extra->recvcounts[i] = recvcounts[i] * dt_size_recv;
749       if (((recvdisps[i] + recvcounts[i]) * dt_size_recv) > maxsize)
750         maxsize = (recvdisps[i] + recvcounts[i]) * dt_size_recv;
751     }
752
753     if (sendbuf == MPI_IN_PLACE) {
754       sendtmpbuf = static_cast<void*>(xbt_malloc(maxsize));
755       memcpy(sendtmpbuf, recvbuf, maxsize);
756       sendtmpcounts = static_cast<int*>(xbt_malloc(size * sizeof(int)));
757       memcpy(sendtmpcounts, recvcounts, size * sizeof(int));
758       sendtmpdisps = static_cast<int*>(xbt_malloc(size * sizeof(int)));
759       memcpy(sendtmpdisps, recvdisps, size * sizeof(int));
760       sendtmptype = recvtype;
761     }
762
763     extra->datatype1 = encode_datatype(sendtmptype, &known);
764     int dt_size_send = sendtmptype->size();
765
766     for (i = 0; i < size; i++) { // copy data to avoid bad free
767       extra->send_size += sendtmpcounts[i] * dt_size_send;
768       extra->sendcounts[i] = sendtmpcounts[i] * dt_size_send;
769     }
770     extra->num_processes = size;
771     TRACE_smpi_collective_in(rank, __FUNCTION__, extra);
772     retval = simgrid::smpi::Colls::alltoallv(sendtmpbuf, sendtmpcounts, sendtmpdisps, sendtmptype, recvbuf, recvcounts,
773                                     recvdisps, recvtype, comm);
774     TRACE_smpi_collective_out(rank, __FUNCTION__);
775
776     if (sendbuf == MPI_IN_PLACE) {
777       xbt_free(sendtmpbuf);
778       xbt_free(sendtmpcounts);
779       xbt_free(sendtmpdisps);
780     }
781   }
782
783   smpi_bench_begin();
784   return retval;
785 }
786
787 }