Logo AND Algorithmique Numérique Distribuée

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