Logo AND Algorithmique Numérique Distribuée

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