Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
[SMPI] Remove now unused datastructure in PMPI_Waitany
[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_replayable() ? 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_replayable() ? 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_replayable() ? 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_replayable() ? 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_replayable() ? 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_replayable() ? 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     std::vector<int>* dst_hack = new std::vector<int>;
412     std::vector<int>* src_hack = new std::vector<int>;
413     dst_hack->push_back(dst_traced);
414     src_hack->push_back(src_traced);
415     TRACE_smpi_comm_in(rank, __FUNCTION__,
416                        new simgrid::instr::VarCollTIData(
417                            "sendRecv", -1, sendtype->is_replayable() ? sendcount : sendcount * sendtype->size(), dst_hack,
418                            recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(), src_hack,
419                            encode_datatype(sendtype), encode_datatype(recvtype)));
420
421     TRACE_smpi_send(rank, rank, dst_traced, sendtag, sendcount * sendtype->size());
422
423     simgrid::smpi::Request::sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf, recvcount, recvtype, src,
424                                      recvtag, comm, status);
425     retval = MPI_SUCCESS;
426
427     TRACE_smpi_recv(src_traced, rank, recvtag);
428     TRACE_smpi_comm_out(rank);
429   }
430
431   smpi_bench_begin();
432   return retval;
433 }
434
435 int PMPI_Sendrecv_replace(void* buf, int count, MPI_Datatype datatype, int dst, int sendtag, int src, int recvtag,
436                           MPI_Comm comm, MPI_Status* status)
437 {
438   int retval = 0;
439   if (not datatype->is_valid()) {
440     return MPI_ERR_TYPE;
441   } else if (count < 0) {
442     return MPI_ERR_COUNT;
443   } else {
444     int size = datatype->get_extent() * count;
445     void* recvbuf = xbt_new0(char, size);
446     retval = MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count, datatype, src, recvtag, comm, status);
447     if(retval==MPI_SUCCESS){
448       simgrid::smpi::Datatype::copy(recvbuf, count, datatype, buf, count, datatype);
449     }
450     xbt_free(recvbuf);
451
452   }
453   return retval;
454 }
455
456 int PMPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
457 {
458   int retval = 0;
459   smpi_bench_end();
460   if (request == nullptr || flag == nullptr) {
461     retval = MPI_ERR_ARG;
462   } else if (*request == MPI_REQUEST_NULL) {
463     *flag= true;
464     simgrid::smpi::Status::empty(status);
465     retval = MPI_SUCCESS;
466   } else {
467     int rank = ((*request)->comm() != MPI_COMM_NULL) ? smpi_process()->index() : -1;
468
469     TRACE_smpi_testing_in(rank);
470
471     *flag = simgrid::smpi::Request::test(request,status);
472
473     TRACE_smpi_testing_out(rank);
474     retval = MPI_SUCCESS;
475   }
476   smpi_bench_begin();
477   return retval;
478 }
479
480 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag, MPI_Status * status)
481 {
482   int retval = 0;
483
484   smpi_bench_end();
485   if (index == nullptr || flag == nullptr) {
486     retval = MPI_ERR_ARG;
487   } else {
488     *flag = simgrid::smpi::Request::testany(count, requests, index, status);
489     retval = MPI_SUCCESS;
490   }
491   smpi_bench_begin();
492   return retval;
493 }
494
495 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses)
496 {
497   int retval = 0;
498
499   smpi_bench_end();
500   if (flag == nullptr) {
501     retval = MPI_ERR_ARG;
502   } else {
503     *flag = simgrid::smpi::Request::testall(count, requests, statuses);
504     retval = MPI_SUCCESS;
505   }
506   smpi_bench_begin();
507   return retval;
508 }
509
510 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
511   int retval = 0;
512   smpi_bench_end();
513
514   if (status == nullptr) {
515     retval = MPI_ERR_ARG;
516   } else if (comm == MPI_COMM_NULL) {
517     retval = MPI_ERR_COMM;
518   } else if (source == MPI_PROC_NULL) {
519     simgrid::smpi::Status::empty(status);
520     status->MPI_SOURCE = MPI_PROC_NULL;
521     retval = MPI_SUCCESS;
522   } else {
523     simgrid::smpi::Request::probe(source, tag, comm, status);
524     retval = MPI_SUCCESS;
525   }
526   smpi_bench_begin();
527   return retval;
528 }
529
530 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
531   int retval = 0;
532   smpi_bench_end();
533
534   if (flag == nullptr) {
535     retval = MPI_ERR_ARG;
536   } else if (comm == MPI_COMM_NULL) {
537     retval = MPI_ERR_COMM;
538   } else if (source == MPI_PROC_NULL) {
539     *flag=true;
540     simgrid::smpi::Status::empty(status);
541     status->MPI_SOURCE = MPI_PROC_NULL;
542     retval = MPI_SUCCESS;
543   } else {
544     simgrid::smpi::Request::iprobe(source, tag, comm, flag, status);
545     retval = MPI_SUCCESS;
546   }
547   smpi_bench_begin();
548   return retval;
549 }
550
551 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
552 {
553   int retval = 0;
554
555   smpi_bench_end();
556
557   simgrid::smpi::Status::empty(status);
558
559   if (request == nullptr) {
560     retval = MPI_ERR_ARG;
561   } else if (*request == MPI_REQUEST_NULL) {
562     retval = MPI_SUCCESS;
563   } else {
564     int rank = (*request)->comm() != MPI_COMM_NULL ? smpi_process()->index() : -1;
565
566     int src_traced = (*request)->src();
567     int dst_traced = (*request)->dst();
568     int tag_traced= (*request)->tag();
569     MPI_Comm comm = (*request)->comm();
570     int is_wait_for_receive = ((*request)->flags() & RECV);
571     TRACE_smpi_comm_in(rank, __FUNCTION__, new simgrid::instr::NoOpTIData("wait"));
572
573     simgrid::smpi::Request::wait(request, status);
574     retval = MPI_SUCCESS;
575
576     //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
577     TRACE_smpi_comm_out(rank);
578     if (is_wait_for_receive) {
579       if(src_traced==MPI_ANY_SOURCE)
580         src_traced = (status != MPI_STATUS_IGNORE) ? comm->group()->rank(status->MPI_SOURCE) : src_traced;
581       TRACE_smpi_recv(src_traced, dst_traced, tag_traced);
582     }
583   }
584
585   smpi_bench_begin();
586   return retval;
587 }
588
589 int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * status)
590 {
591   if (index == nullptr)
592     return MPI_ERR_ARG;
593
594   if (count <= 0)
595     return MPI_SUCCESS;
596
597   smpi_bench_end();
598
599   int rank_traced = smpi_process()->index();
600   TRACE_smpi_comm_in(rank_traced, __FUNCTION__, new simgrid::instr::CpuTIData("waitAny", static_cast<double>(count)));
601
602   *index = simgrid::smpi::Request::waitany(count, requests, status);
603
604   if(*index!=MPI_UNDEFINED){
605     MPI_Request req = requests[*index];
606     if (req != nullptr) { // Requests that were already received will be a nullptr
607       int src_traced = req->src();
608       // the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
609       int dst_traced          = req->dst();
610       int is_wait_for_receive = req->flags() & RECV;
611       if (is_wait_for_receive) {
612         if (req->src() == MPI_ANY_SOURCE)
613           src_traced = (status != MPI_STATUSES_IGNORE) ? req->comm()->group()->rank(status->MPI_SOURCE) : req->src();
614         TRACE_smpi_recv(src_traced, dst_traced, req->tag());
615       }
616     }
617     TRACE_smpi_comm_out(rank_traced);
618   }
619
620   smpi_bench_begin();
621   return MPI_SUCCESS;
622 }
623
624 int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
625 {
626   smpi_bench_end();
627   //save information from requests
628   struct savedvalstype {
629     int src;
630     int dst;
631     int recv;
632     int tag;
633     int valid;
634     MPI_Comm comm;
635   };
636   savedvalstype* savedvals=xbt_new0(savedvalstype, count);
637
638   for (int i = 0; i < count; i++) {
639     MPI_Request req = requests[i];
640     if(req!=MPI_REQUEST_NULL){
641       savedvals[i]=(savedvalstype){req->src(), req->dst(), (req->flags() & RECV), req->tag(), 1, req->comm()};
642     }else{
643       savedvals[i].valid=0;
644     }
645   }
646   int rank_traced = smpi_process()->index();
647   TRACE_smpi_comm_in(rank_traced, __FUNCTION__, new simgrid::instr::CpuTIData("waitAll", static_cast<double>(count)));
648
649   int retval = simgrid::smpi::Request::waitall(count, requests, status);
650
651   for (int i = 0; i < count; i++) {
652     if(savedvals[i].valid){
653       // the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
654       int src_traced = savedvals[i].src;
655       int dst_traced = savedvals[i].dst;
656       int is_wait_for_receive = savedvals[i].recv;
657       if (is_wait_for_receive) {
658         if(src_traced==MPI_ANY_SOURCE)
659           src_traced = (status != MPI_STATUSES_IGNORE) ? savedvals[i].comm->group()->rank(status[i].MPI_SOURCE)
660                                                        : savedvals[i].src;
661         TRACE_smpi_recv(src_traced, dst_traced,savedvals[i].tag);
662       }
663     }
664   }
665   TRACE_smpi_comm_out(rank_traced);
666   xbt_free(savedvals);
667
668   smpi_bench_begin();
669   return retval;
670 }
671
672 int PMPI_Waitsome(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::waitsome(incount, requests, indices, status);
681     retval = MPI_SUCCESS;
682   }
683   smpi_bench_begin();
684   return retval;
685 }
686
687 int PMPI_Testsome(int incount, MPI_Request requests[], int* outcount, int* indices, MPI_Status status[])
688 {
689   int retval = 0;
690
691   smpi_bench_end();
692   if (outcount == nullptr) {
693     retval = MPI_ERR_ARG;
694   } else {
695     *outcount = simgrid::smpi::Request::testsome(incount, requests, indices, status);
696     retval    = MPI_SUCCESS;
697   }
698   smpi_bench_begin();
699   return retval;
700 }
701
702 MPI_Request PMPI_Request_f2c(MPI_Fint request){
703   return static_cast<MPI_Request>(simgrid::smpi::Request::f2c(request));
704 }
705
706 MPI_Fint PMPI_Request_c2f(MPI_Request request) {
707   return request->c2f();
708 }
709
710 }