Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
obey our coding standards and cosmetics
[simgrid.git] / src / smpi / bindings / smpi_pmpi_request.cpp
1 /* Copyright (c) 2007-2018. The SimGrid Team. All rights reserved.          */
2
3 /* This program is free software; you can redistribute it and/or modify it
4  * under the terms of the license (GNU LGPL) which comes with this package. */
5
6 #include "private.hpp"
7 #include "smpi_comm.hpp"
8 #include "smpi_datatype.hpp"
9 #include "smpi_process.hpp"
10 #include "smpi_request.hpp"
11
12 XBT_LOG_EXTERNAL_DEFAULT_CATEGORY(smpi_pmpi);
13
14 static int getPid(MPI_Comm, int);
15 static int getPid(MPI_Comm comm, int id)
16 {
17   simgrid::s4u::ActorPtr actor = comm->group()->actor(id);
18   return (actor == nullptr) ? MPI_UNDEFINED : actor->get_pid();
19 }
20
21 /* PMPI User level calls */
22
23 int PMPI_Send_init(void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request * request)
24 {
25   int retval = 0;
26
27   smpi_bench_end();
28   if (request == nullptr) {
29     retval = MPI_ERR_ARG;
30   } else if (comm == MPI_COMM_NULL) {
31     retval = MPI_ERR_COMM;
32   } else if (not datatype->is_valid()) {
33     retval = MPI_ERR_TYPE;
34   } else if (dst == MPI_PROC_NULL) {
35     retval = MPI_SUCCESS;
36   } else {
37     *request = simgrid::smpi::Request::send_init(buf, count, datatype, dst, tag, comm);
38     retval   = MPI_SUCCESS;
39   }
40   smpi_bench_begin();
41   if (retval != MPI_SUCCESS && request != nullptr)
42     *request = MPI_REQUEST_NULL;
43   return retval;
44 }
45
46 int PMPI_Recv_init(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Request * request)
47 {
48   int retval = 0;
49
50   smpi_bench_end();
51   if (request == nullptr) {
52     retval = MPI_ERR_ARG;
53   } else if (comm == MPI_COMM_NULL) {
54     retval = MPI_ERR_COMM;
55   } else if (not datatype->is_valid()) {
56     retval = MPI_ERR_TYPE;
57   } else if (src == MPI_PROC_NULL) {
58     retval = MPI_SUCCESS;
59   } else {
60     *request = simgrid::smpi::Request::recv_init(buf, count, datatype, src, tag, comm);
61     retval = MPI_SUCCESS;
62   }
63   smpi_bench_begin();
64   if (retval != MPI_SUCCESS && request != nullptr)
65     *request = MPI_REQUEST_NULL;
66   return retval;
67 }
68
69 int PMPI_Ssend_init(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request)
70 {
71   int retval = 0;
72
73   smpi_bench_end();
74   if (request == nullptr) {
75     retval = MPI_ERR_ARG;
76   } else if (comm == MPI_COMM_NULL) {
77     retval = MPI_ERR_COMM;
78   } else if (not datatype->is_valid()) {
79     retval = MPI_ERR_TYPE;
80   } else if (dst == MPI_PROC_NULL) {
81     retval = MPI_SUCCESS;
82   } else {
83     *request = simgrid::smpi::Request::ssend_init(buf, count, datatype, dst, tag, comm);
84     retval = MPI_SUCCESS;
85   }
86   smpi_bench_begin();
87   if (retval != MPI_SUCCESS && request != nullptr)
88     *request = MPI_REQUEST_NULL;
89   return retval;
90 }
91
92 /*
93  * This function starts a request returned by init functions such as
94  * MPI_Send_init(), MPI_Ssend_init (see above), and friends.
95  * They should already have performed sanity checks.
96  */
97 int PMPI_Start(MPI_Request * request)
98 {
99   int retval = 0;
100
101   smpi_bench_end();
102   if (request == nullptr || *request == MPI_REQUEST_NULL) {
103     retval = MPI_ERR_REQUEST;
104   } else {
105     MPI_Request req = *request;
106     int my_proc_id = (req->comm() != MPI_COMM_NULL) ? simgrid::s4u::this_actor::get_pid() : -1;
107     TRACE_smpi_comm_in(my_proc_id, __func__,
108                        new simgrid::instr::Pt2PtTIData("Start", req->dst(),
109                                                        req->size(),
110                                                        req->tag(), 
111                                                        simgrid::smpi::Datatype::encode(req->type())));
112     if (not TRACE_smpi_view_internals() && req->flags() & MPI_REQ_SEND)
113       TRACE_smpi_send(my_proc_id, my_proc_id, getPid(req->comm(), req->dst()), req->tag(), req->size());
114
115     req->start();
116
117     if (not TRACE_smpi_view_internals() && req->flags() & MPI_REQ_RECV)
118       TRACE_smpi_recv(getPid(req->comm(), req->src()), my_proc_id, req->tag());
119     retval = MPI_SUCCESS;
120     TRACE_smpi_comm_out(my_proc_id);
121   }
122   smpi_bench_begin();
123   return retval;
124 }
125
126 int PMPI_Startall(int count, MPI_Request * requests)
127 {
128   int retval;
129   smpi_bench_end();
130   if (requests == nullptr) {
131     retval = MPI_ERR_ARG;
132   } else {
133     retval = MPI_SUCCESS;
134     for (int i = 0; i < count; i++) {
135       if(requests[i] == MPI_REQUEST_NULL) {
136         retval = MPI_ERR_REQUEST;
137       }
138     }
139     if(retval != MPI_ERR_REQUEST) {
140       int my_proc_id = simgrid::s4u::this_actor::get_pid();
141       TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("Startall"));
142       MPI_Request req = MPI_REQUEST_NULL;
143       if (not TRACE_smpi_view_internals())
144         for (int i = 0; i < count; i++) {
145           req = requests[i];
146           if (req->flags() & MPI_REQ_SEND)
147             TRACE_smpi_send(my_proc_id, my_proc_id, getPid(req->comm(), req->dst()), req->tag(), req->size());
148         }
149
150       simgrid::smpi::Request::startall(count, requests);
151
152       if (not TRACE_smpi_view_internals())
153         for (int i = 0; i < count; i++) {
154           req = requests[i];
155           if (req->flags() & MPI_REQ_RECV)
156             TRACE_smpi_recv(getPid(req->comm(), req->src()), my_proc_id, req->tag());
157         }
158       TRACE_smpi_comm_out(my_proc_id);
159     }
160   }
161   smpi_bench_begin();
162   return retval;
163 }
164
165 int PMPI_Request_free(MPI_Request * request)
166 {
167   int retval = 0;
168
169   smpi_bench_end();
170   if (*request == MPI_REQUEST_NULL) {
171     retval = MPI_ERR_ARG;
172   } else {
173     simgrid::smpi::Request::unref(request);
174     retval = MPI_SUCCESS;
175   }
176   smpi_bench_begin();
177   return retval;
178 }
179
180 int PMPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Request * request)
181 {
182   int retval = 0;
183
184   smpi_bench_end();
185
186   if (request == nullptr) {
187     retval = MPI_ERR_ARG;
188   } else if (comm == MPI_COMM_NULL) {
189     retval = MPI_ERR_COMM;
190   } else if (src == MPI_PROC_NULL) {
191     *request = MPI_REQUEST_NULL;
192     retval = MPI_SUCCESS;
193   } else if (src!=MPI_ANY_SOURCE && (src >= comm->group()->size() || src <0)){
194     retval = MPI_ERR_RANK;
195   } else if ((count < 0) || (buf==nullptr && count > 0)) {
196     retval = MPI_ERR_COUNT;
197   } else if (not datatype->is_valid()) {
198     retval = MPI_ERR_TYPE;
199   } else if(tag<0 && tag !=  MPI_ANY_TAG){
200     retval = MPI_ERR_TAG;
201   } else {
202
203     int my_proc_id = simgrid::s4u::this_actor::get_pid();
204
205     TRACE_smpi_comm_in(my_proc_id, __func__,
206                        new simgrid::instr::Pt2PtTIData("irecv", src,
207                                                        datatype->is_replayable() ? count : count * datatype->size(),
208                                                        tag, simgrid::smpi::Datatype::encode(datatype)));
209
210     *request = simgrid::smpi::Request::irecv(buf, count, datatype, src, tag, comm);
211     retval = MPI_SUCCESS;
212
213     TRACE_smpi_comm_out(my_proc_id);
214   }
215
216   smpi_bench_begin();
217   if (retval != MPI_SUCCESS && request != nullptr)
218     *request = MPI_REQUEST_NULL;
219   return retval;
220 }
221
222
223 int PMPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request * request)
224 {
225   int retval = 0;
226
227   smpi_bench_end();
228   if (request == nullptr) {
229     retval = MPI_ERR_ARG;
230   } else if (comm == MPI_COMM_NULL) {
231     retval = MPI_ERR_COMM;
232   } else if (dst == MPI_PROC_NULL) {
233     *request = MPI_REQUEST_NULL;
234     retval = MPI_SUCCESS;
235   } else if (dst >= comm->group()->size() || dst <0){
236     retval = MPI_ERR_RANK;
237   } else if ((count < 0) || (buf==nullptr && count > 0)) {
238     retval = MPI_ERR_COUNT;
239   } else if (not datatype->is_valid()) {
240     retval = MPI_ERR_TYPE;
241   } else if(tag<0 && tag !=  MPI_ANY_TAG){
242     retval = MPI_ERR_TAG;
243   } else {
244     int my_proc_id = simgrid::s4u::this_actor::get_pid();
245     int trace_dst = getPid(comm, dst);
246     TRACE_smpi_comm_in(my_proc_id, __func__,
247                        new simgrid::instr::Pt2PtTIData("isend", dst,
248                                                        datatype->is_replayable() ? count : count * datatype->size(),
249                                                        tag, simgrid::smpi::Datatype::encode(datatype)));
250
251     TRACE_smpi_send(my_proc_id, my_proc_id, trace_dst, tag, count * datatype->size());
252
253     *request = simgrid::smpi::Request::isend(buf, count, datatype, dst, tag, comm);
254     retval = MPI_SUCCESS;
255
256     TRACE_smpi_comm_out(my_proc_id);
257   }
258
259   smpi_bench_begin();
260   if (retval != MPI_SUCCESS && request!=nullptr)
261     *request = MPI_REQUEST_NULL;
262   return retval;
263 }
264
265 int PMPI_Issend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request)
266 {
267   int retval = 0;
268
269   smpi_bench_end();
270   if (request == nullptr) {
271     retval = MPI_ERR_ARG;
272   } else if (comm == MPI_COMM_NULL) {
273     retval = MPI_ERR_COMM;
274   } else if (dst == MPI_PROC_NULL) {
275     *request = MPI_REQUEST_NULL;
276     retval = MPI_SUCCESS;
277   } else if (dst >= comm->group()->size() || dst <0){
278     retval = MPI_ERR_RANK;
279   } else if ((count < 0)|| (buf==nullptr && count > 0)) {
280     retval = MPI_ERR_COUNT;
281   } else if (not datatype->is_valid()) {
282     retval = MPI_ERR_TYPE;
283   } else if(tag<0 && tag !=  MPI_ANY_TAG){
284     retval = MPI_ERR_TAG;
285   } else {
286     int my_proc_id = simgrid::s4u::this_actor::get_pid();
287     int trace_dst = getPid(comm, dst);
288     TRACE_smpi_comm_in(my_proc_id, __func__,
289                        new simgrid::instr::Pt2PtTIData("ISsend", dst,
290                                                        datatype->is_replayable() ? count : count * datatype->size(),
291                                                        tag, simgrid::smpi::Datatype::encode(datatype)));
292     TRACE_smpi_send(my_proc_id, my_proc_id, trace_dst, tag, count * datatype->size());
293
294     *request = simgrid::smpi::Request::issend(buf, count, datatype, dst, tag, comm);
295     retval = MPI_SUCCESS;
296
297     TRACE_smpi_comm_out(my_proc_id);
298   }
299
300   smpi_bench_begin();
301   if (retval != MPI_SUCCESS && request!=nullptr)
302     *request = MPI_REQUEST_NULL;
303   return retval;
304 }
305
306 int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Status * status)
307 {
308   int retval = 0;
309
310   smpi_bench_end();
311   if (comm == MPI_COMM_NULL) {
312     retval = MPI_ERR_COMM;
313   } else if (src == MPI_PROC_NULL) {
314     if(status != MPI_STATUS_IGNORE){
315       simgrid::smpi::Status::empty(status);
316       status->MPI_SOURCE = MPI_PROC_NULL;
317     }
318     retval = MPI_SUCCESS;
319   } else if (src!=MPI_ANY_SOURCE && (src >= comm->group()->size() || src <0)){
320     retval = MPI_ERR_RANK;
321   } else if ((count < 0) || (buf==nullptr && count > 0)) {
322     retval = MPI_ERR_COUNT;
323   } else if (not datatype->is_valid()) {
324     retval = MPI_ERR_TYPE;
325   } else if(tag<0 && tag !=  MPI_ANY_TAG){
326     retval = MPI_ERR_TAG;
327   } else {
328     int my_proc_id = simgrid::s4u::this_actor::get_pid();
329     TRACE_smpi_comm_in(my_proc_id, __func__,
330                        new simgrid::instr::Pt2PtTIData("recv", src,
331                                                        datatype->is_replayable() ? count : count * datatype->size(),
332                                                        tag, simgrid::smpi::Datatype::encode(datatype)));
333
334     simgrid::smpi::Request::recv(buf, count, datatype, src, tag, comm, status);
335     retval = MPI_SUCCESS;
336
337     // the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
338     int src_traced=0;
339     if (status != MPI_STATUS_IGNORE) 
340       src_traced = getPid(comm, status->MPI_SOURCE);
341     else
342       src_traced = getPid(comm, src);
343     if (not TRACE_smpi_view_internals()) {
344       TRACE_smpi_recv(src_traced, my_proc_id, tag);
345     }
346     
347     TRACE_smpi_comm_out(my_proc_id);
348   }
349
350   smpi_bench_begin();
351   return retval;
352 }
353
354 int PMPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm)
355 {
356   int retval = 0;
357
358   smpi_bench_end();
359
360   if (comm == MPI_COMM_NULL) {
361     retval = MPI_ERR_COMM;
362   } else if (dst == MPI_PROC_NULL) {
363     retval = MPI_SUCCESS;
364   } else if (dst >= comm->group()->size() || dst <0){
365     retval = MPI_ERR_RANK;
366   } else if ((count < 0) || (buf == nullptr && count > 0)) {
367     retval = MPI_ERR_COUNT;
368   } else if (not datatype->is_valid()) {
369     retval = MPI_ERR_TYPE;
370   } else if(tag < 0 && tag !=  MPI_ANY_TAG){
371     retval = MPI_ERR_TAG;
372   } else {
373     int my_proc_id         = simgrid::s4u::this_actor::get_pid();
374     int dst_traced         = getPid(comm, dst);
375     TRACE_smpi_comm_in(my_proc_id, __func__,
376                        new simgrid::instr::Pt2PtTIData("send", dst,
377                                                        datatype->is_replayable() ? count : count * datatype->size(),
378                                                        tag, simgrid::smpi::Datatype::encode(datatype)));
379     if (not TRACE_smpi_view_internals()) {
380       TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, tag, count * datatype->size());
381     }
382
383     simgrid::smpi::Request::send(buf, count, datatype, dst, tag, comm);
384     retval = MPI_SUCCESS;
385
386     TRACE_smpi_comm_out(my_proc_id);
387   }
388
389   smpi_bench_begin();
390   return retval;
391 }
392
393 int PMPI_Ssend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm) {
394   int retval = 0;
395
396   smpi_bench_end();
397
398   if (comm == MPI_COMM_NULL) {
399     retval = MPI_ERR_COMM;
400   } else if (dst == MPI_PROC_NULL) {
401     retval = MPI_SUCCESS;
402   } else if (dst >= comm->group()->size() || dst <0){
403     retval = MPI_ERR_RANK;
404   } else if ((count < 0) || (buf==nullptr && count > 0)) {
405     retval = MPI_ERR_COUNT;
406   } else if (not datatype->is_valid()) {
407     retval = MPI_ERR_TYPE;
408   } else if(tag<0 && tag !=  MPI_ANY_TAG){
409     retval = MPI_ERR_TAG;
410   } else {
411     int my_proc_id         = simgrid::s4u::this_actor::get_pid();
412     int dst_traced         = getPid(comm, dst);
413     TRACE_smpi_comm_in(my_proc_id, __func__,
414                        new simgrid::instr::Pt2PtTIData("Ssend", dst,
415                                                        datatype->is_replayable() ? count : count * datatype->size(),
416                                                        tag, simgrid::smpi::Datatype::encode(datatype)));
417     TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, tag, count * datatype->size());
418
419     simgrid::smpi::Request::ssend(buf, count, datatype, dst, tag, comm);
420     retval = MPI_SUCCESS;
421
422     TRACE_smpi_comm_out(my_proc_id);
423   }
424
425   smpi_bench_begin();
426   return retval;
427 }
428
429 int PMPI_Sendrecv(void* sendbuf, int sendcount, MPI_Datatype sendtype, int dst, int sendtag, void* recvbuf,
430                   int recvcount, MPI_Datatype recvtype, int src, int recvtag, MPI_Comm comm, MPI_Status* status)
431 {
432   int retval = 0;
433
434   smpi_bench_end();
435
436   if (comm == MPI_COMM_NULL) {
437     retval = MPI_ERR_COMM;
438   } else if (not sendtype->is_valid() || not recvtype->is_valid()) {
439     retval = MPI_ERR_TYPE;
440   } else if (src == MPI_PROC_NULL) {
441     if(status!=MPI_STATUS_IGNORE){
442       simgrid::smpi::Status::empty(status);
443       status->MPI_SOURCE = MPI_PROC_NULL;
444     }
445     if(dst != MPI_PROC_NULL)
446       simgrid::smpi::Request::send(sendbuf, sendcount, sendtype, dst, sendtag, comm);
447     retval = MPI_SUCCESS;
448   }else if (dst == MPI_PROC_NULL){
449     simgrid::smpi::Request::recv(recvbuf, recvcount, recvtype, src, recvtag, comm, status);
450     retval = MPI_SUCCESS;
451   }else if (dst >= comm->group()->size() || dst <0 ||
452       (src!=MPI_ANY_SOURCE && (src >= comm->group()->size() || src <0))){
453     retval = MPI_ERR_RANK;
454   } else if ((sendcount < 0 || recvcount<0) ||
455       (sendbuf==nullptr && sendcount > 0) || (recvbuf==nullptr && recvcount>0)) {
456     retval = MPI_ERR_COUNT;
457   } else if((sendtag<0 && sendtag !=  MPI_ANY_TAG)||(recvtag<0 && recvtag != MPI_ANY_TAG)){
458     retval = MPI_ERR_TAG;
459   } else {
460     int my_proc_id         = simgrid::s4u::this_actor::get_pid();
461     int dst_traced         = getPid(comm, dst);
462     int src_traced         = getPid(comm, src);
463
464     // FIXME: Hack the way to trace this one
465     std::vector<int>* dst_hack = new std::vector<int>;
466     std::vector<int>* src_hack = new std::vector<int>;
467     dst_hack->push_back(dst_traced);
468     src_hack->push_back(src_traced);
469     TRACE_smpi_comm_in(my_proc_id, __func__,
470                        new simgrid::instr::VarCollTIData(
471                            "sendRecv", -1, sendtype->is_replayable() ? sendcount : sendcount * sendtype->size(),
472                            dst_hack, recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(), src_hack,
473                            simgrid::smpi::Datatype::encode(sendtype), simgrid::smpi::Datatype::encode(recvtype)));
474
475     TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, sendtag, sendcount * sendtype->size());
476
477     simgrid::smpi::Request::sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf, recvcount, recvtype, src,
478                                      recvtag, comm, status);
479     retval = MPI_SUCCESS;
480
481     TRACE_smpi_recv(src_traced, my_proc_id, recvtag);
482     TRACE_smpi_comm_out(my_proc_id);
483   }
484
485   smpi_bench_begin();
486   return retval;
487 }
488
489 int PMPI_Sendrecv_replace(void* buf, int count, MPI_Datatype datatype, int dst, int sendtag, int src, int recvtag,
490                           MPI_Comm comm, MPI_Status* status)
491 {
492   int retval = 0;
493   if (not datatype->is_valid()) {
494     return MPI_ERR_TYPE;
495   } else if (count < 0) {
496     return MPI_ERR_COUNT;
497   } else {
498     int size = datatype->get_extent() * count;
499     void* recvbuf = xbt_new0(char, size);
500     retval = MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count, datatype, src, recvtag, comm, status);
501     if(retval==MPI_SUCCESS){
502       simgrid::smpi::Datatype::copy(recvbuf, count, datatype, buf, count, datatype);
503     }
504     xbt_free(recvbuf);
505
506   }
507   return retval;
508 }
509
510 int PMPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
511 {
512   int retval = 0;
513   smpi_bench_end();
514   if (request == nullptr || flag == nullptr) {
515     retval = MPI_ERR_ARG;
516   } else if (*request == MPI_REQUEST_NULL) {
517     if (status != MPI_STATUS_IGNORE){
518       *flag= true;
519       simgrid::smpi::Status::empty(status);
520     }
521     retval = MPI_SUCCESS;
522   } else {
523     int my_proc_id = ((*request)->comm() != MPI_COMM_NULL) ? simgrid::s4u::this_actor::get_pid() : -1;
524
525     TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("test"));
526     
527     *flag = simgrid::smpi::Request::test(request,status);
528
529     TRACE_smpi_comm_out(my_proc_id);
530     retval = MPI_SUCCESS;
531   }
532   smpi_bench_begin();
533   return retval;
534 }
535
536 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag, MPI_Status * status)
537 {
538   int retval = 0;
539
540   smpi_bench_end();
541   if (index == nullptr || flag == nullptr) {
542     retval = MPI_ERR_ARG;
543   } else {
544     int my_proc_id = simgrid::s4u::this_actor::get_pid();
545     TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("testany"));
546     *flag = simgrid::smpi::Request::testany(count, requests, index, status);
547     TRACE_smpi_comm_out(my_proc_id);
548     retval = MPI_SUCCESS;
549   }
550   smpi_bench_begin();
551   return retval;
552 }
553
554 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses)
555 {
556   int retval = 0;
557
558   smpi_bench_end();
559   if (flag == nullptr) {
560     retval = MPI_ERR_ARG;
561   } else {
562     int my_proc_id = simgrid::s4u::this_actor::get_pid();
563     TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("testall"));
564     *flag = simgrid::smpi::Request::testall(count, requests, statuses);
565     TRACE_smpi_comm_out(my_proc_id);
566     retval = MPI_SUCCESS;
567   }
568   smpi_bench_begin();
569   return retval;
570 }
571
572 int PMPI_Testsome(int incount, MPI_Request requests[], int* outcount, int* indices, MPI_Status status[])
573 {
574   int retval = 0;
575
576   smpi_bench_end();
577   if (outcount == nullptr) {
578     retval = MPI_ERR_ARG;
579   } else {
580     int my_proc_id = simgrid::s4u::this_actor::get_pid();
581     TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("testsome"));
582     *outcount = simgrid::smpi::Request::testsome(incount, requests, indices, status);
583     TRACE_smpi_comm_out(my_proc_id);
584     retval    = MPI_SUCCESS;
585   }
586   smpi_bench_begin();
587   return retval;
588 }
589
590 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
591   int retval = 0;
592   smpi_bench_end();
593
594   if (comm == MPI_COMM_NULL) {
595     retval = MPI_ERR_COMM;
596   } else if (source == MPI_PROC_NULL) {
597     if (status != MPI_STATUS_IGNORE){
598       simgrid::smpi::Status::empty(status);
599       status->MPI_SOURCE = MPI_PROC_NULL;
600     }
601     retval = MPI_SUCCESS;
602   } else {
603     simgrid::smpi::Request::probe(source, tag, comm, status);
604     retval = MPI_SUCCESS;
605   }
606   smpi_bench_begin();
607   return retval;
608 }
609
610 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
611   int retval = 0;
612   smpi_bench_end();
613
614   if (flag == nullptr) {
615     retval = MPI_ERR_ARG;
616   } else if (comm == MPI_COMM_NULL) {
617     retval = MPI_ERR_COMM;
618   } else if (source == MPI_PROC_NULL) {
619     *flag=true;
620     if (status != MPI_STATUS_IGNORE){
621       simgrid::smpi::Status::empty(status);
622       status->MPI_SOURCE = MPI_PROC_NULL;
623     }
624     retval = MPI_SUCCESS;
625   } else {
626     simgrid::smpi::Request::iprobe(source, tag, comm, flag, status);
627     retval = MPI_SUCCESS;
628   }
629   smpi_bench_begin();
630   return retval;
631 }
632
633 // TODO: cheinrich: Move declaration to other file? Rename this function - it's used for PMPI_Wait*?
634 static void trace_smpi_recv_helper(MPI_Request* request, MPI_Status* status);
635 static void trace_smpi_recv_helper(MPI_Request* request, MPI_Status* status)
636 {
637   MPI_Request req = *request;
638   if (req != MPI_REQUEST_NULL) { // Received requests become null
639     int src_traced = req->src();
640     // the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
641     int dst_traced = req->dst();
642     if (req->flags() & MPI_REQ_RECV) { // Is this request a wait for RECV?
643       if (src_traced == MPI_ANY_SOURCE)
644         src_traced = (status != MPI_STATUS_IGNORE) ? req->comm()->group()->rank(status->MPI_SOURCE) : req->src();
645       TRACE_smpi_recv(src_traced, dst_traced, req->tag());
646     }
647   }
648 }
649
650 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
651 {
652   int retval = 0;
653
654   smpi_bench_end();
655
656   simgrid::smpi::Status::empty(status);
657
658   if (request == nullptr) {
659     retval = MPI_ERR_ARG;
660   } else if (*request == MPI_REQUEST_NULL) {
661     retval = MPI_SUCCESS;
662   } else {
663     //for tracing, save the handle which might get overriden before we can use the helper on it
664     MPI_Request savedreq = *request;
665     if (savedreq != MPI_REQUEST_NULL && not(savedreq->flags() & MPI_REQ_FINISHED))
666       savedreq->ref();//don't erase te handle in Request::wait, we'll need it later
667     else
668       savedreq = MPI_REQUEST_NULL;
669
670     int my_proc_id = (*request)->comm() != MPI_COMM_NULL
671                          ? simgrid::s4u::this_actor::get_pid()
672                          : -1; // TODO: cheinrich: Check if this correct or if it should be MPI_UNDEFINED
673     TRACE_smpi_comm_in(my_proc_id, __func__,
674                        new simgrid::instr::WaitTIData((*request)->src(), (*request)->dst(), (*request)->tag()));
675
676     simgrid::smpi::Request::wait(request, status);
677     retval = MPI_SUCCESS;
678
679     //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
680     TRACE_smpi_comm_out(my_proc_id);
681     trace_smpi_recv_helper(&savedreq, status);
682     if (savedreq != MPI_REQUEST_NULL)
683       simgrid::smpi::Request::unref(&savedreq);
684   }
685
686   smpi_bench_begin();
687   return retval;
688 }
689
690 int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * status)
691 {
692   if (index == nullptr)
693     return MPI_ERR_ARG;
694
695   if (count <= 0)
696     return MPI_SUCCESS;
697
698   smpi_bench_end();
699   //for tracing, save the handles which might get overriden before we can use the helper on it
700   std::vector<MPI_Request> savedreqs(requests, requests + count);
701   for (MPI_Request& req : savedreqs) {
702     if (req != MPI_REQUEST_NULL && not(req->flags() & MPI_REQ_FINISHED))
703       req->ref();
704     else
705       req = MPI_REQUEST_NULL;
706   }
707
708   int rank_traced = simgrid::s4u::this_actor::get_pid(); // FIXME: In PMPI_Wait, we check if the comm is null?
709   TRACE_smpi_comm_in(rank_traced, __func__, new simgrid::instr::CpuTIData("waitAny", static_cast<double>(count)));
710
711   *index = simgrid::smpi::Request::waitany(count, requests, status);
712
713   if(*index!=MPI_UNDEFINED){
714     trace_smpi_recv_helper(&savedreqs[*index], status);
715     TRACE_smpi_comm_out(rank_traced);
716   }
717
718   for (MPI_Request& req : savedreqs)
719     if (req != MPI_REQUEST_NULL)
720       simgrid::smpi::Request::unref(&req);
721
722   smpi_bench_begin();
723   return MPI_SUCCESS;
724 }
725
726 int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
727 {
728   smpi_bench_end();
729
730   //for tracing, save the handles which might get overriden before we can use the helper on it
731   std::vector<MPI_Request> savedreqs(requests, requests + count);
732   for (MPI_Request& req : savedreqs) {
733     if (req != MPI_REQUEST_NULL && not(req->flags() & MPI_REQ_FINISHED))
734       req->ref();
735     else
736       req = MPI_REQUEST_NULL;
737   }
738
739   int rank_traced = simgrid::s4u::this_actor::get_pid(); // FIXME: In PMPI_Wait, we check if the comm is null?
740   TRACE_smpi_comm_in(rank_traced, __func__, new simgrid::instr::CpuTIData("waitall", static_cast<double>(count)));
741
742   int retval = simgrid::smpi::Request::waitall(count, requests, status);
743
744   for (int i = 0; i < count; i++) {
745     trace_smpi_recv_helper(&savedreqs[i], status!=MPI_STATUSES_IGNORE ? &status[i]: MPI_STATUS_IGNORE);
746   }
747   TRACE_smpi_comm_out(rank_traced);
748
749   for (MPI_Request& req : savedreqs)
750     if (req != MPI_REQUEST_NULL)
751       simgrid::smpi::Request::unref(&req);
752
753   smpi_bench_begin();
754   return retval;
755 }
756
757 int PMPI_Waitsome(int incount, MPI_Request requests[], int *outcount, int *indices, MPI_Status status[])
758 {
759   int retval = 0;
760
761   smpi_bench_end();
762   if (outcount == nullptr) {
763     retval = MPI_ERR_ARG;
764   } else {
765     *outcount = simgrid::smpi::Request::waitsome(incount, requests, indices, status);
766     retval = MPI_SUCCESS;
767   }
768   smpi_bench_begin();
769   return retval;
770 }
771
772 int PMPI_Cancel(MPI_Request* request)
773 {
774   int retval = 0;
775
776   smpi_bench_end();
777   if (*request == MPI_REQUEST_NULL) {
778     retval = MPI_ERR_REQUEST;
779   } else {
780     (*request)->cancel();
781     retval = MPI_SUCCESS;
782   }
783   smpi_bench_begin();
784   return retval;
785 }
786
787 int PMPI_Test_cancelled(MPI_Status* status, int* flag){
788   if(status==MPI_STATUS_IGNORE){
789     *flag=0;
790     return MPI_ERR_ARG;
791   }
792   *flag=simgrid::smpi::Status::cancelled(status);
793   return MPI_SUCCESS;  
794 }
795
796 MPI_Request PMPI_Request_f2c(MPI_Fint request){
797   return static_cast<MPI_Request>(simgrid::smpi::Request::f2c(request));
798 }
799
800 MPI_Fint PMPI_Request_c2f(MPI_Request request) {
801   return request->c2f();
802 }