Logo AND Algorithmique Numérique Distribuée

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