Logo AND Algorithmique Numérique Distribuée

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