Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Update copyright lines for 2022.
[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   int size = datatype->get_extent() * count;
301   if (bsend_buf == nullptr || bsend_buf_size < size + MPI_BSEND_OVERHEAD)
302     return MPI_ERR_BUFFER;
303   TRACE_smpi_comm_in(my_proc_id, __func__,
304                      new simgrid::instr::Pt2PtTIData("bsend", MPI_COMM_WORLD->group()->rank(dst_traced),
305                                                      count,
306                                                      tag, simgrid::smpi::Datatype::encode(datatype)));
307   if (not TRACE_smpi_view_internals()) {
308     TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, tag, count * datatype->size());
309   }
310   simgrid::smpi::Request::bsend(buf, count, datatype, dst, tag, comm);
311   TRACE_smpi_comm_out(my_proc_id);
312   return MPI_SUCCESS;
313 }
314
315 int PMPI_Ibsend(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request)
316 {
317   CHECK_ISEND_INPUTS
318
319   const SmpiBenchGuard suspend_bench;
320   aid_t my_proc_id   = simgrid::s4u::this_actor::get_pid();
321   aid_t trace_dst    = getPid(comm, dst);
322   int bsend_buf_size = 0;
323   void* bsend_buf = nullptr;
324   smpi_process()->bsend_buffer(&bsend_buf, &bsend_buf_size);
325   int size = datatype->get_extent() * count;
326   if (bsend_buf == nullptr || bsend_buf_size < size + MPI_BSEND_OVERHEAD)
327     return MPI_ERR_BUFFER;
328   TRACE_smpi_comm_in(my_proc_id, __func__,
329                      new simgrid::instr::Pt2PtTIData("ibsend", MPI_COMM_WORLD->group()->rank(trace_dst),
330                                                      count,
331                                                      tag, simgrid::smpi::Datatype::encode(datatype)));
332   TRACE_smpi_send(my_proc_id, my_proc_id, trace_dst, tag, count * datatype->size());
333   *request = simgrid::smpi::Request::ibsend(buf, count, datatype, dst, tag, comm);
334   TRACE_smpi_comm_out(my_proc_id);
335   return MPI_SUCCESS;
336 }
337
338 int PMPI_Bsend_init(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request)
339 {
340   CHECK_ISEND_INPUTS
341
342   const SmpiBenchGuard suspend_bench;
343   int retval = 0;
344   int bsend_buf_size = 0;
345   void* bsend_buf = nullptr;
346   smpi_process()->bsend_buffer(&bsend_buf, &bsend_buf_size);
347   if( bsend_buf==nullptr || bsend_buf_size < datatype->get_extent() * count + MPI_BSEND_OVERHEAD ) {
348     retval = MPI_ERR_BUFFER;
349   } else {
350     *request = simgrid::smpi::Request::bsend_init(buf, count, datatype, dst, tag, comm);
351     retval   = MPI_SUCCESS;
352   }
353   return retval;
354 }
355
356 int PMPI_Ssend(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm)
357 {
358   CHECK_SEND_INPUTS
359
360   const SmpiBenchGuard suspend_bench;
361   aid_t my_proc_id = simgrid::s4u::this_actor::get_pid();
362   aid_t dst_traced = getPid(comm, dst);
363   TRACE_smpi_comm_in(my_proc_id, __func__,
364                      new simgrid::instr::Pt2PtTIData("Ssend", MPI_COMM_WORLD->group()->rank(dst_traced),
365                                                      count,
366                                                      tag, simgrid::smpi::Datatype::encode(datatype)));
367   TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, tag, count * datatype->size());
368   simgrid::smpi::Request::ssend(buf, count, datatype, dst, tag, comm);
369   TRACE_smpi_comm_out(my_proc_id);
370   return MPI_SUCCESS;
371 }
372
373 int PMPI_Sendrecv(const void* sendbuf, int sendcount, MPI_Datatype sendtype, int dst, int sendtag, void* recvbuf,
374                   int recvcount, MPI_Datatype recvtype, int src, int recvtag, MPI_Comm comm, MPI_Status* status)
375 {
376   int retval = 0;
377   SET_BUF1(sendbuf)
378   SET_BUF2(recvbuf)
379   CHECK_COUNT(2, sendcount)
380   CHECK_TYPE(3, sendtype)
381   CHECK_TAG(5, sendtag)
382   CHECK_COUNT(7, recvcount)
383   CHECK_TYPE(8, recvtype)
384   CHECK_BUFFER(1, sendbuf, sendcount, sendtype)
385   CHECK_BUFFER(6, recvbuf, recvcount, recvtype)
386   CHECK_TAG(10, recvtag)
387   CHECK_COMM(11)
388   const SmpiBenchGuard suspend_bench;
389
390   if (src == MPI_PROC_NULL) {
391     if(status!=MPI_STATUS_IGNORE){
392       simgrid::smpi::Status::empty(status);
393       status->MPI_SOURCE = MPI_PROC_NULL;
394     }
395     if(dst != MPI_PROC_NULL)
396       simgrid::smpi::Request::send(sendbuf, sendcount, sendtype, dst, sendtag, comm);
397     retval = MPI_SUCCESS;
398   } else if (dst == MPI_PROC_NULL){
399     retval = simgrid::smpi::Request::recv(recvbuf, recvcount, recvtype, src, recvtag, comm, status);
400   } else if (dst >= comm->group()->size() || dst <0 ||
401       (src!=MPI_ANY_SOURCE && (src >= comm->group()->size() || src <0))){
402     retval = MPI_ERR_RANK;
403   } else {
404     aid_t my_proc_id = simgrid::s4u::this_actor::get_pid();
405     aid_t dst_traced = MPI_COMM_WORLD->group()->rank(getPid(comm, dst));
406     aid_t src_traced = MPI_COMM_WORLD->group()->rank(getPid(comm, src));
407
408     // FIXME: Hack the way to trace this one
409     auto dst_hack = std::make_shared<std::vector<int>>();
410     auto src_hack = std::make_shared<std::vector<int>>();
411     dst_hack->push_back(dst_traced);
412     src_hack->push_back(src_traced);
413     TRACE_smpi_comm_in(my_proc_id, __func__,
414                        new simgrid::instr::VarCollTIData(
415                            "sendRecv", -1, sendcount,
416                            dst_hack, recvcount, src_hack,
417                            simgrid::smpi::Datatype::encode(sendtype), simgrid::smpi::Datatype::encode(recvtype)));
418
419     TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, sendtag, sendcount * sendtype->size());
420
421     simgrid::smpi::Request::sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf, recvcount, recvtype, src,
422                                      recvtag, comm, status);
423     retval = MPI_SUCCESS;
424
425     TRACE_smpi_recv(src_traced, my_proc_id, recvtag);
426     TRACE_smpi_comm_out(my_proc_id);
427   }
428
429   return retval;
430 }
431
432 int PMPI_Sendrecv_replace(void* buf, int count, MPI_Datatype datatype, int dst, int sendtag, int src, int recvtag,
433                           MPI_Comm comm, MPI_Status* status)
434 {
435   int retval = 0;
436   SET_BUF1(buf)
437   CHECK_COUNT(2, count)
438   CHECK_TYPE(3, datatype)
439   CHECK_BUFFER(1, buf, count, datatype)
440
441   int size = datatype->get_extent() * count;
442   xbt_assert(size > 0);
443   std::vector<char> recvbuf(size);
444   retval =
445       MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf.data(), count, datatype, src, recvtag, comm, status);
446   if(retval==MPI_SUCCESS){
447     simgrid::smpi::Datatype::copy(recvbuf.data(), count, datatype, buf, count, datatype);
448   }
449   return retval;
450 }
451
452 int PMPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
453 {
454   int retval = 0;
455   const SmpiBenchGuard suspend_bench;
456   if (request == nullptr || flag == nullptr) {
457     retval = MPI_ERR_ARG;
458   } else if (*request == MPI_REQUEST_NULL) {
459     if (status != MPI_STATUS_IGNORE){
460       *flag= true;
461       simgrid::smpi::Status::empty(status);
462     }
463     retval = MPI_SUCCESS;
464   } else {
465     aid_t my_proc_id = ((*request)->comm() != MPI_COMM_NULL) ? simgrid::s4u::this_actor::get_pid() : -1;
466
467     TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("test"));
468     retval = simgrid::smpi::Request::test(request,status, flag);
469
470     TRACE_smpi_comm_out(my_proc_id);
471   }
472   return retval;
473 }
474
475 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag, MPI_Status * status)
476 {
477   int retval = 0;
478   CHECK_COUNT(1, count)
479   const SmpiBenchGuard suspend_bench;
480   if (index == nullptr || flag == nullptr) {
481     retval = MPI_ERR_ARG;
482   } else {
483     aid_t my_proc_id = simgrid::s4u::this_actor::get_pid();
484     TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("testany"));
485     retval = simgrid::smpi::Request::testany(count, requests, index, flag, status);
486     TRACE_smpi_comm_out(my_proc_id);
487   }
488   return retval;
489 }
490
491 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses)
492 {
493   int retval = 0;
494   CHECK_COUNT(1, count)
495   const SmpiBenchGuard suspend_bench;
496   if (flag == nullptr) {
497     retval = MPI_ERR_ARG;
498   } else {
499     aid_t my_proc_id = simgrid::s4u::this_actor::get_pid();
500     TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("testall"));
501     retval = simgrid::smpi::Request::testall(count, requests, flag, statuses);
502     TRACE_smpi_comm_out(my_proc_id);
503   }
504   return retval;
505 }
506
507 int PMPI_Testsome(int incount, MPI_Request requests[], int* outcount, int* indices, MPI_Status status[])
508 {
509   int retval = 0;
510   CHECK_COUNT(1, incount)
511   const SmpiBenchGuard suspend_bench;
512   if (outcount == nullptr) {
513     retval = MPI_ERR_ARG;
514   } else {
515     aid_t my_proc_id = simgrid::s4u::this_actor::get_pid();
516     TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("testsome"));
517     retval = simgrid::smpi::Request::testsome(incount, requests, outcount, indices, status);
518     TRACE_smpi_comm_out(my_proc_id);
519   }
520   return retval;
521 }
522
523 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
524   int retval = 0;
525   const SmpiBenchGuard suspend_bench;
526
527   CHECK_COMM(6)
528   if(source!=MPI_ANY_SOURCE && source!=MPI_PROC_NULL)\
529     CHECK_RANK(1, source, comm)
530   CHECK_TAG(2, tag)
531   if (source == MPI_PROC_NULL) {
532     if (status != MPI_STATUS_IGNORE){
533       simgrid::smpi::Status::empty(status);
534       status->MPI_SOURCE = MPI_PROC_NULL;
535     }
536     retval = MPI_SUCCESS;
537   } else {
538     simgrid::smpi::Request::probe(source, tag, comm, status);
539     retval = MPI_SUCCESS;
540   }
541   return retval;
542 }
543
544 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
545   int retval = 0;
546   const SmpiBenchGuard suspend_bench;
547   CHECK_COMM(6)
548   if(source!=MPI_ANY_SOURCE && source!=MPI_PROC_NULL)\
549     CHECK_RANK(1, source, comm)
550   CHECK_TAG(2, tag)
551   if (flag == nullptr) {
552     retval = MPI_ERR_ARG;
553   } else if (source == MPI_PROC_NULL) {
554     *flag=true;
555     if (status != MPI_STATUS_IGNORE){
556       simgrid::smpi::Status::empty(status);
557       status->MPI_SOURCE = MPI_PROC_NULL;
558     }
559     retval = MPI_SUCCESS;
560   } else {
561     simgrid::smpi::Request::iprobe(source, tag, comm, flag, status);
562     retval = MPI_SUCCESS;
563   }
564   return retval;
565 }
566
567 static void trace_smpi_wait_recv_helper(MPI_Request* request, MPI_Status* status)
568 {
569   const simgrid::smpi::Request* req = *request;
570   // Requests already received are null. Is this request a wait for RECV?
571   if (req != MPI_REQUEST_NULL && (req->flags() & MPI_REQ_RECV)) {
572     aid_t src_traced = req->src();
573     aid_t dst_traced = req->dst();
574     // the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
575     if (src_traced == MPI_ANY_SOURCE && status != MPI_STATUS_IGNORE)
576       src_traced = req->comm()->group()->actor(status->MPI_SOURCE);
577     TRACE_smpi_recv(src_traced, dst_traced, req->tag());
578   }
579 }
580
581 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
582 {
583   int retval = 0;
584
585   const SmpiBenchGuard suspend_bench;
586
587   simgrid::smpi::Status::empty(status);
588
589   CHECK_NULL(1, MPI_ERR_ARG, request)
590   if (*request == MPI_REQUEST_NULL) {
591     retval = MPI_SUCCESS;
592   } else {
593     // for tracing, save the handle which might get overridden before we can use the helper on it
594     MPI_Request savedreq = *request;
595     if (not(savedreq->flags() & (MPI_REQ_FINISHED | MPI_REQ_GENERALIZED | MPI_REQ_NBC)))
596       savedreq->ref();//don't erase the handle in Request::wait, we'll need it later
597     else
598       savedreq = MPI_REQUEST_NULL;
599
600     aid_t my_proc_id = (*request)->comm() != MPI_COMM_NULL ? simgrid::s4u::this_actor::get_pid() : -1;
601     TRACE_smpi_comm_in(my_proc_id, __func__,
602                        new simgrid::instr::WaitTIData(MPI_COMM_WORLD->group()->rank((*request)->src()),
603                                                       MPI_COMM_WORLD->group()->rank((*request)->dst()),
604                                                       (*request)->tag()));
605
606     retval = simgrid::smpi::Request::wait(request, status);
607
608     //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
609     TRACE_smpi_comm_out(my_proc_id);
610     trace_smpi_wait_recv_helper(&savedreq, status);
611     if (savedreq != MPI_REQUEST_NULL)
612       simgrid::smpi::Request::unref(&savedreq);
613   }
614
615   return retval;
616 }
617
618 int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * status)
619 {
620   if (index == nullptr)
621     return MPI_ERR_ARG;
622
623   if (count <= 0)
624     return MPI_SUCCESS;
625
626   const SmpiBenchGuard suspend_bench;
627   // for tracing, save the handles which might get overridden before we can use the helper on it
628   std::vector<MPI_Request> savedreqs(requests, requests + count);
629   for (MPI_Request& req : savedreqs) {
630     if (req != MPI_REQUEST_NULL && not(req->flags() & (MPI_REQ_FINISHED | MPI_REQ_NBC)))
631       req->ref();
632     else
633       req = MPI_REQUEST_NULL;
634   }
635
636   aid_t rank_traced = simgrid::s4u::this_actor::get_pid(); // FIXME: In PMPI_Wait, we check if the comm is null?
637   TRACE_smpi_comm_in(rank_traced, __func__, new simgrid::instr::CpuTIData("waitAny", count));
638
639   *index = simgrid::smpi::Request::waitany(count, requests, status);
640
641   if(*index!=MPI_UNDEFINED){
642     trace_smpi_wait_recv_helper(&savedreqs[*index], status);
643     TRACE_smpi_comm_out(rank_traced);
644   }
645
646   for (MPI_Request& req : savedreqs)
647     if (req != MPI_REQUEST_NULL)
648       simgrid::smpi::Request::unref(&req);
649
650   return MPI_SUCCESS;
651 }
652
653 int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
654 {
655   const SmpiBenchGuard suspend_bench;
656   CHECK_COUNT(1, count)
657   // for tracing, save the handles which might get overridden before we can use the helper on it
658   std::vector<MPI_Request> savedreqs(requests, requests + count);
659   for (MPI_Request& req : savedreqs) {
660     if (req != MPI_REQUEST_NULL && not(req->flags() & (MPI_REQ_FINISHED | MPI_REQ_NBC)))
661       req->ref();
662     else
663       req = MPI_REQUEST_NULL;
664   }
665
666   aid_t rank_traced = simgrid::s4u::this_actor::get_pid(); // FIXME: In PMPI_Wait, we check if the comm is null?
667   TRACE_smpi_comm_in(rank_traced, __func__, new simgrid::instr::CpuTIData("waitall", count));
668
669   int retval = simgrid::smpi::Request::waitall(count, requests, status);
670
671   for (int i = 0; i < count; i++) {
672     trace_smpi_wait_recv_helper(&savedreqs[i], status != MPI_STATUSES_IGNORE ? &status[i] : MPI_STATUS_IGNORE);
673   }
674   TRACE_smpi_comm_out(rank_traced);
675
676   for (MPI_Request& req : savedreqs)
677     if (req != MPI_REQUEST_NULL)
678       simgrid::smpi::Request::unref(&req);
679
680   return retval;
681 }
682
683 int PMPI_Waitsome(int incount, MPI_Request requests[], int *outcount, int *indices, MPI_Status status[])
684 {
685   int retval = 0;
686   CHECK_COUNT(1, incount)
687   const SmpiBenchGuard suspend_bench;
688   if (outcount == nullptr) {
689     retval = MPI_ERR_ARG;
690   } else {
691     *outcount = simgrid::smpi::Request::waitsome(incount, requests, indices, status);
692     retval = MPI_SUCCESS;
693   }
694   return retval;
695 }
696
697 int PMPI_Cancel(MPI_Request* request)
698 {
699   int retval = 0;
700
701   const SmpiBenchGuard suspend_bench;
702   CHECK_REQUEST_VALID(1)
703   if (*request == MPI_REQUEST_NULL) {
704     retval = MPI_ERR_REQUEST;
705   } else {
706     (*request)->cancel();
707     retval = MPI_SUCCESS;
708   }
709   return retval;
710 }
711
712 int PMPI_Test_cancelled(const MPI_Status* status, int* flag){
713   if(status==MPI_STATUS_IGNORE){
714     *flag=0;
715     return MPI_ERR_ARG;
716   }
717   *flag=simgrid::smpi::Status::cancelled(status);
718   return MPI_SUCCESS;
719 }
720
721 int PMPI_Status_set_cancelled(MPI_Status* status, int flag){
722   if(status==MPI_STATUS_IGNORE){
723     return MPI_ERR_ARG;
724   }
725   simgrid::smpi::Status::set_cancelled(status,flag);
726   return MPI_SUCCESS;
727 }
728
729 int PMPI_Status_set_elements(MPI_Status* status, MPI_Datatype datatype, int count){
730   if(status==MPI_STATUS_IGNORE){
731     return MPI_ERR_ARG;
732   }
733   simgrid::smpi::Status::set_elements(status,datatype, count);
734   return MPI_SUCCESS;
735 }
736
737 int PMPI_Status_set_elements_x(MPI_Status* status, MPI_Datatype datatype, MPI_Count count){
738   if(status==MPI_STATUS_IGNORE){
739     return MPI_ERR_ARG;
740   }
741   simgrid::smpi::Status::set_elements(status,datatype, static_cast<int>(count));
742   return MPI_SUCCESS;
743 }
744
745 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){
746   return simgrid::smpi::Request::grequest_start(query_fn, free_fn,cancel_fn, extra_state, request);
747 }
748
749 int PMPI_Grequest_complete( MPI_Request request){
750   return simgrid::smpi::Request::grequest_complete(request);
751 }
752
753 int PMPI_Request_get_status( MPI_Request request, int *flag, MPI_Status *status){
754   if(request==MPI_REQUEST_NULL){
755     *flag=1;
756     simgrid::smpi::Status::empty(status);
757     return MPI_SUCCESS;
758   } else if (flag == nullptr) {
759     return MPI_ERR_ARG;
760   }
761   return simgrid::smpi::Request::get_status(request,flag,status);
762 }
763
764 MPI_Request PMPI_Request_f2c(MPI_Fint request){
765   if(request==-1)
766     return MPI_REQUEST_NULL;
767   return simgrid::smpi::Request::f2c(request);
768 }
769
770 MPI_Fint PMPI_Request_c2f(MPI_Request request) {
771   if(request==MPI_REQUEST_NULL)
772     return -1;
773   return request->c2f();
774 }