1 /* Copyright (c) 2007-2018. 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_process.hpp"
10 #include "smpi_request.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(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 (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 (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_Ssend_init(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request)
74 if (request == nullptr) {
76 } else if (comm == MPI_COMM_NULL) {
77 retval = MPI_ERR_COMM;
78 } else if (not datatype->is_valid()) {
79 retval = MPI_ERR_TYPE;
80 } else if (dst == MPI_PROC_NULL) {
83 *request = simgrid::smpi::Request::ssend_init(buf, count, datatype, dst, tag, comm);
87 if (retval != MPI_SUCCESS && request != nullptr)
88 *request = MPI_REQUEST_NULL;
93 * This function starts a request returned by init functions such as
94 * MPI_Send_init(), MPI_Ssend_init (see above), and friends.
95 * They should already have performed sanity checks.
97 int PMPI_Start(MPI_Request * request)
102 if (request == nullptr || *request == MPI_REQUEST_NULL) {
103 retval = MPI_ERR_REQUEST;
105 MPI_Request req = *request;
106 int my_proc_id = (req->comm() != MPI_COMM_NULL) ? simgrid::s4u::this_actor::get_pid() : -1;
107 TRACE_smpi_comm_in(my_proc_id, __func__,
108 new simgrid::instr::Pt2PtTIData("Start", req->dst(),
111 simgrid::smpi::Datatype::encode(req->type())));
112 if (not TRACE_smpi_view_internals() && req->flags() & MPI_REQ_SEND)
113 TRACE_smpi_send(my_proc_id, my_proc_id, getPid(req->comm(), req->dst()), req->tag(), req->size());
117 if (not TRACE_smpi_view_internals() && req->flags() & MPI_REQ_RECV)
118 TRACE_smpi_recv(getPid(req->comm(), req->src()), my_proc_id, req->tag());
119 retval = MPI_SUCCESS;
120 TRACE_smpi_comm_out(my_proc_id);
126 int PMPI_Startall(int count, MPI_Request * requests)
130 if (requests == nullptr) {
131 retval = MPI_ERR_ARG;
133 retval = MPI_SUCCESS;
134 for (int i = 0; i < count; i++) {
135 if(requests[i] == MPI_REQUEST_NULL) {
136 retval = MPI_ERR_REQUEST;
139 if(retval != MPI_ERR_REQUEST) {
140 int my_proc_id = simgrid::s4u::this_actor::get_pid();
141 TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("Startall"));
142 MPI_Request req = MPI_REQUEST_NULL;
143 if (not TRACE_smpi_view_internals())
144 for (int i = 0; i < count; i++) {
146 if (req->flags() & MPI_REQ_SEND)
147 TRACE_smpi_send(my_proc_id, my_proc_id, getPid(req->comm(), req->dst()), req->tag(), req->size());
150 simgrid::smpi::Request::startall(count, requests);
152 if (not TRACE_smpi_view_internals())
153 for (int i = 0; i < count; i++) {
155 if (req->flags() & MPI_REQ_RECV)
156 TRACE_smpi_recv(getPid(req->comm(), req->src()), my_proc_id, req->tag());
158 TRACE_smpi_comm_out(my_proc_id);
165 int PMPI_Request_free(MPI_Request * request)
170 if (*request == MPI_REQUEST_NULL) {
171 retval = MPI_ERR_ARG;
173 simgrid::smpi::Request::unref(request);
174 retval = MPI_SUCCESS;
180 int PMPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Request * request)
186 if (request == nullptr) {
187 retval = MPI_ERR_ARG;
188 } else if (comm == MPI_COMM_NULL) {
189 retval = MPI_ERR_COMM;
190 } else if (src == MPI_PROC_NULL) {
191 *request = MPI_REQUEST_NULL;
192 retval = MPI_SUCCESS;
193 } else if (src!=MPI_ANY_SOURCE && (src >= comm->group()->size() || src <0)){
194 retval = MPI_ERR_RANK;
195 } else if ((count < 0) || (buf==nullptr && count > 0)) {
196 retval = MPI_ERR_COUNT;
197 } else if (not datatype->is_valid()) {
198 retval = MPI_ERR_TYPE;
199 } else if(tag<0 && tag != MPI_ANY_TAG){
200 retval = MPI_ERR_TAG;
203 int my_proc_id = simgrid::s4u::this_actor::get_pid();
205 TRACE_smpi_comm_in(my_proc_id, __func__,
206 new simgrid::instr::Pt2PtTIData("Irecv", src,
207 datatype->is_replayable() ? count : count * datatype->size(),
208 tag, simgrid::smpi::Datatype::encode(datatype)));
210 *request = simgrid::smpi::Request::irecv(buf, count, datatype, src, tag, comm);
211 retval = MPI_SUCCESS;
213 TRACE_smpi_comm_out(my_proc_id);
217 if (retval != MPI_SUCCESS && request != nullptr)
218 *request = MPI_REQUEST_NULL;
223 int PMPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request * request)
228 if (request == nullptr) {
229 retval = MPI_ERR_ARG;
230 } else if (comm == MPI_COMM_NULL) {
231 retval = MPI_ERR_COMM;
232 } else if (dst == MPI_PROC_NULL) {
233 *request = MPI_REQUEST_NULL;
234 retval = MPI_SUCCESS;
235 } else if (dst >= comm->group()->size() || dst <0){
236 retval = MPI_ERR_RANK;
237 } else if ((count < 0) || (buf==nullptr && count > 0)) {
238 retval = MPI_ERR_COUNT;
239 } else if (not datatype->is_valid()) {
240 retval = MPI_ERR_TYPE;
241 } else if(tag<0 && tag != MPI_ANY_TAG){
242 retval = MPI_ERR_TAG;
244 int my_proc_id = simgrid::s4u::this_actor::get_pid();
245 int trace_dst = getPid(comm, dst);
246 TRACE_smpi_comm_in(my_proc_id, __func__,
247 new simgrid::instr::Pt2PtTIData("Isend", dst,
248 datatype->is_replayable() ? count : count * datatype->size(),
249 tag, simgrid::smpi::Datatype::encode(datatype)));
251 TRACE_smpi_send(my_proc_id, my_proc_id, trace_dst, tag, count * datatype->size());
253 *request = simgrid::smpi::Request::isend(buf, count, datatype, dst, tag, comm);
254 retval = MPI_SUCCESS;
256 TRACE_smpi_comm_out(my_proc_id);
260 if (retval != MPI_SUCCESS && request!=nullptr)
261 *request = MPI_REQUEST_NULL;
265 int PMPI_Issend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request)
270 if (request == nullptr) {
271 retval = MPI_ERR_ARG;
272 } else if (comm == MPI_COMM_NULL) {
273 retval = MPI_ERR_COMM;
274 } else if (dst == MPI_PROC_NULL) {
275 *request = MPI_REQUEST_NULL;
276 retval = MPI_SUCCESS;
277 } else if (dst >= comm->group()->size() || dst <0){
278 retval = MPI_ERR_RANK;
279 } else if ((count < 0)|| (buf==nullptr && count > 0)) {
280 retval = MPI_ERR_COUNT;
281 } else if (not datatype->is_valid()) {
282 retval = MPI_ERR_TYPE;
283 } else if(tag<0 && tag != MPI_ANY_TAG){
284 retval = MPI_ERR_TAG;
286 int my_proc_id = simgrid::s4u::this_actor::get_pid();
287 int trace_dst = getPid(comm, dst);
288 TRACE_smpi_comm_in(my_proc_id, __func__,
289 new simgrid::instr::Pt2PtTIData("ISsend", dst,
290 datatype->is_replayable() ? count : count * datatype->size(),
291 tag, simgrid::smpi::Datatype::encode(datatype)));
292 TRACE_smpi_send(my_proc_id, my_proc_id, trace_dst, tag, count * datatype->size());
294 *request = simgrid::smpi::Request::issend(buf, count, datatype, dst, tag, comm);
295 retval = MPI_SUCCESS;
297 TRACE_smpi_comm_out(my_proc_id);
301 if (retval != MPI_SUCCESS && request!=nullptr)
302 *request = MPI_REQUEST_NULL;
306 int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Status * status)
311 if (comm == MPI_COMM_NULL) {
312 retval = MPI_ERR_COMM;
313 } else if (src == MPI_PROC_NULL) {
314 if(status != MPI_STATUS_IGNORE){
315 simgrid::smpi::Status::empty(status);
316 status->MPI_SOURCE = MPI_PROC_NULL;
318 retval = MPI_SUCCESS;
319 } else if (src!=MPI_ANY_SOURCE && (src >= comm->group()->size() || src <0)){
320 retval = MPI_ERR_RANK;
321 } else if ((count < 0) || (buf==nullptr && count > 0)) {
322 retval = MPI_ERR_COUNT;
323 } else if (not datatype->is_valid()) {
324 retval = MPI_ERR_TYPE;
325 } else if(tag<0 && tag != MPI_ANY_TAG){
326 retval = MPI_ERR_TAG;
328 int my_proc_id = simgrid::s4u::this_actor::get_pid();
329 TRACE_smpi_comm_in(my_proc_id, __func__,
330 new simgrid::instr::Pt2PtTIData("recv", src,
331 datatype->is_replayable() ? count : count * datatype->size(),
332 tag, simgrid::smpi::Datatype::encode(datatype)));
334 simgrid::smpi::Request::recv(buf, count, datatype, src, tag, comm, status);
335 retval = MPI_SUCCESS;
337 // the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
339 if (status != MPI_STATUS_IGNORE)
340 src_traced = getPid(comm, status->MPI_SOURCE);
342 src_traced = getPid(comm, src);
343 if (not TRACE_smpi_view_internals()) {
344 TRACE_smpi_recv(src_traced, my_proc_id, tag);
347 TRACE_smpi_comm_out(my_proc_id);
354 int PMPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm)
360 if (comm == MPI_COMM_NULL) {
361 retval = MPI_ERR_COMM;
362 } else if (dst == MPI_PROC_NULL) {
363 retval = MPI_SUCCESS;
364 } else if (dst >= comm->group()->size() || dst <0){
365 retval = MPI_ERR_RANK;
366 } else if ((count < 0) || (buf == nullptr && count > 0)) {
367 retval = MPI_ERR_COUNT;
368 } else if (not datatype->is_valid()) {
369 retval = MPI_ERR_TYPE;
370 } else if(tag < 0 && tag != MPI_ANY_TAG){
371 retval = MPI_ERR_TAG;
373 int my_proc_id = simgrid::s4u::this_actor::get_pid();
374 int dst_traced = getPid(comm, dst);
375 TRACE_smpi_comm_in(my_proc_id, __func__,
376 new simgrid::instr::Pt2PtTIData("send", dst,
377 datatype->is_replayable() ? count : count * datatype->size(),
378 tag, simgrid::smpi::Datatype::encode(datatype)));
379 if (not TRACE_smpi_view_internals()) {
380 TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, tag, count * datatype->size());
383 simgrid::smpi::Request::send(buf, count, datatype, dst, tag, comm);
384 retval = MPI_SUCCESS;
386 TRACE_smpi_comm_out(my_proc_id);
393 int PMPI_Ssend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm) {
398 if (comm == MPI_COMM_NULL) {
399 retval = MPI_ERR_COMM;
400 } else if (dst == MPI_PROC_NULL) {
401 retval = MPI_SUCCESS;
402 } else if (dst >= comm->group()->size() || dst <0){
403 retval = MPI_ERR_RANK;
404 } else if ((count < 0) || (buf==nullptr && count > 0)) {
405 retval = MPI_ERR_COUNT;
406 } else if (not datatype->is_valid()) {
407 retval = MPI_ERR_TYPE;
408 } else if(tag<0 && tag != MPI_ANY_TAG){
409 retval = MPI_ERR_TAG;
411 int my_proc_id = simgrid::s4u::this_actor::get_pid();
412 int dst_traced = getPid(comm, dst);
413 TRACE_smpi_comm_in(my_proc_id, __func__,
414 new simgrid::instr::Pt2PtTIData("Ssend", dst,
415 datatype->is_replayable() ? count : count * datatype->size(),
416 tag, simgrid::smpi::Datatype::encode(datatype)));
417 TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, tag, count * datatype->size());
419 simgrid::smpi::Request::ssend(buf, count, datatype, dst, tag, comm);
420 retval = MPI_SUCCESS;
422 TRACE_smpi_comm_out(my_proc_id);
429 int PMPI_Sendrecv(void* sendbuf, int sendcount, MPI_Datatype sendtype, int dst, int sendtag, void* recvbuf,
430 int recvcount, MPI_Datatype recvtype, int src, int recvtag, MPI_Comm comm, MPI_Status* status)
436 if (comm == MPI_COMM_NULL) {
437 retval = MPI_ERR_COMM;
438 } else if (not sendtype->is_valid() || not recvtype->is_valid()) {
439 retval = MPI_ERR_TYPE;
440 } else if (src == MPI_PROC_NULL) {
441 if(status!=MPI_STATUS_IGNORE){
442 simgrid::smpi::Status::empty(status);
443 status->MPI_SOURCE = MPI_PROC_NULL;
445 if(dst != MPI_PROC_NULL)
446 simgrid::smpi::Request::send(sendbuf, sendcount, sendtype, dst, sendtag, comm);
447 retval = MPI_SUCCESS;
448 }else if (dst == MPI_PROC_NULL){
449 simgrid::smpi::Request::recv(recvbuf, recvcount, recvtype, src, recvtag, comm, status);
450 retval = MPI_SUCCESS;
451 }else if (dst >= comm->group()->size() || dst <0 ||
452 (src!=MPI_ANY_SOURCE && (src >= comm->group()->size() || src <0))){
453 retval = MPI_ERR_RANK;
454 } else if ((sendcount < 0 || recvcount<0) ||
455 (sendbuf==nullptr && sendcount > 0) || (recvbuf==nullptr && recvcount>0)) {
456 retval = MPI_ERR_COUNT;
457 } else if((sendtag<0 && sendtag != MPI_ANY_TAG)||(recvtag<0 && recvtag != MPI_ANY_TAG)){
458 retval = MPI_ERR_TAG;
460 int my_proc_id = simgrid::s4u::this_actor::get_pid();
461 int dst_traced = getPid(comm, dst);
462 int src_traced = getPid(comm, src);
464 // FIXME: Hack the way to trace this one
465 std::vector<int>* dst_hack = new std::vector<int>;
466 std::vector<int>* src_hack = new std::vector<int>;
467 dst_hack->push_back(dst_traced);
468 src_hack->push_back(src_traced);
469 TRACE_smpi_comm_in(my_proc_id, __func__,
470 new simgrid::instr::VarCollTIData(
471 "sendRecv", -1, sendtype->is_replayable() ? sendcount : sendcount * sendtype->size(),
472 dst_hack, recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(), src_hack,
473 simgrid::smpi::Datatype::encode(sendtype), simgrid::smpi::Datatype::encode(recvtype)));
475 TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, sendtag, sendcount * sendtype->size());
477 simgrid::smpi::Request::sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf, recvcount, recvtype, src,
478 recvtag, comm, status);
479 retval = MPI_SUCCESS;
481 TRACE_smpi_recv(src_traced, my_proc_id, recvtag);
482 TRACE_smpi_comm_out(my_proc_id);
489 int PMPI_Sendrecv_replace(void* buf, int count, MPI_Datatype datatype, int dst, int sendtag, int src, int recvtag,
490 MPI_Comm comm, MPI_Status* status)
493 if (not datatype->is_valid()) {
495 } else if (count < 0) {
496 return MPI_ERR_COUNT;
498 int size = datatype->get_extent() * count;
499 void* recvbuf = xbt_new0(char, size);
500 retval = MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count, datatype, src, recvtag, comm, status);
501 if(retval==MPI_SUCCESS){
502 simgrid::smpi::Datatype::copy(recvbuf, count, datatype, buf, count, datatype);
510 int PMPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
514 if (request == nullptr || flag == nullptr) {
515 retval = MPI_ERR_ARG;
516 } else if (*request == MPI_REQUEST_NULL) {
517 if (status != MPI_STATUS_IGNORE){
519 simgrid::smpi::Status::empty(status);
521 retval = MPI_SUCCESS;
523 int my_proc_id = ((*request)->comm() != MPI_COMM_NULL) ? simgrid::s4u::this_actor::get_pid() : -1;
525 TRACE_smpi_testing_in(my_proc_id);
527 *flag = simgrid::smpi::Request::test(request,status);
529 TRACE_smpi_testing_out(my_proc_id);
530 retval = MPI_SUCCESS;
536 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag, MPI_Status * status)
541 if (index == nullptr || flag == nullptr) {
542 retval = MPI_ERR_ARG;
544 *flag = simgrid::smpi::Request::testany(count, requests, index, status);
545 retval = MPI_SUCCESS;
551 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses)
556 if (flag == nullptr) {
557 retval = MPI_ERR_ARG;
559 *flag = simgrid::smpi::Request::testall(count, requests, statuses);
560 retval = MPI_SUCCESS;
566 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
570 if (comm == MPI_COMM_NULL) {
571 retval = MPI_ERR_COMM;
572 } else if (source == MPI_PROC_NULL) {
573 if (status != MPI_STATUS_IGNORE){
574 simgrid::smpi::Status::empty(status);
575 status->MPI_SOURCE = MPI_PROC_NULL;
577 retval = MPI_SUCCESS;
579 simgrid::smpi::Request::probe(source, tag, comm, status);
580 retval = MPI_SUCCESS;
586 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
590 if (flag == nullptr) {
591 retval = MPI_ERR_ARG;
592 } else if (comm == MPI_COMM_NULL) {
593 retval = MPI_ERR_COMM;
594 } else if (source == MPI_PROC_NULL) {
596 if (status != MPI_STATUS_IGNORE){
597 simgrid::smpi::Status::empty(status);
598 status->MPI_SOURCE = MPI_PROC_NULL;
600 retval = MPI_SUCCESS;
602 simgrid::smpi::Request::iprobe(source, tag, comm, flag, status);
603 retval = MPI_SUCCESS;
609 // TODO: cheinrich: Move declaration to other file? Rename this function - it's used for PMPI_Wait*?
610 static void trace_smpi_recv_helper(MPI_Request* request, MPI_Status* status);
611 static void trace_smpi_recv_helper(MPI_Request* request, MPI_Status* status)
613 MPI_Request req = *request;
614 if (req != MPI_REQUEST_NULL) { // Received requests become null
615 int src_traced = req->src();
616 // the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
617 int dst_traced = req->dst();
618 if (req->flags() & MPI_REQ_RECV) { // Is this request a wait for RECV?
619 if (src_traced == MPI_ANY_SOURCE)
620 src_traced = (status != MPI_STATUS_IGNORE) ? req->comm()->group()->rank(status->MPI_SOURCE) : req->src();
621 TRACE_smpi_recv(src_traced, dst_traced, req->tag());
626 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
632 simgrid::smpi::Status::empty(status);
634 if (request == nullptr) {
635 retval = MPI_ERR_ARG;
636 } else if (*request == MPI_REQUEST_NULL) {
637 retval = MPI_SUCCESS;
639 //for tracing, save the handle which might get overriden before we can use the helper on it
640 MPI_Request savedreq = *request;
641 if (savedreq != MPI_REQUEST_NULL && not(savedreq->flags() & MPI_REQ_FINISHED))
642 savedreq->ref();//don't erase te handle in Request::wait, we'll need it later
644 savedreq = MPI_REQUEST_NULL;
646 int my_proc_id = (*request)->comm() != MPI_COMM_NULL
647 ? simgrid::s4u::this_actor::get_pid()
648 : -1; // TODO: cheinrich: Check if this correct or if it should be MPI_UNDEFINED
649 TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("wait"));
651 simgrid::smpi::Request::wait(request, status);
652 retval = MPI_SUCCESS;
654 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
655 TRACE_smpi_comm_out(my_proc_id);
656 trace_smpi_recv_helper(&savedreq, status);
657 if (savedreq != MPI_REQUEST_NULL)
658 simgrid::smpi::Request::unref(&savedreq);
665 int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * status)
667 if (index == nullptr)
674 //for tracing, save the handles which might get overriden before we can use the helper on it
675 std::vector<MPI_Request> savedreqs(requests, requests + count);
676 for (MPI_Request& req : savedreqs) {
677 if (req != MPI_REQUEST_NULL && not(req->flags() & MPI_REQ_FINISHED))
680 req = MPI_REQUEST_NULL;
683 int rank_traced = simgrid::s4u::this_actor::get_pid(); // FIXME: In PMPI_Wait, we check if the comm is null?
684 TRACE_smpi_comm_in(rank_traced, __func__, new simgrid::instr::CpuTIData("waitAny", static_cast<double>(count)));
686 *index = simgrid::smpi::Request::waitany(count, requests, status);
688 if(*index!=MPI_UNDEFINED){
689 trace_smpi_recv_helper(&savedreqs[*index], status);
690 TRACE_smpi_comm_out(rank_traced);
693 for (MPI_Request& req : savedreqs)
694 if (req != MPI_REQUEST_NULL)
695 simgrid::smpi::Request::unref(&req);
701 int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
705 //for tracing, save the handles which might get overriden before we can use the helper on it
706 std::vector<MPI_Request> savedreqs(requests, requests + count);
707 for (MPI_Request& req : savedreqs) {
708 if (req != MPI_REQUEST_NULL && not(req->flags() & MPI_REQ_FINISHED))
711 req = MPI_REQUEST_NULL;
714 int rank_traced = simgrid::s4u::this_actor::get_pid(); // FIXME: In PMPI_Wait, we check if the comm is null?
715 TRACE_smpi_comm_in(rank_traced, __func__, new simgrid::instr::CpuTIData("waitAll", static_cast<double>(count)));
717 int retval = simgrid::smpi::Request::waitall(count, requests, status);
719 for (int i = 0; i < count; i++) {
720 trace_smpi_recv_helper(&savedreqs[i], status!=MPI_STATUSES_IGNORE ? &status[i]: MPI_STATUS_IGNORE);
722 TRACE_smpi_comm_out(rank_traced);
724 for (MPI_Request& req : savedreqs)
725 if (req != MPI_REQUEST_NULL)
726 simgrid::smpi::Request::unref(&req);
732 int PMPI_Waitsome(int incount, MPI_Request requests[], int *outcount, int *indices, MPI_Status status[])
737 if (outcount == nullptr) {
738 retval = MPI_ERR_ARG;
740 *outcount = simgrid::smpi::Request::waitsome(incount, requests, indices, status);
741 retval = MPI_SUCCESS;
747 int PMPI_Testsome(int incount, MPI_Request requests[], int* outcount, int* indices, MPI_Status status[])
752 if (outcount == nullptr) {
753 retval = MPI_ERR_ARG;
755 *outcount = simgrid::smpi::Request::testsome(incount, requests, indices, status);
756 retval = MPI_SUCCESS;
762 int PMPI_Cancel(MPI_Request* request)
767 if (*request == MPI_REQUEST_NULL) {
768 retval = MPI_ERR_REQUEST;
770 (*request)->cancel();
771 retval = MPI_SUCCESS;
777 int PMPI_Test_cancelled(MPI_Status* status, int* flag){
778 if(status==MPI_STATUS_IGNORE){
782 *flag=simgrid::smpi::Status::cancelled(status);
786 MPI_Request PMPI_Request_f2c(MPI_Fint request){
787 return static_cast<MPI_Request>(simgrid::smpi::Request::f2c(request));
790 MPI_Fint PMPI_Request_c2f(MPI_Request request) {
791 return request->c2f();