Logo AND Algorithmique Numérique Distribuée

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