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 *request = MPI_REQUEST_NULL;
179 retval = MPI_SUCCESS;
185 int PMPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Request * request)
191 if (request == nullptr) {
192 retval = MPI_ERR_ARG;
193 } else if (comm == MPI_COMM_NULL) {
194 retval = MPI_ERR_COMM;
195 } else if (src == MPI_PROC_NULL) {
196 *request = MPI_REQUEST_NULL;
197 retval = MPI_SUCCESS;
198 } else if (src!=MPI_ANY_SOURCE && (src >= comm->group()->size() || src <0)){
199 retval = MPI_ERR_RANK;
200 } else if ((count < 0) || (buf==nullptr && count > 0)) {
201 retval = MPI_ERR_COUNT;
202 } else if (datatype==MPI_DATATYPE_NULL || not datatype->is_valid()) {
203 retval = MPI_ERR_TYPE;
204 } else if(tag<0 && tag != MPI_ANY_TAG){
205 retval = MPI_ERR_TAG;
208 int my_proc_id = simgrid::s4u::this_actor::get_pid();
210 TRACE_smpi_comm_in(my_proc_id, __func__,
211 new simgrid::instr::Pt2PtTIData("irecv", src,
212 datatype->is_replayable() ? count : count * datatype->size(),
213 tag, simgrid::smpi::Datatype::encode(datatype)));
215 *request = simgrid::smpi::Request::irecv(buf, count, datatype, src, tag, comm);
216 retval = MPI_SUCCESS;
218 TRACE_smpi_comm_out(my_proc_id);
222 if (retval != MPI_SUCCESS && request != nullptr)
223 *request = MPI_REQUEST_NULL;
228 int PMPI_Isend(const void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request * request)
233 if (request == nullptr) {
234 retval = MPI_ERR_ARG;
235 } else if (comm == MPI_COMM_NULL) {
236 retval = MPI_ERR_COMM;
237 } else if (dst == MPI_PROC_NULL) {
238 *request = MPI_REQUEST_NULL;
239 retval = MPI_SUCCESS;
240 } else if (dst >= comm->group()->size() || dst <0){
241 retval = MPI_ERR_RANK;
242 } else if ((count < 0) || (buf==nullptr && count > 0)) {
243 retval = MPI_ERR_COUNT;
244 } else if (datatype==MPI_DATATYPE_NULL || not datatype->is_valid()) {
245 retval = MPI_ERR_TYPE;
246 } else if(tag<0 && tag != MPI_ANY_TAG){
247 retval = MPI_ERR_TAG;
249 int my_proc_id = simgrid::s4u::this_actor::get_pid();
250 int trace_dst = getPid(comm, dst);
251 TRACE_smpi_comm_in(my_proc_id, __func__,
252 new simgrid::instr::Pt2PtTIData("isend", dst,
253 datatype->is_replayable() ? count : count * datatype->size(),
254 tag, simgrid::smpi::Datatype::encode(datatype)));
256 TRACE_smpi_send(my_proc_id, my_proc_id, trace_dst, tag, count * datatype->size());
258 *request = simgrid::smpi::Request::isend(buf, count, datatype, dst, tag, comm);
259 retval = MPI_SUCCESS;
261 TRACE_smpi_comm_out(my_proc_id);
265 if (retval != MPI_SUCCESS && request!=nullptr)
266 *request = MPI_REQUEST_NULL;
270 int PMPI_Irsend(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm,
271 MPI_Request* request)
273 return PMPI_Isend(buf, count, datatype, dst, tag, comm, request);
276 int PMPI_Issend(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request)
281 if (request == nullptr) {
282 retval = MPI_ERR_ARG;
283 } else if (comm == MPI_COMM_NULL) {
284 retval = MPI_ERR_COMM;
285 } else if (dst == MPI_PROC_NULL) {
286 *request = MPI_REQUEST_NULL;
287 retval = MPI_SUCCESS;
288 } else if (dst >= comm->group()->size() || dst <0){
289 retval = MPI_ERR_RANK;
290 } else if ((count < 0)|| (buf==nullptr && count > 0)) {
291 retval = MPI_ERR_COUNT;
292 } else if (datatype==MPI_DATATYPE_NULL || not datatype->is_valid()) {
293 retval = MPI_ERR_TYPE;
294 } else if(tag<0 && tag != MPI_ANY_TAG){
295 retval = MPI_ERR_TAG;
297 int my_proc_id = simgrid::s4u::this_actor::get_pid();
298 int trace_dst = getPid(comm, dst);
299 TRACE_smpi_comm_in(my_proc_id, __func__,
300 new simgrid::instr::Pt2PtTIData("ISsend", dst,
301 datatype->is_replayable() ? count : count * datatype->size(),
302 tag, simgrid::smpi::Datatype::encode(datatype)));
303 TRACE_smpi_send(my_proc_id, my_proc_id, trace_dst, tag, count * datatype->size());
305 *request = simgrid::smpi::Request::issend(buf, count, datatype, dst, tag, comm);
306 retval = MPI_SUCCESS;
308 TRACE_smpi_comm_out(my_proc_id);
312 if (retval != MPI_SUCCESS && request!=nullptr)
313 *request = MPI_REQUEST_NULL;
317 int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Status * status)
322 if (comm == MPI_COMM_NULL) {
323 retval = MPI_ERR_COMM;
324 } else if (src == MPI_PROC_NULL) {
325 if(status != MPI_STATUS_IGNORE){
326 simgrid::smpi::Status::empty(status);
327 status->MPI_SOURCE = MPI_PROC_NULL;
329 retval = MPI_SUCCESS;
330 } else if (src!=MPI_ANY_SOURCE && (src >= comm->group()->size() || src <0)){
331 retval = MPI_ERR_RANK;
332 } else if ((count < 0) || (buf==nullptr && count > 0)) {
333 retval = MPI_ERR_COUNT;
334 } else if (datatype==MPI_DATATYPE_NULL || not datatype->is_valid()) {
335 retval = MPI_ERR_TYPE;
336 } else if(tag<0 && tag != MPI_ANY_TAG){
337 retval = MPI_ERR_TAG;
339 int my_proc_id = simgrid::s4u::this_actor::get_pid();
340 TRACE_smpi_comm_in(my_proc_id, __func__,
341 new simgrid::instr::Pt2PtTIData("recv", src,
342 datatype->is_replayable() ? count : count * datatype->size(),
343 tag, simgrid::smpi::Datatype::encode(datatype)));
345 simgrid::smpi::Request::recv(buf, count, datatype, src, tag, comm, status);
346 retval = MPI_SUCCESS;
348 // the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
350 if (status != MPI_STATUS_IGNORE)
351 src_traced = getPid(comm, status->MPI_SOURCE);
353 src_traced = getPid(comm, src);
354 if (not TRACE_smpi_view_internals()) {
355 TRACE_smpi_recv(src_traced, my_proc_id, tag);
358 TRACE_smpi_comm_out(my_proc_id);
365 int PMPI_Send(const void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm)
371 if (comm == MPI_COMM_NULL) {
372 retval = MPI_ERR_COMM;
373 } else if (dst == MPI_PROC_NULL) {
374 retval = MPI_SUCCESS;
375 } else if (dst >= comm->group()->size() || dst <0){
376 retval = MPI_ERR_RANK;
377 } else if ((count < 0) || (buf == nullptr && count > 0)) {
378 retval = MPI_ERR_COUNT;
379 } else if (datatype==MPI_DATATYPE_NULL || not datatype->is_valid()) {
380 retval = MPI_ERR_TYPE;
381 } else if(tag < 0 && tag != MPI_ANY_TAG){
382 retval = MPI_ERR_TAG;
384 int my_proc_id = simgrid::s4u::this_actor::get_pid();
385 int dst_traced = getPid(comm, dst);
386 TRACE_smpi_comm_in(my_proc_id, __func__,
387 new simgrid::instr::Pt2PtTIData("send", dst,
388 datatype->is_replayable() ? count : count * datatype->size(),
389 tag, simgrid::smpi::Datatype::encode(datatype)));
390 if (not TRACE_smpi_view_internals()) {
391 TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, tag, count * datatype->size());
394 simgrid::smpi::Request::send(buf, count, datatype, dst, tag, comm);
395 retval = MPI_SUCCESS;
397 TRACE_smpi_comm_out(my_proc_id);
404 int PMPI_Rsend(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm)
406 return PMPI_Send(buf, count, datatype, dst, tag, comm);
409 int PMPI_Bsend(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 int bsend_buf_size = 0;
431 void* bsend_buf = nullptr;
432 smpi_process()->bsend_buffer(&bsend_buf, &bsend_buf_size);
433 int size = datatype->get_extent() * count;
434 if(bsend_buf==nullptr || bsend_buf_size < size + MPI_BSEND_OVERHEAD )
435 return MPI_ERR_BUFFER;
436 TRACE_smpi_comm_in(my_proc_id, __func__,
437 new simgrid::instr::Pt2PtTIData("bsend", dst,
438 datatype->is_replayable() ? count : count * datatype->size(),
439 tag, simgrid::smpi::Datatype::encode(datatype)));
440 if (not TRACE_smpi_view_internals()) {
441 TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, tag, count * datatype->size());
444 simgrid::smpi::Request::bsend(buf, count, datatype, dst, tag, comm);
445 retval = MPI_SUCCESS;
447 TRACE_smpi_comm_out(my_proc_id);
454 int PMPI_Ibsend(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request)
459 if (request == nullptr) {
460 retval = MPI_ERR_ARG;
461 } else if (comm == MPI_COMM_NULL) {
462 retval = MPI_ERR_COMM;
463 } else if (dst == MPI_PROC_NULL) {
464 *request = MPI_REQUEST_NULL;
465 retval = MPI_SUCCESS;
466 } else if (dst >= comm->group()->size() || dst <0){
467 retval = MPI_ERR_RANK;
468 } else if ((count < 0) || (buf==nullptr && count > 0)) {
469 retval = MPI_ERR_COUNT;
470 } else if (datatype==MPI_DATATYPE_NULL || not datatype->is_valid()) {
471 retval = MPI_ERR_TYPE;
472 } else if(tag<0 && tag != MPI_ANY_TAG){
473 retval = MPI_ERR_TAG;
475 int my_proc_id = simgrid::s4u::this_actor::get_pid();
476 int trace_dst = getPid(comm, dst);
477 int bsend_buf_size = 0;
478 void* bsend_buf = nullptr;
479 smpi_process()->bsend_buffer(&bsend_buf, &bsend_buf_size);
480 int size = datatype->get_extent() * count;
481 if(bsend_buf==nullptr || bsend_buf_size < size + MPI_BSEND_OVERHEAD )
482 return MPI_ERR_BUFFER;
483 TRACE_smpi_comm_in(my_proc_id, __func__,
484 new simgrid::instr::Pt2PtTIData("ibsend", dst,
485 datatype->is_replayable() ? count : count * datatype->size(),
486 tag, simgrid::smpi::Datatype::encode(datatype)));
488 TRACE_smpi_send(my_proc_id, my_proc_id, trace_dst, tag, count * datatype->size());
490 *request = simgrid::smpi::Request::ibsend(buf, count, datatype, dst, tag, comm);
491 retval = MPI_SUCCESS;
493 TRACE_smpi_comm_out(my_proc_id);
497 if (retval != MPI_SUCCESS && request!=nullptr)
498 *request = MPI_REQUEST_NULL;
502 int PMPI_Bsend_init(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request)
508 if (request == nullptr) {
509 retval = MPI_ERR_ARG;
510 } else if (comm == MPI_COMM_NULL) {
511 retval = MPI_ERR_COMM;
512 } else if (datatype==MPI_DATATYPE_NULL || not datatype->is_valid()) {
513 retval = MPI_ERR_TYPE;
514 } else if (dst == MPI_PROC_NULL) {
515 retval = MPI_SUCCESS;
517 int bsend_buf_size = 0;
518 void* bsend_buf = nullptr;
519 smpi_process()->bsend_buffer(&bsend_buf, &bsend_buf_size);
520 if( bsend_buf==nullptr || bsend_buf_size < datatype->get_extent() * count + MPI_BSEND_OVERHEAD ) {
521 retval = MPI_ERR_BUFFER;
523 *request = simgrid::smpi::Request::bsend_init(buf, count, datatype, dst, tag, comm);
524 retval = MPI_SUCCESS;
528 if (retval != MPI_SUCCESS && request != nullptr)
529 *request = MPI_REQUEST_NULL;
533 int PMPI_Ssend(const void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm) {
538 if (comm == MPI_COMM_NULL) {
539 retval = MPI_ERR_COMM;
540 } else if (dst == MPI_PROC_NULL) {
541 retval = MPI_SUCCESS;
542 } else if (dst >= comm->group()->size() || dst <0){
543 retval = MPI_ERR_RANK;
544 } else if ((count < 0) || (buf==nullptr && count > 0)) {
545 retval = MPI_ERR_COUNT;
546 } else if (datatype==MPI_DATATYPE_NULL || not datatype->is_valid()) {
547 retval = MPI_ERR_TYPE;
548 } else if(tag<0 && tag != MPI_ANY_TAG){
549 retval = MPI_ERR_TAG;
551 int my_proc_id = simgrid::s4u::this_actor::get_pid();
552 int dst_traced = getPid(comm, dst);
553 TRACE_smpi_comm_in(my_proc_id, __func__,
554 new simgrid::instr::Pt2PtTIData("Ssend", dst,
555 datatype->is_replayable() ? count : count * datatype->size(),
556 tag, simgrid::smpi::Datatype::encode(datatype)));
557 TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, tag, count * datatype->size());
559 simgrid::smpi::Request::ssend(buf, count, datatype, dst, tag, comm);
560 retval = MPI_SUCCESS;
562 TRACE_smpi_comm_out(my_proc_id);
569 int PMPI_Sendrecv(const void* sendbuf, int sendcount, MPI_Datatype sendtype, int dst, int sendtag, void* recvbuf,
570 int recvcount, MPI_Datatype recvtype, int src, int recvtag, MPI_Comm comm, MPI_Status* status)
576 if (comm == MPI_COMM_NULL) {
577 retval = MPI_ERR_COMM;
578 } else if (not sendtype->is_valid() || not recvtype->is_valid()) {
579 retval = MPI_ERR_TYPE;
580 } else if (src == MPI_PROC_NULL) {
581 if(status!=MPI_STATUS_IGNORE){
582 simgrid::smpi::Status::empty(status);
583 status->MPI_SOURCE = MPI_PROC_NULL;
585 if(dst != MPI_PROC_NULL)
586 simgrid::smpi::Request::send(sendbuf, sendcount, sendtype, dst, sendtag, comm);
587 retval = MPI_SUCCESS;
588 }else if (dst == MPI_PROC_NULL){
589 simgrid::smpi::Request::recv(recvbuf, recvcount, recvtype, src, recvtag, comm, status);
590 retval = MPI_SUCCESS;
591 }else if (dst >= comm->group()->size() || dst <0 ||
592 (src!=MPI_ANY_SOURCE && (src >= comm->group()->size() || src <0))){
593 retval = MPI_ERR_RANK;
594 } else if ((sendcount < 0 || recvcount<0) ||
595 (sendbuf==nullptr && sendcount > 0) || (recvbuf==nullptr && recvcount>0)) {
596 retval = MPI_ERR_COUNT;
597 } else if((sendtag<0 && sendtag != MPI_ANY_TAG)||(recvtag<0 && recvtag != MPI_ANY_TAG)){
598 retval = MPI_ERR_TAG;
600 int my_proc_id = simgrid::s4u::this_actor::get_pid();
601 int dst_traced = getPid(comm, dst);
602 int src_traced = getPid(comm, src);
604 // FIXME: Hack the way to trace this one
605 std::vector<int>* dst_hack = new std::vector<int>;
606 std::vector<int>* src_hack = new std::vector<int>;
607 dst_hack->push_back(dst_traced);
608 src_hack->push_back(src_traced);
609 TRACE_smpi_comm_in(my_proc_id, __func__,
610 new simgrid::instr::VarCollTIData(
611 "sendRecv", -1, sendtype->is_replayable() ? sendcount : sendcount * sendtype->size(),
612 dst_hack, recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(), src_hack,
613 simgrid::smpi::Datatype::encode(sendtype), simgrid::smpi::Datatype::encode(recvtype)));
615 TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, sendtag, sendcount * sendtype->size());
617 simgrid::smpi::Request::sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf, recvcount, recvtype, src,
618 recvtag, comm, status);
619 retval = MPI_SUCCESS;
621 TRACE_smpi_recv(src_traced, my_proc_id, recvtag);
622 TRACE_smpi_comm_out(my_proc_id);
629 int PMPI_Sendrecv_replace(void* buf, int count, MPI_Datatype datatype, int dst, int sendtag, int src, int recvtag,
630 MPI_Comm comm, MPI_Status* status)
633 if (datatype==MPI_DATATYPE_NULL || not datatype->is_valid()) {
635 } else if (count < 0) {
636 return MPI_ERR_COUNT;
638 int size = datatype->get_extent() * count;
639 void* recvbuf = xbt_new0(char, size);
640 retval = MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count, datatype, src, recvtag, comm, status);
641 if(retval==MPI_SUCCESS){
642 simgrid::smpi::Datatype::copy(recvbuf, count, datatype, buf, count, datatype);
650 int PMPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
654 if (request == nullptr || flag == nullptr) {
655 retval = MPI_ERR_ARG;
656 } else if (*request == MPI_REQUEST_NULL) {
657 if (status != MPI_STATUS_IGNORE){
659 simgrid::smpi::Status::empty(status);
661 retval = MPI_SUCCESS;
663 int my_proc_id = ((*request)->comm() != MPI_COMM_NULL) ? simgrid::s4u::this_actor::get_pid() : -1;
665 TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("test"));
667 retval = simgrid::smpi::Request::test(request,status, flag);
669 TRACE_smpi_comm_out(my_proc_id);
675 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag, MPI_Status * status)
680 if (index == nullptr || flag == nullptr) {
681 retval = MPI_ERR_ARG;
683 int my_proc_id = simgrid::s4u::this_actor::get_pid();
684 TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("testany"));
685 retval = simgrid::smpi::Request::testany(count, requests, index, flag, status);
686 TRACE_smpi_comm_out(my_proc_id);
692 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses)
697 if (flag == nullptr) {
698 retval = MPI_ERR_ARG;
700 int my_proc_id = simgrid::s4u::this_actor::get_pid();
701 TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("testall"));
702 retval = simgrid::smpi::Request::testall(count, requests, flag, statuses);
703 TRACE_smpi_comm_out(my_proc_id);
709 int PMPI_Testsome(int incount, MPI_Request requests[], int* outcount, int* indices, MPI_Status status[])
714 if (outcount == nullptr) {
715 retval = MPI_ERR_ARG;
717 int my_proc_id = simgrid::s4u::this_actor::get_pid();
718 TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("testsome"));
719 retval = simgrid::smpi::Request::testsome(incount, requests, outcount, indices, status);
720 TRACE_smpi_comm_out(my_proc_id);
726 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
730 if (comm == MPI_COMM_NULL) {
731 retval = MPI_ERR_COMM;
732 } else if (source == MPI_PROC_NULL) {
733 if (status != MPI_STATUS_IGNORE){
734 simgrid::smpi::Status::empty(status);
735 status->MPI_SOURCE = MPI_PROC_NULL;
737 retval = MPI_SUCCESS;
739 simgrid::smpi::Request::probe(source, tag, comm, status);
740 retval = MPI_SUCCESS;
746 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
750 if (flag == nullptr) {
751 retval = MPI_ERR_ARG;
752 } else if (comm == MPI_COMM_NULL) {
753 retval = MPI_ERR_COMM;
754 } else if (source == MPI_PROC_NULL) {
756 if (status != MPI_STATUS_IGNORE){
757 simgrid::smpi::Status::empty(status);
758 status->MPI_SOURCE = MPI_PROC_NULL;
760 retval = MPI_SUCCESS;
762 simgrid::smpi::Request::iprobe(source, tag, comm, flag, status);
763 retval = MPI_SUCCESS;
769 // TODO: cheinrich: Move declaration to other file? Rename this function - it's used for PMPI_Wait*?
770 static void trace_smpi_recv_helper(MPI_Request* request, MPI_Status* status);
771 static void trace_smpi_recv_helper(MPI_Request* request, MPI_Status* status)
773 MPI_Request req = *request;
774 if (req != MPI_REQUEST_NULL) { // Received requests become null
775 int src_traced = req->src();
776 // the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
777 int dst_traced = req->dst();
778 if (req->flags() & MPI_REQ_RECV) { // Is this request a wait for RECV?
779 if (src_traced == MPI_ANY_SOURCE)
780 src_traced = (status != MPI_STATUS_IGNORE) ? req->comm()->group()->rank(status->MPI_SOURCE) : req->src();
781 TRACE_smpi_recv(src_traced, dst_traced, req->tag());
786 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
792 simgrid::smpi::Status::empty(status);
794 if (request == nullptr) {
795 retval = MPI_ERR_ARG;
796 } else if (*request == MPI_REQUEST_NULL) {
797 retval = MPI_SUCCESS;
799 //for tracing, save the handle which might get overriden before we can use the helper on it
800 MPI_Request savedreq = *request;
801 if (savedreq != MPI_REQUEST_NULL && not(savedreq->flags() & MPI_REQ_FINISHED)
802 && not(savedreq->flags() & MPI_REQ_GENERALIZED))
803 savedreq->ref();//don't erase te handle in Request::wait, we'll need it later
805 savedreq = MPI_REQUEST_NULL;
807 int my_proc_id = (*request)->comm() != MPI_COMM_NULL
808 ? simgrid::s4u::this_actor::get_pid()
809 : -1; // TODO: cheinrich: Check if this correct or if it should be MPI_UNDEFINED
810 TRACE_smpi_comm_in(my_proc_id, __func__,
811 new simgrid::instr::WaitTIData((*request)->src(), (*request)->dst(), (*request)->tag()));
813 retval = simgrid::smpi::Request::wait(request, status);
815 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
816 TRACE_smpi_comm_out(my_proc_id);
817 trace_smpi_recv_helper(&savedreq, status);
818 if (savedreq != MPI_REQUEST_NULL)
819 simgrid::smpi::Request::unref(&savedreq);
826 int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * status)
828 if (index == nullptr)
835 //for tracing, save the handles which might get overriden before we can use the helper on it
836 std::vector<MPI_Request> savedreqs(requests, requests + count);
837 for (MPI_Request& req : savedreqs) {
838 if (req != MPI_REQUEST_NULL && not(req->flags() & MPI_REQ_FINISHED))
841 req = MPI_REQUEST_NULL;
844 int rank_traced = simgrid::s4u::this_actor::get_pid(); // FIXME: In PMPI_Wait, we check if the comm is null?
845 TRACE_smpi_comm_in(rank_traced, __func__, new simgrid::instr::CpuTIData("waitAny", count));
847 *index = simgrid::smpi::Request::waitany(count, requests, status);
849 if(*index!=MPI_UNDEFINED){
850 trace_smpi_recv_helper(&savedreqs[*index], status);
851 TRACE_smpi_comm_out(rank_traced);
854 for (MPI_Request& req : savedreqs)
855 if (req != MPI_REQUEST_NULL)
856 simgrid::smpi::Request::unref(&req);
862 int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
866 //for tracing, save the handles which might get overriden before we can use the helper on it
867 std::vector<MPI_Request> savedreqs(requests, requests + count);
868 for (MPI_Request& req : savedreqs) {
869 if (req != MPI_REQUEST_NULL && not(req->flags() & MPI_REQ_FINISHED))
872 req = MPI_REQUEST_NULL;
875 int rank_traced = simgrid::s4u::this_actor::get_pid(); // FIXME: In PMPI_Wait, we check if the comm is null?
876 TRACE_smpi_comm_in(rank_traced, __func__, new simgrid::instr::CpuTIData("waitall", count));
878 int retval = simgrid::smpi::Request::waitall(count, requests, status);
880 for (int i = 0; i < count; i++) {
881 trace_smpi_recv_helper(&savedreqs[i], status!=MPI_STATUSES_IGNORE ? &status[i]: MPI_STATUS_IGNORE);
883 TRACE_smpi_comm_out(rank_traced);
885 for (MPI_Request& req : savedreqs)
886 if (req != MPI_REQUEST_NULL)
887 simgrid::smpi::Request::unref(&req);
893 int PMPI_Waitsome(int incount, MPI_Request requests[], int *outcount, int *indices, MPI_Status status[])
898 if (outcount == nullptr) {
899 retval = MPI_ERR_ARG;
901 *outcount = simgrid::smpi::Request::waitsome(incount, requests, indices, status);
902 retval = MPI_SUCCESS;
908 int PMPI_Cancel(MPI_Request* request)
913 if (*request == MPI_REQUEST_NULL) {
914 retval = MPI_ERR_REQUEST;
916 (*request)->cancel();
917 retval = MPI_SUCCESS;
923 int PMPI_Test_cancelled(const MPI_Status* status, int* flag){
924 if(status==MPI_STATUS_IGNORE){
928 *flag=simgrid::smpi::Status::cancelled(status);
932 int PMPI_Status_set_cancelled(MPI_Status* status, int flag){
933 if(status==MPI_STATUS_IGNORE){
936 simgrid::smpi::Status::set_cancelled(status,flag);
940 int PMPI_Status_set_elements(MPI_Status* status, MPI_Datatype datatype, int count){
941 if(status==MPI_STATUS_IGNORE){
944 simgrid::smpi::Status::set_elements(status,datatype, count);
948 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){
949 return simgrid::smpi::Request::grequest_start(query_fn, free_fn,cancel_fn, extra_state, request);
952 int PMPI_Grequest_complete( MPI_Request request){
953 return simgrid::smpi::Request::grequest_complete(request);
956 int PMPI_Request_get_status( MPI_Request request, int *flag, MPI_Status *status){
957 if(request==MPI_REQUEST_NULL){
959 simgrid::smpi::Status::empty(status);
961 } else if (flag==NULL || status ==NULL){
964 return simgrid::smpi::Request::get_status(request,flag,status);
967 MPI_Request PMPI_Request_f2c(MPI_Fint request){
969 return MPI_REQUEST_NULL;
970 return static_cast<MPI_Request>(simgrid::smpi::Request::f2c(request));
973 MPI_Fint PMPI_Request_c2f(MPI_Request request) {
974 if(request==MPI_REQUEST_NULL)
976 return request->c2f();