MFEM v4.8.0
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-2025, Lawrence Livermore National Security, LLC. Produced
2// at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3// LICENSE and NOTICE for details. LLNL-CODE-806117.
4//
5// This file is part of the MFEM library. For more information and source code
6// availability visit https://mfem.org.
7//
8// MFEM is free software; you can redistribute it and/or modify it under the
9// terms of the BSD-3 license. We welcome feedback and contributions, see file
10// CONTRIBUTING.md for details.
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 "globals.hpp"
23#include <mpi.h>
24
25// can't directly use MPI_CXX_BOOL because Microsoft's MPI implementation
26// doesn't include MPI_CXX_BOOL. Fallback to MPI_C_BOOL if unavailable.
27#ifdef MPI_CXX_BOOL
28#define MFEM_MPI_CXX_BOOL MPI_CXX_BOOL
29#else
30#define MFEM_MPI_CXX_BOOL MPI_C_BOOL
31#endif
32
33namespace mfem
34{
35
36/** @brief A simple singleton class that calls MPI_Init() at construction and
37 MPI_Finalize() at destruction. It also provides easy access to
38 MPI_COMM_WORLD's rank and size. */
39class Mpi
40{
41public:
42 /// Singleton creation with Mpi::Init(argc, argv).
43 static void Init(int &argc, char **&argv,
44 int required = default_thread_required,
45 int *provided = nullptr)
46 { Init(&argc, &argv, required, provided); }
47 /// Singleton creation with Mpi::Init().
48 static void Init(int *argc = nullptr, char ***argv = nullptr,
49 int required = default_thread_required,
50 int *provided = nullptr)
51 {
52 MFEM_VERIFY(!IsInitialized(), "MPI already initialized!");
53 if (required == MPI_THREAD_SINGLE)
54 {
55 int mpi_err = MPI_Init(argc, argv);
56 MFEM_VERIFY(!mpi_err, "error in MPI_Init()!");
57 if (provided) { *provided = MPI_THREAD_SINGLE; }
58 }
59 else
60 {
61 int mpi_provided;
62 int mpi_err = MPI_Init_thread(argc, argv, required, &mpi_provided);
63 MFEM_VERIFY(!mpi_err, "error in MPI_Init()!");
64 if (provided) { *provided = mpi_provided; }
65 }
66 // The Mpi singleton object below needs to be created after MPI_Init() for
67 // some MPI implementations.
68 Singleton();
69 }
70 /// Finalize MPI (if it has been initialized and not yet already finalized).
71 static void Finalize()
72 {
73 if (IsInitialized() && !IsFinalized()) { MPI_Finalize(); }
74 }
75 /// Return true if MPI has been initialized.
76 static bool IsInitialized()
77 {
78 int mpi_is_initialized;
79 int mpi_err = MPI_Initialized(&mpi_is_initialized);
80 return (mpi_err == MPI_SUCCESS) && mpi_is_initialized;
81 }
82 /// Return true if MPI has been finalized.
83 static bool IsFinalized()
84 {
85 int mpi_is_finalized;
86 int mpi_err = MPI_Finalized(&mpi_is_finalized);
87 return (mpi_err == MPI_SUCCESS) && mpi_is_finalized;
88 }
89 /// Return the MPI rank in MPI_COMM_WORLD.
90 static int WorldRank()
91 {
92 int world_rank;
93 MPI_Comm_rank(MPI_COMM_WORLD, &world_rank);
94 return world_rank;
95 }
96 /// Return the size of MPI_COMM_WORLD.
97 static int WorldSize()
98 {
99 int world_size;
100 MPI_Comm_size(MPI_COMM_WORLD, &world_size);
101 return world_size;
102 }
103 /// Return true if the rank in MPI_COMM_WORLD is zero.
104 static bool Root() { return WorldRank() == 0; }
105 /// Default level of thread support for MPI_Init_thread.
106 static MFEM_EXPORT int default_thread_required;
107private:
108 /// Initialize the Mpi singleton.
109 static Mpi &Singleton()
110 {
111 static Mpi mpi;
112 return mpi;
113 }
114 /// Finalize MPI.
115 ~Mpi() { Finalize(); }
116 /// Prevent direct construction of objects of this class.
117 Mpi() {}
118};
119
120/** @brief A simple convenience class based on the Mpi singleton class above.
121 Preserved for backward compatibility. New code should use Mpi::Init() and
122 other Mpi methods instead. */
124{
125public:
127 MPI_Session(int &argc, char **&argv) { Mpi::Init(argc, argv); }
128 /// Return MPI_COMM_WORLD's rank.
129 int WorldRank() const { return Mpi::WorldRank(); }
130 /// Return MPI_COMM_WORLD's size.
131 int WorldSize() const { return Mpi::WorldSize(); }
132 /// Return true if WorldRank() == 0.
133 bool Root() const { return Mpi::Root(); }
134};
135
136
137/** The shared entities (e.g. vertices, faces and edges) are split into groups,
138 each group determined by the set of participating processors. They are
139 numbered locally in lproc. Assumptions:
140 - group 0 is the 'local' group
141 - groupmaster_lproc[0] = 0
142 - lproc_proc[0] = MyRank */
144{
145private:
146 MPI_Comm MyComm;
147
148 /// Neighbor ids (lproc) in each group.
149 Table group_lproc;
150 /// Master neighbor id for each group.
151 Array<int> groupmaster_lproc;
152 /// MPI rank of each neighbor.
153 Array<int> lproc_proc;
154 /// Group --> Group number in the master.
155 Array<int> group_mgroup;
156
157 void ProcToLProc();
158
159public:
160 /// Constructor with the MPI communicator = 0.
161 GroupTopology() : MyComm(0) {}
162
163 /// Constructor given the MPI communicator 'comm'.
164 GroupTopology(MPI_Comm comm) { MyComm = comm; }
165
166 /// Copy constructor
167 GroupTopology(const GroupTopology &gt);
168
169 /// Set the MPI communicator to 'comm'.
170 void SetComm(MPI_Comm comm) { MyComm = comm; }
171
172 /// Return the MPI communicator.
173 MPI_Comm GetComm() const { return MyComm; }
174
175 /// Return the MPI rank within this object's communicator.
176 int MyRank() const { int r; MPI_Comm_rank(MyComm, &r); return r; }
177
178 /// Return the number of MPI ranks within this object's communicator.
179 int NRanks() const { int s; MPI_Comm_size(MyComm, &s); return s; }
180
181 /// Set up the group topology given the list of sets of shared entities.
182 void Create(ListOfIntegerSets &groups, int mpitag);
183
184 /// Return the number of groups.
185 int NGroups() const { return group_lproc.Size(); }
186
187 /// Return the number of neighbors including the local processor.
188 int GetNumNeighbors() const { return lproc_proc.Size(); }
189
190 /// Return the MPI rank of neighbor 'i'.
191 int GetNeighborRank(int i) const { return lproc_proc[i]; }
192
193 /// Return true if I am master for group 'g'.
194 bool IAmMaster(int g) const { return (groupmaster_lproc[g] == 0); }
195
196 /** @brief Return the neighbor index of the group master for a given group.
197 Neighbor 0 is the local processor. */
198 int GetGroupMaster(int g) const { return groupmaster_lproc[g]; }
199
200 /// Return the rank of the group master for group 'g'.
201 int GetGroupMasterRank(int g) const
202 { return lproc_proc[groupmaster_lproc[g]]; }
203
204 /// Return the group number in the master for group 'g'.
205 int GetGroupMasterGroup(int g) const { return group_mgroup[g]; }
206
207 /// Get the number of processors in a group
208 int GetGroupSize(int g) const { return group_lproc.RowSize(g); }
209
210 /** @brief Return a pointer to a list of neighbors for a given group.
211 Neighbor 0 is the local processor */
212 const int *GetGroup(int g) const { return group_lproc.GetRow(g); }
213
214 /// Save the data in a stream.
215 void Save(std::ostream &out) const;
216
217 /// Load the data from a stream.
218 void Load(std::istream &in);
219
220 /// Copy the internal data to the external 'copy'.
221 void Copy(GroupTopology & copy) const;
222
223 /// Swap the internal data with another @a GroupTopology object.
224 void Swap(GroupTopology &other);
225
226 virtual ~GroupTopology() {}
227};
228
229/** @brief Communicator performing operations within groups defined by a
230 GroupTopology with arbitrary-size data associated with each group. */
232{
233public:
234 /// Communication mode.
235 enum Mode
236 {
237 byGroup, ///< Communications are performed one group at a time.
238 byNeighbor /**< Communications are performed one neighbor at a time,
239 aggregating over groups. */
240 };
241
242protected:
246 Table group_ltdof; // only for groups for which this processor is master.
249 MPI_Request *requests;
250 // MPI_Status *statuses;
251 // comm_lock: 0 - no lock, 1 - locked for Bcast, 2 - locked for Reduce
252 mutable int comm_lock;
253 mutable int num_requests;
255 int *buf_offsets; // size = max(number of groups, number of neighbors)
257
258public:
259 /// Construct a GroupCommunicator object.
260 /** The object must be initialized before it can be used to perform any
261 operations. To initialize the object, either
262 - call Create() or
263 - initialize the Table reference returned by GroupLDofTable() and then
264 call Finalize().
265 */
267
268 /** @brief Initialize the communicator from a local-dof to group map.
269 Finalize() is called internally. */
270 void Create(const Array<int> &ldof_group);
271
272 /** @brief Fill-in the returned Table reference to initialize the
273 GroupCommunicator then call Finalize(). */
275
276 /// Read-only access to group-ldof Table.
277 const Table &GroupLDofTable() const { return group_ldof; }
278
279 /// Allocate internal buffers after the GroupLDofTable is defined
280 void Finalize();
281
282 /// Initialize the internal group_ltdof Table.
283 /** This method must be called before performing operations that use local
284 data layout 2, see CopyGroupToBuffer() for layout descriptions. */
285 void SetLTDofTable(const Array<int> &ldof_ltdof);
286
287 /// Get a reference to the associated GroupTopology object
289
290 /// Get a const reference to the associated GroupTopology object
291 const GroupTopology &GetGroupTopology() const { return gtopo; }
292
293 /// Dofs to be sent to communication neighbors
294 void GetNeighborLTDofTable(Table &nbr_ltdof) const;
295
296 /// Dofs to be received from communication neighbors
297 void GetNeighborLDofTable(Table &nbr_ldof) const;
298
299 /** @brief Data structure on which we define reduce operations.
300 The data is associated with (and the operation is performed on) one
301 group at a time. */
302 template <class T> struct OpData
303 {
304 int nldofs, nb;
305 const int *ldofs;
306 T *ldata, *buf;
307 };
308
309 /** @brief Copy the entries corresponding to the group @a group from the
310 local array @a ldata to the buffer @a buf. */
311 /** The @a layout of the local array can be:
312 - 0 - @a ldata is an array on all ldofs: copied indices:
313 `{ J[j] : I[group] <= j < I[group+1] }` where `I,J=group_ldof.{I,J}`
314 - 1 - @a ldata is an array on the shared ldofs: copied indices:
315 `{ j : I[group] <= j < I[group+1] }` where `I,J=group_ldof.{I,J}`
316 - 2 - @a ldata is an array on the true ldofs, ltdofs: copied indices:
317 `{ J[j] : I[group] <= j < I[group+1] }` where `I,J=group_ltdof.{I,J}`.
318 @returns The pointer @a buf plus the number of elements in the group. */
319 template <class T>
320 T *CopyGroupToBuffer(const T *ldata, T *buf, int group, int layout) const;
321
322 /** @brief Copy the entries corresponding to the group @a group from the
323 buffer @a buf to the local array @a ldata. */
324 /** For a description of @a layout, see CopyGroupToBuffer().
325 @returns The pointer @a buf plus the number of elements in the group. */
326 template <class T>
327 const T *CopyGroupFromBuffer(const T *buf, T *ldata, int group,
328 int layout) const;
329
330 /** @brief Perform the reduction operation @a Op to the entries of group
331 @a group using the values from the buffer @a buf and the values from the
332 local array @a ldata, saving the result in the latter. */
333 /** For a description of @a layout, see CopyGroupToBuffer().
334 @returns The pointer @a buf plus the number of elements in the group. */
335 template <class T>
336 const T *ReduceGroupFromBuffer(const T *buf, T *ldata, int group,
337 int layout, void (*Op)(OpData<T>)) const;
338
339 /// Begin a broadcast within each group where the master is the root.
340 /** For a description of @a layout, see CopyGroupToBuffer(). */
341 template <class T> void BcastBegin(T *ldata, int layout) const;
342
343 /** @brief Finalize a broadcast started with BcastBegin().
344
345 The output data @a layout can be:
346 - 0 - @a ldata is an array on all ldofs; the input layout should be
347 either 0 or 2
348 - 1 - @a ldata is the same array as given to BcastBegin(); the input
349 layout should be 1.
350
351 For more details about @a layout, see CopyGroupToBuffer(). */
352 template <class T> void BcastEnd(T *ldata, int layout) const;
353
354 /** @brief Broadcast within each group where the master is the root.
355
356 The data @a layout can be either 0 or 1.
357
358 For a description of @a layout, see CopyGroupToBuffer(). */
359 template <class T> void Bcast(T *ldata, int layout) const
360 {
361 BcastBegin(ldata, layout);
362 BcastEnd(ldata, layout);
363 }
364
365 /// Broadcast within each group where the master is the root.
366 template <class T> void Bcast(T *ldata) const { Bcast<T>(ldata, 0); }
367 /// Broadcast within each group where the master is the root.
368 template <class T> void Bcast(Array<T> &ldata) const
369 { Bcast<T>((T *)ldata); }
370
371 /** @brief Begin reduction operation within each group where the master is
372 the root. */
373 /** The input data layout is an array on all ldofs, i.e. layout 0, see
374 CopyGroupToBuffer().
375
376 The reduce operation will be specified when calling ReduceEnd(). This
377 method is instantiated for int and double. */
378 template <class T> void ReduceBegin(const T *ldata) const;
379
380 /** @brief Finalize reduction operation started with ReduceBegin().
381
382 The output data @a layout can be either 0 or 2, see CopyGroupToBuffer().
383
384 The reduce operation is given by the third argument (see below for list
385 of the supported operations.) This method is instantiated for int and
386 double.
387
388 @note If the output data layout is 2, then the data from the @a ldata
389 array passed to this call is used in the reduction operation, instead of
390 the data from the @a ldata array passed to ReduceBegin(). Therefore, the
391 data for master-groups has to be identical in both arrays.
392 */
393 template <class T> void ReduceEnd(T *ldata, int layout,
394 void (*Op)(OpData<T>)) const;
395
396 /** @brief Reduce within each group where the master is the root.
397
398 The reduce operation is given by the second argument (see below for list
399 of the supported operations.) */
400 template <class T> void Reduce(T *ldata, void (*Op)(OpData<T>)) const
401 {
402 ReduceBegin(ldata);
403 ReduceEnd(ldata, 0, Op);
404 }
405
406 /// Reduce within each group where the master is the root.
407 template <class T> void Reduce(Array<T> &ldata, void (*Op)(OpData<T>)) const
408 { Reduce<T>((T *)ldata, Op); }
409
410 /// Reduce operation Sum, instantiated for int and double
411 template <class T> static void Sum(OpData<T>);
412 /// Reduce operation Min, instantiated for int and double
413 template <class T> static void Min(OpData<T>);
414 /// Reduce operation Max, instantiated for int and double
415 template <class T> static void Max(OpData<T>);
416 /// Reduce operation bitwise OR, instantiated for int only
417 template <class T> static void BitOR(OpData<T>);
418
419 /// Print information about the GroupCommunicator from all MPI ranks.
420 void PrintInfo(std::ostream &out = mfem::out) const;
421
422 /** @brief Destroy a GroupCommunicator object, deallocating internal data
423 structures and buffers. */
425};
426
428{
429 NEIGHBOR_ELEMENT_RANK_VM, ///< NeighborElementRankMessage
430 NEIGHBOR_ORDER_VM, ///< NeighborOrderMessage
431 NEIGHBOR_DEREFINEMENT_VM, ///< NeighborDerefinementMessage
432 NEIGHBOR_REFINEMENT_VM, ///< NeighborRefinementMessage
433 NEIGHBOR_PREFINEMENT_VM, ///< NeighborPRefinementMessage
434 NEIGHBOR_ROW_VM, ///< NeighborRowMessage
435 REBALANCE_VM, ///< RebalanceMessage
436 REBALANCE_DOF_VM ///< RebalanceDofMessage
438
439/// \brief Variable-length MPI message containing unspecific binary data.
440template<int Tag>
442{
443 std::string data;
444 MPI_Request send_request;
445
446 /** @brief Non-blocking send to processor 'rank'.
447 Returns immediately. Completion (as tested by MPI_Wait/Test) does not
448 mean the message was received -- it may be on its way or just buffered
449 locally. */
450 void Isend(int rank, MPI_Comm comm)
451 {
452 Encode(rank);
453 MPI_Isend((void*) data.data(), static_cast<int>(data.length()), MPI_BYTE, rank,
454 Tag, comm, &send_request);
455 }
456
457 /** @brief Non-blocking synchronous send to processor 'rank'.
458 Returns immediately. Completion (MPI_Wait/Test) means that the message
459 was received. */
460 void Issend(int rank, MPI_Comm comm)
461 {
462 Encode(rank);
463 MPI_Issend((void*) data.data(), static_cast<int>(data.length()), MPI_BYTE, rank,
464 Tag, comm, &send_request);
465 }
466
467 /// Helper to send all messages in a rank-to-message map container.
468 template<typename MapT>
469 static void IsendAll(MapT& rank_msg, MPI_Comm comm)
470 {
471 for (auto it = rank_msg.begin(); it != rank_msg.end(); ++it)
472 {
473 it->second.Isend(it->first, comm);
474 }
475 }
476
477 /// Helper to wait for all messages in a map container to be sent.
478 template<typename MapT>
479 static void WaitAllSent(MapT& rank_msg)
480 {
481 for (auto it = rank_msg.begin(); it != rank_msg.end(); ++it)
482 {
483 MPI_Wait(&it->second.send_request, MPI_STATUS_IGNORE);
484 it->second.Clear();
485 }
486 }
487
488 /** @brief Return true if all messages in the map container were sent,
489 otherwise return false, without waiting. */
490 template<typename MapT>
491 static bool TestAllSent(MapT& rank_msg)
492 {
493 for (auto it = rank_msg.begin(); it != rank_msg.end(); ++it)
494 {
495 VarMessage &msg = it->second;
496 if (msg.send_request != MPI_REQUEST_NULL)
497 {
498 int sent;
499 MPI_Test(&msg.send_request, &sent, MPI_STATUS_IGNORE);
500 if (!sent) { return false; }
501 msg.Clear();
502 }
503 }
504 return true;
505 }
506
507 /** @brief Blocking probe for incoming message of this type from any rank.
508 Returns the rank and message size. */
509 static void Probe(int &rank, int &size, MPI_Comm comm)
510 {
511 MPI_Status status;
512 MPI_Probe(MPI_ANY_SOURCE, Tag, comm, &status);
513 rank = status.MPI_SOURCE;
514 MPI_Get_count(&status, MPI_BYTE, &size);
515 }
516
517 /** @brief Non-blocking probe for incoming message of this type from any
518 rank. If there is an incoming message, returns true and sets 'rank' and
519 'size'. Otherwise returns false. */
520 static bool IProbe(int &rank, int &size, MPI_Comm comm)
521 {
522 int flag;
523 MPI_Status status;
524 MPI_Iprobe(MPI_ANY_SOURCE, Tag, comm, &flag, &status);
525 if (!flag) { return false; }
526
527 rank = status.MPI_SOURCE;
528 MPI_Get_count(&status, MPI_BYTE, &size);
529 return true;
530 }
531
532 /// Post-probe receive from processor 'rank' of message size 'size'.
533 void Recv(int rank, int size, MPI_Comm comm)
534 {
535 MFEM_ASSERT(size >= 0, "");
536 data.resize(size);
537 MPI_Status status;
538 MPI_Recv((void*) data.data(), size, MPI_BYTE, rank, Tag, comm, &status);
539#ifdef MFEM_DEBUG
540 int count;
541 MPI_Get_count(&status, MPI_BYTE, &count);
542 MFEM_VERIFY(count == size, "");
543#endif
544 Decode(rank);
545 }
546
547 /// Like Recv(), but throw away the message.
548 void RecvDrop(int rank, int size, MPI_Comm comm)
549 {
550 data.resize(size);
551 MPI_Status status;
552 MPI_Recv((void*) data.data(), size, MPI_BYTE, rank, Tag, comm, &status);
553 data.resize(0); // don't decode
554 }
555
556 /// Helper to receive all messages in a rank-to-message map container.
557 template<typename MapT>
558 static void RecvAll(MapT& rank_msg, MPI_Comm comm)
559 {
560 int recv_left = static_cast<int>(rank_msg.size());
561 while (recv_left > 0)
562 {
563 int rank, size;
564 Probe(rank, size, comm);
565 MFEM_ASSERT(rank_msg.find(rank) != rank_msg.end(), "Unexpected message"
566 " (tag " << Tag << ") from rank " << rank);
567 // NOTE: no guard against receiving two messages from the same rank
568 rank_msg[rank].Recv(rank, size, comm);
569 --recv_left;
570 }
571 }
572
573 VarMessage() : send_request(MPI_REQUEST_NULL) {}
574
575 /// Clear the message and associated request.
576 void Clear() { data.clear(); send_request = MPI_REQUEST_NULL; }
577
578 virtual ~VarMessage()
579 {
580 MFEM_ASSERT(send_request == MPI_REQUEST_NULL,
581 "WaitAllSent was not called after Isend");
582 }
583
584 VarMessage(const VarMessage &other)
585 : data(other.data), send_request(other.send_request)
586 {
587 MFEM_ASSERT(send_request == MPI_REQUEST_NULL,
588 "Cannot copy message with a pending send.");
589 }
590
591protected:
592 virtual void Encode(int rank) = 0;
593 virtual void Decode(int rank) = 0;
594};
595
596
597/// Helper struct to convert a C++ type to an MPI type
598template <typename Type> struct MPITypeMap;
599
600// Specializations of MPITypeMap; mpi_type initialized in communication.cpp:
601template<> struct MPITypeMap<int>
602{
603 static MFEM_EXPORT const MPI_Datatype mpi_type;
604};
605template<> struct MPITypeMap<double>
606{
607 static MFEM_EXPORT const MPI_Datatype mpi_type;
608};
609template<> struct MPITypeMap<float>
610{
611 static MFEM_EXPORT const MPI_Datatype mpi_type;
612};
613
614
615/** Reorder MPI ranks to follow the Z-curve within the physical machine topology
616 (provided that functions to query physical node coordinates are available).
617 Returns a new communicator with reordered ranks. */
618MPI_Comm ReorderRanksZCurve(MPI_Comm comm);
619
620
621} // namespace mfem
622
623#endif
624
625#endif
int Size() const
Return the logical size of the array.
Definition array.hpp:147
Communicator performing operations within groups defined by a GroupTopology with arbitrary-size data ...
const Table & GroupLDofTable() const
Read-only access to group-ldof Table.
GroupCommunicator(const GroupTopology &gt, Mode m=byNeighbor)
Construct a GroupCommunicator object.
void GetNeighborLTDofTable(Table &nbr_ltdof) const
Dofs to be sent to communication neighbors.
Table & GroupLDofTable()
Fill-in the returned Table reference to initialize the GroupCommunicator then call Finalize().
Mode
Communication mode.
@ byGroup
Communications are performed one group at a time.
void GetNeighborLDofTable(Table &nbr_ldof) const
Dofs to be received from communication neighbors.
void BcastBegin(T *ldata, int layout) const
Begin a broadcast within each group where the master is the root.
const T * ReduceGroupFromBuffer(const T *buf, T *ldata, int group, int layout, void(*Op)(OpData< T >)) const
Perform the reduction operation Op to the entries of group group using the values from the buffer buf...
void BcastEnd(T *ldata, int layout) const
Finalize a broadcast started with BcastBegin().
void ReduceBegin(const T *ldata) const
Begin reduction operation within each group where the master is the root.
void ReduceEnd(T *ldata, int layout, void(*Op)(OpData< T >)) const
Finalize reduction operation started with ReduceBegin().
void Reduce(Array< T > &ldata, void(*Op)(OpData< T >)) const
Reduce within each group where the master is the root.
const GroupTopology & GetGroupTopology() const
Get a const reference to the associated GroupTopology object.
const GroupTopology & gtopo
void Reduce(T *ldata, void(*Op)(OpData< T >)) const
Reduce within each group where the master is the root.
T * CopyGroupToBuffer(const T *ldata, T *buf, int group, int layout) const
Copy the entries corresponding to the group group from the local array ldata to the buffer buf.
const T * CopyGroupFromBuffer(const T *buf, T *ldata, int group, int layout) const
Copy the entries corresponding to the group group from the buffer buf to the local array ldata.
const GroupTopology & GetGroupTopology()
Get a reference to the associated GroupTopology object.
void Create(const Array< int > &ldof_group)
Initialize the communicator from a local-dof to group map. Finalize() is called internally.
~GroupCommunicator()
Destroy a GroupCommunicator object, deallocating internal data structures and buffers.
static void Sum(OpData< T >)
Reduce operation Sum, instantiated for int and double.
void Bcast(T *ldata, int layout) const
Broadcast within each group where the master is the root.
void SetLTDofTable(const Array< int > &ldof_ltdof)
Initialize the internal group_ltdof Table.
static void BitOR(OpData< T >)
Reduce operation bitwise OR, instantiated for int only.
static void Max(OpData< T >)
Reduce operation Max, instantiated for int and double.
static void Min(OpData< T >)
Reduce operation Min, instantiated for int and double.
void Bcast(Array< T > &ldata) const
Broadcast within each group where the master is the root.
void Bcast(T *ldata) const
Broadcast within each group where the master is the root.
void Finalize()
Allocate internal buffers after the GroupLDofTable is defined.
void PrintInfo(std::ostream &out=mfem::out) const
Print information about the GroupCommunicator from all MPI ranks.
int GetNeighborRank(int i) const
Return the MPI rank of neighbor 'i'.
void SetComm(MPI_Comm comm)
Set the MPI communicator to 'comm'.
int NRanks() const
Return the number of MPI ranks within this object's communicator.
bool IAmMaster(int g) const
Return true if I am master for group 'g'.
void Swap(GroupTopology &other)
Swap the internal data with another GroupTopology object.
void Save(std::ostream &out) const
Save the data in a stream.
const int * GetGroup(int g) const
Return a pointer to a list of neighbors for a given group. Neighbor 0 is the local processor.
int MyRank() const
Return the MPI rank within this object's communicator.
int GetGroupSize(int g) const
Get the number of processors in a group.
int GetGroupMaster(int g) const
Return the neighbor index of the group master for a given group. Neighbor 0 is the local processor.
void Load(std::istream &in)
Load the data from a stream.
GroupTopology()
Constructor with the MPI communicator = 0.
GroupTopology(MPI_Comm comm)
Constructor given the MPI communicator 'comm'.
int GetGroupMasterRank(int g) const
Return the rank of the group master for group 'g'.
void Create(ListOfIntegerSets &groups, int mpitag)
Set up the group topology given the list of sets of shared entities.
MPI_Comm GetComm() const
Return the MPI communicator.
void Copy(GroupTopology &copy) const
Copy the internal data to the external 'copy'.
int GetNumNeighbors() const
Return the number of neighbors including the local processor.
int NGroups() const
Return the number of groups.
int GetGroupMasterGroup(int g) const
Return the group number in the master for group 'g'.
List of integer sets.
Definition sets.hpp:51
A simple convenience class based on the Mpi singleton class above. Preserved for backward compatibili...
bool Root() const
Return true if WorldRank() == 0.
int WorldSize() const
Return MPI_COMM_WORLD's size.
MPI_Session(int &argc, char **&argv)
int WorldRank() const
Return MPI_COMM_WORLD's rank.
A simple singleton class that calls MPI_Init() at construction and MPI_Finalize() at destruction....
static MFEM_EXPORT int default_thread_required
Default level of thread support for MPI_Init_thread.
static bool IsFinalized()
Return true if MPI has been finalized.
static bool Root()
Return true if the rank in MPI_COMM_WORLD is zero.
static void Finalize()
Finalize MPI (if it has been initialized and not yet already finalized).
static void Init(int *argc=nullptr, char ***argv=nullptr, int required=default_thread_required, int *provided=nullptr)
Singleton creation with Mpi::Init().
static bool IsInitialized()
Return true if MPI has been initialized.
static int WorldRank()
Return the MPI rank in MPI_COMM_WORLD.
static int WorldSize()
Return the size of MPI_COMM_WORLD.
static void Init(int &argc, char **&argv, int required=default_thread_required, int *provided=nullptr)
Singleton creation with Mpi::Init(argc, argv).
int RowSize(int i) const
Definition table.hpp:116
void GetRow(int i, Array< int > &row) const
Return row i in array row (the Table must be finalized)
Definition table.cpp:259
int Size() const
Returns the number of TYPE I elements.
Definition table.hpp:100
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition globals.hpp:66
MPI_Comm ReorderRanksZCurve(MPI_Comm comm)
@ REBALANCE_DOF_VM
RebalanceDofMessage.
@ NEIGHBOR_REFINEMENT_VM
NeighborRefinementMessage.
@ NEIGHBOR_ROW_VM
NeighborRowMessage.
@ NEIGHBOR_PREFINEMENT_VM
NeighborPRefinementMessage.
@ NEIGHBOR_ELEMENT_RANK_VM
NeighborElementRankMessage.
@ NEIGHBOR_ORDER_VM
NeighborOrderMessage.
@ REBALANCE_VM
RebalanceMessage.
@ NEIGHBOR_DEREFINEMENT_VM
NeighborDerefinementMessage.
Data structure on which we define reduce operations. The data is associated with (and the operation i...
static MFEM_EXPORT const MPI_Datatype mpi_type
static MFEM_EXPORT const MPI_Datatype mpi_type
static MFEM_EXPORT const MPI_Datatype mpi_type
Helper struct to convert a C++ type to an MPI type.
Variable-length MPI message containing unspecific binary data.
void Issend(int rank, MPI_Comm comm)
Non-blocking synchronous send to processor 'rank'. Returns immediately. Completion (MPI_Wait/Test) me...
static void WaitAllSent(MapT &rank_msg)
Helper to wait for all messages in a map container to be sent.
virtual void Encode(int rank)=0
static bool TestAllSent(MapT &rank_msg)
Return true if all messages in the map container were sent, otherwise return false,...
void Isend(int rank, MPI_Comm comm)
Non-blocking send to processor 'rank'. Returns immediately. Completion (as tested by MPI_Wait/Test) d...
MPI_Request send_request
static void IsendAll(MapT &rank_msg, MPI_Comm comm)
Helper to send all messages in a rank-to-message map container.
virtual void Decode(int rank)=0
static void RecvAll(MapT &rank_msg, MPI_Comm comm)
Helper to receive all messages in a rank-to-message map container.
VarMessage(const VarMessage &other)
void Clear()
Clear the message and associated request.
void RecvDrop(int rank, int size, MPI_Comm comm)
Like Recv(), but throw away the message.
static bool IProbe(int &rank, int &size, MPI_Comm comm)
Non-blocking probe for incoming message of this type from any rank. If there is an incoming message,...
void Recv(int rank, int size, MPI_Comm comm)
Post-probe receive from processor 'rank' of message size 'size'.
static void Probe(int &rank, int &size, MPI_Comm comm)
Blocking probe for incoming message of this type from any rank. Returns the rank and message size.