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,
62 int root, MPI_Comm comm, MPI_Request* request)
67 CHECK_TYPE(3, datatype)
68 CHECK_BUFFER(1, buf, count, datatype)
72 const SmpiBenchGuard suspend_bench;
73 aid_t pid = simgrid::s4u::this_actor::get_pid();
74 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Bcast" : "PMPI_Ibcast",
75 new simgrid::instr::CollTIData(request == MPI_REQUEST_IGNORED ? "bcast" : "ibcast", root, -1.0,
76 datatype->is_replayable() ? count : count * datatype->size(), 0,
77 simgrid::smpi::Datatype::encode(datatype), ""));
78 if (comm->size() > 1) {
79 if (request == MPI_REQUEST_IGNORED)
80 simgrid::smpi::colls::bcast(buf, count, datatype, root, comm);
82 simgrid::smpi::colls::ibcast(buf, count, datatype, root, comm, request);
84 if (request != MPI_REQUEST_IGNORED)
85 *request = MPI_REQUEST_NULL;
88 TRACE_smpi_comm_out(pid);
92 int PMPI_Gather(const void *sendbuf, int sendcount, MPI_Datatype sendtype,void *recvbuf, int recvcount, MPI_Datatype recvtype,
93 int root, MPI_Comm comm){
94 return PMPI_Igather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm, MPI_REQUEST_IGNORED);
97 int PMPI_Igather(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
98 MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request* request)
103 int rank = comm->rank();
104 if(sendbuf != MPI_IN_PLACE){
105 CHECK_COUNT(2, sendcount)
106 CHECK_TYPE(3, sendtype)
107 CHECK_BUFFER(1,sendbuf, sendcount, sendtype)
110 CHECK_NOT_IN_PLACE_ROOT(4, recvbuf)
111 CHECK_TYPE(6, recvtype)
112 CHECK_COUNT(5, recvcount)
113 CHECK_BUFFER(4, recvbuf, recvcount, recvtype)
115 CHECK_NOT_IN_PLACE_ROOT(1, sendbuf)
120 const void* real_sendbuf = sendbuf;
121 int real_sendcount = sendcount;
122 MPI_Datatype real_sendtype = sendtype;
124 if (sendbuf == MPI_IN_PLACE) {
126 real_sendtype = recvtype;
127 } else if(recvtype->size() * recvcount != sendtype->size() * sendcount){
128 XBT_WARN("MPI_(I)Gather : received size at root differs from sent size : %zu vs %zu", recvtype->size() * recvcount , sendtype->size() * sendcount);
129 return MPI_ERR_TRUNCATE;
133 const SmpiBenchGuard suspend_bench;
135 aid_t pid = simgrid::s4u::this_actor::get_pid();
137 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Gather" : "PMPI_Igather",
138 new simgrid::instr::CollTIData(
139 request == MPI_REQUEST_IGNORED ? "gather" : "igather", root, -1.0,
140 real_sendtype->is_replayable() ? real_sendcount : real_sendcount * real_sendtype->size(),
141 (comm->rank() != root || recvtype->is_replayable()) ? recvcount : recvcount * recvtype->size(),
142 simgrid::smpi::Datatype::encode(real_sendtype), simgrid::smpi::Datatype::encode(recvtype)));
143 if (request == MPI_REQUEST_IGNORED)
144 simgrid::smpi::colls::gather(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcount, recvtype, root, comm);
146 simgrid::smpi::colls::igather(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcount, recvtype, root, comm,
149 TRACE_smpi_comm_out(pid);
153 int PMPI_Gatherv(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, const int *recvcounts, const int *displs,
154 MPI_Datatype recvtype, int root, MPI_Comm comm){
155 return PMPI_Igatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, root, comm, MPI_REQUEST_IGNORED);
158 int PMPI_Igatherv(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, const int* recvcounts, const int* displs,
159 MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request* request)
164 int rank = comm->rank();
165 if(sendbuf != MPI_IN_PLACE){
166 CHECK_TYPE(3, sendtype)
167 CHECK_COUNT(2, sendcount)
169 CHECK_BUFFER(1, sendbuf, sendcount, sendtype)
171 CHECK_NOT_IN_PLACE_ROOT(4, recvbuf)
172 CHECK_TYPE(6, recvtype)
173 CHECK_NULL(5, MPI_ERR_COUNT, recvcounts)
174 CHECK_NULL(6, MPI_ERR_ARG, displs)
176 CHECK_NOT_IN_PLACE_ROOT(1, sendbuf)
182 for (int i = 0; i < comm->size(); i++) {
183 CHECK_COUNT(5, recvcounts[i])
184 CHECK_BUFFER(4,recvbuf,recvcounts[i], recvtype)
188 const SmpiBenchGuard suspend_bench;
189 const void* real_sendbuf = sendbuf;
190 int real_sendcount = sendcount;
191 MPI_Datatype real_sendtype = sendtype;
192 if ((rank == root) && (sendbuf == MPI_IN_PLACE)) {
194 real_sendtype = recvtype;
197 aid_t pid = simgrid::s4u::this_actor::get_pid();
198 int dt_size_recv = recvtype->is_replayable() ? 1 : recvtype->size();
200 auto trace_recvcounts = std::make_shared<std::vector<int>>();
202 for (int i = 0; i < comm->size(); i++) // copy data to avoid bad free
203 trace_recvcounts->push_back(recvcounts[i] * dt_size_recv);
206 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Gatherv" : "PMPI_Igatherv",
207 new simgrid::instr::VarCollTIData(
208 request == MPI_REQUEST_IGNORED ? "gatherv" : "igatherv", root,
209 real_sendtype->is_replayable() ? real_sendcount : real_sendcount * real_sendtype->size(),
210 nullptr, dt_size_recv, trace_recvcounts, simgrid::smpi::Datatype::encode(real_sendtype),
211 simgrid::smpi::Datatype::encode(recvtype)));
212 if (request == MPI_REQUEST_IGNORED)
213 simgrid::smpi::colls::gatherv(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcounts, displs, recvtype,
216 simgrid::smpi::colls::igatherv(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcounts, displs, recvtype,
217 root, comm, request);
219 TRACE_smpi_comm_out(pid);
223 int PMPI_Allgather(const void *sendbuf, int sendcount, MPI_Datatype sendtype,
224 void *recvbuf, int recvcount, MPI_Datatype recvtype, MPI_Comm comm){
225 return PMPI_Iallgather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm, MPI_REQUEST_IGNORED);
228 int PMPI_Iallgather(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
229 MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request)
234 int rank = comm->rank();
235 CHECK_NOT_IN_PLACE(4, recvbuf)
236 if(sendbuf != MPI_IN_PLACE){
237 CHECK_COUNT(2, sendcount)
238 CHECK_TYPE(3, sendtype)
240 CHECK_TYPE(6, recvtype)
241 CHECK_COUNT(5, recvcount)
242 CHECK_BUFFER(1, sendbuf, sendcount, sendtype)
243 CHECK_BUFFER(4, recvbuf, recvcount, recvtype)
246 if (sendbuf == MPI_IN_PLACE) {
247 sendbuf = static_cast<char*>(recvbuf) + recvtype->get_extent() * recvcount * comm->rank();
248 sendcount = recvcount;
252 if(recvtype->size() * recvcount != sendtype->size() * sendcount){
253 XBT_WARN("MPI_(I)Allgather : received size from each process differs from sent size : %zu vs %zu", recvtype->size() * recvcount, sendtype->size() * sendcount);
254 return MPI_ERR_TRUNCATE;
257 const SmpiBenchGuard suspend_bench;
259 aid_t pid = simgrid::s4u::this_actor::get_pid();
261 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Allgather" : "PMPI_Iallggather",
262 new simgrid::instr::CollTIData(
263 request == MPI_REQUEST_IGNORED ? "allgather" : "iallgather", -1, -1.0,
264 sendtype->is_replayable() ? sendcount : sendcount * sendtype->size(),
265 recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(),
266 simgrid::smpi::Datatype::encode(sendtype), simgrid::smpi::Datatype::encode(recvtype)));
267 if (request == MPI_REQUEST_IGNORED)
268 simgrid::smpi::colls::allgather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
270 simgrid::smpi::colls::iallgather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm, request);
272 TRACE_smpi_comm_out(pid);
276 int PMPI_Allgatherv(const void *sendbuf, int sendcount, MPI_Datatype sendtype,
277 void *recvbuf, const int *recvcounts, const int *displs, MPI_Datatype recvtype, MPI_Comm comm){
278 return PMPI_Iallgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm, MPI_REQUEST_IGNORED);
281 int PMPI_Iallgatherv(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, const int* recvcounts, const int* displs,
282 MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request)
287 int rank = comm->rank();
288 if(sendbuf != MPI_IN_PLACE)
289 CHECK_TYPE(3, sendtype)
290 CHECK_TYPE(6, recvtype)
291 CHECK_NULL(5, MPI_ERR_COUNT, recvcounts)
292 CHECK_NULL(6, MPI_ERR_ARG, displs)
293 if(sendbuf != MPI_IN_PLACE){
294 CHECK_COUNT(2, sendcount)
295 CHECK_BUFFER(1, sendbuf, sendcount, sendtype)
298 CHECK_NOT_IN_PLACE(4, recvbuf)
299 for (int i = 0; i < comm->size(); i++) {
300 CHECK_COUNT(5, recvcounts[i])
301 CHECK_BUFFER(4, recvbuf, recvcounts[i], recvtype)
304 const SmpiBenchGuard suspend_bench;
305 if (sendbuf == MPI_IN_PLACE) {
306 sendbuf = static_cast<char*>(recvbuf) + recvtype->get_extent() * displs[comm->rank()];
307 sendcount = recvcounts[comm->rank()];
310 aid_t pid = simgrid::s4u::this_actor::get_pid();
311 int dt_size_recv = recvtype->is_replayable() ? 1 : recvtype->size();
313 auto trace_recvcounts = std::make_shared<std::vector<int>>();
314 for (int i = 0; i < comm->size(); i++) { // copy data to avoid bad free
315 trace_recvcounts->push_back(recvcounts[i] * dt_size_recv);
319 pid, request == MPI_REQUEST_IGNORED ? "PMPI_Allgatherv" : "PMPI_Iallgatherv",
320 new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "allgatherv" : "iallgatherv", -1,
321 sendtype->is_replayable() ? sendcount : sendcount * sendtype->size(), nullptr,
322 dt_size_recv, trace_recvcounts, simgrid::smpi::Datatype::encode(sendtype),
323 simgrid::smpi::Datatype::encode(recvtype)));
324 if (request == MPI_REQUEST_IGNORED)
325 simgrid::smpi::colls::allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm);
327 simgrid::smpi::colls::iallgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm,
330 TRACE_smpi_comm_out(pid);
334 int PMPI_Scatter(const void *sendbuf, int sendcount, MPI_Datatype sendtype,
335 void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm){
336 return PMPI_Iscatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm, MPI_REQUEST_IGNORED);
339 int PMPI_Iscatter(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
340 MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request* request)
345 int rank = comm->rank();
347 CHECK_NOT_IN_PLACE_ROOT(1, sendbuf)
348 CHECK_COUNT(2, sendcount)
349 CHECK_TYPE(3, sendtype)
350 CHECK_BUFFER(1, sendbuf, sendcount, sendtype)
352 CHECK_NOT_IN_PLACE_ROOT(4, recvbuf)
354 if(recvbuf != MPI_IN_PLACE){
355 CHECK_COUNT(5, recvcount)
356 CHECK_TYPE(6, recvtype)
357 CHECK_BUFFER(4, recvbuf, recvcount, recvtype)
362 if (recvbuf == MPI_IN_PLACE) {
364 recvcount = sendcount;
367 if((rank == root) && (recvtype->size() * recvcount != sendtype->size() * sendcount)){
368 XBT_WARN("MPI_(I)Scatter : sent size to each process differs from receive size");
369 return MPI_ERR_TRUNCATE;
372 const SmpiBenchGuard suspend_bench;
374 aid_t pid = simgrid::s4u::this_actor::get_pid();
376 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Scatter" : "PMPI_Iscatter",
377 new simgrid::instr::CollTIData(
378 request == MPI_REQUEST_IGNORED ? "scatter" : "iscatter", root, -1.0,
379 (rank != root || sendtype->is_replayable()) ? sendcount : sendcount * sendtype->size(),
380 recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(),
381 simgrid::smpi::Datatype::encode(sendtype), simgrid::smpi::Datatype::encode(recvtype)));
382 if (request == MPI_REQUEST_IGNORED)
383 simgrid::smpi::colls::scatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm);
385 simgrid::smpi::colls::iscatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm, request);
387 TRACE_smpi_comm_out(pid);
391 int PMPI_Scatterv(const void *sendbuf, const int *sendcounts, const int *displs,
392 MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm){
393 return PMPI_Iscatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm, MPI_REQUEST_IGNORED);
396 int PMPI_Iscatterv(const void* sendbuf, const int* sendcounts, const int* displs, MPI_Datatype sendtype, void* recvbuf, int recvcount,
397 MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request* request)
402 int rank = comm->rank();
403 if(recvbuf != MPI_IN_PLACE){
404 CHECK_NOT_IN_PLACE_ROOT(1, sendbuf)
405 CHECK_COUNT(5, recvcount)
406 CHECK_TYPE(7, recvtype)
407 CHECK_BUFFER(4, recvbuf, recvcount, recvtype)
412 CHECK_NULL(2, MPI_ERR_COUNT, sendcounts)
413 CHECK_NULL(3, MPI_ERR_ARG, displs)
414 CHECK_TYPE(4, sendtype)
415 for (int i = 0; i < comm->size(); i++){
416 CHECK_COUNT(2, sendcounts[i])
417 CHECK_BUFFER(1, sendbuf, sendcounts[i], sendtype)
419 if (recvbuf == MPI_IN_PLACE) {
421 recvcount = sendcounts[rank];
424 CHECK_NOT_IN_PLACE_ROOT(4, recvbuf)
427 const SmpiBenchGuard suspend_bench;
429 aid_t pid = simgrid::s4u::this_actor::get_pid();
430 int dt_size_send = sendtype->is_replayable() ? 1 : sendtype->size();
432 auto trace_sendcounts = std::make_shared<std::vector<int>>();
434 for (int i = 0; i < comm->size(); i++) { // copy data to avoid bad free
435 trace_sendcounts->push_back(sendcounts[i] * dt_size_send);
439 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Scatterv" : "PMPI_Iscatterv",
440 new simgrid::instr::VarCollTIData(
441 request == MPI_REQUEST_IGNORED ? "scatterv" : "iscatterv", root, dt_size_send,
442 trace_sendcounts, recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(),
443 nullptr, simgrid::smpi::Datatype::encode(sendtype),
444 simgrid::smpi::Datatype::encode(recvtype)));
445 if (request == MPI_REQUEST_IGNORED)
446 simgrid::smpi::colls::scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm);
448 simgrid::smpi::colls::iscatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm,
451 TRACE_smpi_comm_out(pid);
455 int PMPI_Reduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
457 return PMPI_Ireduce(sendbuf, recvbuf, count, datatype, op, root, comm, MPI_REQUEST_IGNORED);
460 int PMPI_Ireduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm, MPI_Request* request)
465 int rank = comm->rank();
466 CHECK_TYPE(4, datatype)
467 CHECK_COUNT(3, count)
468 CHECK_BUFFER(1, sendbuf, count, datatype)
470 CHECK_NOT_IN_PLACE(2, recvbuf)
471 CHECK_BUFFER(5, recvbuf, count, datatype)
473 CHECK_OP(5, op, datatype)
477 const SmpiBenchGuard suspend_bench;
478 aid_t pid = simgrid::s4u::this_actor::get_pid();
480 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Reduce" : "PMPI_Ireduce",
481 new simgrid::instr::CollTIData(request == MPI_REQUEST_IGNORED ? "reduce" : "ireduce", root, 0,
482 datatype->is_replayable() ? count : count * datatype->size(), 0,
483 simgrid::smpi::Datatype::encode(datatype), ""));
484 if (request == MPI_REQUEST_IGNORED)
485 simgrid::smpi::colls::reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
487 simgrid::smpi::colls::ireduce(sendbuf, recvbuf, count, datatype, op, root, comm, request);
489 TRACE_smpi_comm_out(pid);
493 int PMPI_Reduce_local(const void* inbuf, void* inoutbuf, int count, MPI_Datatype datatype, MPI_Op op)
497 CHECK_TYPE(4, datatype)
498 CHECK_COUNT(3, count)
499 CHECK_BUFFER(1, inbuf, count, datatype)
500 CHECK_BUFFER(2, inoutbuf, count, datatype)
501 CHECK_OP(5, op, datatype)
503 const SmpiBenchGuard suspend_bench;
504 op->apply(inbuf, inoutbuf, &count, datatype);
508 int PMPI_Allreduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
510 return PMPI_Iallreduce(sendbuf, recvbuf, count, datatype, op, comm, MPI_REQUEST_IGNORED);
513 int PMPI_Iallreduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request *request)
518 int rank = comm->rank();
519 CHECK_NOT_IN_PLACE(2, recvbuf)
520 CHECK_TYPE(4, datatype)
521 CHECK_OP(5, op, datatype)
522 CHECK_COUNT(3, count)
523 CHECK_BUFFER(1, sendbuf, count, datatype)
524 CHECK_BUFFER(2, recvbuf, count, datatype)
527 const SmpiBenchGuard suspend_bench;
528 std::vector<unsigned char> tmp_sendbuf;
529 const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, count, datatype);
531 aid_t pid = simgrid::s4u::this_actor::get_pid();
533 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Allreduce" : "PMPI_Iallreduce",
534 new simgrid::instr::CollTIData(request == MPI_REQUEST_IGNORED ? "allreduce" : "iallreduce", -1, 0,
535 datatype->is_replayable() ? count : count * datatype->size(), 0,
536 simgrid::smpi::Datatype::encode(datatype), ""));
538 if (request == MPI_REQUEST_IGNORED)
539 simgrid::smpi::colls::allreduce(real_sendbuf, recvbuf, count, datatype, op, comm);
541 simgrid::smpi::colls::iallreduce(real_sendbuf, recvbuf, count, datatype, op, comm, request);
543 TRACE_smpi_comm_out(pid);
547 int PMPI_Scan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
549 return PMPI_Iscan(sendbuf, recvbuf, count, datatype, op, comm, MPI_REQUEST_IGNORED);
552 int PMPI_Iscan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request* request)
557 CHECK_TYPE(4, datatype)
558 CHECK_COUNT(3, count)
559 CHECK_BUFFER(1,sendbuf,count, datatype)
560 CHECK_BUFFER(2,recvbuf,count, datatype)
562 CHECK_OP(5, op, datatype)
564 const SmpiBenchGuard suspend_bench;
565 aid_t pid = simgrid::s4u::this_actor::get_pid();
566 std::vector<unsigned char> tmp_sendbuf;
567 const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, count, datatype);
569 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Scan" : "PMPI_Iscan",
570 new simgrid::instr::Pt2PtTIData(request == MPI_REQUEST_IGNORED ? "scan" : "iscan", -1,
571 datatype->is_replayable() ? count : count * datatype->size(),
572 simgrid::smpi::Datatype::encode(datatype)));
575 if (request == MPI_REQUEST_IGNORED)
576 retval = simgrid::smpi::colls::scan(real_sendbuf, recvbuf, count, datatype, op, comm);
578 retval = simgrid::smpi::colls::iscan(real_sendbuf, recvbuf, count, datatype, op, comm, request);
580 TRACE_smpi_comm_out(pid);
584 int PMPI_Exscan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
586 return PMPI_Iexscan(sendbuf, recvbuf, count, datatype, op, comm, MPI_REQUEST_IGNORED);
589 int PMPI_Iexscan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request* request){
593 CHECK_TYPE(4, datatype)
594 CHECK_COUNT(3, count)
595 CHECK_BUFFER(1, sendbuf, count, datatype)
596 CHECK_BUFFER(2, recvbuf, count, datatype)
598 CHECK_OP(5, op, datatype)
600 const SmpiBenchGuard suspend_bench;
601 aid_t pid = simgrid::s4u::this_actor::get_pid();
602 std::vector<unsigned char> tmp_sendbuf;
603 const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, count, datatype);
605 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Exscan" : "PMPI_Iexscan",
606 new simgrid::instr::Pt2PtTIData(request == MPI_REQUEST_IGNORED ? "exscan" : "iexscan", -1,
607 datatype->is_replayable() ? count : count * datatype->size(),
608 simgrid::smpi::Datatype::encode(datatype)));
611 if (request == MPI_REQUEST_IGNORED)
612 retval = simgrid::smpi::colls::exscan(real_sendbuf, recvbuf, count, datatype, op, comm);
614 retval = simgrid::smpi::colls::iexscan(real_sendbuf, recvbuf, count, datatype, op, comm, request);
616 TRACE_smpi_comm_out(pid);
620 int PMPI_Reduce_scatter(const void *sendbuf, void *recvbuf, const int *recvcounts, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
622 return PMPI_Ireduce_scatter(sendbuf, recvbuf, recvcounts, datatype, op, comm, MPI_REQUEST_IGNORED);
625 int PMPI_Ireduce_scatter(const void *sendbuf, void *recvbuf, const int *recvcounts, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request *request)
630 int rank = comm->rank();
631 CHECK_NOT_IN_PLACE(2, recvbuf)
632 CHECK_TYPE(4, datatype)
633 CHECK_NULL(3, MPI_ERR_COUNT, recvcounts)
635 CHECK_OP(5, op, datatype)
636 for (int i = 0; i < comm->size(); i++) {
637 CHECK_COUNT(3, recvcounts[i])
638 CHECK_BUFFER(1, sendbuf, recvcounts[i], datatype)
639 CHECK_BUFFER(2, recvbuf, recvcounts[i], datatype)
642 const SmpiBenchGuard suspend_bench;
643 aid_t pid = simgrid::s4u::this_actor::get_pid();
644 auto trace_recvcounts = std::make_shared<std::vector<int>>();
645 int dt_send_size = datatype->is_replayable() ? 1 : datatype->size();
648 for (int i = 0; i < comm->size(); i++) { // copy data to avoid bad free
649 trace_recvcounts->push_back(recvcounts[i] * dt_send_size);
650 totalcount += recvcounts[i];
652 std::vector<unsigned char> tmp_sendbuf;
653 const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, totalcount, datatype);
655 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Reduce_scatter" : "PMPI_Ireduce_scatter",
656 new simgrid::instr::VarCollTIData(
657 request == MPI_REQUEST_IGNORED ? "reducescatter" : "ireducescatter", -1, dt_send_size, nullptr,
658 0, trace_recvcounts, simgrid::smpi::Datatype::encode(datatype), ""));
660 if (request == MPI_REQUEST_IGNORED)
661 simgrid::smpi::colls::reduce_scatter(real_sendbuf, recvbuf, recvcounts, datatype, op, comm);
663 simgrid::smpi::colls::ireduce_scatter(real_sendbuf, recvbuf, recvcounts, datatype, op, comm, request);
665 TRACE_smpi_comm_out(pid);
669 int PMPI_Reduce_scatter_block(const void *sendbuf, void *recvbuf, int recvcount,
670 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
672 return PMPI_Ireduce_scatter_block(sendbuf, recvbuf, recvcount, datatype, op, comm, MPI_REQUEST_IGNORED);
675 int PMPI_Ireduce_scatter_block(const void* sendbuf, void* recvbuf, int recvcount, MPI_Datatype datatype, MPI_Op op,
676 MPI_Comm comm, MPI_Request* request)
681 CHECK_TYPE(4, datatype)
682 CHECK_COUNT(3, recvcount)
683 CHECK_BUFFER(1, sendbuf, recvcount, datatype)
684 CHECK_BUFFER(2, recvbuf, recvcount, datatype)
686 CHECK_OP(5, op, datatype)
688 const SmpiBenchGuard suspend_bench;
689 int count = comm->size();
691 aid_t pid = simgrid::s4u::this_actor::get_pid();
692 int dt_send_size = datatype->is_replayable() ? 1 : datatype->size();
693 auto trace_recvcounts = std::make_shared<std::vector<int>>(recvcount * dt_send_size); // copy data to avoid bad free
694 std::vector<unsigned char> tmp_sendbuf;
695 const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, recvcount * count, datatype);
698 pid, request == MPI_REQUEST_IGNORED ? "PMPI_Reduce_scatter_block" : "PMPI_Ireduce_scatter_block",
699 new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "reducescatter" : "ireducescatter", -1, 0,
700 nullptr, 0, trace_recvcounts, simgrid::smpi::Datatype::encode(datatype), ""));
702 std::vector<int> recvcounts(count);
703 for (int i = 0; i < count; i++)
704 recvcounts[i] = recvcount;
705 if (request == MPI_REQUEST_IGNORED)
706 simgrid::smpi::colls::reduce_scatter(real_sendbuf, recvbuf, recvcounts.data(), datatype, op, comm);
708 simgrid::smpi::colls::ireduce_scatter(real_sendbuf, recvbuf, recvcounts.data(), datatype, op, comm, request);
710 TRACE_smpi_comm_out(pid);
714 int PMPI_Alltoall(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
715 MPI_Datatype recvtype, MPI_Comm comm){
716 return PMPI_Ialltoall(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm, MPI_REQUEST_IGNORED);
719 int PMPI_Ialltoall(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
720 MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request)
725 if(sendbuf != MPI_IN_PLACE){
726 CHECK_TYPE(3, sendtype)
727 CHECK_COUNT(2, sendcount)
728 CHECK_BUFFER(1, sendbuf, sendcount, sendtype)
730 CHECK_TYPE(6, recvtype)
731 CHECK_COUNT(5, recvcount)
732 CHECK_COUNT(5, recvcount)
733 CHECK_BUFFER(4, recvbuf, recvcount, recvtype)
736 aid_t pid = simgrid::s4u::this_actor::get_pid();
737 int real_sendcount = sendcount;
738 MPI_Datatype real_sendtype = sendtype;
740 std::vector<unsigned char> tmp_sendbuf;
741 const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, recvcount * comm->size(), recvtype);
743 if (sendbuf == MPI_IN_PLACE) {
744 real_sendcount = recvcount;
745 real_sendtype = recvtype;
748 if(recvtype->size() * recvcount != real_sendtype->size() * real_sendcount){
749 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);
750 return MPI_ERR_TRUNCATE;
753 const SmpiBenchGuard suspend_bench;
755 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Alltoall" : "PMPI_Ialltoall",
756 new simgrid::instr::CollTIData(
757 request == MPI_REQUEST_IGNORED ? "alltoall" : "ialltoall", -1, -1.0,
758 real_sendtype->is_replayable() ? real_sendcount : real_sendcount * real_sendtype->size(),
759 recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(),
760 simgrid::smpi::Datatype::encode(real_sendtype), simgrid::smpi::Datatype::encode(recvtype)));
762 if (request == MPI_REQUEST_IGNORED)
764 simgrid::smpi::colls::alltoall(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcount, recvtype, comm);
766 retval = simgrid::smpi::colls::ialltoall(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcount, recvtype,
769 TRACE_smpi_comm_out(pid);
773 int PMPI_Alltoallv(const void* sendbuf, const int* sendcounts, const int* senddispls, MPI_Datatype sendtype, void* recvbuf,
774 const int* recvcounts, const int* recvdispls, MPI_Datatype recvtype, MPI_Comm comm)
776 return PMPI_Ialltoallv(sendbuf, sendcounts, senddispls, sendtype, recvbuf, recvcounts, recvdispls, recvtype, comm, MPI_REQUEST_IGNORED);
779 int PMPI_Ialltoallv(const void* sendbuf, const int* sendcounts, const int* senddispls, MPI_Datatype sendtype, void* recvbuf,
780 const int* recvcounts, const int* recvdispls, MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request)
785 if(sendbuf != MPI_IN_PLACE){
786 CHECK_NULL(2, MPI_ERR_COUNT, sendcounts)
787 CHECK_NULL(3, MPI_ERR_ARG, senddispls)
788 CHECK_TYPE(4, sendtype)
790 CHECK_TYPE(8, recvtype)
791 CHECK_NULL(6, MPI_ERR_COUNT, recvcounts)
792 CHECK_NULL(7, MPI_ERR_ARG, recvdispls)
795 aid_t pid = simgrid::s4u::this_actor::get_pid();
796 int size = comm->size();
797 for (int i = 0; i < size; i++) {
798 if(sendbuf != MPI_IN_PLACE){
799 CHECK_BUFFER(1, sendbuf, sendcounts[i], sendtype)
800 CHECK_COUNT(2, sendcounts[i])
802 CHECK_BUFFER(5, recvbuf, recvcounts[i], recvtype)
803 CHECK_COUNT(6, recvcounts[i])
806 const SmpiBenchGuard suspend_bench;
809 auto trace_sendcounts = std::make_shared<std::vector<int>>();
810 auto trace_recvcounts = std::make_shared<std::vector<int>>();
811 int dt_size_recv = recvtype->size();
813 const int* real_sendcounts = sendcounts;
814 const int* real_senddispls = senddispls;
815 MPI_Datatype real_sendtype = sendtype;
817 for (int i = 0; i < size; i++) { // copy data to avoid bad free
818 recv_size += recvcounts[i] * dt_size_recv;
819 trace_recvcounts->push_back(recvcounts[i] * dt_size_recv);
820 if (((recvdispls[i] + recvcounts[i]) * dt_size_recv) > maxsize)
821 maxsize = (recvdispls[i] + recvcounts[i]) * dt_size_recv;
824 std::vector<unsigned char> tmp_sendbuf;
825 std::vector<int> tmp_sendcounts;
826 std::vector<int> tmp_senddispls;
827 const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, maxsize, MPI_CHAR);
828 if (sendbuf == MPI_IN_PLACE) {
829 tmp_sendcounts.assign(recvcounts, recvcounts + size);
830 real_sendcounts = tmp_sendcounts.data();
831 tmp_senddispls.assign(recvdispls, recvdispls + size);
832 real_senddispls = tmp_senddispls.data();
833 real_sendtype = recvtype;
836 if(recvtype->size() * recvcounts[comm->rank()] != real_sendtype->size() * real_sendcounts[comm->rank()]){
837 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()]);
838 return MPI_ERR_TRUNCATE;
841 int dt_size_send = real_sendtype->size();
843 for (int i = 0; i < size; i++) { // copy data to avoid bad free
844 send_size += real_sendcounts[i] * dt_size_send;
845 trace_sendcounts->push_back(real_sendcounts[i] * dt_size_send);
848 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Alltoallv" : "PMPI_Ialltoallv",
849 new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "alltoallv" : "ialltoallv", -1,
850 send_size, trace_sendcounts, recv_size, trace_recvcounts,
851 simgrid::smpi::Datatype::encode(real_sendtype),
852 simgrid::smpi::Datatype::encode(recvtype)));
855 if (request == MPI_REQUEST_IGNORED)
856 retval = simgrid::smpi::colls::alltoallv(real_sendbuf, real_sendcounts, real_senddispls, real_sendtype, recvbuf,
857 recvcounts, recvdispls, recvtype, comm);
859 retval = simgrid::smpi::colls::ialltoallv(real_sendbuf, real_sendcounts, real_senddispls, real_sendtype, recvbuf,
860 recvcounts, recvdispls, recvtype, comm, request);
862 TRACE_smpi_comm_out(pid);
866 int PMPI_Alltoallw(const void* sendbuf, const int* sendcounts, const int* senddispls, const MPI_Datatype* sendtypes, void* recvbuf,
867 const int* recvcounts, const int* recvdispls, const MPI_Datatype* recvtypes, MPI_Comm comm)
869 return PMPI_Ialltoallw(sendbuf, sendcounts, senddispls, sendtypes, recvbuf, recvcounts, recvdispls, recvtypes, comm, MPI_REQUEST_IGNORED);
872 int PMPI_Ialltoallw(const void* sendbuf, const int* sendcounts, const int* senddispls, const MPI_Datatype* sendtypes, void* recvbuf,
873 const int* recvcounts, const int* recvdispls, const MPI_Datatype* recvtypes, MPI_Comm comm, MPI_Request* request)
878 if(sendbuf != MPI_IN_PLACE){
879 CHECK_NULL(2, MPI_ERR_COUNT, sendcounts)
880 CHECK_NULL(3, MPI_ERR_ARG, senddispls)
881 CHECK_NULL(4, MPI_ERR_TYPE, sendtypes)
883 CHECK_NULL(6, MPI_ERR_COUNT, recvcounts)
884 CHECK_NULL(7, MPI_ERR_ARG, recvdispls)
885 CHECK_NULL(8, MPI_ERR_TYPE, recvtypes)
887 aid_t pid = simgrid::s4u::this_actor::get_pid();
888 int size = comm->size();
889 for (int i = 0; i < size; i++) {
890 if(sendbuf != MPI_IN_PLACE){
891 CHECK_COUNT(2, sendcounts[i])
892 CHECK_TYPE(4, sendtypes[i])
893 CHECK_BUFFER(1, sendbuf, sendcounts[i], sendtypes[i])
895 CHECK_COUNT(6, recvcounts[i])
896 CHECK_TYPE(8, recvtypes[i])
897 CHECK_BUFFER(5, recvbuf, recvcounts[i], recvtypes[i])
900 const SmpiBenchGuard suspend_bench;
904 auto trace_sendcounts = std::make_shared<std::vector<int>>();
905 auto trace_recvcounts = std::make_shared<std::vector<int>>();
907 const int* real_sendcounts = sendcounts;
908 const int* real_senddispls = senddispls;
909 const MPI_Datatype* real_sendtypes = sendtypes;
910 unsigned long maxsize = 0;
911 for (int i = 0; i < size; i++) { // copy data to avoid bad free
912 if (recvtypes[i] == MPI_DATATYPE_NULL)
914 recv_size += recvcounts[i] * recvtypes[i]->size();
915 trace_recvcounts->push_back(recvcounts[i] * recvtypes[i]->size());
916 if ((recvdispls[i] + (recvcounts[i] * recvtypes[i]->size())) > maxsize)
917 maxsize = recvdispls[i] + (recvcounts[i] * recvtypes[i]->size());
920 std::vector<unsigned char> tmp_sendbuf;
921 std::vector<int> tmp_sendcounts;
922 std::vector<int> tmp_senddispls;
923 std::vector<MPI_Datatype> tmp_sendtypes;
924 const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, maxsize, MPI_CHAR);
925 if (sendbuf == MPI_IN_PLACE) {
926 tmp_sendcounts.assign(recvcounts, recvcounts + size);
927 real_sendcounts = tmp_sendcounts.data();
928 tmp_senddispls.assign(recvdispls, recvdispls + size);
929 real_senddispls = tmp_senddispls.data();
930 tmp_sendtypes.assign(recvtypes, recvtypes + size);
931 real_sendtypes = tmp_sendtypes.data();
935 if(recvtypes[comm->rank()]->size() * recvcounts[comm->rank()] != real_sendtypes[comm->rank()]->size() * real_sendcounts[comm->rank()]){
936 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()]);
937 return MPI_ERR_TRUNCATE;
940 for (int i = 0; i < size; i++) { // copy data to avoid bad free
941 send_size += real_sendcounts[i] * real_sendtypes[i]->size();
942 trace_sendcounts->push_back(real_sendcounts[i] * real_sendtypes[i]->size());
945 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Alltoallw" : "PMPI_Ialltoallw",
946 new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "alltoallv" : "ialltoallv", -1,
947 send_size, trace_sendcounts, recv_size, trace_recvcounts,
948 simgrid::smpi::Datatype::encode(real_sendtypes[0]),
949 simgrid::smpi::Datatype::encode(recvtypes[0])));
952 if (request == MPI_REQUEST_IGNORED)
953 retval = simgrid::smpi::colls::alltoallw(real_sendbuf, real_sendcounts, real_senddispls, real_sendtypes, recvbuf,
954 recvcounts, recvdispls, recvtypes, comm);
956 retval = simgrid::smpi::colls::ialltoallw(real_sendbuf, real_sendcounts, real_senddispls, real_sendtypes, recvbuf,
957 recvcounts, recvdispls, recvtypes, comm, request);
959 TRACE_smpi_comm_out(pid);