Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
snake_case Engine signals
[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     simgrid::smpi::Status::empty(status);
315     status->MPI_SOURCE = MPI_PROC_NULL;
316     retval = MPI_SUCCESS;
317   } else if (src!=MPI_ANY_SOURCE && (src >= comm->group()->size() || src <0)){
318     retval = MPI_ERR_RANK;
319   } else if ((count < 0) || (buf==nullptr && count > 0)) {
320     retval = MPI_ERR_COUNT;
321   } else if (not datatype->is_valid()) {
322     retval = MPI_ERR_TYPE;
323   } else if(tag<0 && tag !=  MPI_ANY_TAG){
324     retval = MPI_ERR_TAG;
325   } else {
326     int my_proc_id = simgrid::s4u::this_actor::get_pid();
327     TRACE_smpi_comm_in(my_proc_id, __func__,
328                        new simgrid::instr::Pt2PtTIData("recv", src,
329                                                        datatype->is_replayable() ? count : count * datatype->size(),
330                                                        tag, simgrid::smpi::Datatype::encode(datatype)));
331
332     simgrid::smpi::Request::recv(buf, count, datatype, src, tag, comm, status);
333     retval = MPI_SUCCESS;
334
335     // the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
336     if (status != MPI_STATUS_IGNORE) {
337       int src_traced = getPid(comm, status->MPI_SOURCE);
338       if (not TRACE_smpi_view_internals()) {
339         TRACE_smpi_recv(src_traced, my_proc_id, tag);
340       }
341     }
342     TRACE_smpi_comm_out(my_proc_id);
343   }
344
345   smpi_bench_begin();
346   return retval;
347 }
348
349 int PMPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm)
350 {
351   int retval = 0;
352
353   smpi_bench_end();
354
355   if (comm == MPI_COMM_NULL) {
356     retval = MPI_ERR_COMM;
357   } else if (dst == MPI_PROC_NULL) {
358     retval = MPI_SUCCESS;
359   } else if (dst >= comm->group()->size() || dst <0){
360     retval = MPI_ERR_RANK;
361   } else if ((count < 0) || (buf == nullptr && count > 0)) {
362     retval = MPI_ERR_COUNT;
363   } else if (not datatype->is_valid()) {
364     retval = MPI_ERR_TYPE;
365   } else if(tag < 0 && tag !=  MPI_ANY_TAG){
366     retval = MPI_ERR_TAG;
367   } else {
368     int my_proc_id         = simgrid::s4u::this_actor::get_pid();
369     int dst_traced         = getPid(comm, dst);
370     TRACE_smpi_comm_in(my_proc_id, __func__,
371                        new simgrid::instr::Pt2PtTIData("send", dst,
372                                                        datatype->is_replayable() ? count : count * datatype->size(),
373                                                        tag, simgrid::smpi::Datatype::encode(datatype)));
374     if (not TRACE_smpi_view_internals()) {
375       TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, tag, count * datatype->size());
376     }
377
378     simgrid::smpi::Request::send(buf, count, datatype, dst, tag, comm);
379     retval = MPI_SUCCESS;
380
381     TRACE_smpi_comm_out(my_proc_id);
382   }
383
384   smpi_bench_begin();
385   return retval;
386 }
387
388 int PMPI_Ssend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm) {
389   int retval = 0;
390
391   smpi_bench_end();
392
393   if (comm == MPI_COMM_NULL) {
394     retval = MPI_ERR_COMM;
395   } else if (dst == MPI_PROC_NULL) {
396     retval = MPI_SUCCESS;
397   } else if (dst >= comm->group()->size() || dst <0){
398     retval = MPI_ERR_RANK;
399   } else if ((count < 0) || (buf==nullptr && count > 0)) {
400     retval = MPI_ERR_COUNT;
401   } else if (not datatype->is_valid()) {
402     retval = MPI_ERR_TYPE;
403   } else if(tag<0 && tag !=  MPI_ANY_TAG){
404     retval = MPI_ERR_TAG;
405   } else {
406     int my_proc_id         = simgrid::s4u::this_actor::get_pid();
407     int dst_traced         = getPid(comm, dst);
408     TRACE_smpi_comm_in(my_proc_id, __func__,
409                        new simgrid::instr::Pt2PtTIData("Ssend", dst,
410                                                        datatype->is_replayable() ? count : count * datatype->size(),
411                                                        tag, simgrid::smpi::Datatype::encode(datatype)));
412     TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, tag, count * datatype->size());
413
414     simgrid::smpi::Request::ssend(buf, count, datatype, dst, tag, comm);
415     retval = MPI_SUCCESS;
416
417     TRACE_smpi_comm_out(my_proc_id);
418   }
419
420   smpi_bench_begin();
421   return retval;
422 }
423
424 int PMPI_Sendrecv(void* sendbuf, int sendcount, MPI_Datatype sendtype, int dst, int sendtag, void* recvbuf,
425                   int recvcount, MPI_Datatype recvtype, int src, int recvtag, MPI_Comm comm, MPI_Status* status)
426 {
427   int retval = 0;
428
429   smpi_bench_end();
430
431   if (comm == MPI_COMM_NULL) {
432     retval = MPI_ERR_COMM;
433   } else if (not sendtype->is_valid() || not recvtype->is_valid()) {
434     retval = MPI_ERR_TYPE;
435   } else if (src == MPI_PROC_NULL) {
436     if(status!=MPI_STATUS_IGNORE){
437       simgrid::smpi::Status::empty(status);
438       status->MPI_SOURCE = MPI_PROC_NULL;
439     }
440     if(dst != MPI_PROC_NULL)
441       simgrid::smpi::Request::send(sendbuf, sendcount, sendtype, dst, sendtag, comm);
442     retval = MPI_SUCCESS;
443   }else if (dst == MPI_PROC_NULL){
444     simgrid::smpi::Request::recv(recvbuf, recvcount, recvtype, src, recvtag, comm, status);
445     retval = MPI_SUCCESS;
446   }else if (dst >= comm->group()->size() || dst <0 ||
447       (src!=MPI_ANY_SOURCE && (src >= comm->group()->size() || src <0))){
448     retval = MPI_ERR_RANK;
449   } else if ((sendcount < 0 || recvcount<0) ||
450       (sendbuf==nullptr && sendcount > 0) || (recvbuf==nullptr && recvcount>0)) {
451     retval = MPI_ERR_COUNT;
452   } else if((sendtag<0 && sendtag !=  MPI_ANY_TAG)||(recvtag<0 && recvtag != MPI_ANY_TAG)){
453     retval = MPI_ERR_TAG;
454   } else {
455     int my_proc_id         = simgrid::s4u::this_actor::get_pid();
456     int dst_traced         = getPid(comm, dst);
457     int src_traced         = getPid(comm, src);
458
459     // FIXME: Hack the way to trace this one
460     std::vector<int>* dst_hack = new std::vector<int>;
461     std::vector<int>* src_hack = new std::vector<int>;
462     dst_hack->push_back(dst_traced);
463     src_hack->push_back(src_traced);
464     TRACE_smpi_comm_in(my_proc_id, __func__,
465                        new simgrid::instr::VarCollTIData(
466                            "sendRecv", -1, sendtype->is_replayable() ? sendcount : sendcount * sendtype->size(),
467                            dst_hack, recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(), src_hack,
468                            simgrid::smpi::Datatype::encode(sendtype), simgrid::smpi::Datatype::encode(recvtype)));
469
470     TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, sendtag, sendcount * sendtype->size());
471
472     simgrid::smpi::Request::sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf, recvcount, recvtype, src,
473                                      recvtag, comm, status);
474     retval = MPI_SUCCESS;
475
476     TRACE_smpi_recv(src_traced, my_proc_id, recvtag);
477     TRACE_smpi_comm_out(my_proc_id);
478   }
479
480   smpi_bench_begin();
481   return retval;
482 }
483
484 int PMPI_Sendrecv_replace(void* buf, int count, MPI_Datatype datatype, int dst, int sendtag, int src, int recvtag,
485                           MPI_Comm comm, MPI_Status* status)
486 {
487   int retval = 0;
488   if (not datatype->is_valid()) {
489     return MPI_ERR_TYPE;
490   } else if (count < 0) {
491     return MPI_ERR_COUNT;
492   } else {
493     int size = datatype->get_extent() * count;
494     void* recvbuf = xbt_new0(char, size);
495     retval = MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count, datatype, src, recvtag, comm, status);
496     if(retval==MPI_SUCCESS){
497       simgrid::smpi::Datatype::copy(recvbuf, count, datatype, buf, count, datatype);
498     }
499     xbt_free(recvbuf);
500
501   }
502   return retval;
503 }
504
505 int PMPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
506 {
507   int retval = 0;
508   smpi_bench_end();
509   if (request == nullptr || flag == nullptr) {
510     retval = MPI_ERR_ARG;
511   } else if (*request == MPI_REQUEST_NULL) {
512     *flag= true;
513     simgrid::smpi::Status::empty(status);
514     retval = MPI_SUCCESS;
515   } else {
516     int my_proc_id = ((*request)->comm() != MPI_COMM_NULL) ? simgrid::s4u::this_actor::get_pid() : -1;
517
518     TRACE_smpi_testing_in(my_proc_id);
519
520     *flag = simgrid::smpi::Request::test(request,status);
521
522     TRACE_smpi_testing_out(my_proc_id);
523     retval = MPI_SUCCESS;
524   }
525   smpi_bench_begin();
526   return retval;
527 }
528
529 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag, MPI_Status * status)
530 {
531   int retval = 0;
532
533   smpi_bench_end();
534   if (index == nullptr || flag == nullptr) {
535     retval = MPI_ERR_ARG;
536   } else {
537     *flag = simgrid::smpi::Request::testany(count, requests, index, status);
538     retval = MPI_SUCCESS;
539   }
540   smpi_bench_begin();
541   return retval;
542 }
543
544 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses)
545 {
546   int retval = 0;
547
548   smpi_bench_end();
549   if (flag == nullptr) {
550     retval = MPI_ERR_ARG;
551   } else {
552     *flag = simgrid::smpi::Request::testall(count, requests, statuses);
553     retval = MPI_SUCCESS;
554   }
555   smpi_bench_begin();
556   return retval;
557 }
558
559 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
560   int retval = 0;
561   smpi_bench_end();
562
563   if (status == nullptr) {
564     retval = MPI_ERR_ARG;
565   } else if (comm == MPI_COMM_NULL) {
566     retval = MPI_ERR_COMM;
567   } else if (source == MPI_PROC_NULL) {
568     simgrid::smpi::Status::empty(status);
569     status->MPI_SOURCE = MPI_PROC_NULL;
570     retval = MPI_SUCCESS;
571   } else {
572     simgrid::smpi::Request::probe(source, tag, comm, status);
573     retval = MPI_SUCCESS;
574   }
575   smpi_bench_begin();
576   return retval;
577 }
578
579 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
580   int retval = 0;
581   smpi_bench_end();
582
583   if (flag == nullptr) {
584     retval = MPI_ERR_ARG;
585   } else if (comm == MPI_COMM_NULL) {
586     retval = MPI_ERR_COMM;
587   } else if (source == MPI_PROC_NULL) {
588     *flag=true;
589     simgrid::smpi::Status::empty(status);
590     status->MPI_SOURCE = MPI_PROC_NULL;
591     retval = MPI_SUCCESS;
592   } else {
593     simgrid::smpi::Request::iprobe(source, tag, comm, flag, status);
594     retval = MPI_SUCCESS;
595   }
596   smpi_bench_begin();
597   return retval;
598 }
599
600 // TODO: cheinrich: Move declaration to other file? Rename this function - it's used for PMPI_Wait*?
601 static void trace_smpi_recv_helper(MPI_Request* request, MPI_Status* status);
602 static void trace_smpi_recv_helper(MPI_Request* request, MPI_Status* status)
603 {
604   MPI_Request req = *request;
605   if (req != MPI_REQUEST_NULL) { // Received requests become null
606     int src_traced = req->src();
607     // the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
608     int dst_traced = req->dst();
609     if (req->flags() & RECV) { // Is this request a wait for RECV?
610       if (src_traced == MPI_ANY_SOURCE)
611         src_traced = (status != MPI_STATUSES_IGNORE) ? req->comm()->group()->rank(status->MPI_SOURCE) : req->src();
612       TRACE_smpi_recv(src_traced, dst_traced, req->tag());
613     }
614   }
615 }
616
617 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
618 {
619   int retval = 0;
620
621   smpi_bench_end();
622
623   simgrid::smpi::Status::empty(status);
624
625   if (request == nullptr) {
626     retval = MPI_ERR_ARG;
627   } else if (*request == MPI_REQUEST_NULL) {
628     retval = MPI_SUCCESS;
629   } else {
630     int my_proc_id = (*request)->comm() != MPI_COMM_NULL
631                          ? simgrid::s4u::this_actor::get_pid()
632                          : -1; // TODO: cheinrich: Check if this correct or if it should be MPI_UNDEFINED
633
634     TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("wait"));
635
636     simgrid::smpi::Request::wait(request, status);
637     retval = MPI_SUCCESS;
638
639     //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
640     TRACE_smpi_comm_out(my_proc_id);
641     trace_smpi_recv_helper(request, status);
642   }
643
644   smpi_bench_begin();
645   return retval;
646 }
647
648 int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * status)
649 {
650   if (index == nullptr)
651     return MPI_ERR_ARG;
652
653   if (count <= 0)
654     return MPI_SUCCESS;
655
656   smpi_bench_end();
657
658   int rank_traced = simgrid::s4u::this_actor::get_pid(); // FIXME: In PMPI_Wait, we check if the comm is null?
659   TRACE_smpi_comm_in(rank_traced, __func__, new simgrid::instr::CpuTIData("waitAny", static_cast<double>(count)));
660
661   *index = simgrid::smpi::Request::waitany(count, requests, status);
662
663   if(*index!=MPI_UNDEFINED){
664     trace_smpi_recv_helper(&requests[*index], status);
665     TRACE_smpi_comm_out(rank_traced);
666   }
667
668   smpi_bench_begin();
669   return MPI_SUCCESS;
670 }
671
672 int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
673 {
674   smpi_bench_end();
675
676   int rank_traced = simgrid::s4u::this_actor::get_pid(); // FIXME: In PMPI_Wait, we check if the comm is null?
677   TRACE_smpi_comm_in(rank_traced, __func__, new simgrid::instr::CpuTIData("waitAll", static_cast<double>(count)));
678
679   int retval = simgrid::smpi::Request::waitall(count, requests, status);
680
681   for (int i = 0; i < count; i++) {
682     trace_smpi_recv_helper(&requests[i], &status[i]);
683   }
684   TRACE_smpi_comm_out(rank_traced);
685
686   smpi_bench_begin();
687   return retval;
688 }
689
690 int PMPI_Waitsome(int incount, MPI_Request requests[], int *outcount, int *indices, MPI_Status status[])
691 {
692   int retval = 0;
693
694   smpi_bench_end();
695   if (outcount == nullptr) {
696     retval = MPI_ERR_ARG;
697   } else {
698     *outcount = simgrid::smpi::Request::waitsome(incount, requests, indices, status);
699     retval = MPI_SUCCESS;
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     *outcount = simgrid::smpi::Request::testsome(incount, requests, indices, status);
714     retval    = MPI_SUCCESS;
715   }
716   smpi_bench_begin();
717   return retval;
718 }
719
720 int PMPI_Cancel(MPI_Request* request)
721 {
722   int retval = 0;
723
724   smpi_bench_end();
725   if (*request == MPI_REQUEST_NULL) {
726     retval = MPI_ERR_REQUEST;
727   } else {
728     (*request)->cancel();
729     retval = MPI_SUCCESS;
730   }
731   smpi_bench_begin();
732   return retval;
733 }
734
735 int PMPI_Test_cancelled(MPI_Status* status, int* flag){
736   if(status==MPI_STATUS_IGNORE){
737     *flag=0;
738     return MPI_ERR_ARG;
739   }
740   *flag=simgrid::smpi::Status::cancelled(status);
741   return MPI_SUCCESS;  
742 }
743
744 MPI_Request PMPI_Request_f2c(MPI_Fint request){
745   return static_cast<MPI_Request>(simgrid::smpi::Request::f2c(request));
746 }
747
748 MPI_Fint PMPI_Request_c2f(MPI_Request request) {
749   return request->c2f();
750 }