Logo AND Algorithmique Numérique Distribuée

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