Logo AND Algorithmique Numérique Distribuée

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