Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
remove useless parameters in (simplified) instr
[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     int src_traced = comm->group()->index(src);
162
163     instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
164     extra->type = TRACING_IRECV;
165     extra->src = src_traced;
166     extra->dst = rank;
167     int known=0;
168     extra->datatype1 = encode_datatype(datatype, &known);
169     int dt_size_send = 1;
170     if(known==0)
171       dt_size_send = datatype->size();
172     extra->send_size = count*dt_size_send;
173     TRACE_smpi_ptp_in(rank, __FUNCTION__, extra);
174
175     *request = simgrid::smpi::Request::irecv(buf, count, datatype, src, tag, comm);
176     retval = MPI_SUCCESS;
177
178     TRACE_smpi_ptp_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 dst_traced = comm->group()->index(dst);
211     instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
212     extra->type = TRACING_ISEND;
213     extra->src = rank;
214     extra->dst = dst_traced;
215     int known=0;
216     extra->datatype1 = encode_datatype(datatype, &known);
217     int dt_size_send = 1;
218     if(known==0)
219       dt_size_send = datatype->size();
220     extra->send_size = count*dt_size_send;
221     TRACE_smpi_ptp_in(rank, __FUNCTION__, extra);
222     TRACE_smpi_send(rank, rank, dst_traced, tag, count*datatype->size());
223
224     *request = simgrid::smpi::Request::isend(buf, count, datatype, dst, tag, comm);
225     retval = MPI_SUCCESS;
226
227     TRACE_smpi_ptp_out(rank);
228   }
229
230   smpi_bench_begin();
231   if (retval != MPI_SUCCESS && request!=nullptr)
232     *request = MPI_REQUEST_NULL;
233   return retval;
234 }
235
236 int PMPI_Issend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request)
237 {
238   int retval = 0;
239
240   smpi_bench_end();
241   if (request == nullptr) {
242     retval = MPI_ERR_ARG;
243   } else if (comm == MPI_COMM_NULL) {
244     retval = MPI_ERR_COMM;
245   } else if (dst == MPI_PROC_NULL) {
246     *request = MPI_REQUEST_NULL;
247     retval = MPI_SUCCESS;
248   } else if (dst >= comm->group()->size() || dst <0){
249     retval = MPI_ERR_RANK;
250   } else if ((count < 0)|| (buf==nullptr && count > 0)) {
251     retval = MPI_ERR_COUNT;
252   } else if (not datatype->is_valid()) {
253     retval = MPI_ERR_TYPE;
254   } else if(tag<0 && tag !=  MPI_ANY_TAG){
255     retval = MPI_ERR_TAG;
256   } else {
257     int rank               = smpi_process()->index();
258     int dst_traced = comm->group()->index(dst);
259     instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
260     extra->type = TRACING_ISSEND;
261     extra->src = rank;
262     extra->dst = dst_traced;
263     int known=0;
264     extra->datatype1 = encode_datatype(datatype, &known);
265     int dt_size_send = 1;
266     if(known==0)
267       dt_size_send = datatype->size();
268     extra->send_size = count*dt_size_send;
269     TRACE_smpi_ptp_in(rank, __FUNCTION__, extra);
270     TRACE_smpi_send(rank, rank, dst_traced, tag, count*datatype->size());
271
272     *request = simgrid::smpi::Request::issend(buf, count, datatype, dst, tag, comm);
273     retval = MPI_SUCCESS;
274
275     TRACE_smpi_ptp_out(rank);
276   }
277
278   smpi_bench_begin();
279   if (retval != MPI_SUCCESS && request!=nullptr)
280     *request = MPI_REQUEST_NULL;
281   return retval;
282 }
283
284 int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Status * status)
285 {
286   int retval = 0;
287
288   smpi_bench_end();
289   if (comm == MPI_COMM_NULL) {
290     retval = MPI_ERR_COMM;
291   } else if (src == MPI_PROC_NULL) {
292     simgrid::smpi::Status::empty(status);
293     status->MPI_SOURCE = MPI_PROC_NULL;
294     retval = MPI_SUCCESS;
295   } else if (src!=MPI_ANY_SOURCE && (src >= comm->group()->size() || src <0)){
296     retval = MPI_ERR_RANK;
297   } else if ((count < 0) || (buf==nullptr && count > 0)) {
298     retval = MPI_ERR_COUNT;
299   } else if (not datatype->is_valid()) {
300     retval = MPI_ERR_TYPE;
301   } else if(tag<0 && tag !=  MPI_ANY_TAG){
302     retval = MPI_ERR_TAG;
303   } else {
304     int rank               = smpi_process()->index();
305     int src_traced         = comm->group()->index(src);
306     instr_extra_data extra = xbt_new0(s_instr_extra_data_t, 1);
307     extra->type            = TRACING_RECV;
308     extra->src             = src_traced;
309     extra->dst             = rank;
310     int known              = 0;
311     extra->datatype1       = encode_datatype(datatype, &known);
312     int dt_size_send       = 1;
313     if (known == 0)
314       dt_size_send   = datatype->size();
315     extra->send_size = count * dt_size_send;
316     TRACE_smpi_ptp_in(rank, __FUNCTION__, extra);
317
318     simgrid::smpi::Request::recv(buf, count, datatype, src, tag, comm, status);
319     retval = MPI_SUCCESS;
320
321     // the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
322     if (status != MPI_STATUS_IGNORE) {
323       src_traced = comm->group()->index(status->MPI_SOURCE);
324       if (not TRACE_smpi_view_internals()) {
325         TRACE_smpi_recv(src_traced, rank, tag);
326       }
327     }
328     TRACE_smpi_ptp_out(rank);
329   }
330
331   smpi_bench_begin();
332   return retval;
333 }
334
335 int PMPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm)
336 {
337   int retval = 0;
338
339   smpi_bench_end();
340
341   if (comm == MPI_COMM_NULL) {
342     retval = MPI_ERR_COMM;
343   } else if (dst == MPI_PROC_NULL) {
344     retval = MPI_SUCCESS;
345   } else if (dst >= comm->group()->size() || dst <0){
346     retval = MPI_ERR_RANK;
347   } else if ((count < 0) || (buf == nullptr && count > 0)) {
348     retval = MPI_ERR_COUNT;
349   } else if (not datatype->is_valid()) {
350     retval = MPI_ERR_TYPE;
351   } else if(tag < 0 && tag !=  MPI_ANY_TAG){
352     retval = MPI_ERR_TAG;
353   } else {
354     int rank               = smpi_process()->index();
355     int dst_traced         = comm->group()->index(dst);
356     instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
357     extra->type            = TRACING_SEND;
358     extra->src             = rank;
359     extra->dst             = dst_traced;
360     int known              = 0;
361     extra->datatype1       = encode_datatype(datatype, &known);
362     int dt_size_send       = 1;
363     if (known == 0) {
364       dt_size_send = datatype->size();
365     }
366     extra->send_size = count*dt_size_send;
367     TRACE_smpi_ptp_in(rank, __FUNCTION__, extra);
368     if (not TRACE_smpi_view_internals()) {
369       TRACE_smpi_send(rank, rank, dst_traced, tag,count*datatype->size());
370     }
371
372     simgrid::smpi::Request::send(buf, count, datatype, dst, tag, comm);
373     retval = MPI_SUCCESS;
374
375     TRACE_smpi_ptp_out(rank);
376   }
377
378   smpi_bench_begin();
379   return retval;
380 }
381
382 int PMPI_Ssend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm) {
383   int retval = 0;
384
385   smpi_bench_end();
386
387   if (comm == MPI_COMM_NULL) {
388     retval = MPI_ERR_COMM;
389   } else if (dst == MPI_PROC_NULL) {
390     retval = MPI_SUCCESS;
391   } else if (dst >= comm->group()->size() || dst <0){
392     retval = MPI_ERR_RANK;
393   } else if ((count < 0) || (buf==nullptr && count > 0)) {
394     retval = MPI_ERR_COUNT;
395   } else if (not datatype->is_valid()) {
396     retval = MPI_ERR_TYPE;
397   } else if(tag<0 && tag !=  MPI_ANY_TAG){
398     retval = MPI_ERR_TAG;
399   } else {
400     int rank               = smpi_process()->index();
401     int dst_traced         = comm->group()->index(dst);
402     instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
403     extra->type            = TRACING_SSEND;
404     extra->src             = rank;
405     extra->dst             = dst_traced;
406     int known              = 0;
407     extra->datatype1       = encode_datatype(datatype, &known);
408     int dt_size_send       = 1;
409     if(known == 0) {
410       dt_size_send = datatype->size();
411     }
412     extra->send_size = count*dt_size_send;
413     TRACE_smpi_ptp_in(rank, __FUNCTION__, extra);
414     TRACE_smpi_send(rank, rank, dst_traced, tag,count*datatype->size());
415
416     simgrid::smpi::Request::ssend(buf, count, datatype, dst, tag, comm);
417     retval = MPI_SUCCESS;
418
419     TRACE_smpi_ptp_out(rank);
420   }
421
422   smpi_bench_begin();
423   return retval;
424 }
425
426 int PMPI_Sendrecv(void* sendbuf, int sendcount, MPI_Datatype sendtype, int dst, int sendtag, void* recvbuf,
427                   int recvcount, MPI_Datatype recvtype, int src, int recvtag, MPI_Comm comm, MPI_Status* status)
428 {
429   int retval = 0;
430
431   smpi_bench_end();
432
433   if (comm == MPI_COMM_NULL) {
434     retval = MPI_ERR_COMM;
435   } else if (not sendtype->is_valid() || not recvtype->is_valid()) {
436     retval = MPI_ERR_TYPE;
437   } else if (src == MPI_PROC_NULL || dst == MPI_PROC_NULL) {
438     simgrid::smpi::Status::empty(status);
439     status->MPI_SOURCE = MPI_PROC_NULL;
440     retval             = MPI_SUCCESS;
441   }else if (dst >= comm->group()->size() || dst <0 ||
442       (src!=MPI_ANY_SOURCE && (src >= comm->group()->size() || src <0))){
443     retval = MPI_ERR_RANK;
444   } else if ((sendcount < 0 || recvcount<0) ||
445       (sendbuf==nullptr && sendcount > 0) || (recvbuf==nullptr && recvcount>0)) {
446     retval = MPI_ERR_COUNT;
447   } else if((sendtag<0 && sendtag !=  MPI_ANY_TAG)||(recvtag<0 && recvtag != MPI_ANY_TAG)){
448     retval = MPI_ERR_TAG;
449   } else {
450
451     int rank               = smpi_process()->index();
452     int dst_traced         = comm->group()->index(dst);
453     int src_traced         = comm->group()->index(src);
454     instr_extra_data extra = xbt_new0(s_instr_extra_data_t, 1);
455     extra->type            = TRACING_SENDRECV;
456     extra->src             = src_traced;
457     extra->dst             = dst_traced;
458     int known              = 0;
459     extra->datatype1       = encode_datatype(sendtype, &known);
460     int dt_size_send       = 1;
461     if (known == 0)
462       dt_size_send   = sendtype->size();
463     extra->send_size = sendcount * dt_size_send;
464     extra->datatype2 = encode_datatype(recvtype, &known);
465     int dt_size_recv = 1;
466     if (known == 0)
467       dt_size_recv   = recvtype->size();
468     extra->recv_size = recvcount * dt_size_recv;
469
470     TRACE_smpi_ptp_in(rank, __FUNCTION__, extra);
471     TRACE_smpi_send(rank, rank, dst_traced, sendtag, sendcount * sendtype->size());
472
473     simgrid::smpi::Request::sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf, recvcount, recvtype, src,
474                                      recvtag, comm, status);
475     retval = MPI_SUCCESS;
476
477     TRACE_smpi_ptp_out(rank);
478     TRACE_smpi_recv(src_traced, rank, recvtag);
479   }
480
481   smpi_bench_begin();
482   return retval;
483 }
484
485 int PMPI_Sendrecv_replace(void* buf, int count, MPI_Datatype datatype, int dst, int sendtag, int src, int recvtag,
486                           MPI_Comm comm, MPI_Status* status)
487 {
488   int retval = 0;
489   if (not datatype->is_valid()) {
490     return MPI_ERR_TYPE;
491   } else if (count < 0) {
492     return MPI_ERR_COUNT;
493   } else {
494     int size = datatype->get_extent() * count;
495     void* recvbuf = xbt_new0(char, size);
496     retval = MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count, datatype, src, recvtag, comm, status);
497     if(retval==MPI_SUCCESS){
498       simgrid::smpi::Datatype::copy(recvbuf, count, datatype, buf, count, datatype);
499     }
500     xbt_free(recvbuf);
501
502   }
503   return retval;
504 }
505
506 int PMPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
507 {
508   int retval = 0;
509   smpi_bench_end();
510   if (request == nullptr || flag == nullptr) {
511     retval = MPI_ERR_ARG;
512   } else if (*request == MPI_REQUEST_NULL) {
513     *flag= true;
514     simgrid::smpi::Status::empty(status);
515     retval = MPI_SUCCESS;
516   } else {
517     int rank = ((*request)->comm() != MPI_COMM_NULL) ? smpi_process()->index() : -1;
518
519     instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
520     extra->type = TRACING_TEST;
521     TRACE_smpi_testing_in(rank, extra);
522
523     *flag = simgrid::smpi::Request::test(request,status);
524
525     TRACE_smpi_testing_out(rank);
526     retval = MPI_SUCCESS;
527   }
528   smpi_bench_begin();
529   return retval;
530 }
531
532 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag, MPI_Status * status)
533 {
534   int retval = 0;
535
536   smpi_bench_end();
537   if (index == nullptr || flag == nullptr) {
538     retval = MPI_ERR_ARG;
539   } else {
540     *flag = simgrid::smpi::Request::testany(count, requests, index, status);
541     retval = MPI_SUCCESS;
542   }
543   smpi_bench_begin();
544   return retval;
545 }
546
547 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses)
548 {
549   int retval = 0;
550
551   smpi_bench_end();
552   if (flag == nullptr) {
553     retval = MPI_ERR_ARG;
554   } else {
555     *flag = simgrid::smpi::Request::testall(count, requests, statuses);
556     retval = MPI_SUCCESS;
557   }
558   smpi_bench_begin();
559   return retval;
560 }
561
562 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
563   int retval = 0;
564   smpi_bench_end();
565
566   if (status == nullptr) {
567     retval = MPI_ERR_ARG;
568   } else if (comm == MPI_COMM_NULL) {
569     retval = MPI_ERR_COMM;
570   } else if (source == MPI_PROC_NULL) {
571     simgrid::smpi::Status::empty(status);
572     status->MPI_SOURCE = MPI_PROC_NULL;
573     retval = MPI_SUCCESS;
574   } else {
575     simgrid::smpi::Request::probe(source, tag, comm, status);
576     retval = MPI_SUCCESS;
577   }
578   smpi_bench_begin();
579   return retval;
580 }
581
582 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
583   int retval = 0;
584   smpi_bench_end();
585
586   if (flag == nullptr) {
587     retval = MPI_ERR_ARG;
588   } else if (comm == MPI_COMM_NULL) {
589     retval = MPI_ERR_COMM;
590   } else if (source == MPI_PROC_NULL) {
591     *flag=true;
592     simgrid::smpi::Status::empty(status);
593     status->MPI_SOURCE = MPI_PROC_NULL;
594     retval = MPI_SUCCESS;
595   } else {
596     simgrid::smpi::Request::iprobe(source, tag, comm, flag, status);
597     retval = MPI_SUCCESS;
598   }
599   smpi_bench_begin();
600   return retval;
601 }
602
603 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
604 {
605   int retval = 0;
606
607   smpi_bench_end();
608
609   simgrid::smpi::Status::empty(status);
610
611   if (request == nullptr) {
612     retval = MPI_ERR_ARG;
613   } else if (*request == MPI_REQUEST_NULL) {
614     retval = MPI_SUCCESS;
615   } else {
616     int rank = (*request)->comm() != MPI_COMM_NULL ? smpi_process()->index() : -1;
617
618     int src_traced = (*request)->src();
619     int dst_traced = (*request)->dst();
620     int tag_traced= (*request)->tag();
621     MPI_Comm comm = (*request)->comm();
622     int is_wait_for_receive = ((*request)->flags() & RECV);
623     instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
624     extra->type = TRACING_WAIT;
625     TRACE_smpi_ptp_in(rank, __FUNCTION__, extra);
626
627     simgrid::smpi::Request::wait(request, status);
628     retval = MPI_SUCCESS;
629
630     //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
631     TRACE_smpi_ptp_out(rank);
632     if (is_wait_for_receive) {
633       if(src_traced==MPI_ANY_SOURCE)
634         src_traced = (status != MPI_STATUS_IGNORE) ? comm->group()->rank(status->MPI_SOURCE) : src_traced;
635       TRACE_smpi_recv(src_traced, dst_traced, tag_traced);
636     }
637   }
638
639   smpi_bench_begin();
640   return retval;
641 }
642
643 int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * status)
644 {
645   if (index == nullptr)
646     return MPI_ERR_ARG;
647
648   if (count <= 0)
649     return MPI_SUCCESS;
650
651   smpi_bench_end();
652   //save requests information for tracing
653   struct savedvalstype {
654     int src;
655     int dst;
656     int recv;
657     int tag;
658     MPI_Comm comm;
659   };
660   savedvalstype* savedvals = xbt_new0(savedvalstype, count);
661
662   for (int i = 0; i < count; i++) {
663     MPI_Request req = requests[i];      //already received requests are no longer valid
664     if (req) {
665       savedvals[i]=(savedvalstype){req->src(), req->dst(), (req->flags() & RECV), req->tag(), req->comm()};
666     }
667   }
668   int rank_traced = smpi_process()->index();
669   instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
670   extra->type = TRACING_WAITANY;
671   extra->send_size=count;
672   TRACE_smpi_ptp_in(rank_traced, __FUNCTION__,extra);
673
674   *index = simgrid::smpi::Request::waitany(count, requests, status);
675
676   if(*index!=MPI_UNDEFINED){
677     int src_traced = savedvals[*index].src;
678     //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
679     int dst_traced = savedvals[*index].dst;
680     int is_wait_for_receive = savedvals[*index].recv;
681     if (is_wait_for_receive) {
682       if(savedvals[*index].src==MPI_ANY_SOURCE)
683         src_traced = (status != MPI_STATUSES_IGNORE) ? savedvals[*index].comm->group()->rank(status->MPI_SOURCE)
684                                                      : savedvals[*index].src;
685       TRACE_smpi_recv(src_traced, dst_traced, savedvals[*index].tag);
686     }
687     TRACE_smpi_ptp_out(rank_traced);
688   }
689   xbt_free(savedvals);
690
691   smpi_bench_begin();
692   return MPI_SUCCESS;
693 }
694
695 int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
696 {
697   smpi_bench_end();
698   //save information from requests
699   struct savedvalstype {
700     int src;
701     int dst;
702     int recv;
703     int tag;
704     int valid;
705     MPI_Comm comm;
706   };
707   savedvalstype* savedvals=xbt_new0(savedvalstype, count);
708
709   for (int i = 0; i < count; i++) {
710     MPI_Request req = requests[i];
711     if(req!=MPI_REQUEST_NULL){
712       savedvals[i]=(savedvalstype){req->src(), req->dst(), (req->flags() & RECV), req->tag(), 1, req->comm()};
713     }else{
714       savedvals[i].valid=0;
715     }
716   }
717   int rank_traced = smpi_process()->index();
718   instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
719   extra->type = TRACING_WAITALL;
720   extra->send_size=count;
721   TRACE_smpi_ptp_in(rank_traced, __FUNCTION__,extra);
722
723   int retval = simgrid::smpi::Request::waitall(count, requests, status);
724
725   for (int i = 0; i < count; i++) {
726     if(savedvals[i].valid){
727       // the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
728       int src_traced = savedvals[i].src;
729       int dst_traced = savedvals[i].dst;
730       int is_wait_for_receive = savedvals[i].recv;
731       if (is_wait_for_receive) {
732         if(src_traced==MPI_ANY_SOURCE)
733           src_traced = (status != MPI_STATUSES_IGNORE) ? savedvals[i].comm->group()->rank(status[i].MPI_SOURCE)
734                                                        : savedvals[i].src;
735         TRACE_smpi_recv(src_traced, dst_traced,savedvals[i].tag);
736       }
737     }
738   }
739   TRACE_smpi_ptp_out(rank_traced);
740   xbt_free(savedvals);
741
742   smpi_bench_begin();
743   return retval;
744 }
745
746 int PMPI_Waitsome(int incount, MPI_Request requests[], int *outcount, int *indices, MPI_Status status[])
747 {
748   int retval = 0;
749
750   smpi_bench_end();
751   if (outcount == nullptr) {
752     retval = MPI_ERR_ARG;
753   } else {
754     *outcount = simgrid::smpi::Request::waitsome(incount, requests, indices, status);
755     retval = MPI_SUCCESS;
756   }
757   smpi_bench_begin();
758   return retval;
759 }
760
761 int PMPI_Testsome(int incount, MPI_Request requests[], int* outcount, int* indices, MPI_Status status[])
762 {
763   int retval = 0;
764
765   smpi_bench_end();
766   if (outcount == nullptr) {
767     retval = MPI_ERR_ARG;
768   } else {
769     *outcount = simgrid::smpi::Request::testsome(incount, requests, indices, status);
770     retval    = MPI_SUCCESS;
771   }
772   smpi_bench_begin();
773   return retval;
774 }
775
776 MPI_Request PMPI_Request_f2c(MPI_Fint request){
777   return static_cast<MPI_Request>(simgrid::smpi::Request::f2c(request));
778 }
779
780 MPI_Fint PMPI_Request_c2f(MPI_Request request) {
781   return request->c2f();
782 }
783
784 }