Logo AND Algorithmique Numérique Distribuée

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