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)
42 int rank = simgrid::s4u::this_actor::get_pid();
43 TRACE_smpi_comm_in(rank, 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(rank);
57 int PMPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
59 return PMPI_Ibcast(buf, count, datatype, root, comm, MPI_REQUEST_IGNORED);
62 int PMPI_Ibcast(void *buf, int count, MPI_Datatype datatype,
63 int root, MPI_Comm comm, MPI_Request* request)
66 CHECK_BUFFER(1, buf, count)
68 CHECK_TYPE(3, datatype)
73 int 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(), -1,
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);
93 int PMPI_Gather(const void *sendbuf, int sendcount, MPI_Datatype sendtype,void *recvbuf, int recvcount, MPI_Datatype recvtype,
94 int root, MPI_Comm comm){
95 return PMPI_Igather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm, MPI_REQUEST_IGNORED);
98 int PMPI_Igather(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
99 MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request* request)
102 int rank = comm->rank();
103 if(sendbuf != MPI_IN_PLACE){
104 CHECK_BUFFER(1,sendbuf, sendcount)
105 CHECK_COUNT(2, sendcount)
106 CHECK_TYPE(3, 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)
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;
134 int 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_sendtype->is_replayable() ? real_sendcount : real_sendcount * real_sendtype->size(),
140 (comm->rank() != root || recvtype->is_replayable()) ? recvcount : recvcount * recvtype->size(),
141 simgrid::smpi::Datatype::encode(real_sendtype), simgrid::smpi::Datatype::encode(recvtype)));
142 if (request == MPI_REQUEST_IGNORED)
143 simgrid::smpi::colls::gather(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcount, recvtype, root, comm);
145 simgrid::smpi::colls::igather(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcount, recvtype, root, comm,
148 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)
162 int rank = comm->rank();
163 CHECK_BUFFER(1, sendbuf, sendcount)
164 if(sendbuf != MPI_IN_PLACE){
165 CHECK_TYPE(3, sendtype)
166 CHECK_COUNT(2, sendcount)
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])
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 int pid = simgrid::s4u::this_actor::get_pid();
196 int dt_size_recv = recvtype->is_replayable() ? 1 : recvtype->size();
198 auto trace_recvcounts = std::make_shared<std::vector<int>>();
200 for (int i = 0; i < comm->size(); i++) // copy data to avoid bad free
201 trace_recvcounts->push_back(recvcounts[i] * dt_size_recv);
204 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Gatherv" : "PMPI_Igatherv",
205 new simgrid::instr::VarCollTIData(
206 request == MPI_REQUEST_IGNORED ? "gatherv" : "igatherv", root,
207 real_sendtype->is_replayable() ? real_sendcount : real_sendcount * real_sendtype->size(),
208 nullptr, dt_size_recv, trace_recvcounts, simgrid::smpi::Datatype::encode(real_sendtype),
209 simgrid::smpi::Datatype::encode(recvtype)));
210 if (request == MPI_REQUEST_IGNORED)
211 simgrid::smpi::colls::gatherv(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcounts, displs, recvtype,
214 simgrid::smpi::colls::igatherv(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcounts, displs, recvtype,
215 root, comm, request);
217 TRACE_smpi_comm_out(pid);
222 int PMPI_Allgather(const void *sendbuf, int sendcount, MPI_Datatype sendtype,
223 void *recvbuf, int recvcount, MPI_Datatype recvtype, MPI_Comm comm){
224 return PMPI_Iallgather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm, MPI_REQUEST_IGNORED);
227 int PMPI_Iallgather(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
228 MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request)
231 int rank = comm->rank();
232 CHECK_BUFFER(1, sendbuf, sendcount)
233 CHECK_BUFFER(4, recvbuf, recvcount)
234 CHECK_NOT_IN_PLACE(4, recvbuf)
235 if(sendbuf != MPI_IN_PLACE){
236 CHECK_COUNT(2, sendcount)
237 CHECK_TYPE(3, sendtype)
239 CHECK_TYPE(6, recvtype)
240 CHECK_COUNT(5, recvcount)
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;
256 int 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 sendtype->is_replayable() ? sendcount : sendcount * sendtype->size(),
262 recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(),
263 simgrid::smpi::Datatype::encode(sendtype), simgrid::smpi::Datatype::encode(recvtype)));
264 if (request == MPI_REQUEST_IGNORED)
265 simgrid::smpi::colls::allgather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
267 simgrid::smpi::colls::iallgather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm, request);
269 TRACE_smpi_comm_out(pid);
274 int PMPI_Allgatherv(const void *sendbuf, int sendcount, MPI_Datatype sendtype,
275 void *recvbuf, const int *recvcounts, const int *displs, MPI_Datatype recvtype, MPI_Comm comm){
276 return PMPI_Iallgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm, MPI_REQUEST_IGNORED);
279 int PMPI_Iallgatherv(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, const int* recvcounts, const int* displs,
280 MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request)
283 int rank = comm->rank();
284 CHECK_BUFFER(1, sendbuf, sendcount)
285 if(sendbuf != MPI_IN_PLACE)
286 CHECK_TYPE(3, sendtype)
287 CHECK_TYPE(6, recvtype)
288 CHECK_NULL(5, MPI_ERR_COUNT, recvcounts)
289 CHECK_NULL(6, MPI_ERR_ARG, displs)
290 if(sendbuf != MPI_IN_PLACE)
291 CHECK_COUNT(2, sendcount)
293 CHECK_NOT_IN_PLACE(4, recvbuf)
294 for (int i = 0; i < comm->size(); i++) {
295 CHECK_COUNT(5, recvcounts[i])
296 CHECK_BUFFER(4, recvbuf, recvcounts[i])
300 if (sendbuf == MPI_IN_PLACE) {
301 sendbuf = static_cast<char*>(recvbuf) + recvtype->get_extent() * displs[comm->rank()];
302 sendcount = recvcounts[comm->rank()];
305 int pid = simgrid::s4u::this_actor::get_pid();
306 int dt_size_recv = recvtype->is_replayable() ? 1 : recvtype->size();
308 auto trace_recvcounts = std::make_shared<std::vector<int>>();
309 for (int i = 0; i < comm->size(); i++) { // copy data to avoid bad free
310 trace_recvcounts->push_back(recvcounts[i] * dt_size_recv);
314 pid, request == MPI_REQUEST_IGNORED ? "PMPI_Allgatherv" : "PMPI_Iallgatherv",
315 new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "allgatherv" : "iallgatherv", -1,
316 sendtype->is_replayable() ? sendcount : sendcount * sendtype->size(), nullptr,
317 dt_size_recv, trace_recvcounts, simgrid::smpi::Datatype::encode(sendtype),
318 simgrid::smpi::Datatype::encode(recvtype)));
319 if (request == MPI_REQUEST_IGNORED)
320 simgrid::smpi::colls::allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm);
322 simgrid::smpi::colls::iallgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm,
325 TRACE_smpi_comm_out(pid);
330 int PMPI_Scatter(const void *sendbuf, int sendcount, MPI_Datatype sendtype,
331 void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm){
332 return PMPI_Iscatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm, MPI_REQUEST_IGNORED);
335 int PMPI_Iscatter(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
336 MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request* request)
339 int rank = comm->rank();
341 CHECK_NOT_IN_PLACE_ROOT(1, sendbuf)
342 CHECK_BUFFER(1, sendbuf, sendcount)
343 CHECK_COUNT(2, sendcount)
344 CHECK_TYPE(3, sendtype)
346 CHECK_NOT_IN_PLACE_ROOT(4, recvbuf)
348 if(recvbuf != MPI_IN_PLACE){
349 CHECK_BUFFER(4, recvbuf, recvcount)
350 CHECK_COUNT(5, recvcount)
351 CHECK_TYPE(6, recvtype)
356 if (recvbuf == MPI_IN_PLACE) {
358 recvcount = sendcount;
361 if((rank == root) && (recvtype->size() * recvcount != sendtype->size() * sendcount)){
362 XBT_WARN("MPI_(I)Scatter : sent size to each process differs from receive size");
363 return MPI_ERR_TRUNCATE;
368 int pid = simgrid::s4u::this_actor::get_pid();
370 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Scatter" : "PMPI_Iscatter",
371 new simgrid::instr::CollTIData(
372 request == MPI_REQUEST_IGNORED ? "scatter" : "iscatter", root, -1.0,
373 (rank != root || sendtype->is_replayable()) ? sendcount : sendcount * sendtype->size(),
374 recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(),
375 simgrid::smpi::Datatype::encode(sendtype), simgrid::smpi::Datatype::encode(recvtype)));
376 if (request == MPI_REQUEST_IGNORED)
377 simgrid::smpi::colls::scatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm);
379 simgrid::smpi::colls::iscatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm, request);
381 TRACE_smpi_comm_out(pid);
386 int PMPI_Scatterv(const void *sendbuf, const int *sendcounts, const int *displs,
387 MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm){
388 return PMPI_Iscatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm, MPI_REQUEST_IGNORED);
391 int PMPI_Iscatterv(const void* sendbuf, const int* sendcounts, const int* displs, MPI_Datatype sendtype, void* recvbuf, int recvcount,
392 MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request* request)
395 int rank = comm->rank();
396 if(recvbuf != MPI_IN_PLACE){
397 CHECK_NOT_IN_PLACE_ROOT(1, sendbuf)
398 CHECK_BUFFER(4, recvbuf, recvcount)
399 CHECK_COUNT(5, recvcount)
400 CHECK_TYPE(7, recvtype)
405 CHECK_NULL(2, MPI_ERR_COUNT, sendcounts)
406 CHECK_NULL(3, MPI_ERR_ARG, displs)
407 CHECK_TYPE(4, sendtype)
408 for (int i = 0; i < comm->size(); i++){
409 CHECK_BUFFER(1, sendbuf, sendcounts[i])
410 CHECK_COUNT(2, sendcounts[i])
412 if (recvbuf == MPI_IN_PLACE) {
414 recvcount = sendcounts[rank];
417 CHECK_NOT_IN_PLACE_ROOT(4, recvbuf)
422 int pid = simgrid::s4u::this_actor::get_pid();
423 int dt_size_send = sendtype->is_replayable() ? 1 : sendtype->size();
425 auto trace_sendcounts = std::make_shared<std::vector<int>>();
427 for (int i = 0; i < comm->size(); i++) { // copy data to avoid bad free
428 trace_sendcounts->push_back(sendcounts[i] * dt_size_send);
432 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Scatterv" : "PMPI_Iscatterv",
433 new simgrid::instr::VarCollTIData(
434 request == MPI_REQUEST_IGNORED ? "scatterv" : "iscatterv", root, dt_size_send,
435 trace_sendcounts, recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(),
436 nullptr, simgrid::smpi::Datatype::encode(sendtype),
437 simgrid::smpi::Datatype::encode(recvtype)));
438 if (request == MPI_REQUEST_IGNORED)
439 simgrid::smpi::colls::scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm);
441 simgrid::smpi::colls::iscatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm,
444 TRACE_smpi_comm_out(pid);
449 int PMPI_Reduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
451 return PMPI_Ireduce(sendbuf, recvbuf, count, datatype, op, root, comm, MPI_REQUEST_IGNORED);
454 int PMPI_Ireduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm, MPI_Request* request)
457 int rank = comm->rank();
458 CHECK_BUFFER(1, sendbuf, count)
460 CHECK_NOT_IN_PLACE(2, recvbuf)
461 CHECK_BUFFER(5, recvbuf, count)
463 CHECK_TYPE(4, datatype)
464 CHECK_COUNT(3, count)
465 CHECK_OP(5, op, datatype)
470 int pid = simgrid::s4u::this_actor::get_pid();
472 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Reduce" : "PMPI_Ireduce",
473 new simgrid::instr::CollTIData(request == MPI_REQUEST_IGNORED ? "reduce" : "ireduce", root, 0,
474 datatype->is_replayable() ? count : count * datatype->size(), -1,
475 simgrid::smpi::Datatype::encode(datatype), ""));
476 if (request == MPI_REQUEST_IGNORED)
477 simgrid::smpi::colls::reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
479 simgrid::smpi::colls::ireduce(sendbuf, recvbuf, count, datatype, op, root, comm, request);
481 TRACE_smpi_comm_out(pid);
486 int PMPI_Reduce_local(const void* inbuf, void* inoutbuf, int count, MPI_Datatype datatype, MPI_Op op)
488 CHECK_BUFFER(1, inbuf, count)
489 CHECK_BUFFER(2, inoutbuf, count)
490 CHECK_TYPE(4, datatype)
491 CHECK_COUNT(3, count)
492 CHECK_OP(5, op, datatype)
495 op->apply(inbuf, inoutbuf, &count, datatype);
500 int PMPI_Allreduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
502 return PMPI_Iallreduce(sendbuf, recvbuf, count, datatype, op, comm, MPI_REQUEST_IGNORED);
505 int PMPI_Iallreduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request *request)
508 int rank = comm->rank();
509 CHECK_BUFFER(1, sendbuf, count)
510 CHECK_BUFFER(2, recvbuf, count)
511 CHECK_NOT_IN_PLACE(2, recvbuf)
512 CHECK_TYPE(4, datatype)
513 CHECK_COUNT(3, count)
515 CHECK_OP(5, op, datatype)
518 std::vector<unsigned char> tmp_sendbuf;
519 const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, count, datatype);
521 int pid = simgrid::s4u::this_actor::get_pid();
523 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Allreduce" : "PMPI_Iallreduce",
524 new simgrid::instr::CollTIData(request == MPI_REQUEST_IGNORED ? "allreduce" : "iallreduce", -1, 0,
525 datatype->is_replayable() ? count : count * datatype->size(), -1,
526 simgrid::smpi::Datatype::encode(datatype), ""));
528 if (request == MPI_REQUEST_IGNORED)
529 simgrid::smpi::colls::allreduce(real_sendbuf, recvbuf, count, datatype, op, comm);
531 simgrid::smpi::colls::iallreduce(real_sendbuf, recvbuf, count, datatype, op, comm, request);
533 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)
546 CHECK_BUFFER(1,sendbuf,count)
547 CHECK_BUFFER(2,recvbuf,count)
548 CHECK_TYPE(4, datatype)
549 CHECK_COUNT(3, count)
551 CHECK_OP(5, op, datatype)
554 int pid = simgrid::s4u::this_actor::get_pid();
555 std::vector<unsigned char> tmp_sendbuf;
556 const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, count, datatype);
558 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Scan" : "PMPI_Iscan",
559 new simgrid::instr::Pt2PtTIData(request == MPI_REQUEST_IGNORED ? "scan" : "iscan", -1,
560 datatype->is_replayable() ? count : count * datatype->size(),
561 simgrid::smpi::Datatype::encode(datatype)));
564 if (request == MPI_REQUEST_IGNORED)
565 retval = simgrid::smpi::colls::scan(real_sendbuf, recvbuf, count, datatype, op, comm);
567 retval = simgrid::smpi::colls::iscan(real_sendbuf, recvbuf, count, datatype, op, comm, request);
569 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){
581 CHECK_BUFFER(1, sendbuf, count)
582 CHECK_BUFFER(2, recvbuf, count)
583 CHECK_TYPE(4, datatype)
584 CHECK_COUNT(3, count)
586 CHECK_OP(5, op, datatype)
589 int pid = simgrid::s4u::this_actor::get_pid();
590 std::vector<unsigned char> tmp_sendbuf;
591 const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, count, datatype);
593 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Exscan" : "PMPI_Iexscan",
594 new simgrid::instr::Pt2PtTIData(request == MPI_REQUEST_IGNORED ? "exscan" : "iexscan", -1,
595 datatype->is_replayable() ? count : count * datatype->size(),
596 simgrid::smpi::Datatype::encode(datatype)));
599 if (request == MPI_REQUEST_IGNORED)
600 retval = simgrid::smpi::colls::exscan(real_sendbuf, recvbuf, count, datatype, op, comm);
602 retval = simgrid::smpi::colls::iexscan(real_sendbuf, recvbuf, count, datatype, op, comm, request);
604 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)
617 int rank = comm->rank();
618 CHECK_NOT_IN_PLACE(2, recvbuf)
619 CHECK_TYPE(4, datatype)
620 CHECK_NULL(3, MPI_ERR_COUNT, recvcounts)
622 CHECK_OP(5, op, datatype)
623 for (int i = 0; i < comm->size(); i++) {
624 CHECK_COUNT(3, recvcounts[i])
625 CHECK_BUFFER(1, sendbuf, recvcounts[i])
626 CHECK_BUFFER(2, recvbuf, recvcounts[i])
630 int pid = simgrid::s4u::this_actor::get_pid();
631 auto trace_recvcounts = std::make_shared<std::vector<int>>();
632 int dt_send_size = datatype->is_replayable() ? 1 : datatype->size();
635 for (int i = 0; i < comm->size(); i++) { // copy data to avoid bad free
636 trace_recvcounts->push_back(recvcounts[i] * dt_send_size);
637 totalcount += recvcounts[i];
639 std::vector<unsigned char> tmp_sendbuf;
640 const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, totalcount, datatype);
642 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Reduce_scatter" : "PMPI_Ireduce_scatter",
643 new simgrid::instr::VarCollTIData(
644 request == MPI_REQUEST_IGNORED ? "reducescatter" : "ireducescatter", -1, dt_send_size, nullptr,
645 -1, trace_recvcounts, simgrid::smpi::Datatype::encode(datatype), ""));
647 if (request == MPI_REQUEST_IGNORED)
648 simgrid::smpi::colls::reduce_scatter(real_sendbuf, recvbuf, recvcounts, datatype, op, comm);
650 simgrid::smpi::colls::ireduce_scatter(real_sendbuf, recvbuf, recvcounts, datatype, op, comm, request);
652 TRACE_smpi_comm_out(pid);
657 int PMPI_Reduce_scatter_block(const void *sendbuf, void *recvbuf, int recvcount,
658 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
660 return PMPI_Ireduce_scatter_block(sendbuf, recvbuf, recvcount, datatype, op, comm, MPI_REQUEST_IGNORED);
663 int PMPI_Ireduce_scatter_block(const void* sendbuf, void* recvbuf, int recvcount, MPI_Datatype datatype, MPI_Op op,
664 MPI_Comm comm, MPI_Request* request)
667 CHECK_BUFFER(1, sendbuf, recvcount)
668 CHECK_BUFFER(2, recvbuf, recvcount)
669 CHECK_TYPE(4, datatype)
670 CHECK_COUNT(3, recvcount)
672 CHECK_OP(5, op, datatype)
675 int count = comm->size();
677 int pid = simgrid::s4u::this_actor::get_pid();
678 int dt_send_size = datatype->is_replayable() ? 1 : datatype->size();
679 auto trace_recvcounts = std::make_shared<std::vector<int>>(recvcount * dt_send_size); // copy data to avoid bad free
680 std::vector<unsigned char> tmp_sendbuf;
681 const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, recvcount * count, datatype);
684 pid, request == MPI_REQUEST_IGNORED ? "PMPI_Reduce_scatter_block" : "PMPI_Ireduce_scatter_block",
685 new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "reducescatter" : "ireducescatter", -1, 0,
686 nullptr, -1, trace_recvcounts, simgrid::smpi::Datatype::encode(datatype), ""));
688 std::vector<int> recvcounts(count);
689 for (int i = 0; i < count; i++)
690 recvcounts[i] = recvcount;
691 if (request == MPI_REQUEST_IGNORED)
692 simgrid::smpi::colls::reduce_scatter(real_sendbuf, recvbuf, recvcounts.data(), datatype, op, comm);
694 simgrid::smpi::colls::ireduce_scatter(real_sendbuf, recvbuf, recvcounts.data(), datatype, op, comm, request);
696 TRACE_smpi_comm_out(pid);
701 int PMPI_Alltoall(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
702 MPI_Datatype recvtype, MPI_Comm comm){
703 return PMPI_Ialltoall(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm, MPI_REQUEST_IGNORED);
706 int PMPI_Ialltoall(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
707 MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request)
710 CHECK_BUFFER(1, sendbuf, sendcount)
711 CHECK_BUFFER(4, recvbuf, recvcount)
712 if(sendbuf != MPI_IN_PLACE)
713 CHECK_TYPE(3, sendtype)
714 CHECK_TYPE(6, recvtype)
715 CHECK_COUNT(5, recvcount)
716 if(sendbuf != MPI_IN_PLACE)
717 CHECK_COUNT(2, sendcount)
718 CHECK_COUNT(5, recvcount)
721 int pid = simgrid::s4u::this_actor::get_pid();
722 int real_sendcount = sendcount;
723 MPI_Datatype real_sendtype = sendtype;
725 std::vector<unsigned char> tmp_sendbuf;
726 const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, recvcount * comm->size(), recvtype);
728 if (sendbuf == MPI_IN_PLACE) {
729 real_sendcount = recvcount;
730 real_sendtype = recvtype;
733 if(recvtype->size() * recvcount != real_sendtype->size() * real_sendcount){
734 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);
735 return MPI_ERR_TRUNCATE;
740 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Alltoall" : "PMPI_Ialltoall",
741 new simgrid::instr::CollTIData(
742 request == MPI_REQUEST_IGNORED ? "alltoall" : "ialltoall", -1, -1.0,
743 real_sendtype->is_replayable() ? real_sendcount : real_sendcount * real_sendtype->size(),
744 recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(),
745 simgrid::smpi::Datatype::encode(real_sendtype), simgrid::smpi::Datatype::encode(recvtype)));
747 if (request == MPI_REQUEST_IGNORED)
749 simgrid::smpi::colls::alltoall(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcount, recvtype, comm);
751 retval = simgrid::smpi::colls::ialltoall(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcount, recvtype,
754 TRACE_smpi_comm_out(pid);
759 int PMPI_Alltoallv(const void* sendbuf, const int* sendcounts, const int* senddispls, MPI_Datatype sendtype, void* recvbuf,
760 const int* recvcounts, const int* recvdispls, MPI_Datatype recvtype, MPI_Comm comm)
762 return PMPI_Ialltoallv(sendbuf, sendcounts, senddispls, sendtype, recvbuf, recvcounts, recvdispls, recvtype, comm, MPI_REQUEST_IGNORED);
765 int PMPI_Ialltoallv(const void* sendbuf, const int* sendcounts, const int* senddispls, MPI_Datatype sendtype, void* recvbuf,
766 const int* recvcounts, const int* recvdispls, MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request)
769 if(sendbuf != MPI_IN_PLACE){
770 CHECK_NULL(2, MPI_ERR_COUNT, sendcounts)
771 CHECK_NULL(3, MPI_ERR_ARG, senddispls)
772 CHECK_TYPE(4, sendtype)
774 CHECK_TYPE(8, recvtype)
775 CHECK_NULL(6, MPI_ERR_COUNT, recvcounts)
776 CHECK_NULL(7, MPI_ERR_ARG, recvdispls)
779 int pid = simgrid::s4u::this_actor::get_pid();
780 int size = comm->size();
781 for (int i = 0; i < size; i++) {
782 if(sendbuf != MPI_IN_PLACE){
783 CHECK_BUFFER(1, sendbuf, sendcounts[i])
784 CHECK_COUNT(2, sendcounts[i])
786 CHECK_BUFFER(5, recvbuf, recvcounts[i])
787 CHECK_COUNT(6, recvcounts[i])
793 auto trace_sendcounts = std::make_shared<std::vector<int>>();
794 auto trace_recvcounts = std::make_shared<std::vector<int>>();
795 int dt_size_recv = recvtype->size();
797 const int* real_sendcounts = sendcounts;
798 const int* real_senddispls = senddispls;
799 MPI_Datatype real_sendtype = sendtype;
801 for (int i = 0; i < size; i++) { // copy data to avoid bad free
802 recv_size += recvcounts[i] * dt_size_recv;
803 trace_recvcounts->push_back(recvcounts[i] * dt_size_recv);
804 if (((recvdispls[i] + recvcounts[i]) * dt_size_recv) > maxsize)
805 maxsize = (recvdispls[i] + recvcounts[i]) * dt_size_recv;
808 std::vector<unsigned char> tmp_sendbuf;
809 std::vector<int> tmp_sendcounts;
810 std::vector<int> tmp_senddispls;
811 const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, maxsize, MPI_CHAR);
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 if(recvtype->size() * recvcounts[comm->rank()] != real_sendtype->size() * real_sendcounts[comm->rank()]){
821 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()]);
823 return MPI_ERR_TRUNCATE;
826 int dt_size_send = real_sendtype->size();
828 for (int i = 0; i < size; i++) { // copy data to avoid bad free
829 send_size += real_sendcounts[i] * dt_size_send;
830 trace_sendcounts->push_back(real_sendcounts[i] * dt_size_send);
833 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Alltoallv" : "PMPI_Ialltoallv",
834 new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "alltoallv" : "ialltoallv", -1,
835 send_size, trace_sendcounts, recv_size, trace_recvcounts,
836 simgrid::smpi::Datatype::encode(real_sendtype),
837 simgrid::smpi::Datatype::encode(recvtype)));
840 if (request == MPI_REQUEST_IGNORED)
841 retval = simgrid::smpi::colls::alltoallv(real_sendbuf, real_sendcounts, real_senddispls, real_sendtype, recvbuf,
842 recvcounts, recvdispls, recvtype, comm);
844 retval = simgrid::smpi::colls::ialltoallv(real_sendbuf, real_sendcounts, real_senddispls, real_sendtype, recvbuf,
845 recvcounts, recvdispls, recvtype, comm, request);
847 TRACE_smpi_comm_out(pid);
852 int PMPI_Alltoallw(const void* sendbuf, const int* sendcounts, const int* senddispls, const MPI_Datatype* sendtypes, void* recvbuf,
853 const int* recvcounts, const int* recvdispls, const MPI_Datatype* recvtypes, MPI_Comm comm)
855 return PMPI_Ialltoallw(sendbuf, sendcounts, senddispls, sendtypes, recvbuf, recvcounts, recvdispls, recvtypes, comm, MPI_REQUEST_IGNORED);
858 int PMPI_Ialltoallw(const void* sendbuf, const int* sendcounts, const int* senddispls, const MPI_Datatype* sendtypes, void* recvbuf,
859 const int* recvcounts, const int* recvdispls, const MPI_Datatype* recvtypes, MPI_Comm comm, MPI_Request* request)
862 if(sendbuf != MPI_IN_PLACE){
863 CHECK_NULL(2, MPI_ERR_COUNT, sendcounts)
864 CHECK_NULL(3, MPI_ERR_ARG, senddispls)
865 CHECK_NULL(4, MPI_ERR_TYPE, sendtypes)
867 CHECK_NULL(6, MPI_ERR_COUNT, recvcounts)
868 CHECK_NULL(7, MPI_ERR_ARG, recvdispls)
869 CHECK_NULL(8, MPI_ERR_TYPE, recvtypes)
871 int pid = simgrid::s4u::this_actor::get_pid();
872 int size = comm->size();
873 for (int i = 0; i < size; i++) {
874 if(sendbuf != MPI_IN_PLACE){
875 CHECK_BUFFER(1, sendbuf, sendcounts[i])
876 CHECK_COUNT(2, sendcounts[i])
877 CHECK_TYPE(4, sendtypes[i])
879 CHECK_BUFFER(5, recvbuf, recvcounts[i])
880 CHECK_COUNT(6, recvcounts[i])
881 CHECK_TYPE(8, recvtypes[i])
888 auto trace_sendcounts = std::make_shared<std::vector<int>>();
889 auto trace_recvcounts = std::make_shared<std::vector<int>>();
891 const int* real_sendcounts = sendcounts;
892 const int* real_senddispls = senddispls;
893 const MPI_Datatype* real_sendtypes = sendtypes;
894 unsigned long maxsize = 0;
895 for (int i = 0; i < size; i++) { // copy data to avoid bad free
896 if (recvtypes[i] == MPI_DATATYPE_NULL)
898 recv_size += recvcounts[i] * recvtypes[i]->size();
899 trace_recvcounts->push_back(recvcounts[i] * recvtypes[i]->size());
900 if ((recvdispls[i] + (recvcounts[i] * recvtypes[i]->size())) > maxsize)
901 maxsize = recvdispls[i] + (recvcounts[i] * recvtypes[i]->size());
904 std::vector<unsigned char> tmp_sendbuf;
905 std::vector<int> tmp_sendcounts;
906 std::vector<int> tmp_senddispls;
907 std::vector<MPI_Datatype> tmp_sendtypes;
908 const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, maxsize, MPI_CHAR);
909 if (sendbuf == MPI_IN_PLACE) {
910 tmp_sendcounts.assign(recvcounts, recvcounts + size);
911 real_sendcounts = tmp_sendcounts.data();
912 tmp_senddispls.assign(recvdispls, recvdispls + size);
913 real_senddispls = tmp_senddispls.data();
914 tmp_sendtypes.assign(recvtypes, recvtypes + size);
915 real_sendtypes = tmp_sendtypes.data();
919 if(recvtypes[comm->rank()]->size() * recvcounts[comm->rank()] != real_sendtypes[comm->rank()]->size() * real_sendcounts[comm->rank()]){
920 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()]);
922 return MPI_ERR_TRUNCATE;
925 for (int i = 0; i < size; i++) { // copy data to avoid bad free
926 send_size += real_sendcounts[i] * real_sendtypes[i]->size();
927 trace_sendcounts->push_back(real_sendcounts[i] * real_sendtypes[i]->size());
930 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Alltoallw" : "PMPI_Ialltoallw",
931 new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "alltoallv" : "ialltoallv", -1,
932 send_size, trace_sendcounts, recv_size, trace_recvcounts,
933 simgrid::smpi::Datatype::encode(real_sendtypes[0]),
934 simgrid::smpi::Datatype::encode(recvtypes[0])));
937 if (request == MPI_REQUEST_IGNORED)
938 retval = simgrid::smpi::colls::alltoallw(real_sendbuf, real_sendcounts, real_senddispls, real_sendtypes, recvbuf,
939 recvcounts, recvdispls, recvtypes, comm);
941 retval = simgrid::smpi::colls::ialltoallw(real_sendbuf, real_sendcounts, real_senddispls, real_sendtypes, recvbuf,
942 recvcounts, recvdispls, recvtypes, comm, request);
944 TRACE_smpi_comm_out(pid);