Logo AND Algorithmique Numérique Distribuée

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