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 retval = MPI_ERR_ARG;
179 simgrid::smpi::Request::unref(request);
180 retval = MPI_SUCCESS;
186 int PMPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Request * request)
192 if (request == nullptr) {
193 retval = MPI_ERR_ARG;
194 } else if (comm == MPI_COMM_NULL) {
195 retval = MPI_ERR_COMM;
196 } else if (src == MPI_PROC_NULL) {
197 *request = MPI_REQUEST_NULL;
198 retval = MPI_SUCCESS;
199 } else if (src!=MPI_ANY_SOURCE && (src >= comm->group()->size() || src <0)){
200 retval = MPI_ERR_RANK;
201 } else if ((count < 0) || (buf==nullptr && count > 0)) {
202 retval = MPI_ERR_COUNT;
203 } else if (datatype==MPI_DATATYPE_NULL || not datatype->is_valid()) {
204 retval = MPI_ERR_TYPE;
205 } else if(tag<0 && tag != MPI_ANY_TAG){
206 retval = MPI_ERR_TAG;
209 int my_proc_id = simgrid::s4u::this_actor::get_pid();
211 TRACE_smpi_comm_in(my_proc_id, __func__,
212 new simgrid::instr::Pt2PtTIData("irecv", src,
213 datatype->is_replayable() ? count : count * datatype->size(),
214 tag, simgrid::smpi::Datatype::encode(datatype)));
216 *request = simgrid::smpi::Request::irecv(buf, count, datatype, src, tag, comm);
217 retval = MPI_SUCCESS;
219 TRACE_smpi_comm_out(my_proc_id);
223 if (retval != MPI_SUCCESS && request != nullptr)
224 *request = MPI_REQUEST_NULL;
229 int PMPI_Isend(const void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request * request)
234 if (request == nullptr) {
235 retval = MPI_ERR_ARG;
236 } else if (comm == MPI_COMM_NULL) {
237 retval = MPI_ERR_COMM;
238 } else if (dst == MPI_PROC_NULL) {
239 *request = MPI_REQUEST_NULL;
240 retval = MPI_SUCCESS;
241 } else if (dst >= comm->group()->size() || dst <0){
242 retval = MPI_ERR_RANK;
243 } else if ((count < 0) || (buf==nullptr && count > 0)) {
244 retval = MPI_ERR_COUNT;
245 } else if (datatype==MPI_DATATYPE_NULL || not datatype->is_valid()) {
246 retval = MPI_ERR_TYPE;
247 } else if(tag<0 && tag != MPI_ANY_TAG){
248 retval = MPI_ERR_TAG;
250 int my_proc_id = simgrid::s4u::this_actor::get_pid();
251 int trace_dst = getPid(comm, dst);
252 TRACE_smpi_comm_in(my_proc_id, __func__,
253 new simgrid::instr::Pt2PtTIData("isend", dst,
254 datatype->is_replayable() ? count : count * datatype->size(),
255 tag, simgrid::smpi::Datatype::encode(datatype)));
257 TRACE_smpi_send(my_proc_id, my_proc_id, trace_dst, tag, count * datatype->size());
259 *request = simgrid::smpi::Request::isend(buf, count, datatype, dst, tag, comm);
260 retval = MPI_SUCCESS;
262 TRACE_smpi_comm_out(my_proc_id);
266 if (retval != MPI_SUCCESS && request!=nullptr)
267 *request = MPI_REQUEST_NULL;
271 int PMPI_Irsend(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm,
272 MPI_Request* request)
274 return PMPI_Isend(buf, count, datatype, dst, tag, comm, request);
277 int PMPI_Issend(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request)
282 if (request == nullptr) {
283 retval = MPI_ERR_ARG;
284 } else if (comm == MPI_COMM_NULL) {
285 retval = MPI_ERR_COMM;
286 } else if (dst == MPI_PROC_NULL) {
287 *request = MPI_REQUEST_NULL;
288 retval = MPI_SUCCESS;
289 } else if (dst >= comm->group()->size() || dst <0){
290 retval = MPI_ERR_RANK;
291 } else if ((count < 0)|| (buf==nullptr && count > 0)) {
292 retval = MPI_ERR_COUNT;
293 } else if (datatype==MPI_DATATYPE_NULL || not datatype->is_valid()) {
294 retval = MPI_ERR_TYPE;
295 } else if(tag<0 && tag != MPI_ANY_TAG){
296 retval = MPI_ERR_TAG;
298 int my_proc_id = simgrid::s4u::this_actor::get_pid();
299 int trace_dst = getPid(comm, dst);
300 TRACE_smpi_comm_in(my_proc_id, __func__,
301 new simgrid::instr::Pt2PtTIData("ISsend", dst,
302 datatype->is_replayable() ? count : count * datatype->size(),
303 tag, simgrid::smpi::Datatype::encode(datatype)));
304 TRACE_smpi_send(my_proc_id, my_proc_id, trace_dst, tag, count * datatype->size());
306 *request = simgrid::smpi::Request::issend(buf, count, datatype, dst, tag, comm);
307 retval = MPI_SUCCESS;
309 TRACE_smpi_comm_out(my_proc_id);
313 if (retval != MPI_SUCCESS && request!=nullptr)
314 *request = MPI_REQUEST_NULL;
318 int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Status * status)
323 if (comm == MPI_COMM_NULL) {
324 retval = MPI_ERR_COMM;
325 } else if (src == MPI_PROC_NULL) {
326 if(status != MPI_STATUS_IGNORE){
327 simgrid::smpi::Status::empty(status);
328 status->MPI_SOURCE = MPI_PROC_NULL;
330 retval = MPI_SUCCESS;
331 } else if (src!=MPI_ANY_SOURCE && (src >= comm->group()->size() || src <0)){
332 retval = MPI_ERR_RANK;
333 } else if ((count < 0) || (buf==nullptr && count > 0)) {
334 retval = MPI_ERR_COUNT;
335 } else if (datatype==MPI_DATATYPE_NULL || not datatype->is_valid()) {
336 retval = MPI_ERR_TYPE;
337 } else if(tag<0 && tag != MPI_ANY_TAG){
338 retval = MPI_ERR_TAG;
340 int my_proc_id = simgrid::s4u::this_actor::get_pid();
341 TRACE_smpi_comm_in(my_proc_id, __func__,
342 new simgrid::instr::Pt2PtTIData("recv", src,
343 datatype->is_replayable() ? count : count * datatype->size(),
344 tag, simgrid::smpi::Datatype::encode(datatype)));
346 simgrid::smpi::Request::recv(buf, count, datatype, src, tag, comm, status);
347 retval = MPI_SUCCESS;
349 // the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
351 if (status != MPI_STATUS_IGNORE)
352 src_traced = getPid(comm, status->MPI_SOURCE);
354 src_traced = getPid(comm, src);
355 if (not TRACE_smpi_view_internals()) {
356 TRACE_smpi_recv(src_traced, my_proc_id, tag);
359 TRACE_smpi_comm_out(my_proc_id);
366 int PMPI_Send(const void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm)
372 if (comm == MPI_COMM_NULL) {
373 retval = MPI_ERR_COMM;
374 } else if (dst == MPI_PROC_NULL) {
375 retval = MPI_SUCCESS;
376 } else if (dst >= comm->group()->size() || dst <0){
377 retval = MPI_ERR_RANK;
378 } else if ((count < 0) || (buf == nullptr && count > 0)) {
379 retval = MPI_ERR_COUNT;
380 } else if (datatype==MPI_DATATYPE_NULL || not datatype->is_valid()) {
381 retval = MPI_ERR_TYPE;
382 } else if(tag < 0 && tag != MPI_ANY_TAG){
383 retval = MPI_ERR_TAG;
385 int my_proc_id = simgrid::s4u::this_actor::get_pid();
386 int dst_traced = getPid(comm, dst);
387 TRACE_smpi_comm_in(my_proc_id, __func__,
388 new simgrid::instr::Pt2PtTIData("send", dst,
389 datatype->is_replayable() ? count : count * datatype->size(),
390 tag, simgrid::smpi::Datatype::encode(datatype)));
391 if (not TRACE_smpi_view_internals()) {
392 TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, tag, count * datatype->size());
395 simgrid::smpi::Request::send(buf, count, datatype, dst, tag, comm);
396 retval = MPI_SUCCESS;
398 TRACE_smpi_comm_out(my_proc_id);
405 int PMPI_Rsend(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm)
407 return PMPI_Send(buf, count, datatype, dst, tag, comm);
410 int PMPI_Ssend(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm) {
415 if (comm == MPI_COMM_NULL) {
416 retval = MPI_ERR_COMM;
417 } else if (dst == MPI_PROC_NULL) {
418 retval = MPI_SUCCESS;
419 } else if (dst >= comm->group()->size() || dst <0){
420 retval = MPI_ERR_RANK;
421 } else if ((count < 0) || (buf==nullptr && count > 0)) {
422 retval = MPI_ERR_COUNT;
423 } else if (datatype==MPI_DATATYPE_NULL || not datatype->is_valid()) {
424 retval = MPI_ERR_TYPE;
425 } else if(tag<0 && tag != MPI_ANY_TAG){
426 retval = MPI_ERR_TAG;
428 int my_proc_id = simgrid::s4u::this_actor::get_pid();
429 int dst_traced = getPid(comm, dst);
430 TRACE_smpi_comm_in(my_proc_id, __func__,
431 new simgrid::instr::Pt2PtTIData("Ssend", dst,
432 datatype->is_replayable() ? count : count * datatype->size(),
433 tag, simgrid::smpi::Datatype::encode(datatype)));
434 TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, tag, count * datatype->size());
436 simgrid::smpi::Request::ssend(buf, count, datatype, dst, tag, comm);
437 retval = MPI_SUCCESS;
439 TRACE_smpi_comm_out(my_proc_id);
446 int PMPI_Sendrecv(const void* sendbuf, int sendcount, MPI_Datatype sendtype, int dst, int sendtag, void* recvbuf,
447 int recvcount, MPI_Datatype recvtype, int src, int recvtag, MPI_Comm comm, MPI_Status* status)
453 if (comm == MPI_COMM_NULL) {
454 retval = MPI_ERR_COMM;
455 } else if (not sendtype->is_valid() || not recvtype->is_valid()) {
456 retval = MPI_ERR_TYPE;
457 } else if (src == MPI_PROC_NULL) {
458 if(status!=MPI_STATUS_IGNORE){
459 simgrid::smpi::Status::empty(status);
460 status->MPI_SOURCE = MPI_PROC_NULL;
462 if(dst != MPI_PROC_NULL)
463 simgrid::smpi::Request::send(sendbuf, sendcount, sendtype, dst, sendtag, comm);
464 retval = MPI_SUCCESS;
465 }else if (dst == MPI_PROC_NULL){
466 simgrid::smpi::Request::recv(recvbuf, recvcount, recvtype, src, recvtag, comm, status);
467 retval = MPI_SUCCESS;
468 }else if (dst >= comm->group()->size() || dst <0 ||
469 (src!=MPI_ANY_SOURCE && (src >= comm->group()->size() || src <0))){
470 retval = MPI_ERR_RANK;
471 } else if ((sendcount < 0 || recvcount<0) ||
472 (sendbuf==nullptr && sendcount > 0) || (recvbuf==nullptr && recvcount>0)) {
473 retval = MPI_ERR_COUNT;
474 } else if((sendtag<0 && sendtag != MPI_ANY_TAG)||(recvtag<0 && recvtag != MPI_ANY_TAG)){
475 retval = MPI_ERR_TAG;
477 int my_proc_id = simgrid::s4u::this_actor::get_pid();
478 int dst_traced = getPid(comm, dst);
479 int src_traced = getPid(comm, src);
481 // FIXME: Hack the way to trace this one
482 std::vector<int>* dst_hack = new std::vector<int>;
483 std::vector<int>* src_hack = new std::vector<int>;
484 dst_hack->push_back(dst_traced);
485 src_hack->push_back(src_traced);
486 TRACE_smpi_comm_in(my_proc_id, __func__,
487 new simgrid::instr::VarCollTIData(
488 "sendRecv", -1, sendtype->is_replayable() ? sendcount : sendcount * sendtype->size(),
489 dst_hack, recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(), src_hack,
490 simgrid::smpi::Datatype::encode(sendtype), simgrid::smpi::Datatype::encode(recvtype)));
492 TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, sendtag, sendcount * sendtype->size());
494 simgrid::smpi::Request::sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf, recvcount, recvtype, src,
495 recvtag, comm, status);
496 retval = MPI_SUCCESS;
498 TRACE_smpi_recv(src_traced, my_proc_id, recvtag);
499 TRACE_smpi_comm_out(my_proc_id);
506 int PMPI_Sendrecv_replace(void* buf, int count, MPI_Datatype datatype, int dst, int sendtag, int src, int recvtag,
507 MPI_Comm comm, MPI_Status* status)
510 if (datatype==MPI_DATATYPE_NULL || not datatype->is_valid()) {
512 } else if (count < 0) {
513 return MPI_ERR_COUNT;
515 int size = datatype->get_extent() * count;
516 void* recvbuf = xbt_new0(char, size);
517 retval = MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count, datatype, src, recvtag, comm, status);
518 if(retval==MPI_SUCCESS){
519 simgrid::smpi::Datatype::copy(recvbuf, count, datatype, buf, count, datatype);
527 int PMPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
531 if (request == nullptr || flag == nullptr) {
532 retval = MPI_ERR_ARG;
533 } else if (*request == MPI_REQUEST_NULL) {
534 if (status != MPI_STATUS_IGNORE){
536 simgrid::smpi::Status::empty(status);
538 retval = MPI_SUCCESS;
540 int my_proc_id = ((*request)->comm() != MPI_COMM_NULL) ? simgrid::s4u::this_actor::get_pid() : -1;
542 TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("test"));
544 retval = simgrid::smpi::Request::test(request,status, flag);
546 TRACE_smpi_comm_out(my_proc_id);
552 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag, MPI_Status * status)
557 if (index == nullptr || flag == nullptr) {
558 retval = MPI_ERR_ARG;
560 int my_proc_id = simgrid::s4u::this_actor::get_pid();
561 TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("testany"));
562 retval = simgrid::smpi::Request::testany(count, requests, index, flag, status);
563 TRACE_smpi_comm_out(my_proc_id);
569 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses)
574 if (flag == nullptr) {
575 retval = MPI_ERR_ARG;
577 int my_proc_id = simgrid::s4u::this_actor::get_pid();
578 TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("testall"));
579 retval = simgrid::smpi::Request::testall(count, requests, flag, statuses);
580 TRACE_smpi_comm_out(my_proc_id);
586 int PMPI_Testsome(int incount, MPI_Request requests[], int* outcount, int* indices, MPI_Status status[])
591 if (outcount == nullptr) {
592 retval = MPI_ERR_ARG;
594 int my_proc_id = simgrid::s4u::this_actor::get_pid();
595 TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("testsome"));
596 retval = simgrid::smpi::Request::testsome(incount, requests, outcount, indices, status);
597 TRACE_smpi_comm_out(my_proc_id);
603 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
607 if (comm == MPI_COMM_NULL) {
608 retval = MPI_ERR_COMM;
609 } else if (source == MPI_PROC_NULL) {
610 if (status != MPI_STATUS_IGNORE){
611 simgrid::smpi::Status::empty(status);
612 status->MPI_SOURCE = MPI_PROC_NULL;
614 retval = MPI_SUCCESS;
616 simgrid::smpi::Request::probe(source, tag, comm, status);
617 retval = MPI_SUCCESS;
623 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
627 if (flag == nullptr) {
628 retval = MPI_ERR_ARG;
629 } else if (comm == MPI_COMM_NULL) {
630 retval = MPI_ERR_COMM;
631 } else if (source == MPI_PROC_NULL) {
633 if (status != MPI_STATUS_IGNORE){
634 simgrid::smpi::Status::empty(status);
635 status->MPI_SOURCE = MPI_PROC_NULL;
637 retval = MPI_SUCCESS;
639 simgrid::smpi::Request::iprobe(source, tag, comm, flag, status);
640 retval = MPI_SUCCESS;
646 // TODO: cheinrich: Move declaration to other file? Rename this function - it's used for PMPI_Wait*?
647 static void trace_smpi_recv_helper(MPI_Request* request, MPI_Status* status);
648 static void trace_smpi_recv_helper(MPI_Request* request, MPI_Status* status)
650 MPI_Request req = *request;
651 if (req != MPI_REQUEST_NULL) { // Received requests become null
652 int src_traced = req->src();
653 // the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
654 int dst_traced = req->dst();
655 if (req->flags() & MPI_REQ_RECV) { // Is this request a wait for RECV?
656 if (src_traced == MPI_ANY_SOURCE)
657 src_traced = (status != MPI_STATUS_IGNORE) ? req->comm()->group()->rank(status->MPI_SOURCE) : req->src();
658 TRACE_smpi_recv(src_traced, dst_traced, req->tag());
663 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
669 simgrid::smpi::Status::empty(status);
671 if (request == nullptr) {
672 retval = MPI_ERR_ARG;
673 } else if (*request == MPI_REQUEST_NULL) {
674 retval = MPI_SUCCESS;
676 //for tracing, save the handle which might get overriden before we can use the helper on it
677 MPI_Request savedreq = *request;
678 if (savedreq != MPI_REQUEST_NULL && not(savedreq->flags() & MPI_REQ_FINISHED)
679 && not(savedreq->flags() & MPI_REQ_GENERALIZED))
680 savedreq->ref();//don't erase te handle in Request::wait, we'll need it later
682 savedreq = MPI_REQUEST_NULL;
684 int my_proc_id = (*request)->comm() != MPI_COMM_NULL
685 ? simgrid::s4u::this_actor::get_pid()
686 : -1; // TODO: cheinrich: Check if this correct or if it should be MPI_UNDEFINED
687 TRACE_smpi_comm_in(my_proc_id, __func__,
688 new simgrid::instr::WaitTIData((*request)->src(), (*request)->dst(), (*request)->tag()));
690 retval = simgrid::smpi::Request::wait(request, status);
692 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
693 TRACE_smpi_comm_out(my_proc_id);
694 trace_smpi_recv_helper(&savedreq, status);
695 if (savedreq != MPI_REQUEST_NULL)
696 simgrid::smpi::Request::unref(&savedreq);
703 int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * status)
705 if (index == nullptr)
712 //for tracing, save the handles which might get overriden before we can use the helper on it
713 std::vector<MPI_Request> savedreqs(requests, requests + count);
714 for (MPI_Request& req : savedreqs) {
715 if (req != MPI_REQUEST_NULL && not(req->flags() & MPI_REQ_FINISHED))
718 req = MPI_REQUEST_NULL;
721 int rank_traced = simgrid::s4u::this_actor::get_pid(); // FIXME: In PMPI_Wait, we check if the comm is null?
722 TRACE_smpi_comm_in(rank_traced, __func__, new simgrid::instr::CpuTIData("waitAny", count));
724 *index = simgrid::smpi::Request::waitany(count, requests, status);
726 if(*index!=MPI_UNDEFINED){
727 trace_smpi_recv_helper(&savedreqs[*index], status);
728 TRACE_smpi_comm_out(rank_traced);
731 for (MPI_Request& req : savedreqs)
732 if (req != MPI_REQUEST_NULL)
733 simgrid::smpi::Request::unref(&req);
739 int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
743 //for tracing, save the handles which might get overriden before we can use the helper on it
744 std::vector<MPI_Request> savedreqs(requests, requests + count);
745 for (MPI_Request& req : savedreqs) {
746 if (req != MPI_REQUEST_NULL && not(req->flags() & MPI_REQ_FINISHED))
749 req = MPI_REQUEST_NULL;
752 int rank_traced = simgrid::s4u::this_actor::get_pid(); // FIXME: In PMPI_Wait, we check if the comm is null?
753 TRACE_smpi_comm_in(rank_traced, __func__, new simgrid::instr::CpuTIData("waitall", count));
755 int retval = simgrid::smpi::Request::waitall(count, requests, status);
757 for (int i = 0; i < count; i++) {
758 trace_smpi_recv_helper(&savedreqs[i], status!=MPI_STATUSES_IGNORE ? &status[i]: MPI_STATUS_IGNORE);
760 TRACE_smpi_comm_out(rank_traced);
762 for (MPI_Request& req : savedreqs)
763 if (req != MPI_REQUEST_NULL)
764 simgrid::smpi::Request::unref(&req);
770 int PMPI_Waitsome(int incount, MPI_Request requests[], int *outcount, int *indices, MPI_Status status[])
775 if (outcount == nullptr) {
776 retval = MPI_ERR_ARG;
778 *outcount = simgrid::smpi::Request::waitsome(incount, requests, indices, status);
779 retval = MPI_SUCCESS;
785 int PMPI_Cancel(MPI_Request* request)
790 if (*request == MPI_REQUEST_NULL) {
791 retval = MPI_ERR_REQUEST;
793 (*request)->cancel();
794 retval = MPI_SUCCESS;
800 int PMPI_Test_cancelled(const MPI_Status* status, int* flag){
801 if(status==MPI_STATUS_IGNORE){
805 *flag=simgrid::smpi::Status::cancelled(status);
809 int PMPI_Status_set_cancelled(MPI_Status* status, int flag){
810 if(status==MPI_STATUS_IGNORE){
813 simgrid::smpi::Status::set_cancelled(status,flag);
817 int PMPI_Status_set_elements(MPI_Status* status, MPI_Datatype datatype, int count){
818 if(status==MPI_STATUS_IGNORE){
821 simgrid::smpi::Status::set_elements(status,datatype, count);
825 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){
826 return simgrid::smpi::Request::grequest_start(query_fn, free_fn,cancel_fn, extra_state, request);
829 int PMPI_Grequest_complete( MPI_Request request){
830 return simgrid::smpi::Request::grequest_complete(request);
834 int PMPI_Request_get_status( MPI_Request request, int *flag, MPI_Status *status){
835 if(request==MPI_REQUEST_NULL){
837 simgrid::smpi::Status::empty(status);
838 return MPI_ERR_REQUEST;
839 } else if (flag==NULL || status ==NULL){
842 return simgrid::smpi::Request::get_status(request,flag,status);
845 MPI_Request PMPI_Request_f2c(MPI_Fint request){
847 return MPI_REQUEST_NULL;
848 return static_cast<MPI_Request>(simgrid::smpi::Request::f2c(request));
851 MPI_Fint PMPI_Request_c2f(MPI_Request request) {
852 if(request==MPI_REQUEST_NULL)
854 return request->c2f();