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() & 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() & 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() & 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() & 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 simgrid::smpi::Status::empty(status);
315 status->MPI_SOURCE = MPI_PROC_NULL;
316 retval = MPI_SUCCESS;
317 } else if (src!=MPI_ANY_SOURCE && (src >= comm->group()->size() || src <0)){
318 retval = MPI_ERR_RANK;
319 } else if ((count < 0) || (buf==nullptr && count > 0)) {
320 retval = MPI_ERR_COUNT;
321 } else if (not datatype->is_valid()) {
322 retval = MPI_ERR_TYPE;
323 } else if(tag<0 && tag != MPI_ANY_TAG){
324 retval = MPI_ERR_TAG;
326 int my_proc_id = simgrid::s4u::this_actor::get_pid();
327 TRACE_smpi_comm_in(my_proc_id, __func__,
328 new simgrid::instr::Pt2PtTIData("recv", src,
329 datatype->is_replayable() ? count : count * datatype->size(),
330 tag, simgrid::smpi::Datatype::encode(datatype)));
332 simgrid::smpi::Request::recv(buf, count, datatype, src, tag, comm, status);
333 retval = MPI_SUCCESS;
335 // the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
336 if (status != MPI_STATUS_IGNORE) {
337 int src_traced = getPid(comm, status->MPI_SOURCE);
338 if (not TRACE_smpi_view_internals()) {
339 TRACE_smpi_recv(src_traced, my_proc_id, tag);
342 TRACE_smpi_comm_out(my_proc_id);
349 int PMPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm)
355 if (comm == MPI_COMM_NULL) {
356 retval = MPI_ERR_COMM;
357 } else if (dst == MPI_PROC_NULL) {
358 retval = MPI_SUCCESS;
359 } else if (dst >= comm->group()->size() || dst <0){
360 retval = MPI_ERR_RANK;
361 } else if ((count < 0) || (buf == nullptr && count > 0)) {
362 retval = MPI_ERR_COUNT;
363 } else if (not datatype->is_valid()) {
364 retval = MPI_ERR_TYPE;
365 } else if(tag < 0 && tag != MPI_ANY_TAG){
366 retval = MPI_ERR_TAG;
368 int my_proc_id = simgrid::s4u::this_actor::get_pid();
369 int dst_traced = getPid(comm, dst);
370 TRACE_smpi_comm_in(my_proc_id, __func__,
371 new simgrid::instr::Pt2PtTIData("send", dst,
372 datatype->is_replayable() ? count : count * datatype->size(),
373 tag, simgrid::smpi::Datatype::encode(datatype)));
374 if (not TRACE_smpi_view_internals()) {
375 TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, tag, count * datatype->size());
378 simgrid::smpi::Request::send(buf, count, datatype, dst, tag, comm);
379 retval = MPI_SUCCESS;
381 TRACE_smpi_comm_out(my_proc_id);
388 int PMPI_Ssend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm) {
393 if (comm == MPI_COMM_NULL) {
394 retval = MPI_ERR_COMM;
395 } else if (dst == MPI_PROC_NULL) {
396 retval = MPI_SUCCESS;
397 } else if (dst >= comm->group()->size() || dst <0){
398 retval = MPI_ERR_RANK;
399 } else if ((count < 0) || (buf==nullptr && count > 0)) {
400 retval = MPI_ERR_COUNT;
401 } else if (not datatype->is_valid()) {
402 retval = MPI_ERR_TYPE;
403 } else if(tag<0 && tag != MPI_ANY_TAG){
404 retval = MPI_ERR_TAG;
406 int my_proc_id = simgrid::s4u::this_actor::get_pid();
407 int dst_traced = getPid(comm, dst);
408 TRACE_smpi_comm_in(my_proc_id, __func__,
409 new simgrid::instr::Pt2PtTIData("Ssend", dst,
410 datatype->is_replayable() ? count : count * datatype->size(),
411 tag, simgrid::smpi::Datatype::encode(datatype)));
412 TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, tag, count * datatype->size());
414 simgrid::smpi::Request::ssend(buf, count, datatype, dst, tag, comm);
415 retval = MPI_SUCCESS;
417 TRACE_smpi_comm_out(my_proc_id);
424 int PMPI_Sendrecv(void* sendbuf, int sendcount, MPI_Datatype sendtype, int dst, int sendtag, void* recvbuf,
425 int recvcount, MPI_Datatype recvtype, int src, int recvtag, MPI_Comm comm, MPI_Status* status)
431 if (comm == MPI_COMM_NULL) {
432 retval = MPI_ERR_COMM;
433 } else if (not sendtype->is_valid() || not recvtype->is_valid()) {
434 retval = MPI_ERR_TYPE;
435 } else if (src == MPI_PROC_NULL) {
436 if(status!=MPI_STATUS_IGNORE){
437 simgrid::smpi::Status::empty(status);
438 status->MPI_SOURCE = MPI_PROC_NULL;
440 if(dst != MPI_PROC_NULL)
441 simgrid::smpi::Request::send(sendbuf, sendcount, sendtype, dst, sendtag, comm);
442 retval = MPI_SUCCESS;
443 }else if (dst == MPI_PROC_NULL){
444 simgrid::smpi::Request::recv(recvbuf, recvcount, recvtype, src, recvtag, comm, status);
445 retval = MPI_SUCCESS;
446 }else if (dst >= comm->group()->size() || dst <0 ||
447 (src!=MPI_ANY_SOURCE && (src >= comm->group()->size() || src <0))){
448 retval = MPI_ERR_RANK;
449 } else if ((sendcount < 0 || recvcount<0) ||
450 (sendbuf==nullptr && sendcount > 0) || (recvbuf==nullptr && recvcount>0)) {
451 retval = MPI_ERR_COUNT;
452 } else if((sendtag<0 && sendtag != MPI_ANY_TAG)||(recvtag<0 && recvtag != MPI_ANY_TAG)){
453 retval = MPI_ERR_TAG;
455 int my_proc_id = simgrid::s4u::this_actor::get_pid();
456 int dst_traced = getPid(comm, dst);
457 int src_traced = getPid(comm, src);
459 // FIXME: Hack the way to trace this one
460 std::vector<int>* dst_hack = new std::vector<int>;
461 std::vector<int>* src_hack = new std::vector<int>;
462 dst_hack->push_back(dst_traced);
463 src_hack->push_back(src_traced);
464 TRACE_smpi_comm_in(my_proc_id, __func__,
465 new simgrid::instr::VarCollTIData(
466 "sendRecv", -1, sendtype->is_replayable() ? sendcount : sendcount * sendtype->size(),
467 dst_hack, recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(), src_hack,
468 simgrid::smpi::Datatype::encode(sendtype), simgrid::smpi::Datatype::encode(recvtype)));
470 TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, sendtag, sendcount * sendtype->size());
472 simgrid::smpi::Request::sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf, recvcount, recvtype, src,
473 recvtag, comm, status);
474 retval = MPI_SUCCESS;
476 TRACE_smpi_recv(src_traced, my_proc_id, recvtag);
477 TRACE_smpi_comm_out(my_proc_id);
484 int PMPI_Sendrecv_replace(void* buf, int count, MPI_Datatype datatype, int dst, int sendtag, int src, int recvtag,
485 MPI_Comm comm, MPI_Status* status)
488 if (not datatype->is_valid()) {
490 } else if (count < 0) {
491 return MPI_ERR_COUNT;
493 int size = datatype->get_extent() * count;
494 void* recvbuf = xbt_new0(char, size);
495 retval = MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count, datatype, src, recvtag, comm, status);
496 if(retval==MPI_SUCCESS){
497 simgrid::smpi::Datatype::copy(recvbuf, count, datatype, buf, count, datatype);
505 int PMPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
509 if (request == nullptr || flag == nullptr) {
510 retval = MPI_ERR_ARG;
511 } else if (*request == MPI_REQUEST_NULL) {
513 simgrid::smpi::Status::empty(status);
514 retval = MPI_SUCCESS;
516 int my_proc_id = ((*request)->comm() != MPI_COMM_NULL) ? simgrid::s4u::this_actor::get_pid() : -1;
518 TRACE_smpi_testing_in(my_proc_id);
520 *flag = simgrid::smpi::Request::test(request,status);
522 TRACE_smpi_testing_out(my_proc_id);
523 retval = MPI_SUCCESS;
529 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag, MPI_Status * status)
534 if (index == nullptr || flag == nullptr) {
535 retval = MPI_ERR_ARG;
537 *flag = simgrid::smpi::Request::testany(count, requests, index, status);
538 retval = MPI_SUCCESS;
544 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses)
549 if (flag == nullptr) {
550 retval = MPI_ERR_ARG;
552 *flag = simgrid::smpi::Request::testall(count, requests, statuses);
553 retval = MPI_SUCCESS;
559 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
563 if (status == nullptr) {
564 retval = MPI_ERR_ARG;
565 } else if (comm == MPI_COMM_NULL) {
566 retval = MPI_ERR_COMM;
567 } else if (source == MPI_PROC_NULL) {
568 simgrid::smpi::Status::empty(status);
569 status->MPI_SOURCE = MPI_PROC_NULL;
570 retval = MPI_SUCCESS;
572 simgrid::smpi::Request::probe(source, tag, comm, status);
573 retval = MPI_SUCCESS;
579 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
583 if (flag == nullptr) {
584 retval = MPI_ERR_ARG;
585 } else if (comm == MPI_COMM_NULL) {
586 retval = MPI_ERR_COMM;
587 } else if (source == MPI_PROC_NULL) {
589 simgrid::smpi::Status::empty(status);
590 status->MPI_SOURCE = MPI_PROC_NULL;
591 retval = MPI_SUCCESS;
593 simgrid::smpi::Request::iprobe(source, tag, comm, flag, status);
594 retval = MPI_SUCCESS;
600 // TODO: cheinrich: Move declaration to other file? Rename this function - it's used for PMPI_Wait*?
601 static void trace_smpi_recv_helper(MPI_Request* request, MPI_Status* status);
602 static void trace_smpi_recv_helper(MPI_Request* request, MPI_Status* status)
604 MPI_Request req = *request;
605 if (req != MPI_REQUEST_NULL) { // Received requests become null
606 int src_traced = req->src();
607 // the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
608 int dst_traced = req->dst();
609 if (req->flags() & RECV) { // Is this request a wait for RECV?
610 if (src_traced == MPI_ANY_SOURCE)
611 src_traced = (status != MPI_STATUSES_IGNORE) ? req->comm()->group()->rank(status->MPI_SOURCE) : req->src();
612 TRACE_smpi_recv(src_traced, dst_traced, req->tag());
617 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
623 simgrid::smpi::Status::empty(status);
625 if (request == nullptr) {
626 retval = MPI_ERR_ARG;
627 } else if (*request == MPI_REQUEST_NULL) {
628 retval = MPI_SUCCESS;
630 int my_proc_id = (*request)->comm() != MPI_COMM_NULL
631 ? simgrid::s4u::this_actor::get_pid()
632 : -1; // TODO: cheinrich: Check if this correct or if it should be MPI_UNDEFINED
634 TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("wait"));
636 simgrid::smpi::Request::wait(request, status);
637 retval = MPI_SUCCESS;
639 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
640 TRACE_smpi_comm_out(my_proc_id);
641 trace_smpi_recv_helper(request, status);
648 int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * status)
650 if (index == nullptr)
658 int rank_traced = simgrid::s4u::this_actor::get_pid(); // FIXME: In PMPI_Wait, we check if the comm is null?
659 TRACE_smpi_comm_in(rank_traced, __func__, new simgrid::instr::CpuTIData("waitAny", static_cast<double>(count)));
661 *index = simgrid::smpi::Request::waitany(count, requests, status);
663 if(*index!=MPI_UNDEFINED){
664 trace_smpi_recv_helper(&requests[*index], status);
665 TRACE_smpi_comm_out(rank_traced);
672 int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
676 int rank_traced = simgrid::s4u::this_actor::get_pid(); // FIXME: In PMPI_Wait, we check if the comm is null?
677 TRACE_smpi_comm_in(rank_traced, __func__, new simgrid::instr::CpuTIData("waitAll", static_cast<double>(count)));
679 int retval = simgrid::smpi::Request::waitall(count, requests, status);
681 for (int i = 0; i < count; i++) {
682 trace_smpi_recv_helper(&requests[i], &status[i]);
684 TRACE_smpi_comm_out(rank_traced);
690 int PMPI_Waitsome(int incount, MPI_Request requests[], int *outcount, int *indices, MPI_Status status[])
695 if (outcount == nullptr) {
696 retval = MPI_ERR_ARG;
698 *outcount = simgrid::smpi::Request::waitsome(incount, requests, indices, status);
699 retval = MPI_SUCCESS;
705 int PMPI_Testsome(int incount, MPI_Request requests[], int* outcount, int* indices, MPI_Status status[])
710 if (outcount == nullptr) {
711 retval = MPI_ERR_ARG;
713 *outcount = simgrid::smpi::Request::testsome(incount, requests, indices, status);
714 retval = MPI_SUCCESS;
720 int PMPI_Cancel(MPI_Request* request)
725 if (*request == MPI_REQUEST_NULL) {
726 retval = MPI_ERR_REQUEST;
728 (*request)->cancel();
729 retval = MPI_SUCCESS;
735 int PMPI_Test_cancelled(MPI_Status* status, int* flag){
736 if(status==MPI_STATUS_IGNORE){
740 *flag=simgrid::smpi::Status::cancelled(status);
744 MPI_Request PMPI_Request_f2c(MPI_Fint request){
745 return static_cast<MPI_Request>(simgrid::smpi::Request::f2c(request));
748 MPI_Fint PMPI_Request_c2f(MPI_Request request) {
749 return request->c2f();