High Performance Plasticity  0.5.0
mpiUtils.h
Go to the documentation of this file.
1 /* @file mpiUtils.h
2 * @author Michael Malahe
3 * @brief Header file for MPI helper functions
4 * @details
5 */
6 
7 #ifndef HPP_MPIUTILS_H
8 #define HPP_MPIUTILS_H
9 
10 #include "mpi.h"
11 #include <cstddef>
12 #include <stdexcept>
13 #include <iostream>
14 #include <vector>
15 #include <hpp/config.h>
16 
17 namespace hpp
18 {
19 
24 template <typename T>
25 MPI_Datatype MPIType() {
26  MPI_Datatype datatype;
27  if (std::is_same<T, double>::value) {
28  datatype = MPI_DOUBLE;
29  }
30  else if (std::is_same<T, float>::value) {
31  datatype = MPI_FLOAT;
32  }
33  else if (std::is_same<T, int>::value) {
34  datatype = MPI_INT;
35  }
36  else if (std::is_same<T, unsigned int>::value) {
37  datatype = MPI_UNSIGNED;
38  }
39  else if (std::is_same<T, long>::value) {
40  datatype = MPI_LONG;
41  }
42  else if (std::is_same<T, unsigned long>::value) {
43  datatype = MPI_UNSIGNED_LONG;
44  }
45  else {
46  throw std::runtime_error("Type not recognised.");
47  }
48  return datatype;
49 }
50 
51 template <typename T>
52 T MPIMin(const T& localVal, MPI_Comm comm) {
53  T globalVal;
54  MPI_Allreduce(&localVal, &globalVal, 1, MPIType<T>(), MPI_MIN, comm);
55  return globalVal;
56 }
57 
58 template <typename T>
59 T MPIMax(const T& localVal, MPI_Comm comm) {
60  T globalVal;
61  MPI_Allreduce(&localVal, &globalVal, 1, MPIType<T>(), MPI_MAX, comm);
62  return globalVal;
63 }
64 
65 template <typename T>
66 T MPISum(const T& localVal, MPI_Comm comm) {
67  T globalVal;
68  MPI_Allreduce(&localVal, &globalVal, 1, MPIType<T>(), MPI_SUM, comm);
69  return globalVal;
70 }
71 
72 template <typename T>
73 std::vector<T> MPIConcatOnRoot(T localVal, MPI_Comm comm) {
74  // Comm
75  int comm_size, comm_rank;
76  MPI_Comm_size(comm, &comm_size);
77  MPI_Comm_rank(comm, &comm_rank);
78 
79  // Gather
80  std::vector<T> globalVec;
81  if (comm_rank == 0) {
82  globalVec.resize(comm_size);
83  }
84  MPI_Datatype dataType = MPIType<T>();
85  MPI_Gather(&localVal, 1, dataType, globalVec.data(), 1, dataType, 0, comm);
86 
87  // Return
88  return globalVec;
89 }
90 
91 template <typename T>
92 std::vector<T> MPIConcatOnRoot(std::vector<T> localVec, MPI_Comm comm) {
93  // Comm
94  int comm_size, comm_rank;
95  MPI_Comm_size(comm, &comm_size);
96  MPI_Comm_rank(comm, &comm_rank);
97 
98  // Let root know how many values are coming from each process
99  int globalRecvCountSize = 0;
100  if (comm_rank == 0) {
101  globalRecvCountSize = comm_size;
102  }
103  std::vector<int> recvcounts(globalRecvCountSize);
104  int nLocal = localVec.size();
105  MPI_Gather(&nLocal, 1, MPI_INT, recvcounts.data(), 1, MPI_INT, 0, comm);
106 
107  // Get displacements and total size
108  int globalVecSize = 0;
109  std::vector<int> displs;
110  if (comm_rank == 0) {
111  displs.push_back(0);
112  for (int i=0; i<comm_size-1; i++) {
113  globalVecSize += recvcounts[i];
114  displs.push_back(globalVecSize);
115  }
116  globalVecSize += recvcounts[comm_size-1];
117  }
118 
119  // Gather the global vector on root
120  std::vector<T> globalVec(globalVecSize);
121  MPI_Datatype dataType = MPIType<T>();
122  MPI_Gatherv(localVec.data(), nLocal, dataType, globalVec.data(), recvcounts.data(),
123  displs.data(), dataType, 0, comm);
124 
125  // Return
126  return globalVec;
127 }
128 
129 template <typename T>
130 T MPIBroadcastFromRoot(T rootVal, MPI_Comm comm) {
131  // Comm
132  int comm_size, comm_rank;
133  MPI_Comm_size(comm, &comm_size);
134  MPI_Comm_rank(comm, &comm_rank);
135 
136  // Scatter value
137  T localVal = rootVal;
138  MPI_Datatype dataType = MPIType<T>();
139  MPI_Bcast(&localVal, 1, dataType, 0, comm);
140 
141  // Return
142  return localVal;
143 }
144 
145 template <typename T>
146 std::vector<T> MPIBroadcastFromRoot(std::vector<T> rootVec, MPI_Comm comm) {
147  // Comm
148  int comm_size, comm_rank;
149  MPI_Comm_size(comm, &comm_size);
150  MPI_Comm_rank(comm, &comm_rank);
151 
152  // Get size of root vector
153  // Size can be anything on non-root vectors
154  unsigned int size = rootVec.size();
155  size = MPIBroadcastFromRoot(size, comm);
156 
157  // Create local vector
158  std::vector<T> localVec(size);
159  if (comm_rank == 0) {
160  localVec = rootVec;
161  }
162 
163  // Scatter values
164  MPI_Datatype dataType = MPIType<T>();
165  MPI_Bcast(localVec.data(), size, dataType, 0, comm);
166 
167  // Return
168  return localVec;
169 }
170 
171 template <typename T>
172 std::vector<T> MPISplitVectorEvenly(const std::vector<T>& rootVec, MPI_Comm comm) {
173  // Comm
174  int comm_size, comm_rank;
175  MPI_Comm_size(comm, &comm_size);
176  MPI_Comm_rank(comm, &comm_rank);
177 
178  // Sizing of local vectors and MPI parameters
179  int rootVecSize = -1;
180  int localVecSize = -1;
181  int localVecSizeFinalProc = -1;
182  std::vector<int> sendcounts(comm_size);
183  std::vector<int> displs(comm_size);
184  if (comm_rank == 0) {
185  // Divide up vector
186  rootVecSize = rootVec.size();
187  localVecSize = rootVecSize/comm_size;
188  localVecSizeFinalProc = rootVecSize-localVecSize*(comm_size-1);
189 
190  // Create MPI parameters
191  displs[0] = 0;
192  for (int i=0; i<comm_size-1; i++) {
193  sendcounts[i] = localVecSize;
194  displs[i+1] = displs[i] + localVecSize;
195  }
196  sendcounts.back() = localVecSizeFinalProc;
197  }
198 
199  // Broadcast parameters from root
200  rootVecSize = MPIBroadcastFromRoot(rootVecSize, comm);
201  localVecSize = MPIBroadcastFromRoot(localVecSize, comm);
202  localVecSizeFinalProc = MPIBroadcastFromRoot(localVecSizeFinalProc, comm);
203  sendcounts = MPIBroadcastFromRoot(sendcounts, comm);
204  displs = MPIBroadcastFromRoot(displs, comm);
205 
206  // Final processor may have a different size
207  if (comm_rank == comm_size-1) {
208  localVecSize = localVecSizeFinalProc;
209  }
210 
211  // Scatter to local vectors
212  std::vector<T> localVec(localVecSize);
213  MPI_Scatterv(rootVec.data(), sendcounts.data(), displs.data(), MPIType<T>(),
214  localVec.data(), localVecSize, MPIType<T>(), 0, comm);
215 
216  // Return
217  return localVec;
218 }
219 
220 template <typename T>
221 std::vector<T> MPISplitVectorEvenly(const std::vector<T>& rootVec, MPI_Comm comm, MPI_Datatype dtype) {
222  // Comm
223  int comm_size, comm_rank;
224  MPI_Comm_size(comm, &comm_size);
225  MPI_Comm_rank(comm, &comm_rank);
226 
227  // Sizing of local vectors and MPI parameters
228  int rootVecSize = -1;
229  int localVecSize = -1;
230  int localVecSizeFinalProc = -1;
231  std::vector<int> sendcounts(comm_size);
232  std::vector<int> displs(comm_size);
233  if (comm_rank == 0) {
234  // Divide up vector
235  rootVecSize = rootVec.size();
236  localVecSize = rootVecSize/comm_size;
237  localVecSizeFinalProc = rootVecSize-localVecSize*(comm_size-1);
238 
239  // Create MPI parameters
240  displs[0] = 0;
241  for (int i=0; i<comm_size-1; i++) {
242  sendcounts[i] = localVecSize;
243  displs[i+1] = displs[i] + localVecSize;
244  }
245  sendcounts.back() = localVecSizeFinalProc;
246  }
247 
248  // Broadcast parameters from root
249  rootVecSize = MPIBroadcastFromRoot(rootVecSize, comm);
250  localVecSize = MPIBroadcastFromRoot(localVecSize, comm);
251  localVecSizeFinalProc = MPIBroadcastFromRoot(localVecSizeFinalProc, comm);
252  sendcounts = MPIBroadcastFromRoot(sendcounts, comm);
253  displs = MPIBroadcastFromRoot(displs, comm);
254 
255  // Final processor may have a different size
256  if (comm_rank == comm_size-1) {
257  localVecSize = localVecSizeFinalProc;
258  }
259 
260  // Scatter to local vectors
261  std::vector<T> localVec(localVecSize);
262  MPI_Scatterv(rootVec.data(), sendcounts.data(), displs.data(), dtype,
263  localVec.data(), localVecSize, dtype, 0, comm);
264 
265  // Return
266  return localVec;
267 }
268 
275 bool MPIAllTrue(bool condition, MPI_Comm comm);
276 
277 }//END NAMESPACE HPP
278 
279 #endif /* HPP_MPIUTILS_H */
Definition: casesUtils.cpp:4
T MPISum(const T &localVal, MPI_Comm comm)
Definition: mpiUtils.h:66
bool MPIAllTrue(bool condition, MPI_Comm comm)
Determine if condition is true for all processes.
Definition: mpiUtils.cpp:6
T MPIMin(const T &localVal, MPI_Comm comm)
Definition: mpiUtils.h:52
std::vector< T > MPIConcatOnRoot(T localVal, MPI_Comm comm)
Definition: mpiUtils.h:73
MPI_Datatype MPIType()
Get the MPI datatype from the C++ type.
Definition: mpiUtils.h:25
std::vector< T > MPISplitVectorEvenly(const std::vector< T > &rootVec, MPI_Comm comm)
Definition: mpiUtils.h:172
T MPIBroadcastFromRoot(T rootVal, MPI_Comm comm)
Definition: mpiUtils.h:130
T MPIMax(const T &localVal, MPI_Comm comm)
Definition: mpiUtils.h:59