MFEM  v3.3
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
communication.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
11 
12 #ifndef MFEM_COMMUNICATION
13 #define MFEM_COMMUNICATION
14 
15 #include "../config/config.hpp"
16 
17 #ifdef MFEM_USE_MPI
18 
19 #include "array.hpp"
20 #include "table.hpp"
21 #include "sets.hpp"
22 #include <mpi.h>
23 
24 namespace mfem
25 {
26 
27 /** @brief A simple convenience class that calls MPI_Init() at construction and
28  MPI_Finalize() at destruction. It also provides easy access to
29  MPI_COMM_WORLD's rank and size. */
31 {
32 protected:
35  {
36  MPI_Comm_rank(MPI_COMM_WORLD, &world_rank);
37  MPI_Comm_size(MPI_COMM_WORLD, &world_size);
38  }
39 public:
40  MPI_Session() { MPI_Init(NULL, NULL); GetRankAndSize(); }
41  MPI_Session(int &argc, char **&argv)
42  { MPI_Init(&argc, &argv); GetRankAndSize(); }
43  ~MPI_Session() { MPI_Finalize(); }
44  /// Return MPI_COMM_WORLD's rank.
45  int WorldRank() const { return world_rank; }
46  /// Return MPI_COMM_WORLD's size.
47  int WorldSize() const { return world_size; }
48  /// Return true if WorldRank() == 0.
49  bool Root() const { return world_rank == 0; }
50 };
51 
53 {
54 private:
55  MPI_Comm MyComm;
56 
57  /* The shared entities (e.g. vertices, faces and edges) are split into
58  groups, each group determined by the set of participating processors.
59  They are numbered locally in lproc. Assumptions:
60  - group 0 is the 'local' group
61  - groupmaster_lproc[0] = 0
62  - lproc_proc[0] = MyRank */
63 
64  // Neighbor ids (lproc) in each group.
65  Table group_lproc;
66  // Master neighbor id for each group.
67  Array<int> groupmaster_lproc;
68  // MPI rank of each neighbor.
69  Array<int> lproc_proc;
70  // Group --> Group number in the master.
71  Array<int> group_mgroup;
72 
73  void ProcToLProc();
74 
75 public:
76  GroupTopology() : MyComm(0) {}
77  GroupTopology(MPI_Comm comm) { MyComm = comm; }
78 
79  /// Copy constructor
80  GroupTopology(const GroupTopology &gt);
81  void SetComm(MPI_Comm comm) { MyComm = comm; }
82 
83  MPI_Comm GetComm() { return MyComm; }
84  int MyRank() { int r; MPI_Comm_rank(MyComm, &r); return r; }
85  int NRanks() { int s; MPI_Comm_size(MyComm, &s); return s; }
86 
87  void Create(ListOfIntegerSets &groups, int mpitag);
88 
89  int NGroups() const { return group_lproc.Size(); }
90  // return the number of neighbors including the local processor
91  int GetNumNeighbors() const { return lproc_proc.Size(); }
92  int GetNeighborRank(int i) const { return lproc_proc[i]; }
93  // am I master for group 'g'?
94  bool IAmMaster(int g) const { return (groupmaster_lproc[g] == 0); }
95  // return the neighbor index of the group master for a given group.
96  // neighbor 0 is the local processor
97  int GetGroupMaster(int g) const { return groupmaster_lproc[g]; }
98  // return the rank of the group master for a given group
99  int GetGroupMasterRank(int g) const
100  { return lproc_proc[groupmaster_lproc[g]]; }
101  // for a given group return the group number in the master
102  int GetGroupMasterGroup(int g) const { return group_mgroup[g]; }
103  // get the number of processors in a group
104  int GetGroupSize(int g) const { return group_lproc.RowSize(g); }
105  // return a pointer to a list of neighbors for a given group.
106  // neighbor 0 is the local processor
107  const int *GetGroup(int g) const { return group_lproc.GetRow(g); }
108 
109  /// Save the data in a stream.
110  void Save(std::ostream &out) const;
111  /// Load the data from a stream.
112  void Load(std::istream &in);
113 
114  virtual ~GroupTopology() {}
115 };
116 
118 {
119 private:
120  GroupTopology &gtopo;
121  Table group_ldof;
122  int group_buf_size;
123  Array<char> group_buf;
124  MPI_Request *requests;
125  MPI_Status *statuses;
126 
127 public:
129  /** Initialize the communicator from a local-dof to group map.
130  Finalize is called internally. */
131  void Create(Array<int> &ldof_group);
132  /** Fill-in the returned Table reference to initialize the communicator
133  then call Finalize. */
134  Table &GroupLDofTable() { return group_ldof; }
135  /// Allocate internal buffers after the GroupLDofTable is defined
136  void Finalize();
137 
138  /// Get a reference to the group topology object
139  GroupTopology & GetGroupTopology() { return gtopo; }
140 
141  /** @brief Broadcast within each group where the master is the root.
142  The data @a layout can be:
143 
144  0 - data is an array on all ldofs
145  1 - data is an array on the shared ldofs as given by group_ldof
146  */
147  template <class T> void Bcast(T *data, int layout);
148 
149  /** @brief Broadcast within each group where the master is the root.
150  This method is instantiated for int and double. */
151  template <class T> void Bcast(T *ldata) { Bcast<T>(ldata, 0); }
152  template <class T> void Bcast(Array<T> &ldata) { Bcast<T>((T *)ldata); }
153 
154  /** @brief Data structure on which we define reduce operations.
155 
156  The data is associated with (and the operation is performed on) one group
157  at a time. */
158  template <class T> struct OpData
159  {
160  int nldofs, nb, *ldofs;
161  T *ldata, *buf;
162  };
163 
164  /** @brief Reduce within each group where the master is the root.
165 
166  The reduce operation is given by the second argument (see below for list
167  of the supported operations.) This method is instantiated for int and
168  double. */
169  template <class T> void Reduce(T *ldata, void (*Op)(OpData<T>));
170  template <class T> void Reduce(Array<T> &ldata, void (*Op)(OpData<T>))
171  { Reduce<T>((T *)ldata, Op); }
172 
173  /// Reduce operation Sum, instantiated for int and double
174  template <class T> static void Sum(OpData<T>);
175  /// Reduce operation Min, instantiated for int and double
176  template <class T> static void Min(OpData<T>);
177  /// Reduce operation Max, instantiated for int and double
178  template <class T> static void Max(OpData<T>);
179  /// Reduce operation bitwise OR, instantiated for int only
180  template <class T> static void BitOR(OpData<T>);
181 
183 };
184 
185 
186 /// \brief Variable-length MPI message containing unspecific binary data.
187 template<int Tag>
189 {
190  std::string data;
191  MPI_Request send_request;
192 
193  /// Non-blocking send to processor 'rank'.
194  void Isend(int rank, MPI_Comm comm)
195  {
196  Encode();
197  MPI_Isend((void*) data.data(), data.length(), MPI_BYTE, rank, Tag, comm,
198  &send_request);
199  }
200 
201  /// Helper to send all messages in a rank-to-message map container.
202  template<typename MapT>
203  static void IsendAll(MapT& rank_msg, MPI_Comm comm)
204  {
205  typename MapT::iterator it;
206  for (it = rank_msg.begin(); it != rank_msg.end(); ++it)
207  {
208  it->second.Isend(it->first, comm);
209  }
210  }
211 
212  /// Helper to wait for all messages in a map container to be sent.
213  template<typename MapT>
214  static void WaitAllSent(MapT& rank_msg)
215  {
216  typename MapT::iterator it;
217  for (it = rank_msg.begin(); it != rank_msg.end(); ++it)
218  {
219  MPI_Wait(&it->second.send_request, MPI_STATUS_IGNORE);
220  it->second.Clear();
221  }
222  }
223 
224  /** Blocking probe for incoming message of this type from any rank.
225  Returns the rank and message size. */
226  static void Probe(int &rank, int &size, MPI_Comm comm)
227  {
228  MPI_Status status;
229  MPI_Probe(MPI_ANY_SOURCE, Tag, comm, &status);
230  rank = status.MPI_SOURCE;
231  MPI_Get_count(&status, MPI_BYTE, &size);
232  }
233 
234  /** Non-blocking probe for incoming message of this type from any rank.
235  If there is an incoming message, returns true and sets 'rank' and 'size'.
236  Otherwise returns false. */
237  static bool IProbe(int &rank, int &size, MPI_Comm comm)
238  {
239  int flag;
240  MPI_Status status;
241  MPI_Iprobe(MPI_ANY_SOURCE, Tag, comm, &flag, &status);
242  if (!flag) { return false; }
243 
244  rank = status.MPI_SOURCE;
245  MPI_Get_count(&status, MPI_BYTE, &size);
246  return true;
247  }
248 
249  /// Post-probe receive from processor 'rank' of message size 'size'.
250  void Recv(int rank, int size, MPI_Comm comm)
251  {
252  MFEM_ASSERT(size >= 0, "");
253  data.resize(size);
254  MPI_Status status;
255  MPI_Recv((void*) data.data(), size, MPI_BYTE, rank, Tag, comm, &status);
256 #ifdef MFEM_DEBUG
257  int count;
258  MPI_Get_count(&status, MPI_BYTE, &count);
259  MFEM_VERIFY(count == size, "");
260 #endif
261  Decode();
262  }
263 
264  /// Helper to receive all messages in a rank-to-message map container.
265  template<typename MapT>
266  static void RecvAll(MapT& rank_msg, MPI_Comm comm)
267  {
268  int recv_left = rank_msg.size();
269  while (recv_left > 0)
270  {
271  int rank, size;
272  Probe(rank, size, comm);
273  MFEM_ASSERT(rank_msg.find(rank) != rank_msg.end(), "Unexpected message"
274  " (tag " << Tag << ") from rank " << rank);
275  // NOTE: no guard against receiving two messages from the same rank
276  rank_msg[rank].Recv(rank, size, comm);
277  --recv_left;
278  }
279  }
280 
281  VarMessage() : send_request(MPI_REQUEST_NULL) {}
282  void Clear() { data.clear(); send_request = MPI_REQUEST_NULL; }
283 
284  virtual ~VarMessage()
285  {
286  MFEM_ASSERT(send_request == MPI_REQUEST_NULL,
287  "WaitAllSent was not called after Isend");
288  }
289 
290  VarMessage(const VarMessage &other)
291  : data(other.data), send_request(other.send_request)
292  {
293  MFEM_ASSERT(send_request == MPI_REQUEST_NULL,
294  "Cannot copy message with a pending send.");
295  }
296 
297 protected:
298  virtual void Encode() {}
299  virtual void Decode() {}
300 };
301 
302 
303 /// Helper struct to convert a C++ type to an MPI type
304 template <typename Type>
305 struct MPITypeMap { static const MPI_Datatype mpi_type; };
306 
307 
308 /** Reorder MPI ranks to follow the Z-curve within the physical machine topology
309  (provided that functions to query physical node coordinates are available).
310  Returns a new communicator with reordered ranks. */
311 MPI_Comm ReorderRanksZCurve(MPI_Comm comm);
312 
313 
314 } // namespace mfem
315 
316 #endif
317 
318 #endif
int WorldSize() const
Return MPI_COMM_WORLD&#39;s size.
int GetGroupMasterRank(int g) const
void Create(ListOfIntegerSets &groups, int mpitag)
int Size() const
Logical size of the array.
Definition: array.hpp:109
void Bcast(Array< T > &ldata)
int GetGroupMasterGroup(int g) const
GroupCommunicator(GroupTopology &gt)
Helper struct to convert a C++ type to an MPI type.
MPI_Comm ReorderRanksZCurve(MPI_Comm comm)
GroupTopology(MPI_Comm comm)
int GetGroupSize(int g) const
virtual void Decode()
const int * GetGroup(int g) const
void GetRow(int i, Array< int > &row) const
Return row i in array row (the Table must be finalized)
Definition: table.cpp:179
void Bcast(T *ldata)
Broadcast within each group where the master is the root. This method is instantiated for int and dou...
void Reduce(Array< T > &ldata, void(*Op)(OpData< T >))
void Recv(int rank, int size, MPI_Comm comm)
Post-probe receive from processor &#39;rank&#39; of message size &#39;size&#39;.
bool IAmMaster(int g) const
VarMessage(const VarMessage &other)
A simple convenience class that calls MPI_Init() at construction and MPI_Finalize() at destruction...
static void Sum(OpData< T >)
Reduce operation Sum, instantiated for int and double.
int GetNumNeighbors() const
void Finalize()
Allocate internal buffers after the GroupLDofTable is defined.
bool Root() const
Return true if WorldRank() == 0.
MPI_Request send_request
static void Min(OpData< T >)
Reduce operation Min, instantiated for int and double.
static void IsendAll(MapT &rank_msg, MPI_Comm comm)
Helper to send all messages in a rank-to-message map container.
int Size() const
Returns the number of TYPE I elements.
Definition: table.hpp:86
GroupTopology & GetGroupTopology()
Get a reference to the group topology object.
int GetGroupMaster(int g) const
virtual void Encode()
void Save(std::ostream &out) const
Save the data in a stream.
int WorldRank() const
Return MPI_COMM_WORLD&#39;s rank.
int GetNeighborRank(int i) const
static void Max(OpData< T >)
Reduce operation Max, instantiated for int and double.
static void RecvAll(MapT &rank_msg, MPI_Comm comm)
Helper to receive all messages in a rank-to-message map container.
Data structure on which we define reduce operations.
MPI_Session(int &argc, char **&argv)
void Reduce(T *ldata, void(*Op)(OpData< T >))
Reduce within each group where the master is the root.
void Isend(int rank, MPI_Comm comm)
Non-blocking send to processor &#39;rank&#39;.
int NGroups() const
void Create(Array< int > &ldof_group)
void SetComm(MPI_Comm comm)
static void WaitAllSent(MapT &rank_msg)
Helper to wait for all messages in a map container to be sent.
int RowSize(int i) const
Definition: table.hpp:102
static bool IProbe(int &rank, int &size, MPI_Comm comm)
List of integer sets.
Definition: sets.hpp:54
static const MPI_Datatype mpi_type
Variable-length MPI message containing unspecific binary data.
void Bcast(T *data, int layout)
Broadcast within each group where the master is the root. The data layout can be: ...
void Load(std::istream &in)
Load the data from a stream.
static void Probe(int &rank, int &size, MPI_Comm comm)
static void BitOR(OpData< T >)
Reduce operation bitwise OR, instantiated for int only.