Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
bd0ac325a3a18dd0f70027db2ddceadad0762cb1
[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   return PMPI_Ireduce(sendbuf, recvbuf, count, datatype, op, root, comm, MPI_REQUEST_IGNORED);
394 }
395
396 int PMPI_Ireduce(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm, MPI_Request* request)
397 {
398   int retval = 0;
399
400   smpi_bench_end();
401
402   if (comm == MPI_COMM_NULL) {
403     retval = MPI_ERR_COMM;
404   } else if (not datatype->is_valid() || op == MPI_OP_NULL) {
405     retval = MPI_ERR_ARG;
406   } else {
407     int rank = simgrid::s4u::this_actor::get_pid();
408
409     TRACE_smpi_comm_in(rank, request==MPI_REQUEST_IGNORED ? "PMPI_Reduce":"PMPI_Ireduce",
410                        new simgrid::instr::CollTIData(request==MPI_REQUEST_IGNORED ? "reduce":"ireduce", root, 0,
411                                                       datatype->is_replayable() ? count : count * datatype->size(), -1,
412                                                       simgrid::smpi::Datatype::encode(datatype), ""));
413     if(request == MPI_REQUEST_IGNORED)
414       simgrid::smpi::Colls::reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
415     else
416       simgrid::smpi::Colls::ireduce(sendbuf, recvbuf, count, datatype, op, root, comm, request);
417
418
419     retval = MPI_SUCCESS;
420     TRACE_smpi_comm_out(rank);
421   }
422
423   smpi_bench_begin();
424   return retval;
425 }
426
427 int PMPI_Reduce_local(void *inbuf, void *inoutbuf, int count, MPI_Datatype datatype, MPI_Op op){
428   int retval = 0;
429
430   smpi_bench_end();
431   if (not datatype->is_valid() || op == MPI_OP_NULL) {
432     retval = MPI_ERR_ARG;
433   } else {
434     op->apply(inbuf, inoutbuf, &count, datatype);
435     retval = MPI_SUCCESS;
436   }
437   smpi_bench_begin();
438   return retval;
439 }
440
441 int PMPI_Allreduce(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
442 {
443   return PMPI_Iallreduce(sendbuf, recvbuf, count, datatype, op, comm, MPI_REQUEST_IGNORED);
444 }
445
446 int PMPI_Iallreduce(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request *request)
447 {
448   int retval = 0;
449
450   smpi_bench_end();
451
452   if (comm == MPI_COMM_NULL) {
453     retval = MPI_ERR_COMM;
454   } else if (not datatype->is_valid()) {
455     retval = MPI_ERR_TYPE;
456   } else if (op == MPI_OP_NULL) {
457     retval = MPI_ERR_OP;
458   } else if (request != MPI_REQUEST_IGNORED) {
459     xbt_die("Iallreduce is not yet implemented. WIP");
460     retval = MPI_ERR_ARG;
461   } else {
462
463     char* sendtmpbuf = static_cast<char*>(sendbuf);
464     if( sendbuf == MPI_IN_PLACE ) {
465       sendtmpbuf = static_cast<char*>(xbt_malloc(count*datatype->get_extent()));
466       simgrid::smpi::Datatype::copy(recvbuf, count, datatype,sendtmpbuf, count, datatype);
467     }
468     int rank = simgrid::s4u::this_actor::get_pid();
469
470     TRACE_smpi_comm_in(rank, __func__,
471                        new simgrid::instr::CollTIData(request==MPI_REQUEST_IGNORED ? "allreduce":"iallreduce", -1, 0,
472                                                       datatype->is_replayable() ? count : count * datatype->size(), -1,
473                                                       simgrid::smpi::Datatype::encode(datatype), ""));
474
475 //    if(request == MPI_REQUEST_IGNORED)
476       simgrid::smpi::Colls::allreduce(sendtmpbuf, recvbuf, count, datatype, op, comm);
477 //    else
478 //      simgrid::smpi::Colls::iallreduce(sendtmpbuf, recvbuf, count, datatype, op, comm, request);
479
480     if( sendbuf == MPI_IN_PLACE )
481       xbt_free(sendtmpbuf);
482
483     retval = MPI_SUCCESS;
484     TRACE_smpi_comm_out(rank);
485   }
486
487   smpi_bench_begin();
488   return retval;
489 }
490
491 int PMPI_Scan(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
492 {
493   int retval = 0;
494
495   smpi_bench_end();
496
497   if (comm == MPI_COMM_NULL) {
498     retval = MPI_ERR_COMM;
499   } else if (not datatype->is_valid()) {
500     retval = MPI_ERR_TYPE;
501   } else if (op == MPI_OP_NULL) {
502     retval = MPI_ERR_OP;
503   } else {
504     int rank = simgrid::s4u::this_actor::get_pid();
505     void* sendtmpbuf = sendbuf;
506     if (sendbuf == MPI_IN_PLACE) {
507       sendtmpbuf = static_cast<void*>(xbt_malloc(count * datatype->size()));
508       memcpy(sendtmpbuf, recvbuf, count * datatype->size());
509     }
510     TRACE_smpi_comm_in(rank, __func__, new simgrid::instr::Pt2PtTIData(
511                                            "scan", -1, datatype->is_replayable() ? count : count * datatype->size(),
512                                            simgrid::smpi::Datatype::encode(datatype)));
513
514     retval = simgrid::smpi::Colls::scan(sendtmpbuf, recvbuf, count, datatype, op, comm);
515
516     TRACE_smpi_comm_out(rank);
517     if (sendbuf == MPI_IN_PLACE)
518       xbt_free(sendtmpbuf);
519   }
520
521   smpi_bench_begin();
522   return retval;
523 }
524
525 int PMPI_Exscan(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm){
526   int retval = 0;
527
528   smpi_bench_end();
529
530   if (comm == MPI_COMM_NULL) {
531     retval = MPI_ERR_COMM;
532   } else if (not datatype->is_valid()) {
533     retval = MPI_ERR_TYPE;
534   } else if (op == MPI_OP_NULL) {
535     retval = MPI_ERR_OP;
536   } else {
537     int rank         = simgrid::s4u::this_actor::get_pid();
538     void* sendtmpbuf = sendbuf;
539     if (sendbuf == MPI_IN_PLACE) {
540       sendtmpbuf = static_cast<void*>(xbt_malloc(count * datatype->size()));
541       memcpy(sendtmpbuf, recvbuf, count * datatype->size());
542     }
543
544     TRACE_smpi_comm_in(rank, __func__, new simgrid::instr::Pt2PtTIData(
545                                            "exscan", -1, datatype->is_replayable() ? count : count * datatype->size(),
546                                            simgrid::smpi::Datatype::encode(datatype)));
547
548     retval = simgrid::smpi::Colls::exscan(sendtmpbuf, recvbuf, count, datatype, op, comm);
549
550     TRACE_smpi_comm_out(rank);
551     if (sendbuf == MPI_IN_PLACE)
552       xbt_free(sendtmpbuf);
553   }
554
555   smpi_bench_begin();
556   return retval;
557 }
558
559 int PMPI_Reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
560 {
561   int retval = 0;
562   smpi_bench_end();
563
564   if (comm == MPI_COMM_NULL) {
565     retval = MPI_ERR_COMM;
566   } else if (not datatype->is_valid()) {
567     retval = MPI_ERR_TYPE;
568   } else if (op == MPI_OP_NULL) {
569     retval = MPI_ERR_OP;
570   } else if (recvcounts == nullptr) {
571     retval = MPI_ERR_ARG;
572   } else {
573     int rank                           = simgrid::s4u::this_actor::get_pid();
574     std::vector<int>* trace_recvcounts = new std::vector<int>;
575     int dt_send_size                   = datatype->is_replayable() ? 1 : datatype->size();
576     int totalcount    = 0;
577
578     for (int i = 0; i < comm->size(); i++) { // copy data to avoid bad free
579       trace_recvcounts->push_back(recvcounts[i] * dt_send_size);
580       totalcount += recvcounts[i];
581     }
582
583     void* sendtmpbuf = sendbuf;
584     if (sendbuf == MPI_IN_PLACE) {
585       sendtmpbuf = static_cast<void*>(xbt_malloc(totalcount * datatype->size()));
586       memcpy(sendtmpbuf, recvbuf, totalcount * datatype->size());
587     }
588
589     TRACE_smpi_comm_in(rank, __func__, new simgrid::instr::VarCollTIData(
590                                            "reducescatter", -1, dt_send_size, nullptr, -1, trace_recvcounts,
591                                            simgrid::smpi::Datatype::encode(datatype), ""));
592
593     simgrid::smpi::Colls::reduce_scatter(sendtmpbuf, recvbuf, recvcounts, datatype, op, comm);
594     retval = MPI_SUCCESS;
595     TRACE_smpi_comm_out(rank);
596
597     if (sendbuf == MPI_IN_PLACE)
598       xbt_free(sendtmpbuf);
599   }
600
601   smpi_bench_begin();
602   return retval;
603 }
604
605 int PMPI_Reduce_scatter_block(void *sendbuf, void *recvbuf, int recvcount,
606                               MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
607 {
608   int retval;
609   smpi_bench_end();
610
611   if (comm == MPI_COMM_NULL) {
612     retval = MPI_ERR_COMM;
613   } else if (not datatype->is_valid()) {
614     retval = MPI_ERR_TYPE;
615   } else if (op == MPI_OP_NULL) {
616     retval = MPI_ERR_OP;
617   } else if (recvcount < 0) {
618     retval = MPI_ERR_ARG;
619   } else {
620     int count = comm->size();
621
622     int rank                           = simgrid::s4u::this_actor::get_pid();
623     int dt_send_size                   = datatype->is_replayable() ? 1 : datatype->size();
624     std::vector<int>* trace_recvcounts = new std::vector<int>(recvcount * dt_send_size); // copy data to avoid bad free
625
626     void* sendtmpbuf       = sendbuf;
627     if (sendbuf == MPI_IN_PLACE) {
628       sendtmpbuf = static_cast<void*>(xbt_malloc(recvcount * count * datatype->size()));
629       memcpy(sendtmpbuf, recvbuf, recvcount * count * datatype->size());
630     }
631
632     TRACE_smpi_comm_in(rank, __func__,
633                        new simgrid::instr::VarCollTIData("reducescatter", -1, 0, nullptr, -1, trace_recvcounts,
634                                                          simgrid::smpi::Datatype::encode(datatype), ""));
635
636     int* recvcounts = new int[count];
637     for (int i      = 0; i < count; i++)
638       recvcounts[i] = recvcount;
639     simgrid::smpi::Colls::reduce_scatter(sendtmpbuf, recvbuf, recvcounts, datatype, op, comm);
640     delete[] recvcounts;
641     retval = MPI_SUCCESS;
642
643     TRACE_smpi_comm_out(rank);
644
645     if (sendbuf == MPI_IN_PLACE)
646       xbt_free(sendtmpbuf);
647   }
648
649   smpi_bench_begin();
650   return retval;
651 }
652 int PMPI_Alltoall(void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
653                   MPI_Datatype recvtype, MPI_Comm comm){
654   return PMPI_Ialltoall(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm, MPI_REQUEST_IGNORED);
655 }
656                   
657 int PMPI_Ialltoall(void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
658                   MPI_Datatype recvtype, MPI_Comm comm, MPI_Request *request)
659 {
660   int retval = 0;
661   smpi_bench_end();
662
663   if (comm == MPI_COMM_NULL) {
664     retval = MPI_ERR_COMM;
665   } else if ((sendbuf != MPI_IN_PLACE && sendtype == MPI_DATATYPE_NULL) || recvtype == MPI_DATATYPE_NULL) {
666     retval = MPI_ERR_TYPE;
667   } else if (request == nullptr){
668     retval = MPI_ERR_ARG;
669   } else {
670     int rank                 = simgrid::s4u::this_actor::get_pid();
671     void* sendtmpbuf         = static_cast<char*>(sendbuf);
672     int sendtmpcount         = sendcount;
673     MPI_Datatype sendtmptype = sendtype;
674     if (sendbuf == MPI_IN_PLACE) {
675       sendtmpbuf = static_cast<void*>(xbt_malloc(recvcount * comm->size() * recvtype->size()));
676       memcpy(sendtmpbuf, recvbuf, recvcount * comm->size() * recvtype->size());
677       sendtmpcount = recvcount;
678       sendtmptype  = recvtype;
679     }
680
681     TRACE_smpi_comm_in(rank, request==MPI_REQUEST_IGNORED?"PMPI_Alltoall":"PMPI_Ialltoall",
682                        new simgrid::instr::CollTIData(
683                            request==MPI_REQUEST_IGNORED ? "alltoall" : "ialltoall", -1, -1.0,
684                            sendtmptype->is_replayable() ? sendtmpcount : sendtmpcount * sendtmptype->size(),
685                            recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(),
686                            simgrid::smpi::Datatype::encode(sendtmptype), simgrid::smpi::Datatype::encode(recvtype)));
687     if(request == MPI_REQUEST_IGNORED)
688       retval = simgrid::smpi::Colls::alltoall(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcount, recvtype, comm);
689     else
690       retval = simgrid::smpi::Colls::ialltoall(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcount, recvtype, comm, request);
691
692     TRACE_smpi_comm_out(rank);
693
694     if (sendbuf == MPI_IN_PLACE)
695       xbt_free(sendtmpbuf);
696   }
697
698   smpi_bench_begin();
699   return retval;
700 }
701
702 int PMPI_Alltoallv(void* sendbuf, int* sendcounts, int* senddisps, MPI_Datatype sendtype, void* recvbuf,
703                    int* recvcounts, int* recvdisps, MPI_Datatype recvtype, MPI_Comm comm)
704 {
705   return PMPI_Ialltoallv(sendbuf, sendcounts, senddisps, sendtype, recvbuf, recvcounts, recvdisps, recvtype, comm, MPI_REQUEST_IGNORED);
706 }
707
708 int PMPI_Ialltoallv(void* sendbuf, int* sendcounts, int* senddisps, MPI_Datatype sendtype, void* recvbuf,
709                    int* recvcounts, int* recvdisps, MPI_Datatype recvtype, MPI_Comm comm, MPI_Request *request)
710 {
711   int retval = 0;
712
713   smpi_bench_end();
714
715   if (comm == MPI_COMM_NULL) {
716     retval = MPI_ERR_COMM;
717   } else if ((sendbuf != MPI_IN_PLACE && sendtype == MPI_DATATYPE_NULL)  || recvtype == MPI_DATATYPE_NULL) {
718     retval = MPI_ERR_TYPE;
719   } else if ((sendbuf != MPI_IN_PLACE && (sendcounts == nullptr || senddisps == nullptr)) || recvcounts == nullptr ||
720              recvdisps == nullptr) {
721     retval = MPI_ERR_ARG;
722   } else if (request == nullptr){
723     retval = MPI_ERR_ARG;
724   }  else {
725     int rank                           = simgrid::s4u::this_actor::get_pid();
726     int size               = comm->size();
727     int send_size                      = 0;
728     int recv_size                      = 0;
729     std::vector<int>* trace_sendcounts = new std::vector<int>;
730     std::vector<int>* trace_recvcounts = new std::vector<int>;
731     int dt_size_recv       = recvtype->size();
732
733     void* sendtmpbuf         = static_cast<char*>(sendbuf);
734     int* sendtmpcounts       = sendcounts;
735     int* sendtmpdisps        = senddisps;
736     MPI_Datatype sendtmptype = sendtype;
737     int maxsize              = 0;
738     for (int i = 0; i < size; i++) { // copy data to avoid bad free
739       recv_size += recvcounts[i] * dt_size_recv;
740       trace_recvcounts->push_back(recvcounts[i] * dt_size_recv);
741       if (((recvdisps[i] + recvcounts[i]) * dt_size_recv) > maxsize)
742         maxsize = (recvdisps[i] + recvcounts[i]) * dt_size_recv;
743     }
744
745     if (sendbuf == MPI_IN_PLACE) {
746       sendtmpbuf = static_cast<void*>(xbt_malloc(maxsize));
747       memcpy(sendtmpbuf, recvbuf, maxsize);
748       sendtmpcounts = static_cast<int*>(xbt_malloc(size * sizeof(int)));
749       memcpy(sendtmpcounts, recvcounts, size * sizeof(int));
750       sendtmpdisps = static_cast<int*>(xbt_malloc(size * sizeof(int)));
751       memcpy(sendtmpdisps, recvdisps, size * sizeof(int));
752       sendtmptype = recvtype;
753     }
754
755     int dt_size_send = sendtmptype->size();
756
757     for (int i = 0; i < size; i++) { // copy data to avoid bad free
758       send_size += sendtmpcounts[i] * dt_size_send;
759       trace_sendcounts->push_back(sendtmpcounts[i] * dt_size_send);
760     }
761
762     TRACE_smpi_comm_in(rank, request==MPI_REQUEST_IGNORED?"PMPI_Alltoallv":"PMPI_Ialltoallv",
763                        new simgrid::instr::VarCollTIData(request==MPI_REQUEST_IGNORED ? "alltoallv":"ialltoallv", -1, send_size, trace_sendcounts, recv_size,
764                                                          trace_recvcounts, simgrid::smpi::Datatype::encode(sendtype),
765                                                          simgrid::smpi::Datatype::encode(recvtype)));
766
767     if(request == MPI_REQUEST_IGNORED)
768       retval = simgrid::smpi::Colls::alltoallv(sendtmpbuf, sendtmpcounts, sendtmpdisps, sendtmptype, recvbuf, recvcounts,
769                                     recvdisps, recvtype, comm);
770     else
771       retval = simgrid::smpi::Colls::ialltoallv(sendtmpbuf, sendtmpcounts, sendtmpdisps, sendtmptype, recvbuf, recvcounts,
772                                     recvdisps, recvtype, comm, request);
773     TRACE_smpi_comm_out(rank);
774
775     if (sendbuf == MPI_IN_PLACE) {
776       xbt_free(sendtmpbuf);
777       xbt_free(sendtmpcounts);
778       xbt_free(sendtmpdisps);
779     }
780   }
781
782   smpi_bench_begin();
783   return retval;
784 }
785
786 int PMPI_Alltoallw(void* sendbuf, int* sendcounts, int* senddisps, MPI_Datatype* sendtypes, void* recvbuf,
787                    int* recvcounts, int* recvdisps, MPI_Datatype* recvtypes, MPI_Comm comm)
788 {
789   return PMPI_Ialltoallw(sendbuf, sendcounts, senddisps, sendtypes, recvbuf, recvcounts, recvdisps, recvtypes, comm, MPI_REQUEST_IGNORED);
790 }
791
792 int PMPI_Ialltoallw(void* sendbuf, int* sendcounts, int* senddisps, MPI_Datatype* sendtypes, void* recvbuf,
793                    int* recvcounts, int* recvdisps, MPI_Datatype* recvtypes, MPI_Comm comm, MPI_Request *request)
794 {
795   int retval = 0;
796
797   smpi_bench_end();
798
799   if (comm == MPI_COMM_NULL) {
800     retval = MPI_ERR_COMM;
801   } else if ((sendbuf != MPI_IN_PLACE && sendtypes == nullptr)  || recvtypes == nullptr) {
802     retval = MPI_ERR_TYPE;
803   } else if ((sendbuf != MPI_IN_PLACE && (sendcounts == nullptr || senddisps == nullptr)) || recvcounts == nullptr ||
804              recvdisps == nullptr) {
805     retval = MPI_ERR_ARG;
806   } else if (request == nullptr){
807     retval = MPI_ERR_ARG;
808   }  else {
809     int rank                           = simgrid::s4u::this_actor::get_pid();
810     int size                           = comm->size();
811     int send_size                      = 0;
812     int recv_size                      = 0;
813     std::vector<int>* trace_sendcounts = new std::vector<int>;
814     std::vector<int>* trace_recvcounts = new std::vector<int>;
815
816     void* sendtmpbuf           = static_cast<char*>(sendbuf);
817     int* sendtmpcounts         = sendcounts;
818     int* sendtmpdisps          = senddisps;
819     MPI_Datatype* sendtmptypes = sendtypes;
820     unsigned long maxsize                = 0;
821     for (int i = 0; i < size; i++) { // copy data to avoid bad free
822       if(recvtypes[i]==MPI_DATATYPE_NULL){
823         delete trace_recvcounts;
824         delete trace_sendcounts;
825         return MPI_ERR_TYPE;
826       }
827       recv_size += recvcounts[i] * recvtypes[i]->size();
828       trace_recvcounts->push_back(recvcounts[i] * recvtypes[i]->size());
829       if ((recvdisps[i] + (recvcounts[i] * recvtypes[i]->size())) > maxsize)
830         maxsize = recvdisps[i] + (recvcounts[i] * recvtypes[i]->size());
831     }
832
833     if (sendbuf == MPI_IN_PLACE) {
834       sendtmpbuf = static_cast<void*>(xbt_malloc(maxsize));
835       memcpy(sendtmpbuf, recvbuf, maxsize);
836       sendtmpcounts = static_cast<int*>(xbt_malloc(size * sizeof(int)));
837       memcpy(sendtmpcounts, recvcounts, size * sizeof(int));
838       sendtmpdisps = static_cast<int*>(xbt_malloc(size * sizeof(int)));
839       memcpy(sendtmpdisps, recvdisps, size * sizeof(int));
840       sendtmptypes = static_cast<MPI_Datatype*>(xbt_malloc(size * sizeof(MPI_Datatype)));
841       memcpy(sendtmptypes, recvtypes, size * sizeof(MPI_Datatype));
842     }
843
844     for (int i = 0; i < size; i++) { // copy data to avoid bad free
845       send_size += sendtmpcounts[i] * sendtmptypes[i]->size();
846       trace_sendcounts->push_back(sendtmpcounts[i] * sendtmptypes[i]->size());
847     }
848
849     TRACE_smpi_comm_in(rank, request==MPI_REQUEST_IGNORED?"PMPI_Alltoallw":"PMPI_Ialltoallw",
850                        new simgrid::instr::VarCollTIData(request==MPI_REQUEST_IGNORED ? "alltoallv":"ialltoallv", -1, send_size, trace_sendcounts, recv_size,
851                                                          trace_recvcounts, simgrid::smpi::Datatype::encode(sendtmptypes[0]),
852                                                          simgrid::smpi::Datatype::encode(recvtypes[0])));
853
854     if(request == MPI_REQUEST_IGNORED)
855       retval = simgrid::smpi::Colls::alltoallw(sendtmpbuf, sendtmpcounts, sendtmpdisps, sendtmptypes, recvbuf, recvcounts,
856                                     recvdisps, recvtypes, comm);
857     else
858       retval = simgrid::smpi::Colls::ialltoallw(sendtmpbuf, sendtmpcounts, sendtmpdisps, sendtmptypes, recvbuf, recvcounts,
859                                     recvdisps, recvtypes, comm, request);
860     TRACE_smpi_comm_out(rank);
861
862     if (sendbuf == MPI_IN_PLACE) {
863       xbt_free(sendtmpbuf);
864       xbt_free(sendtmpcounts);
865       xbt_free(sendtmpdisps);
866       xbt_free(sendtmptypes);
867     }
868   }
869
870   smpi_bench_begin();
871   return retval;
872 }