1 /* Copyright (c) 2007-2021. 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 #define CHECK_SEND_INPUTS\
23 CHECK_COUNT(2, count)\
24 CHECK_TYPE(3, datatype)\
25 CHECK_BUFFER(1, buf, count, datatype)\
27 if(dst!= MPI_PROC_NULL)\
28 CHECK_RANK(4, dst, comm)\
31 #define CHECK_ISEND_INPUTS\
33 *request = MPI_REQUEST_NULL;\
36 #define CHECK_IRECV_INPUTS\
39 *request = MPI_REQUEST_NULL;\
40 CHECK_COUNT(2, count)\
41 CHECK_TYPE(3, datatype)\
42 CHECK_BUFFER(1, buf, count, datatype)\
44 if(src!=MPI_ANY_SOURCE && src!=MPI_PROC_NULL)\
45 CHECK_RANK(4, src, comm)\
47 /* PMPI User level calls */
49 int PMPI_Send_init(const void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request * request)
54 *request = simgrid::smpi::Request::send_init(buf, count, datatype, dst, tag, comm);
60 int PMPI_Recv_init(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Request * request)
65 *request = simgrid::smpi::Request::recv_init(buf, count, datatype, src, tag, comm);
70 int PMPI_Rsend_init(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm,
73 return PMPI_Send_init(buf, count, datatype, dst, tag, comm, request);
76 int PMPI_Ssend_init(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request)
82 *request = simgrid::smpi::Request::ssend_init(buf, count, datatype, dst, tag, comm);
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.
94 int PMPI_Start(MPI_Request * request)
99 CHECK_REQUEST_VALID(1)
100 if ( *request == MPI_REQUEST_NULL) {
101 retval = MPI_ERR_REQUEST;
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(),
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());
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);
124 int PMPI_Startall(int count, MPI_Request * requests)
128 if (requests == nullptr) {
129 retval = MPI_ERR_ARG;
131 retval = MPI_SUCCESS;
132 for (int i = 0; i < count; i++) {
133 if(requests[i] == MPI_REQUEST_NULL) {
134 retval = MPI_ERR_REQUEST;
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());
147 simgrid::smpi::Request::startall(count, requests);
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());
155 TRACE_smpi_comm_out(my_proc_id);
162 int PMPI_Request_free(MPI_Request * request)
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;
177 int PMPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Request * request)
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);
194 int PMPI_Isend(const void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request * request)
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);
215 int PMPI_Irsend(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm,
216 MPI_Request* request)
218 return PMPI_Isend(buf, count, datatype, dst, tag, comm, request);
221 int PMPI_Issend(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request)
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);
239 int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Status * status)
243 CHECK_COUNT(2, count)
244 CHECK_TYPE(3, datatype)
245 CHECK_BUFFER(1, buf, count, datatype)
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;
255 retval = MPI_SUCCESS;
256 } else if (src!=MPI_ANY_SOURCE && (src >= comm->group()->size() || src <0)){
257 retval = MPI_ERR_RANK;
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)));
265 retval = simgrid::smpi::Request::recv(buf, count, datatype, src, tag, comm, status);
267 // the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
269 if (status != MPI_STATUS_IGNORE)
270 src_traced = getPid(comm, status->MPI_SOURCE);
272 src_traced = getPid(comm, src);
273 if (not TRACE_smpi_view_internals()) {
274 TRACE_smpi_recv(src_traced, my_proc_id, tag);
277 TRACE_smpi_comm_out(my_proc_id);
284 int PMPI_Send(const void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm)
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());
298 simgrid::smpi::Request::send(buf, count, datatype, dst, tag, comm);
299 TRACE_smpi_comm_out(my_proc_id);
304 int PMPI_Rsend(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm)
306 return PMPI_Send(buf, count, datatype, dst, tag, comm);
309 int PMPI_Bsend(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm)
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());
329 simgrid::smpi::Request::bsend(buf, count, datatype, dst, tag, comm);
330 TRACE_smpi_comm_out(my_proc_id);
335 int PMPI_Ibsend(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request)
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);
359 int PMPI_Bsend_init(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request)
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;
371 *request = simgrid::smpi::Request::bsend_init(buf, count, datatype, dst, tag, comm);
372 retval = MPI_SUCCESS;
378 int PMPI_Ssend(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm)
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);
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)
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)
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;
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;
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);
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)));
442 TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, sendtag, sendcount * sendtype->size());
444 simgrid::smpi::Request::sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf, recvcount, recvtype, src,
445 recvtag, comm, status);
446 retval = MPI_SUCCESS;
448 TRACE_smpi_recv(src_traced, my_proc_id, recvtag);
449 TRACE_smpi_comm_out(my_proc_id);
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)
461 CHECK_COUNT(2, count)
462 CHECK_TYPE(3, datatype)
463 CHECK_BUFFER(1, buf, count, datatype)
465 int size = datatype->get_extent() * count;
466 xbt_assert(size > 0);
467 std::vector<char> recvbuf(size);
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);
476 int PMPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
480 if (request == nullptr || flag == nullptr) {
481 retval = MPI_ERR_ARG;
482 } else if (*request == MPI_REQUEST_NULL) {
483 if (status != MPI_STATUS_IGNORE){
485 simgrid::smpi::Status::empty(status);
487 retval = MPI_SUCCESS;
489 int my_proc_id = ((*request)->comm() != MPI_COMM_NULL) ? simgrid::s4u::this_actor::get_pid() : -1;
491 TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("test"));
492 retval = simgrid::smpi::Request::test(request,status, flag);
494 TRACE_smpi_comm_out(my_proc_id);
500 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag, MPI_Status * status)
503 CHECK_COUNT(1, count)
505 if (index == nullptr || flag == nullptr) {
506 retval = MPI_ERR_ARG;
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);
517 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses)
520 CHECK_COUNT(1, count)
522 if (flag == nullptr) {
523 retval = MPI_ERR_ARG;
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);
534 int PMPI_Testsome(int incount, MPI_Request requests[], int* outcount, int* indices, MPI_Status status[])
537 CHECK_COUNT(1, incount)
539 if (outcount == nullptr) {
540 retval = MPI_ERR_ARG;
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);
551 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
556 if(source!=MPI_ANY_SOURCE && source!=MPI_PROC_NULL)\
557 CHECK_RANK(1, source, comm)
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;
564 retval = MPI_SUCCESS;
566 simgrid::smpi::Request::probe(source, tag, comm, status);
567 retval = MPI_SUCCESS;
573 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
577 if(source!=MPI_ANY_SOURCE && source!=MPI_PROC_NULL)\
578 CHECK_RANK(1, source, comm)
580 if (flag == nullptr) {
581 retval = MPI_ERR_ARG;
582 } else if (source == MPI_PROC_NULL) {
584 if (status != MPI_STATUS_IGNORE){
585 simgrid::smpi::Status::empty(status);
586 status->MPI_SOURCE = MPI_PROC_NULL;
588 retval = MPI_SUCCESS;
590 simgrid::smpi::Request::iprobe(source, tag, comm, flag, status);
591 retval = MPI_SUCCESS;
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)
600 const simgrid::smpi::Request* req = *request;
601 // Requests already received are null. Is this request a wait for RECV?
602 if (req != MPI_REQUEST_NULL && (req->flags() & MPI_REQ_RECV)) {
603 int src_traced = req->src();
604 int dst_traced = req->dst();
605 // the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
606 if (src_traced == MPI_ANY_SOURCE && status != MPI_STATUS_IGNORE)
607 src_traced = req->comm()->group()->actor(status->MPI_SOURCE)->get_pid();
608 TRACE_smpi_recv(src_traced, dst_traced, req->tag());
612 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
618 simgrid::smpi::Status::empty(status);
620 CHECK_NULL(1, MPI_ERR_ARG, request)
621 if (*request == MPI_REQUEST_NULL) {
622 retval = MPI_SUCCESS;
624 // for tracing, save the handle which might get overridden before we can use the helper on it
625 MPI_Request savedreq = *request;
626 if (savedreq != MPI_REQUEST_NULL && not(savedreq->flags() & MPI_REQ_FINISHED)
627 && not(savedreq->flags() & MPI_REQ_GENERALIZED))
628 savedreq->ref();//don't erase the handle in Request::wait, we'll need it later
630 savedreq = MPI_REQUEST_NULL;
632 int my_proc_id = (*request)->comm() != MPI_COMM_NULL ? simgrid::s4u::this_actor::get_pid() : -1;
633 TRACE_smpi_comm_in(my_proc_id, __func__,
634 new simgrid::instr::WaitTIData((*request)->src(), (*request)->dst(), (*request)->tag()));
636 retval = simgrid::smpi::Request::wait(request, status);
638 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
639 TRACE_smpi_comm_out(my_proc_id);
640 trace_smpi_recv_helper(&savedreq, status);
641 if (savedreq != MPI_REQUEST_NULL)
642 simgrid::smpi::Request::unref(&savedreq);
649 int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * status)
651 if (index == nullptr)
658 // for tracing, save the handles which might get overridden before we can use the helper on it
659 std::vector<MPI_Request> savedreqs(requests, requests + count);
660 for (MPI_Request& req : savedreqs) {
661 if (req != MPI_REQUEST_NULL && not(req->flags() & MPI_REQ_FINISHED))
664 req = MPI_REQUEST_NULL;
667 int rank_traced = simgrid::s4u::this_actor::get_pid(); // FIXME: In PMPI_Wait, we check if the comm is null?
668 TRACE_smpi_comm_in(rank_traced, __func__, new simgrid::instr::CpuTIData("waitAny", count));
670 *index = simgrid::smpi::Request::waitany(count, requests, status);
672 if(*index!=MPI_UNDEFINED){
673 trace_smpi_recv_helper(&savedreqs[*index], status);
674 TRACE_smpi_comm_out(rank_traced);
677 for (MPI_Request& req : savedreqs)
678 if (req != MPI_REQUEST_NULL)
679 simgrid::smpi::Request::unref(&req);
685 int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
688 CHECK_COUNT(1, count)
689 // for tracing, save the handles which might get overridden before we can use the helper on it
690 std::vector<MPI_Request> savedreqs(requests, requests + count);
691 for (MPI_Request& req : savedreqs) {
692 if (req != MPI_REQUEST_NULL && not(req->flags() & MPI_REQ_FINISHED))
695 req = MPI_REQUEST_NULL;
698 int rank_traced = simgrid::s4u::this_actor::get_pid(); // FIXME: In PMPI_Wait, we check if the comm is null?
699 TRACE_smpi_comm_in(rank_traced, __func__, new simgrid::instr::CpuTIData("waitall", count));
701 int retval = simgrid::smpi::Request::waitall(count, requests, status);
703 for (int i = 0; i < count; i++) {
704 trace_smpi_recv_helper(&savedreqs[i], status!=MPI_STATUSES_IGNORE ? &status[i]: MPI_STATUS_IGNORE);
706 TRACE_smpi_comm_out(rank_traced);
708 for (MPI_Request& req : savedreqs)
709 if (req != MPI_REQUEST_NULL)
710 simgrid::smpi::Request::unref(&req);
716 int PMPI_Waitsome(int incount, MPI_Request requests[], int *outcount, int *indices, MPI_Status status[])
719 CHECK_COUNT(1, incount)
721 if (outcount == nullptr) {
722 retval = MPI_ERR_ARG;
724 *outcount = simgrid::smpi::Request::waitsome(incount, requests, indices, status);
725 retval = MPI_SUCCESS;
731 int PMPI_Cancel(MPI_Request* request)
736 CHECK_REQUEST_VALID(1)
737 if (*request == MPI_REQUEST_NULL) {
738 retval = MPI_ERR_REQUEST;
740 (*request)->cancel();
741 retval = MPI_SUCCESS;
747 int PMPI_Test_cancelled(const MPI_Status* status, int* flag){
748 if(status==MPI_STATUS_IGNORE){
752 *flag=simgrid::smpi::Status::cancelled(status);
756 int PMPI_Status_set_cancelled(MPI_Status* status, int flag){
757 if(status==MPI_STATUS_IGNORE){
760 simgrid::smpi::Status::set_cancelled(status,flag);
764 int PMPI_Status_set_elements(MPI_Status* status, MPI_Datatype datatype, int count){
765 if(status==MPI_STATUS_IGNORE){
768 simgrid::smpi::Status::set_elements(status,datatype, count);
772 int PMPI_Status_set_elements_x(MPI_Status* status, MPI_Datatype datatype, MPI_Count count){
773 if(status==MPI_STATUS_IGNORE){
776 simgrid::smpi::Status::set_elements(status,datatype, static_cast<int>(count));
780 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){
781 return simgrid::smpi::Request::grequest_start(query_fn, free_fn,cancel_fn, extra_state, request);
784 int PMPI_Grequest_complete( MPI_Request request){
785 return simgrid::smpi::Request::grequest_complete(request);
788 int PMPI_Request_get_status( MPI_Request request, int *flag, MPI_Status *status){
789 if(request==MPI_REQUEST_NULL){
791 simgrid::smpi::Status::empty(status);
793 } else if (flag == nullptr) {
796 return simgrid::smpi::Request::get_status(request,flag,status);
799 MPI_Request PMPI_Request_f2c(MPI_Fint request){
801 return MPI_REQUEST_NULL;
802 return simgrid::smpi::Request::f2c(request);
805 MPI_Fint PMPI_Request_c2f(MPI_Request request) {
806 if(request==MPI_REQUEST_NULL)
808 return request->c2f();