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 && sendcount > 0) || ((comm->rank() == root) && recvbuf == nullptr && recvcount > 0))
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 CHECK_ARGS(comm == MPI_COMM_NULL, MPI_ERR_COMM, "Iallgather: the communicator cannot be MPI_COMM_NULL");
217 CHECK_ARGS(recvbuf == nullptr && recvcount > 0, MPI_ERR_BUFFER, "Iallgather: param 4 recvbuf cannot be NULL");
218 CHECK_ARGS(sendbuf == nullptr && sendcount > 0, MPI_ERR_BUFFER,
219 "Iallgather: param 1 sendbuf cannot be NULL when sendcound > 0");
220 CHECK_ARGS((sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL), MPI_ERR_TYPE,
221 "Iallgather: param 3 sendtype cannot be MPI_DATATYPE_NULL when sendbuff is not MPI_IN_PLACE");
222 CHECK_ARGS(recvtype == MPI_DATATYPE_NULL, MPI_ERR_TYPE, "Iallgather: param 6 recvtype cannot be MPI_DATATYPE_NULL");
223 CHECK_ARGS(recvcount < 0, MPI_ERR_COUNT, "Iallgather: param 5 recvcount cannot be negative");
224 CHECK_ARGS((sendbuf != MPI_IN_PLACE) && (sendcount < 0), MPI_ERR_COUNT,
225 "Iallgather: param 2 sendcount cannot be negative when sendbuf is not MPI_IN_PLACE");
226 CHECK_ARGS(request == nullptr, MPI_ERR_ARG, "Iallgather: param 8 request cannot be NULL");
229 if (sendbuf == MPI_IN_PLACE) {
230 sendbuf = static_cast<char*>(recvbuf) + recvtype->get_extent() * recvcount * comm->rank();
231 sendcount = recvcount;
234 int rank = simgrid::s4u::this_actor::get_pid();
236 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Allgather" : "PMPI_Iallggather",
237 new simgrid::instr::CollTIData(
238 request == MPI_REQUEST_IGNORED ? "allgather" : "iallgather", -1, -1.0,
239 sendtype->is_replayable() ? sendcount : sendcount * sendtype->size(),
240 recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(),
241 simgrid::smpi::Datatype::encode(sendtype), simgrid::smpi::Datatype::encode(recvtype)));
242 if (request == MPI_REQUEST_IGNORED)
243 simgrid::smpi::colls::allgather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
245 simgrid::smpi::colls::iallgather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm, request);
247 TRACE_smpi_comm_out(rank);
252 int PMPI_Allgatherv(const void *sendbuf, int sendcount, MPI_Datatype sendtype,
253 void *recvbuf, const int *recvcounts, const int *displs, MPI_Datatype recvtype, MPI_Comm comm){
254 return PMPI_Iallgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm, MPI_REQUEST_IGNORED);
257 int PMPI_Iallgatherv(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, const int* recvcounts, const int* displs,
258 MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request)
260 if (comm == MPI_COMM_NULL)
262 if (sendbuf == nullptr && sendcount > 0)
263 return MPI_ERR_BUFFER;
264 if (((sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) || (recvtype == MPI_DATATYPE_NULL))
266 if ((sendbuf != MPI_IN_PLACE) && (sendcount < 0))
267 return MPI_ERR_COUNT;
268 if (recvcounts == nullptr || displs == nullptr)
270 if (request == nullptr)
273 for (int i = 0; i < comm->size(); i++) {
274 if (recvcounts[i] < 0)
275 return MPI_ERR_COUNT;
276 else if (recvcounts[i] > 0 && recvbuf == nullptr)
277 return MPI_ERR_BUFFER;
281 if (sendbuf == MPI_IN_PLACE) {
282 sendbuf = static_cast<char*>(recvbuf) + recvtype->get_extent() * displs[comm->rank()];
283 sendcount = recvcounts[comm->rank()];
286 int rank = simgrid::s4u::this_actor::get_pid();
287 int dt_size_recv = recvtype->is_replayable() ? 1 : recvtype->size();
289 std::vector<int>* trace_recvcounts = new std::vector<int>;
290 for (int i = 0; i < comm->size(); i++) { // copy data to avoid bad free
291 trace_recvcounts->push_back(recvcounts[i] * dt_size_recv);
295 rank, request == MPI_REQUEST_IGNORED ? "PMPI_Allgatherv" : "PMPI_Iallgatherv",
296 new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "allgatherv" : "iallgatherv", -1,
297 sendtype->is_replayable() ? sendcount : sendcount * sendtype->size(), nullptr,
298 dt_size_recv, trace_recvcounts, simgrid::smpi::Datatype::encode(sendtype),
299 simgrid::smpi::Datatype::encode(recvtype)));
300 if (request == MPI_REQUEST_IGNORED)
301 simgrid::smpi::colls::allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm);
303 simgrid::smpi::colls::iallgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm,
306 TRACE_smpi_comm_out(rank);
311 int PMPI_Scatter(const void *sendbuf, int sendcount, MPI_Datatype sendtype,
312 void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm){
313 return PMPI_Iscatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm, MPI_REQUEST_IGNORED);
316 int PMPI_Iscatter(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
317 MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request* request)
319 if (comm == MPI_COMM_NULL)
321 if (((comm->rank() == root) && (sendtype == MPI_DATATYPE_NULL || not sendtype->is_valid())) ||
322 ((recvbuf != MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL || not recvtype->is_valid())))
324 if (((comm->rank() == root) && (sendcount < 0)) || ((recvbuf != MPI_IN_PLACE) && (recvcount < 0)))
325 return MPI_ERR_COUNT;
326 if ((sendbuf == recvbuf) || ((comm->rank() == root) && sendcount > 0 && (sendbuf == nullptr)) ||
327 (recvcount > 0 && recvbuf == nullptr))
328 return MPI_ERR_BUFFER;
329 if (root < 0 || root >= comm->size())
331 if (request == nullptr)
335 if (recvbuf == MPI_IN_PLACE) {
337 recvcount = sendcount;
339 int rank = simgrid::s4u::this_actor::get_pid();
341 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Scatter" : "PMPI_Iscatter",
342 new simgrid::instr::CollTIData(
343 request == MPI_REQUEST_IGNORED ? "scatter" : "iscatter", root, -1.0,
344 (comm->rank() != root || sendtype->is_replayable()) ? sendcount : sendcount * sendtype->size(),
345 recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(),
346 simgrid::smpi::Datatype::encode(sendtype), simgrid::smpi::Datatype::encode(recvtype)));
347 if (request == MPI_REQUEST_IGNORED)
348 simgrid::smpi::colls::scatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm);
350 simgrid::smpi::colls::iscatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm, request);
352 TRACE_smpi_comm_out(rank);
357 int PMPI_Scatterv(const void *sendbuf, const int *sendcounts, const int *displs,
358 MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm){
359 return PMPI_Iscatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm, MPI_REQUEST_IGNORED);
362 int PMPI_Iscatterv(const void* sendbuf, const int* sendcounts, const int* displs, MPI_Datatype sendtype, void* recvbuf, int recvcount,
363 MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request* request)
365 CHECK_ARGS(comm == MPI_COMM_NULL, MPI_ERR_COMM, "Iscatterv: the communicator cannot be MPI_COMM_NULL");
366 CHECK_ARGS((comm->rank() == root) && (sendcounts == nullptr), MPI_ERR_ARG,
367 "Iscatterv: param 2 sendcounts cannot be NULL on the root rank");
368 CHECK_ARGS((comm->rank() == root) && (displs == nullptr), MPI_ERR_ARG,
369 "Iscatterv: param 3 displs cannot be NULL on the root rank");
370 CHECK_ARGS((comm->rank() == root) && (sendtype == MPI_DATATYPE_NULL), MPI_ERR_TYPE,
371 "Iscatterv: The sendtype cannot be NULL on the root rank");
372 CHECK_ARGS((recvbuf != MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL), MPI_ERR_TYPE,
373 "Iscatterv: the recvtype cannot be NULL when not receiving in place");
374 CHECK_ARGS(request == nullptr, MPI_ERR_ARG, "Iscatterv: param 10 request cannot be NULL");
375 CHECK_ARGS(recvbuf != MPI_IN_PLACE && recvcount < 0, MPI_ERR_COUNT,
376 "Iscatterv: When not receiving in place, the recvcound cannot be negative");
377 CHECK_ARGS(root < 0, MPI_ERR_ROOT, "Iscatterv: root cannot be negative");
378 CHECK_ARGS(root >= comm->size(), MPI_ERR_ROOT, "Iscatterv: root (=%d) is larger than communicator size (=%d)", root,
381 if (comm->rank() == root) {
382 if (recvbuf == MPI_IN_PLACE) {
384 recvcount = sendcounts[comm->rank()];
386 for (int i = 0; i < comm->size(); i++)
387 CHECK_ARGS(sendcounts[i] < 0, MPI_ERR_COUNT, "Iscatterv: sendcounts[%d]=%d but this cannot be negative", i,
393 int rank = simgrid::s4u::this_actor::get_pid();
394 int dt_size_send = sendtype->is_replayable() ? 1 : sendtype->size();
396 std::vector<int>* trace_sendcounts = new std::vector<int>;
397 if (comm->rank() == root) {
398 for (int i = 0; i < comm->size(); i++) { // copy data to avoid bad free
399 trace_sendcounts->push_back(sendcounts[i] * dt_size_send);
403 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Scatterv" : "PMPI_Iscatterv",
404 new simgrid::instr::VarCollTIData(
405 request == MPI_REQUEST_IGNORED ? "scatterv" : "iscatterv", root, dt_size_send,
406 trace_sendcounts, recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(),
407 nullptr, simgrid::smpi::Datatype::encode(sendtype),
408 simgrid::smpi::Datatype::encode(recvtype)));
409 if (request == MPI_REQUEST_IGNORED)
410 simgrid::smpi::colls::scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm);
412 simgrid::smpi::colls::iscatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm,
415 TRACE_smpi_comm_out(rank);
420 int PMPI_Reduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
422 return PMPI_Ireduce(sendbuf, recvbuf, count, datatype, op, root, comm, MPI_REQUEST_IGNORED);
425 int PMPI_Ireduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm, MPI_Request* request)
427 if (comm == MPI_COMM_NULL)
429 if ((sendbuf == nullptr && count > 0) || ((comm->rank() == root) && recvbuf == nullptr))
430 return MPI_ERR_BUFFER;
431 if (datatype == MPI_DATATYPE_NULL || not datatype->is_valid())
433 if (op == MPI_OP_NULL)
435 if (request == nullptr)
437 if (root < 0 || root >= comm->size())
440 return MPI_ERR_COUNT;
443 int rank = simgrid::s4u::this_actor::get_pid();
445 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Reduce" : "PMPI_Ireduce",
446 new simgrid::instr::CollTIData(request == MPI_REQUEST_IGNORED ? "reduce" : "ireduce", root, 0,
447 datatype->is_replayable() ? count : count * datatype->size(), -1,
448 simgrid::smpi::Datatype::encode(datatype), ""));
449 if (request == MPI_REQUEST_IGNORED)
450 simgrid::smpi::colls::reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
452 simgrid::smpi::colls::ireduce(sendbuf, recvbuf, count, datatype, op, root, comm, request);
454 TRACE_smpi_comm_out(rank);
459 int PMPI_Reduce_local(const void* inbuf, void* inoutbuf, int count, MPI_Datatype datatype, MPI_Op op)
461 if (datatype == MPI_DATATYPE_NULL || not datatype->is_valid())
463 if (op == MPI_OP_NULL)
466 return MPI_ERR_COUNT;
469 op->apply(inbuf, inoutbuf, &count, datatype);
474 int PMPI_Allreduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
476 return PMPI_Iallreduce(sendbuf, recvbuf, count, datatype, op, comm, MPI_REQUEST_IGNORED);
479 int PMPI_Iallreduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request *request)
481 if (comm == MPI_COMM_NULL)
483 if ((sendbuf == nullptr && count > 0) || (recvbuf == nullptr))
484 return MPI_ERR_BUFFER;
485 if (datatype == MPI_DATATYPE_NULL || not datatype->is_valid())
488 return MPI_ERR_COUNT;
489 if (op == MPI_OP_NULL)
491 if (request == nullptr)
495 const void* real_sendbuf = sendbuf;
496 std::unique_ptr<unsigned char[]> tmp_sendbuf;
497 if (sendbuf == MPI_IN_PLACE) {
498 tmp_sendbuf.reset(new unsigned char[count * datatype->get_extent()]);
499 simgrid::smpi::Datatype::copy(recvbuf, count, datatype, tmp_sendbuf.get(), count, datatype);
500 real_sendbuf = tmp_sendbuf.get();
502 int rank = simgrid::s4u::this_actor::get_pid();
504 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Allreduce" : "PMPI_Iallreduce",
505 new simgrid::instr::CollTIData(request == MPI_REQUEST_IGNORED ? "allreduce" : "iallreduce", -1, 0,
506 datatype->is_replayable() ? count : count * datatype->size(), -1,
507 simgrid::smpi::Datatype::encode(datatype), ""));
509 if (request == MPI_REQUEST_IGNORED)
510 simgrid::smpi::colls::allreduce(real_sendbuf, recvbuf, count, datatype, op, comm);
512 simgrid::smpi::colls::iallreduce(real_sendbuf, recvbuf, count, datatype, op, comm, request);
514 TRACE_smpi_comm_out(rank);
519 int PMPI_Scan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
521 return PMPI_Iscan(sendbuf, recvbuf, count, datatype, op, comm, MPI_REQUEST_IGNORED);
524 int PMPI_Iscan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request* request)
526 if (comm == MPI_COMM_NULL)
528 if (datatype == MPI_DATATYPE_NULL || not datatype->is_valid())
530 if (op == MPI_OP_NULL)
532 if (request == nullptr)
535 return MPI_ERR_COUNT;
536 if (sendbuf == nullptr || recvbuf == nullptr)
537 return MPI_ERR_BUFFER;
540 int rank = simgrid::s4u::this_actor::get_pid();
541 const void* real_sendbuf = sendbuf;
542 std::unique_ptr<unsigned char[]> tmp_sendbuf;
543 if (sendbuf == MPI_IN_PLACE) {
544 tmp_sendbuf.reset(new unsigned char[count * datatype->size()]);
545 real_sendbuf = memcpy(tmp_sendbuf.get(), recvbuf, count * datatype->size());
547 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Scan" : "PMPI_Iscan",
548 new simgrid::instr::Pt2PtTIData(request == MPI_REQUEST_IGNORED ? "scan" : "iscan", -1,
549 datatype->is_replayable() ? count : count * datatype->size(),
550 simgrid::smpi::Datatype::encode(datatype)));
553 if (request == MPI_REQUEST_IGNORED)
554 retval = simgrid::smpi::colls::scan(real_sendbuf, recvbuf, count, datatype, op, comm);
556 retval = simgrid::smpi::colls::iscan(real_sendbuf, recvbuf, count, datatype, op, comm, request);
558 TRACE_smpi_comm_out(rank);
563 int PMPI_Exscan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
565 return PMPI_Iexscan(sendbuf, recvbuf, count, datatype, op, comm, MPI_REQUEST_IGNORED);
568 int PMPI_Iexscan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request* request){
569 if (comm == MPI_COMM_NULL)
571 if (not datatype->is_valid())
573 if (op == MPI_OP_NULL)
575 if (request == nullptr)
578 return MPI_ERR_COUNT;
579 if (sendbuf == nullptr || recvbuf == nullptr)
580 return MPI_ERR_BUFFER;
583 int rank = simgrid::s4u::this_actor::get_pid();
584 const void* real_sendbuf = sendbuf;
585 std::unique_ptr<unsigned char[]> tmp_sendbuf;
586 if (sendbuf == MPI_IN_PLACE) {
587 tmp_sendbuf.reset(new unsigned char[count * datatype->size()]);
588 real_sendbuf = memcpy(tmp_sendbuf.get(), recvbuf, count * datatype->size());
591 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Exscan" : "PMPI_Iexscan",
592 new simgrid::instr::Pt2PtTIData(request == MPI_REQUEST_IGNORED ? "exscan" : "iexscan", -1,
593 datatype->is_replayable() ? count : count * datatype->size(),
594 simgrid::smpi::Datatype::encode(datatype)));
597 if (request == MPI_REQUEST_IGNORED)
598 retval = simgrid::smpi::colls::exscan(real_sendbuf, recvbuf, count, datatype, op, comm);
600 retval = simgrid::smpi::colls::iexscan(real_sendbuf, recvbuf, count, datatype, op, comm, request);
602 TRACE_smpi_comm_out(rank);
607 int PMPI_Reduce_scatter(const void *sendbuf, void *recvbuf, const int *recvcounts, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
609 return PMPI_Ireduce_scatter(sendbuf, recvbuf, recvcounts, datatype, op, comm, MPI_REQUEST_IGNORED);
612 int PMPI_Ireduce_scatter(const void *sendbuf, void *recvbuf, const int *recvcounts, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request *request)
614 if (comm == MPI_COMM_NULL)
616 if ((sendbuf == nullptr) || (recvbuf == nullptr))
617 return MPI_ERR_BUFFER;
618 if (datatype == MPI_DATATYPE_NULL || not datatype->is_valid())
620 if (op == MPI_OP_NULL)
622 if (recvcounts == nullptr)
624 if (request == nullptr)
627 for (int i = 0; i < comm->size(); i++) { // copy data to avoid bad free
628 if (recvcounts[i] < 0)
629 return MPI_ERR_COUNT;
633 int rank = simgrid::s4u::this_actor::get_pid();
634 std::vector<int>* trace_recvcounts = new std::vector<int>;
635 int dt_send_size = datatype->is_replayable() ? 1 : datatype->size();
638 for (int i = 0; i < comm->size(); i++) { // copy data to avoid bad free
639 trace_recvcounts->push_back(recvcounts[i] * dt_send_size);
640 totalcount += recvcounts[i];
643 const void* real_sendbuf = sendbuf;
644 std::unique_ptr<unsigned char[]> tmp_sendbuf;
645 if (sendbuf == MPI_IN_PLACE) {
646 tmp_sendbuf.reset(new unsigned char[totalcount * datatype->size()]);
647 real_sendbuf = memcpy(tmp_sendbuf.get(), recvbuf, totalcount * datatype->size());
650 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Reduce_scatter" : "PMPI_Ireduce_scatter",
651 new simgrid::instr::VarCollTIData(
652 request == MPI_REQUEST_IGNORED ? "reducescatter" : "ireducescatter", -1, dt_send_size, nullptr,
653 -1, trace_recvcounts, simgrid::smpi::Datatype::encode(datatype), ""));
655 if (request == MPI_REQUEST_IGNORED)
656 simgrid::smpi::colls::reduce_scatter(real_sendbuf, recvbuf, recvcounts, datatype, op, comm);
658 simgrid::smpi::colls::ireduce_scatter(real_sendbuf, recvbuf, recvcounts, datatype, op, comm, request);
660 TRACE_smpi_comm_out(rank);
665 int PMPI_Reduce_scatter_block(const void *sendbuf, void *recvbuf, int recvcount,
666 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
668 return PMPI_Ireduce_scatter_block(sendbuf, recvbuf, recvcount, datatype, op, comm, MPI_REQUEST_IGNORED);
671 int PMPI_Ireduce_scatter_block(const void* sendbuf, void* recvbuf, int recvcount, MPI_Datatype datatype, MPI_Op op,
672 MPI_Comm comm, MPI_Request* request)
674 if (comm == MPI_COMM_NULL)
676 if (not datatype->is_valid())
678 if (op == MPI_OP_NULL)
682 if (request == nullptr)
686 int count = comm->size();
688 int rank = simgrid::s4u::this_actor::get_pid();
689 int dt_send_size = datatype->is_replayable() ? 1 : datatype->size();
690 std::vector<int>* trace_recvcounts = new std::vector<int>(recvcount * dt_send_size); // copy data to avoid bad free
692 const void* real_sendbuf = sendbuf;
693 std::unique_ptr<unsigned char[]> tmp_sendbuf;
694 if (sendbuf == MPI_IN_PLACE) {
695 tmp_sendbuf.reset(new unsigned char[recvcount * count * datatype->size()]);
696 real_sendbuf = memcpy(tmp_sendbuf.get(), recvbuf, recvcount * count * datatype->size());
700 rank, request == MPI_REQUEST_IGNORED ? "PMPI_Reduce_scatter_block" : "PMPI_Ireduce_scatter_block",
701 new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "reducescatter" : "ireducescatter", -1, 0,
702 nullptr, -1, trace_recvcounts, simgrid::smpi::Datatype::encode(datatype), ""));
704 int* recvcounts = new int[count];
705 for (int i = 0; i < count; i++)
706 recvcounts[i] = recvcount;
707 if (request == MPI_REQUEST_IGNORED)
708 simgrid::smpi::colls::reduce_scatter(real_sendbuf, recvbuf, recvcounts, datatype, op, comm);
710 simgrid::smpi::colls::ireduce_scatter(real_sendbuf, recvbuf, recvcounts, datatype, op, comm, request);
713 TRACE_smpi_comm_out(rank);
718 int PMPI_Alltoall(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
719 MPI_Datatype recvtype, MPI_Comm comm){
720 return PMPI_Ialltoall(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm, MPI_REQUEST_IGNORED);
723 int PMPI_Ialltoall(const void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
724 MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request)
726 if (comm == MPI_COMM_NULL)
728 if ((sendbuf == nullptr && sendcount > 0) || (recvbuf == nullptr && recvcount > 0))
729 return MPI_ERR_BUFFER;
730 if ((sendbuf != MPI_IN_PLACE && sendtype == MPI_DATATYPE_NULL) || recvtype == MPI_DATATYPE_NULL)
732 if ((sendbuf != MPI_IN_PLACE && sendcount < 0) || recvcount < 0)
733 return MPI_ERR_COUNT;
734 if (request == nullptr)
738 int rank = simgrid::s4u::this_actor::get_pid();
739 const void* real_sendbuf = sendbuf;
740 int real_sendcount = sendcount;
741 MPI_Datatype real_sendtype = sendtype;
742 std::unique_ptr<unsigned char[]> tmp_sendbuf;
743 if (sendbuf == MPI_IN_PLACE) {
744 tmp_sendbuf.reset(new unsigned char[recvcount * comm->size() * recvtype->size()]);
745 // memcpy(??,nullptr,0) is actually undefined behavior, even if harmless.
746 if (recvbuf != nullptr)
747 memcpy(tmp_sendbuf.get(), recvbuf, recvcount * comm->size() * recvtype->size());
748 real_sendbuf = tmp_sendbuf.get();
749 real_sendcount = recvcount;
750 real_sendtype = recvtype;
753 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Alltoall" : "PMPI_Ialltoall",
754 new simgrid::instr::CollTIData(
755 request == MPI_REQUEST_IGNORED ? "alltoall" : "ialltoall", -1, -1.0,
756 real_sendtype->is_replayable() ? real_sendcount : real_sendcount * real_sendtype->size(),
757 recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(),
758 simgrid::smpi::Datatype::encode(real_sendtype), simgrid::smpi::Datatype::encode(recvtype)));
760 if (request == MPI_REQUEST_IGNORED)
762 simgrid::smpi::colls::alltoall(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcount, recvtype, comm);
764 retval = simgrid::smpi::colls::ialltoall(real_sendbuf, real_sendcount, real_sendtype, recvbuf, recvcount, recvtype,
767 TRACE_smpi_comm_out(rank);
772 int PMPI_Alltoallv(const void* sendbuf, const int* sendcounts, const int* senddisps, MPI_Datatype sendtype, void* recvbuf,
773 const int* recvcounts, const int* recvdisps, MPI_Datatype recvtype, MPI_Comm comm)
775 return PMPI_Ialltoallv(sendbuf, sendcounts, senddisps, sendtype, recvbuf, recvcounts, recvdisps, recvtype, comm, MPI_REQUEST_IGNORED);
778 int PMPI_Ialltoallv(const void* sendbuf, const int* sendcounts, const int* senddisps, MPI_Datatype sendtype, void* recvbuf,
779 const int* recvcounts, const int* recvdisps, MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request)
781 if (comm == MPI_COMM_NULL)
783 if (sendbuf == nullptr || recvbuf == nullptr)
784 return MPI_ERR_BUFFER;
785 if ((sendbuf != MPI_IN_PLACE && sendtype == MPI_DATATYPE_NULL) || recvtype == MPI_DATATYPE_NULL)
787 if ((sendbuf != MPI_IN_PLACE && (sendcounts == nullptr || senddisps == nullptr)) || recvcounts == nullptr ||
788 recvdisps == nullptr)
790 if (request == nullptr)
793 int rank = simgrid::s4u::this_actor::get_pid();
794 int size = comm->size();
795 for (int i = 0; i < size; i++) {
796 if (recvcounts[i] < 0 || (sendbuf != MPI_IN_PLACE && sendcounts[i] < 0))
797 return MPI_ERR_COUNT;
803 std::vector<int>* trace_sendcounts = new std::vector<int>;
804 std::vector<int>* trace_recvcounts = new std::vector<int>;
805 int dt_size_recv = recvtype->size();
807 const void* real_sendbuf = sendbuf;
808 const int* real_sendcounts = sendcounts;
809 const int* real_senddisps = senddisps;
810 MPI_Datatype real_sendtype = sendtype;
812 for (int i = 0; i < size; i++) { // copy data to avoid bad free
813 recv_size += recvcounts[i] * dt_size_recv;
814 trace_recvcounts->push_back(recvcounts[i] * dt_size_recv);
815 if (((recvdisps[i] + recvcounts[i]) * dt_size_recv) > maxsize)
816 maxsize = (recvdisps[i] + recvcounts[i]) * dt_size_recv;
819 std::unique_ptr<unsigned char[]> tmp_sendbuf;
820 std::unique_ptr<int[]> tmp_sendcounts;
821 std::unique_ptr<int[]> tmp_senddisps;
822 if (sendbuf == MPI_IN_PLACE) {
823 tmp_sendbuf.reset(new unsigned char[maxsize]);
824 real_sendbuf = memcpy(tmp_sendbuf.get(), recvbuf, maxsize);
825 tmp_sendcounts.reset(new int[size]);
826 std::copy(recvcounts, recvcounts + size, tmp_sendcounts.get());
827 real_sendcounts = tmp_sendcounts.get();
828 tmp_senddisps.reset(new int[size]);
829 std::copy(recvdisps, recvdisps + size, tmp_senddisps.get());
830 real_senddisps = tmp_senddisps.get();
831 real_sendtype = recvtype;
834 int dt_size_send = real_sendtype->size();
836 for (int i = 0; i < size; i++) { // copy data to avoid bad free
837 send_size += real_sendcounts[i] * dt_size_send;
838 trace_sendcounts->push_back(real_sendcounts[i] * dt_size_send);
841 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Alltoallv" : "PMPI_Ialltoallv",
842 new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "alltoallv" : "ialltoallv", -1,
843 send_size, trace_sendcounts, recv_size, trace_recvcounts,
844 simgrid::smpi::Datatype::encode(real_sendtype),
845 simgrid::smpi::Datatype::encode(recvtype)));
848 if (request == MPI_REQUEST_IGNORED)
849 retval = simgrid::smpi::colls::alltoallv(real_sendbuf, real_sendcounts, real_senddisps, real_sendtype, recvbuf,
850 recvcounts, recvdisps, recvtype, comm);
852 retval = simgrid::smpi::colls::ialltoallv(real_sendbuf, real_sendcounts, real_senddisps, real_sendtype, recvbuf,
853 recvcounts, recvdisps, recvtype, comm, request);
855 TRACE_smpi_comm_out(rank);
860 int PMPI_Alltoallw(const void* sendbuf, const int* sendcounts, const int* senddisps, const MPI_Datatype* sendtypes, void* recvbuf,
861 const int* recvcounts, const int* recvdisps, const MPI_Datatype* recvtypes, MPI_Comm comm)
863 return PMPI_Ialltoallw(sendbuf, sendcounts, senddisps, sendtypes, recvbuf, recvcounts, recvdisps, recvtypes, comm, MPI_REQUEST_IGNORED);
866 int PMPI_Ialltoallw(const void* sendbuf, const int* sendcounts, const int* senddisps, const MPI_Datatype* sendtypes, void* recvbuf,
867 const int* recvcounts, const int* recvdisps, const MPI_Datatype* recvtypes, MPI_Comm comm, MPI_Request* request)
869 if (comm == MPI_COMM_NULL)
871 if (sendbuf == nullptr || recvbuf == nullptr)
872 return MPI_ERR_BUFFER;
873 if ((sendbuf != MPI_IN_PLACE && sendtypes == nullptr) || recvtypes == nullptr)
875 if ((sendbuf != MPI_IN_PLACE && (sendcounts == nullptr || senddisps == nullptr)) || recvcounts == nullptr ||
876 recvdisps == nullptr)
878 if (request == nullptr)
882 int rank = simgrid::s4u::this_actor::get_pid();
883 int size = comm->size();
884 for (int i = 0; i < size; i++) {
885 if (recvcounts[i] < 0 || (sendbuf != MPI_IN_PLACE && sendcounts[i] < 0))
886 return MPI_ERR_COUNT;
890 std::vector<int>* trace_sendcounts = new std::vector<int>;
891 std::vector<int>* trace_recvcounts = new std::vector<int>;
893 const void* real_sendbuf = sendbuf;
894 const int* real_sendcounts = sendcounts;
895 const int* real_senddisps = senddisps;
896 const MPI_Datatype* real_sendtypes = sendtypes;
897 unsigned long maxsize = 0;
898 for (int i = 0; i < size; i++) { // copy data to avoid bad free
899 if (recvtypes[i] == MPI_DATATYPE_NULL) {
900 delete trace_recvcounts;
901 delete trace_sendcounts;
904 recv_size += recvcounts[i] * recvtypes[i]->size();
905 trace_recvcounts->push_back(recvcounts[i] * recvtypes[i]->size());
906 if ((recvdisps[i] + (recvcounts[i] * recvtypes[i]->size())) > maxsize)
907 maxsize = recvdisps[i] + (recvcounts[i] * recvtypes[i]->size());
910 std::unique_ptr<unsigned char[]> tmp_sendbuf;
911 std::unique_ptr<int[]> tmp_sendcounts;
912 std::unique_ptr<int[]> tmp_senddisps;
913 std::unique_ptr<MPI_Datatype[]> tmp_sendtypes;
914 if (sendbuf == MPI_IN_PLACE) {
915 tmp_sendbuf.reset(new unsigned char[maxsize]);
916 real_sendbuf = memcpy(tmp_sendbuf.get(), recvbuf, maxsize);
917 tmp_sendcounts.reset(new int[size]);
918 std::copy(recvcounts, recvcounts + size, tmp_sendcounts.get());
919 real_sendcounts = tmp_sendcounts.get();
920 tmp_senddisps.reset(new int[size]);
921 std::copy(recvdisps, recvdisps + size, tmp_senddisps.get());
922 real_senddisps = tmp_senddisps.get();
923 tmp_sendtypes.reset(new MPI_Datatype[size]);
924 std::copy(recvtypes, recvtypes + size, tmp_sendtypes.get());
925 real_sendtypes = tmp_sendtypes.get();
928 for (int i = 0; i < size; i++) { // copy data to avoid bad free
929 send_size += real_sendcounts[i] * real_sendtypes[i]->size();
930 trace_sendcounts->push_back(real_sendcounts[i] * real_sendtypes[i]->size());
933 TRACE_smpi_comm_in(rank, request == MPI_REQUEST_IGNORED ? "PMPI_Alltoallw" : "PMPI_Ialltoallw",
934 new simgrid::instr::VarCollTIData(request == MPI_REQUEST_IGNORED ? "alltoallv" : "ialltoallv", -1,
935 send_size, trace_sendcounts, recv_size, trace_recvcounts,
936 simgrid::smpi::Datatype::encode(real_sendtypes[0]),
937 simgrid::smpi::Datatype::encode(recvtypes[0])));
940 if (request == MPI_REQUEST_IGNORED)
941 retval = simgrid::smpi::colls::alltoallw(real_sendbuf, real_sendcounts, real_senddisps, real_sendtypes, recvbuf,
942 recvcounts, recvdisps, recvtypes, comm);
944 retval = simgrid::smpi::colls::ialltoallw(real_sendbuf, real_sendcounts, real_senddisps, real_sendtypes, recvbuf,
945 recvcounts, recvdisps, recvtypes, comm, request);
947 TRACE_smpi_comm_out(rank);