SOPT
Sparse OPTimisation
communicator.cc
Go to the documentation of this file.
2 #include <exception>
3 #include <mpi.h>
4 #include "sopt/mpi/session.h"
5 
6 namespace sopt::mpi {
7 
8 void Communicator::delete_comm(Communicator::Impl *const impl) {
9  if (impl->comm != MPI_COMM_WORLD and impl->comm != MPI_COMM_SELF and impl->comm != MPI_COMM_NULL)
10  MPI_Comm_free(&impl->comm);
11  delete impl;
12 }
13 
14 Communicator::Communicator(MPI_Comm const &comm) : impl(nullptr), session(session_singleton()) {
15  if (comm == MPI_COMM_NULL) return;
16  int size;
17  int rank;
18  MPI_Comm_size(comm, &size);
19  MPI_Comm_rank(comm, &rank);
20 
21  Impl const data{comm, static_cast<t_uint>(size), static_cast<t_uint>(rank)};
22  impl = std::shared_ptr<Impl const>(new Impl(data), &delete_comm);
23 }
24 
25 void Communicator::abort(const std::string &reason) const {
26  fprintf(stderr, "MPI Error on Rank %lu: %s\n", rank(), reason.c_str());
27  MPI_Abort(**this, MPI_ERR_OTHER);
28 }
29 
30 Communicator Communicator::duplicate() const {
31  if (not impl) return Communicator(MPI_COMM_NULL);
32  MPI_Comm comm;
33  MPI_Comm_dup(**this, &comm);
34  return comm;
35 }
36 
37 std::string Communicator::broadcast(std::string const &input, t_uint const root) const {
38  if (not impl) return input;
39  if (rank() == root) {
40  auto const N = broadcast(input.size(), root);
41  MPI_Bcast(const_cast<std::string::pointer>(input.c_str()), N,
42  Type<std::string::value_type>::value, root, **this);
43  return input;
44  }
45  auto const N = broadcast(input.size(), root);
46  std::string result(N, ' ');
47  MPI_Bcast(const_cast<std::string::pointer>(result.c_str()), N,
48  Type<std::string::value_type>::value, root, **this);
49  return result;
50 }
51 
52 } // namespace sopt::mpi
constexpr auto N
Definition: wavelets.cc:57
std::shared_ptr< details::initializer > session_singleton()
Definition: session.cc:50
size_t t_uint
Root of the type hierarchy for unsigned integers.
Definition: types.h:15