12 #ifndef MFEM_COMMUNICATION
13 #define MFEM_COMMUNICATION
15 #include "../config/config.hpp"
76 void SetComm(MPI_Comm comm) { MyComm = comm; }
79 int MyRank() {
int r; MPI_Comm_rank(MyComm, &r);
return r; }
80 int NRanks() {
int s; MPI_Comm_size(MyComm, &s);
return s; }
89 bool IAmMaster(
int g)
const {
return (groupmaster_lproc[g] == 0); }
95 {
return lproc_proc[groupmaster_lproc[g]]; }
114 MPI_Request *requests;
115 MPI_Status *statuses;
133 template <
class T>
void Bcast(T *ldata);
150 { Reduce<T>((T *)ldata, Op); }
153 template <
class T>
static void Sum(OpData<T>);
155 template <
class T>
static void Min(OpData<T>);
157 template <
class T>
static void Max(OpData<T>);
159 template <
class T>
static void BitOR(OpData<T>);
176 MPI_Isend((
void*)
data.data(),
data.length(), MPI_BYTE, rank, Tag, comm,
181 template<
typename MapT>
182 static void IsendAll(MapT& rank_msg, MPI_Comm comm)
184 typename MapT::iterator it;
185 for (it = rank_msg.begin(); it != rank_msg.end(); ++it)
187 it->second.Isend(it->first, comm);
192 template<
typename MapT>
195 typename MapT::iterator it;
196 for (it = rank_msg.begin(); it != rank_msg.end(); ++it)
198 MPI_Wait(&it->second.send_request, MPI_STATUS_IGNORE);
205 static void Probe(
int &rank,
int &size, MPI_Comm comm)
208 MPI_Probe(MPI_ANY_SOURCE, Tag, comm, &status);
209 rank = status.MPI_SOURCE;
210 MPI_Get_count(&status, MPI_BYTE, &size);
216 static bool IProbe(
int &rank,
int &size, MPI_Comm comm)
220 MPI_Iprobe(MPI_ANY_SOURCE, Tag, comm, &flag, &status);
221 if (!flag) {
return false; }
223 rank = status.MPI_SOURCE;
224 MPI_Get_count(&status, MPI_BYTE, &size);
229 void Recv(
int rank,
int size, MPI_Comm comm)
231 MFEM_ASSERT(size >= 0,
"");
234 MPI_Recv((
void*)
data.data(), size, MPI_BYTE, rank, Tag, comm, &status);
237 MPI_Get_count(&status, MPI_BYTE, &count);
238 MFEM_VERIFY(count == size,
"");
244 template<
typename MapT>
245 static void RecvAll(MapT& rank_msg, MPI_Comm comm)
247 int recv_left = rank_msg.size();
248 while (recv_left > 0)
251 Probe(rank, size, comm);
252 MFEM_ASSERT(rank_msg.find(rank) != rank_msg.end(),
"Unexpected message"
253 " (tag " << Tag <<
") from rank " << rank);
255 rank_msg[rank].Recv(rank, size, comm);
266 "WaitAllSent was not called after Isend");
273 "Cannot copy message with a pending send.");
283 template <
typename Type>
int WorldSize() const
Return MPI_COMM_WORLD's size.
int GetGroupMasterRank(int g) const
void Create(ListOfIntegerSets &groups, int mpitag)
int Size() const
Logical size of the array.
void Bcast(Array< T > &ldata)
int GetGroupMasterGroup(int g) const
GroupCommunicator(GroupTopology >)
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
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)
void Reduce(Array< T > &ldata, void(*Op)(OpData< T >))
void Recv(int rank, int size, MPI_Comm comm)
Post-probe receive from processor 'rank' of message size 'size'.
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.
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.
GroupTopology & GetGroupTopology()
Get a reference to the group topology object.
int GetGroupMaster(int g) const
int WorldRank() const
Return MPI_COMM_WORLD'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.
MPI_Session(int &argc, char **&argv)
void Reduce(T *ldata, void(*Op)(OpData< T >))
void Isend(int rank, MPI_Comm comm)
Non-blocking send to processor 'rank'.
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.
static bool IProbe(int &rank, int &size, MPI_Comm comm)
static const MPI_Datatype mpi_type
Variable-length MPI message containing unspecific binary data.
static void Probe(int &rank, int &size, MPI_Comm comm)
static void BitOR(OpData< T >)
Reduce operation bitwise OR, instantiated for int only.