Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Finally rename smpi::Group::actor_pid() back to actor().
[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 comm, int id)
15 {
16   return comm->group()->actor(id);
17 }
18
19 #define CHECK_SEND_INPUTS\
20   SET_BUF1(buf)\
21   CHECK_COUNT(2, count)\
22   CHECK_TYPE(3, datatype)\
23   CHECK_BUFFER(1, buf, count, datatype)\
24   CHECK_COMM(6)\
25   if(dst!= MPI_PROC_NULL)\
26     CHECK_RANK(4, dst, comm)\
27   CHECK_TAG(5, tag)
28
29 #define CHECK_ISEND_INPUTS\
30   CHECK_REQUEST(7)\
31   *request = MPI_REQUEST_NULL;\
32   CHECK_SEND_INPUTS
33   
34 #define CHECK_IRECV_INPUTS\
35   SET_BUF1(buf)\
36   CHECK_REQUEST(7)\
37   *request = MPI_REQUEST_NULL;\
38   CHECK_COUNT(2, count)\
39   CHECK_TYPE(3, datatype)\
40   CHECK_BUFFER(1, buf, count, datatype)\
41   CHECK_COMM(6)\
42   if(src!=MPI_ANY_SOURCE && src!=MPI_PROC_NULL)\
43     CHECK_RANK(4, src, comm)\
44   CHECK_TAG(5, tag)
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_VALID(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     (*request)->mark_as_deleted();
167     simgrid::smpi::Request::unref(request);
168     *request = MPI_REQUEST_NULL;
169     retval = MPI_SUCCESS;
170   }
171   smpi_bench_begin();
172   return retval;
173 }
174
175 int PMPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Request * request)
176 {
177   CHECK_IRECV_INPUTS
178
179   smpi_bench_end();
180   int my_proc_id = simgrid::s4u::this_actor::get_pid();
181   TRACE_smpi_comm_in(my_proc_id, __func__,
182                      new simgrid::instr::Pt2PtTIData("irecv", src,
183                                                      datatype->is_replayable() ? count : count * datatype->size(),
184                                                      tag, simgrid::smpi::Datatype::encode(datatype)));
185   *request = simgrid::smpi::Request::irecv(buf, count, datatype, src, tag, comm);
186   TRACE_smpi_comm_out(my_proc_id);
187   smpi_bench_begin();
188   return MPI_SUCCESS;
189 }
190
191
192 int PMPI_Isend(const void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request * request)
193 {
194   CHECK_ISEND_INPUTS
195
196   smpi_bench_end();
197   int retval = 0;
198   int my_proc_id = simgrid::s4u::this_actor::get_pid();
199   int trace_dst = getPid(comm, dst);
200   TRACE_smpi_comm_in(my_proc_id, __func__,
201                      new simgrid::instr::Pt2PtTIData("isend", dst,
202                                                      datatype->is_replayable() ? count : count * datatype->size(),
203                                                      tag, simgrid::smpi::Datatype::encode(datatype)));
204   TRACE_smpi_send(my_proc_id, my_proc_id, trace_dst, tag, count * datatype->size());
205   *request = simgrid::smpi::Request::isend(buf, count, datatype, dst, tag, comm);
206   retval = MPI_SUCCESS;
207   TRACE_smpi_comm_out(my_proc_id);
208   smpi_bench_begin();
209
210   return retval;
211 }
212
213 int PMPI_Irsend(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm,
214                 MPI_Request* request)
215 {
216   return PMPI_Isend(buf, count, datatype, dst, tag, comm, request);
217 }
218
219 int PMPI_Issend(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request)
220 {
221   CHECK_ISEND_INPUTS
222
223   smpi_bench_end();
224   int my_proc_id = simgrid::s4u::this_actor::get_pid();
225   int trace_dst = getPid(comm, dst);
226   TRACE_smpi_comm_in(my_proc_id, __func__,
227                      new simgrid::instr::Pt2PtTIData("ISsend", dst,
228                                                      datatype->is_replayable() ? count : count * datatype->size(),
229                                                      tag, simgrid::smpi::Datatype::encode(datatype)));
230   TRACE_smpi_send(my_proc_id, my_proc_id, trace_dst, tag, count * datatype->size());
231   *request = simgrid::smpi::Request::issend(buf, count, datatype, dst, tag, comm);
232   TRACE_smpi_comm_out(my_proc_id);
233   smpi_bench_begin();
234   return MPI_SUCCESS;
235 }
236
237 int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Status * status)
238 {
239   int retval = 0;
240   SET_BUF1(buf)
241   CHECK_COUNT(2, count)
242   CHECK_TYPE(3, datatype)
243   CHECK_BUFFER(1, buf, count, datatype)
244   CHECK_TAG(5, tag)
245   CHECK_COMM(6)
246
247   smpi_bench_end();
248   if (src == MPI_PROC_NULL) {
249     if(status != MPI_STATUS_IGNORE){
250       simgrid::smpi::Status::empty(status);
251       status->MPI_SOURCE = MPI_PROC_NULL;
252     }
253     retval = MPI_SUCCESS;
254   } else if (src!=MPI_ANY_SOURCE && (src >= comm->group()->size() || src <0)){
255     retval = MPI_ERR_RANK;
256   } else {
257     int my_proc_id = simgrid::s4u::this_actor::get_pid();
258     TRACE_smpi_comm_in(my_proc_id, __func__,
259                        new simgrid::instr::Pt2PtTIData("recv", src,
260                                                        datatype->is_replayable() ? count : count * datatype->size(),
261                                                        tag, simgrid::smpi::Datatype::encode(datatype)));
262
263     retval = simgrid::smpi::Request::recv(buf, count, datatype, src, tag, comm, status);
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   SET_BUF1(sendbuf)
399   SET_BUF2(recvbuf)
400   CHECK_COUNT(2, sendcount)
401   CHECK_TYPE(3, sendtype)
402   CHECK_TAG(5, sendtag)
403   CHECK_COUNT(7, recvcount)
404   CHECK_TYPE(8, recvtype)
405   CHECK_BUFFER(1, sendbuf, sendcount, sendtype)
406   CHECK_BUFFER(6, recvbuf, recvcount, recvtype)
407   CHECK_TAG(10, recvtag)
408   CHECK_COMM(11)
409   smpi_bench_end();
410
411   if (src == MPI_PROC_NULL) {
412     if(status!=MPI_STATUS_IGNORE){
413       simgrid::smpi::Status::empty(status);
414       status->MPI_SOURCE = MPI_PROC_NULL;
415     }
416     if(dst != MPI_PROC_NULL)
417       simgrid::smpi::Request::send(sendbuf, sendcount, sendtype, dst, sendtag, comm);
418     retval = MPI_SUCCESS;
419   } else if (dst == MPI_PROC_NULL){
420     retval = simgrid::smpi::Request::recv(recvbuf, recvcount, recvtype, src, recvtag, comm, status);
421   } else if (dst >= comm->group()->size() || dst <0 ||
422       (src!=MPI_ANY_SOURCE && (src >= comm->group()->size() || src <0))){
423     retval = MPI_ERR_RANK;
424   } else {
425     int my_proc_id         = simgrid::s4u::this_actor::get_pid();
426     int dst_traced         = getPid(comm, dst);
427     int src_traced         = getPid(comm, src);
428
429     // FIXME: Hack the way to trace this one
430     auto dst_hack = std::make_shared<std::vector<int>>();
431     auto src_hack = std::make_shared<std::vector<int>>();
432     dst_hack->push_back(dst_traced);
433     src_hack->push_back(src_traced);
434     TRACE_smpi_comm_in(my_proc_id, __func__,
435                        new simgrid::instr::VarCollTIData(
436                            "sendRecv", -1, sendtype->is_replayable() ? sendcount : sendcount * sendtype->size(),
437                            dst_hack, recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(), src_hack,
438                            simgrid::smpi::Datatype::encode(sendtype), simgrid::smpi::Datatype::encode(recvtype)));
439
440     TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, sendtag, sendcount * sendtype->size());
441
442     simgrid::smpi::Request::sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf, recvcount, recvtype, src,
443                                      recvtag, comm, status);
444     retval = MPI_SUCCESS;
445
446     TRACE_smpi_recv(src_traced, my_proc_id, recvtag);
447     TRACE_smpi_comm_out(my_proc_id);
448   }
449
450   smpi_bench_begin();
451   return retval;
452 }
453
454 int PMPI_Sendrecv_replace(void* buf, int count, MPI_Datatype datatype, int dst, int sendtag, int src, int recvtag,
455                           MPI_Comm comm, MPI_Status* status)
456 {
457   int retval = 0;
458   SET_BUF1(buf)
459   CHECK_COUNT(2, count)
460   CHECK_TYPE(3, datatype)
461   CHECK_BUFFER(1, buf, count, datatype)
462
463   int size = datatype->get_extent() * count;
464   xbt_assert(size > 0);
465   std::vector<char> recvbuf(size);
466   retval =
467       MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf.data(), count, datatype, src, recvtag, comm, status);
468   if(retval==MPI_SUCCESS){
469     simgrid::smpi::Datatype::copy(recvbuf.data(), count, datatype, buf, count, datatype);
470   }
471   return retval;
472 }
473
474 int PMPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
475 {
476   int retval = 0;
477   smpi_bench_end();
478   if (request == nullptr || flag == nullptr) {
479     retval = MPI_ERR_ARG;
480   } else if (*request == MPI_REQUEST_NULL) {
481     if (status != MPI_STATUS_IGNORE){
482       *flag= true;
483       simgrid::smpi::Status::empty(status);
484     }
485     retval = MPI_SUCCESS;
486   } else {
487     int my_proc_id = ((*request)->comm() != MPI_COMM_NULL) ? simgrid::s4u::this_actor::get_pid() : -1;
488
489     TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("test"));
490     retval = simgrid::smpi::Request::test(request,status, flag);
491
492     TRACE_smpi_comm_out(my_proc_id);
493   }
494   smpi_bench_begin();
495   return retval;
496 }
497
498 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag, MPI_Status * status)
499 {
500   int retval = 0;
501   CHECK_COUNT(1, count)
502   smpi_bench_end();
503   if (index == nullptr || flag == nullptr) {
504     retval = MPI_ERR_ARG;
505   } else {
506     int my_proc_id = simgrid::s4u::this_actor::get_pid();
507     TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("testany"));
508     retval = simgrid::smpi::Request::testany(count, requests, index, flag, status);
509     TRACE_smpi_comm_out(my_proc_id);
510   }
511   smpi_bench_begin();
512   return retval;
513 }
514
515 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses)
516 {
517   int retval = 0;
518   CHECK_COUNT(1, count)
519   smpi_bench_end();
520   if (flag == nullptr) {
521     retval = MPI_ERR_ARG;
522   } else {
523     int my_proc_id = simgrid::s4u::this_actor::get_pid();
524     TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("testall"));
525     retval = simgrid::smpi::Request::testall(count, requests, flag, statuses);
526     TRACE_smpi_comm_out(my_proc_id);
527   }
528   smpi_bench_begin();
529   return retval;
530 }
531
532 int PMPI_Testsome(int incount, MPI_Request requests[], int* outcount, int* indices, MPI_Status status[])
533 {
534   int retval = 0;
535   CHECK_COUNT(1, incount)
536   smpi_bench_end();
537   if (outcount == nullptr) {
538     retval = MPI_ERR_ARG;
539   } else {
540     int my_proc_id = simgrid::s4u::this_actor::get_pid();
541     TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("testsome"));
542     retval = simgrid::smpi::Request::testsome(incount, requests, outcount, indices, status);
543     TRACE_smpi_comm_out(my_proc_id);
544   }
545   smpi_bench_begin();
546   return retval;
547 }
548
549 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
550   int retval = 0;
551   smpi_bench_end();
552
553   CHECK_COMM(6)
554   if(source!=MPI_ANY_SOURCE && source!=MPI_PROC_NULL)\
555     CHECK_RANK(1, source, comm)
556   CHECK_TAG(2, tag)
557   if (source == MPI_PROC_NULL) {
558     if (status != MPI_STATUS_IGNORE){
559       simgrid::smpi::Status::empty(status);
560       status->MPI_SOURCE = MPI_PROC_NULL;
561     }
562     retval = MPI_SUCCESS;
563   } else {
564     simgrid::smpi::Request::probe(source, tag, comm, status);
565     retval = MPI_SUCCESS;
566   }
567   smpi_bench_begin();
568   return retval;
569 }
570
571 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
572   int retval = 0;
573   smpi_bench_end();
574   CHECK_COMM(6)
575   if(source!=MPI_ANY_SOURCE && source!=MPI_PROC_NULL)\
576     CHECK_RANK(1, source, comm)
577   CHECK_TAG(2, tag)
578   if (flag == nullptr) {
579     retval = MPI_ERR_ARG;
580   } else if (source == MPI_PROC_NULL) {
581     *flag=true;
582     if (status != MPI_STATUS_IGNORE){
583       simgrid::smpi::Status::empty(status);
584       status->MPI_SOURCE = MPI_PROC_NULL;
585     }
586     retval = MPI_SUCCESS;
587   } else {
588     simgrid::smpi::Request::iprobe(source, tag, comm, flag, status);
589     retval = MPI_SUCCESS;
590   }
591   smpi_bench_begin();
592   return retval;
593 }
594
595 // TODO: cheinrich: Move declaration to other file? Rename this function - it's used for PMPI_Wait*?
596 static void trace_smpi_recv_helper(MPI_Request* request, MPI_Status* status)
597 {
598   const simgrid::smpi::Request* req = *request;
599   // Requests already received are null. Is this request a wait for RECV?
600   if (req != MPI_REQUEST_NULL && (req->flags() & MPI_REQ_RECV)) {
601     int src_traced = req->src();
602     int dst_traced = req->dst();
603     // the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
604     if (src_traced == MPI_ANY_SOURCE && status != MPI_STATUS_IGNORE)
605       src_traced = req->comm()->group()->actor(status->MPI_SOURCE);
606     TRACE_smpi_recv(src_traced, dst_traced, req->tag());
607   }
608 }
609
610 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
611 {
612   int retval = 0;
613
614   smpi_bench_end();
615
616   simgrid::smpi::Status::empty(status);
617
618   CHECK_NULL(1, MPI_ERR_ARG, request) 
619   if (*request == MPI_REQUEST_NULL) {
620     retval = MPI_SUCCESS;
621   } else {
622     // for tracing, save the handle which might get overridden before we can use the helper on it
623     MPI_Request savedreq = *request;
624     if (savedreq != MPI_REQUEST_NULL && not(savedreq->flags() & MPI_REQ_FINISHED)
625     && not(savedreq->flags() & MPI_REQ_GENERALIZED))
626       savedreq->ref();//don't erase the handle in Request::wait, we'll need it later
627     else
628       savedreq = MPI_REQUEST_NULL;
629
630     int my_proc_id = (*request)->comm() != MPI_COMM_NULL ? simgrid::s4u::this_actor::get_pid() : -1;
631     TRACE_smpi_comm_in(my_proc_id, __func__,
632                        new simgrid::instr::WaitTIData((*request)->src(), (*request)->dst(), (*request)->tag()));
633
634     retval = simgrid::smpi::Request::wait(request, status);
635
636     //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
637     TRACE_smpi_comm_out(my_proc_id);
638     trace_smpi_recv_helper(&savedreq, status);
639     if (savedreq != MPI_REQUEST_NULL)
640       simgrid::smpi::Request::unref(&savedreq);
641   }
642
643   smpi_bench_begin();
644   return retval;
645 }
646
647 int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * status)
648 {
649   if (index == nullptr)
650     return MPI_ERR_ARG;
651
652   if (count <= 0)
653     return MPI_SUCCESS;
654
655   smpi_bench_end();
656   // for tracing, save the handles which might get overridden before we can use the helper on it
657   std::vector<MPI_Request> savedreqs(requests, requests + count);
658   for (MPI_Request& req : savedreqs) {
659     if (req != MPI_REQUEST_NULL && not(req->flags() & MPI_REQ_FINISHED))
660       req->ref();
661     else
662       req = MPI_REQUEST_NULL;
663   }
664
665   int rank_traced = simgrid::s4u::this_actor::get_pid(); // FIXME: In PMPI_Wait, we check if the comm is null?
666   TRACE_smpi_comm_in(rank_traced, __func__, new simgrid::instr::CpuTIData("waitAny", count));
667
668   *index = simgrid::smpi::Request::waitany(count, requests, status);
669
670   if(*index!=MPI_UNDEFINED){
671     trace_smpi_recv_helper(&savedreqs[*index], status);
672     TRACE_smpi_comm_out(rank_traced);
673   }
674
675   for (MPI_Request& req : savedreqs)
676     if (req != MPI_REQUEST_NULL)
677       simgrid::smpi::Request::unref(&req);
678
679   smpi_bench_begin();
680   return MPI_SUCCESS;
681 }
682
683 int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
684 {
685   smpi_bench_end();
686   CHECK_COUNT(1, count)
687   // for tracing, save the handles which might get overridden before we can use the helper on it
688   std::vector<MPI_Request> savedreqs(requests, requests + count);
689   for (MPI_Request& req : savedreqs) {
690     if (req != MPI_REQUEST_NULL && not(req->flags() & MPI_REQ_FINISHED))
691       req->ref();
692     else
693       req = MPI_REQUEST_NULL;
694   }
695
696   int rank_traced = simgrid::s4u::this_actor::get_pid(); // FIXME: In PMPI_Wait, we check if the comm is null?
697   TRACE_smpi_comm_in(rank_traced, __func__, new simgrid::instr::CpuTIData("waitall", count));
698
699   int retval = simgrid::smpi::Request::waitall(count, requests, status);
700
701   for (int i = 0; i < count; i++) {
702     trace_smpi_recv_helper(&savedreqs[i], status!=MPI_STATUSES_IGNORE ? &status[i]: MPI_STATUS_IGNORE);
703   }
704   TRACE_smpi_comm_out(rank_traced);
705
706   for (MPI_Request& req : savedreqs)
707     if (req != MPI_REQUEST_NULL)
708       simgrid::smpi::Request::unref(&req);
709
710   smpi_bench_begin();
711   return retval;
712 }
713
714 int PMPI_Waitsome(int incount, MPI_Request requests[], int *outcount, int *indices, MPI_Status status[])
715 {
716   int retval = 0;
717   CHECK_COUNT(1, incount)
718   smpi_bench_end();
719   if (outcount == nullptr) {
720     retval = MPI_ERR_ARG;
721   } else {
722     *outcount = simgrid::smpi::Request::waitsome(incount, requests, indices, status);
723     retval = MPI_SUCCESS;
724   }
725   smpi_bench_begin();
726   return retval;
727 }
728
729 int PMPI_Cancel(MPI_Request* request)
730 {
731   int retval = 0;
732
733   smpi_bench_end();
734   CHECK_REQUEST_VALID(1)
735   if (*request == MPI_REQUEST_NULL) {
736     retval = MPI_ERR_REQUEST;
737   } else {
738     (*request)->cancel();
739     retval = MPI_SUCCESS;
740   }
741   smpi_bench_begin();
742   return retval;
743 }
744
745 int PMPI_Test_cancelled(const MPI_Status* status, int* flag){
746   if(status==MPI_STATUS_IGNORE){
747     *flag=0;
748     return MPI_ERR_ARG;
749   }
750   *flag=simgrid::smpi::Status::cancelled(status);
751   return MPI_SUCCESS;  
752 }
753
754 int PMPI_Status_set_cancelled(MPI_Status* status, int flag){
755   if(status==MPI_STATUS_IGNORE){
756     return MPI_ERR_ARG;
757   }
758   simgrid::smpi::Status::set_cancelled(status,flag);
759   return MPI_SUCCESS;  
760 }
761
762 int PMPI_Status_set_elements(MPI_Status* status, MPI_Datatype datatype, int count){
763   if(status==MPI_STATUS_IGNORE){
764     return MPI_ERR_ARG;
765   }
766   simgrid::smpi::Status::set_elements(status,datatype, count);
767   return MPI_SUCCESS;
768 }
769
770 int PMPI_Status_set_elements_x(MPI_Status* status, MPI_Datatype datatype, MPI_Count count){
771   if(status==MPI_STATUS_IGNORE){
772     return MPI_ERR_ARG;
773   }
774   simgrid::smpi::Status::set_elements(status,datatype, static_cast<int>(count));
775   return MPI_SUCCESS;
776 }
777
778 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){
779   return simgrid::smpi::Request::grequest_start(query_fn, free_fn,cancel_fn, extra_state, request);
780 }
781
782 int PMPI_Grequest_complete( MPI_Request request){
783   return simgrid::smpi::Request::grequest_complete(request);
784 }
785
786 int PMPI_Request_get_status( MPI_Request request, int *flag, MPI_Status *status){
787   if(request==MPI_REQUEST_NULL){
788     *flag=1;
789     simgrid::smpi::Status::empty(status);
790     return MPI_SUCCESS;
791   } else if (flag == nullptr) {
792     return MPI_ERR_ARG;
793   }
794   return simgrid::smpi::Request::get_status(request,flag,status);
795 }
796
797 MPI_Request PMPI_Request_f2c(MPI_Fint request){
798   if(request==-1)
799     return MPI_REQUEST_NULL;
800   return simgrid::smpi::Request::f2c(request);
801 }
802
803 MPI_Fint PMPI_Request_c2f(MPI_Request request) {
804   if(request==MPI_REQUEST_NULL)
805     return -1;
806   return request->c2f();
807 }