Logo AND Algorithmique Numérique Distribuée

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