Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Useless test; TODO--.
[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 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", req->dst(),
102                                                        req->size(),
103                                                        req->tag(), 
104                                                        simgrid::smpi::Datatype::encode(req->type())));
105     if (not TRACE_smpi_view_internals() && req->flags() & MPI_REQ_SEND)
106       TRACE_smpi_send(my_proc_id, my_proc_id, getPid(req->comm(), req->dst()), req->tag(), req->size());
107
108     req->start();
109
110     if (not TRACE_smpi_view_internals() && req->flags() & MPI_REQ_RECV)
111       TRACE_smpi_recv(getPid(req->comm(), req->src()), my_proc_id, req->tag());
112     retval = MPI_SUCCESS;
113     TRACE_smpi_comm_out(my_proc_id);
114   }
115   return retval;
116 }
117
118 int PMPI_Startall(int count, MPI_Request * requests)
119 {
120   int retval;
121   const SmpiBenchGuard suspend_bench;
122   if (requests == nullptr) {
123     retval = MPI_ERR_ARG;
124   } else {
125     retval = MPI_SUCCESS;
126     for (int i = 0; i < count; i++) {
127       if(requests[i] == MPI_REQUEST_NULL) {
128         retval = MPI_ERR_REQUEST;
129       }
130     }
131     if(retval != MPI_ERR_REQUEST) {
132       aid_t my_proc_id = simgrid::s4u::this_actor::get_pid();
133       TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("Startall"));
134       if (not TRACE_smpi_view_internals())
135         for (int i = 0; i < count; i++) {
136           const simgrid::smpi::Request* req = requests[i];
137           if (req->flags() & MPI_REQ_SEND)
138             TRACE_smpi_send(my_proc_id, my_proc_id, getPid(req->comm(), req->dst()), req->tag(), req->size());
139         }
140
141       simgrid::smpi::Request::startall(count, requests);
142
143       if (not TRACE_smpi_view_internals())
144         for (int i = 0; i < count; i++) {
145           const simgrid::smpi::Request* req = requests[i];
146           if (req->flags() & MPI_REQ_RECV)
147             TRACE_smpi_recv(getPid(req->comm(), req->src()), my_proc_id, req->tag());
148         }
149       TRACE_smpi_comm_out(my_proc_id);
150     }
151   }
152   return retval;
153 }
154
155 int PMPI_Request_free(MPI_Request * request)
156 {
157   int retval = 0;
158
159   const SmpiBenchGuard suspend_bench;
160   if (*request != MPI_REQUEST_NULL) {
161     (*request)->mark_as_deleted();
162     simgrid::smpi::Request::unref(request);
163     *request = MPI_REQUEST_NULL;
164     retval = MPI_SUCCESS;
165   }
166   return retval;
167 }
168
169 int PMPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Request * request)
170 {
171   CHECK_IRECV_INPUTS
172
173   const SmpiBenchGuard suspend_bench;
174   aid_t my_proc_id = simgrid::s4u::this_actor::get_pid();
175   TRACE_smpi_comm_in(my_proc_id, __func__,
176                      new simgrid::instr::Pt2PtTIData("irecv", src,
177                                                      datatype->is_replayable() ? count : count * datatype->size(),
178                                                      tag, simgrid::smpi::Datatype::encode(datatype)));
179   *request = simgrid::smpi::Request::irecv(buf, count, datatype, src, tag, comm);
180   TRACE_smpi_comm_out(my_proc_id);
181   return MPI_SUCCESS;
182 }
183
184
185 int PMPI_Isend(const void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request * request)
186 {
187   CHECK_ISEND_INPUTS
188
189   const SmpiBenchGuard suspend_bench;
190   int retval = 0;
191   aid_t my_proc_id = simgrid::s4u::this_actor::get_pid();
192   aid_t trace_dst  = getPid(comm, dst);
193   TRACE_smpi_comm_in(my_proc_id, __func__,
194                      new simgrid::instr::Pt2PtTIData("isend", dst,
195                                                      datatype->is_replayable() ? count : count * datatype->size(),
196                                                      tag, simgrid::smpi::Datatype::encode(datatype)));
197   TRACE_smpi_send(my_proc_id, my_proc_id, trace_dst, tag, count * datatype->size());
198   *request = simgrid::smpi::Request::isend(buf, count, datatype, dst, tag, comm);
199   retval = MPI_SUCCESS;
200   TRACE_smpi_comm_out(my_proc_id);
201
202   return retval;
203 }
204
205 int PMPI_Irsend(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm,
206                 MPI_Request* request)
207 {
208   return PMPI_Isend(buf, count, datatype, dst, tag, comm, request);
209 }
210
211 int PMPI_Issend(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request)
212 {
213   CHECK_ISEND_INPUTS
214
215   const SmpiBenchGuard suspend_bench;
216   aid_t my_proc_id = simgrid::s4u::this_actor::get_pid();
217   aid_t trace_dst  = getPid(comm, dst);
218   TRACE_smpi_comm_in(my_proc_id, __func__,
219                      new simgrid::instr::Pt2PtTIData("ISsend", dst,
220                                                      datatype->is_replayable() ? count : count * datatype->size(),
221                                                      tag, simgrid::smpi::Datatype::encode(datatype)));
222   TRACE_smpi_send(my_proc_id, my_proc_id, trace_dst, tag, count * datatype->size());
223   *request = simgrid::smpi::Request::issend(buf, count, datatype, dst, tag, comm);
224   TRACE_smpi_comm_out(my_proc_id);
225   return MPI_SUCCESS;
226 }
227
228 int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Status * status)
229 {
230   int retval = 0;
231   SET_BUF1(buf)
232   CHECK_COUNT(2, count)
233   CHECK_TYPE(3, datatype)
234   CHECK_BUFFER(1, buf, count, datatype)
235   CHECK_TAG(5, tag)
236   CHECK_COMM(6)
237
238   const SmpiBenchGuard suspend_bench;
239   if (src == MPI_PROC_NULL) {
240     if(status != MPI_STATUS_IGNORE){
241       simgrid::smpi::Status::empty(status);
242       status->MPI_SOURCE = MPI_PROC_NULL;
243     }
244     retval = MPI_SUCCESS;
245   } else if (src!=MPI_ANY_SOURCE && (src >= comm->group()->size() || src <0)){
246     retval = MPI_ERR_RANK;
247   } else {
248     aid_t my_proc_id = simgrid::s4u::this_actor::get_pid();
249     TRACE_smpi_comm_in(my_proc_id, __func__,
250                        new simgrid::instr::Pt2PtTIData("recv", src,
251                                                        datatype->is_replayable() ? count : count * datatype->size(),
252                                                        tag, simgrid::smpi::Datatype::encode(datatype)));
253
254     retval = simgrid::smpi::Request::recv(buf, count, datatype, src, tag, comm, status);
255
256     // the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
257     if (not TRACE_smpi_view_internals()) {
258       aid_t src_traced = (status != MPI_STATUS_IGNORE) ? getPid(comm, status->MPI_SOURCE) : getPid(comm, src);
259       TRACE_smpi_recv(src_traced, my_proc_id, tag);
260     }
261     
262     TRACE_smpi_comm_out(my_proc_id);
263   }
264
265   return retval;
266 }
267
268 int PMPI_Send(const void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm)
269 {
270   CHECK_SEND_INPUTS
271
272   const SmpiBenchGuard suspend_bench;
273   aid_t my_proc_id = simgrid::s4u::this_actor::get_pid();
274   aid_t dst_traced = getPid(comm, dst);
275   TRACE_smpi_comm_in(my_proc_id, __func__,
276                      new simgrid::instr::Pt2PtTIData("send", dst,
277                                                      datatype->is_replayable() ? count : count * datatype->size(),
278                                                      tag, simgrid::smpi::Datatype::encode(datatype)));
279   if (not TRACE_smpi_view_internals()) {
280     TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, tag, count * datatype->size());
281   }
282   simgrid::smpi::Request::send(buf, count, datatype, dst, tag, comm);
283   TRACE_smpi_comm_out(my_proc_id);
284   return MPI_SUCCESS;
285 }
286
287 int PMPI_Rsend(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm)
288 {
289   return PMPI_Send(buf, count, datatype, dst, tag, comm);
290 }
291
292 int PMPI_Bsend(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm)
293 {
294   CHECK_SEND_INPUTS
295
296   const SmpiBenchGuard suspend_bench;
297   aid_t my_proc_id   = simgrid::s4u::this_actor::get_pid();
298   aid_t dst_traced   = getPid(comm, dst);
299   int bsend_buf_size = 0;
300   void* bsend_buf = nullptr;
301   smpi_process()->bsend_buffer(&bsend_buf, &bsend_buf_size);
302   int size = datatype->get_extent() * count;
303   if (bsend_buf == nullptr || bsend_buf_size < size + MPI_BSEND_OVERHEAD)
304     return MPI_ERR_BUFFER;
305   TRACE_smpi_comm_in(my_proc_id, __func__,
306                      new simgrid::instr::Pt2PtTIData("bsend", dst,
307                                                      datatype->is_replayable() ? count : count * datatype->size(),
308                                                      tag, simgrid::smpi::Datatype::encode(datatype)));
309   if (not TRACE_smpi_view_internals()) {
310     TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, tag, count * datatype->size());
311   }
312   simgrid::smpi::Request::bsend(buf, count, datatype, dst, tag, comm);
313   TRACE_smpi_comm_out(my_proc_id);
314   return MPI_SUCCESS;
315 }
316
317 int PMPI_Ibsend(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request)
318 {
319   CHECK_ISEND_INPUTS
320
321   const SmpiBenchGuard suspend_bench;
322   aid_t my_proc_id   = simgrid::s4u::this_actor::get_pid();
323   aid_t trace_dst    = getPid(comm, dst);
324   int bsend_buf_size = 0;
325   void* bsend_buf = nullptr;
326   smpi_process()->bsend_buffer(&bsend_buf, &bsend_buf_size);
327   int size = datatype->get_extent() * count;
328   if (bsend_buf == nullptr || bsend_buf_size < size + MPI_BSEND_OVERHEAD)
329     return MPI_ERR_BUFFER;
330   TRACE_smpi_comm_in(my_proc_id, __func__,
331                      new simgrid::instr::Pt2PtTIData("ibsend", dst,
332                                                      datatype->is_replayable() ? count : count * datatype->size(),
333                                                      tag, simgrid::smpi::Datatype::encode(datatype)));
334   TRACE_smpi_send(my_proc_id, my_proc_id, trace_dst, tag, count * datatype->size());
335   *request = simgrid::smpi::Request::ibsend(buf, count, datatype, dst, tag, comm);
336   TRACE_smpi_comm_out(my_proc_id);
337   return MPI_SUCCESS;
338 }
339
340 int PMPI_Bsend_init(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request)
341 {
342   CHECK_ISEND_INPUTS
343
344   const SmpiBenchGuard suspend_bench;
345   int retval = 0;
346   int bsend_buf_size = 0;
347   void* bsend_buf = nullptr;
348   smpi_process()->bsend_buffer(&bsend_buf, &bsend_buf_size);
349   if( bsend_buf==nullptr || bsend_buf_size < datatype->get_extent() * count + MPI_BSEND_OVERHEAD ) {
350     retval = MPI_ERR_BUFFER;
351   } else {
352     *request = simgrid::smpi::Request::bsend_init(buf, count, datatype, dst, tag, comm);
353     retval   = MPI_SUCCESS;
354   }
355   return retval;
356 }
357
358 int PMPI_Ssend(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm)
359 {
360   CHECK_SEND_INPUTS
361
362   const SmpiBenchGuard suspend_bench;
363   aid_t my_proc_id = simgrid::s4u::this_actor::get_pid();
364   aid_t dst_traced = getPid(comm, dst);
365   TRACE_smpi_comm_in(my_proc_id, __func__,
366                      new simgrid::instr::Pt2PtTIData("Ssend", dst,
367                                                      datatype->is_replayable() ? count : count * datatype->size(),
368                                                      tag, simgrid::smpi::Datatype::encode(datatype)));
369   TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, tag, count * datatype->size());
370   simgrid::smpi::Request::ssend(buf, count, datatype, dst, tag, comm);
371   TRACE_smpi_comm_out(my_proc_id);
372   return MPI_SUCCESS;
373 }
374
375 int PMPI_Sendrecv(const void* sendbuf, int sendcount, MPI_Datatype sendtype, int dst, int sendtag, void* recvbuf,
376                   int recvcount, MPI_Datatype recvtype, int src, int recvtag, MPI_Comm comm, MPI_Status* status)
377 {
378   int retval = 0;
379   SET_BUF1(sendbuf)
380   SET_BUF2(recvbuf)
381   CHECK_COUNT(2, sendcount)
382   CHECK_TYPE(3, sendtype)
383   CHECK_TAG(5, sendtag)
384   CHECK_COUNT(7, recvcount)
385   CHECK_TYPE(8, recvtype)
386   CHECK_BUFFER(1, sendbuf, sendcount, sendtype)
387   CHECK_BUFFER(6, recvbuf, recvcount, recvtype)
388   CHECK_TAG(10, recvtag)
389   CHECK_COMM(11)
390   const SmpiBenchGuard suspend_bench;
391
392   if (src == MPI_PROC_NULL) {
393     if(status!=MPI_STATUS_IGNORE){
394       simgrid::smpi::Status::empty(status);
395       status->MPI_SOURCE = MPI_PROC_NULL;
396     }
397     if(dst != MPI_PROC_NULL)
398       simgrid::smpi::Request::send(sendbuf, sendcount, sendtype, dst, sendtag, comm);
399     retval = MPI_SUCCESS;
400   } else if (dst == MPI_PROC_NULL){
401     retval = simgrid::smpi::Request::recv(recvbuf, recvcount, recvtype, src, recvtag, comm, status);
402   } else if (dst >= comm->group()->size() || dst <0 ||
403       (src!=MPI_ANY_SOURCE && (src >= comm->group()->size() || src <0))){
404     retval = MPI_ERR_RANK;
405   } else {
406     aid_t my_proc_id = simgrid::s4u::this_actor::get_pid();
407     aid_t dst_traced = getPid(comm, dst);
408     aid_t src_traced = getPid(comm, src);
409
410     // FIXME: Hack the way to trace this one
411     auto dst_hack = std::make_shared<std::vector<int>>();
412     auto src_hack = std::make_shared<std::vector<int>>();
413     dst_hack->push_back(dst_traced);
414     src_hack->push_back(src_traced);
415     TRACE_smpi_comm_in(my_proc_id, __func__,
416                        new simgrid::instr::VarCollTIData(
417                            "sendRecv", -1, sendtype->is_replayable() ? sendcount : sendcount * sendtype->size(),
418                            dst_hack, recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(), src_hack,
419                            simgrid::smpi::Datatype::encode(sendtype), simgrid::smpi::Datatype::encode(recvtype)));
420
421     TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, sendtag, sendcount * sendtype->size());
422
423     simgrid::smpi::Request::sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf, recvcount, recvtype, src,
424                                      recvtag, comm, status);
425     retval = MPI_SUCCESS;
426
427     TRACE_smpi_recv(src_traced, my_proc_id, recvtag);
428     TRACE_smpi_comm_out(my_proc_id);
429   }
430
431   return retval;
432 }
433
434 int PMPI_Sendrecv_replace(void* buf, int count, MPI_Datatype datatype, int dst, int sendtag, int src, int recvtag,
435                           MPI_Comm comm, MPI_Status* status)
436 {
437   int retval = 0;
438   SET_BUF1(buf)
439   CHECK_COUNT(2, count)
440   CHECK_TYPE(3, datatype)
441   CHECK_BUFFER(1, buf, count, datatype)
442
443   int size = datatype->get_extent() * count;
444   xbt_assert(size > 0);
445   std::vector<char> recvbuf(size);
446   retval =
447       MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf.data(), count, datatype, src, recvtag, comm, status);
448   if(retval==MPI_SUCCESS){
449     simgrid::smpi::Datatype::copy(recvbuf.data(), count, datatype, buf, count, datatype);
450   }
451   return retval;
452 }
453
454 int PMPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
455 {
456   int retval = 0;
457   const SmpiBenchGuard suspend_bench;
458   if (request == nullptr || flag == nullptr) {
459     retval = MPI_ERR_ARG;
460   } else if (*request == MPI_REQUEST_NULL) {
461     if (status != MPI_STATUS_IGNORE){
462       *flag= true;
463       simgrid::smpi::Status::empty(status);
464     }
465     retval = MPI_SUCCESS;
466   } else {
467     aid_t my_proc_id = ((*request)->comm() != MPI_COMM_NULL) ? simgrid::s4u::this_actor::get_pid() : -1;
468
469     TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("test"));
470     retval = simgrid::smpi::Request::test(request,status, flag);
471
472     TRACE_smpi_comm_out(my_proc_id);
473   }
474   return retval;
475 }
476
477 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag, MPI_Status * status)
478 {
479   int retval = 0;
480   CHECK_COUNT(1, count)
481   const SmpiBenchGuard suspend_bench;
482   if (index == nullptr || flag == nullptr) {
483     retval = MPI_ERR_ARG;
484   } else {
485     aid_t my_proc_id = simgrid::s4u::this_actor::get_pid();
486     TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("testany"));
487     retval = simgrid::smpi::Request::testany(count, requests, index, flag, status);
488     TRACE_smpi_comm_out(my_proc_id);
489   }
490   return retval;
491 }
492
493 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses)
494 {
495   int retval = 0;
496   CHECK_COUNT(1, count)
497   const SmpiBenchGuard suspend_bench;
498   if (flag == nullptr) {
499     retval = MPI_ERR_ARG;
500   } else {
501     aid_t my_proc_id = simgrid::s4u::this_actor::get_pid();
502     TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("testall"));
503     retval = simgrid::smpi::Request::testall(count, requests, flag, statuses);
504     TRACE_smpi_comm_out(my_proc_id);
505   }
506   return retval;
507 }
508
509 int PMPI_Testsome(int incount, MPI_Request requests[], int* outcount, int* indices, MPI_Status status[])
510 {
511   int retval = 0;
512   CHECK_COUNT(1, incount)
513   const SmpiBenchGuard suspend_bench;
514   if (outcount == nullptr) {
515     retval = MPI_ERR_ARG;
516   } else {
517     aid_t my_proc_id = simgrid::s4u::this_actor::get_pid();
518     TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("testsome"));
519     retval = simgrid::smpi::Request::testsome(incount, requests, outcount, indices, status);
520     TRACE_smpi_comm_out(my_proc_id);
521   }
522   return retval;
523 }
524
525 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
526   int retval = 0;
527   const SmpiBenchGuard suspend_bench;
528
529   CHECK_COMM(6)
530   if(source!=MPI_ANY_SOURCE && source!=MPI_PROC_NULL)\
531     CHECK_RANK(1, source, comm)
532   CHECK_TAG(2, tag)
533   if (source == MPI_PROC_NULL) {
534     if (status != MPI_STATUS_IGNORE){
535       simgrid::smpi::Status::empty(status);
536       status->MPI_SOURCE = MPI_PROC_NULL;
537     }
538     retval = MPI_SUCCESS;
539   } else {
540     simgrid::smpi::Request::probe(source, tag, comm, status);
541     retval = MPI_SUCCESS;
542   }
543   return retval;
544 }
545
546 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
547   int retval = 0;
548   const SmpiBenchGuard suspend_bench;
549   CHECK_COMM(6)
550   if(source!=MPI_ANY_SOURCE && source!=MPI_PROC_NULL)\
551     CHECK_RANK(1, source, comm)
552   CHECK_TAG(2, tag)
553   if (flag == nullptr) {
554     retval = MPI_ERR_ARG;
555   } else if (source == MPI_PROC_NULL) {
556     *flag=true;
557     if (status != MPI_STATUS_IGNORE){
558       simgrid::smpi::Status::empty(status);
559       status->MPI_SOURCE = MPI_PROC_NULL;
560     }
561     retval = MPI_SUCCESS;
562   } else {
563     simgrid::smpi::Request::iprobe(source, tag, comm, flag, status);
564     retval = MPI_SUCCESS;
565   }
566   return retval;
567 }
568
569 static void trace_smpi_wait_recv_helper(MPI_Request* request, MPI_Status* status)
570 {
571   const simgrid::smpi::Request* req = *request;
572   // Requests already received are null. Is this request a wait for RECV?
573   if (req != MPI_REQUEST_NULL && (req->flags() & MPI_REQ_RECV)) {
574     aid_t src_traced = req->src();
575     aid_t dst_traced = req->dst();
576     // the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
577     if (src_traced == MPI_ANY_SOURCE && status != MPI_STATUS_IGNORE)
578       src_traced = req->comm()->group()->actor(status->MPI_SOURCE);
579     TRACE_smpi_recv(src_traced, dst_traced, req->tag());
580   }
581 }
582
583 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
584 {
585   int retval = 0;
586
587   const SmpiBenchGuard suspend_bench;
588
589   simgrid::smpi::Status::empty(status);
590
591   CHECK_NULL(1, MPI_ERR_ARG, request) 
592   if (*request == MPI_REQUEST_NULL) {
593     retval = MPI_SUCCESS;
594   } else {
595     // for tracing, save the handle which might get overridden before we can use the helper on it
596     MPI_Request savedreq = *request;
597     if (not(savedreq->flags() & (MPI_REQ_FINISHED | MPI_REQ_GENERALIZED | MPI_REQ_NBC)))
598       savedreq->ref();//don't erase the handle in Request::wait, we'll need it later
599     else
600       savedreq = MPI_REQUEST_NULL;
601
602     aid_t my_proc_id = (*request)->comm() != MPI_COMM_NULL ? simgrid::s4u::this_actor::get_pid() : -1;
603     TRACE_smpi_comm_in(my_proc_id, __func__,
604                        new simgrid::instr::WaitTIData((*request)->src(), (*request)->dst(), (*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 }