Logo AND Algorithmique Numérique Distribuée

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