29#ifdef MFEM_USE_STRUMPACK
30#include <StrumpackConfig.hpp>
41#if defined(MFEM_USE_STRUMPACK) && \
42 (defined(STRUMPACK_USE_PTSCOTCH) || defined(STRUMPACK_USE_SLATE_SCALAPACK))
51 group_lproc(gt.group_lproc)
53 gt.groupmaster_lproc.
Copy(groupmaster_lproc);
54 gt.lproc_proc.
Copy(lproc_proc);
55 gt.group_mgroup.
Copy(group_mgroup);
58void GroupTopology::ProcToLProc()
61 MPI_Comm_size(MyComm, &
NRanks);
63 map<int, int> proc_lproc;
68 int lproc_counter = 0;
71 const pair<const int, int>
p(group_lproc.
GetJ()[i], lproc_counter);
72 if (proc_lproc.insert(
p).second)
79 lproc_proc.
SetSize(lproc_counter);
80 for (map<int, int>::iterator it = proc_lproc.begin();
81 it != proc_lproc.end(); ++it)
83 lproc_proc[it->second] = it->first;
88 group_lproc.
GetJ()[i] = proc_lproc[group_lproc.
GetJ()[i]];
91 for (
int i = 0; i <
NGroups(); i++)
93 groupmaster_lproc[i] = proc_lproc[groupmaster_lproc[i]];
101 Table group_mgroupandproc;
104 for (
int i = 0; i <
NGroups(); i++)
106 int j = group_mgroupandproc.
GetI()[i];
107 group_mgroupandproc.
GetI()[i+1] = j + group_lproc.
RowSize(i) + 1;
108 group_mgroupandproc.
GetJ()[j] = i;
110 for (
int k = group_lproc.
GetI()[i];
111 j < group_mgroupandproc.
GetI()[i+1]; j++, k++)
113 group_mgroupandproc.
GetJ()[j] = group_lproc.
GetJ()[k];
121 for (
int i = 0; i <
NGroups(); i++)
139 MFEM_DEBUG_DO(group_mgroup = -1);
140 for (
int g = 0; g <
NGroups(); g++)
142 if (
IAmMaster(g)) { group_mgroup[g] = g; }
150 for (
int g = 1; g <
NGroups(); g++)
156 for (
int i = 0; i < gs; i++)
169 lproc_cgroup_list.
Sort();
170 lproc_cgroup_list.
Unique();
182 for (
int nbr = 1; nbr < lproc_cgroup.
Size(); nbr++)
184 const int send_row = 2*(nbr-1);
185 const int recv_row = send_row+1;
186 const int ng = lproc_cgroup.
RowSize(nbr);
187 const int *g = lproc_cgroup.
GetRow(nbr);
188 for (
int j = 0; j < ng; j++)
190 const int gs = group_lproc.
RowSize(g[j]);
203 for (
int nbr = 1; nbr < lproc_cgroup.
Size(); nbr++)
205 const int send_row = 2*(nbr-1);
206 const int recv_row = send_row+1;
207 const int ng = lproc_cgroup.
RowSize(nbr);
208 const int *g = lproc_cgroup.
GetRow(nbr);
209 for (
int j = 0; j < ng; j++)
211 const int gs = group_lproc.
RowSize(g[j]);
216 send_row, group_mgroupandproc.
GetRow(g[j]), gs+1);
227 send_requests = MPI_REQUEST_NULL;
228 recv_requests = MPI_REQUEST_NULL;
229 for (
int nbr = 1; nbr < lproc_cgroup.
Size(); nbr++)
231 const int send_row = 2*(nbr-1);
232 const int recv_row = send_row+1;
233 const int send_size = buffer.
RowSize(send_row);
234 const int recv_size = buffer.
RowSize(recv_row);
237 MPI_Isend(buffer.
GetRow(send_row), send_size, MPI_INT, lproc_proc[nbr],
238 mpitag, MyComm, &send_requests[nbr-1]);
242 MPI_Irecv(buffer.
GetRow(recv_row), recv_size, MPI_INT, lproc_proc[nbr],
243 mpitag, MyComm, &recv_requests[nbr-1]);
247 if (recv_requests.
Size() > 0)
251 while (MPI_Waitany(recv_requests.
Size(), recv_requests.
GetData(), &idx,
253 idx != MPI_UNDEFINED)
255 const int recv_size = buffer.
RowSize(2*idx+1);
256 const int *recv_buf = buffer.
GetRow(2*idx+1);
257 for (
int s = 0; s < recv_size; s += recv_buf[s]+2)
259 group.
Recreate(recv_buf[s], recv_buf+s+2);
260 const int g = groups.
Lookup(group);
261 MFEM_ASSERT(group_mgroup[g] == -1,
"communication error");
262 group_mgroup[g] = recv_buf[s+1];
267 MPI_Waitall(send_requests.
Size(), send_requests.
GetData(),
268 MPI_STATUSES_IGNORE);
275 os <<
"\ncommunication_groups\n";
276 os <<
"number_of_groups " <<
NGroups() <<
"\n\n";
278 os <<
"# number of entities in each group, followed by ranks in group\n";
279 for (
int group_id = 0; group_id <
NGroups(); ++group_id)
282 const int * group_ptr =
GetGroup(group_id);
284 for (
int group_member_index = 0; group_member_index < group_size;
285 ++group_member_index)
303 int number_of_groups = -1;
305 MFEM_VERIFY(ident ==
"number_of_groups",
306 "GroupTopology::Load - expected 'number_of_groups' entry.");
307 in >> number_of_groups;
313 for (
int group_id = 0; group_id < number_of_groups; ++group_id)
319 array.Reserve(group_size);
326 integer_sets.
Insert(integer_set);
329 Create(integer_sets, 823);
334 copy.SetComm(MyComm);
335 group_lproc.
Copy(copy.group_lproc);
336 groupmaster_lproc.
Copy(copy.groupmaster_lproc);
337 lproc_proc.
Copy(copy.lproc_proc);
338 group_mgroup.
Copy(copy.group_mgroup);
345 mfem::Swap(groupmaster_lproc, other.groupmaster_lproc);
363 MPI_UNSIGNED_LONG_LONG;
383 for (
int i = 0; i < ldof_group.
Size(); i++)
385 int group = ldof_group[i];
393 for (
int i = 0; i < ldof_group.
Size(); i++)
395 int group = ldof_group[i];
408 int request_counter = 0;
427 request_counter += gr_requests;
432 requests =
new MPI_Request[request_counter];
442 if (nldofs == 0) {
continue; }
452 for (
int i = 0; i < grp_size; i++)
454 if (grp_nbr_list[i] != 0)
466 if (nldofs == 0) {
continue; }
476 for (
int i = 0; i < grp_size; i++)
478 if (grp_nbr_list[i] != 0)
497 if (num_recv_groups > 0)
500 group_ids.
SetSize(num_recv_groups);
501 for (
int i = 0; i < num_recv_groups; i++)
504 group_ids[i].two = grp_list[i];
507 for (
int i = 0; i < num_recv_groups; i++)
509 grp_list[i] = group_ids[i].two;
534 for (
int i = 0; i < nldofs; i++)
549 if (num_send_groups > 0)
552 for (
int i = 0; i < num_send_groups; i++)
554 const int group = grp_list[i];
564 if (num_send_groups > 0)
567 for (
int i = 0; i < num_send_groups; i++)
569 const int group = grp_list[i];
585 if (num_recv_groups > 0)
588 for (
int i = 0; i < num_recv_groups; i++)
590 const int group = grp_list[i];
600 if (num_recv_groups > 0)
603 for (
int i = 0; i < num_recv_groups; i++)
605 const int group = grp_list[i];
631 for (
int j = 0; j < nltdofs; j++)
633 buf[j] = ldata[ltdofs[j]];
635 return buf + nltdofs;
641 for (
int j = 0; j < nldofs; j++)
643 buf[j] = ldata[ldofs[j]];
652 int group,
int layout)
const
665 for (
int j = 0; j < nldofs; j++)
667 ldata[ltdofs[j]] = buf[j];
674 for (
int j = 0; j < nldofs; j++)
676 ldata[ldofs[j]] = buf[j];
686 int group,
int layout,
693 opd.
buf =
const_cast<T*
>(buf);
699 MFEM_ABORT(
"layout 1 is not supported");
701 for (
int j = 0; j < opd.
nldofs; j++)
726 MFEM_VERIFY(
comm_lock == 0,
"object is already in use");
730 int request_counter = 0;
741 "'group_ltdof' is not set, use SetLTDofTable()");
753 if (nldofs == 0) {
continue; }
775 for (
int i = 0; i < gs; i++)
803 if (num_send_groups > 0)
810 for (
int i = 0; i < num_send_groups; i++)
826 if (num_recv_groups > 0)
834 for (
int i = 0; i < num_recv_groups; i++)
865 MFEM_VERIFY(
comm_lock == 1,
"object is NOT locked for Bcast");
875 else if (layout == 0)
880 idx != MPI_UNDEFINED)
883 if (gr == -1) {
continue; }
898 idx != MPI_UNDEFINED)
901 if (nbr == -1) {
continue; }
904 if (num_recv_groups > 0)
908 for (
int i = 0; i < num_recv_groups; i++)
925 MFEM_VERIFY(
comm_lock == 0,
"object is already in use");
929 int request_counter = 0;
940 if (nldofs == 0) {
continue; }
944 const int layout = 0;
962 for (
int i = 0; i < gs; i++)
989 if (num_send_groups > 0)
993 for (
int i = 0; i < num_send_groups; i++)
995 const int layout = 0;
1011 if (num_recv_groups > 0)
1015 for (
int i = 0; i < num_recv_groups; i++)
1047 MFEM_VERIFY(
comm_lock == 2,
"object is NOT locked for Reduce");
1063 idx != MPI_UNDEFINED)
1066 if (gr == -1) {
continue; }
1070 if ((--group_num_req[gr]) != 0) {
continue; }
1076 opd.
ldofs = (layout == 0) ?
1092 if (num_recv_groups > 0)
1096 for (
int i = 0; i < num_recv_groups; i++)
1116 for (
int i = 0; i < opd.
nldofs; i++)
1123 for (
int i = 0; i < opd.
nldofs; i++)
1126 for (
int j = 0; j < opd.
nb; j++)
1138 for (
int i = 0; i < opd.
nldofs; i++)
1141 for (
int j = 0; j < opd.
nb; j++)
1156 for (
int i = 0; i < opd.
nldofs; i++)
1159 for (
int j = 0; j < opd.
nb; j++)
1174 for (
int i = 0; i < opd.
nldofs; i++)
1177 for (
int j = 0; j < opd.
nb; j++)
1188 const int tag = 46800;
1191 int num_sends = 0, num_recvs = 0;
1192 size_t mem_sends = 0, mem_recvs = 0;
1193 int num_master_groups = 0, num_empty_groups = 0;
1194 int num_active_neighbors = 0;
1210 num_master_groups++;
1215 mem_recvs +=
sizeof(double)*nldofs;
1231 num_master_groups++;
1237 if (num_send_groups > 0)
1240 for (
int i = 0; i < num_send_groups; i++)
1248 if (num_recv_groups > 0)
1251 for (
int i = 0; i < num_recv_groups; i++)
1257 if (num_send_groups > 0 || num_recv_groups > 0)
1259 num_active_neighbors++;
1266 MPI_Recv(&c, 1, MPI_CHAR, myid-1, tag,
gtopo.
GetComm(),
1271 os <<
"\nGroupCommunicator:\n";
1273 os <<
"Rank " << myid <<
":\n"
1275 (
mode ==
byGroup ?
"byGroup" :
"byNeighbor") <<
"\n"
1276 " number of sends = " << num_sends <<
1277 " (" << mem_sends <<
" bytes)\n"
1278 " number of recvs = " << num_recvs <<
1279 " (" << mem_recvs <<
" bytes)\n";
1282 num_master_groups <<
" + " <<
1284 num_empty_groups <<
" (master + slave + empty)\n";
1289 num_active_neighbors <<
" + " <<
1291 " (active + inactive)\n";
1296 MPI_Send(&c, 1, MPI_CHAR, myid+1, tag,
gtopo.
GetComm());
1320 int *,
int,
void (*)(OpData<int>))
const;
1326 double *,
int,
void (*)(OpData<double>))
const;
1332 float *,
int,
void (*)(OpData<float>))
const;
1352static void DebugRankCoords(
int** coords,
int dim,
int size)
1354 for (
int i = 0; i < size; i++)
1356 mfem::out <<
"Rank " << i <<
" coords: ";
1357 for (
int j = 0; j <
dim; j++)
1367 CompareCoords(
int coord) : coord(coord) {}
1370 bool operator()(
int*
const &
a,
int*
const &
b)
const
1371 {
return a[coord] <
b[coord]; }
1378 bool all_same =
true;
1379 for (
int i = 1; i < size && all_same; i++)
1381 for (
int j = 0; j <
dim; j++)
1383 if (coords[i][j] != coords[0][j]) { all_same =
false;
break; }
1386 if (all_same) {
return; }
1389 std::sort(coords, coords + size, CompareCoords(d));
1390 int next = (d + 1) %
dim;
1392 if (coords[0][d] < coords[size-1][d])
1410 MPI_Comm_rank(comm, &rank);
1411 MPI_Comm_size(comm, &size);
1414 MPIX_Torus_ndims(&
dim);
1416 int* mycoords =
new int[
dim + 1];
1417 MPIX_Rank2torus(rank, mycoords);
1419 MPI_Send(mycoords,
dim, MPI_INT, 0, 111, comm);
1424 int** coords =
new int*[size];
1425 for (
int i = 0; i < size; i++)
1427 coords[i] =
new int[
dim + 1];
1429 MPI_Recv(coords[i],
dim, MPI_INT, i, 111, comm, &status);
1436 for (
int i = 0; i < size; i++)
1438 MPI_Send(&coords[i][
dim], 1, MPI_INT, i, 112, comm);
1439 delete [] coords[i];
1445 MPI_Recv(&new_rank, 1, MPI_INT, 0, 112, comm, &status);
1448 MPI_Comm_split(comm, 0, new_rank, &new_comm);
void Sort()
Sorts the array in ascending order. This requires operator< to be defined for T.
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
int Size() const
Return the logical size of the array.
int Append(const T &el)
Append element 'el' to array, resize if necessary.
T * GetData()
Returns the data.
void Copy(Array ©) const
Create a copy of the internal array to the provided copy.
void Unique()
Removes duplicities from a sorted array. This requires operator== to be defined for T.
GroupCommunicator(const GroupTopology >, Mode m=byNeighbor)
Construct a GroupCommunicator object.
void GetNeighborLTDofTable(Table &nbr_ltdof) const
Dofs to be sent to communication neighbors.
@ 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().
const GroupTopology & gtopo
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.
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 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 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'.
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.
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 ©) 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'.
void Recreate(const int n, const int *p)
Create an integer set from C-array 'p' of 'n' integers. Overwrites any existing set data.
int Insert(const IntegerSet &s)
Check to see if set 's' is in the list. If not append it to the end of the list. Returns the index of...
void AsTable(Table &t) const
Write the list of sets into table 't'.
int PickElementInSet(int i) const
Return the value of the first element of the ith set.
int Lookup(const IntegerSet &s) const
static MFEM_EXPORT int default_thread_required
Default level of thread support for MPI_Init_thread.
Table stores the connectivity of elements of TYPE I to elements of TYPE II. For example,...
void AddConnections(int r, const int *c, int nc)
void GetRow(int i, Array< int > &row) const
Return row i in array row (the Table must be finalized)
void AddConnection(int r, int c)
int Size() const
Returns the number of TYPE I elements.
int Size_of_connections() const
Returns the number of connections in the table.
void AddColumnsInRow(int r, int ncol)
void MakeFromList(int nrows, const Array< Connection > &list)
Create the table from a list of connections {(from, to)}, where 'from' is a TYPE I index and 'to' is ...
void Copy(Table ©) const
void AddAColumnInRow(int r)
void SetDims(int rows, int nnz)
Set the rows and the number of all connections for the table.
int index(int i, int j, int nx, int ny)
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
void Swap(T &a, T &b)
Swap objects of type T. The operation is performed using the most specialized swap function from the ...
void KdTreeSort(int **coords, int d, int dim, int size)
MPI_Comm ReorderRanksZCurve(MPI_Comm comm)
void skip_comment_lines(std::istream &is, const char comment_char)
Check if the stream starts with comment_char. If so skip it.
real_t p(const Vector &x, real_t t)
Helper struct for defining a connectivity table, see Table::MakeFromList.
Data structure on which we define reduce operations. The data is associated with (and the operation i...
Helper struct to convert a C++ type to an MPI type.