Logo AND Algorithmique Numérique Distribuée

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