Logo AND Algorithmique Numérique Distribuée

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