1 /* Copyright (c) 2007-2021. The SimGrid Team. All rights reserved. */
3 /* This program is free software; you can redistribute it and/or modify it
4 * under the terms of the license (GNU LGPL) which comes with this package. */
7 #include "smpi_coll.hpp"
8 #include "smpi_comm.hpp"
9 #include "smpi_request.hpp"
10 #include "smpi_datatype_derived.hpp"
11 #include "smpi_op.hpp"
12 #include "src/smpi/include/smpi_actor.hpp"
16 XBT_LOG_EXTERNAL_DEFAULT_CATEGORY(smpi_pmpi);
18 static const void* smpi_get_in_place_buf(const void* inplacebuf, const void* otherbuf,
19 std::vector<unsigned char>& tmp_sendbuf, int count, MPI_Datatype datatype)
21 if (inplacebuf == MPI_IN_PLACE) {
22 tmp_sendbuf.resize(count * datatype->get_extent());
23 simgrid::smpi::Datatype::copy(otherbuf, count, datatype, tmp_sendbuf.data(), count, datatype);
24 return tmp_sendbuf.data();
29 /* PMPI User level calls */
31 int PMPI_Barrier(MPI_Comm comm)
33 return PMPI_Ibarrier(comm, MPI_REQUEST_IGNORED);
36 int PMPI_Ibarrier(MPI_Comm comm, MPI_Request *request)
41 const SmpiBenchGuard suspend_bench;
42 aid_t pid = simgrid::s4u::this_actor::get_pid();
43 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Barrier" : "PMPI_Ibarrier",
44 new simgrid::instr::NoOpTIData(request == MPI_REQUEST_IGNORED ? "barrier" : "ibarrier"));
45 if (request == MPI_REQUEST_IGNORED) {
46 simgrid::smpi::colls::barrier(comm);
47 // Barrier can be used to synchronize RMA calls. Finish all requests from comm before.
48 comm->finish_rma_calls();
50 simgrid::smpi::colls::ibarrier(comm, request);
52 TRACE_smpi_comm_out(pid);
56 int PMPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
58 return PMPI_Ibcast(buf, count, datatype, root, comm, MPI_REQUEST_IGNORED);
61 int PMPI_Ibcast(void* buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm, MPI_Request* request)
66 CHECK_TYPE(3, datatype)
67 CHECK_BUFFER(1, buf, count, datatype)
71 const SmpiBenchGuard suspend_bench;
72 aid_t pid = simgrid::s4u::this_actor::get_pid();
73 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Bcast" : "PMPI_Ibcast",
74 new simgrid::instr::CollTIData(request == MPI_REQUEST_IGNORED ? "bcast" : "ibcast", root, -1.0,
76 simgrid::smpi::Datatype::encode(datatype), ""));
77 if (comm->size() > 1) {
78 if (request == MPI_REQUEST_IGNORED)
79 simgrid::smpi::colls::bcast(buf, count, datatype, root, comm);
81 simgrid::smpi::colls::ibcast(buf, count, datatype, root, comm, request);
83 if (request != MPI_REQUEST_IGNORED)
84 *request = MPI_REQUEST_NULL;
87 TRACE_smpi_comm_out(pid);
91 int PMPI_Gather(const void *sendbuf, int sendcount, MPI_Datatype sendtype,void *recvbuf, int recvcount, MPI_Datatype recvtype,
92 int root, MPI_Comm comm){
93 return PMPI_Igather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm, MPI_REQUEST_IGNORED);
96 int PMPI_Igather(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
97 MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request* request)
101 int rank = comm->rank();
102 if(sendbuf != MPI_IN_PLACE){
103 CHECK_COUNT(2, sendcount)
104 CHECK_TYPE(3, sendtype)
105 CHECK_BUFFER(1,sendbuf, sendcount, sendtype)
109 CHECK_NOT_IN_PLACE_ROOT(4, recvbuf)
110 CHECK_TYPE(6, recvtype)
111 CHECK_COUNT(5, recvcount)
112 CHECK_BUFFER(4, recvbuf, recvcount, recvtype)
114 CHECK_NOT_IN_PLACE_ROOT(1, sendbuf)
119 const void* real_sendbuf = sendbuf;
120 int real_sendcount = sendcount;
121 MPI_Datatype real_sendtype = sendtype;
123 if (sendbuf == MPI_IN_PLACE) {
125 real_sendtype = recvtype;
126 } else if(recvtype->size() * recvcount != sendtype->size() * sendcount){
127 XBT_WARN("MPI_(I)Gather : received size at root differs from sent size : %zu vs %zu", recvtype->size() * recvcount , sendtype->size() * sendcount);
128 return MPI_ERR_TRUNCATE;
132 const SmpiBenchGuard suspend_bench;
134 aid_t pid = simgrid::s4u::this_actor::get_pid();
136 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Gather" : "PMPI_Igather",
137 new simgrid::instr::CollTIData(
138 request == MPI_REQUEST_IGNORED ? "gather" : "igather", root, -1.0,
139 real_sendcount, recvcount,
140 simgrid::smpi::Datatype::encode(real_sendtype), simgrid::smpi::Datatype::encode(recvtype)));
141 if (request == MPI_REQUEST_IGNORED)
142 simgrid::smpi::colls::gather(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcount, recvtype, root, comm);
144 simgrid::smpi::colls::igather(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcount, recvtype, root, comm,
147 TRACE_smpi_comm_out(pid);
151 int PMPI_Gatherv(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, const int *recvcounts, const int *displs,
152 MPI_Datatype recvtype, int root, MPI_Comm comm){
153 return PMPI_Igatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, root, comm, MPI_REQUEST_IGNORED);
156 int PMPI_Igatherv(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, const int* recvcounts, const int* displs,
157 MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request* request)
161 int rank = comm->rank();
162 if(sendbuf != MPI_IN_PLACE){
163 CHECK_TYPE(3, sendtype)
164 CHECK_COUNT(2, sendcount)
166 CHECK_BUFFER(1, sendbuf, sendcount, sendtype)
169 CHECK_NOT_IN_PLACE_ROOT(4, recvbuf)
170 CHECK_TYPE(6, recvtype)
171 CHECK_NULL(5, MPI_ERR_COUNT, recvcounts)
172 CHECK_NULL(6, MPI_ERR_ARG, displs)
174 CHECK_NOT_IN_PLACE_ROOT(1, sendbuf)
180 for (int i = 0; i < comm->size(); i++) {
181 CHECK_COUNT(5, recvcounts[i])
182 CHECK_BUFFER(4,recvbuf,recvcounts[i], recvtype)
186 const SmpiBenchGuard suspend_bench;
187 const void* real_sendbuf = sendbuf;
188 int real_sendcount = sendcount;
189 MPI_Datatype real_sendtype = sendtype;
190 if ((rank == root) && (sendbuf == MPI_IN_PLACE)) {
192 real_sendtype = recvtype;
195 aid_t pid = simgrid::s4u::this_actor::get_pid();
197 auto trace_recvcounts = std::make_shared<std::vector<int>>();
199 trace_recvcounts->insert(trace_recvcounts->end(), &recvcounts[0], &recvcounts[comm->size()]);
200 else //this is not significant outside of root, put 0 as we don't know if recvcounts is initialized
201 trace_recvcounts->insert(trace_recvcounts->end(), comm->size(), 0);
203 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Gatherv" : "PMPI_Igatherv",
204 new simgrid::instr::VarCollTIData(
205 request == MPI_REQUEST_IGNORED ? "gatherv" : "igatherv", root,
207 nullptr, -1, trace_recvcounts, simgrid::smpi::Datatype::encode(real_sendtype),
208 simgrid::smpi::Datatype::encode(recvtype)));
209 if (request == MPI_REQUEST_IGNORED)
210 simgrid::smpi::colls::gatherv(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcounts, displs, recvtype,
213 simgrid::smpi::colls::igatherv(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcounts, displs, recvtype,
214 root, comm, request);
216 TRACE_smpi_comm_out(pid);
220 int PMPI_Allgather(const void *sendbuf, int sendcount, MPI_Datatype sendtype,
221 void *recvbuf, int recvcount, MPI_Datatype recvtype, MPI_Comm comm){
222 return PMPI_Iallgather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm, MPI_REQUEST_IGNORED);
225 int PMPI_Iallgather(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
226 MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request)
231 int rank = comm->rank();
232 CHECK_NOT_IN_PLACE(4, recvbuf)
233 if(sendbuf != MPI_IN_PLACE){
234 CHECK_COUNT(2, sendcount)
235 CHECK_TYPE(3, sendtype)
237 CHECK_TYPE(6, recvtype)
238 CHECK_COUNT(5, recvcount)
239 CHECK_BUFFER(1, sendbuf, sendcount, sendtype)
240 CHECK_BUFFER(4, recvbuf, recvcount, recvtype)
243 if (sendbuf == MPI_IN_PLACE) {
244 sendbuf = static_cast<char*>(recvbuf) + recvtype->get_extent() * recvcount * comm->rank();
245 sendcount = recvcount;
249 if(recvtype->size() * recvcount != sendtype->size() * sendcount){
250 XBT_WARN("MPI_(I)Allgather : received size from each process differs from sent size : %zu vs %zu", recvtype->size() * recvcount, sendtype->size() * sendcount);
251 return MPI_ERR_TRUNCATE;
254 const SmpiBenchGuard suspend_bench;
256 aid_t pid = simgrid::s4u::this_actor::get_pid();
258 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Allgather" : "PMPI_Iallggather",
259 new simgrid::instr::CollTIData(
260 request == MPI_REQUEST_IGNORED ? "allgather" : "iallgather", -1, -1.0,
261 sendcount, recvcount,
262 simgrid::smpi::Datatype::encode(sendtype), simgrid::smpi::Datatype::encode(recvtype)));
263 if (request == MPI_REQUEST_IGNORED)
264 simgrid::smpi::colls::allgather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
266 simgrid::smpi::colls::iallgather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm, request);
268 TRACE_smpi_comm_out(pid);
272 int PMPI_Allgatherv(const void *sendbuf, int sendcount, MPI_Datatype sendtype,
273 void *recvbuf, const int *recvcounts, const int *displs, MPI_Datatype recvtype, MPI_Comm comm){
274 return PMPI_Iallgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm, MPI_REQUEST_IGNORED);
277 int PMPI_Iallgatherv(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, const int* recvcounts, const int* displs,
278 MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request)
283 int rank = comm->rank();
284 if(sendbuf != MPI_IN_PLACE)
285 CHECK_TYPE(3, sendtype)
286 CHECK_TYPE(6, recvtype)
287 CHECK_NULL(5, MPI_ERR_COUNT, recvcounts)
288 CHECK_NULL(6, MPI_ERR_ARG, displs)
289 if(sendbuf != MPI_IN_PLACE){
290 CHECK_COUNT(2, sendcount)
291 CHECK_BUFFER(1, sendbuf, sendcount, sendtype)
294 CHECK_NOT_IN_PLACE(4, recvbuf)
295 for (int i = 0; i < comm->size(); i++) {
296 CHECK_COUNT(5, recvcounts[i])
297 CHECK_BUFFER(4, recvbuf, recvcounts[i], recvtype)
300 const SmpiBenchGuard suspend_bench;
301 if (sendbuf == MPI_IN_PLACE) {
302 sendbuf = static_cast<char*>(recvbuf) + recvtype->get_extent() * displs[comm->rank()];
303 sendcount = recvcounts[comm->rank()];
306 aid_t pid = simgrid::s4u::this_actor::get_pid();
308 auto trace_recvcounts = std::make_shared<std::vector<int>>();
309 trace_recvcounts->insert(trace_recvcounts->end(), &recvcounts[0], &recvcounts[comm->size()]);
312 pid, request == MPI_REQUEST_IGNORED ? "PMPI_Allgatherv" : "PMPI_Iallgatherv",
313 new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "allgatherv" : "iallgatherv", -1,
315 -1, trace_recvcounts, simgrid::smpi::Datatype::encode(sendtype),
316 simgrid::smpi::Datatype::encode(recvtype)));
317 if (request == MPI_REQUEST_IGNORED)
318 simgrid::smpi::colls::allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm);
320 simgrid::smpi::colls::iallgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm,
323 TRACE_smpi_comm_out(pid);
327 int PMPI_Scatter(const void *sendbuf, int sendcount, MPI_Datatype sendtype,
328 void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm){
329 return PMPI_Iscatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm, MPI_REQUEST_IGNORED);
332 int PMPI_Iscatter(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
333 MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request* request)
337 int rank = comm->rank();
340 CHECK_NOT_IN_PLACE_ROOT(1, sendbuf)
341 CHECK_COUNT(2, sendcount)
342 CHECK_TYPE(3, sendtype)
343 CHECK_BUFFER(1, sendbuf, sendcount, sendtype)
345 CHECK_NOT_IN_PLACE_ROOT(4, recvbuf)
347 if(recvbuf != MPI_IN_PLACE){
348 CHECK_COUNT(5, recvcount)
349 CHECK_TYPE(6, recvtype)
350 CHECK_BUFFER(4, recvbuf, recvcount, recvtype)
355 if (recvbuf == MPI_IN_PLACE) {
357 recvcount = sendcount;
360 if((rank == root) && (recvtype->size() * recvcount != sendtype->size() * sendcount)){
361 XBT_WARN("MPI_(I)Scatter : sent size to each process differs from receive size");
362 return MPI_ERR_TRUNCATE;
365 const SmpiBenchGuard suspend_bench;
367 aid_t pid = simgrid::s4u::this_actor::get_pid();
369 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Scatter" : "PMPI_Iscatter",
370 new simgrid::instr::CollTIData(
371 request == MPI_REQUEST_IGNORED ? "scatter" : "iscatter", root, -1.0,
372 sendcount, recvcount,
373 simgrid::smpi::Datatype::encode(sendtype), simgrid::smpi::Datatype::encode(recvtype)));
374 if (request == MPI_REQUEST_IGNORED)
375 simgrid::smpi::colls::scatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm);
377 simgrid::smpi::colls::iscatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm, request);
379 TRACE_smpi_comm_out(pid);
383 int PMPI_Scatterv(const void *sendbuf, const int *sendcounts, const int *displs,
384 MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm){
385 return PMPI_Iscatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm, MPI_REQUEST_IGNORED);
388 int PMPI_Iscatterv(const void* sendbuf, const int* sendcounts, const int* displs, MPI_Datatype sendtype, void* recvbuf, int recvcount,
389 MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request* request)
393 int rank = comm->rank();
394 if(recvbuf != MPI_IN_PLACE){
395 CHECK_COUNT(5, recvcount)
396 CHECK_TYPE(7, recvtype)
397 CHECK_BUFFER(4, recvbuf, recvcount, recvtype)
403 CHECK_NOT_IN_PLACE_ROOT(1, sendbuf)
404 CHECK_NULL(2, MPI_ERR_COUNT, sendcounts)
405 CHECK_NULL(3, MPI_ERR_ARG, displs)
406 CHECK_TYPE(4, sendtype)
407 for (int i = 0; i < comm->size(); i++){
408 CHECK_COUNT(2, sendcounts[i])
409 CHECK_BUFFER(1, sendbuf, sendcounts[i], sendtype)
411 if (recvbuf == MPI_IN_PLACE) {
413 recvcount = sendcounts[rank];
416 CHECK_NOT_IN_PLACE_ROOT(4, recvbuf)
419 const SmpiBenchGuard suspend_bench;
421 aid_t pid = simgrid::s4u::this_actor::get_pid();
423 auto trace_sendcounts = std::make_shared<std::vector<int>>();
425 trace_sendcounts->insert(trace_sendcounts->end(), &sendcounts[0], &sendcounts[comm->size()]);
426 else //this is not significant outside of root, put 0 as we don't know if sendcounts is initialized
427 trace_sendcounts->insert(trace_sendcounts->end(), comm->size(), 0);
430 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Scatterv" : "PMPI_Iscatterv",
431 new simgrid::instr::VarCollTIData(
432 request == MPI_REQUEST_IGNORED ? "scatterv" : "iscatterv", root, -1,
433 trace_sendcounts, recvcount,
434 nullptr, simgrid::smpi::Datatype::encode(sendtype),
435 simgrid::smpi::Datatype::encode(recvtype)));
436 if (request == MPI_REQUEST_IGNORED)
437 simgrid::smpi::colls::scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm);
439 simgrid::smpi::colls::iscatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm,
442 TRACE_smpi_comm_out(pid);
446 int PMPI_Reduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
448 return PMPI_Ireduce(sendbuf, recvbuf, count, datatype, op, root, comm, MPI_REQUEST_IGNORED);
451 int PMPI_Ireduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm, MPI_Request* request)
455 int rank = comm->rank();
456 CHECK_TYPE(4, datatype)
457 CHECK_COUNT(3, count)
458 CHECK_BUFFER(1, sendbuf, count, datatype)
461 CHECK_NOT_IN_PLACE(2, recvbuf)
462 CHECK_BUFFER(5, recvbuf, count, datatype)
464 CHECK_OP(5, op, datatype)
468 const SmpiBenchGuard suspend_bench;
469 aid_t pid = simgrid::s4u::this_actor::get_pid();
471 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Reduce" : "PMPI_Ireduce",
472 new simgrid::instr::CollTIData(request == MPI_REQUEST_IGNORED ? "reduce" : "ireduce", root, 0,
474 simgrid::smpi::Datatype::encode(datatype), ""));
475 if (request == MPI_REQUEST_IGNORED)
476 simgrid::smpi::colls::reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
478 simgrid::smpi::colls::ireduce(sendbuf, recvbuf, count, datatype, op, root, comm, request);
480 TRACE_smpi_comm_out(pid);
484 int PMPI_Reduce_local(const void* inbuf, void* inoutbuf, int count, MPI_Datatype datatype, MPI_Op op)
488 CHECK_TYPE(4, datatype)
489 CHECK_COUNT(3, count)
490 CHECK_BUFFER(1, inbuf, count, datatype)
491 CHECK_BUFFER(2, inoutbuf, count, datatype)
492 CHECK_OP(5, op, datatype)
494 const SmpiBenchGuard suspend_bench;
495 op->apply(inbuf, inoutbuf, &count, datatype);
499 int PMPI_Allreduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
501 return PMPI_Iallreduce(sendbuf, recvbuf, count, datatype, op, comm, MPI_REQUEST_IGNORED);
504 int PMPI_Iallreduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request *request)
509 int rank = comm->rank();
510 CHECK_NOT_IN_PLACE(2, recvbuf)
511 CHECK_TYPE(4, datatype)
512 CHECK_OP(5, op, datatype)
513 CHECK_COUNT(3, count)
514 CHECK_BUFFER(1, sendbuf, count, datatype)
515 CHECK_BUFFER(2, recvbuf, count, datatype)
518 const SmpiBenchGuard suspend_bench;
519 std::vector<unsigned char> tmp_sendbuf;
520 const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, count, datatype);
522 aid_t pid = simgrid::s4u::this_actor::get_pid();
524 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Allreduce" : "PMPI_Iallreduce",
525 new simgrid::instr::CollTIData(request == MPI_REQUEST_IGNORED ? "allreduce" : "iallreduce", -1, 0,
527 simgrid::smpi::Datatype::encode(datatype), ""));
529 if (request == MPI_REQUEST_IGNORED)
530 simgrid::smpi::colls::allreduce(real_sendbuf, recvbuf, count, datatype, op, comm);
532 simgrid::smpi::colls::iallreduce(real_sendbuf, recvbuf, count, datatype, op, comm, request);
534 TRACE_smpi_comm_out(pid);
538 int PMPI_Scan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
540 return PMPI_Iscan(sendbuf, recvbuf, count, datatype, op, comm, MPI_REQUEST_IGNORED);
543 int PMPI_Iscan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request* request)
548 CHECK_TYPE(4, datatype)
549 CHECK_COUNT(3, count)
550 CHECK_BUFFER(1,sendbuf,count, datatype)
551 CHECK_BUFFER(2,recvbuf,count, datatype)
553 CHECK_OP(5, op, datatype)
555 const SmpiBenchGuard suspend_bench;
556 aid_t pid = simgrid::s4u::this_actor::get_pid();
557 std::vector<unsigned char> tmp_sendbuf;
558 const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, count, datatype);
560 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Scan" : "PMPI_Iscan",
561 new simgrid::instr::CollTIData(request == MPI_REQUEST_IGNORED ? "scan" : "iscan", -1, 0.0,
562 count, 0, simgrid::smpi::Datatype::encode(datatype), ""));
565 if (request == MPI_REQUEST_IGNORED)
566 retval = simgrid::smpi::colls::scan(real_sendbuf, recvbuf, count, datatype, op, comm);
568 retval = simgrid::smpi::colls::iscan(real_sendbuf, recvbuf, count, datatype, op, comm, request);
570 TRACE_smpi_comm_out(pid);
574 int PMPI_Exscan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
576 return PMPI_Iexscan(sendbuf, recvbuf, count, datatype, op, comm, MPI_REQUEST_IGNORED);
579 int PMPI_Iexscan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request* request){
583 CHECK_TYPE(4, datatype)
584 CHECK_COUNT(3, count)
585 CHECK_BUFFER(1, sendbuf, count, datatype)
586 CHECK_BUFFER(2, recvbuf, count, datatype)
588 CHECK_OP(5, op, datatype)
590 const SmpiBenchGuard suspend_bench;
591 aid_t pid = simgrid::s4u::this_actor::get_pid();
592 std::vector<unsigned char> tmp_sendbuf;
593 const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, count, datatype);
595 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Exscan" : "PMPI_Iexscan",
596 new simgrid::instr::CollTIData(request == MPI_REQUEST_IGNORED ? "exscan" : "iexscan", -1, 0.0,
597 count, 0, simgrid::smpi::Datatype::encode(datatype), ""));
600 if (request == MPI_REQUEST_IGNORED)
601 retval = simgrid::smpi::colls::exscan(real_sendbuf, recvbuf, count, datatype, op, comm);
603 retval = simgrid::smpi::colls::iexscan(real_sendbuf, recvbuf, count, datatype, op, comm, request);
605 TRACE_smpi_comm_out(pid);
609 int PMPI_Reduce_scatter(const void *sendbuf, void *recvbuf, const int *recvcounts, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
611 return PMPI_Ireduce_scatter(sendbuf, recvbuf, recvcounts, datatype, op, comm, MPI_REQUEST_IGNORED);
614 int PMPI_Ireduce_scatter(const void *sendbuf, void *recvbuf, const int *recvcounts, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request *request)
619 int rank = comm->rank();
620 CHECK_NOT_IN_PLACE(2, recvbuf)
621 CHECK_TYPE(4, datatype)
622 CHECK_NULL(3, MPI_ERR_COUNT, recvcounts)
624 CHECK_OP(5, op, datatype)
625 for (int i = 0; i < comm->size(); i++) {
626 CHECK_COUNT(3, recvcounts[i])
627 CHECK_BUFFER(1, sendbuf, recvcounts[i], datatype)
628 CHECK_BUFFER(2, recvbuf, recvcounts[i], datatype)
631 const SmpiBenchGuard suspend_bench;
632 aid_t pid = simgrid::s4u::this_actor::get_pid();
633 auto trace_recvcounts = std::make_shared<std::vector<int>>();
634 trace_recvcounts->insert(trace_recvcounts->end(), &recvcounts[0], &recvcounts[comm->size()]);
638 for (int i = 0; i < comm->size(); i++) { // copy data to avoid bad free
639 totalcount += recvcounts[i];
641 std::vector<unsigned char> tmp_sendbuf;
642 const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, totalcount, datatype);
644 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Reduce_scatter" : "PMPI_Ireduce_scatter",
645 new simgrid::instr::VarCollTIData(
646 request == MPI_REQUEST_IGNORED ? "reducescatter" : "ireducescatter", -1, -1, nullptr,
647 -1 , trace_recvcounts, std::to_string(0), simgrid::smpi::Datatype::encode(datatype)));
649 if (request == MPI_REQUEST_IGNORED)
650 simgrid::smpi::colls::reduce_scatter(real_sendbuf, recvbuf, recvcounts, datatype, op, comm);
652 simgrid::smpi::colls::ireduce_scatter(real_sendbuf, recvbuf, recvcounts, datatype, op, comm, request);
654 TRACE_smpi_comm_out(pid);
658 int PMPI_Reduce_scatter_block(const void *sendbuf, void *recvbuf, int recvcount,
659 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
661 return PMPI_Ireduce_scatter_block(sendbuf, recvbuf, recvcount, datatype, op, comm, MPI_REQUEST_IGNORED);
664 int PMPI_Ireduce_scatter_block(const void* sendbuf, void* recvbuf, int recvcount, MPI_Datatype datatype, MPI_Op op,
665 MPI_Comm comm, MPI_Request* request)
670 CHECK_TYPE(4, datatype)
671 CHECK_COUNT(3, recvcount)
672 CHECK_BUFFER(1, sendbuf, recvcount, datatype)
673 CHECK_BUFFER(2, recvbuf, recvcount, datatype)
675 CHECK_OP(5, op, datatype)
677 const SmpiBenchGuard suspend_bench;
678 int count = comm->size();
680 aid_t pid = simgrid::s4u::this_actor::get_pid();
681 auto trace_recvcounts = std::make_shared<std::vector<int>>(recvcount);
683 std::vector<unsigned char> tmp_sendbuf;
684 const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, recvcount * count, datatype);
687 pid, request == MPI_REQUEST_IGNORED ? "PMPI_Reduce_scatter_block" : "PMPI_Ireduce_scatter_block",
688 new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "reducescatter" : "ireducescatter", -1, -1,
689 nullptr, -1, trace_recvcounts, simgrid::smpi::Datatype::encode(datatype), ""));
691 std::vector<int> recvcounts(count);
692 for (int i = 0; i < count; i++)
693 recvcounts[i] = recvcount;
694 if (request == MPI_REQUEST_IGNORED)
695 simgrid::smpi::colls::reduce_scatter(real_sendbuf, recvbuf, recvcounts.data(), datatype, op, comm);
697 simgrid::smpi::colls::ireduce_scatter(real_sendbuf, recvbuf, recvcounts.data(), datatype, op, comm, request);
699 TRACE_smpi_comm_out(pid);
703 int PMPI_Alltoall(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
704 MPI_Datatype recvtype, MPI_Comm comm){
705 return PMPI_Ialltoall(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm, MPI_REQUEST_IGNORED);
708 int PMPI_Ialltoall(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
709 MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request)
714 if(sendbuf != MPI_IN_PLACE){
715 CHECK_TYPE(3, sendtype)
716 CHECK_COUNT(2, sendcount)
717 CHECK_BUFFER(1, sendbuf, sendcount, sendtype)
719 CHECK_TYPE(6, recvtype)
720 CHECK_COUNT(5, recvcount)
721 CHECK_COUNT(5, recvcount)
722 CHECK_BUFFER(4, recvbuf, recvcount, recvtype)
725 aid_t pid = simgrid::s4u::this_actor::get_pid();
726 int real_sendcount = sendcount;
727 MPI_Datatype real_sendtype = sendtype;
729 std::vector<unsigned char> tmp_sendbuf;
730 const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, recvcount * comm->size(), recvtype);
732 if (sendbuf == MPI_IN_PLACE) {
733 real_sendcount = recvcount;
734 real_sendtype = recvtype;
737 if(recvtype->size() * recvcount != real_sendtype->size() * real_sendcount){
738 XBT_WARN("MPI_(I)Alltoall : receive size from each process differs from sent size : %zu vs %zu", recvtype->size() * recvcount, real_sendtype->size() * real_sendcount);
739 return MPI_ERR_TRUNCATE;
742 const SmpiBenchGuard suspend_bench;
744 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Alltoall" : "PMPI_Ialltoall",
745 new simgrid::instr::CollTIData(
746 request == MPI_REQUEST_IGNORED ? "alltoall" : "ialltoall", -1, -1.0,
747 real_sendcount, recvcount,
748 simgrid::smpi::Datatype::encode(real_sendtype), simgrid::smpi::Datatype::encode(recvtype)));
750 if (request == MPI_REQUEST_IGNORED)
752 simgrid::smpi::colls::alltoall(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcount, recvtype, comm);
754 retval = simgrid::smpi::colls::ialltoall(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcount, recvtype,
757 TRACE_smpi_comm_out(pid);
761 int PMPI_Alltoallv(const void* sendbuf, const int* sendcounts, const int* senddispls, MPI_Datatype sendtype, void* recvbuf,
762 const int* recvcounts, const int* recvdispls, MPI_Datatype recvtype, MPI_Comm comm)
764 return PMPI_Ialltoallv(sendbuf, sendcounts, senddispls, sendtype, recvbuf, recvcounts, recvdispls, recvtype, comm, MPI_REQUEST_IGNORED);
767 int PMPI_Ialltoallv(const void* sendbuf, const int* sendcounts, const int* senddispls, MPI_Datatype sendtype, void* recvbuf,
768 const int* recvcounts, const int* recvdispls, MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request)
773 if(sendbuf != MPI_IN_PLACE){
774 CHECK_NULL(2, MPI_ERR_COUNT, sendcounts)
775 CHECK_NULL(3, MPI_ERR_ARG, senddispls)
776 CHECK_TYPE(4, sendtype)
778 CHECK_TYPE(8, recvtype)
779 CHECK_NULL(6, MPI_ERR_COUNT, recvcounts)
780 CHECK_NULL(7, MPI_ERR_ARG, recvdispls)
783 aid_t pid = simgrid::s4u::this_actor::get_pid();
784 int size = comm->size();
785 for (int i = 0; i < size; i++) {
786 if(sendbuf != MPI_IN_PLACE){
787 CHECK_BUFFER(1, sendbuf, sendcounts[i], sendtype)
788 CHECK_COUNT(2, sendcounts[i])
790 CHECK_BUFFER(5, recvbuf, recvcounts[i], recvtype)
791 CHECK_COUNT(6, recvcounts[i])
794 const SmpiBenchGuard suspend_bench;
797 auto trace_sendcounts = std::make_shared<std::vector<int>>();
798 auto trace_recvcounts = std::make_shared<std::vector<int>>();
799 trace_recvcounts->insert(trace_recvcounts->end(), &recvcounts[0], &recvcounts[size]);
801 int dt_size_recv = recvtype->size();
803 const int* real_sendcounts = sendcounts;
804 const int* real_senddispls = senddispls;
805 MPI_Datatype real_sendtype = sendtype;
807 std::vector<unsigned char> tmp_sendbuf;
808 std::vector<int> tmp_sendcounts;
809 std::vector<int> tmp_senddispls;
810 const void* real_sendbuf;
812 if (sendbuf == MPI_IN_PLACE) {
813 tmp_sendcounts.assign(recvcounts, recvcounts + size);
814 real_sendcounts = tmp_sendcounts.data();
815 tmp_senddispls.assign(recvdispls, recvdispls + size);
816 real_senddispls = tmp_senddispls.data();
817 real_sendtype = recvtype;
820 for (int i = 0; i < size; i++) { // copy data to avoid bad free
821 send_size += real_sendcounts[i] ;
822 recv_size += recvcounts[i];
823 if (((recvdispls[i] + recvcounts[i]) * dt_size_recv) > maxsize)
824 maxsize = (recvdispls[i] + recvcounts[i]) * dt_size_recv;
826 real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, maxsize, MPI_CHAR);
828 if(recvtype->size() * recvcounts[comm->rank()] != real_sendtype->size() * real_sendcounts[comm->rank()]){
829 XBT_WARN("MPI_(I)Alltoallv : receive size from me differs from sent size to me : %zu vs %zu", recvtype->size() * recvcounts[comm->rank()], real_sendtype->size() * real_sendcounts[comm->rank()]);
830 return MPI_ERR_TRUNCATE;
833 trace_sendcounts->insert(trace_sendcounts->end(), &real_sendcounts[0], &real_sendcounts[size]);
835 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Alltoallv" : "PMPI_Ialltoallv",
836 new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "alltoallv" : "ialltoallv", -1,
837 send_size, trace_sendcounts, recv_size, trace_recvcounts,
838 simgrid::smpi::Datatype::encode(real_sendtype),
839 simgrid::smpi::Datatype::encode(recvtype)));
842 if (request == MPI_REQUEST_IGNORED)
843 retval = simgrid::smpi::colls::alltoallv(real_sendbuf, real_sendcounts, real_senddispls, real_sendtype, recvbuf,
844 recvcounts, recvdispls, recvtype, comm);
846 retval = simgrid::smpi::colls::ialltoallv(real_sendbuf, real_sendcounts, real_senddispls, real_sendtype, recvbuf,
847 recvcounts, recvdispls, recvtype, comm, request);
849 TRACE_smpi_comm_out(pid);
853 int PMPI_Alltoallw(const void* sendbuf, const int* sendcounts, const int* senddispls, const MPI_Datatype* sendtypes, void* recvbuf,
854 const int* recvcounts, const int* recvdispls, const MPI_Datatype* recvtypes, MPI_Comm comm)
856 return PMPI_Ialltoallw(sendbuf, sendcounts, senddispls, sendtypes, recvbuf, recvcounts, recvdispls, recvtypes, comm, MPI_REQUEST_IGNORED);
859 int PMPI_Ialltoallw(const void* sendbuf, const int* sendcounts, const int* senddispls, const MPI_Datatype* sendtypes, void* recvbuf,
860 const int* recvcounts, const int* recvdispls, const MPI_Datatype* recvtypes, MPI_Comm comm, MPI_Request* request)
865 if(sendbuf != MPI_IN_PLACE){
866 CHECK_NULL(2, MPI_ERR_COUNT, sendcounts)
867 CHECK_NULL(3, MPI_ERR_ARG, senddispls)
868 CHECK_NULL(4, MPI_ERR_TYPE, sendtypes)
870 CHECK_NULL(6, MPI_ERR_COUNT, recvcounts)
871 CHECK_NULL(7, MPI_ERR_ARG, recvdispls)
872 CHECK_NULL(8, MPI_ERR_TYPE, recvtypes)
874 aid_t pid = simgrid::s4u::this_actor::get_pid();
875 int size = comm->size();
876 for (int i = 0; i < size; i++) {
877 if(sendbuf != MPI_IN_PLACE){
878 CHECK_COUNT(2, sendcounts[i])
879 CHECK_TYPE(4, sendtypes[i])
880 CHECK_BUFFER(1, sendbuf, sendcounts[i], sendtypes[i])
882 CHECK_COUNT(6, recvcounts[i])
883 CHECK_TYPE(8, recvtypes[i])
884 CHECK_BUFFER(5, recvbuf, recvcounts[i], recvtypes[i])
887 const SmpiBenchGuard suspend_bench;
891 auto trace_sendcounts = std::make_shared<std::vector<int>>();
892 auto trace_recvcounts = std::make_shared<std::vector<int>>();
893 trace_recvcounts->insert(trace_recvcounts->end(), &recvcounts[0], &recvcounts[size]);
895 const int* real_sendcounts = sendcounts;
896 const int* real_senddispls = senddispls;
897 const MPI_Datatype* real_sendtypes = sendtypes;
899 unsigned long maxsize = 0;
900 for (int i = 0; i < size; i++) { // copy data to avoid bad free
901 if (recvtypes[i] == MPI_DATATYPE_NULL)
903 recv_size += recvcounts[i] * recvtypes[i]->size();
904 if ((recvdispls[i] + (recvcounts[i] * recvtypes[i]->size())) > maxsize)
905 maxsize = recvdispls[i] + (recvcounts[i] * recvtypes[i]->size());
908 std::vector<unsigned char> tmp_sendbuf;
909 std::vector<int> tmp_sendcounts;
910 std::vector<int> tmp_senddispls;
911 std::vector<MPI_Datatype> tmp_sendtypes;
912 const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, maxsize, MPI_CHAR);
913 if (sendbuf == MPI_IN_PLACE) {
914 tmp_sendcounts.assign(recvcounts, recvcounts + size);
915 real_sendcounts = tmp_sendcounts.data();
916 tmp_senddispls.assign(recvdispls, recvdispls + size);
917 real_senddispls = tmp_senddispls.data();
918 tmp_sendtypes.assign(recvtypes, recvtypes + size);
919 real_sendtypes = tmp_sendtypes.data();
923 if(recvtypes[comm->rank()]->size() * recvcounts[comm->rank()] != real_sendtypes[comm->rank()]->size() * real_sendcounts[comm->rank()]){
924 XBT_WARN("MPI_(I)Alltoallw : receive size from me differs from sent size to me : %zu vs %zu", recvtypes[comm->rank()]->size() * recvcounts[comm->rank()], real_sendtypes[comm->rank()]->size() * real_sendcounts[comm->rank()]);
925 return MPI_ERR_TRUNCATE;
928 trace_sendcounts->insert(trace_sendcounts->end(), &real_sendcounts[0], &real_sendcounts[size]);
929 for (int i = 0; i < size; i++) { // copy data to avoid bad free
930 send_size += real_sendcounts[i] * real_sendtypes[i]->size();
933 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Alltoallw" : "PMPI_Ialltoallw",
934 new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "alltoallv" : "ialltoallv", -1,
935 send_size, trace_sendcounts, recv_size, trace_recvcounts,
936 simgrid::smpi::Datatype::encode(real_sendtypes[0]),
937 simgrid::smpi::Datatype::encode(recvtypes[0])));
940 if (request == MPI_REQUEST_IGNORED)
941 retval = simgrid::smpi::colls::alltoallw(real_sendbuf, real_sendcounts, real_senddispls, real_sendtypes, recvbuf,
942 recvcounts, recvdispls, recvtypes, comm);
944 retval = simgrid::smpi::colls::ialltoallw(real_sendbuf, real_sendcounts, real_senddispls, real_sendtypes, recvbuf,
945 recvcounts, recvdispls, recvtypes, comm, request);
947 TRACE_smpi_comm_out(pid);