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,
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)
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_sendcount, recvcount,
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);
152 int PMPI_Gatherv(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, const int *recvcounts, const int *displs,
153 MPI_Datatype recvtype, int root, MPI_Comm comm){
154 return PMPI_Igatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, root, comm, MPI_REQUEST_IGNORED);
157 int PMPI_Igatherv(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, const int* recvcounts, const int* displs,
158 MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request* request)
162 int rank = comm->rank();
163 if(sendbuf != MPI_IN_PLACE){
164 CHECK_TYPE(3, sendtype)
165 CHECK_COUNT(2, sendcount)
167 CHECK_BUFFER(1, sendbuf, sendcount, sendtype)
170 CHECK_NOT_IN_PLACE_ROOT(4, recvbuf)
171 CHECK_TYPE(6, recvtype)
172 CHECK_NULL(5, MPI_ERR_COUNT, recvcounts)
173 CHECK_NULL(6, MPI_ERR_ARG, displs)
175 CHECK_NOT_IN_PLACE_ROOT(1, sendbuf)
181 for (int i = 0; i < comm->size(); i++) {
182 CHECK_COUNT(5, recvcounts[i])
183 CHECK_BUFFER(4,recvbuf,recvcounts[i], recvtype)
187 const SmpiBenchGuard suspend_bench;
188 const void* real_sendbuf = sendbuf;
189 int real_sendcount = sendcount;
190 MPI_Datatype real_sendtype = sendtype;
191 if ((rank == root) && (sendbuf == MPI_IN_PLACE)) {
193 real_sendtype = recvtype;
196 aid_t pid = simgrid::s4u::this_actor::get_pid();
198 auto trace_recvcounts = std::make_shared<std::vector<int>>();
200 trace_recvcounts->insert(trace_recvcounts->end(), &recvcounts[0], &recvcounts[comm->size()]);
201 else //this is not significant outside of root, put 0 as we don't know if recvcounts is initialized
202 trace_recvcounts->insert(trace_recvcounts->end(), comm->size(), 0);
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,
208 nullptr, -1, 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);
221 int PMPI_Allgather(const void *sendbuf, int sendcount, MPI_Datatype sendtype,
222 void *recvbuf, int recvcount, MPI_Datatype recvtype, MPI_Comm comm){
223 return PMPI_Iallgather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm, MPI_REQUEST_IGNORED);
226 int PMPI_Iallgather(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
227 MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request)
232 int rank = comm->rank();
233 CHECK_NOT_IN_PLACE(4, recvbuf)
234 if(sendbuf != MPI_IN_PLACE){
235 CHECK_COUNT(2, sendcount)
236 CHECK_TYPE(3, sendtype)
238 CHECK_TYPE(6, recvtype)
239 CHECK_COUNT(5, recvcount)
240 CHECK_BUFFER(1, sendbuf, sendcount, sendtype)
241 CHECK_BUFFER(4, recvbuf, recvcount, recvtype)
244 if (sendbuf == MPI_IN_PLACE) {
245 sendbuf = static_cast<char*>(recvbuf) + recvtype->get_extent() * recvcount * comm->rank();
246 sendcount = recvcount;
250 if(recvtype->size() * recvcount != sendtype->size() * sendcount){
251 XBT_WARN("MPI_(I)Allgather : received size from each process differs from sent size : %zu vs %zu", recvtype->size() * recvcount, sendtype->size() * sendcount);
252 return MPI_ERR_TRUNCATE;
255 const SmpiBenchGuard suspend_bench;
257 aid_t pid = simgrid::s4u::this_actor::get_pid();
259 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Allgather" : "PMPI_Iallggather",
260 new simgrid::instr::CollTIData(
261 request == MPI_REQUEST_IGNORED ? "allgather" : "iallgather", -1, -1.0,
262 sendcount, recvcount,
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);
273 int PMPI_Allgatherv(const void *sendbuf, int sendcount, MPI_Datatype sendtype,
274 void *recvbuf, const int *recvcounts, const int *displs, MPI_Datatype recvtype, MPI_Comm comm){
275 return PMPI_Iallgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm, MPI_REQUEST_IGNORED);
278 int PMPI_Iallgatherv(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, const int* recvcounts, const int* displs,
279 MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request)
284 int rank = comm->rank();
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)
292 CHECK_BUFFER(1, sendbuf, sendcount, sendtype)
295 CHECK_NOT_IN_PLACE(4, recvbuf)
296 for (int i = 0; i < comm->size(); i++) {
297 CHECK_COUNT(5, recvcounts[i])
298 CHECK_BUFFER(4, recvbuf, recvcounts[i], recvtype)
301 const SmpiBenchGuard suspend_bench;
302 if (sendbuf == MPI_IN_PLACE) {
303 sendbuf = static_cast<char*>(recvbuf) + recvtype->get_extent() * displs[comm->rank()];
304 sendcount = recvcounts[comm->rank()];
307 aid_t pid = simgrid::s4u::this_actor::get_pid();
309 auto trace_recvcounts = std::make_shared<std::vector<int>>();
310 trace_recvcounts->insert(trace_recvcounts->end(), &recvcounts[0], &recvcounts[comm->size()]);
313 pid, request == MPI_REQUEST_IGNORED ? "PMPI_Allgatherv" : "PMPI_Iallgatherv",
314 new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "allgatherv" : "iallgatherv", -1,
316 -1, trace_recvcounts, simgrid::smpi::Datatype::encode(sendtype),
317 simgrid::smpi::Datatype::encode(recvtype)));
318 if (request == MPI_REQUEST_IGNORED)
319 simgrid::smpi::colls::allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm);
321 simgrid::smpi::colls::iallgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm,
324 TRACE_smpi_comm_out(pid);
328 int PMPI_Scatter(const void *sendbuf, int sendcount, MPI_Datatype sendtype,
329 void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm){
330 return PMPI_Iscatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm, MPI_REQUEST_IGNORED);
333 int PMPI_Iscatter(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
334 MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request* request)
338 int rank = comm->rank();
341 CHECK_NOT_IN_PLACE_ROOT(1, sendbuf)
342 CHECK_COUNT(2, sendcount)
343 CHECK_TYPE(3, sendtype)
344 CHECK_BUFFER(1, sendbuf, sendcount, sendtype)
346 CHECK_NOT_IN_PLACE_ROOT(4, recvbuf)
348 if(recvbuf != MPI_IN_PLACE){
349 CHECK_COUNT(5, recvcount)
350 CHECK_TYPE(6, recvtype)
351 CHECK_BUFFER(4, recvbuf, recvcount, 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;
366 const SmpiBenchGuard suspend_bench;
368 aid_t 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 sendcount, recvcount,
374 simgrid::smpi::Datatype::encode(sendtype), simgrid::smpi::Datatype::encode(recvtype)));
375 if (request == MPI_REQUEST_IGNORED)
376 simgrid::smpi::colls::scatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm);
378 simgrid::smpi::colls::iscatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm, request);
380 TRACE_smpi_comm_out(pid);
384 int PMPI_Scatterv(const void *sendbuf, const int *sendcounts, const int *displs,
385 MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm){
386 return PMPI_Iscatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm, MPI_REQUEST_IGNORED);
389 int PMPI_Iscatterv(const void* sendbuf, const int* sendcounts, const int* displs, MPI_Datatype sendtype, void* recvbuf, int recvcount,
390 MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request* request)
394 int rank = comm->rank();
395 if(recvbuf != MPI_IN_PLACE){
396 CHECK_COUNT(5, recvcount)
397 CHECK_TYPE(7, recvtype)
398 CHECK_BUFFER(4, recvbuf, recvcount, recvtype)
404 CHECK_NOT_IN_PLACE_ROOT(1, sendbuf)
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_COUNT(2, sendcounts[i])
410 CHECK_BUFFER(1, sendbuf, sendcounts[i], sendtype)
412 if (recvbuf == MPI_IN_PLACE) {
414 recvcount = sendcounts[rank];
417 CHECK_NOT_IN_PLACE_ROOT(4, recvbuf)
420 const SmpiBenchGuard suspend_bench;
422 aid_t pid = simgrid::s4u::this_actor::get_pid();
424 auto trace_sendcounts = std::make_shared<std::vector<int>>();
426 trace_sendcounts->insert(trace_sendcounts->end(), &sendcounts[0], &sendcounts[comm->size()]);
427 else //this is not significant outside of root, put 0 as we don't know if sendcounts is initialized
428 trace_sendcounts->insert(trace_sendcounts->end(), comm->size(), 0);
431 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Scatterv" : "PMPI_Iscatterv",
432 new simgrid::instr::VarCollTIData(
433 request == MPI_REQUEST_IGNORED ? "scatterv" : "iscatterv", root, -1,
434 trace_sendcounts, recvcount,
435 nullptr, simgrid::smpi::Datatype::encode(sendtype),
436 simgrid::smpi::Datatype::encode(recvtype)));
437 if (request == MPI_REQUEST_IGNORED)
438 simgrid::smpi::colls::scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm);
440 simgrid::smpi::colls::iscatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm,
443 TRACE_smpi_comm_out(pid);
447 int PMPI_Reduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
449 return PMPI_Ireduce(sendbuf, recvbuf, count, datatype, op, root, comm, MPI_REQUEST_IGNORED);
452 int PMPI_Ireduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm, MPI_Request* request)
456 int rank = comm->rank();
457 CHECK_TYPE(4, datatype)
458 CHECK_COUNT(3, count)
459 CHECK_BUFFER(1, sendbuf, count, datatype)
462 CHECK_NOT_IN_PLACE(2, recvbuf)
463 CHECK_BUFFER(5, recvbuf, count, datatype)
465 CHECK_OP(5, op, datatype)
469 const SmpiBenchGuard suspend_bench;
470 aid_t 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,
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);
485 int PMPI_Reduce_local(const void* inbuf, void* inoutbuf, int count, MPI_Datatype datatype, MPI_Op op)
489 CHECK_TYPE(4, datatype)
490 CHECK_COUNT(3, count)
491 CHECK_BUFFER(1, inbuf, count, datatype)
492 CHECK_BUFFER(2, inoutbuf, count, datatype)
493 CHECK_OP(5, op, datatype)
495 const SmpiBenchGuard suspend_bench;
496 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)
510 int rank = comm->rank();
511 CHECK_NOT_IN_PLACE(2, recvbuf)
512 CHECK_TYPE(4, datatype)
513 CHECK_OP(5, op, datatype)
514 CHECK_COUNT(3, count)
515 CHECK_BUFFER(1, sendbuf, count, datatype)
516 CHECK_BUFFER(2, recvbuf, count, datatype)
519 const SmpiBenchGuard suspend_bench;
520 std::vector<unsigned char> tmp_sendbuf;
521 const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, count, datatype);
523 aid_t pid = simgrid::s4u::this_actor::get_pid();
525 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Allreduce" : "PMPI_Iallreduce",
526 new simgrid::instr::CollTIData(request == MPI_REQUEST_IGNORED ? "allreduce" : "iallreduce", -1, 0,
528 simgrid::smpi::Datatype::encode(datatype), ""));
530 if (request == MPI_REQUEST_IGNORED)
531 simgrid::smpi::colls::allreduce(real_sendbuf, recvbuf, count, datatype, op, comm);
533 simgrid::smpi::colls::iallreduce(real_sendbuf, recvbuf, count, datatype, op, comm, request);
535 TRACE_smpi_comm_out(pid);
539 int PMPI_Scan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
541 return PMPI_Iscan(sendbuf, recvbuf, count, datatype, op, comm, MPI_REQUEST_IGNORED);
544 int PMPI_Iscan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request* request)
549 CHECK_TYPE(4, datatype)
550 CHECK_COUNT(3, count)
551 CHECK_BUFFER(1,sendbuf,count, datatype)
552 CHECK_BUFFER(2,recvbuf,count, datatype)
554 CHECK_OP(5, op, datatype)
556 const SmpiBenchGuard suspend_bench;
557 aid_t pid = simgrid::s4u::this_actor::get_pid();
558 std::vector<unsigned char> tmp_sendbuf;
559 const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, count, datatype);
561 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Scan" : "PMPI_Iscan",
562 new simgrid::instr::CollTIData(request == MPI_REQUEST_IGNORED ? "scan" : "iscan", -1, 0.0,
563 count, 0, simgrid::smpi::Datatype::encode(datatype), ""));
566 if (request == MPI_REQUEST_IGNORED)
567 retval = simgrid::smpi::colls::scan(real_sendbuf, recvbuf, count, datatype, op, comm);
569 retval = simgrid::smpi::colls::iscan(real_sendbuf, recvbuf, count, datatype, op, comm, request);
571 TRACE_smpi_comm_out(pid);
575 int PMPI_Exscan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
577 return PMPI_Iexscan(sendbuf, recvbuf, count, datatype, op, comm, MPI_REQUEST_IGNORED);
580 int PMPI_Iexscan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request* request){
584 CHECK_TYPE(4, datatype)
585 CHECK_COUNT(3, count)
586 CHECK_BUFFER(1, sendbuf, count, datatype)
587 CHECK_BUFFER(2, recvbuf, count, datatype)
589 CHECK_OP(5, op, datatype)
591 const SmpiBenchGuard suspend_bench;
592 aid_t pid = simgrid::s4u::this_actor::get_pid();
593 std::vector<unsigned char> tmp_sendbuf;
594 const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, count, datatype);
596 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Exscan" : "PMPI_Iexscan",
597 new simgrid::instr::CollTIData(request == MPI_REQUEST_IGNORED ? "exscan" : "iexscan", -1, 0.0,
598 count, 0, simgrid::smpi::Datatype::encode(datatype), ""));
601 if (request == MPI_REQUEST_IGNORED)
602 retval = simgrid::smpi::colls::exscan(real_sendbuf, recvbuf, count, datatype, op, comm);
604 retval = simgrid::smpi::colls::iexscan(real_sendbuf, recvbuf, count, datatype, op, comm, request);
606 TRACE_smpi_comm_out(pid);
610 int PMPI_Reduce_scatter(const void *sendbuf, void *recvbuf, const int *recvcounts, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
612 return PMPI_Ireduce_scatter(sendbuf, recvbuf, recvcounts, datatype, op, comm, MPI_REQUEST_IGNORED);
615 int PMPI_Ireduce_scatter(const void *sendbuf, void *recvbuf, const int *recvcounts, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request *request)
620 int rank = comm->rank();
621 CHECK_NOT_IN_PLACE(2, recvbuf)
622 CHECK_TYPE(4, datatype)
623 CHECK_NULL(3, MPI_ERR_COUNT, recvcounts)
625 CHECK_OP(5, op, datatype)
626 for (int i = 0; i < comm->size(); i++) {
627 CHECK_COUNT(3, recvcounts[i])
628 CHECK_BUFFER(1, sendbuf, recvcounts[i], datatype)
629 CHECK_BUFFER(2, recvbuf, recvcounts[i], datatype)
632 const SmpiBenchGuard suspend_bench;
633 aid_t pid = simgrid::s4u::this_actor::get_pid();
634 auto trace_recvcounts = std::make_shared<std::vector<int>>();
635 trace_recvcounts->insert(trace_recvcounts->end(), &recvcounts[0], &recvcounts[comm->size()]);
639 for (int i = 0; i < comm->size(); i++) { // copy data to avoid bad free
640 totalcount += recvcounts[i];
642 std::vector<unsigned char> tmp_sendbuf;
643 const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, totalcount, datatype);
645 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Reduce_scatter" : "PMPI_Ireduce_scatter",
646 new simgrid::instr::VarCollTIData(
647 request == MPI_REQUEST_IGNORED ? "reducescatter" : "ireducescatter", -1, -1, nullptr,
648 0, trace_recvcounts, simgrid::smpi::Datatype::encode(datatype), ""));
650 if (request == MPI_REQUEST_IGNORED)
651 simgrid::smpi::colls::reduce_scatter(real_sendbuf, recvbuf, recvcounts, datatype, op, comm);
653 simgrid::smpi::colls::ireduce_scatter(real_sendbuf, recvbuf, recvcounts, datatype, op, comm, request);
655 TRACE_smpi_comm_out(pid);
659 int PMPI_Reduce_scatter_block(const void *sendbuf, void *recvbuf, int recvcount,
660 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
662 return PMPI_Ireduce_scatter_block(sendbuf, recvbuf, recvcount, datatype, op, comm, MPI_REQUEST_IGNORED);
665 int PMPI_Ireduce_scatter_block(const void* sendbuf, void* recvbuf, int recvcount, MPI_Datatype datatype, MPI_Op op,
666 MPI_Comm comm, MPI_Request* request)
671 CHECK_TYPE(4, datatype)
672 CHECK_COUNT(3, recvcount)
673 CHECK_BUFFER(1, sendbuf, recvcount, datatype)
674 CHECK_BUFFER(2, recvbuf, recvcount, datatype)
676 CHECK_OP(5, op, datatype)
678 const SmpiBenchGuard suspend_bench;
679 int count = comm->size();
681 aid_t pid = simgrid::s4u::this_actor::get_pid();
682 auto trace_recvcounts = std::make_shared<std::vector<int>>(recvcount);
684 std::vector<unsigned char> tmp_sendbuf;
685 const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, recvcount * count, datatype);
688 pid, request == MPI_REQUEST_IGNORED ? "PMPI_Reduce_scatter_block" : "PMPI_Ireduce_scatter_block",
689 new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "reducescatter" : "ireducescatter", -1, -1,
690 nullptr, -1, trace_recvcounts, simgrid::smpi::Datatype::encode(datatype), ""));
692 std::vector<int> recvcounts(count);
693 for (int i = 0; i < count; i++)
694 recvcounts[i] = recvcount;
695 if (request == MPI_REQUEST_IGNORED)
696 simgrid::smpi::colls::reduce_scatter(real_sendbuf, recvbuf, recvcounts.data(), datatype, op, comm);
698 simgrid::smpi::colls::ireduce_scatter(real_sendbuf, recvbuf, recvcounts.data(), datatype, op, comm, request);
700 TRACE_smpi_comm_out(pid);
704 int PMPI_Alltoall(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
705 MPI_Datatype recvtype, MPI_Comm comm){
706 return PMPI_Ialltoall(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm, MPI_REQUEST_IGNORED);
709 int PMPI_Ialltoall(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
710 MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request)
715 if(sendbuf != MPI_IN_PLACE){
716 CHECK_TYPE(3, sendtype)
717 CHECK_COUNT(2, sendcount)
718 CHECK_BUFFER(1, sendbuf, sendcount, sendtype)
720 CHECK_TYPE(6, recvtype)
721 CHECK_COUNT(5, recvcount)
722 CHECK_COUNT(5, recvcount)
723 CHECK_BUFFER(4, recvbuf, recvcount, recvtype)
726 aid_t pid = simgrid::s4u::this_actor::get_pid();
727 int real_sendcount = sendcount;
728 MPI_Datatype real_sendtype = sendtype;
730 std::vector<unsigned char> tmp_sendbuf;
731 const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, recvcount * comm->size(), recvtype);
733 if (sendbuf == MPI_IN_PLACE) {
734 real_sendcount = recvcount;
735 real_sendtype = recvtype;
738 if(recvtype->size() * recvcount != real_sendtype->size() * real_sendcount){
739 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);
740 return MPI_ERR_TRUNCATE;
743 const SmpiBenchGuard suspend_bench;
745 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Alltoall" : "PMPI_Ialltoall",
746 new simgrid::instr::CollTIData(
747 request == MPI_REQUEST_IGNORED ? "alltoall" : "ialltoall", -1, -1.0,
748 real_sendcount, recvcount,
749 simgrid::smpi::Datatype::encode(real_sendtype), simgrid::smpi::Datatype::encode(recvtype)));
751 if (request == MPI_REQUEST_IGNORED)
753 simgrid::smpi::colls::alltoall(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcount, recvtype, comm);
755 retval = simgrid::smpi::colls::ialltoall(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcount, recvtype,
758 TRACE_smpi_comm_out(pid);
762 int PMPI_Alltoallv(const void* sendbuf, const int* sendcounts, const int* senddispls, MPI_Datatype sendtype, void* recvbuf,
763 const int* recvcounts, const int* recvdispls, MPI_Datatype recvtype, MPI_Comm comm)
765 return PMPI_Ialltoallv(sendbuf, sendcounts, senddispls, sendtype, recvbuf, recvcounts, recvdispls, recvtype, comm, MPI_REQUEST_IGNORED);
768 int PMPI_Ialltoallv(const void* sendbuf, const int* sendcounts, const int* senddispls, MPI_Datatype sendtype, void* recvbuf,
769 const int* recvcounts, const int* recvdispls, MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request)
774 if(sendbuf != MPI_IN_PLACE){
775 CHECK_NULL(2, MPI_ERR_COUNT, sendcounts)
776 CHECK_NULL(3, MPI_ERR_ARG, senddispls)
777 CHECK_TYPE(4, sendtype)
779 CHECK_TYPE(8, recvtype)
780 CHECK_NULL(6, MPI_ERR_COUNT, recvcounts)
781 CHECK_NULL(7, MPI_ERR_ARG, recvdispls)
784 aid_t pid = simgrid::s4u::this_actor::get_pid();
785 int size = comm->size();
786 for (int i = 0; i < size; i++) {
787 if(sendbuf != MPI_IN_PLACE){
788 CHECK_BUFFER(1, sendbuf, sendcounts[i], sendtype)
789 CHECK_COUNT(2, sendcounts[i])
791 CHECK_BUFFER(5, recvbuf, recvcounts[i], recvtype)
792 CHECK_COUNT(6, recvcounts[i])
795 const SmpiBenchGuard suspend_bench;
798 auto trace_sendcounts = std::make_shared<std::vector<int>>();
799 auto trace_recvcounts = std::make_shared<std::vector<int>>();
800 trace_recvcounts->insert(trace_recvcounts->end(), &recvcounts[0], &recvcounts[size]);
802 int dt_size_recv = recvtype->size();
804 const int* real_sendcounts = sendcounts;
805 const int* real_senddispls = senddispls;
806 MPI_Datatype real_sendtype = sendtype;
808 std::vector<unsigned char> tmp_sendbuf;
809 std::vector<int> tmp_sendcounts;
810 std::vector<int> tmp_senddispls;
811 const void* real_sendbuf;
813 if (sendbuf == MPI_IN_PLACE) {
814 tmp_sendcounts.assign(recvcounts, recvcounts + size);
815 real_sendcounts = tmp_sendcounts.data();
816 tmp_senddispls.assign(recvdispls, recvdispls + size);
817 real_senddispls = tmp_senddispls.data();
818 real_sendtype = recvtype;
821 for (int i = 0; i < size; i++) { // copy data to avoid bad free
822 send_size += real_sendcounts[i] ;
823 recv_size += recvcounts[i];
824 if (((recvdispls[i] + recvcounts[i]) * dt_size_recv) > maxsize)
825 maxsize = (recvdispls[i] + recvcounts[i]) * dt_size_recv;
827 real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, maxsize, MPI_CHAR);
829 if(recvtype->size() * recvcounts[comm->rank()] != real_sendtype->size() * real_sendcounts[comm->rank()]){
830 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()]);
831 return MPI_ERR_TRUNCATE;
834 trace_sendcounts->insert(trace_sendcounts->end(), &real_sendcounts[0], &real_sendcounts[size]);
836 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Alltoallv" : "PMPI_Ialltoallv",
837 new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "alltoallv" : "ialltoallv", -1,
838 send_size, trace_sendcounts, recv_size, trace_recvcounts,
839 simgrid::smpi::Datatype::encode(real_sendtype),
840 simgrid::smpi::Datatype::encode(recvtype)));
843 if (request == MPI_REQUEST_IGNORED)
844 retval = simgrid::smpi::colls::alltoallv(real_sendbuf, real_sendcounts, real_senddispls, real_sendtype, recvbuf,
845 recvcounts, recvdispls, recvtype, comm);
847 retval = simgrid::smpi::colls::ialltoallv(real_sendbuf, real_sendcounts, real_senddispls, real_sendtype, recvbuf,
848 recvcounts, recvdispls, recvtype, comm, request);
850 TRACE_smpi_comm_out(pid);
854 int PMPI_Alltoallw(const void* sendbuf, const int* sendcounts, const int* senddispls, const MPI_Datatype* sendtypes, void* recvbuf,
855 const int* recvcounts, const int* recvdispls, const MPI_Datatype* recvtypes, MPI_Comm comm)
857 return PMPI_Ialltoallw(sendbuf, sendcounts, senddispls, sendtypes, recvbuf, recvcounts, recvdispls, recvtypes, comm, MPI_REQUEST_IGNORED);
860 int PMPI_Ialltoallw(const void* sendbuf, const int* sendcounts, const int* senddispls, const MPI_Datatype* sendtypes, void* recvbuf,
861 const int* recvcounts, const int* recvdispls, const MPI_Datatype* recvtypes, MPI_Comm comm, MPI_Request* request)
866 if(sendbuf != MPI_IN_PLACE){
867 CHECK_NULL(2, MPI_ERR_COUNT, sendcounts)
868 CHECK_NULL(3, MPI_ERR_ARG, senddispls)
869 CHECK_NULL(4, MPI_ERR_TYPE, sendtypes)
871 CHECK_NULL(6, MPI_ERR_COUNT, recvcounts)
872 CHECK_NULL(7, MPI_ERR_ARG, recvdispls)
873 CHECK_NULL(8, MPI_ERR_TYPE, recvtypes)
875 aid_t pid = simgrid::s4u::this_actor::get_pid();
876 int size = comm->size();
877 for (int i = 0; i < size; i++) {
878 if(sendbuf != MPI_IN_PLACE){
879 CHECK_COUNT(2, sendcounts[i])
880 CHECK_TYPE(4, sendtypes[i])
881 CHECK_BUFFER(1, sendbuf, sendcounts[i], sendtypes[i])
883 CHECK_COUNT(6, recvcounts[i])
884 CHECK_TYPE(8, recvtypes[i])
885 CHECK_BUFFER(5, recvbuf, recvcounts[i], recvtypes[i])
888 const SmpiBenchGuard suspend_bench;
892 auto trace_sendcounts = std::make_shared<std::vector<int>>();
893 auto trace_recvcounts = std::make_shared<std::vector<int>>();
894 trace_recvcounts->insert(trace_recvcounts->end(), &recvcounts[0], &recvcounts[size]);
896 const int* real_sendcounts = sendcounts;
897 const int* real_senddispls = senddispls;
898 const MPI_Datatype* real_sendtypes = sendtypes;
900 unsigned long maxsize = 0;
901 for (int i = 0; i < size; i++) { // copy data to avoid bad free
902 if (recvtypes[i] == MPI_DATATYPE_NULL)
904 recv_size += recvcounts[i] * recvtypes[i]->size();
905 if ((recvdispls[i] + (recvcounts[i] * recvtypes[i]->size())) > maxsize)
906 maxsize = recvdispls[i] + (recvcounts[i] * recvtypes[i]->size());
909 std::vector<unsigned char> tmp_sendbuf;
910 std::vector<int> tmp_sendcounts;
911 std::vector<int> tmp_senddispls;
912 std::vector<MPI_Datatype> tmp_sendtypes;
913 const void* real_sendbuf = smpi_get_in_place_buf(sendbuf, recvbuf, tmp_sendbuf, maxsize, MPI_CHAR);
914 if (sendbuf == MPI_IN_PLACE) {
915 tmp_sendcounts.assign(recvcounts, recvcounts + size);
916 real_sendcounts = tmp_sendcounts.data();
917 tmp_senddispls.assign(recvdispls, recvdispls + size);
918 real_senddispls = tmp_senddispls.data();
919 tmp_sendtypes.assign(recvtypes, recvtypes + size);
920 real_sendtypes = tmp_sendtypes.data();
924 if(recvtypes[comm->rank()]->size() * recvcounts[comm->rank()] != real_sendtypes[comm->rank()]->size() * real_sendcounts[comm->rank()]){
925 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()]);
926 return MPI_ERR_TRUNCATE;
929 trace_sendcounts->insert(trace_sendcounts->end(), &real_sendcounts[0], &real_sendcounts[size]);
930 for (int i = 0; i < size; i++) { // copy data to avoid bad free
931 send_size += real_sendcounts[i] * real_sendtypes[i]->size();
934 TRACE_smpi_comm_in(pid, request == MPI_REQUEST_IGNORED ? "PMPI_Alltoallw" : "PMPI_Ialltoallw",
935 new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "alltoallv" : "ialltoallv", -1,
936 send_size, trace_sendcounts, recv_size, trace_recvcounts,
937 simgrid::smpi::Datatype::encode(real_sendtypes[0]),
938 simgrid::smpi::Datatype::encode(recvtypes[0])));
941 if (request == MPI_REQUEST_IGNORED)
942 retval = simgrid::smpi::colls::alltoallw(real_sendbuf, real_sendcounts, real_senddispls, real_sendtypes, recvbuf,
943 recvcounts, recvdispls, recvtypes, comm);
945 retval = simgrid::smpi::colls::ialltoallw(real_sendbuf, real_sendcounts, real_senddispls, real_sendtypes, recvbuf,
946 recvcounts, recvdispls, recvtypes, comm, request);
948 TRACE_smpi_comm_out(pid);