1 /* Copyright (c) 2007-2019. The SimGrid Team. All rights reserved. */
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. */
7 #include "smpi_comm.hpp"
8 #include "smpi_datatype.hpp"
9 #include "smpi_request.hpp"
10 #include "src/smpi/include/smpi_actor.hpp"
12 XBT_LOG_EXTERNAL_DEFAULT_CATEGORY(smpi_pmpi);
14 static int getPid(MPI_Comm, int);
15 static int getPid(MPI_Comm comm, int id)
17 simgrid::s4u::ActorPtr actor = comm->group()->actor(id);
18 return (actor == nullptr) ? MPI_UNDEFINED : actor->get_pid();
21 /* PMPI User level calls */
23 int PMPI_Send_init(const void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request * request)
28 if (request == nullptr) {
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) {
37 *request = simgrid::smpi::Request::send_init(buf, count, datatype, dst, tag, comm);
41 if (retval != MPI_SUCCESS && request != nullptr)
42 *request = MPI_REQUEST_NULL;
46 int PMPI_Recv_init(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Request * request)
51 if (request == nullptr) {
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) {
60 *request = simgrid::smpi::Request::recv_init(buf, count, datatype, src, tag, comm);
64 if (retval != MPI_SUCCESS && request != nullptr)
65 *request = MPI_REQUEST_NULL;
69 int PMPI_Rsend_init(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm,
72 return PMPI_Send_init(buf, count, datatype, dst, tag, comm, request);
75 int PMPI_Ssend_init(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request)
80 if (request == nullptr) {
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) {
89 *request = simgrid::smpi::Request::ssend_init(buf, count, datatype, dst, tag, comm);
93 if (retval != MPI_SUCCESS && request != nullptr)
94 *request = MPI_REQUEST_NULL;
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.
103 int PMPI_Start(MPI_Request * request)
108 if (request == nullptr || *request == MPI_REQUEST_NULL) {
109 retval = MPI_ERR_REQUEST;
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(),
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());
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);
132 int PMPI_Startall(int count, MPI_Request * requests)
136 if (requests == nullptr) {
137 retval = MPI_ERR_ARG;
139 retval = MPI_SUCCESS;
140 for (int i = 0; i < count; i++) {
141 if(requests[i] == MPI_REQUEST_NULL) {
142 retval = MPI_ERR_REQUEST;
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 MPI_Request req = MPI_REQUEST_NULL;
149 if (not TRACE_smpi_view_internals())
150 for (int i = 0; i < count; i++) {
152 if (req->flags() & MPI_REQ_SEND)
153 TRACE_smpi_send(my_proc_id, my_proc_id, getPid(req->comm(), req->dst()), req->tag(), req->size());
156 simgrid::smpi::Request::startall(count, requests);
158 if (not TRACE_smpi_view_internals())
159 for (int i = 0; i < count; i++) {
161 if (req->flags() & MPI_REQ_RECV)
162 TRACE_smpi_recv(getPid(req->comm(), req->src()), my_proc_id, req->tag());
164 TRACE_smpi_comm_out(my_proc_id);
171 int PMPI_Request_free(MPI_Request * request)
176 if (*request != MPI_REQUEST_NULL) {
177 simgrid::smpi::Request::unref(request);
178 retval = MPI_SUCCESS;
184 int PMPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Request * request)
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;
207 int my_proc_id = simgrid::s4u::this_actor::get_pid();
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)));
214 *request = simgrid::smpi::Request::irecv(buf, count, datatype, src, tag, comm);
215 retval = MPI_SUCCESS;
217 TRACE_smpi_comm_out(my_proc_id);
221 if (retval != MPI_SUCCESS && request != nullptr)
222 *request = MPI_REQUEST_NULL;
227 int PMPI_Isend(const void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request * request)
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;
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)));
255 TRACE_smpi_send(my_proc_id, my_proc_id, trace_dst, tag, count * datatype->size());
257 *request = simgrid::smpi::Request::isend(buf, count, datatype, dst, tag, comm);
258 retval = MPI_SUCCESS;
260 TRACE_smpi_comm_out(my_proc_id);
264 if (retval != MPI_SUCCESS && request!=nullptr)
265 *request = MPI_REQUEST_NULL;
269 int PMPI_Irsend(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm,
270 MPI_Request* request)
272 return PMPI_Isend(buf, count, datatype, dst, tag, comm, request);
275 int PMPI_Issend(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request)
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;
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());
304 *request = simgrid::smpi::Request::issend(buf, count, datatype, dst, tag, comm);
305 retval = MPI_SUCCESS;
307 TRACE_smpi_comm_out(my_proc_id);
311 if (retval != MPI_SUCCESS && request!=nullptr)
312 *request = MPI_REQUEST_NULL;
316 int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Status * status)
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;
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;
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)));
344 simgrid::smpi::Request::recv(buf, count, datatype, src, tag, comm, status);
345 retval = MPI_SUCCESS;
347 // the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
349 if (status != MPI_STATUS_IGNORE)
350 src_traced = getPid(comm, status->MPI_SOURCE);
352 src_traced = getPid(comm, src);
353 if (not TRACE_smpi_view_internals()) {
354 TRACE_smpi_recv(src_traced, my_proc_id, tag);
357 TRACE_smpi_comm_out(my_proc_id);
364 int PMPI_Send(const void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm)
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;
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());
393 simgrid::smpi::Request::send(buf, count, datatype, dst, tag, comm);
394 retval = MPI_SUCCESS;
396 TRACE_smpi_comm_out(my_proc_id);
403 int PMPI_Rsend(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm)
405 return PMPI_Send(buf, count, datatype, dst, tag, comm);
408 int PMPI_Ssend(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm) {
413 if (comm == MPI_COMM_NULL) {
414 retval = MPI_ERR_COMM;
415 } else if (dst == MPI_PROC_NULL) {
416 retval = MPI_SUCCESS;
417 } else if (dst >= comm->group()->size() || dst <0){
418 retval = MPI_ERR_RANK;
419 } else if ((count < 0) || (buf==nullptr && count > 0)) {
420 retval = MPI_ERR_COUNT;
421 } else if (datatype==MPI_DATATYPE_NULL || not datatype->is_valid()) {
422 retval = MPI_ERR_TYPE;
423 } else if(tag<0 && tag != MPI_ANY_TAG){
424 retval = MPI_ERR_TAG;
426 int my_proc_id = simgrid::s4u::this_actor::get_pid();
427 int dst_traced = getPid(comm, dst);
428 TRACE_smpi_comm_in(my_proc_id, __func__,
429 new simgrid::instr::Pt2PtTIData("Ssend", dst,
430 datatype->is_replayable() ? count : count * datatype->size(),
431 tag, simgrid::smpi::Datatype::encode(datatype)));
432 TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, tag, count * datatype->size());
434 simgrid::smpi::Request::ssend(buf, count, datatype, dst, tag, comm);
435 retval = MPI_SUCCESS;
437 TRACE_smpi_comm_out(my_proc_id);
444 int PMPI_Sendrecv(const void* sendbuf, int sendcount, MPI_Datatype sendtype, int dst, int sendtag, void* recvbuf,
445 int recvcount, MPI_Datatype recvtype, int src, int recvtag, MPI_Comm comm, MPI_Status* status)
451 if (comm == MPI_COMM_NULL) {
452 retval = MPI_ERR_COMM;
453 } else if (not sendtype->is_valid() || not recvtype->is_valid()) {
454 retval = MPI_ERR_TYPE;
455 } else if (src == MPI_PROC_NULL) {
456 if(status!=MPI_STATUS_IGNORE){
457 simgrid::smpi::Status::empty(status);
458 status->MPI_SOURCE = MPI_PROC_NULL;
460 if(dst != MPI_PROC_NULL)
461 simgrid::smpi::Request::send(sendbuf, sendcount, sendtype, dst, sendtag, comm);
462 retval = MPI_SUCCESS;
463 }else if (dst == MPI_PROC_NULL){
464 simgrid::smpi::Request::recv(recvbuf, recvcount, recvtype, src, recvtag, comm, status);
465 retval = MPI_SUCCESS;
466 }else if (dst >= comm->group()->size() || dst <0 ||
467 (src!=MPI_ANY_SOURCE && (src >= comm->group()->size() || src <0))){
468 retval = MPI_ERR_RANK;
469 } else if ((sendcount < 0 || recvcount<0) ||
470 (sendbuf==nullptr && sendcount > 0) || (recvbuf==nullptr && recvcount>0)) {
471 retval = MPI_ERR_COUNT;
472 } else if((sendtag<0 && sendtag != MPI_ANY_TAG)||(recvtag<0 && recvtag != MPI_ANY_TAG)){
473 retval = MPI_ERR_TAG;
475 int my_proc_id = simgrid::s4u::this_actor::get_pid();
476 int dst_traced = getPid(comm, dst);
477 int src_traced = getPid(comm, src);
479 // FIXME: Hack the way to trace this one
480 std::vector<int>* dst_hack = new std::vector<int>;
481 std::vector<int>* src_hack = new std::vector<int>;
482 dst_hack->push_back(dst_traced);
483 src_hack->push_back(src_traced);
484 TRACE_smpi_comm_in(my_proc_id, __func__,
485 new simgrid::instr::VarCollTIData(
486 "sendRecv", -1, sendtype->is_replayable() ? sendcount : sendcount * sendtype->size(),
487 dst_hack, recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(), src_hack,
488 simgrid::smpi::Datatype::encode(sendtype), simgrid::smpi::Datatype::encode(recvtype)));
490 TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, sendtag, sendcount * sendtype->size());
492 simgrid::smpi::Request::sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf, recvcount, recvtype, src,
493 recvtag, comm, status);
494 retval = MPI_SUCCESS;
496 TRACE_smpi_recv(src_traced, my_proc_id, recvtag);
497 TRACE_smpi_comm_out(my_proc_id);
504 int PMPI_Sendrecv_replace(void* buf, int count, MPI_Datatype datatype, int dst, int sendtag, int src, int recvtag,
505 MPI_Comm comm, MPI_Status* status)
508 if (datatype==MPI_DATATYPE_NULL || not datatype->is_valid()) {
510 } else if (count < 0) {
511 return MPI_ERR_COUNT;
513 int size = datatype->get_extent() * count;
514 void* recvbuf = xbt_new0(char, size);
515 retval = MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count, datatype, src, recvtag, comm, status);
516 if(retval==MPI_SUCCESS){
517 simgrid::smpi::Datatype::copy(recvbuf, count, datatype, buf, count, datatype);
525 int PMPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
529 if (request == nullptr || flag == nullptr) {
530 retval = MPI_ERR_ARG;
531 } else if (*request == MPI_REQUEST_NULL) {
532 if (status != MPI_STATUS_IGNORE){
534 simgrid::smpi::Status::empty(status);
536 retval = MPI_SUCCESS;
538 int my_proc_id = ((*request)->comm() != MPI_COMM_NULL) ? simgrid::s4u::this_actor::get_pid() : -1;
540 TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("test"));
542 retval = simgrid::smpi::Request::test(request,status, flag);
544 TRACE_smpi_comm_out(my_proc_id);
550 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag, MPI_Status * status)
555 if (index == nullptr || flag == nullptr) {
556 retval = MPI_ERR_ARG;
558 int my_proc_id = simgrid::s4u::this_actor::get_pid();
559 TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("testany"));
560 retval = simgrid::smpi::Request::testany(count, requests, index, flag, status);
561 TRACE_smpi_comm_out(my_proc_id);
567 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses)
572 if (flag == nullptr) {
573 retval = MPI_ERR_ARG;
575 int my_proc_id = simgrid::s4u::this_actor::get_pid();
576 TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("testall"));
577 retval = simgrid::smpi::Request::testall(count, requests, flag, statuses);
578 TRACE_smpi_comm_out(my_proc_id);
584 int PMPI_Testsome(int incount, MPI_Request requests[], int* outcount, int* indices, MPI_Status status[])
589 if (outcount == nullptr) {
590 retval = MPI_ERR_ARG;
592 int my_proc_id = simgrid::s4u::this_actor::get_pid();
593 TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("testsome"));
594 retval = simgrid::smpi::Request::testsome(incount, requests, outcount, indices, status);
595 TRACE_smpi_comm_out(my_proc_id);
601 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
605 if (comm == MPI_COMM_NULL) {
606 retval = MPI_ERR_COMM;
607 } else if (source == MPI_PROC_NULL) {
608 if (status != MPI_STATUS_IGNORE){
609 simgrid::smpi::Status::empty(status);
610 status->MPI_SOURCE = MPI_PROC_NULL;
612 retval = MPI_SUCCESS;
614 simgrid::smpi::Request::probe(source, tag, comm, status);
615 retval = MPI_SUCCESS;
621 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
625 if (flag == nullptr) {
626 retval = MPI_ERR_ARG;
627 } else if (comm == MPI_COMM_NULL) {
628 retval = MPI_ERR_COMM;
629 } else if (source == MPI_PROC_NULL) {
631 if (status != MPI_STATUS_IGNORE){
632 simgrid::smpi::Status::empty(status);
633 status->MPI_SOURCE = MPI_PROC_NULL;
635 retval = MPI_SUCCESS;
637 simgrid::smpi::Request::iprobe(source, tag, comm, flag, status);
638 retval = MPI_SUCCESS;
644 // TODO: cheinrich: Move declaration to other file? Rename this function - it's used for PMPI_Wait*?
645 static void trace_smpi_recv_helper(MPI_Request* request, MPI_Status* status);
646 static void trace_smpi_recv_helper(MPI_Request* request, MPI_Status* status)
648 MPI_Request req = *request;
649 if (req != MPI_REQUEST_NULL) { // Received requests become null
650 int src_traced = req->src();
651 // the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
652 int dst_traced = req->dst();
653 if (req->flags() & MPI_REQ_RECV) { // Is this request a wait for RECV?
654 if (src_traced == MPI_ANY_SOURCE)
655 src_traced = (status != MPI_STATUS_IGNORE) ? req->comm()->group()->rank(status->MPI_SOURCE) : req->src();
656 TRACE_smpi_recv(src_traced, dst_traced, req->tag());
661 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
667 simgrid::smpi::Status::empty(status);
669 if (request == nullptr) {
670 retval = MPI_ERR_ARG;
671 } else if (*request == MPI_REQUEST_NULL) {
672 retval = MPI_SUCCESS;
674 //for tracing, save the handle which might get overriden before we can use the helper on it
675 MPI_Request savedreq = *request;
676 if (savedreq != MPI_REQUEST_NULL && not(savedreq->flags() & MPI_REQ_FINISHED)
677 && not(savedreq->flags() & MPI_REQ_GENERALIZED))
678 savedreq->ref();//don't erase te handle in Request::wait, we'll need it later
680 savedreq = MPI_REQUEST_NULL;
682 int my_proc_id = (*request)->comm() != MPI_COMM_NULL
683 ? simgrid::s4u::this_actor::get_pid()
684 : -1; // TODO: cheinrich: Check if this correct or if it should be MPI_UNDEFINED
685 TRACE_smpi_comm_in(my_proc_id, __func__,
686 new simgrid::instr::WaitTIData((*request)->src(), (*request)->dst(), (*request)->tag()));
688 retval = simgrid::smpi::Request::wait(request, status);
690 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
691 TRACE_smpi_comm_out(my_proc_id);
692 trace_smpi_recv_helper(&savedreq, status);
693 if (savedreq != MPI_REQUEST_NULL)
694 simgrid::smpi::Request::unref(&savedreq);
701 int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * status)
703 if (index == nullptr)
710 //for tracing, save the handles which might get overriden before we can use the helper on it
711 std::vector<MPI_Request> savedreqs(requests, requests + count);
712 for (MPI_Request& req : savedreqs) {
713 if (req != MPI_REQUEST_NULL && not(req->flags() & MPI_REQ_FINISHED))
716 req = MPI_REQUEST_NULL;
719 int rank_traced = simgrid::s4u::this_actor::get_pid(); // FIXME: In PMPI_Wait, we check if the comm is null?
720 TRACE_smpi_comm_in(rank_traced, __func__, new simgrid::instr::CpuTIData("waitAny", count));
722 *index = simgrid::smpi::Request::waitany(count, requests, status);
724 if(*index!=MPI_UNDEFINED){
725 trace_smpi_recv_helper(&savedreqs[*index], status);
726 TRACE_smpi_comm_out(rank_traced);
729 for (MPI_Request& req : savedreqs)
730 if (req != MPI_REQUEST_NULL)
731 simgrid::smpi::Request::unref(&req);
737 int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
741 //for tracing, save the handles which might get overriden before we can use the helper on it
742 std::vector<MPI_Request> savedreqs(requests, requests + count);
743 for (MPI_Request& req : savedreqs) {
744 if (req != MPI_REQUEST_NULL && not(req->flags() & MPI_REQ_FINISHED))
747 req = MPI_REQUEST_NULL;
750 int rank_traced = simgrid::s4u::this_actor::get_pid(); // FIXME: In PMPI_Wait, we check if the comm is null?
751 TRACE_smpi_comm_in(rank_traced, __func__, new simgrid::instr::CpuTIData("waitall", count));
753 int retval = simgrid::smpi::Request::waitall(count, requests, status);
755 for (int i = 0; i < count; i++) {
756 trace_smpi_recv_helper(&savedreqs[i], status!=MPI_STATUSES_IGNORE ? &status[i]: MPI_STATUS_IGNORE);
758 TRACE_smpi_comm_out(rank_traced);
760 for (MPI_Request& req : savedreqs)
761 if (req != MPI_REQUEST_NULL)
762 simgrid::smpi::Request::unref(&req);
768 int PMPI_Waitsome(int incount, MPI_Request requests[], int *outcount, int *indices, MPI_Status status[])
773 if (outcount == nullptr) {
774 retval = MPI_ERR_ARG;
776 *outcount = simgrid::smpi::Request::waitsome(incount, requests, indices, status);
777 retval = MPI_SUCCESS;
783 int PMPI_Cancel(MPI_Request* request)
788 if (*request == MPI_REQUEST_NULL) {
789 retval = MPI_ERR_REQUEST;
791 (*request)->cancel();
792 retval = MPI_SUCCESS;
798 int PMPI_Test_cancelled(const MPI_Status* status, int* flag){
799 if(status==MPI_STATUS_IGNORE){
803 *flag=simgrid::smpi::Status::cancelled(status);
807 int PMPI_Status_set_cancelled(MPI_Status* status, int flag){
808 if(status==MPI_STATUS_IGNORE){
811 simgrid::smpi::Status::set_cancelled(status,flag);
815 int PMPI_Status_set_elements(MPI_Status* status, MPI_Datatype datatype, int count){
816 if(status==MPI_STATUS_IGNORE){
819 simgrid::smpi::Status::set_elements(status,datatype, count);
823 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){
824 return simgrid::smpi::Request::grequest_start(query_fn, free_fn,cancel_fn, extra_state, request);
827 int PMPI_Grequest_complete( MPI_Request request){
828 return simgrid::smpi::Request::grequest_complete(request);
831 int PMPI_Request_get_status( MPI_Request request, int *flag, MPI_Status *status){
832 if(request==MPI_REQUEST_NULL){
834 simgrid::smpi::Status::empty(status);
836 } else if (flag==NULL || status ==NULL){
839 return simgrid::smpi::Request::get_status(request,flag,status);
842 MPI_Request PMPI_Request_f2c(MPI_Fint request){
844 return MPI_REQUEST_NULL;
845 return static_cast<MPI_Request>(simgrid::smpi::Request::f2c(request));
848 MPI_Fint PMPI_Request_c2f(MPI_Request request) {
849 if(request==MPI_REQUEST_NULL)
851 return request->c2f();