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 #define CHECK_ARGS(test, errcode, ...) \
18 XBT_WARN(__VA_ARGS__); \
22 /* PMPI User level calls */
24 int PMPI_Barrier(MPI_Comm comm)
26 return PMPI_Ibarrier(comm, MPI_REQUEST_IGNORED);
29 int PMPI_Ibarrier(MPI_Comm comm, MPI_Request *request)
31 if (comm == MPI_COMM_NULL)
33 if (request == nullptr)
37 int rank = simgrid::s4u::this_actor::get_pid();
38 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Barrier" : "PMPI_Ibarrier",
39 new simgrid::instr::NoOpTIData(request == MPI_REQUEST_IGNORED ? "barrier" : "ibarrier"));
40 if (request == MPI_REQUEST_IGNORED) {
41 simgrid::smpi::Colls::barrier(comm);
42 // Barrier can be used to synchronize RMA calls. Finish all requests from comm before.
43 comm->finish_rma_calls();
45 simgrid::smpi::Colls::ibarrier(comm, request);
47 TRACE_smpi_comm_out(rank);
52 int PMPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
54 return PMPI_Ibcast(buf, count, datatype, root, comm, MPI_REQUEST_IGNORED);
57 int PMPI_Ibcast(void *buf, int count, MPI_Datatype datatype,
58 int root, MPI_Comm comm, MPI_Request* request)
60 if (comm == MPI_COMM_NULL)
62 if (buf == nullptr && count > 0)
63 return MPI_ERR_BUFFER;
64 if (datatype == MPI_DATATYPE_NULL || not datatype->is_valid())
68 if (root < 0 || root >= comm->size())
70 if (request == nullptr)
74 int rank = simgrid::s4u::this_actor::get_pid();
75 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Bcast" : "PMPI_Ibcast",
76 new simgrid::instr::CollTIData(request == MPI_REQUEST_IGNORED ? "bcast" : "ibcast", root, -1.0,
77 datatype->is_replayable() ? count : count * datatype->size(), -1,
78 simgrid::smpi::Datatype::encode(datatype), ""));
79 if (comm->size() > 1) {
80 if (request == MPI_REQUEST_IGNORED)
81 simgrid::smpi::Colls::bcast(buf, count, datatype, root, comm);
83 simgrid::smpi::Colls::ibcast(buf, count, datatype, root, comm, request);
85 if (request != MPI_REQUEST_IGNORED)
86 *request = MPI_REQUEST_NULL;
89 TRACE_smpi_comm_out(rank);
94 int PMPI_Gather(const void *sendbuf, int sendcount, MPI_Datatype sendtype,void *recvbuf, int recvcount, MPI_Datatype recvtype,
95 int root, MPI_Comm comm){
96 return PMPI_Igather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm, MPI_REQUEST_IGNORED);
99 int PMPI_Igather(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
100 MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request* request)
102 if (comm == MPI_COMM_NULL)
104 if ((sendbuf == nullptr) || ((comm->rank() == root) && recvbuf == nullptr))
105 return MPI_ERR_BUFFER;
106 if (((sendbuf != MPI_IN_PLACE && sendcount > 0) && (sendtype == MPI_DATATYPE_NULL)) ||
107 ((comm->rank() == root) && (recvtype == MPI_DATATYPE_NULL)))
109 if (((sendbuf != MPI_IN_PLACE) && (sendcount < 0)) || ((comm->rank() == root) && (recvcount < 0)))
110 return MPI_ERR_COUNT;
111 if (root < 0 || root >= comm->size())
113 if (request == nullptr)
117 const void* real_sendbuf = sendbuf;
118 int real_sendcount = sendcount;
119 MPI_Datatype real_sendtype = sendtype;
120 if ((comm->rank() == root) && (sendbuf == MPI_IN_PLACE)) {
122 real_sendtype = recvtype;
124 int rank = simgrid::s4u::this_actor::get_pid();
126 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Gather" : "PMPI_Igather",
127 new simgrid::instr::CollTIData(
128 request == MPI_REQUEST_IGNORED ? "gather" : "igather", root, -1.0,
129 real_sendtype->is_replayable() ? real_sendcount : real_sendcount * real_sendtype->size(),
130 (comm->rank() != root || recvtype->is_replayable()) ? recvcount : recvcount * recvtype->size(),
131 simgrid::smpi::Datatype::encode(real_sendtype), simgrid::smpi::Datatype::encode(recvtype)));
132 if (request == MPI_REQUEST_IGNORED)
133 simgrid::smpi::Colls::gather(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcount, recvtype, root, comm);
135 simgrid::smpi::Colls::igather(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcount, recvtype, root, comm,
138 TRACE_smpi_comm_out(rank);
143 int PMPI_Gatherv(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, const int *recvcounts, const int *displs,
144 MPI_Datatype recvtype, int root, MPI_Comm comm){
145 return PMPI_Igatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, root, comm, MPI_REQUEST_IGNORED);
148 int PMPI_Igatherv(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, const int* recvcounts, const int* displs,
149 MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request* request)
151 if (comm == MPI_COMM_NULL)
153 if ((sendbuf == nullptr && sendcount > 0) || ((comm->rank() == root) && recvbuf == nullptr))
154 return MPI_ERR_BUFFER;
155 if (((sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
156 ((comm->rank() == root) && (recvtype == MPI_DATATYPE_NULL)))
158 if ((sendbuf != MPI_IN_PLACE) && (sendcount < 0))
159 return MPI_ERR_COUNT;
160 if ((comm->rank() == root) && (recvcounts == nullptr || displs == nullptr))
162 if (root < 0 || root >= comm->size())
164 if (request == nullptr)
167 for (int i = 0; i < comm->size(); i++) {
168 if ((comm->rank() == root) && (recvcounts[i] < 0))
169 return MPI_ERR_COUNT;
173 const void* real_sendbuf = sendbuf;
174 int real_sendcount = sendcount;
175 MPI_Datatype real_sendtype = sendtype;
176 if ((comm->rank() == root) && (sendbuf == MPI_IN_PLACE)) {
178 real_sendtype = recvtype;
181 int rank = simgrid::s4u::this_actor::get_pid();
182 int dt_size_recv = recvtype->is_replayable() ? 1 : recvtype->size();
184 std::vector<int>* trace_recvcounts = new std::vector<int>;
185 if (comm->rank() == root) {
186 for (int i = 0; i < comm->size(); i++) // copy data to avoid bad free
187 trace_recvcounts->push_back(recvcounts[i] * dt_size_recv);
190 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Gatherv" : "PMPI_Igatherv",
191 new simgrid::instr::VarCollTIData(
192 request == MPI_REQUEST_IGNORED ? "gatherv" : "igatherv", root,
193 real_sendtype->is_replayable() ? real_sendcount : real_sendcount * real_sendtype->size(),
194 nullptr, dt_size_recv, trace_recvcounts, simgrid::smpi::Datatype::encode(real_sendtype),
195 simgrid::smpi::Datatype::encode(recvtype)));
196 if (request == MPI_REQUEST_IGNORED)
197 simgrid::smpi::Colls::gatherv(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcounts, displs, recvtype,
200 simgrid::smpi::Colls::igatherv(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcounts, displs, recvtype,
201 root, comm, request);
203 TRACE_smpi_comm_out(rank);
208 int PMPI_Allgather(const void *sendbuf, int sendcount, MPI_Datatype sendtype,
209 void *recvbuf, int recvcount, MPI_Datatype recvtype, MPI_Comm comm){
210 return PMPI_Iallgather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm, MPI_REQUEST_IGNORED);
213 int PMPI_Iallgather(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
214 MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request)
216 if (comm == MPI_COMM_NULL)
218 if ((sendbuf == nullptr && sendcount > 0) || (recvbuf == nullptr))
219 return MPI_ERR_BUFFER;
220 if (((sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) || (recvtype == MPI_DATATYPE_NULL))
222 if (((sendbuf != MPI_IN_PLACE) && (sendcount < 0)) || (recvcount < 0))
223 return MPI_ERR_COUNT;
224 if (request == nullptr)
228 if (sendbuf == MPI_IN_PLACE) {
229 sendbuf = static_cast<char*>(recvbuf) + recvtype->get_extent() * recvcount * comm->rank();
230 sendcount = recvcount;
233 int rank = simgrid::s4u::this_actor::get_pid();
235 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Allgather" : "PMPI_Iallggather",
236 new simgrid::instr::CollTIData(
237 request == MPI_REQUEST_IGNORED ? "allgather" : "iallgather", -1, -1.0,
238 sendtype->is_replayable() ? sendcount : sendcount * sendtype->size(),
239 recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(),
240 simgrid::smpi::Datatype::encode(sendtype), simgrid::smpi::Datatype::encode(recvtype)));
241 if (request == MPI_REQUEST_IGNORED)
242 simgrid::smpi::Colls::allgather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
244 simgrid::smpi::Colls::iallgather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm, request);
246 TRACE_smpi_comm_out(rank);
251 int PMPI_Allgatherv(const void *sendbuf, int sendcount, MPI_Datatype sendtype,
252 void *recvbuf, const int *recvcounts, const int *displs, MPI_Datatype recvtype, MPI_Comm comm){
253 return PMPI_Iallgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm, MPI_REQUEST_IGNORED);
256 int PMPI_Iallgatherv(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, const int* recvcounts, const int* displs,
257 MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request)
259 if (comm == MPI_COMM_NULL)
261 if ((sendbuf == nullptr && sendcount > 0) || (recvbuf == nullptr))
262 return MPI_ERR_BUFFER;
263 if (((sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) || (recvtype == MPI_DATATYPE_NULL))
265 if ((sendbuf != MPI_IN_PLACE) && (sendcount < 0))
266 return MPI_ERR_COUNT;
267 if (recvcounts == nullptr || displs == nullptr)
269 if (request == nullptr)
272 for (int i = 0; i < comm->size(); i++) { // copy data to avoid bad free
273 if (recvcounts[i] < 0)
274 return MPI_ERR_COUNT;
278 if (sendbuf == MPI_IN_PLACE) {
279 sendbuf = static_cast<char*>(recvbuf) + recvtype->get_extent() * displs[comm->rank()];
280 sendcount = recvcounts[comm->rank()];
283 int rank = simgrid::s4u::this_actor::get_pid();
284 int dt_size_recv = recvtype->is_replayable() ? 1 : recvtype->size();
286 std::vector<int>* trace_recvcounts = new std::vector<int>;
287 for (int i = 0; i < comm->size(); i++) { // copy data to avoid bad free
288 trace_recvcounts->push_back(recvcounts[i] * dt_size_recv);
292 rank, request == MPI_REQUEST_IGNORED ? "PMPI_Allgatherv" : "PMPI_Iallgatherv",
293 new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "allgatherv" : "iallgatherv", -1,
294 sendtype->is_replayable() ? sendcount : sendcount * sendtype->size(), nullptr,
295 dt_size_recv, trace_recvcounts, simgrid::smpi::Datatype::encode(sendtype),
296 simgrid::smpi::Datatype::encode(recvtype)));
297 if (request == MPI_REQUEST_IGNORED)
298 simgrid::smpi::Colls::allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm);
300 simgrid::smpi::Colls::iallgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm,
303 TRACE_smpi_comm_out(rank);
308 int PMPI_Scatter(const void *sendbuf, int sendcount, MPI_Datatype sendtype,
309 void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm){
310 return PMPI_Iscatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm, MPI_REQUEST_IGNORED);
313 int PMPI_Iscatter(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
314 MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request* request)
316 if (comm == MPI_COMM_NULL)
318 if (((comm->rank() == root) && (sendtype == MPI_DATATYPE_NULL || not sendtype->is_valid())) ||
319 ((recvbuf != MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL || not recvtype->is_valid())))
321 if (((comm->rank() == root) && (sendcount < 0)) || ((recvbuf != MPI_IN_PLACE) && (recvcount < 0)))
322 return MPI_ERR_COUNT;
323 if ((sendbuf == recvbuf) || ((comm->rank() == root) && sendcount > 0 && (sendbuf == nullptr)) ||
324 (recvcount > 0 && recvbuf == nullptr))
325 return MPI_ERR_BUFFER;
326 if (root < 0 || root >= comm->size())
328 if (request == nullptr)
332 if (recvbuf == MPI_IN_PLACE) {
334 recvcount = sendcount;
336 int rank = simgrid::s4u::this_actor::get_pid();
338 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Scatter" : "PMPI_Iscatter",
339 new simgrid::instr::CollTIData(
340 request == MPI_REQUEST_IGNORED ? "scatter" : "iscatter", root, -1.0,
341 (comm->rank() != root || sendtype->is_replayable()) ? sendcount : sendcount * sendtype->size(),
342 recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(),
343 simgrid::smpi::Datatype::encode(sendtype), simgrid::smpi::Datatype::encode(recvtype)));
344 if (request == MPI_REQUEST_IGNORED)
345 simgrid::smpi::Colls::scatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm);
347 simgrid::smpi::Colls::iscatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm, request);
349 TRACE_smpi_comm_out(rank);
354 int PMPI_Scatterv(const void *sendbuf, const int *sendcounts, const int *displs,
355 MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm){
356 return PMPI_Iscatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm, MPI_REQUEST_IGNORED);
359 int PMPI_Iscatterv(const void* sendbuf, const int* sendcounts, const int* displs, MPI_Datatype sendtype, void* recvbuf, int recvcount,
360 MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request* request)
362 CHECK_ARGS(comm == MPI_COMM_NULL, MPI_ERR_COMM, "Iscatterv: the communicator cannot be MPI_COMM_NULL");
363 CHECK_ARGS((comm->rank() == root) && (sendcounts == nullptr), MPI_ERR_ARG,
364 "Iscatterv: param 2 sendcounts cannot be NULL on the root rank");
365 CHECK_ARGS((comm->rank() == root) && (displs == nullptr), MPI_ERR_ARG,
366 "Iscatterv: param 3 displs cannot be NULL on the root rank");
367 CHECK_ARGS((comm->rank() == root) && (sendtype == MPI_DATATYPE_NULL), MPI_ERR_TYPE,
368 "Iscatterv: The sendtype cannot be NULL on the root rank");
369 CHECK_ARGS((recvbuf != MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL), MPI_ERR_TYPE,
370 "Iscatterv: the recvtype cannot be NULL when not receiving in place");
371 CHECK_ARGS(request == nullptr, MPI_ERR_ARG, "Iscatterv: param 10 request cannot be NULL");
372 CHECK_ARGS(recvbuf != MPI_IN_PLACE && recvcount < 0, MPI_ERR_COUNT,
373 "Iscatterv: When not receiving in place, the recvcound cannot be negative");
374 CHECK_ARGS(root < 0, MPI_ERR_ROOT, "Iscatterv: root cannot be negative");
375 CHECK_ARGS(root >= comm->size(), MPI_ERR_ROOT, "Iscatterv: root (=%d) is larger than communicator size (=%d)", root,
378 if (comm->rank() == root) {
379 if (recvbuf == MPI_IN_PLACE) {
381 recvcount = sendcounts[comm->rank()];
383 for (int i = 0; i < comm->size(); i++)
384 CHECK_ARGS(sendcounts[i] < 0, MPI_ERR_COUNT, "Iscatterv: sendcounts[%d]=%d but this cannot be negative", i,
390 int rank = simgrid::s4u::this_actor::get_pid();
391 int dt_size_send = sendtype->is_replayable() ? 1 : sendtype->size();
393 std::vector<int>* trace_sendcounts = new std::vector<int>;
394 if (comm->rank() == root) {
395 for (int i = 0; i < comm->size(); i++) { // copy data to avoid bad free
396 trace_sendcounts->push_back(sendcounts[i] * dt_size_send);
400 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Scatterv" : "PMPI_Iscatterv",
401 new simgrid::instr::VarCollTIData(
402 request == MPI_REQUEST_IGNORED ? "scatterv" : "iscatterv", root, dt_size_send,
403 trace_sendcounts, recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(),
404 nullptr, simgrid::smpi::Datatype::encode(sendtype),
405 simgrid::smpi::Datatype::encode(recvtype)));
406 if (request == MPI_REQUEST_IGNORED)
407 simgrid::smpi::Colls::scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm);
409 simgrid::smpi::Colls::iscatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm,
412 TRACE_smpi_comm_out(rank);
417 int PMPI_Reduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
419 return PMPI_Ireduce(sendbuf, recvbuf, count, datatype, op, root, comm, MPI_REQUEST_IGNORED);
422 int PMPI_Ireduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm, MPI_Request* request)
424 if (comm == MPI_COMM_NULL)
426 if ((sendbuf == nullptr && count > 0) || ((comm->rank() == root) && recvbuf == nullptr))
427 return MPI_ERR_BUFFER;
428 if (datatype == MPI_DATATYPE_NULL || not datatype->is_valid())
430 if (op == MPI_OP_NULL)
432 if (request == nullptr)
434 if (root < 0 || root >= comm->size())
437 return MPI_ERR_COUNT;
440 int rank = simgrid::s4u::this_actor::get_pid();
442 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Reduce" : "PMPI_Ireduce",
443 new simgrid::instr::CollTIData(request == MPI_REQUEST_IGNORED ? "reduce" : "ireduce", root, 0,
444 datatype->is_replayable() ? count : count * datatype->size(), -1,
445 simgrid::smpi::Datatype::encode(datatype), ""));
446 if (request == MPI_REQUEST_IGNORED)
447 simgrid::smpi::Colls::reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
449 simgrid::smpi::Colls::ireduce(sendbuf, recvbuf, count, datatype, op, root, comm, request);
451 TRACE_smpi_comm_out(rank);
456 int PMPI_Reduce_local(const void* inbuf, void* inoutbuf, int count, MPI_Datatype datatype, MPI_Op op)
458 if (datatype == MPI_DATATYPE_NULL || not datatype->is_valid())
460 if (op == MPI_OP_NULL)
463 return MPI_ERR_COUNT;
466 op->apply(inbuf, inoutbuf, &count, datatype);
471 int PMPI_Allreduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
473 return PMPI_Iallreduce(sendbuf, recvbuf, count, datatype, op, comm, MPI_REQUEST_IGNORED);
476 int PMPI_Iallreduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request *request)
478 if (comm == MPI_COMM_NULL)
480 if ((sendbuf == nullptr && count > 0) || (recvbuf == nullptr))
481 return MPI_ERR_BUFFER;
482 if (datatype == MPI_DATATYPE_NULL || not datatype->is_valid())
485 return MPI_ERR_COUNT;
486 if (op == MPI_OP_NULL)
488 if (request == nullptr)
492 const void* real_sendbuf = sendbuf;
493 std::unique_ptr<unsigned char[]> tmp_sendbuf;
494 if (sendbuf == MPI_IN_PLACE) {
495 tmp_sendbuf.reset(new unsigned char[count * datatype->get_extent()]);
496 simgrid::smpi::Datatype::copy(recvbuf, count, datatype, tmp_sendbuf.get(), count, datatype);
497 real_sendbuf = tmp_sendbuf.get();
499 int rank = simgrid::s4u::this_actor::get_pid();
501 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Allreduce" : "PMPI_Iallreduce",
502 new simgrid::instr::CollTIData(request == MPI_REQUEST_IGNORED ? "allreduce" : "iallreduce", -1, 0,
503 datatype->is_replayable() ? count : count * datatype->size(), -1,
504 simgrid::smpi::Datatype::encode(datatype), ""));
506 if (request == MPI_REQUEST_IGNORED)
507 simgrid::smpi::Colls::allreduce(real_sendbuf, recvbuf, count, datatype, op, comm);
509 simgrid::smpi::Colls::iallreduce(real_sendbuf, recvbuf, count, datatype, op, comm, request);
511 TRACE_smpi_comm_out(rank);
516 int PMPI_Scan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
518 return PMPI_Iscan(sendbuf, recvbuf, count, datatype, op, comm, MPI_REQUEST_IGNORED);
521 int PMPI_Iscan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request* request)
523 if (comm == MPI_COMM_NULL)
525 if (datatype == MPI_DATATYPE_NULL || not datatype->is_valid())
527 if (op == MPI_OP_NULL)
529 if (request == nullptr)
532 return MPI_ERR_COUNT;
533 if (sendbuf == nullptr || recvbuf == nullptr)
534 return MPI_ERR_BUFFER;
537 int rank = simgrid::s4u::this_actor::get_pid();
538 const void* real_sendbuf = sendbuf;
539 std::unique_ptr<unsigned char[]> tmp_sendbuf;
540 if (sendbuf == MPI_IN_PLACE) {
541 tmp_sendbuf.reset(new unsigned char[count * datatype->size()]);
542 real_sendbuf = memcpy(tmp_sendbuf.get(), recvbuf, count * datatype->size());
544 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Scan" : "PMPI_Iscan",
545 new simgrid::instr::Pt2PtTIData(request == MPI_REQUEST_IGNORED ? "scan" : "iscan", -1,
546 datatype->is_replayable() ? count : count * datatype->size(),
547 simgrid::smpi::Datatype::encode(datatype)));
550 if (request == MPI_REQUEST_IGNORED)
551 retval = simgrid::smpi::Colls::scan(real_sendbuf, recvbuf, count, datatype, op, comm);
553 retval = simgrid::smpi::Colls::iscan(real_sendbuf, recvbuf, count, datatype, op, comm, request);
555 TRACE_smpi_comm_out(rank);
560 int PMPI_Exscan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
562 return PMPI_Iexscan(sendbuf, recvbuf, count, datatype, op, comm, MPI_REQUEST_IGNORED);
565 int PMPI_Iexscan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request* request){
566 if (comm == MPI_COMM_NULL)
568 if (not datatype->is_valid())
570 if (op == MPI_OP_NULL)
572 if (request == nullptr)
575 return MPI_ERR_COUNT;
576 if (sendbuf == nullptr || recvbuf == nullptr)
577 return MPI_ERR_BUFFER;
580 int rank = simgrid::s4u::this_actor::get_pid();
581 const void* real_sendbuf = sendbuf;
582 std::unique_ptr<unsigned char[]> tmp_sendbuf;
583 if (sendbuf == MPI_IN_PLACE) {
584 tmp_sendbuf.reset(new unsigned char[count * datatype->size()]);
585 real_sendbuf = memcpy(tmp_sendbuf.get(), recvbuf, count * datatype->size());
588 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Exscan" : "PMPI_Iexscan",
589 new simgrid::instr::Pt2PtTIData(request == MPI_REQUEST_IGNORED ? "exscan" : "iexscan", -1,
590 datatype->is_replayable() ? count : count * datatype->size(),
591 simgrid::smpi::Datatype::encode(datatype)));
594 if (request == MPI_REQUEST_IGNORED)
595 retval = simgrid::smpi::Colls::exscan(real_sendbuf, recvbuf, count, datatype, op, comm);
597 retval = simgrid::smpi::Colls::iexscan(real_sendbuf, recvbuf, count, datatype, op, comm, request);
599 TRACE_smpi_comm_out(rank);
604 int PMPI_Reduce_scatter(const void *sendbuf, void *recvbuf, const int *recvcounts, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
606 return PMPI_Ireduce_scatter(sendbuf, recvbuf, recvcounts, datatype, op, comm, MPI_REQUEST_IGNORED);
609 int PMPI_Ireduce_scatter(const void *sendbuf, void *recvbuf, const int *recvcounts, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request *request)
611 if (comm == MPI_COMM_NULL)
613 if ((sendbuf == nullptr) || (recvbuf == nullptr))
614 return MPI_ERR_BUFFER;
615 if (datatype == MPI_DATATYPE_NULL || not datatype->is_valid())
617 if (op == MPI_OP_NULL)
619 if (recvcounts == nullptr)
621 if (request == nullptr)
624 for (int i = 0; i < comm->size(); i++) { // copy data to avoid bad free
625 if (recvcounts[i] < 0)
626 return MPI_ERR_COUNT;
630 int rank = simgrid::s4u::this_actor::get_pid();
631 std::vector<int>* trace_recvcounts = new std::vector<int>;
632 int dt_send_size = datatype->is_replayable() ? 1 : datatype->size();
635 for (int i = 0; i < comm->size(); i++) { // copy data to avoid bad free
636 trace_recvcounts->push_back(recvcounts[i] * dt_send_size);
637 totalcount += recvcounts[i];
640 const void* real_sendbuf = sendbuf;
641 std::unique_ptr<unsigned char[]> tmp_sendbuf;
642 if (sendbuf == MPI_IN_PLACE) {
643 tmp_sendbuf.reset(new unsigned char[totalcount * datatype->size()]);
644 real_sendbuf = memcpy(tmp_sendbuf.get(), recvbuf, totalcount * datatype->size());
647 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Reduce_scatter" : "PMPI_Ireduce_scatter",
648 new simgrid::instr::VarCollTIData(
649 request == MPI_REQUEST_IGNORED ? "reducescatter" : "ireducescatter", -1, dt_send_size, nullptr,
650 -1, trace_recvcounts, simgrid::smpi::Datatype::encode(datatype), ""));
652 if (request == MPI_REQUEST_IGNORED)
653 simgrid::smpi::Colls::reduce_scatter(real_sendbuf, recvbuf, recvcounts, datatype, op, comm);
655 simgrid::smpi::Colls::ireduce_scatter(real_sendbuf, recvbuf, recvcounts, datatype, op, comm, request);
657 TRACE_smpi_comm_out(rank);
662 int PMPI_Reduce_scatter_block(const void *sendbuf, void *recvbuf, int recvcount,
663 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
665 return PMPI_Ireduce_scatter_block(sendbuf, recvbuf, recvcount, datatype, op, comm, MPI_REQUEST_IGNORED);
668 int PMPI_Ireduce_scatter_block(const void* sendbuf, void* recvbuf, int recvcount, MPI_Datatype datatype, MPI_Op op,
669 MPI_Comm comm, MPI_Request* request)
671 if (comm == MPI_COMM_NULL)
673 if (not datatype->is_valid())
675 if (op == MPI_OP_NULL)
679 if (request == nullptr)
683 int count = comm->size();
685 int rank = simgrid::s4u::this_actor::get_pid();
686 int dt_send_size = datatype->is_replayable() ? 1 : datatype->size();
687 std::vector<int>* trace_recvcounts = new std::vector<int>(recvcount * dt_send_size); // copy data to avoid bad free
689 const void* real_sendbuf = sendbuf;
690 std::unique_ptr<unsigned char[]> tmp_sendbuf;
691 if (sendbuf == MPI_IN_PLACE) {
692 tmp_sendbuf.reset(new unsigned char[recvcount * count * datatype->size()]);
693 real_sendbuf = memcpy(tmp_sendbuf.get(), recvbuf, recvcount * count * datatype->size());
697 rank, request == MPI_REQUEST_IGNORED ? "PMPI_Reduce_scatter_block" : "PMPI_Ireduce_scatter_block",
698 new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "reducescatter" : "ireducescatter", -1, 0,
699 nullptr, -1, trace_recvcounts, simgrid::smpi::Datatype::encode(datatype), ""));
701 int* recvcounts = new int[count];
702 for (int i = 0; i < count; i++)
703 recvcounts[i] = recvcount;
704 if (request == MPI_REQUEST_IGNORED)
705 simgrid::smpi::Colls::reduce_scatter(real_sendbuf, recvbuf, recvcounts, datatype, op, comm);
707 simgrid::smpi::Colls::ireduce_scatter(real_sendbuf, recvbuf, recvcounts, datatype, op, comm, request);
710 TRACE_smpi_comm_out(rank);
715 int PMPI_Alltoall(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
716 MPI_Datatype recvtype, MPI_Comm comm){
717 return PMPI_Ialltoall(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm, MPI_REQUEST_IGNORED);
720 int PMPI_Ialltoall(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
721 MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request)
723 if (comm == MPI_COMM_NULL)
725 if ((sendbuf == nullptr && sendcount > 0) || (recvbuf == nullptr && recvcount > 0))
726 return MPI_ERR_BUFFER;
727 if ((sendbuf != MPI_IN_PLACE && sendtype == MPI_DATATYPE_NULL) || recvtype == MPI_DATATYPE_NULL)
729 if ((sendbuf != MPI_IN_PLACE && sendcount < 0) || recvcount < 0)
730 return MPI_ERR_COUNT;
731 if (request == nullptr)
735 int rank = simgrid::s4u::this_actor::get_pid();
736 const void* real_sendbuf = sendbuf;
737 int real_sendcount = sendcount;
738 MPI_Datatype real_sendtype = sendtype;
739 std::unique_ptr<unsigned char[]> tmp_sendbuf;
740 if (sendbuf == MPI_IN_PLACE) {
741 tmp_sendbuf.reset(new unsigned char[recvcount * comm->size() * recvtype->size()]);
742 // memcpy(??,nullptr,0) is actually undefined behavor, even if harmless.
743 if (recvbuf != nullptr)
744 memcpy(tmp_sendbuf.get(), recvbuf, recvcount * comm->size() * recvtype->size());
745 real_sendbuf = tmp_sendbuf.get();
746 real_sendcount = recvcount;
747 real_sendtype = recvtype;
750 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Alltoall" : "PMPI_Ialltoall",
751 new simgrid::instr::CollTIData(
752 request == MPI_REQUEST_IGNORED ? "alltoall" : "ialltoall", -1, -1.0,
753 real_sendtype->is_replayable() ? real_sendcount : real_sendcount * real_sendtype->size(),
754 recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(),
755 simgrid::smpi::Datatype::encode(real_sendtype), simgrid::smpi::Datatype::encode(recvtype)));
757 if (request == MPI_REQUEST_IGNORED)
759 simgrid::smpi::Colls::alltoall(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcount, recvtype, comm);
761 retval = simgrid::smpi::Colls::ialltoall(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcount, recvtype,
764 TRACE_smpi_comm_out(rank);
769 int PMPI_Alltoallv(const void* sendbuf, const int* sendcounts, const int* senddisps, MPI_Datatype sendtype, void* recvbuf,
770 const int* recvcounts, const int* recvdisps, MPI_Datatype recvtype, MPI_Comm comm)
772 return PMPI_Ialltoallv(sendbuf, sendcounts, senddisps, sendtype, recvbuf, recvcounts, recvdisps, recvtype, comm, MPI_REQUEST_IGNORED);
775 int PMPI_Ialltoallv(const void* sendbuf, const int* sendcounts, const int* senddisps, MPI_Datatype sendtype, void* recvbuf,
776 const int* recvcounts, const int* recvdisps, MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request)
778 if (comm == MPI_COMM_NULL)
780 if (sendbuf == nullptr || recvbuf == nullptr)
781 return MPI_ERR_BUFFER;
782 if ((sendbuf != MPI_IN_PLACE && sendtype == MPI_DATATYPE_NULL) || recvtype == MPI_DATATYPE_NULL)
784 if ((sendbuf != MPI_IN_PLACE && (sendcounts == nullptr || senddisps == nullptr)) || recvcounts == nullptr ||
785 recvdisps == nullptr)
787 if (request == nullptr)
790 int rank = simgrid::s4u::this_actor::get_pid();
791 int size = comm->size();
792 for (int i = 0; i < size; i++) {
793 if (recvcounts[i] < 0 || (sendbuf != MPI_IN_PLACE && sendcounts[i] < 0))
794 return MPI_ERR_COUNT;
800 std::vector<int>* trace_sendcounts = new std::vector<int>;
801 std::vector<int>* trace_recvcounts = new std::vector<int>;
802 int dt_size_recv = recvtype->size();
804 const void* real_sendbuf = sendbuf;
805 const int* real_sendcounts = sendcounts;
806 const int* real_senddisps = senddisps;
807 MPI_Datatype real_sendtype = sendtype;
809 for (int i = 0; i < size; i++) { // copy data to avoid bad free
810 recv_size += recvcounts[i] * dt_size_recv;
811 trace_recvcounts->push_back(recvcounts[i] * dt_size_recv);
812 if (((recvdisps[i] + recvcounts[i]) * dt_size_recv) > maxsize)
813 maxsize = (recvdisps[i] + recvcounts[i]) * dt_size_recv;
816 std::unique_ptr<unsigned char[]> tmp_sendbuf;
817 std::unique_ptr<int[]> tmp_sendcounts;
818 std::unique_ptr<int[]> tmp_senddisps;
819 if (sendbuf == MPI_IN_PLACE) {
820 tmp_sendbuf.reset(new unsigned char[maxsize]);
821 real_sendbuf = memcpy(tmp_sendbuf.get(), recvbuf, maxsize);
822 tmp_sendcounts.reset(new int[size]);
823 std::copy(recvcounts, recvcounts + size, tmp_sendcounts.get());
824 real_sendcounts = tmp_sendcounts.get();
825 tmp_senddisps.reset(new int[size]);
826 std::copy(recvdisps, recvdisps + size, tmp_senddisps.get());
827 real_senddisps = tmp_senddisps.get();
828 real_sendtype = recvtype;
831 int dt_size_send = real_sendtype->size();
833 for (int i = 0; i < size; i++) { // copy data to avoid bad free
834 send_size += real_sendcounts[i] * dt_size_send;
835 trace_sendcounts->push_back(real_sendcounts[i] * dt_size_send);
838 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Alltoallv" : "PMPI_Ialltoallv",
839 new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "alltoallv" : "ialltoallv", -1,
840 send_size, trace_sendcounts, recv_size, trace_recvcounts,
841 simgrid::smpi::Datatype::encode(real_sendtype),
842 simgrid::smpi::Datatype::encode(recvtype)));
845 if (request == MPI_REQUEST_IGNORED)
846 retval = simgrid::smpi::Colls::alltoallv(real_sendbuf, real_sendcounts, real_senddisps, real_sendtype, recvbuf,
847 recvcounts, recvdisps, recvtype, comm);
849 retval = simgrid::smpi::Colls::ialltoallv(real_sendbuf, real_sendcounts, real_senddisps, real_sendtype, recvbuf,
850 recvcounts, recvdisps, recvtype, comm, request);
852 TRACE_smpi_comm_out(rank);
857 int PMPI_Alltoallw(const void* sendbuf, const int* sendcounts, const int* senddisps, const MPI_Datatype* sendtypes, void* recvbuf,
858 const int* recvcounts, const int* recvdisps, const MPI_Datatype* recvtypes, MPI_Comm comm)
860 return PMPI_Ialltoallw(sendbuf, sendcounts, senddisps, sendtypes, recvbuf, recvcounts, recvdisps, recvtypes, comm, MPI_REQUEST_IGNORED);
863 int PMPI_Ialltoallw(const void* sendbuf, const int* sendcounts, const int* senddisps, const MPI_Datatype* sendtypes, void* recvbuf,
864 const int* recvcounts, const int* recvdisps, const MPI_Datatype* recvtypes, MPI_Comm comm, MPI_Request* request)
866 if (comm == MPI_COMM_NULL)
868 if (sendbuf == nullptr || recvbuf == nullptr)
869 return MPI_ERR_BUFFER;
870 if ((sendbuf != MPI_IN_PLACE && sendtypes == nullptr) || recvtypes == nullptr)
872 if ((sendbuf != MPI_IN_PLACE && (sendcounts == nullptr || senddisps == nullptr)) || recvcounts == nullptr ||
873 recvdisps == nullptr)
875 if (request == nullptr)
879 int rank = simgrid::s4u::this_actor::get_pid();
880 int size = comm->size();
881 for (int i = 0; i < size; i++) { // copy data to avoid bad free
882 if (recvcounts[i] < 0 || (sendbuf != MPI_IN_PLACE && sendcounts[i] < 0))
883 return MPI_ERR_COUNT;
887 std::vector<int>* trace_sendcounts = new std::vector<int>;
888 std::vector<int>* trace_recvcounts = new std::vector<int>;
890 const void* real_sendbuf = sendbuf;
891 const int* real_sendcounts = sendcounts;
892 const int* real_senddisps = senddisps;
893 const MPI_Datatype* real_sendtypes = sendtypes;
894 unsigned long maxsize = 0;
895 for (int i = 0; i < size; i++) { // copy data to avoid bad free
896 if (recvtypes[i] == MPI_DATATYPE_NULL) {
897 delete trace_recvcounts;
898 delete trace_sendcounts;
901 recv_size += recvcounts[i] * recvtypes[i]->size();
902 trace_recvcounts->push_back(recvcounts[i] * recvtypes[i]->size());
903 if ((recvdisps[i] + (recvcounts[i] * recvtypes[i]->size())) > maxsize)
904 maxsize = recvdisps[i] + (recvcounts[i] * recvtypes[i]->size());
907 std::unique_ptr<unsigned char[]> tmp_sendbuf;
908 std::unique_ptr<int[]> tmp_sendcounts;
909 std::unique_ptr<int[]> tmp_senddisps;
910 std::unique_ptr<MPI_Datatype[]> tmp_sendtypes;
911 if (sendbuf == MPI_IN_PLACE) {
912 tmp_sendbuf.reset(new unsigned char[maxsize]);
913 real_sendbuf = memcpy(tmp_sendbuf.get(), recvbuf, maxsize);
914 tmp_sendcounts.reset(new int[size]);
915 std::copy(recvcounts, recvcounts + size, tmp_sendcounts.get());
916 real_sendcounts = tmp_sendcounts.get();
917 tmp_senddisps.reset(new int[size]);
918 std::copy(recvdisps, recvdisps + size, tmp_senddisps.get());
919 real_senddisps = tmp_senddisps.get();
920 tmp_sendtypes.reset(new MPI_Datatype[size]);
921 std::copy(recvtypes, recvtypes + size, tmp_sendtypes.get());
922 real_sendtypes = tmp_sendtypes.get();
925 for (int i = 0; i < size; i++) { // copy data to avoid bad free
926 send_size += real_sendcounts[i] * real_sendtypes[i]->size();
927 trace_sendcounts->push_back(real_sendcounts[i] * real_sendtypes[i]->size());
930 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Alltoallw" : "PMPI_Ialltoallw",
931 new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "alltoallv" : "ialltoallv", -1,
932 send_size, trace_sendcounts, recv_size, trace_recvcounts,
933 simgrid::smpi::Datatype::encode(real_sendtypes[0]),
934 simgrid::smpi::Datatype::encode(recvtypes[0])));
937 if (request == MPI_REQUEST_IGNORED)
938 retval = simgrid::smpi::Colls::alltoallw(real_sendbuf, real_sendcounts, real_senddisps, real_sendtypes, recvbuf,
939 recvcounts, recvdisps, recvtypes, comm);
941 retval = simgrid::smpi::Colls::ialltoallw(real_sendbuf, real_sendcounts, real_senddisps, real_sendtypes, recvbuf,
942 recvcounts, recvdisps, recvtypes, comm, request);
944 TRACE_smpi_comm_out(rank);