Logo AND Algorithmique Numérique Distribuée

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