1 /* Copyright (c) 2007-2022. 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)
40 CHECK_COLLECTIVE(comm, request == MPI_REQUEST_IGNORED ? "PMPI_Barrier" : "PMPI_Ibarrier")
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)
70 CHECK_COLLECTIVE(comm, request == MPI_REQUEST_IGNORED ? "PMPI_Bcast" : "PMPI_Ibcast")
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,
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)
102 int rank = comm->rank();
103 if(sendbuf != MPI_IN_PLACE){
104 CHECK_COUNT(2, sendcount)
105 CHECK_TYPE(3, sendtype)
106 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)
119 CHECK_COLLECTIVE(comm, request == MPI_REQUEST_IGNORED ? "PMPI_Gather" : "PMPI_Igather")
121 const void* real_sendbuf = sendbuf;
122 int real_sendcount = sendcount;
123 MPI_Datatype real_sendtype = sendtype;
125 if (sendbuf == MPI_IN_PLACE) {
127 real_sendtype = recvtype;
128 } else if(recvtype->size() * recvcount != sendtype->size() * sendcount){
129 XBT_WARN("MPI_(I)Gather : received size at root differs from sent size : %zu vs %zu", recvtype->size() * recvcount , sendtype->size() * sendcount);
130 return MPI_ERR_TRUNCATE;
134 const SmpiBenchGuard suspend_bench;
136 aid_t pid = simgrid::s4u::this_actor::get_pid();
138 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Gather" : "PMPI_Igather",
139 new simgrid::instr::CollTIData(
140 request == MPI_REQUEST_IGNORED ? "gather" : "igather", root, -1.0,
141 real_sendcount, recvcount,
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)
163 int rank = comm->rank();
164 if(sendbuf != MPI_IN_PLACE){
165 CHECK_TYPE(3, sendtype)
166 CHECK_COUNT(2, sendcount)
168 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)
180 CHECK_COLLECTIVE(comm, request == MPI_REQUEST_IGNORED ? "PMPI_Gatherv" : "PMPI_Igatherv")
183 for (int i = 0; i < comm->size(); i++) {
184 CHECK_COUNT(5, recvcounts[i])
185 CHECK_BUFFER(4,recvbuf,recvcounts[i], recvtype)
189 const SmpiBenchGuard suspend_bench;
190 const void* real_sendbuf = sendbuf;
191 int real_sendcount = sendcount;
192 MPI_Datatype real_sendtype = sendtype;
193 if ((rank == root) && (sendbuf == MPI_IN_PLACE)) {
195 real_sendtype = recvtype;
198 aid_t pid = simgrid::s4u::this_actor::get_pid();
200 auto trace_recvcounts = std::make_shared<std::vector<int>>();
202 trace_recvcounts->insert(trace_recvcounts->end(), &recvcounts[0], &recvcounts[comm->size()]);
203 else //this is not significant outside of root, put 0 as we don't know if recvcounts is initialized
204 trace_recvcounts->insert(trace_recvcounts->end(), comm->size(), 0);
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,
210 nullptr, -1, 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)
245 CHECK_COLLECTIVE(comm, request == MPI_REQUEST_IGNORED ? "PMPI_Allgather" : "PMPI_Iallggather")
247 if (sendbuf == MPI_IN_PLACE) {
248 sendbuf = static_cast<char*>(recvbuf) + recvtype->get_extent() * recvcount * comm->rank();
249 sendcount = recvcount;
253 if(recvtype->size() * recvcount != sendtype->size() * sendcount){
254 XBT_WARN("MPI_(I)Allgather : received size from each process differs from sent size : %zu vs %zu", recvtype->size() * recvcount, sendtype->size() * sendcount);
255 return MPI_ERR_TRUNCATE;
258 const SmpiBenchGuard suspend_bench;
260 aid_t pid = simgrid::s4u::this_actor::get_pid();
262 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Allgather" : "PMPI_Iallggather",
263 new simgrid::instr::CollTIData(
264 request == MPI_REQUEST_IGNORED ? "allgather" : "iallgather", -1, -1.0,
265 sendcount, recvcount,
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)
303 CHECK_COLLECTIVE(comm, MPI_REQUEST_IGNORED ? "PMPI_Allgatherv" : "PMPI_Iallgatherv")
305 const SmpiBenchGuard suspend_bench;
306 if (sendbuf == MPI_IN_PLACE) {
307 sendbuf = static_cast<char*>(recvbuf) + recvtype->get_extent() * displs[comm->rank()];
308 sendcount = recvcounts[comm->rank()];
311 aid_t pid = simgrid::s4u::this_actor::get_pid();
313 auto trace_recvcounts = std::make_shared<std::vector<int>>();
314 trace_recvcounts->insert(trace_recvcounts->end(), &recvcounts[0], &recvcounts[comm->size()]);
317 pid, request == MPI_REQUEST_IGNORED ? "PMPI_Allgatherv" : "PMPI_Iallgatherv",
318 new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "allgatherv" : "iallgatherv", -1,
320 -1, trace_recvcounts, simgrid::smpi::Datatype::encode(sendtype),
321 simgrid::smpi::Datatype::encode(recvtype)));
322 if (request == MPI_REQUEST_IGNORED)
323 simgrid::smpi::colls::allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm);
325 simgrid::smpi::colls::iallgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm,
328 TRACE_smpi_comm_out(pid);
332 int PMPI_Scatter(const void *sendbuf, int sendcount, MPI_Datatype sendtype,
333 void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm){
334 return PMPI_Iscatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm, MPI_REQUEST_IGNORED);
337 int PMPI_Iscatter(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
338 MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request* request)
342 int rank = comm->rank();
345 CHECK_NOT_IN_PLACE_ROOT(1, sendbuf)
346 CHECK_COUNT(2, sendcount)
347 CHECK_TYPE(3, sendtype)
348 CHECK_BUFFER(1, sendbuf, sendcount, sendtype)
350 CHECK_NOT_IN_PLACE_ROOT(4, recvbuf)
352 if(recvbuf != MPI_IN_PLACE){
353 CHECK_COUNT(5, recvcount)
354 CHECK_TYPE(6, recvtype)
355 CHECK_BUFFER(4, recvbuf, recvcount, recvtype)
359 CHECK_COLLECTIVE(comm, request == MPI_REQUEST_IGNORED ? "PMPI_Scatter" : "PMPI_Iscatter")
361 if (recvbuf == MPI_IN_PLACE) {
363 recvcount = sendcount;
366 if((rank == root) && (recvtype->size() * recvcount != sendtype->size() * sendcount)){
367 XBT_WARN("MPI_(I)Scatter : sent size to each process differs from receive size");
368 return MPI_ERR_TRUNCATE;
371 const SmpiBenchGuard suspend_bench;
373 aid_t pid = simgrid::s4u::this_actor::get_pid();
375 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Scatter" : "PMPI_Iscatter",
376 new simgrid::instr::CollTIData(
377 request == MPI_REQUEST_IGNORED ? "scatter" : "iscatter", root, -1.0,
378 sendcount, recvcount,
379 simgrid::smpi::Datatype::encode(sendtype), simgrid::smpi::Datatype::encode(recvtype)));
380 if (request == MPI_REQUEST_IGNORED)
381 simgrid::smpi::colls::scatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm);
383 simgrid::smpi::colls::iscatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm, request);
385 TRACE_smpi_comm_out(pid);
389 int PMPI_Scatterv(const void *sendbuf, const int *sendcounts, const int *displs,
390 MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm){
391 return PMPI_Iscatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm, MPI_REQUEST_IGNORED);
394 int PMPI_Iscatterv(const void* sendbuf, const int* sendcounts, const int* displs, MPI_Datatype sendtype, void* recvbuf, int recvcount,
395 MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request* request)
399 int rank = comm->rank();
400 if(recvbuf != MPI_IN_PLACE){
401 CHECK_COUNT(5, recvcount)
402 CHECK_TYPE(7, recvtype)
403 CHECK_BUFFER(4, recvbuf, recvcount, recvtype)
409 CHECK_NOT_IN_PLACE_ROOT(1, sendbuf)
410 CHECK_NULL(2, MPI_ERR_COUNT, sendcounts)
411 CHECK_NULL(3, MPI_ERR_ARG, displs)
412 CHECK_TYPE(4, sendtype)
413 for (int i = 0; i < comm->size(); i++){
414 CHECK_COUNT(2, sendcounts[i])
415 CHECK_BUFFER(1, sendbuf, sendcounts[i], sendtype)
417 if (recvbuf == MPI_IN_PLACE) {
419 recvcount = sendcounts[rank];
422 CHECK_NOT_IN_PLACE_ROOT(4, recvbuf)
424 CHECK_COLLECTIVE(comm, request == MPI_REQUEST_IGNORED ? "PMPI_Scatterv" : "PMPI_Iscatterv")
426 const SmpiBenchGuard suspend_bench;
428 aid_t pid = simgrid::s4u::this_actor::get_pid();
430 auto trace_sendcounts = std::make_shared<std::vector<int>>();
432 trace_sendcounts->insert(trace_sendcounts->end(), &sendcounts[0], &sendcounts[comm->size()]);
433 else //this is not significant outside of root, put 0 as we don't know if sendcounts is initialized
434 trace_sendcounts->insert(trace_sendcounts->end(), comm->size(), 0);
437 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Scatterv" : "PMPI_Iscatterv",
438 new simgrid::instr::VarCollTIData(
439 request == MPI_REQUEST_IGNORED ? "scatterv" : "iscatterv", root, -1,
440 trace_sendcounts, recvcount,
441 nullptr, simgrid::smpi::Datatype::encode(sendtype),
442 simgrid::smpi::Datatype::encode(recvtype)));
443 if (request == MPI_REQUEST_IGNORED)
444 simgrid::smpi::colls::scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm);
446 simgrid::smpi::colls::iscatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm,
449 TRACE_smpi_comm_out(pid);
453 int PMPI_Reduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
455 return PMPI_Ireduce(sendbuf, recvbuf, count, datatype, op, root, comm, MPI_REQUEST_IGNORED);
458 int PMPI_Ireduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm, MPI_Request* request)
462 int rank = comm->rank();
463 CHECK_TYPE(4, datatype)
464 CHECK_COUNT(3, count)
465 CHECK_BUFFER(1, sendbuf, count, datatype)
468 CHECK_NOT_IN_PLACE(2, recvbuf)
469 CHECK_BUFFER(5, recvbuf, count, datatype)
471 CHECK_OP(5, op, datatype)
474 CHECK_COLLECTIVE(comm, request == MPI_REQUEST_IGNORED ? "PMPI_Reduce" : "PMPI_Ireduce")
476 const SmpiBenchGuard suspend_bench;
477 aid_t pid = simgrid::s4u::this_actor::get_pid();
479 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Reduce" : "PMPI_Ireduce",
480 new simgrid::instr::CollTIData(request == MPI_REQUEST_IGNORED ? "reduce" : "ireduce", root, 0,
482 simgrid::smpi::Datatype::encode(datatype), ""));
483 if (request == MPI_REQUEST_IGNORED)
484 simgrid::smpi::colls::reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
486 simgrid::smpi::colls::ireduce(sendbuf, recvbuf, count, datatype, op, root, comm, request);
488 TRACE_smpi_comm_out(pid);
492 int PMPI_Reduce_local(const void* inbuf, void* inoutbuf, int count, MPI_Datatype datatype, MPI_Op op)
496 CHECK_TYPE(4, datatype)
497 CHECK_COUNT(3, count)
498 CHECK_BUFFER(1, inbuf, count, datatype)
499 CHECK_BUFFER(2, inoutbuf, count, datatype)
500 CHECK_OP(5, op, datatype)
502 const SmpiBenchGuard suspend_bench;
503 op->apply(inbuf, inoutbuf, &count, datatype);
507 int PMPI_Allreduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
509 return PMPI_Iallreduce(sendbuf, recvbuf, count, datatype, op, comm, MPI_REQUEST_IGNORED);
512 int PMPI_Iallreduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request *request)
517 int rank = comm->rank();
518 CHECK_NOT_IN_PLACE(2, recvbuf)
519 CHECK_TYPE(4, datatype)
520 CHECK_OP(5, op, datatype)
521 CHECK_COUNT(3, count)
522 CHECK_BUFFER(1, sendbuf, count, datatype)
523 CHECK_BUFFER(2, recvbuf, count, datatype)
525 CHECK_COLLECTIVE(comm, request == MPI_REQUEST_IGNORED ? "PMPI_Allreduce" : "PMPI_Iallreduce")
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,
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)
563 CHECK_COLLECTIVE(comm, request == MPI_REQUEST_IGNORED ? "PMPI_Scan" : "PMPI_Iscan")
565 const SmpiBenchGuard suspend_bench;
566 aid_t pid = simgrid::s4u::this_actor::get_pid();
567 std::vector<unsigned char> tmp_sendbuf;
568 const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, count, datatype);
570 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Scan" : "PMPI_Iscan",
571 new simgrid::instr::CollTIData(request == MPI_REQUEST_IGNORED ? "scan" : "iscan", -1, 0.0,
572 count, 0, 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)
599 CHECK_COLLECTIVE(comm, request == MPI_REQUEST_IGNORED ? "PMPI_Exscan" : "PMPI_Iexscan")
601 const SmpiBenchGuard suspend_bench;
602 aid_t pid = simgrid::s4u::this_actor::get_pid();
603 std::vector<unsigned char> tmp_sendbuf;
604 const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, count, datatype);
606 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Exscan" : "PMPI_Iexscan",
607 new simgrid::instr::CollTIData(request == MPI_REQUEST_IGNORED ? "exscan" : "iexscan", -1, 0.0,
608 count, 0, 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)
641 CHECK_COLLECTIVE(comm, request == MPI_REQUEST_IGNORED ? "PMPI_Reduce_scatter" : "PMPI_Ireduce_scatter")
643 const SmpiBenchGuard suspend_bench;
644 aid_t pid = simgrid::s4u::this_actor::get_pid();
645 auto trace_recvcounts = std::make_shared<std::vector<int>>();
646 trace_recvcounts->insert(trace_recvcounts->end(), &recvcounts[0], &recvcounts[comm->size()]);
650 for (int i = 0; i < comm->size(); i++) { // copy data to avoid bad free
651 totalcount += recvcounts[i];
653 std::vector<unsigned char> tmp_sendbuf;
654 const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, totalcount, datatype);
656 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Reduce_scatter" : "PMPI_Ireduce_scatter",
657 new simgrid::instr::VarCollTIData(
658 request == MPI_REQUEST_IGNORED ? "reducescatter" : "ireducescatter", -1, -1, nullptr,
659 -1 , trace_recvcounts, std::to_string(0), simgrid::smpi::Datatype::encode(datatype)));
661 if (request == MPI_REQUEST_IGNORED)
662 simgrid::smpi::colls::reduce_scatter(real_sendbuf, recvbuf, recvcounts, datatype, op, comm);
664 simgrid::smpi::colls::ireduce_scatter(real_sendbuf, recvbuf, recvcounts, datatype, op, comm, request);
666 TRACE_smpi_comm_out(pid);
670 int PMPI_Reduce_scatter_block(const void *sendbuf, void *recvbuf, int recvcount,
671 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
673 return PMPI_Ireduce_scatter_block(sendbuf, recvbuf, recvcount, datatype, op, comm, MPI_REQUEST_IGNORED);
676 int PMPI_Ireduce_scatter_block(const void* sendbuf, void* recvbuf, int recvcount, MPI_Datatype datatype, MPI_Op op,
677 MPI_Comm comm, MPI_Request* request)
682 CHECK_TYPE(4, datatype)
683 CHECK_COUNT(3, recvcount)
684 CHECK_BUFFER(1, sendbuf, recvcount, datatype)
685 CHECK_BUFFER(2, recvbuf, recvcount, datatype)
687 CHECK_OP(5, op, datatype)
688 CHECK_COLLECTIVE(comm, request == MPI_REQUEST_IGNORED ? "PMPI_Reduce_scatter_block" : "PMPI_Ireduce_scatter_block")
690 const SmpiBenchGuard suspend_bench;
691 int count = comm->size();
693 aid_t pid = simgrid::s4u::this_actor::get_pid();
694 auto trace_recvcounts = std::make_shared<std::vector<int>>(recvcount);
696 std::vector<unsigned char> tmp_sendbuf;
697 const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, recvcount * count, datatype);
700 pid, request == MPI_REQUEST_IGNORED ? "PMPI_Reduce_scatter_block" : "PMPI_Ireduce_scatter_block",
701 new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "reducescatter" : "ireducescatter", -1, -1,
702 nullptr, -1, trace_recvcounts, simgrid::smpi::Datatype::encode(datatype), ""));
704 std::vector<int> recvcounts(count);
705 for (int i = 0; i < count; i++)
706 recvcounts[i] = recvcount;
707 if (request == MPI_REQUEST_IGNORED)
708 simgrid::smpi::colls::reduce_scatter(real_sendbuf, recvbuf, recvcounts.data(), datatype, op, comm);
710 simgrid::smpi::colls::ireduce_scatter(real_sendbuf, recvbuf, recvcounts.data(), datatype, op, comm, request);
712 TRACE_smpi_comm_out(pid);
716 int PMPI_Alltoall(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
717 MPI_Datatype recvtype, MPI_Comm comm){
718 return PMPI_Ialltoall(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm, MPI_REQUEST_IGNORED);
721 int PMPI_Ialltoall(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
722 MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request)
727 if(sendbuf != MPI_IN_PLACE){
728 CHECK_TYPE(3, sendtype)
729 CHECK_COUNT(2, sendcount)
730 CHECK_BUFFER(1, sendbuf, sendcount, sendtype)
732 CHECK_TYPE(6, recvtype)
733 CHECK_COUNT(5, recvcount)
734 CHECK_COUNT(5, recvcount)
735 CHECK_BUFFER(4, recvbuf, recvcount, recvtype)
737 CHECK_COLLECTIVE(comm, request == MPI_REQUEST_IGNORED ? "PMPI_Alltoall" : "PMPI_Ialltoall")
739 aid_t pid = simgrid::s4u::this_actor::get_pid();
740 int real_sendcount = sendcount;
741 MPI_Datatype real_sendtype = sendtype;
743 std::vector<unsigned char> tmp_sendbuf;
744 const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, recvcount * comm->size(), recvtype);
746 if (sendbuf == MPI_IN_PLACE) {
747 real_sendcount = recvcount;
748 real_sendtype = recvtype;
751 if(recvtype->size() * recvcount != real_sendtype->size() * real_sendcount){
752 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);
753 return MPI_ERR_TRUNCATE;
756 const SmpiBenchGuard suspend_bench;
758 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Alltoall" : "PMPI_Ialltoall",
759 new simgrid::instr::CollTIData(
760 request == MPI_REQUEST_IGNORED ? "alltoall" : "ialltoall", -1, -1.0,
761 real_sendcount, recvcount,
762 simgrid::smpi::Datatype::encode(real_sendtype), simgrid::smpi::Datatype::encode(recvtype)));
764 if (request == MPI_REQUEST_IGNORED)
766 simgrid::smpi::colls::alltoall(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcount, recvtype, comm);
768 retval = simgrid::smpi::colls::ialltoall(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcount, recvtype,
771 TRACE_smpi_comm_out(pid);
775 int PMPI_Alltoallv(const void* sendbuf, const int* sendcounts, const int* senddispls, MPI_Datatype sendtype, void* recvbuf,
776 const int* recvcounts, const int* recvdispls, MPI_Datatype recvtype, MPI_Comm comm)
778 return PMPI_Ialltoallv(sendbuf, sendcounts, senddispls, sendtype, recvbuf, recvcounts, recvdispls, recvtype, comm, MPI_REQUEST_IGNORED);
781 int PMPI_Ialltoallv(const void* sendbuf, const int* sendcounts, const int* senddispls, MPI_Datatype sendtype, void* recvbuf,
782 const int* recvcounts, const int* recvdispls, MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request)
787 if(sendbuf != MPI_IN_PLACE){
788 CHECK_NULL(2, MPI_ERR_COUNT, sendcounts)
789 CHECK_NULL(3, MPI_ERR_ARG, senddispls)
790 CHECK_TYPE(4, sendtype)
792 CHECK_TYPE(8, recvtype)
793 CHECK_NULL(6, MPI_ERR_COUNT, recvcounts)
794 CHECK_NULL(7, MPI_ERR_ARG, recvdispls)
796 CHECK_COLLECTIVE(comm, request == MPI_REQUEST_IGNORED ? "PMPI_Alltoallv" : "PMPI_Ialltoallv")
798 aid_t pid = simgrid::s4u::this_actor::get_pid();
799 int size = comm->size();
800 for (int i = 0; i < size; i++) {
801 if(sendbuf != MPI_IN_PLACE){
802 CHECK_BUFFER(1, sendbuf, sendcounts[i], sendtype)
803 CHECK_COUNT(2, sendcounts[i])
805 CHECK_BUFFER(5, recvbuf, recvcounts[i], recvtype)
806 CHECK_COUNT(6, recvcounts[i])
809 const SmpiBenchGuard suspend_bench;
812 auto trace_sendcounts = std::make_shared<std::vector<int>>();
813 auto trace_recvcounts = std::make_shared<std::vector<int>>();
814 trace_recvcounts->insert(trace_recvcounts->end(), &recvcounts[0], &recvcounts[size]);
816 int dt_size_recv = recvtype->size();
818 const int* real_sendcounts = sendcounts;
819 const int* real_senddispls = senddispls;
820 MPI_Datatype real_sendtype = sendtype;
822 std::vector<unsigned char> tmp_sendbuf;
823 std::vector<int> tmp_sendcounts;
824 std::vector<int> tmp_senddispls;
825 const void* real_sendbuf;
827 if (sendbuf == MPI_IN_PLACE) {
828 tmp_sendcounts.assign(recvcounts, recvcounts + size);
829 real_sendcounts = tmp_sendcounts.data();
830 tmp_senddispls.assign(recvdispls, recvdispls + size);
831 real_senddispls = tmp_senddispls.data();
832 real_sendtype = recvtype;
835 for (int i = 0; i < size; i++) { // copy data to avoid bad free
836 send_size += real_sendcounts[i] ;
837 recv_size += recvcounts[i];
838 if (((recvdispls[i] + recvcounts[i]) * dt_size_recv) > maxsize)
839 maxsize = (recvdispls[i] + recvcounts[i]) * dt_size_recv;
841 real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, maxsize, MPI_CHAR);
843 if(recvtype->size() * recvcounts[comm->rank()] != real_sendtype->size() * real_sendcounts[comm->rank()]){
844 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()]);
845 return MPI_ERR_TRUNCATE;
848 trace_sendcounts->insert(trace_sendcounts->end(), &real_sendcounts[0], &real_sendcounts[size]);
850 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Alltoallv" : "PMPI_Ialltoallv",
851 new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "alltoallv" : "ialltoallv", -1,
852 send_size, trace_sendcounts, recv_size, trace_recvcounts,
853 simgrid::smpi::Datatype::encode(real_sendtype),
854 simgrid::smpi::Datatype::encode(recvtype)));
857 if (request == MPI_REQUEST_IGNORED)
858 retval = simgrid::smpi::colls::alltoallv(real_sendbuf, real_sendcounts, real_senddispls, real_sendtype, recvbuf,
859 recvcounts, recvdispls, recvtype, comm);
861 retval = simgrid::smpi::colls::ialltoallv(real_sendbuf, real_sendcounts, real_senddispls, real_sendtype, recvbuf,
862 recvcounts, recvdispls, recvtype, comm, request);
864 TRACE_smpi_comm_out(pid);
868 int PMPI_Alltoallw(const void* sendbuf, const int* sendcounts, const int* senddispls, const MPI_Datatype* sendtypes, void* recvbuf,
869 const int* recvcounts, const int* recvdispls, const MPI_Datatype* recvtypes, MPI_Comm comm)
871 return PMPI_Ialltoallw(sendbuf, sendcounts, senddispls, sendtypes, recvbuf, recvcounts, recvdispls, recvtypes, comm, MPI_REQUEST_IGNORED);
874 int PMPI_Ialltoallw(const void* sendbuf, const int* sendcounts, const int* senddispls, const MPI_Datatype* sendtypes, void* recvbuf,
875 const int* recvcounts, const int* recvdispls, const MPI_Datatype* recvtypes, MPI_Comm comm, MPI_Request* request)
880 if(sendbuf != MPI_IN_PLACE){
881 CHECK_NULL(2, MPI_ERR_COUNT, sendcounts)
882 CHECK_NULL(3, MPI_ERR_ARG, senddispls)
883 CHECK_NULL(4, MPI_ERR_TYPE, sendtypes)
885 CHECK_NULL(6, MPI_ERR_COUNT, recvcounts)
886 CHECK_NULL(7, MPI_ERR_ARG, recvdispls)
887 CHECK_NULL(8, MPI_ERR_TYPE, recvtypes)
889 aid_t pid = simgrid::s4u::this_actor::get_pid();
890 int size = comm->size();
891 for (int i = 0; i < size; i++) {
892 if(sendbuf != MPI_IN_PLACE){
893 CHECK_COUNT(2, sendcounts[i])
894 CHECK_TYPE(4, sendtypes[i])
895 CHECK_BUFFER(1, sendbuf, sendcounts[i], sendtypes[i])
897 CHECK_COUNT(6, recvcounts[i])
898 CHECK_TYPE(8, recvtypes[i])
899 CHECK_BUFFER(5, recvbuf, recvcounts[i], recvtypes[i])
901 CHECK_COLLECTIVE(comm, request == MPI_REQUEST_IGNORED ? "PMPI_Alltoallw" : "PMPI_Ialltoallw")
903 const SmpiBenchGuard suspend_bench;
907 auto trace_sendcounts = std::make_shared<std::vector<int>>();
908 auto trace_recvcounts = std::make_shared<std::vector<int>>();
909 trace_recvcounts->insert(trace_recvcounts->end(), &recvcounts[0], &recvcounts[size]);
911 const int* real_sendcounts = sendcounts;
912 const int* real_senddispls = senddispls;
913 const MPI_Datatype* real_sendtypes = sendtypes;
915 unsigned long maxsize = 0;
916 for (int i = 0; i < size; i++) { // copy data to avoid bad free
917 if (recvtypes[i] == MPI_DATATYPE_NULL)
919 recv_size += recvcounts[i] * recvtypes[i]->size();
920 if ((recvdispls[i] + (recvcounts[i] * recvtypes[i]->size())) > maxsize)
921 maxsize = recvdispls[i] + (recvcounts[i] * recvtypes[i]->size());
924 std::vector<unsigned char> tmp_sendbuf;
925 std::vector<int> tmp_sendcounts;
926 std::vector<int> tmp_senddispls;
927 std::vector<MPI_Datatype> tmp_sendtypes;
928 const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, maxsize, MPI_CHAR);
929 if (sendbuf == MPI_IN_PLACE) {
930 tmp_sendcounts.assign(recvcounts, recvcounts + size);
931 real_sendcounts = tmp_sendcounts.data();
932 tmp_senddispls.assign(recvdispls, recvdispls + size);
933 real_senddispls = tmp_senddispls.data();
934 tmp_sendtypes.assign(recvtypes, recvtypes + size);
935 real_sendtypes = tmp_sendtypes.data();
939 if(recvtypes[comm->rank()]->size() * recvcounts[comm->rank()] != real_sendtypes[comm->rank()]->size() * real_sendcounts[comm->rank()]){
940 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()]);
941 return MPI_ERR_TRUNCATE;
944 trace_sendcounts->insert(trace_sendcounts->end(), &real_sendcounts[0], &real_sendcounts[size]);
945 for (int i = 0; i < size; i++) { // copy data to avoid bad free
946 send_size += real_sendcounts[i] * real_sendtypes[i]->size();
949 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Alltoallw" : "PMPI_Ialltoallw",
950 new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "alltoallv" : "ialltoallv", -1,
951 send_size, trace_sendcounts, recv_size, trace_recvcounts,
952 simgrid::smpi::Datatype::encode(real_sendtypes[0]),
953 simgrid::smpi::Datatype::encode(recvtypes[0])));
956 if (request == MPI_REQUEST_IGNORED)
957 retval = simgrid::smpi::colls::alltoallw(real_sendbuf, real_sendcounts, real_senddispls, real_sendtypes, recvbuf,
958 recvcounts, recvdispls, recvtypes, comm);
960 retval = simgrid::smpi::colls::ialltoallw(real_sendbuf, real_sendcounts, real_senddispls, real_sendtypes, recvbuf,
961 recvcounts, recvdispls, recvtypes, comm, request);
963 TRACE_smpi_comm_out(pid);