MPQC  3.0.0-alpha
comm.hpp
1 #ifndef MPQC_MPI_COMM_HPP
2 #define MPQC_MPI_COMM_HPP
3 
4 #include "mpqc/mpi.hpp"
5 #include "mpqc/utility/string.hpp"
6 
7 #include <stdarg.h>
8 
9 namespace mpqc {
10 namespace MPI {
11 
14  struct Comm {
15 
16  protected:
17 #ifdef HAVE_MPI
18  // warning: comm_ must be initialized prior to cout
19  MPI_Comm comm_;
20 #endif
21 
22  public:
23  struct OStream {
24  explicit OStream(const MPI::Comm &comm) {
25  this->rank_ = comm.rank();
26  }
27  std::ostream& stream() const {
28  std::cout << rank_ << ": ";
29  return std::cout;
30  }
31  private:
32  int rank_;
33  };
34 
35  const OStream cout;
36 
37  void printf(std::string fmt, ...) const {
38  fmt = string_cast(this->rank()) + ": " + fmt;
39  va_list args;
40  va_start(args, fmt);
41  vfprintf(stdout, fmt.c_str(), args);
42  va_end(args);
43  }
44 
45  // serial "MPI" stubs, PTP aren't available in serial
46 #ifndef HAVE_MPI
47  public:
48  static MPI::Comm Self() { return Comm(); }
49  static MPI::Comm World() { return Comm(); }
50  static MPI::Comm dup(MPI::Comm comm) { return Comm(); }
51  bool operator==(const Comm &comm) const { return true; }
52  void free() {}
53  int rank() const { return 0; }
54  int size() const { return 1; }
55  void barrier() const {}
56  private:
57  Comm() : cout(*this) {}
58 #endif
59 
60  // following functions are used if MPI is available
61 #ifdef HAVE_MPI
62 
63  public:
64 
65  static MPI::Comm Self() {
66  return Comm(MPI_COMM_SELF);
67  }
68 
69  static MPI::Comm World() {
70  return Comm(MPI_COMM_WORLD);
71  }
72 
74  static MPI::Comm dup(MPI::Comm comm) {
75  MPI_Comm dup;
76  MPI_Comm_dup(comm, &dup);
77  return Comm(dup);
78  }
79 
82  void free() {
83  MPI_Comm_free(&this->comm_);
84  }
85 
86  explicit Comm(MPI_Comm comm)
87  : comm_(comm), cout(*this) {}
88 
89  operator MPI_Comm() const {
90  return comm_;
91  }
92 
93  bool operator==(const Comm &comm) const {
94  return this->comm_ == comm.comm_;
95  }
96 
97  int rank() const {
98  MPQC_MPI_THREADSAFE;
99  int rank;
100  MPI_Comm_rank(comm_, &rank);
101  return rank;
102  }
103 
104  int size() const {
105  MPQC_MPI_THREADSAFE;
106  int size;
107  MPI_Comm_size(comm_, &size);
108  return size;
109  }
110 
111  void barrier() const {
112  //printf("%i: Comm::barrier()\n", this->rank());
113  MPI_Barrier(this->comm_);
114  }
115 
116  MPI_Status recv(void *data, int count, MPI_Datatype type,
117  int src, int tag) const {
118  MPQC_MPI_THREADSAFE;
119  // printf("MPI::recv(data=%p, count=%i, src=%i, tag=%i)\n",
120  // data, count, src, tag);
121  MPI_Status status;
122  MPI_Recv(data, count, type, src, tag, comm_, &status);
123  return status;
124  }
125 
126  template <typename T>
127  T recv(int src, int tag) const {
128  T data;
129  MPI::Comm::recv(&data, sizeof(T), MPI_BYTE, src, tag);
130  return data;
131  }
132 
133  void send(const void *data, int count, MPI_Datatype type,
134  int dst, int tag) const {
135  MPQC_MPI_THREADSAFE;
136  // printf("MPI::send(data=%p, count=%i, dst=%i, tag=%i)\n",
137  // data, count, dst, tag);
138  MPI_Send((void*)data, count, type, dst, tag, comm_);
139  }
140 
141  template <typename T>
142  void send(const T &data, int dst, int tag) const {
143  MPI::Comm::send((void*)&data, sizeof(T), MPI_BYTE, dst, tag);
144  }
145 
146 
147  void ssend(const void *data, int count, MPI_Datatype type,
148  int dst, int tag) const {
149  MPQC_MPI_THREADSAFE;
150  // printf("MPI::ssend(data=%p, count=%i, dst=%i, tag=%i)\n",
151  // data, count, dst, tag);
152  MPI_Ssend((void*)data, count, type, dst, tag, comm_);
153  }
154 
155  template <typename T>
156  void ssend(const T &data, int dst, int tag) const {
157  MPI::Comm::ssend((void*)&data, sizeof(T), MPI_BYTE, dst, tag);
158  }
159 
160  MPI_Request irecv(void *data, int count, MPI_Datatype type,
161  int src, int tag) const {
162  MPQC_MPI_THREADSAFE;
163  //printf("Thread::recv(count=%i, src=%i, tag=%i)\n", count, src, tag);
164  MPI_Request request;
165  MPI_Irecv(data, count, type, src, tag, comm_, &request);
166  return request;
167  }
168 
169  MPI_Request isend(const void *data, int count, MPI_Datatype type,
170  int src, int tag) const {
171  MPQC_MPI_THREADSAFE;
172  //printf("Thread::recv(count=%i, src=%i, tag=%i)\n", count, src, tag);
173  MPI_Request request;
174  MPI_Isend((void*)data, count, type, src, tag, comm_, &request);
175  return request;
176  }
177 
178 #endif // HAVE_MPI
179 
180  // Collective calls, either MPI or serial stubs
181 
182  public:
183 
184  template <typename T>
185  void broadcast(T &value, int root) const {
186 #ifdef HAVE_MPI
187  MPI_Bcast(&value, sizeof(T), MPI_BYTE, root, comm_);
188 #endif
189  }
190 
191  template <typename T>
192  void broadcast(T *data, int count, int root) const {
193 #ifdef HAVE_MPI
194  MPI_Bcast(data, sizeof(T)*count, MPI_BYTE, root, comm_);
195 #endif
196  }
197 
198  bool any(const bool &value) const {
199  int result = value;
200 #ifdef HAVE_MPI
201  MPI_Allreduce(MPI_IN_PLACE, &result, 1, MPI_INT, MPI_LOR, comm_);
202 #endif
203  return result;
204  }
205 
206  template <typename T>
207  void sum(T *value, int count) const {
208 #ifdef HAVE_MPI
209  MPI_Allreduce(MPI_IN_PLACE, value, count, MPI::type<T>(), MPI_SUM, comm_);
210 #endif
211  }
212 
213  template <typename T>
214  void sum(T &value) const {
215  sum(&value, 1);
216  }
217 
218  template <typename T>
219  std::vector<T> allgather(T value) const {
220  int bytes = sizeof(T);
221  std::vector<T> data(this->size(), T());
222 #ifdef HAVE_MPI
223  MPI_Allgather(&value, bytes, MPI_BYTE, &data[0], bytes, MPI_BYTE, comm_);
224 #else
225  data[0] = value;
226 #endif
227  return data;
228  }
229 
230  };
231 
232 
233  template <typename T>
234  std::ostream& operator<<(const MPI::Comm::OStream &s, const T &t) {
235  s.stream() << t;
236  return s.stream();
237  }
238 
239 
240 }
241 }
242 
243 #endif // MPQC_MPI_COMM_HPP
mpqc
Contains new MPQC code since version 3.
Definition: integralenginepool.hpp:37
mpqc::World
World is a wrapper around madness::World.
Definition: world.h:46
mpqc::string_cast
std::string string_cast(const T &value)
cast type T to string
Definition: string.hpp:14
mpqc::MPI::Comm::OStream
Definition: comm.hpp:23
mpqc::operator<<
void operator<<(Array< T > A, const V &v)
Write to Array from a generic vector V.
Definition: array.hpp:191
mpqc::MPI::Comm
MPI_Comm object wrapper/stub.
Definition: comm.hpp:14

Generated at Sun Jan 26 2020 23:24:01 for MPQC 3.0.0-alpha using the documentation package Doxygen 1.8.16.