3 #include "../constants.h"
4 #include "operations.hpp"
13 template<
class T,
class Op>
16 typedef detail::mpi_datatype<T> Datatype;
18 static void broadcast(
const communicator& comm, T& x,
int root)
21 MPI_Bcast(Datatype::address(x),
23 Datatype::datatype(), root, comm);
31 static void broadcast(
const communicator& comm, std::vector<T>& x,
int root)
35 int elem_size = Datatype::count(x[0]);
38 if (comm.rank() != root)
41 MPI_Bcast(Datatype::address(x[0]),
43 Datatype::datatype(), root, comm);
51 static request
ibroadcast(
const communicator& comm, T& x,
int root)
55 MPI_Ibcast(Datatype::address(x),
57 Datatype::datatype(), root, comm, &r.r);
63 DIY_UNSUPPORTED_MPI_CALL(MPI_Ibcast);
67 static void gather(
const communicator& comm,
const T&
in, std::vector<T>& out,
int root)
69 out.resize(comm.size());
71 MPI_Gather(Datatype::address(const_cast<T&>(in)),
74 Datatype::address(out[0]),
85 static void gather(
const communicator& comm,
const std::vector<T>& in, std::vector< std::vector<T> >& out,
int root)
88 std::vector<int> counts(comm.size());
89 int elem_size = Datatype::count(in[0]);
92 std::vector<int> offsets(comm.size(), 0);
93 for (
unsigned i = 1; i < offsets.size(); ++i)
94 offsets[i] = offsets[i-1] + counts[i-1];
96 std::vector<T> buffer((offsets.back() + counts.back()) / elem_size);
97 MPI_Gatherv(Datatype::address(const_cast<T&>(in[0])),
98 elem_size * in.size(),
100 Datatype::address(buffer[0]),
103 Datatype::datatype(),
106 out.resize(comm.size());
108 for (
unsigned i = 0; i < (unsigned)comm.size(); ++i)
110 out[i].reserve(counts[i] / elem_size);
111 for (
unsigned j = 0; j < (unsigned)(counts[i] / elem_size); ++j)
112 out[i].push_back(buffer[cur++]);
122 static void gather(
const communicator& comm,
const T& in,
int root)
125 MPI_Gather(Datatype::address(const_cast<T&>(in)),
127 Datatype::datatype(),
128 Datatype::address(const_cast<T&>(in)),
130 Datatype::datatype(),
136 DIY_UNSUPPORTED_MPI_CALL(
"MPI_Gather");
140 static void gather(
const communicator& comm,
const std::vector<T>& in,
int root)
143 int elem_size = Datatype::count(in[0]);
146 MPI_Gatherv(Datatype::address(const_cast<T&>(in[0])),
147 elem_size * in.size(),
148 Datatype::datatype(),
150 Datatype::datatype(),
156 DIY_UNSUPPORTED_MPI_CALL(
"MPI_Gatherv");
160 static void all_gather(
const communicator& comm,
const T& in, std::vector<T>& out)
162 out.resize(comm.size());
164 MPI_Allgather(Datatype::address(const_cast<T&>(in)),
166 Datatype::datatype(),
167 Datatype::address(out[0]),
169 Datatype::datatype(),
177 static void all_gather(
const communicator& comm,
const std::vector<T>& in, std::vector< std::vector<T> >& out)
180 std::vector<int> counts(comm.size());
181 int elem_size = Datatype::count(in[0]);
184 std::vector<int> offsets(comm.size(), 0);
185 for (
unsigned i = 1; i < offsets.size(); ++i)
186 offsets[i] = offsets[i-1] + counts[i-1];
188 std::vector<T> buffer((offsets.back() + counts.back()) / elem_size);
189 MPI_Allgatherv(Datatype::address(const_cast<T&>(in[0])),
190 elem_size * in.size(),
191 Datatype::datatype(),
192 Datatype::address(buffer[0]),
195 Datatype::datatype(),
198 out.resize(comm.size());
200 for (
int i = 0; i < comm.size(); ++i)
202 out[i].reserve(counts[i] / elem_size);
203 for (
int j = 0; j < (int)(counts[i] / elem_size); ++j)
204 out[i].push_back(buffer[cur++]);
213 static void reduce(
const communicator& comm,
const T& in, T& out,
int root,
const Op&)
216 MPI_Reduce(Datatype::address(const_cast<T&>(in)),
217 Datatype::address(out),
219 Datatype::datatype(),
220 detail::mpi_op<Op>::get(),
229 static void reduce(
const communicator& comm,
const T& in,
int root,
const Op&)
232 MPI_Reduce(Datatype::address(const_cast<T&>(in)),
233 Datatype::address(const_cast<T&>(in)),
235 Datatype::datatype(),
236 detail::mpi_op<Op>::get(),
242 DIY_UNSUPPORTED_MPI_CALL(
"MPI_Reduce");
246 static void all_reduce(
const communicator& comm,
const T& in, T& out,
const Op&)
249 MPI_Allreduce(Datatype::address(const_cast<T&>(in)),
250 Datatype::address(out),
252 Datatype::datatype(),
253 detail::mpi_op<Op>::get(),
261 static void all_reduce(
const communicator& comm,
const std::vector<T>& in, std::vector<T>& out,
const Op&)
264 out.resize(in.size());
265 int elem_size = Datatype::count(in[0]);
266 MPI_Allreduce(Datatype::address(const_cast<T&>(in[0])),
267 Datatype::address(out[0]),
268 elem_size * in.size(),
269 Datatype::datatype(),
270 detail::mpi_op<Op>::get(),
278 static void scan(
const communicator& comm,
const T& in, T& out,
const Op&)
281 MPI_Scan(Datatype::address(const_cast<T&>(in)),
282 Datatype::address(out),
284 Datatype::datatype(),
285 detail::mpi_op<Op>::get(),
293 static void all_to_all(
const communicator& comm,
const std::vector<T>& in, std::vector<T>& out,
int n = 1)
296 int elem_size = Datatype::count(in[0]);
298 MPI_Alltoall(Datatype::address(const_cast<T&>(in[0])),
300 Datatype::datatype(),
301 Datatype::address(out[0]),
303 Datatype::datatype(),
345 void gather(
const communicator& comm,
const std::vector<T>& in, std::vector< std::vector<T> >& out,
int root)
381 template<
class T,
class Op>
388 template<
class T,
class Op>
395 template<
class T,
class Op>
402 template<
class T,
class Op>
409 template<
class T,
class Op>
request ibroadcast(const communicator &comm, T &x, int root)
iBroadcast to all processes in comm.
Definition: collectives.hpp:329
void all_reduce(const communicator &comm, const std::vector< T > &in, std::vector< T > &out, const Op &op)
Same as above, but for vectors.
Definition: collectives.hpp:403
void all_to_all(const communicator &comm, const std::vector< T > &in, std::vector< T > &out, int n=1)
all_to_all
Definition: collectives.hpp:417
void reduce(const communicator &comm, const T &in, T &out, int root, const Op &op)
reduce
Definition: collectives.hpp:382
void broadcast(const communicator &comm, std::vector< T > &x, int root)
Broadcast for vectors.
Definition: collectives.hpp:322
void scan(const communicator &comm, const T &in, T &out, const Op &op)
scan
Definition: collectives.hpp:410
void in(const RegularLink< Bounds > &link, const Point &p, OutIter out, const Bounds &domain)
Finds the neighbor(s) containing the target point.
Definition: pick.hpp:102
Simple wrapper around MPI_Comm.
Definition: communicator.hpp:8
void all_gather(const communicator &comm, const T &in, std::vector< T > &out)
all_gather from all processes in comm. out is resized to comm.size() and filled with elements from th...
Definition: collectives.hpp:368
void reduce(const communicator &comm, const T &in, int root, const Op &op)
Simplified version (without out) for use on non-root processes.
Definition: collectives.hpp:389
void all_reduce(const communicator &comm, const T &in, T &out, const Op &op)
all_reduce
Definition: collectives.hpp:396
Definition: request.hpp:5
void all_gather(const communicator &comm, const std::vector< T > &in, std::vector< std::vector< T > > &out)
Same as above, but for vectors.
Definition: collectives.hpp:375
void broadcast(const communicator &comm, T &x, int root)
Broadcast to all processes in comm.
Definition: collectives.hpp:315
void gather(const communicator &comm, const T &in, std::vector< T > &out, int root)
Gather from all processes in comm. On root process, out is resized to comm.size() and filled with ele...
Definition: collectives.hpp:338
void gather(const communicator &comm, const std::vector< T > &in, int root)
Simplified version (without out) for use on non-root processes.
Definition: collectives.hpp:359