Logo AND Algorithmique Numérique Distribuée

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