Logo AND Algorithmique Numérique Distribuée

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