1 /* Copyright (c) 2007-2019. 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"
14 XBT_LOG_EXTERNAL_DEFAULT_CATEGORY(smpi_pmpi);
16 /* PMPI User level calls */
18 int PMPI_Barrier(MPI_Comm comm)
20 return PMPI_Ibarrier(comm, MPI_REQUEST_IGNORED);
23 int PMPI_Ibarrier(MPI_Comm comm, MPI_Request *request)
25 if (comm == MPI_COMM_NULL)
27 if (request == nullptr)
31 int rank = simgrid::s4u::this_actor::get_pid();
32 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Barrier" : "PMPI_Ibarrier",
33 new simgrid::instr::NoOpTIData(request == MPI_REQUEST_IGNORED ? "barrier" : "ibarrier"));
34 if (request == MPI_REQUEST_IGNORED) {
35 simgrid::smpi::Colls::barrier(comm);
36 // Barrier can be used to synchronize RMA calls. Finish all requests from comm before.
37 comm->finish_rma_calls();
39 simgrid::smpi::Colls::ibarrier(comm, request);
41 TRACE_smpi_comm_out(rank);
46 int PMPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
48 return PMPI_Ibcast(buf, count, datatype, root, comm, MPI_REQUEST_IGNORED);
51 int PMPI_Ibcast(void *buf, int count, MPI_Datatype datatype,
52 int root, MPI_Comm comm, MPI_Request* request)
54 if (comm == MPI_COMM_NULL)
56 if (buf == nullptr && count > 0)
57 return MPI_ERR_BUFFER;
58 if (datatype == MPI_DATATYPE_NULL || not datatype->is_valid())
62 if (root < 0 || root >= comm->size())
64 if (request == nullptr)
68 int rank = simgrid::s4u::this_actor::get_pid();
69 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Bcast" : "PMPI_Ibcast",
70 new simgrid::instr::CollTIData(request == MPI_REQUEST_IGNORED ? "bcast" : "ibcast", root, -1.0,
71 datatype->is_replayable() ? count : count * datatype->size(), -1,
72 simgrid::smpi::Datatype::encode(datatype), ""));
73 if (comm->size() > 1) {
74 if (request == MPI_REQUEST_IGNORED)
75 simgrid::smpi::Colls::bcast(buf, count, datatype, root, comm);
77 simgrid::smpi::Colls::ibcast(buf, count, datatype, root, comm, request);
79 if (request != MPI_REQUEST_IGNORED)
80 *request = MPI_REQUEST_NULL;
83 TRACE_smpi_comm_out(rank);
88 int PMPI_Gather(const void *sendbuf, int sendcount, MPI_Datatype sendtype,void *recvbuf, int recvcount, MPI_Datatype recvtype,
89 int root, MPI_Comm comm){
90 return PMPI_Igather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm, MPI_REQUEST_IGNORED);
93 int PMPI_Igather(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
94 MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request* request)
96 if (comm == MPI_COMM_NULL)
98 if ((sendbuf == nullptr) || ((comm->rank() == root) && recvbuf == nullptr))
99 return MPI_ERR_BUFFER;
100 if (((sendbuf != MPI_IN_PLACE && sendcount > 0) && (sendtype == MPI_DATATYPE_NULL)) ||
101 ((comm->rank() == root) && (recvtype == MPI_DATATYPE_NULL)))
103 if (((sendbuf != MPI_IN_PLACE) && (sendcount < 0)) || ((comm->rank() == root) && (recvcount < 0)))
104 return MPI_ERR_COUNT;
105 if (root < 0 || root >= comm->size())
107 if (request == nullptr)
111 const void* real_sendbuf = sendbuf;
112 int real_sendcount = sendcount;
113 MPI_Datatype real_sendtype = sendtype;
114 if ((comm->rank() == root) && (sendbuf == MPI_IN_PLACE)) {
116 real_sendtype = recvtype;
118 int rank = simgrid::s4u::this_actor::get_pid();
120 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Gather" : "PMPI_Igather",
121 new simgrid::instr::CollTIData(
122 request == MPI_REQUEST_IGNORED ? "gather" : "igather", root, -1.0,
123 real_sendtype->is_replayable() ? real_sendcount : real_sendcount * real_sendtype->size(),
124 (comm->rank() != root || recvtype->is_replayable()) ? recvcount : recvcount * recvtype->size(),
125 simgrid::smpi::Datatype::encode(real_sendtype), simgrid::smpi::Datatype::encode(recvtype)));
126 if (request == MPI_REQUEST_IGNORED)
127 simgrid::smpi::Colls::gather(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcount, recvtype, root, comm);
129 simgrid::smpi::Colls::igather(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcount, recvtype, root, comm,
132 TRACE_smpi_comm_out(rank);
137 int PMPI_Gatherv(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, const int *recvcounts, const int *displs,
138 MPI_Datatype recvtype, int root, MPI_Comm comm){
139 return PMPI_Igatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, root, comm, MPI_REQUEST_IGNORED);
142 int PMPI_Igatherv(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, const int* recvcounts, const int* displs,
143 MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request* request)
145 if (comm == MPI_COMM_NULL)
147 if ((sendbuf == nullptr && sendcount > 0) || ((comm->rank() == root) && recvbuf == nullptr))
148 return MPI_ERR_BUFFER;
149 if (((sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
150 ((comm->rank() == root) && (recvtype == MPI_DATATYPE_NULL)))
152 if ((sendbuf != MPI_IN_PLACE) && (sendcount < 0))
153 return MPI_ERR_COUNT;
154 if ((comm->rank() == root) && (recvcounts == nullptr || displs == nullptr))
156 if (root < 0 || root >= comm->size())
158 if (request == nullptr)
161 for (int i = 0; i < comm->size(); i++) {
162 if ((comm->rank() == root) && (recvcounts[i] < 0))
163 return MPI_ERR_COUNT;
167 const void* real_sendbuf = sendbuf;
168 int real_sendcount = sendcount;
169 MPI_Datatype real_sendtype = sendtype;
170 if ((comm->rank() == root) && (sendbuf == MPI_IN_PLACE)) {
172 real_sendtype = recvtype;
175 int rank = simgrid::s4u::this_actor::get_pid();
176 int dt_size_recv = recvtype->is_replayable() ? 1 : recvtype->size();
178 std::vector<int>* trace_recvcounts = new std::vector<int>;
179 if (comm->rank() == root) {
180 for (int i = 0; i < comm->size(); i++) // copy data to avoid bad free
181 trace_recvcounts->push_back(recvcounts[i] * dt_size_recv);
184 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Gatherv" : "PMPI_Igatherv",
185 new simgrid::instr::VarCollTIData(
186 request == MPI_REQUEST_IGNORED ? "gatherv" : "igatherv", root,
187 real_sendtype->is_replayable() ? real_sendcount : real_sendcount * real_sendtype->size(),
188 nullptr, dt_size_recv, trace_recvcounts, simgrid::smpi::Datatype::encode(real_sendtype),
189 simgrid::smpi::Datatype::encode(recvtype)));
190 if (request == MPI_REQUEST_IGNORED)
191 simgrid::smpi::Colls::gatherv(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcounts, displs, recvtype,
194 simgrid::smpi::Colls::igatherv(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcounts, displs, recvtype,
195 root, comm, request);
197 TRACE_smpi_comm_out(rank);
202 int PMPI_Allgather(const void *sendbuf, int sendcount, MPI_Datatype sendtype,
203 void *recvbuf, int recvcount, MPI_Datatype recvtype, MPI_Comm comm){
204 return PMPI_Iallgather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm, MPI_REQUEST_IGNORED);
207 int PMPI_Iallgather(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
208 MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request)
210 if (comm == MPI_COMM_NULL)
212 if ((sendbuf == nullptr && sendcount > 0) || (recvbuf == nullptr))
213 return MPI_ERR_BUFFER;
214 if (((sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) || (recvtype == MPI_DATATYPE_NULL))
216 if (((sendbuf != MPI_IN_PLACE) && (sendcount < 0)) || (recvcount < 0))
217 return MPI_ERR_COUNT;
218 if (request == nullptr)
222 if (sendbuf == MPI_IN_PLACE) {
223 sendbuf = static_cast<char*>(recvbuf) + recvtype->get_extent() * recvcount * comm->rank();
224 sendcount = recvcount;
227 int rank = simgrid::s4u::this_actor::get_pid();
229 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Allgather" : "PMPI_Iallggather",
230 new simgrid::instr::CollTIData(
231 request == MPI_REQUEST_IGNORED ? "allgather" : "iallgather", -1, -1.0,
232 sendtype->is_replayable() ? sendcount : sendcount * sendtype->size(),
233 recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(),
234 simgrid::smpi::Datatype::encode(sendtype), simgrid::smpi::Datatype::encode(recvtype)));
235 if (request == MPI_REQUEST_IGNORED)
236 simgrid::smpi::Colls::allgather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
238 simgrid::smpi::Colls::iallgather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm, request);
240 TRACE_smpi_comm_out(rank);
245 int PMPI_Allgatherv(const void *sendbuf, int sendcount, MPI_Datatype sendtype,
246 void *recvbuf, const int *recvcounts, const int *displs, MPI_Datatype recvtype, MPI_Comm comm){
247 return PMPI_Iallgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm, MPI_REQUEST_IGNORED);
250 int PMPI_Iallgatherv(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, const int* recvcounts, const int* displs,
251 MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request)
253 if (comm == MPI_COMM_NULL)
255 if ((sendbuf == nullptr && sendcount > 0) || (recvbuf == nullptr))
256 return MPI_ERR_BUFFER;
257 if (((sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) || (recvtype == MPI_DATATYPE_NULL))
259 if ((sendbuf != MPI_IN_PLACE) && (sendcount < 0))
260 return MPI_ERR_COUNT;
261 if (recvcounts == nullptr || displs == nullptr)
263 if (request == nullptr)
266 for (int i = 0; i < comm->size(); i++) { // copy data to avoid bad free
267 if (recvcounts[i] < 0)
268 return MPI_ERR_COUNT;
272 if (sendbuf == MPI_IN_PLACE) {
273 sendbuf = static_cast<char*>(recvbuf) + recvtype->get_extent() * displs[comm->rank()];
274 sendcount = recvcounts[comm->rank()];
277 int rank = simgrid::s4u::this_actor::get_pid();
278 int dt_size_recv = recvtype->is_replayable() ? 1 : recvtype->size();
280 std::vector<int>* trace_recvcounts = new std::vector<int>;
281 for (int i = 0; i < comm->size(); i++) { // copy data to avoid bad free
282 trace_recvcounts->push_back(recvcounts[i] * dt_size_recv);
286 rank, request == MPI_REQUEST_IGNORED ? "PMPI_Allgatherv" : "PMPI_Iallgatherv",
287 new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "allgatherv" : "iallgatherv", -1,
288 sendtype->is_replayable() ? sendcount : sendcount * sendtype->size(), nullptr,
289 dt_size_recv, trace_recvcounts, simgrid::smpi::Datatype::encode(sendtype),
290 simgrid::smpi::Datatype::encode(recvtype)));
291 if (request == MPI_REQUEST_IGNORED)
292 simgrid::smpi::Colls::allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm);
294 simgrid::smpi::Colls::iallgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm,
297 TRACE_smpi_comm_out(rank);
302 int PMPI_Scatter(const void *sendbuf, int sendcount, MPI_Datatype sendtype,
303 void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm){
304 return PMPI_Iscatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm, MPI_REQUEST_IGNORED);
307 int PMPI_Iscatter(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
308 MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request* request)
310 if (comm == MPI_COMM_NULL)
312 if (((comm->rank() == root) && (sendtype == MPI_DATATYPE_NULL || not sendtype->is_valid())) ||
313 ((recvbuf != MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL || not recvtype->is_valid())))
315 if (((comm->rank() == root) && (sendcount < 0)) || ((recvbuf != MPI_IN_PLACE) && (recvcount < 0)))
316 return MPI_ERR_COUNT;
317 if ((sendbuf == recvbuf) || ((comm->rank() == root) && sendcount > 0 && (sendbuf == nullptr)) ||
318 (recvcount > 0 && recvbuf == nullptr))
319 return MPI_ERR_BUFFER;
320 if (root < 0 || root >= comm->size())
322 if (request == nullptr)
326 if (recvbuf == MPI_IN_PLACE) {
328 recvcount = sendcount;
330 int rank = simgrid::s4u::this_actor::get_pid();
332 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Scatter" : "PMPI_Iscatter",
333 new simgrid::instr::CollTIData(
334 request == MPI_REQUEST_IGNORED ? "scatter" : "iscatter", root, -1.0,
335 (comm->rank() != root || sendtype->is_replayable()) ? sendcount : sendcount * sendtype->size(),
336 recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(),
337 simgrid::smpi::Datatype::encode(sendtype), simgrid::smpi::Datatype::encode(recvtype)));
338 if (request == MPI_REQUEST_IGNORED)
339 simgrid::smpi::Colls::scatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm);
341 simgrid::smpi::Colls::iscatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm, request);
343 TRACE_smpi_comm_out(rank);
348 int PMPI_Scatterv(const void *sendbuf, const int *sendcounts, const int *displs,
349 MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm){
350 return PMPI_Iscatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm, MPI_REQUEST_IGNORED);
353 int PMPI_Iscatterv(const void* sendbuf, const int* sendcounts, const int* displs, MPI_Datatype sendtype, void* recvbuf, int recvcount,
354 MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request* request)
356 if (comm == MPI_COMM_NULL)
358 if (sendcounts == nullptr || displs == nullptr)
360 if (((comm->rank() == root) && (sendtype == MPI_DATATYPE_NULL)) ||
361 ((recvbuf != MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL)))
363 if (request == nullptr)
365 if (recvbuf != MPI_IN_PLACE && recvcount < 0)
366 return MPI_ERR_COUNT;
367 if (root < 0 || root >= comm->size())
370 if (comm->rank() == root) {
371 if (recvbuf == MPI_IN_PLACE) {
373 recvcount = sendcounts[comm->rank()];
375 for (int i = 0; i < comm->size(); i++) {
376 if (sendcounts[i] < 0)
377 return MPI_ERR_COUNT;
383 int rank = simgrid::s4u::this_actor::get_pid();
384 int dt_size_send = sendtype->is_replayable() ? 1 : sendtype->size();
386 std::vector<int>* trace_sendcounts = new std::vector<int>;
387 if (comm->rank() == root) {
388 for (int i = 0; i < comm->size(); i++) { // copy data to avoid bad free
389 trace_sendcounts->push_back(sendcounts[i] * dt_size_send);
393 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Scatterv" : "PMPI_Iscatterv",
394 new simgrid::instr::VarCollTIData(
395 request == MPI_REQUEST_IGNORED ? "scatterv" : "iscatterv", root, dt_size_send,
396 trace_sendcounts, recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(),
397 nullptr, simgrid::smpi::Datatype::encode(sendtype),
398 simgrid::smpi::Datatype::encode(recvtype)));
399 if (request == MPI_REQUEST_IGNORED)
400 simgrid::smpi::Colls::scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm);
402 simgrid::smpi::Colls::iscatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm,
405 TRACE_smpi_comm_out(rank);
410 int PMPI_Reduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
412 return PMPI_Ireduce(sendbuf, recvbuf, count, datatype, op, root, comm, MPI_REQUEST_IGNORED);
415 int PMPI_Ireduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm, MPI_Request* request)
417 if (comm == MPI_COMM_NULL)
419 if ((sendbuf == nullptr && count > 0) || ((comm->rank() == root) && recvbuf == nullptr))
420 return MPI_ERR_BUFFER;
421 if (datatype == MPI_DATATYPE_NULL || not datatype->is_valid())
423 if (op == MPI_OP_NULL)
425 if (request == nullptr)
427 if (root < 0 || root >= comm->size())
430 return MPI_ERR_COUNT;
433 int rank = simgrid::s4u::this_actor::get_pid();
435 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Reduce" : "PMPI_Ireduce",
436 new simgrid::instr::CollTIData(request == MPI_REQUEST_IGNORED ? "reduce" : "ireduce", root, 0,
437 datatype->is_replayable() ? count : count * datatype->size(), -1,
438 simgrid::smpi::Datatype::encode(datatype), ""));
439 if (request == MPI_REQUEST_IGNORED)
440 simgrid::smpi::Colls::reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
442 simgrid::smpi::Colls::ireduce(sendbuf, recvbuf, count, datatype, op, root, comm, request);
444 TRACE_smpi_comm_out(rank);
449 int PMPI_Reduce_local(const void* inbuf, void* inoutbuf, int count, MPI_Datatype datatype, MPI_Op op)
451 if (datatype == MPI_DATATYPE_NULL || not datatype->is_valid())
453 if (op == MPI_OP_NULL)
456 return MPI_ERR_COUNT;
459 op->apply(inbuf, inoutbuf, &count, datatype);
464 int PMPI_Allreduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
466 return PMPI_Iallreduce(sendbuf, recvbuf, count, datatype, op, comm, MPI_REQUEST_IGNORED);
469 int PMPI_Iallreduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request *request)
471 if (comm == MPI_COMM_NULL)
473 if ((sendbuf == nullptr && count > 0) || (recvbuf == nullptr))
474 return MPI_ERR_BUFFER;
475 if (datatype == MPI_DATATYPE_NULL || not datatype->is_valid())
478 return MPI_ERR_COUNT;
479 if (op == MPI_OP_NULL)
481 if (request == nullptr)
485 const void* real_sendbuf = sendbuf;
486 std::unique_ptr<unsigned char[]> tmp_sendbuf;
487 if (sendbuf == MPI_IN_PLACE) {
488 tmp_sendbuf.reset(new unsigned char[count * datatype->get_extent()]);
489 simgrid::smpi::Datatype::copy(recvbuf, count, datatype, tmp_sendbuf.get(), count, datatype);
490 real_sendbuf = tmp_sendbuf.get();
492 int rank = simgrid::s4u::this_actor::get_pid();
494 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Allreduce" : "PMPI_Iallreduce",
495 new simgrid::instr::CollTIData(request == MPI_REQUEST_IGNORED ? "allreduce" : "iallreduce", -1, 0,
496 datatype->is_replayable() ? count : count * datatype->size(), -1,
497 simgrid::smpi::Datatype::encode(datatype), ""));
499 if (request == MPI_REQUEST_IGNORED)
500 simgrid::smpi::Colls::allreduce(real_sendbuf, recvbuf, count, datatype, op, comm);
502 simgrid::smpi::Colls::iallreduce(real_sendbuf, recvbuf, count, datatype, op, comm, request);
504 TRACE_smpi_comm_out(rank);
509 int PMPI_Scan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
511 return PMPI_Iscan(sendbuf, recvbuf, count, datatype, op, comm, MPI_REQUEST_IGNORED);
514 int PMPI_Iscan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request* request)
516 if (comm == MPI_COMM_NULL)
518 if (datatype == MPI_DATATYPE_NULL || not datatype->is_valid())
520 if (op == MPI_OP_NULL)
522 if (request == nullptr)
525 return MPI_ERR_COUNT;
526 if (sendbuf == nullptr || recvbuf == nullptr)
527 return MPI_ERR_BUFFER;
530 int rank = simgrid::s4u::this_actor::get_pid();
531 const void* real_sendbuf = sendbuf;
532 std::unique_ptr<unsigned char[]> tmp_sendbuf;
533 if (sendbuf == MPI_IN_PLACE) {
534 tmp_sendbuf.reset(new unsigned char[count * datatype->size()]);
535 real_sendbuf = memcpy(tmp_sendbuf.get(), recvbuf, count * datatype->size());
537 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Scan" : "PMPI_Iscan",
538 new simgrid::instr::Pt2PtTIData(request == MPI_REQUEST_IGNORED ? "scan" : "iscan", -1,
539 datatype->is_replayable() ? count : count * datatype->size(),
540 simgrid::smpi::Datatype::encode(datatype)));
543 if (request == MPI_REQUEST_IGNORED)
544 retval = simgrid::smpi::Colls::scan(real_sendbuf, recvbuf, count, datatype, op, comm);
546 retval = simgrid::smpi::Colls::iscan(real_sendbuf, recvbuf, count, datatype, op, comm, request);
548 TRACE_smpi_comm_out(rank);
553 int PMPI_Exscan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
555 return PMPI_Iexscan(sendbuf, recvbuf, count, datatype, op, comm, MPI_REQUEST_IGNORED);
558 int PMPI_Iexscan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request* request){
559 if (comm == MPI_COMM_NULL)
561 if (not datatype->is_valid())
563 if (op == MPI_OP_NULL)
565 if (request == nullptr)
568 return MPI_ERR_COUNT;
569 if (sendbuf == nullptr || recvbuf == nullptr)
570 return MPI_ERR_BUFFER;
573 int rank = simgrid::s4u::this_actor::get_pid();
574 const void* real_sendbuf = sendbuf;
575 std::unique_ptr<unsigned char[]> tmp_sendbuf;
576 if (sendbuf == MPI_IN_PLACE) {
577 tmp_sendbuf.reset(new unsigned char[count * datatype->size()]);
578 real_sendbuf = memcpy(tmp_sendbuf.get(), recvbuf, count * datatype->size());
581 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Exscan" : "PMPI_Iexscan",
582 new simgrid::instr::Pt2PtTIData(request == MPI_REQUEST_IGNORED ? "exscan" : "iexscan", -1,
583 datatype->is_replayable() ? count : count * datatype->size(),
584 simgrid::smpi::Datatype::encode(datatype)));
587 if (request == MPI_REQUEST_IGNORED)
588 retval = simgrid::smpi::Colls::exscan(real_sendbuf, recvbuf, count, datatype, op, comm);
590 retval = simgrid::smpi::Colls::iexscan(real_sendbuf, recvbuf, count, datatype, op, comm, request);
592 TRACE_smpi_comm_out(rank);
597 int PMPI_Reduce_scatter(const void *sendbuf, void *recvbuf, const int *recvcounts, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
599 return PMPI_Ireduce_scatter(sendbuf, recvbuf, recvcounts, datatype, op, comm, MPI_REQUEST_IGNORED);
602 int PMPI_Ireduce_scatter(const void *sendbuf, void *recvbuf, const int *recvcounts, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request *request)
604 if (comm == MPI_COMM_NULL)
606 if ((sendbuf == nullptr) || (recvbuf == nullptr))
607 return MPI_ERR_BUFFER;
608 if (datatype == MPI_DATATYPE_NULL || not datatype->is_valid())
610 if (op == MPI_OP_NULL)
612 if (recvcounts == nullptr)
614 if (request == nullptr)
617 for (int i = 0; i < comm->size(); i++) { // copy data to avoid bad free
618 if (recvcounts[i] < 0)
619 return MPI_ERR_COUNT;
623 int rank = simgrid::s4u::this_actor::get_pid();
624 std::vector<int>* trace_recvcounts = new std::vector<int>;
625 int dt_send_size = datatype->is_replayable() ? 1 : datatype->size();
628 for (int i = 0; i < comm->size(); i++) { // copy data to avoid bad free
629 trace_recvcounts->push_back(recvcounts[i] * dt_send_size);
630 totalcount += recvcounts[i];
633 const void* real_sendbuf = sendbuf;
634 std::unique_ptr<unsigned char[]> tmp_sendbuf;
635 if (sendbuf == MPI_IN_PLACE) {
636 tmp_sendbuf.reset(new unsigned char[totalcount * datatype->size()]);
637 real_sendbuf = memcpy(tmp_sendbuf.get(), recvbuf, totalcount * datatype->size());
640 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Reduce_scatter" : "PMPI_Ireduce_scatter",
641 new simgrid::instr::VarCollTIData(
642 request == MPI_REQUEST_IGNORED ? "reducescatter" : "ireducescatter", -1, dt_send_size, nullptr,
643 -1, trace_recvcounts, simgrid::smpi::Datatype::encode(datatype), ""));
645 if (request == MPI_REQUEST_IGNORED)
646 simgrid::smpi::Colls::reduce_scatter(real_sendbuf, recvbuf, recvcounts, datatype, op, comm);
648 simgrid::smpi::Colls::ireduce_scatter(real_sendbuf, recvbuf, recvcounts, datatype, op, comm, request);
650 TRACE_smpi_comm_out(rank);
655 int PMPI_Reduce_scatter_block(const void *sendbuf, void *recvbuf, int recvcount,
656 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
658 return PMPI_Ireduce_scatter_block(sendbuf, recvbuf, recvcount, datatype, op, comm, MPI_REQUEST_IGNORED);
661 int PMPI_Ireduce_scatter_block(const void* sendbuf, void* recvbuf, int recvcount, MPI_Datatype datatype, MPI_Op op,
662 MPI_Comm comm, MPI_Request* request)
664 if (comm == MPI_COMM_NULL)
666 if (not datatype->is_valid())
668 if (op == MPI_OP_NULL)
672 if (request == nullptr)
676 int count = comm->size();
678 int rank = simgrid::s4u::this_actor::get_pid();
679 int dt_send_size = datatype->is_replayable() ? 1 : datatype->size();
680 std::vector<int>* trace_recvcounts = new std::vector<int>(recvcount * dt_send_size); // copy data to avoid bad free
682 const void* real_sendbuf = sendbuf;
683 std::unique_ptr<unsigned char[]> tmp_sendbuf;
684 if (sendbuf == MPI_IN_PLACE) {
685 tmp_sendbuf.reset(new unsigned char[recvcount * count * datatype->size()]);
686 real_sendbuf = memcpy(tmp_sendbuf.get(), recvbuf, recvcount * count * datatype->size());
690 rank, request == MPI_REQUEST_IGNORED ? "PMPI_Reduce_scatter_block" : "PMPI_Ireduce_scatter_block",
691 new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "reducescatter" : "ireducescatter", -1, 0,
692 nullptr, -1, trace_recvcounts, simgrid::smpi::Datatype::encode(datatype), ""));
694 int* recvcounts = new int[count];
695 for (int i = 0; i < count; i++)
696 recvcounts[i] = recvcount;
697 if (request == MPI_REQUEST_IGNORED)
698 simgrid::smpi::Colls::reduce_scatter(real_sendbuf, recvbuf, recvcounts, datatype, op, comm);
700 simgrid::smpi::Colls::ireduce_scatter(real_sendbuf, recvbuf, recvcounts, datatype, op, comm, request);
703 TRACE_smpi_comm_out(rank);
708 int PMPI_Alltoall(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
709 MPI_Datatype recvtype, MPI_Comm comm){
710 return PMPI_Ialltoall(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm, MPI_REQUEST_IGNORED);
713 int PMPI_Ialltoall(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
714 MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request)
716 if (comm == MPI_COMM_NULL)
718 if ((sendbuf == nullptr && sendcount > 0) || (recvbuf == nullptr && recvcount > 0))
719 return MPI_ERR_BUFFER;
720 if ((sendbuf != MPI_IN_PLACE && sendtype == MPI_DATATYPE_NULL) || recvtype == MPI_DATATYPE_NULL)
722 if ((sendbuf != MPI_IN_PLACE && sendcount < 0) || recvcount < 0)
723 return MPI_ERR_COUNT;
724 if (request == nullptr)
728 int rank = simgrid::s4u::this_actor::get_pid();
729 const void* real_sendbuf = sendbuf;
730 int real_sendcount = sendcount;
731 MPI_Datatype real_sendtype = sendtype;
732 std::unique_ptr<unsigned char[]> tmp_sendbuf;
733 if (sendbuf == MPI_IN_PLACE) {
734 tmp_sendbuf.reset(new unsigned char[recvcount * comm->size() * recvtype->size()]);
735 // memcpy(??,nullptr,0) is actually undefined behavor, even if harmless.
736 if (recvbuf != nullptr)
737 memcpy(tmp_sendbuf.get(), recvbuf, recvcount * comm->size() * recvtype->size());
738 real_sendbuf = tmp_sendbuf.get();
739 real_sendcount = recvcount;
740 real_sendtype = recvtype;
743 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Alltoall" : "PMPI_Ialltoall",
744 new simgrid::instr::CollTIData(
745 request == MPI_REQUEST_IGNORED ? "alltoall" : "ialltoall", -1, -1.0,
746 real_sendtype->is_replayable() ? real_sendcount : real_sendcount * real_sendtype->size(),
747 recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(),
748 simgrid::smpi::Datatype::encode(real_sendtype), simgrid::smpi::Datatype::encode(recvtype)));
750 if (request == MPI_REQUEST_IGNORED)
752 simgrid::smpi::Colls::alltoall(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcount, recvtype, comm);
754 retval = simgrid::smpi::Colls::ialltoall(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcount, recvtype,
757 TRACE_smpi_comm_out(rank);
762 int PMPI_Alltoallv(const void* sendbuf, const int* sendcounts, const int* senddisps, MPI_Datatype sendtype, void* recvbuf,
763 const int* recvcounts, const int* recvdisps, MPI_Datatype recvtype, MPI_Comm comm)
765 return PMPI_Ialltoallv(sendbuf, sendcounts, senddisps, sendtype, recvbuf, recvcounts, recvdisps, recvtype, comm, MPI_REQUEST_IGNORED);
768 int PMPI_Ialltoallv(const void* sendbuf, const int* sendcounts, const int* senddisps, MPI_Datatype sendtype, void* recvbuf,
769 const int* recvcounts, const int* recvdisps, MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request)
771 if (comm == MPI_COMM_NULL)
773 if (sendbuf == nullptr || recvbuf == nullptr)
774 return MPI_ERR_BUFFER;
775 if ((sendbuf != MPI_IN_PLACE && sendtype == MPI_DATATYPE_NULL) || recvtype == MPI_DATATYPE_NULL)
777 if ((sendbuf != MPI_IN_PLACE && (sendcounts == nullptr || senddisps == nullptr)) || recvcounts == nullptr ||
778 recvdisps == nullptr)
780 if (request == nullptr)
783 int rank = simgrid::s4u::this_actor::get_pid();
784 int size = comm->size();
785 for (int i = 0; i < size; i++) {
786 if (recvcounts[i] < 0 || (sendbuf != MPI_IN_PLACE && sendcounts[i] < 0))
787 return MPI_ERR_COUNT;
793 std::vector<int>* trace_sendcounts = new std::vector<int>;
794 std::vector<int>* trace_recvcounts = new std::vector<int>;
795 int dt_size_recv = recvtype->size();
797 const void* real_sendbuf = sendbuf;
798 const int* real_sendcounts = sendcounts;
799 const int* real_senddisps = senddisps;
800 MPI_Datatype real_sendtype = sendtype;
802 for (int i = 0; i < size; i++) { // copy data to avoid bad free
803 recv_size += recvcounts[i] * dt_size_recv;
804 trace_recvcounts->push_back(recvcounts[i] * dt_size_recv);
805 if (((recvdisps[i] + recvcounts[i]) * dt_size_recv) > maxsize)
806 maxsize = (recvdisps[i] + recvcounts[i]) * dt_size_recv;
809 std::unique_ptr<unsigned char[]> tmp_sendbuf;
810 std::unique_ptr<int[]> tmp_sendcounts;
811 std::unique_ptr<int[]> tmp_senddisps;
812 if (sendbuf == MPI_IN_PLACE) {
813 tmp_sendbuf.reset(new unsigned char[maxsize]);
814 real_sendbuf = memcpy(tmp_sendbuf.get(), recvbuf, maxsize);
815 tmp_sendcounts.reset(new int[size]);
816 std::copy(recvcounts, recvcounts + size, tmp_sendcounts.get());
817 real_sendcounts = tmp_sendcounts.get();
818 tmp_senddisps.reset(new int[size]);
819 std::copy(recvdisps, recvdisps + size, tmp_senddisps.get());
820 real_senddisps = tmp_senddisps.get();
821 real_sendtype = recvtype;
824 int dt_size_send = real_sendtype->size();
826 for (int i = 0; i < size; i++) { // copy data to avoid bad free
827 send_size += real_sendcounts[i] * dt_size_send;
828 trace_sendcounts->push_back(real_sendcounts[i] * dt_size_send);
831 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Alltoallv" : "PMPI_Ialltoallv",
832 new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "alltoallv" : "ialltoallv", -1,
833 send_size, trace_sendcounts, recv_size, trace_recvcounts,
834 simgrid::smpi::Datatype::encode(real_sendtype),
835 simgrid::smpi::Datatype::encode(recvtype)));
838 if (request == MPI_REQUEST_IGNORED)
839 retval = simgrid::smpi::Colls::alltoallv(real_sendbuf, real_sendcounts, real_senddisps, real_sendtype, recvbuf,
840 recvcounts, recvdisps, recvtype, comm);
842 retval = simgrid::smpi::Colls::ialltoallv(real_sendbuf, real_sendcounts, real_senddisps, real_sendtype, recvbuf,
843 recvcounts, recvdisps, recvtype, comm, request);
845 TRACE_smpi_comm_out(rank);
850 int PMPI_Alltoallw(const void* sendbuf, const int* sendcounts, const int* senddisps, const MPI_Datatype* sendtypes, void* recvbuf,
851 const int* recvcounts, const int* recvdisps, const MPI_Datatype* recvtypes, MPI_Comm comm)
853 return PMPI_Ialltoallw(sendbuf, sendcounts, senddisps, sendtypes, recvbuf, recvcounts, recvdisps, recvtypes, comm, MPI_REQUEST_IGNORED);
856 int PMPI_Ialltoallw(const void* sendbuf, const int* sendcounts, const int* senddisps, const MPI_Datatype* sendtypes, void* recvbuf,
857 const int* recvcounts, const int* recvdisps, const MPI_Datatype* recvtypes, MPI_Comm comm, MPI_Request* request)
859 if (comm == MPI_COMM_NULL)
861 if (sendbuf == nullptr || recvbuf == nullptr)
862 return MPI_ERR_BUFFER;
863 if ((sendbuf != MPI_IN_PLACE && sendtypes == nullptr) || recvtypes == nullptr)
865 if ((sendbuf != MPI_IN_PLACE && (sendcounts == nullptr || senddisps == nullptr)) || recvcounts == nullptr ||
866 recvdisps == nullptr)
868 if (request == nullptr)
872 int rank = simgrid::s4u::this_actor::get_pid();
873 int size = comm->size();
874 for (int i = 0; i < size; i++) { // copy data to avoid bad free
875 if (recvcounts[i] < 0 || (sendbuf != MPI_IN_PLACE && sendcounts[i] < 0))
876 return MPI_ERR_COUNT;
880 std::vector<int>* trace_sendcounts = new std::vector<int>;
881 std::vector<int>* trace_recvcounts = new std::vector<int>;
883 const void* real_sendbuf = sendbuf;
884 const int* real_sendcounts = sendcounts;
885 const int* real_senddisps = senddisps;
886 const MPI_Datatype* real_sendtypes = sendtypes;
887 unsigned long maxsize = 0;
888 for (int i = 0; i < size; i++) { // copy data to avoid bad free
889 if (recvtypes[i] == MPI_DATATYPE_NULL) {
890 delete trace_recvcounts;
891 delete trace_sendcounts;
894 recv_size += recvcounts[i] * recvtypes[i]->size();
895 trace_recvcounts->push_back(recvcounts[i] * recvtypes[i]->size());
896 if ((recvdisps[i] + (recvcounts[i] * recvtypes[i]->size())) > maxsize)
897 maxsize = recvdisps[i] + (recvcounts[i] * recvtypes[i]->size());
900 std::unique_ptr<unsigned char[]> tmp_sendbuf;
901 std::unique_ptr<int[]> tmp_sendcounts;
902 std::unique_ptr<int[]> tmp_senddisps;
903 std::unique_ptr<MPI_Datatype[]> tmp_sendtypes;
904 if (sendbuf == MPI_IN_PLACE) {
905 tmp_sendbuf.reset(new unsigned char[maxsize]);
906 real_sendbuf = memcpy(tmp_sendbuf.get(), recvbuf, maxsize);
907 tmp_sendcounts.reset(new int[size]);
908 std::copy(recvcounts, recvcounts + size, tmp_sendcounts.get());
909 real_sendcounts = tmp_sendcounts.get();
910 tmp_senddisps.reset(new int[size]);
911 std::copy(recvdisps, recvdisps + size, tmp_senddisps.get());
912 real_senddisps = tmp_senddisps.get();
913 tmp_sendtypes.reset(new MPI_Datatype[size]);
914 std::copy(recvtypes, recvtypes + size, tmp_sendtypes.get());
915 real_sendtypes = tmp_sendtypes.get();
918 for (int i = 0; i < size; i++) { // copy data to avoid bad free
919 send_size += real_sendcounts[i] * real_sendtypes[i]->size();
920 trace_sendcounts->push_back(real_sendcounts[i] * real_sendtypes[i]->size());
923 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Alltoallw" : "PMPI_Ialltoallw",
924 new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "alltoallv" : "ialltoallv", -1,
925 send_size, trace_sendcounts, recv_size, trace_recvcounts,
926 simgrid::smpi::Datatype::encode(real_sendtypes[0]),
927 simgrid::smpi::Datatype::encode(recvtypes[0])));
930 if (request == MPI_REQUEST_IGNORED)
931 retval = simgrid::smpi::Colls::alltoallw(real_sendbuf, real_sendcounts, real_senddisps, real_sendtypes, recvbuf,
932 recvcounts, recvdisps, recvtypes, comm);
934 retval = simgrid::smpi::Colls::ialltoallw(real_sendbuf, real_sendcounts, real_senddisps, real_sendtypes, recvbuf,
935 recvcounts, recvdisps, recvtypes, comm, request);
937 TRACE_smpi_comm_out(rank);