Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
fix bug on MPI_Recv when MPI_IGNORE_STATUS is used with MPI_PROC_NULL as a sender
[simgrid.git] / src / smpi / bindings / smpi_pmpi_request.cpp
1 /* Copyright (c) 2007-2018. 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_process.hpp"
10 #include "smpi_request.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(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 (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 (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_Ssend_init(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request)
70 {
71   int retval = 0;
72
73   smpi_bench_end();
74   if (request == nullptr) {
75     retval = MPI_ERR_ARG;
76   } else if (comm == MPI_COMM_NULL) {
77     retval = MPI_ERR_COMM;
78   } else if (not datatype->is_valid()) {
79     retval = MPI_ERR_TYPE;
80   } else if (dst == MPI_PROC_NULL) {
81     retval = MPI_SUCCESS;
82   } else {
83     *request = simgrid::smpi::Request::ssend_init(buf, count, datatype, dst, tag, comm);
84     retval = MPI_SUCCESS;
85   }
86   smpi_bench_begin();
87   if (retval != MPI_SUCCESS && request != nullptr)
88     *request = MPI_REQUEST_NULL;
89   return retval;
90 }
91
92 /*
93  * This function starts a request returned by init functions such as
94  * MPI_Send_init(), MPI_Ssend_init (see above), and friends.
95  * They should already have performed sanity checks.
96  */
97 int PMPI_Start(MPI_Request * request)
98 {
99   int retval = 0;
100
101   smpi_bench_end();
102   if (request == nullptr || *request == MPI_REQUEST_NULL) {
103     retval = MPI_ERR_REQUEST;
104   } else {
105     MPI_Request req = *request;
106     int my_proc_id = (req->comm() != MPI_COMM_NULL) ? simgrid::s4u::this_actor::get_pid() : -1;
107     TRACE_smpi_comm_in(my_proc_id, __func__,
108                        new simgrid::instr::Pt2PtTIData("Start", req->dst(),
109                                                        req->size(),
110                                                        req->tag(), 
111                                                        simgrid::smpi::Datatype::encode(req->type())));
112     if (not TRACE_smpi_view_internals() && req->flags() & SEND)
113       TRACE_smpi_send(my_proc_id, my_proc_id, getPid(req->comm(), req->dst()), req->tag(), req->size());
114
115     req->start();
116
117     if (not TRACE_smpi_view_internals() && req->flags() & RECV)
118       TRACE_smpi_recv(getPid(req->comm(), req->src()), my_proc_id, req->tag());
119     retval = MPI_SUCCESS;
120     TRACE_smpi_comm_out(my_proc_id);
121   }
122   smpi_bench_begin();
123   return retval;
124 }
125
126 int PMPI_Startall(int count, MPI_Request * requests)
127 {
128   int retval;
129   smpi_bench_end();
130   if (requests == nullptr) {
131     retval = MPI_ERR_ARG;
132   } else {
133     retval = MPI_SUCCESS;
134     for (int i = 0; i < count; i++) {
135       if(requests[i] == MPI_REQUEST_NULL) {
136         retval = MPI_ERR_REQUEST;
137       }
138     }
139     if(retval != MPI_ERR_REQUEST) {
140       int my_proc_id = simgrid::s4u::this_actor::get_pid();
141       TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("Startall"));
142       MPI_Request req = MPI_REQUEST_NULL;
143       if (not TRACE_smpi_view_internals())
144         for (int i = 0; i < count; i++) {
145           req = requests[i];
146           if (req->flags() & SEND)
147             TRACE_smpi_send(my_proc_id, my_proc_id, getPid(req->comm(), req->dst()), req->tag(), req->size());
148         }
149
150       simgrid::smpi::Request::startall(count, requests);
151
152       if (not TRACE_smpi_view_internals())
153         for (int i = 0; i < count; i++) {
154           req = requests[i];
155           if (req->flags() & RECV)
156             TRACE_smpi_recv(getPid(req->comm(), req->src()), my_proc_id, req->tag());
157         }
158       TRACE_smpi_comm_out(my_proc_id);
159     }
160   }
161   smpi_bench_begin();
162   return retval;
163 }
164
165 int PMPI_Request_free(MPI_Request * request)
166 {
167   int retval = 0;
168
169   smpi_bench_end();
170   if (*request == MPI_REQUEST_NULL) {
171     retval = MPI_ERR_ARG;
172   } else {
173     simgrid::smpi::Request::unref(request);
174     retval = MPI_SUCCESS;
175   }
176   smpi_bench_begin();
177   return retval;
178 }
179
180 int PMPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Request * request)
181 {
182   int retval = 0;
183
184   smpi_bench_end();
185
186   if (request == nullptr) {
187     retval = MPI_ERR_ARG;
188   } else if (comm == MPI_COMM_NULL) {
189     retval = MPI_ERR_COMM;
190   } else if (src == MPI_PROC_NULL) {
191     *request = MPI_REQUEST_NULL;
192     retval = MPI_SUCCESS;
193   } else if (src!=MPI_ANY_SOURCE && (src >= comm->group()->size() || src <0)){
194     retval = MPI_ERR_RANK;
195   } else if ((count < 0) || (buf==nullptr && count > 0)) {
196     retval = MPI_ERR_COUNT;
197   } else if (not datatype->is_valid()) {
198     retval = MPI_ERR_TYPE;
199   } else if(tag<0 && tag !=  MPI_ANY_TAG){
200     retval = MPI_ERR_TAG;
201   } else {
202
203     int my_proc_id = simgrid::s4u::this_actor::get_pid();
204
205     TRACE_smpi_comm_in(my_proc_id, __func__,
206                        new simgrid::instr::Pt2PtTIData("Irecv", src,
207                                                        datatype->is_replayable() ? count : count * datatype->size(),
208                                                        tag, simgrid::smpi::Datatype::encode(datatype)));
209
210     *request = simgrid::smpi::Request::irecv(buf, count, datatype, src, tag, comm);
211     retval = MPI_SUCCESS;
212
213     TRACE_smpi_comm_out(my_proc_id);
214   }
215
216   smpi_bench_begin();
217   if (retval != MPI_SUCCESS && request != nullptr)
218     *request = MPI_REQUEST_NULL;
219   return retval;
220 }
221
222
223 int PMPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request * request)
224 {
225   int retval = 0;
226
227   smpi_bench_end();
228   if (request == nullptr) {
229     retval = MPI_ERR_ARG;
230   } else if (comm == MPI_COMM_NULL) {
231     retval = MPI_ERR_COMM;
232   } else if (dst == MPI_PROC_NULL) {
233     *request = MPI_REQUEST_NULL;
234     retval = MPI_SUCCESS;
235   } else if (dst >= comm->group()->size() || dst <0){
236     retval = MPI_ERR_RANK;
237   } else if ((count < 0) || (buf==nullptr && count > 0)) {
238     retval = MPI_ERR_COUNT;
239   } else if (not datatype->is_valid()) {
240     retval = MPI_ERR_TYPE;
241   } else if(tag<0 && tag !=  MPI_ANY_TAG){
242     retval = MPI_ERR_TAG;
243   } else {
244     int my_proc_id = simgrid::s4u::this_actor::get_pid();
245     int trace_dst = getPid(comm, dst);
246     TRACE_smpi_comm_in(my_proc_id, __func__,
247                        new simgrid::instr::Pt2PtTIData("Isend", dst,
248                                                        datatype->is_replayable() ? count : count * datatype->size(),
249                                                        tag, simgrid::smpi::Datatype::encode(datatype)));
250
251     TRACE_smpi_send(my_proc_id, my_proc_id, trace_dst, tag, count * datatype->size());
252
253     *request = simgrid::smpi::Request::isend(buf, count, datatype, dst, tag, comm);
254     retval = MPI_SUCCESS;
255
256     TRACE_smpi_comm_out(my_proc_id);
257   }
258
259   smpi_bench_begin();
260   if (retval != MPI_SUCCESS && request!=nullptr)
261     *request = MPI_REQUEST_NULL;
262   return retval;
263 }
264
265 int PMPI_Issend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request)
266 {
267   int retval = 0;
268
269   smpi_bench_end();
270   if (request == nullptr) {
271     retval = MPI_ERR_ARG;
272   } else if (comm == MPI_COMM_NULL) {
273     retval = MPI_ERR_COMM;
274   } else if (dst == MPI_PROC_NULL) {
275     *request = MPI_REQUEST_NULL;
276     retval = MPI_SUCCESS;
277   } else if (dst >= comm->group()->size() || dst <0){
278     retval = MPI_ERR_RANK;
279   } else if ((count < 0)|| (buf==nullptr && count > 0)) {
280     retval = MPI_ERR_COUNT;
281   } else if (not datatype->is_valid()) {
282     retval = MPI_ERR_TYPE;
283   } else if(tag<0 && tag !=  MPI_ANY_TAG){
284     retval = MPI_ERR_TAG;
285   } else {
286     int my_proc_id = simgrid::s4u::this_actor::get_pid();
287     int trace_dst = getPid(comm, dst);
288     TRACE_smpi_comm_in(my_proc_id, __func__,
289                        new simgrid::instr::Pt2PtTIData("ISsend", dst,
290                                                        datatype->is_replayable() ? count : count * datatype->size(),
291                                                        tag, simgrid::smpi::Datatype::encode(datatype)));
292     TRACE_smpi_send(my_proc_id, my_proc_id, trace_dst, tag, count * datatype->size());
293
294     *request = simgrid::smpi::Request::issend(buf, count, datatype, dst, tag, comm);
295     retval = MPI_SUCCESS;
296
297     TRACE_smpi_comm_out(my_proc_id);
298   }
299
300   smpi_bench_begin();
301   if (retval != MPI_SUCCESS && request!=nullptr)
302     *request = MPI_REQUEST_NULL;
303   return retval;
304 }
305
306 int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Status * status)
307 {
308   int retval = 0;
309
310   smpi_bench_end();
311   if (comm == MPI_COMM_NULL) {
312     retval = MPI_ERR_COMM;
313   } else if (src == MPI_PROC_NULL) {
314     if(status != MPI_STATUS_IGNORE){
315       simgrid::smpi::Status::empty(status);
316       status->MPI_SOURCE = MPI_PROC_NULL;
317     }
318     retval = MPI_SUCCESS;
319   } else if (src!=MPI_ANY_SOURCE && (src >= comm->group()->size() || src <0)){
320     retval = MPI_ERR_RANK;
321   } else if ((count < 0) || (buf==nullptr && count > 0)) {
322     retval = MPI_ERR_COUNT;
323   } else if (not datatype->is_valid()) {
324     retval = MPI_ERR_TYPE;
325   } else if(tag<0 && tag !=  MPI_ANY_TAG){
326     retval = MPI_ERR_TAG;
327   } else {
328     int my_proc_id = simgrid::s4u::this_actor::get_pid();
329     TRACE_smpi_comm_in(my_proc_id, __func__,
330                        new simgrid::instr::Pt2PtTIData("recv", src,
331                                                        datatype->is_replayable() ? count : count * datatype->size(),
332                                                        tag, simgrid::smpi::Datatype::encode(datatype)));
333
334     simgrid::smpi::Request::recv(buf, count, datatype, src, tag, comm, status);
335     retval = MPI_SUCCESS;
336
337     // the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
338     if (status != MPI_STATUS_IGNORE) {
339       int src_traced = getPid(comm, status->MPI_SOURCE);
340       if (not TRACE_smpi_view_internals()) {
341         TRACE_smpi_recv(src_traced, my_proc_id, tag);
342       }
343     }
344     TRACE_smpi_comm_out(my_proc_id);
345   }
346
347   smpi_bench_begin();
348   return retval;
349 }
350
351 int PMPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm)
352 {
353   int retval = 0;
354
355   smpi_bench_end();
356
357   if (comm == MPI_COMM_NULL) {
358     retval = MPI_ERR_COMM;
359   } else if (dst == MPI_PROC_NULL) {
360     retval = MPI_SUCCESS;
361   } else if (dst >= comm->group()->size() || dst <0){
362     retval = MPI_ERR_RANK;
363   } else if ((count < 0) || (buf == nullptr && count > 0)) {
364     retval = MPI_ERR_COUNT;
365   } else if (not datatype->is_valid()) {
366     retval = MPI_ERR_TYPE;
367   } else if(tag < 0 && tag !=  MPI_ANY_TAG){
368     retval = MPI_ERR_TAG;
369   } else {
370     int my_proc_id         = simgrid::s4u::this_actor::get_pid();
371     int dst_traced         = getPid(comm, dst);
372     TRACE_smpi_comm_in(my_proc_id, __func__,
373                        new simgrid::instr::Pt2PtTIData("send", dst,
374                                                        datatype->is_replayable() ? count : count * datatype->size(),
375                                                        tag, simgrid::smpi::Datatype::encode(datatype)));
376     if (not TRACE_smpi_view_internals()) {
377       TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, tag, count * datatype->size());
378     }
379
380     simgrid::smpi::Request::send(buf, count, datatype, dst, tag, comm);
381     retval = MPI_SUCCESS;
382
383     TRACE_smpi_comm_out(my_proc_id);
384   }
385
386   smpi_bench_begin();
387   return retval;
388 }
389
390 int PMPI_Ssend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm) {
391   int retval = 0;
392
393   smpi_bench_end();
394
395   if (comm == MPI_COMM_NULL) {
396     retval = MPI_ERR_COMM;
397   } else if (dst == MPI_PROC_NULL) {
398     retval = MPI_SUCCESS;
399   } else if (dst >= comm->group()->size() || dst <0){
400     retval = MPI_ERR_RANK;
401   } else if ((count < 0) || (buf==nullptr && count > 0)) {
402     retval = MPI_ERR_COUNT;
403   } else if (not datatype->is_valid()) {
404     retval = MPI_ERR_TYPE;
405   } else if(tag<0 && tag !=  MPI_ANY_TAG){
406     retval = MPI_ERR_TAG;
407   } else {
408     int my_proc_id         = simgrid::s4u::this_actor::get_pid();
409     int dst_traced         = getPid(comm, dst);
410     TRACE_smpi_comm_in(my_proc_id, __func__,
411                        new simgrid::instr::Pt2PtTIData("Ssend", dst,
412                                                        datatype->is_replayable() ? count : count * datatype->size(),
413                                                        tag, simgrid::smpi::Datatype::encode(datatype)));
414     TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, tag, count * datatype->size());
415
416     simgrid::smpi::Request::ssend(buf, count, datatype, dst, tag, comm);
417     retval = MPI_SUCCESS;
418
419     TRACE_smpi_comm_out(my_proc_id);
420   }
421
422   smpi_bench_begin();
423   return retval;
424 }
425
426 int PMPI_Sendrecv(void* sendbuf, int sendcount, MPI_Datatype sendtype, int dst, int sendtag, void* recvbuf,
427                   int recvcount, MPI_Datatype recvtype, int src, int recvtag, MPI_Comm comm, MPI_Status* status)
428 {
429   int retval = 0;
430
431   smpi_bench_end();
432
433   if (comm == MPI_COMM_NULL) {
434     retval = MPI_ERR_COMM;
435   } else if (not sendtype->is_valid() || not recvtype->is_valid()) {
436     retval = MPI_ERR_TYPE;
437   } else if (src == MPI_PROC_NULL) {
438     if(status!=MPI_STATUS_IGNORE){
439       simgrid::smpi::Status::empty(status);
440       status->MPI_SOURCE = MPI_PROC_NULL;
441     }
442     if(dst != MPI_PROC_NULL)
443       simgrid::smpi::Request::send(sendbuf, sendcount, sendtype, dst, sendtag, comm);
444     retval = MPI_SUCCESS;
445   }else if (dst == MPI_PROC_NULL){
446     simgrid::smpi::Request::recv(recvbuf, recvcount, recvtype, src, recvtag, comm, status);
447     retval = MPI_SUCCESS;
448   }else if (dst >= comm->group()->size() || dst <0 ||
449       (src!=MPI_ANY_SOURCE && (src >= comm->group()->size() || src <0))){
450     retval = MPI_ERR_RANK;
451   } else if ((sendcount < 0 || recvcount<0) ||
452       (sendbuf==nullptr && sendcount > 0) || (recvbuf==nullptr && recvcount>0)) {
453     retval = MPI_ERR_COUNT;
454   } else if((sendtag<0 && sendtag !=  MPI_ANY_TAG)||(recvtag<0 && recvtag != MPI_ANY_TAG)){
455     retval = MPI_ERR_TAG;
456   } else {
457     int my_proc_id         = simgrid::s4u::this_actor::get_pid();
458     int dst_traced         = getPid(comm, dst);
459     int src_traced         = getPid(comm, src);
460
461     // FIXME: Hack the way to trace this one
462     std::vector<int>* dst_hack = new std::vector<int>;
463     std::vector<int>* src_hack = new std::vector<int>;
464     dst_hack->push_back(dst_traced);
465     src_hack->push_back(src_traced);
466     TRACE_smpi_comm_in(my_proc_id, __func__,
467                        new simgrid::instr::VarCollTIData(
468                            "sendRecv", -1, sendtype->is_replayable() ? sendcount : sendcount * sendtype->size(),
469                            dst_hack, recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(), src_hack,
470                            simgrid::smpi::Datatype::encode(sendtype), simgrid::smpi::Datatype::encode(recvtype)));
471
472     TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, sendtag, sendcount * sendtype->size());
473
474     simgrid::smpi::Request::sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf, recvcount, recvtype, src,
475                                      recvtag, comm, status);
476     retval = MPI_SUCCESS;
477
478     TRACE_smpi_recv(src_traced, my_proc_id, recvtag);
479     TRACE_smpi_comm_out(my_proc_id);
480   }
481
482   smpi_bench_begin();
483   return retval;
484 }
485
486 int PMPI_Sendrecv_replace(void* buf, int count, MPI_Datatype datatype, int dst, int sendtag, int src, int recvtag,
487                           MPI_Comm comm, MPI_Status* status)
488 {
489   int retval = 0;
490   if (not datatype->is_valid()) {
491     return MPI_ERR_TYPE;
492   } else if (count < 0) {
493     return MPI_ERR_COUNT;
494   } else {
495     int size = datatype->get_extent() * count;
496     void* recvbuf = xbt_new0(char, size);
497     retval = MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count, datatype, src, recvtag, comm, status);
498     if(retval==MPI_SUCCESS){
499       simgrid::smpi::Datatype::copy(recvbuf, count, datatype, buf, count, datatype);
500     }
501     xbt_free(recvbuf);
502
503   }
504   return retval;
505 }
506
507 int PMPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
508 {
509   int retval = 0;
510   smpi_bench_end();
511   if (request == nullptr || flag == nullptr) {
512     retval = MPI_ERR_ARG;
513   } else if (*request == MPI_REQUEST_NULL) {
514     *flag= true;
515     simgrid::smpi::Status::empty(status);
516     retval = MPI_SUCCESS;
517   } else {
518     int my_proc_id = ((*request)->comm() != MPI_COMM_NULL) ? simgrid::s4u::this_actor::get_pid() : -1;
519
520     TRACE_smpi_testing_in(my_proc_id);
521
522     *flag = simgrid::smpi::Request::test(request,status);
523
524     TRACE_smpi_testing_out(my_proc_id);
525     retval = MPI_SUCCESS;
526   }
527   smpi_bench_begin();
528   return retval;
529 }
530
531 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag, MPI_Status * status)
532 {
533   int retval = 0;
534
535   smpi_bench_end();
536   if (index == nullptr || flag == nullptr) {
537     retval = MPI_ERR_ARG;
538   } else {
539     *flag = simgrid::smpi::Request::testany(count, requests, index, status);
540     retval = MPI_SUCCESS;
541   }
542   smpi_bench_begin();
543   return retval;
544 }
545
546 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses)
547 {
548   int retval = 0;
549
550   smpi_bench_end();
551   if (flag == nullptr) {
552     retval = MPI_ERR_ARG;
553   } else {
554     *flag = simgrid::smpi::Request::testall(count, requests, statuses);
555     retval = MPI_SUCCESS;
556   }
557   smpi_bench_begin();
558   return retval;
559 }
560
561 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
562   int retval = 0;
563   smpi_bench_end();
564
565   if (status == nullptr) {
566     retval = MPI_ERR_ARG;
567   } else if (comm == MPI_COMM_NULL) {
568     retval = MPI_ERR_COMM;
569   } else if (source == MPI_PROC_NULL) {
570     simgrid::smpi::Status::empty(status);
571     status->MPI_SOURCE = MPI_PROC_NULL;
572     retval = MPI_SUCCESS;
573   } else {
574     simgrid::smpi::Request::probe(source, tag, comm, status);
575     retval = MPI_SUCCESS;
576   }
577   smpi_bench_begin();
578   return retval;
579 }
580
581 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
582   int retval = 0;
583   smpi_bench_end();
584
585   if (flag == nullptr) {
586     retval = MPI_ERR_ARG;
587   } else if (comm == MPI_COMM_NULL) {
588     retval = MPI_ERR_COMM;
589   } else if (source == MPI_PROC_NULL) {
590     *flag=true;
591     simgrid::smpi::Status::empty(status);
592     status->MPI_SOURCE = MPI_PROC_NULL;
593     retval = MPI_SUCCESS;
594   } else {
595     simgrid::smpi::Request::iprobe(source, tag, comm, flag, status);
596     retval = MPI_SUCCESS;
597   }
598   smpi_bench_begin();
599   return retval;
600 }
601
602 // TODO: cheinrich: Move declaration to other file? Rename this function - it's used for PMPI_Wait*?
603 static void trace_smpi_recv_helper(MPI_Request* request, MPI_Status* status);
604 static void trace_smpi_recv_helper(MPI_Request* request, MPI_Status* status)
605 {
606   MPI_Request req = *request;
607   if (req != MPI_REQUEST_NULL) { // Received requests become null
608     int src_traced = req->src();
609     // the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
610     int dst_traced = req->dst();
611     if (req->flags() & RECV) { // Is this request a wait for RECV?
612       if (src_traced == MPI_ANY_SOURCE)
613         src_traced = (status != MPI_STATUSES_IGNORE) ? req->comm()->group()->rank(status->MPI_SOURCE) : req->src();
614       TRACE_smpi_recv(src_traced, dst_traced, req->tag());
615     }
616   }
617 }
618
619 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
620 {
621   int retval = 0;
622
623   smpi_bench_end();
624
625   simgrid::smpi::Status::empty(status);
626
627   if (request == nullptr) {
628     retval = MPI_ERR_ARG;
629   } else if (*request == MPI_REQUEST_NULL) {
630     retval = MPI_SUCCESS;
631   } else {
632     int my_proc_id = (*request)->comm() != MPI_COMM_NULL
633                          ? simgrid::s4u::this_actor::get_pid()
634                          : -1; // TODO: cheinrich: Check if this correct or if it should be MPI_UNDEFINED
635
636     TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("wait"));
637
638     simgrid::smpi::Request::wait(request, status);
639     retval = MPI_SUCCESS;
640
641     //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
642     TRACE_smpi_comm_out(my_proc_id);
643     trace_smpi_recv_helper(request, status);
644   }
645
646   smpi_bench_begin();
647   return retval;
648 }
649
650 int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * status)
651 {
652   if (index == nullptr)
653     return MPI_ERR_ARG;
654
655   if (count <= 0)
656     return MPI_SUCCESS;
657
658   smpi_bench_end();
659
660   int rank_traced = simgrid::s4u::this_actor::get_pid(); // FIXME: In PMPI_Wait, we check if the comm is null?
661   TRACE_smpi_comm_in(rank_traced, __func__, new simgrid::instr::CpuTIData("waitAny", static_cast<double>(count)));
662
663   *index = simgrid::smpi::Request::waitany(count, requests, status);
664
665   if(*index!=MPI_UNDEFINED){
666     trace_smpi_recv_helper(&requests[*index], status);
667     TRACE_smpi_comm_out(rank_traced);
668   }
669
670   smpi_bench_begin();
671   return MPI_SUCCESS;
672 }
673
674 int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
675 {
676   smpi_bench_end();
677
678   int rank_traced = simgrid::s4u::this_actor::get_pid(); // FIXME: In PMPI_Wait, we check if the comm is null?
679   TRACE_smpi_comm_in(rank_traced, __func__, new simgrid::instr::CpuTIData("waitAll", static_cast<double>(count)));
680
681   int retval = simgrid::smpi::Request::waitall(count, requests, status);
682
683   for (int i = 0; i < count; i++) {
684     trace_smpi_recv_helper(&requests[i], &status[i]);
685   }
686   TRACE_smpi_comm_out(rank_traced);
687
688   smpi_bench_begin();
689   return retval;
690 }
691
692 int PMPI_Waitsome(int incount, MPI_Request requests[], int *outcount, int *indices, MPI_Status status[])
693 {
694   int retval = 0;
695
696   smpi_bench_end();
697   if (outcount == nullptr) {
698     retval = MPI_ERR_ARG;
699   } else {
700     *outcount = simgrid::smpi::Request::waitsome(incount, requests, indices, status);
701     retval = MPI_SUCCESS;
702   }
703   smpi_bench_begin();
704   return retval;
705 }
706
707 int PMPI_Testsome(int incount, MPI_Request requests[], int* outcount, int* indices, MPI_Status status[])
708 {
709   int retval = 0;
710
711   smpi_bench_end();
712   if (outcount == nullptr) {
713     retval = MPI_ERR_ARG;
714   } else {
715     *outcount = simgrid::smpi::Request::testsome(incount, requests, indices, status);
716     retval    = MPI_SUCCESS;
717   }
718   smpi_bench_begin();
719   return retval;
720 }
721
722 int PMPI_Cancel(MPI_Request* request)
723 {
724   int retval = 0;
725
726   smpi_bench_end();
727   if (*request == MPI_REQUEST_NULL) {
728     retval = MPI_ERR_REQUEST;
729   } else {
730     (*request)->cancel();
731     retval = MPI_SUCCESS;
732   }
733   smpi_bench_begin();
734   return retval;
735 }
736
737 int PMPI_Test_cancelled(MPI_Status* status, int* flag){
738   if(status==MPI_STATUS_IGNORE){
739     *flag=0;
740     return MPI_ERR_ARG;
741   }
742   *flag=simgrid::smpi::Status::cancelled(status);
743   return MPI_SUCCESS;  
744 }
745
746 MPI_Request PMPI_Request_f2c(MPI_Fint request){
747   return static_cast<MPI_Request>(simgrid::smpi::Request::f2c(request));
748 }
749
750 MPI_Fint PMPI_Request_c2f(MPI_Request request) {
751   return request->c2f();
752 }