Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
62e4822d5785a3c8b3edd7313997459123403a71
[simgrid.git] / src / smpi / bindings / smpi_pmpi_request.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_comm.hpp"
8 #include "smpi_datatype.hpp"
9 #include "smpi_request.hpp"
10 #include "src/smpi/include/smpi_actor.hpp"
11
12 XBT_LOG_EXTERNAL_DEFAULT_CATEGORY(smpi_pmpi);
13
14 static int getPid(MPI_Comm, int);
15 static int getPid(MPI_Comm comm, int id)
16 {
17   simgrid::s4u::ActorPtr actor = comm->group()->actor(id);
18   return (actor == nullptr) ? MPI_UNDEFINED : actor->get_pid();
19 }
20
21 /* PMPI User level calls */
22
23 int PMPI_Send_init(const void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request * request)
24 {
25   int retval = 0;
26
27   smpi_bench_end();
28   if (request == nullptr) {
29     retval = MPI_ERR_ARG;
30   } else if (comm == MPI_COMM_NULL) {
31     retval = MPI_ERR_COMM;
32   } else if (datatype==MPI_DATATYPE_NULL || not datatype->is_valid()) {
33     retval = MPI_ERR_TYPE;
34   } else if (dst == MPI_PROC_NULL) {
35     retval = MPI_SUCCESS;
36   } else {
37     *request = simgrid::smpi::Request::send_init(buf, count, datatype, dst, tag, comm);
38     retval   = MPI_SUCCESS;
39   }
40   smpi_bench_begin();
41   if (retval != MPI_SUCCESS && request != nullptr)
42     *request = MPI_REQUEST_NULL;
43   return retval;
44 }
45
46 int PMPI_Recv_init(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Request * request)
47 {
48   int retval = 0;
49
50   smpi_bench_end();
51   if (request == nullptr) {
52     retval = MPI_ERR_ARG;
53   } else if (comm == MPI_COMM_NULL) {
54     retval = MPI_ERR_COMM;
55   } else if (datatype==MPI_DATATYPE_NULL || not datatype->is_valid()) {
56     retval = MPI_ERR_TYPE;
57   } else if (src == MPI_PROC_NULL) {
58     retval = MPI_SUCCESS;
59   } else {
60     *request = simgrid::smpi::Request::recv_init(buf, count, datatype, src, tag, comm);
61     retval = MPI_SUCCESS;
62   }
63   smpi_bench_begin();
64   if (retval != MPI_SUCCESS && request != nullptr)
65     *request = MPI_REQUEST_NULL;
66   return retval;
67 }
68
69 int PMPI_Rsend_init(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm,
70                     MPI_Request* request)
71 {
72   return PMPI_Send_init(buf, count, datatype, dst, tag, comm, request);
73 }
74
75 int PMPI_Ssend_init(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request)
76 {
77   int retval = 0;
78
79   smpi_bench_end();
80   if (request == nullptr) {
81     retval = MPI_ERR_ARG;
82   } else if (comm == MPI_COMM_NULL) {
83     retval = MPI_ERR_COMM;
84   } else if (datatype==MPI_DATATYPE_NULL || not datatype->is_valid()) {
85     retval = MPI_ERR_TYPE;
86   } else if (dst == MPI_PROC_NULL) {
87     retval = MPI_SUCCESS;
88   } else {
89     *request = simgrid::smpi::Request::ssend_init(buf, count, datatype, dst, tag, comm);
90     retval = MPI_SUCCESS;
91   }
92   smpi_bench_begin();
93   if (retval != MPI_SUCCESS && request != nullptr)
94     *request = MPI_REQUEST_NULL;
95   return retval;
96 }
97
98 /*
99  * This function starts a request returned by init functions such as
100  * MPI_Send_init(), MPI_Ssend_init (see above), and friends.
101  * They should already have performed sanity checks.
102  */
103 int PMPI_Start(MPI_Request * request)
104 {
105   int retval = 0;
106
107   smpi_bench_end();
108   if (request == nullptr || *request == MPI_REQUEST_NULL) {
109     retval = MPI_ERR_REQUEST;
110   } else {
111     MPI_Request req = *request;
112     int my_proc_id = (req->comm() != MPI_COMM_NULL) ? simgrid::s4u::this_actor::get_pid() : -1;
113     TRACE_smpi_comm_in(my_proc_id, __func__,
114                        new simgrid::instr::Pt2PtTIData("Start", req->dst(),
115                                                        req->size(),
116                                                        req->tag(), 
117                                                        simgrid::smpi::Datatype::encode(req->type())));
118     if (not TRACE_smpi_view_internals() && req->flags() & MPI_REQ_SEND)
119       TRACE_smpi_send(my_proc_id, my_proc_id, getPid(req->comm(), req->dst()), req->tag(), req->size());
120
121     req->start();
122
123     if (not TRACE_smpi_view_internals() && req->flags() & MPI_REQ_RECV)
124       TRACE_smpi_recv(getPid(req->comm(), req->src()), my_proc_id, req->tag());
125     retval = MPI_SUCCESS;
126     TRACE_smpi_comm_out(my_proc_id);
127   }
128   smpi_bench_begin();
129   return retval;
130 }
131
132 int PMPI_Startall(int count, MPI_Request * requests)
133 {
134   int retval;
135   smpi_bench_end();
136   if (requests == nullptr) {
137     retval = MPI_ERR_ARG;
138   } else {
139     retval = MPI_SUCCESS;
140     for (int i = 0; i < count; i++) {
141       if(requests[i] == MPI_REQUEST_NULL) {
142         retval = MPI_ERR_REQUEST;
143       }
144     }
145     if(retval != MPI_ERR_REQUEST) {
146       int my_proc_id = simgrid::s4u::this_actor::get_pid();
147       TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("Startall"));
148       if (not TRACE_smpi_view_internals())
149         for (int i = 0; i < count; i++) {
150           MPI_Request req = requests[i];
151           if (req->flags() & MPI_REQ_SEND)
152             TRACE_smpi_send(my_proc_id, my_proc_id, getPid(req->comm(), req->dst()), req->tag(), req->size());
153         }
154
155       simgrid::smpi::Request::startall(count, requests);
156
157       if (not TRACE_smpi_view_internals())
158         for (int i = 0; i < count; i++) {
159           MPI_Request req = requests[i];
160           if (req->flags() & MPI_REQ_RECV)
161             TRACE_smpi_recv(getPid(req->comm(), req->src()), my_proc_id, req->tag());
162         }
163       TRACE_smpi_comm_out(my_proc_id);
164     }
165   }
166   smpi_bench_begin();
167   return retval;
168 }
169
170 int PMPI_Request_free(MPI_Request * request)
171 {
172   int retval = 0;
173
174   smpi_bench_end();
175   if (*request != MPI_REQUEST_NULL) {
176     simgrid::smpi::Request::unref(request);
177     *request = MPI_REQUEST_NULL;
178     retval = MPI_SUCCESS;
179   }
180   smpi_bench_begin();
181   return retval;
182 }
183
184 int PMPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Request * request)
185 {
186   int retval = 0;
187
188   smpi_bench_end();
189
190   if (request == nullptr) {
191     retval = MPI_ERR_ARG;
192   } else if (comm == MPI_COMM_NULL) {
193     retval = MPI_ERR_COMM;
194   } else if (src == MPI_PROC_NULL) {
195     *request = MPI_REQUEST_NULL;
196     retval = MPI_SUCCESS;
197   } else if (src!=MPI_ANY_SOURCE && (src >= comm->group()->size() || src <0)){
198     retval = MPI_ERR_RANK;
199   } else if ((count < 0) || (buf==nullptr && count > 0)) {
200     retval = MPI_ERR_COUNT;
201   } else if (datatype==MPI_DATATYPE_NULL || not datatype->is_valid()) {
202     retval = MPI_ERR_TYPE;
203   } else if(tag<0 && tag !=  MPI_ANY_TAG){
204     retval = MPI_ERR_TAG;
205   } else {
206     int my_proc_id = simgrid::s4u::this_actor::get_pid();
207
208     TRACE_smpi_comm_in(my_proc_id, __func__,
209                        new simgrid::instr::Pt2PtTIData("irecv", src,
210                                                        datatype->is_replayable() ? count : count * datatype->size(),
211                                                        tag, simgrid::smpi::Datatype::encode(datatype)));
212
213     *request = simgrid::smpi::Request::irecv(buf, count, datatype, src, tag, comm);
214     retval = MPI_SUCCESS;
215
216     TRACE_smpi_comm_out(my_proc_id);
217   }
218
219   smpi_bench_begin();
220   if (retval != MPI_SUCCESS && request != nullptr)
221     *request = MPI_REQUEST_NULL;
222   return retval;
223 }
224
225
226 int PMPI_Isend(const void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request * request)
227 {
228   int retval = 0;
229
230   smpi_bench_end();
231   if (request == nullptr) {
232     retval = MPI_ERR_ARG;
233   } else if (comm == MPI_COMM_NULL) {
234     retval = MPI_ERR_COMM;
235   } else if (dst == MPI_PROC_NULL) {
236     *request = MPI_REQUEST_NULL;
237     retval = MPI_SUCCESS;
238   } else if (dst >= comm->group()->size() || dst <0){
239     retval = MPI_ERR_RANK;
240   } else if ((count < 0) || (buf==nullptr && count > 0)) {
241     retval = MPI_ERR_COUNT;
242   } else if (datatype==MPI_DATATYPE_NULL || not datatype->is_valid()) {
243     retval = MPI_ERR_TYPE;
244   } else if(tag<0 && tag !=  MPI_ANY_TAG){
245     retval = MPI_ERR_TAG;
246   } else {
247     int my_proc_id = simgrid::s4u::this_actor::get_pid();
248     int trace_dst = getPid(comm, dst);
249     TRACE_smpi_comm_in(my_proc_id, __func__,
250                        new simgrid::instr::Pt2PtTIData("isend", dst,
251                                                        datatype->is_replayable() ? count : count * datatype->size(),
252                                                        tag, simgrid::smpi::Datatype::encode(datatype)));
253
254     TRACE_smpi_send(my_proc_id, my_proc_id, trace_dst, tag, count * datatype->size());
255
256     *request = simgrid::smpi::Request::isend(buf, count, datatype, dst, tag, comm);
257     retval = MPI_SUCCESS;
258
259     TRACE_smpi_comm_out(my_proc_id);
260   }
261
262   smpi_bench_begin();
263   if (retval != MPI_SUCCESS && request!=nullptr)
264     *request = MPI_REQUEST_NULL;
265   return retval;
266 }
267
268 int PMPI_Irsend(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm,
269                 MPI_Request* request)
270 {
271   return PMPI_Isend(buf, count, datatype, dst, tag, comm, request);
272 }
273
274 int PMPI_Issend(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request)
275 {
276   int retval = 0;
277
278   smpi_bench_end();
279   if (request == nullptr) {
280     retval = MPI_ERR_ARG;
281   } else if (comm == MPI_COMM_NULL) {
282     retval = MPI_ERR_COMM;
283   } else if (dst == MPI_PROC_NULL) {
284     *request = MPI_REQUEST_NULL;
285     retval = MPI_SUCCESS;
286   } else if (dst >= comm->group()->size() || dst <0){
287     retval = MPI_ERR_RANK;
288   } else if ((count < 0)|| (buf==nullptr && count > 0)) {
289     retval = MPI_ERR_COUNT;
290   } else if (datatype==MPI_DATATYPE_NULL || not datatype->is_valid()) {
291     retval = MPI_ERR_TYPE;
292   } else if(tag<0 && tag !=  MPI_ANY_TAG){
293     retval = MPI_ERR_TAG;
294   } else {
295     int my_proc_id = simgrid::s4u::this_actor::get_pid();
296     int trace_dst = getPid(comm, dst);
297     TRACE_smpi_comm_in(my_proc_id, __func__,
298                        new simgrid::instr::Pt2PtTIData("ISsend", dst,
299                                                        datatype->is_replayable() ? count : count * datatype->size(),
300                                                        tag, simgrid::smpi::Datatype::encode(datatype)));
301     TRACE_smpi_send(my_proc_id, my_proc_id, trace_dst, tag, count * datatype->size());
302
303     *request = simgrid::smpi::Request::issend(buf, count, datatype, dst, tag, comm);
304     retval = MPI_SUCCESS;
305
306     TRACE_smpi_comm_out(my_proc_id);
307   }
308
309   smpi_bench_begin();
310   if (retval != MPI_SUCCESS && request!=nullptr)
311     *request = MPI_REQUEST_NULL;
312   return retval;
313 }
314
315 int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Status * status)
316 {
317   int retval = 0;
318
319   smpi_bench_end();
320   if (comm == MPI_COMM_NULL) {
321     retval = MPI_ERR_COMM;
322   } else if (src == MPI_PROC_NULL) {
323     if(status != MPI_STATUS_IGNORE){
324       simgrid::smpi::Status::empty(status);
325       status->MPI_SOURCE = MPI_PROC_NULL;
326     }
327     retval = MPI_SUCCESS;
328   } else if (src!=MPI_ANY_SOURCE && (src >= comm->group()->size() || src <0)){
329     retval = MPI_ERR_RANK;
330   } else if ((count < 0) || (buf==nullptr && count > 0)) {
331     retval = MPI_ERR_COUNT;
332   } else if (datatype==MPI_DATATYPE_NULL || not datatype->is_valid()) {
333     retval = MPI_ERR_TYPE;
334   } else if(tag<0 && tag !=  MPI_ANY_TAG){
335     retval = MPI_ERR_TAG;
336   } else {
337     int my_proc_id = simgrid::s4u::this_actor::get_pid();
338     TRACE_smpi_comm_in(my_proc_id, __func__,
339                        new simgrid::instr::Pt2PtTIData("recv", src,
340                                                        datatype->is_replayable() ? count : count * datatype->size(),
341                                                        tag, simgrid::smpi::Datatype::encode(datatype)));
342
343     simgrid::smpi::Request::recv(buf, count, datatype, src, tag, comm, status);
344     retval = MPI_SUCCESS;
345
346     // the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
347     int src_traced=0;
348     if (status != MPI_STATUS_IGNORE) 
349       src_traced = getPid(comm, status->MPI_SOURCE);
350     else
351       src_traced = getPid(comm, src);
352     if (not TRACE_smpi_view_internals()) {
353       TRACE_smpi_recv(src_traced, my_proc_id, tag);
354     }
355     
356     TRACE_smpi_comm_out(my_proc_id);
357   }
358
359   smpi_bench_begin();
360   return retval;
361 }
362
363 int PMPI_Send(const void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm)
364 {
365   int retval = 0;
366
367   smpi_bench_end();
368
369   if (comm == MPI_COMM_NULL) {
370     retval = MPI_ERR_COMM;
371   } else if (dst == MPI_PROC_NULL) {
372     retval = MPI_SUCCESS;
373   } else if (dst >= comm->group()->size() || dst <0){
374     retval = MPI_ERR_RANK;
375   } else if ((count < 0) || (buf == nullptr && count > 0)) {
376     retval = MPI_ERR_COUNT;
377   } else if (datatype==MPI_DATATYPE_NULL || not datatype->is_valid()) {
378     retval = MPI_ERR_TYPE;
379   } else if(tag < 0 && tag !=  MPI_ANY_TAG){
380     retval = MPI_ERR_TAG;
381   } else {
382     int my_proc_id         = simgrid::s4u::this_actor::get_pid();
383     int dst_traced         = getPid(comm, dst);
384     TRACE_smpi_comm_in(my_proc_id, __func__,
385                        new simgrid::instr::Pt2PtTIData("send", dst,
386                                                        datatype->is_replayable() ? count : count * datatype->size(),
387                                                        tag, simgrid::smpi::Datatype::encode(datatype)));
388     if (not TRACE_smpi_view_internals()) {
389       TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, tag, count * datatype->size());
390     }
391
392     simgrid::smpi::Request::send(buf, count, datatype, dst, tag, comm);
393     retval = MPI_SUCCESS;
394
395     TRACE_smpi_comm_out(my_proc_id);
396   }
397
398   smpi_bench_begin();
399   return retval;
400 }
401
402 int PMPI_Rsend(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm)
403 {
404   return PMPI_Send(buf, count, datatype, dst, tag, comm);
405 }
406
407 int PMPI_Bsend(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm)
408 {
409   int retval = 0;
410
411   smpi_bench_end();
412
413   if (comm == MPI_COMM_NULL) {
414     retval = MPI_ERR_COMM;
415   } else if (dst == MPI_PROC_NULL) {
416     retval = MPI_SUCCESS;
417   } else if (dst >= comm->group()->size() || dst <0){
418     retval = MPI_ERR_RANK;
419   } else if ((count < 0) || (buf == nullptr && count > 0)) {
420     retval = MPI_ERR_COUNT;
421   } else if (datatype==MPI_DATATYPE_NULL || not datatype->is_valid()) {
422     retval = MPI_ERR_TYPE;
423   } else if(tag < 0 && tag !=  MPI_ANY_TAG){
424     retval = MPI_ERR_TAG;
425   } else {
426     int my_proc_id         = simgrid::s4u::this_actor::get_pid();
427     int dst_traced         = getPid(comm, dst);
428     int bsend_buf_size = 0;
429     void* bsend_buf = nullptr;
430     smpi_process()->bsend_buffer(&bsend_buf, &bsend_buf_size);
431     int size = datatype->get_extent() * count;
432     if(bsend_buf==nullptr || bsend_buf_size < size + MPI_BSEND_OVERHEAD )
433       return MPI_ERR_BUFFER;
434     TRACE_smpi_comm_in(my_proc_id, __func__,
435                        new simgrid::instr::Pt2PtTIData("bsend", dst,
436                                                        datatype->is_replayable() ? count : count * datatype->size(),
437                                                        tag, simgrid::smpi::Datatype::encode(datatype)));
438     if (not TRACE_smpi_view_internals()) {
439       TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, tag, count * datatype->size());
440     }
441
442     simgrid::smpi::Request::bsend(buf, count, datatype, dst, tag, comm);
443     retval = MPI_SUCCESS;
444
445     TRACE_smpi_comm_out(my_proc_id);
446   }
447
448   smpi_bench_begin();
449   return retval;
450 }
451
452 int PMPI_Ibsend(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request)
453 {
454   int retval = 0;
455
456   smpi_bench_end();
457   if (request == nullptr) {
458     retval = MPI_ERR_ARG;
459   } else if (comm == MPI_COMM_NULL) {
460     retval = MPI_ERR_COMM;
461   } else if (dst == MPI_PROC_NULL) {
462     *request = MPI_REQUEST_NULL;
463     retval = MPI_SUCCESS;
464   } else if (dst >= comm->group()->size() || dst <0){
465     retval = MPI_ERR_RANK;
466   } else if ((count < 0) || (buf==nullptr && count > 0)) {
467     retval = MPI_ERR_COUNT;
468   } else if (datatype==MPI_DATATYPE_NULL || not datatype->is_valid()) {
469     retval = MPI_ERR_TYPE;
470   } else if(tag<0 && tag !=  MPI_ANY_TAG){
471     retval = MPI_ERR_TAG;
472   } else {
473     int my_proc_id = simgrid::s4u::this_actor::get_pid();
474     int trace_dst = getPid(comm, dst);
475     int bsend_buf_size = 0;
476     void* bsend_buf = nullptr;
477     smpi_process()->bsend_buffer(&bsend_buf, &bsend_buf_size);
478     int size = datatype->get_extent() * count;
479     if(bsend_buf==nullptr || bsend_buf_size < size + MPI_BSEND_OVERHEAD )
480       return MPI_ERR_BUFFER;
481     TRACE_smpi_comm_in(my_proc_id, __func__,
482                        new simgrid::instr::Pt2PtTIData("ibsend", dst,
483                                                        datatype->is_replayable() ? count : count * datatype->size(),
484                                                        tag, simgrid::smpi::Datatype::encode(datatype)));
485
486     TRACE_smpi_send(my_proc_id, my_proc_id, trace_dst, tag, count * datatype->size());
487
488     *request = simgrid::smpi::Request::ibsend(buf, count, datatype, dst, tag, comm);
489     retval = MPI_SUCCESS;
490
491     TRACE_smpi_comm_out(my_proc_id);
492   }
493
494   smpi_bench_begin();
495   if (retval != MPI_SUCCESS && request!=nullptr)
496     *request = MPI_REQUEST_NULL;
497   return retval;
498 }
499
500 int PMPI_Bsend_init(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request)
501 {
502   int retval = 0;
503
504   smpi_bench_end();
505   if (request == nullptr) {
506     retval = MPI_ERR_ARG;
507   } else if (comm == MPI_COMM_NULL) {
508     retval = MPI_ERR_COMM;
509   } else if (datatype==MPI_DATATYPE_NULL || not datatype->is_valid()) {
510     retval = MPI_ERR_TYPE;
511   } else if (dst == MPI_PROC_NULL) {
512     retval = MPI_SUCCESS;
513   } else {
514     int bsend_buf_size = 0;
515     void* bsend_buf = nullptr;
516     smpi_process()->bsend_buffer(&bsend_buf, &bsend_buf_size);
517     if( bsend_buf==nullptr || bsend_buf_size < datatype->get_extent() * count + MPI_BSEND_OVERHEAD ) {
518       retval = MPI_ERR_BUFFER;
519     } else {
520       *request = simgrid::smpi::Request::bsend_init(buf, count, datatype, dst, tag, comm);
521       retval   = MPI_SUCCESS;
522     }
523   }
524   smpi_bench_begin();
525   if (retval != MPI_SUCCESS && request != nullptr)
526     *request = MPI_REQUEST_NULL;
527   return retval;
528 }
529
530 int PMPI_Ssend(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm) {
531   int retval = 0;
532
533   smpi_bench_end();
534
535   if (comm == MPI_COMM_NULL) {
536     retval = MPI_ERR_COMM;
537   } else if (dst == MPI_PROC_NULL) {
538     retval = MPI_SUCCESS;
539   } else if (dst >= comm->group()->size() || dst <0){
540     retval = MPI_ERR_RANK;
541   } else if ((count < 0) || (buf==nullptr && count > 0)) {
542     retval = MPI_ERR_COUNT;
543   } else if (datatype==MPI_DATATYPE_NULL || not datatype->is_valid()) {
544     retval = MPI_ERR_TYPE;
545   } else if(tag<0 && tag !=  MPI_ANY_TAG){
546     retval = MPI_ERR_TAG;
547   } else {
548     int my_proc_id         = simgrid::s4u::this_actor::get_pid();
549     int dst_traced         = getPid(comm, dst);
550     TRACE_smpi_comm_in(my_proc_id, __func__,
551                        new simgrid::instr::Pt2PtTIData("Ssend", dst,
552                                                        datatype->is_replayable() ? count : count * datatype->size(),
553                                                        tag, simgrid::smpi::Datatype::encode(datatype)));
554     TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, tag, count * datatype->size());
555
556     simgrid::smpi::Request::ssend(buf, count, datatype, dst, tag, comm);
557     retval = MPI_SUCCESS;
558
559     TRACE_smpi_comm_out(my_proc_id);
560   }
561
562   smpi_bench_begin();
563   return retval;
564 }
565
566 int PMPI_Sendrecv(const void* sendbuf, int sendcount, MPI_Datatype sendtype, int dst, int sendtag, void* recvbuf,
567                   int recvcount, MPI_Datatype recvtype, int src, int recvtag, MPI_Comm comm, MPI_Status* status)
568 {
569   int retval = 0;
570
571   smpi_bench_end();
572
573   if (comm == MPI_COMM_NULL) {
574     retval = MPI_ERR_COMM;
575   } else if (not sendtype->is_valid() || not recvtype->is_valid()) {
576     retval = MPI_ERR_TYPE;
577   } else if (src == MPI_PROC_NULL) {
578     if(status!=MPI_STATUS_IGNORE){
579       simgrid::smpi::Status::empty(status);
580       status->MPI_SOURCE = MPI_PROC_NULL;
581     }
582     if(dst != MPI_PROC_NULL)
583       simgrid::smpi::Request::send(sendbuf, sendcount, sendtype, dst, sendtag, comm);
584     retval = MPI_SUCCESS;
585   }else if (dst == MPI_PROC_NULL){
586     simgrid::smpi::Request::recv(recvbuf, recvcount, recvtype, src, recvtag, comm, status);
587     retval = MPI_SUCCESS;
588   }else if (dst >= comm->group()->size() || dst <0 ||
589       (src!=MPI_ANY_SOURCE && (src >= comm->group()->size() || src <0))){
590     retval = MPI_ERR_RANK;
591   } else if ((sendcount < 0 || recvcount<0) ||
592       (sendbuf==nullptr && sendcount > 0) || (recvbuf==nullptr && recvcount>0)) {
593     retval = MPI_ERR_COUNT;
594   } else if((sendtag<0 && sendtag !=  MPI_ANY_TAG)||(recvtag<0 && recvtag != MPI_ANY_TAG)){
595     retval = MPI_ERR_TAG;
596   } else {
597     int my_proc_id         = simgrid::s4u::this_actor::get_pid();
598     int dst_traced         = getPid(comm, dst);
599     int src_traced         = getPid(comm, src);
600
601     // FIXME: Hack the way to trace this one
602     std::vector<int>* dst_hack = new std::vector<int>;
603     std::vector<int>* src_hack = new std::vector<int>;
604     dst_hack->push_back(dst_traced);
605     src_hack->push_back(src_traced);
606     TRACE_smpi_comm_in(my_proc_id, __func__,
607                        new simgrid::instr::VarCollTIData(
608                            "sendRecv", -1, sendtype->is_replayable() ? sendcount : sendcount * sendtype->size(),
609                            dst_hack, recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(), src_hack,
610                            simgrid::smpi::Datatype::encode(sendtype), simgrid::smpi::Datatype::encode(recvtype)));
611
612     TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, sendtag, sendcount * sendtype->size());
613
614     simgrid::smpi::Request::sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf, recvcount, recvtype, src,
615                                      recvtag, comm, status);
616     retval = MPI_SUCCESS;
617
618     TRACE_smpi_recv(src_traced, my_proc_id, recvtag);
619     TRACE_smpi_comm_out(my_proc_id);
620   }
621
622   smpi_bench_begin();
623   return retval;
624 }
625
626 int PMPI_Sendrecv_replace(void* buf, int count, MPI_Datatype datatype, int dst, int sendtag, int src, int recvtag,
627                           MPI_Comm comm, MPI_Status* status)
628 {
629   int retval = 0;
630   if (datatype==MPI_DATATYPE_NULL || not datatype->is_valid()) {
631     return MPI_ERR_TYPE;
632   } else if (count < 0) {
633     return MPI_ERR_COUNT;
634   } else {
635     int size = datatype->get_extent() * count;
636     void* recvbuf = xbt_new0(char, size);
637     retval = MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count, datatype, src, recvtag, comm, status);
638     if(retval==MPI_SUCCESS){
639       simgrid::smpi::Datatype::copy(recvbuf, count, datatype, buf, count, datatype);
640     }
641     xbt_free(recvbuf);
642   }
643   return retval;
644 }
645
646 int PMPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
647 {
648   int retval = 0;
649   smpi_bench_end();
650   if (request == nullptr || flag == nullptr) {
651     retval = MPI_ERR_ARG;
652   } else if (*request == MPI_REQUEST_NULL) {
653     if (status != MPI_STATUS_IGNORE){
654       *flag= true;
655       simgrid::smpi::Status::empty(status);
656     }
657     retval = MPI_SUCCESS;
658   } else {
659     int my_proc_id = ((*request)->comm() != MPI_COMM_NULL) ? simgrid::s4u::this_actor::get_pid() : -1;
660
661     TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("test"));
662     
663     retval = simgrid::smpi::Request::test(request,status, flag);
664
665     TRACE_smpi_comm_out(my_proc_id);
666   }
667   smpi_bench_begin();
668   return retval;
669 }
670
671 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag, MPI_Status * status)
672 {
673   int retval = 0;
674
675   smpi_bench_end();
676   if (index == nullptr || flag == nullptr) {
677     retval = MPI_ERR_ARG;
678   } else {
679     int my_proc_id = simgrid::s4u::this_actor::get_pid();
680     TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("testany"));
681     retval = simgrid::smpi::Request::testany(count, requests, index, flag, status);
682     TRACE_smpi_comm_out(my_proc_id);
683   }
684   smpi_bench_begin();
685   return retval;
686 }
687
688 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses)
689 {
690   int retval = 0;
691
692   smpi_bench_end();
693   if (flag == nullptr) {
694     retval = MPI_ERR_ARG;
695   } else {
696     int my_proc_id = simgrid::s4u::this_actor::get_pid();
697     TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("testall"));
698     retval = simgrid::smpi::Request::testall(count, requests, flag, statuses);
699     TRACE_smpi_comm_out(my_proc_id);
700   }
701   smpi_bench_begin();
702   return retval;
703 }
704
705 int PMPI_Testsome(int incount, MPI_Request requests[], int* outcount, int* indices, MPI_Status status[])
706 {
707   int retval = 0;
708
709   smpi_bench_end();
710   if (outcount == nullptr) {
711     retval = MPI_ERR_ARG;
712   } else {
713     int my_proc_id = simgrid::s4u::this_actor::get_pid();
714     TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("testsome"));
715     retval = simgrid::smpi::Request::testsome(incount, requests, outcount, indices, status);
716     TRACE_smpi_comm_out(my_proc_id);
717   }
718   smpi_bench_begin();
719   return retval;
720 }
721
722 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
723   int retval = 0;
724   smpi_bench_end();
725
726   if (comm == MPI_COMM_NULL) {
727     retval = MPI_ERR_COMM;
728   } else if (source == MPI_PROC_NULL) {
729     if (status != MPI_STATUS_IGNORE){
730       simgrid::smpi::Status::empty(status);
731       status->MPI_SOURCE = MPI_PROC_NULL;
732     }
733     retval = MPI_SUCCESS;
734   } else {
735     simgrid::smpi::Request::probe(source, tag, comm, status);
736     retval = MPI_SUCCESS;
737   }
738   smpi_bench_begin();
739   return retval;
740 }
741
742 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
743   int retval = 0;
744   smpi_bench_end();
745
746   if (flag == nullptr) {
747     retval = MPI_ERR_ARG;
748   } else if (comm == MPI_COMM_NULL) {
749     retval = MPI_ERR_COMM;
750   } else if (source == MPI_PROC_NULL) {
751     *flag=true;
752     if (status != MPI_STATUS_IGNORE){
753       simgrid::smpi::Status::empty(status);
754       status->MPI_SOURCE = MPI_PROC_NULL;
755     }
756     retval = MPI_SUCCESS;
757   } else {
758     simgrid::smpi::Request::iprobe(source, tag, comm, flag, status);
759     retval = MPI_SUCCESS;
760   }
761   smpi_bench_begin();
762   return retval;
763 }
764
765 // TODO: cheinrich: Move declaration to other file? Rename this function - it's used for PMPI_Wait*?
766 static void trace_smpi_recv_helper(MPI_Request* request, MPI_Status* status);
767 static void trace_smpi_recv_helper(MPI_Request* request, MPI_Status* status)
768 {
769   MPI_Request req = *request;
770   if (req != MPI_REQUEST_NULL) { // Received requests become null
771     int src_traced = req->src();
772     // the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
773     int dst_traced = req->dst();
774     if (req->flags() & MPI_REQ_RECV) { // Is this request a wait for RECV?
775       if (src_traced == MPI_ANY_SOURCE)
776         src_traced = (status != MPI_STATUS_IGNORE) ? req->comm()->group()->rank(status->MPI_SOURCE) : req->src();
777       TRACE_smpi_recv(src_traced, dst_traced, req->tag());
778     }
779   }
780 }
781
782 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
783 {
784   int retval = 0;
785
786   smpi_bench_end();
787
788   simgrid::smpi::Status::empty(status);
789
790   if (request == nullptr) {
791     retval = MPI_ERR_ARG;
792   } else if (*request == MPI_REQUEST_NULL) {
793     retval = MPI_SUCCESS;
794   } else {
795     // for tracing, save the handle which might get overridden before we can use the helper on it
796     MPI_Request savedreq = *request;
797     if (savedreq != MPI_REQUEST_NULL && not(savedreq->flags() & MPI_REQ_FINISHED)
798     && not(savedreq->flags() & MPI_REQ_GENERALIZED))
799       savedreq->ref();//don't erase te handle in Request::wait, we'll need it later
800     else
801       savedreq = MPI_REQUEST_NULL;
802
803     int my_proc_id = (*request)->comm() != MPI_COMM_NULL
804                          ? simgrid::s4u::this_actor::get_pid()
805                          : -1; // TODO: cheinrich: Check if this correct or if it should be MPI_UNDEFINED
806     TRACE_smpi_comm_in(my_proc_id, __func__,
807                        new simgrid::instr::WaitTIData((*request)->src(), (*request)->dst(), (*request)->tag()));
808
809     retval = simgrid::smpi::Request::wait(request, status);
810
811     //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
812     TRACE_smpi_comm_out(my_proc_id);
813     trace_smpi_recv_helper(&savedreq, status);
814     if (savedreq != MPI_REQUEST_NULL)
815       simgrid::smpi::Request::unref(&savedreq);
816   }
817
818   smpi_bench_begin();
819   return retval;
820 }
821
822 int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * status)
823 {
824   if (index == nullptr)
825     return MPI_ERR_ARG;
826
827   if (count <= 0)
828     return MPI_SUCCESS;
829
830   smpi_bench_end();
831   // for tracing, save the handles which might get overridden before we can use the helper on it
832   std::vector<MPI_Request> savedreqs(requests, requests + count);
833   for (MPI_Request& req : savedreqs) {
834     if (req != MPI_REQUEST_NULL && not(req->flags() & MPI_REQ_FINISHED))
835       req->ref();
836     else
837       req = MPI_REQUEST_NULL;
838   }
839
840   int rank_traced = simgrid::s4u::this_actor::get_pid(); // FIXME: In PMPI_Wait, we check if the comm is null?
841   TRACE_smpi_comm_in(rank_traced, __func__, new simgrid::instr::CpuTIData("waitAny", count));
842
843   *index = simgrid::smpi::Request::waitany(count, requests, status);
844
845   if(*index!=MPI_UNDEFINED){
846     trace_smpi_recv_helper(&savedreqs[*index], status);
847     TRACE_smpi_comm_out(rank_traced);
848   }
849
850   for (MPI_Request& req : savedreqs)
851     if (req != MPI_REQUEST_NULL)
852       simgrid::smpi::Request::unref(&req);
853
854   smpi_bench_begin();
855   return MPI_SUCCESS;
856 }
857
858 int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
859 {
860   smpi_bench_end();
861
862   // for tracing, save the handles which might get overridden before we can use the helper on it
863   std::vector<MPI_Request> savedreqs(requests, requests + count);
864   for (MPI_Request& req : savedreqs) {
865     if (req != MPI_REQUEST_NULL && not(req->flags() & MPI_REQ_FINISHED))
866       req->ref();
867     else
868       req = MPI_REQUEST_NULL;
869   }
870
871   int rank_traced = simgrid::s4u::this_actor::get_pid(); // FIXME: In PMPI_Wait, we check if the comm is null?
872   TRACE_smpi_comm_in(rank_traced, __func__, new simgrid::instr::CpuTIData("waitall", count));
873
874   int retval = simgrid::smpi::Request::waitall(count, requests, status);
875
876   for (int i = 0; i < count; i++) {
877     trace_smpi_recv_helper(&savedreqs[i], status!=MPI_STATUSES_IGNORE ? &status[i]: MPI_STATUS_IGNORE);
878   }
879   TRACE_smpi_comm_out(rank_traced);
880
881   for (MPI_Request& req : savedreqs)
882     if (req != MPI_REQUEST_NULL)
883       simgrid::smpi::Request::unref(&req);
884
885   smpi_bench_begin();
886   return retval;
887 }
888
889 int PMPI_Waitsome(int incount, MPI_Request requests[], int *outcount, int *indices, MPI_Status status[])
890 {
891   int retval = 0;
892
893   smpi_bench_end();
894   if (outcount == nullptr) {
895     retval = MPI_ERR_ARG;
896   } else {
897     *outcount = simgrid::smpi::Request::waitsome(incount, requests, indices, status);
898     retval = MPI_SUCCESS;
899   }
900   smpi_bench_begin();
901   return retval;
902 }
903
904 int PMPI_Cancel(MPI_Request* request)
905 {
906   int retval = 0;
907
908   smpi_bench_end();
909   if (*request == MPI_REQUEST_NULL) {
910     retval = MPI_ERR_REQUEST;
911   } else {
912     (*request)->cancel();
913     retval = MPI_SUCCESS;
914   }
915   smpi_bench_begin();
916   return retval;
917 }
918
919 int PMPI_Test_cancelled(const MPI_Status* status, int* flag){
920   if(status==MPI_STATUS_IGNORE){
921     *flag=0;
922     return MPI_ERR_ARG;
923   }
924   *flag=simgrid::smpi::Status::cancelled(status);
925   return MPI_SUCCESS;  
926 }
927
928 int PMPI_Status_set_cancelled(MPI_Status* status, int flag){
929   if(status==MPI_STATUS_IGNORE){
930     return MPI_ERR_ARG;
931   }
932   simgrid::smpi::Status::set_cancelled(status,flag);
933   return MPI_SUCCESS;  
934 }
935
936 int PMPI_Status_set_elements(MPI_Status* status, MPI_Datatype datatype, int count){
937   if(status==MPI_STATUS_IGNORE){
938     return MPI_ERR_ARG;
939   }
940   simgrid::smpi::Status::set_elements(status,datatype, count);
941   return MPI_SUCCESS;  
942 }
943
944 int PMPI_Grequest_start( MPI_Grequest_query_function *query_fn, MPI_Grequest_free_function *free_fn, MPI_Grequest_cancel_function *cancel_fn, void *extra_state, MPI_Request *request){
945   return simgrid::smpi::Request::grequest_start(query_fn, free_fn,cancel_fn, extra_state, request);
946 }
947
948 int PMPI_Grequest_complete( MPI_Request request){
949   return simgrid::smpi::Request::grequest_complete(request);
950 }
951
952 int PMPI_Request_get_status( MPI_Request request, int *flag, MPI_Status *status){
953   if(request==MPI_REQUEST_NULL){
954     *flag=1;
955     simgrid::smpi::Status::empty(status);
956     return MPI_SUCCESS;
957   } else if (flag==NULL || status ==NULL){
958     return MPI_ERR_ARG;
959   }
960   return simgrid::smpi::Request::get_status(request,flag,status);
961 }
962
963 MPI_Request PMPI_Request_f2c(MPI_Fint request){
964   if(request==-1)
965     return MPI_REQUEST_NULL;
966   return static_cast<MPI_Request>(simgrid::smpi::Request::f2c(request));
967 }
968
969 MPI_Fint PMPI_Request_c2f(MPI_Request request) {
970   if(request==MPI_REQUEST_NULL)
971     return -1;
972   return request->c2f();
973 }