Logo AND Algorithmique Numérique Distribuée

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