Logo AND Algorithmique Numérique Distribuée

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