Logo AND Algorithmique Numérique Distribuée

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