Hex  1.0
Hydrogen-electron collision solver
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
parallel.h
Go to the documentation of this file.
1 /* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *\
2  * *
3  * / / / / __ \ \ / / *
4  * / /__ / / / _ \ \ \/ / *
5  * / ___ / | |/_/ / /\ \ *
6  * / / / / \_\ / / \ \ *
7  * *
8  * Jakub Benda (c) 2014 *
9  * Charles University in Prague *
10  * *
11 \* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
12 
13 #ifndef HEX_PARALLEL
14 #define HEX_PARALLEL
15 
16 #ifndef NO_MPI
17  #include <mpi.h>
18 #endif
19 
20 #include "arrays.h"
21 
29 class Parallel
30 {
31  public:
32 
43  Parallel (int* argc, char*** argv, bool active)
44  : active_(active), iproc_(0), Nproc_(1)
45  {
46 #ifndef NO_MPI
47  MPI_Init (argc, argv);
48  MPI_Comm_size (MPI_COMM_WORLD, &Nproc_);
49  MPI_Comm_rank (MPI_COMM_WORLD, &iproc_);
50 #else
51  active_ = false;
52 #endif
53  }
54 
56  {
57 #ifndef NO_MPI
58  if (active_)
59  MPI_Finalize();
60 #endif
61  }
62 
63  //
64  // getters
65  //
66 
73  inline bool IamMaster () const
74  {
75  return iproc_ == 0;
76  }
77 
87  inline bool isMyWork (int i) const
88  {
89  return i % Nproc_ == iproc_;
90  }
91 
95  inline bool active () const
96  {
97  return active_;
98  }
99 
103  inline int iproc () const
104  {
105  return iproc_;
106  }
107 
111  inline int Nproc () const
112  {
113  return Nproc_;
114  }
115 
129  template <class T> void sync (ArrayView<T> array, size_t chunksize, size_t Nchunk) const
130  {
131 #ifndef NO_MPI
132  if (active_)
133  {
134  for (unsigned ichunk = 0; ichunk < Nchunk; ichunk++)
135  {
136  // relevant process will broadcast this chunk's data
137  MPI_Bcast (
138  array.data() + ichunk * chunksize,
139  chunksize,
140  MPI_DOUBLE_COMPLEX,
141  ichunk % Nproc_,
142  MPI_COMM_WORLD
143  );
144  }
145  }
146 #endif
147  }
148 
154  void wait () const
155  {
156 #ifndef NO_MPI
157  if (active_)
158  MPI_Barrier(MPI_COMM_WORLD);
159 #endif
160  }
161 
162  private:
163 
164  // whether the MPI is on
165  bool active_;
166 
167  // communicator rank
168  int iproc_;
169 
170  // communicator size
171  int Nproc_;
172 };
173 
174 #endif
Array view.
Definition: arrays.h:186
int iproc() const
Returns the rank of the communicator (process is).
Definition: parallel.h:103
bool active() const
Returns true if the MPI is active.
Definition: parallel.h:95
Parallel(int *argc, char ***argv, bool active)
Constructor.
Definition: parallel.h:43
bool isMyWork(int i) const
Returns true if the work item is assigned to this process.
Definition: parallel.h:87
void sync(ArrayView< T > array, size_t chunksize, size_t Nchunk) const
Synchronize across processes.
Definition: parallel.h:129
MPI info.
Definition: parallel.h:29
~Parallel()
Definition: parallel.h:55
void wait() const
Wait for completition of all running tasks.
Definition: parallel.h:154
virtual T * data()
Pointer to the data.
Definition: arrays.h:280
int Nproc() const
Returns the size of the communicator (process count).
Definition: parallel.h:111
bool IamMaster() const
Returns true if this process is the master process.
Definition: parallel.h:73