MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
communication.cpp
Go to the documentation of this file.
1// Copyright (c) 2010-2024, 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#include "../config/config.hpp"
13
14#ifdef MFEM_USE_MPI
15
16#include <mpi.h>
17#ifdef __bgq__
18#include <mpix.h>
19#endif
20
21#include "array.hpp"
22#include "table.hpp"
23#include "sets.hpp"
24#include "communication.hpp"
25#include "text.hpp"
26#include "sort_pairs.hpp"
27#include "globals.hpp"
28
29#ifdef MFEM_USE_STRUMPACK
30#include <StrumpackConfig.hpp> // STRUMPACK_USE_PTSCOTCH, etc.
31#endif
32
33#include <iostream>
34#include <map>
35
36using namespace std;
37
38namespace mfem
39{
40
41#if defined(MFEM_USE_STRUMPACK) && \
42 (defined(STRUMPACK_USE_PTSCOTCH) || defined(STRUMPACK_USE_SLATE_SCALAPACK))
43int Mpi::default_thread_required = MPI_THREAD_MULTIPLE;
44#else
45int Mpi::default_thread_required = MPI_THREAD_SINGLE;
46#endif
47
48
50 : MyComm(gt.MyComm),
51 group_lproc(gt.group_lproc)
52{
53 gt.groupmaster_lproc.Copy(groupmaster_lproc);
54 gt.lproc_proc.Copy(lproc_proc);
55 gt.group_mgroup.Copy(group_mgroup);
56}
57
58void GroupTopology::ProcToLProc()
59{
60 int NRanks;
61 MPI_Comm_size(MyComm, &NRanks);
62
63 map<int, int> proc_lproc;
64
65 // The local processor ids are assigned following the group order and within
66 // a group following their ordering in the group. In other words, the ids are
67 // assigned based on their order in the J array of group_lproc.
68 int lproc_counter = 0;
69 for (int i = 0; i < group_lproc.Size_of_connections(); i++)
70 {
71 const pair<const int, int> p(group_lproc.GetJ()[i], lproc_counter);
72 if (proc_lproc.insert(p).second)
73 {
74 lproc_counter++;
75 }
76 }
77 // Note: group_lproc.GetJ()[0] == MyRank --> proc_lproc[MyRank] == 0
78
79 lproc_proc.SetSize(lproc_counter);
80 for (map<int, int>::iterator it = proc_lproc.begin();
81 it != proc_lproc.end(); ++it)
82 {
83 lproc_proc[it->second] = it->first;
84 }
85
86 for (int i = 0; i < group_lproc.Size_of_connections(); i++)
87 {
88 group_lproc.GetJ()[i] = proc_lproc[group_lproc.GetJ()[i]];
89 }
90
91 for (int i = 0; i < NGroups(); i++)
92 {
93 groupmaster_lproc[i] = proc_lproc[groupmaster_lproc[i]];
94 }
95}
96
98{
99 groups.AsTable(group_lproc); // group_lproc = group_proc
100
101 Table group_mgroupandproc;
102 group_mgroupandproc.SetDims(NGroups(),
103 group_lproc.Size_of_connections() + NGroups());
104 for (int i = 0; i < NGroups(); i++)
105 {
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;
109 j++;
110 for (int k = group_lproc.GetI()[i];
111 j < group_mgroupandproc.GetI()[i+1]; j++, k++)
112 {
113 group_mgroupandproc.GetJ()[j] = group_lproc.GetJ()[k];
114 }
115 }
116
117 // build groupmaster_lproc with lproc = proc
118 groupmaster_lproc.SetSize(NGroups());
119
120 // simplest choice of the group owner
121 for (int i = 0; i < NGroups(); i++)
122 {
123 groupmaster_lproc[i] = groups.PickElementInSet(i);
124 }
125
126 // load-balanced choice of the group owner, which however can lead to
127 // isolated dofs
128 // for (i = 0; i < NGroups(); i++)
129 // groupmaster_lproc[i] = groups.PickRandomElementInSet(i);
130
131 ProcToLProc();
132
133 // Build 'group_mgroup':
134
135 // Use aggregated neighbor communication: at most one send to and/or one
136 // receive from each neighbor.
137
138 group_mgroup.SetSize(NGroups());
139 MFEM_DEBUG_DO(group_mgroup = -1);
140 for (int g = 0; g < NGroups(); g++)
141 {
142 if (IAmMaster(g)) { group_mgroup[g] = g; }
143 }
144
145 // The Table 'lproc_cgroup': for each lproc, list the groups that are owned
146 // by this rank or by that lproc.
147 Table lproc_cgroup;
148 {
149 Array<Connection> lproc_cgroup_list;
150 for (int g = 1; g < NGroups(); g++)
151 {
152 if (IAmMaster(g))
153 {
154 const int gs = GetGroupSize(g);
155 const int *lprocs = GetGroup(g);
156 for (int i = 0; i < gs; i++)
157 {
158 if (lprocs[i])
159 {
160 lproc_cgroup_list.Append(Connection(lprocs[i],g));
161 }
162 }
163 }
164 else
165 {
166 lproc_cgroup_list.Append(Connection(GetGroupMaster(g),g));
167 }
168 }
169 lproc_cgroup_list.Sort();
170 lproc_cgroup_list.Unique();
171 lproc_cgroup.MakeFromList(GetNumNeighbors(), lproc_cgroup_list);
172 }
173
174 // Determine size of the send-receive buffer. For each neighbor the buffer
175 // contains: <send-part><receive-part> with each part consisting of a list of
176 // groups. Each group, g, has group_lproc.RowSize(g)+2 integers: the first
177 // entry is group_lproc.RowSize(g) - the number of processors in the group,
178 // followed by the group-id in the master processor, followed by the ranks of
179 // the processors in the group.
180 Table buffer;
181 buffer.MakeI(2*lproc_cgroup.Size()-2); // excluding the "local" lproc, 0
182 for (int nbr = 1; nbr < lproc_cgroup.Size(); nbr++)
183 {
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++)
189 {
190 const int gs = group_lproc.RowSize(g[j]);
191 if (IAmMaster(g[j]))
192 {
193 buffer.AddColumnsInRow(send_row, gs+2);
194 }
195 else
196 {
197 MFEM_ASSERT(GetGroupMaster(g[j]) == nbr, "internal error");
198 buffer.AddColumnsInRow(recv_row, gs+2);
199 }
200 }
201 }
202 buffer.MakeJ();
203 for (int nbr = 1; nbr < lproc_cgroup.Size(); nbr++)
204 {
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++)
210 {
211 const int gs = group_lproc.RowSize(g[j]);
212 if (IAmMaster(g[j]))
213 {
214 buffer.AddConnection(send_row, gs);
215 buffer.AddConnections(
216 send_row, group_mgroupandproc.GetRow(g[j]), gs+1);
217 }
218 else
219 {
220 buffer.AddColumnsInRow(recv_row, gs+2);
221 }
222 }
223 }
224 buffer.ShiftUpI();
225 Array<MPI_Request> send_requests(lproc_cgroup.Size()-1);
226 Array<MPI_Request> recv_requests(lproc_cgroup.Size()-1);
227 send_requests = MPI_REQUEST_NULL;
228 recv_requests = MPI_REQUEST_NULL;
229 for (int nbr = 1; nbr < lproc_cgroup.Size(); nbr++)
230 {
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);
235 if (send_size > 0)
236 {
237 MPI_Isend(buffer.GetRow(send_row), send_size, MPI_INT, lproc_proc[nbr],
238 mpitag, MyComm, &send_requests[nbr-1]);
239 }
240 if (recv_size > 0)
241 {
242 MPI_Irecv(buffer.GetRow(recv_row), recv_size, MPI_INT, lproc_proc[nbr],
243 mpitag, MyComm, &recv_requests[nbr-1]);
244 }
245 }
246
247 if (recv_requests.Size() > 0)
248 {
249 int idx;
250 IntegerSet group;
251 while (MPI_Waitany(recv_requests.Size(), recv_requests.GetData(), &idx,
252 MPI_STATUS_IGNORE),
253 idx != MPI_UNDEFINED)
254 {
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)
258 {
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];
263 }
264 }
265 }
266
267 MPI_Waitall(send_requests.Size(), send_requests.GetData(),
268 MPI_STATUSES_IGNORE);
269
270 // debug barrier: MPI_Barrier(MyComm);
271}
272
273void GroupTopology::Save(ostream &os) const
274{
275 os << "\ncommunication_groups\n";
276 os << "number_of_groups " << NGroups() << "\n\n";
277
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)
280 {
281 int group_size = GetGroupSize(group_id);
282 const int * group_ptr = GetGroup(group_id);
283 os << group_size;
284 for ( int group_member_index = 0; group_member_index < group_size;
285 ++group_member_index)
286 {
287 os << " " << GetNeighborRank( group_ptr[group_member_index] );
288 }
289 os << "\n";
290 }
291
292 // For future use, optional ownership strategy.
293 // os << "# ownership";
294}
295
296void GroupTopology::Load(istream &in)
297{
298 // Load in group topology and create list of integer sets. Use constructor
299 // that uses list of integer sets.
300 std::string ident;
301
302 // Read in number of groups
303 int number_of_groups = -1;
304 in >> ident;
305 MFEM_VERIFY(ident == "number_of_groups",
306 "GroupTopology::Load - expected 'number_of_groups' entry.");
307 in >> number_of_groups;
308
309 // Skip number of entries in each group comment.
310 skip_comment_lines(in, '#');
311
312 ListOfIntegerSets integer_sets;
313 for (int group_id = 0; group_id < number_of_groups; ++group_id)
314 {
315 IntegerSet integer_set;
316 Array<int>& array = integer_set;
317 int group_size;
318 in >> group_size;
319 array.Reserve(group_size);
320 for ( int index = 0; index < group_size; ++index )
321 {
322 int value;
323 in >> value;
324 array.Append(value);
325 }
326 integer_sets.Insert(integer_set);
327 }
328
329 Create(integer_sets, 823);
330}
331
333{
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);
339}
340
342{
343 mfem::Swap(MyComm, other.MyComm);
344 mfem::Swap(group_lproc, other.group_lproc);
345 mfem::Swap(groupmaster_lproc, other.groupmaster_lproc);
346 mfem::Swap(lproc_proc, other.lproc_proc);
347 mfem::Swap(group_mgroup, other.group_mgroup);
348}
349
350// Initialize the static mpi_type for the specializations of MPITypeMap:
351const MPI_Datatype MPITypeMap<int>::mpi_type = MPI_INT;
352const MPI_Datatype MPITypeMap<float>::mpi_type = MPI_FLOAT;
353const MPI_Datatype MPITypeMap<double>::mpi_type = MPI_DOUBLE;
354
355
357 : gtopo(gt), mode(m)
358{
359 group_buf_size = 0;
360 requests = NULL;
361 // statuses = NULL;
362 comm_lock = 0;
363 num_requests = 0;
364 request_marker = NULL;
365 buf_offsets = NULL;
366}
367
369{
371 for (int i = 0; i < ldof_group.Size(); i++)
372 {
373 int group = ldof_group[i];
374 if (group != 0)
375 {
377 }
378 }
380
381 for (int i = 0; i < ldof_group.Size(); i++)
382 {
383 int group = ldof_group[i];
384 if (group != 0)
385 {
386 group_ldof.AddConnection(group, i);
387 }
388 }
390
391 Finalize();
392}
393
395{
396 int request_counter = 0;
397
398 // size buf_offsets = max(number of groups, number of neighbors)
399 buf_offsets = new int[max(group_ldof.Size(), gtopo.GetNumNeighbors())];
400 buf_offsets[0] = 0;
401 for (int gr = 1; gr < group_ldof.Size(); gr++)
402 {
403 if (group_ldof.RowSize(gr) != 0)
404 {
405 int gr_requests;
406 if (!gtopo.IAmMaster(gr)) // we are not the master
407 {
408 gr_requests = 1;
409 }
410 else
411 {
412 gr_requests = gtopo.GetGroupSize(gr)-1;
413 }
414
415 request_counter += gr_requests;
416 group_buf_size += gr_requests * group_ldof.RowSize(gr);
417 }
418 }
419
420 requests = new MPI_Request[request_counter];
421 // statuses = new MPI_Status[request_counter];
422 request_marker = new int[request_counter];
423
424 // Construct nbr_send_groups and nbr_recv_groups: (nbr 0 = me)
427 for (int gr = 1; gr < group_ldof.Size(); gr++)
428 {
429 const int nldofs = group_ldof.RowSize(gr);
430 if (nldofs == 0) { continue; }
431
432 if (!gtopo.IAmMaster(gr)) // we are not the master
433 {
435 }
436 else // we are the master
437 {
438 const int grp_size = gtopo.GetGroupSize(gr);
439 const int *grp_nbr_list = gtopo.GetGroup(gr);
440 for (int i = 0; i < grp_size; i++)
441 {
442 if (grp_nbr_list[i] != 0)
443 {
444 nbr_send_groups.AddAColumnInRow(grp_nbr_list[i]);
445 }
446 }
447 }
448 }
451 for (int gr = 1; gr < group_ldof.Size(); gr++)
452 {
453 const int nldofs = group_ldof.RowSize(gr);
454 if (nldofs == 0) { continue; }
455
456 if (!gtopo.IAmMaster(gr)) // we are not the master
457 {
459 }
460 else // we are the master
461 {
462 const int grp_size = gtopo.GetGroupSize(gr);
463 const int *grp_nbr_list = gtopo.GetGroup(gr);
464 for (int i = 0; i < grp_size; i++)
465 {
466 if (grp_nbr_list[i] != 0)
467 {
468 nbr_send_groups.AddConnection(grp_nbr_list[i], gr);
469 }
470 }
471 }
472 }
475 // The above construction creates the Tables with the column indices
476 // sorted, i.e. the group lists are sorted. To coordinate this order between
477 // processors, we will sort the group lists in the nbr_recv_groups Table
478 // according to their indices in the master. This does not require any
479 // communication because we have access to the group indices in the master
480 // by calling: master_group_id = gtopo.GetGroupMasterGroup(my_group_id).
481 Array<Pair<int,int> > group_ids;
482 for (int nbr = 1; nbr < nbr_recv_groups.Size(); nbr++)
483 {
484 const int num_recv_groups = nbr_recv_groups.RowSize(nbr);
485 if (num_recv_groups > 0)
486 {
487 int *grp_list = nbr_recv_groups.GetRow(nbr);
488 group_ids.SetSize(num_recv_groups);
489 for (int i = 0; i < num_recv_groups; i++)
490 {
491 group_ids[i].one = gtopo.GetGroupMasterGroup(grp_list[i]);
492 group_ids[i].two = grp_list[i]; // my_group_id
493 }
494 group_ids.Sort();
495 for (int i = 0; i < num_recv_groups; i++)
496 {
497 grp_list[i] = group_ids[i].two;
498 }
499 }
500 }
501}
502
504{
505 if (group_ltdof.Size() == group_ldof.Size()) { return; }
506
508 for (int gr = 1; gr < group_ldof.Size(); gr++)
509 {
510 if (gtopo.IAmMaster(gr))
511 {
513 }
514 }
516 for (int gr = 1; gr < group_ldof.Size(); gr++)
517 {
518 if (gtopo.IAmMaster(gr))
519 {
520 const int *ldofs = group_ldof.GetRow(gr);
521 const int nldofs = group_ldof.RowSize(gr);
522 for (int i = 0; i < nldofs; i++)
523 {
524 group_ltdof.AddConnection(gr, ldof_ltdof[ldofs[i]]);
525 }
526 }
527 }
529}
530
532{
533 nbr_ltdof.MakeI(nbr_send_groups.Size());
534 for (int nbr = 1; nbr < nbr_send_groups.Size(); nbr++)
535 {
536 const int num_send_groups = nbr_send_groups.RowSize(nbr);
537 if (num_send_groups > 0)
538 {
539 const int *grp_list = nbr_send_groups.GetRow(nbr);
540 for (int i = 0; i < num_send_groups; i++)
541 {
542 const int group = grp_list[i];
543 const int nltdofs = group_ltdof.RowSize(group);
544 nbr_ltdof.AddColumnsInRow(nbr, nltdofs);
545 }
546 }
547 }
548 nbr_ltdof.MakeJ();
549 for (int nbr = 1; nbr < nbr_send_groups.Size(); nbr++)
550 {
551 const int num_send_groups = nbr_send_groups.RowSize(nbr);
552 if (num_send_groups > 0)
553 {
554 const int *grp_list = nbr_send_groups.GetRow(nbr);
555 for (int i = 0; i < num_send_groups; i++)
556 {
557 const int group = grp_list[i];
558 const int nltdofs = group_ltdof.RowSize(group);
559 const int *ltdofs = group_ltdof.GetRow(group);
560 nbr_ltdof.AddConnections(nbr, ltdofs, nltdofs);
561 }
562 }
563 }
564 nbr_ltdof.ShiftUpI();
565}
566
568{
569 nbr_ldof.MakeI(nbr_recv_groups.Size());
570 for (int nbr = 1; nbr < nbr_recv_groups.Size(); nbr++)
571 {
572 const int num_recv_groups = nbr_recv_groups.RowSize(nbr);
573 if (num_recv_groups > 0)
574 {
575 const int *grp_list = nbr_recv_groups.GetRow(nbr);
576 for (int i = 0; i < num_recv_groups; i++)
577 {
578 const int group = grp_list[i];
579 const int nldofs = group_ldof.RowSize(group);
580 nbr_ldof.AddColumnsInRow(nbr, nldofs);
581 }
582 }
583 }
584 nbr_ldof.MakeJ();
585 for (int nbr = 1; nbr < nbr_recv_groups.Size(); nbr++)
586 {
587 const int num_recv_groups = nbr_recv_groups.RowSize(nbr);
588 if (num_recv_groups > 0)
589 {
590 const int *grp_list = nbr_recv_groups.GetRow(nbr);
591 for (int i = 0; i < num_recv_groups; i++)
592 {
593 const int group = grp_list[i];
594 const int nldofs = group_ldof.RowSize(group);
595 const int *ldofs = group_ldof.GetRow(group);
596 nbr_ldof.AddConnections(nbr, ldofs, nldofs);
597 }
598 }
599 }
600 nbr_ldof.ShiftUpI();
601}
602
603template <class T>
604T *GroupCommunicator::CopyGroupToBuffer(const T *ldata, T *buf, int group,
605 int layout) const
606{
607 switch (layout)
608 {
609 case 1:
610 {
611 return std::copy(ldata + group_ldof.GetI()[group],
612 ldata + group_ldof.GetI()[group+1],
613 buf);
614 }
615 case 2:
616 {
617 const int nltdofs = group_ltdof.RowSize(group);
618 const int *ltdofs = group_ltdof.GetRow(group);
619 for (int j = 0; j < nltdofs; j++)
620 {
621 buf[j] = ldata[ltdofs[j]];
622 }
623 return buf + nltdofs;
624 }
625 default:
626 {
627 const int nldofs = group_ldof.RowSize(group);
628 const int *ldofs = group_ldof.GetRow(group);
629 for (int j = 0; j < nldofs; j++)
630 {
631 buf[j] = ldata[ldofs[j]];
632 }
633 return buf + nldofs;
634 }
635 }
636}
637
638template <class T>
639const T *GroupCommunicator::CopyGroupFromBuffer(const T *buf, T *ldata,
640 int group, int layout) const
641{
642 const int nldofs = group_ldof.RowSize(group);
643 switch (layout)
644 {
645 case 1:
646 {
647 std::copy(buf, buf + nldofs, ldata + group_ldof.GetI()[group]);
648 break;
649 }
650 case 2:
651 {
652 const int *ltdofs = group_ltdof.GetRow(group);
653 for (int j = 0; j < nldofs; j++)
654 {
655 ldata[ltdofs[j]] = buf[j];
656 }
657 break;
658 }
659 default:
660 {
661 const int *ldofs = group_ldof.GetRow(group);
662 for (int j = 0; j < nldofs; j++)
663 {
664 ldata[ldofs[j]] = buf[j];
665 }
666 break;
667 }
668 }
669 return buf + nldofs;
670}
671
672template <class T>
673const T *GroupCommunicator::ReduceGroupFromBuffer(const T *buf, T *ldata,
674 int group, int layout,
675 void (*Op)(OpData<T>)) const
676{
677 OpData<T> opd;
678 opd.ldata = ldata;
679 opd.nldofs = group_ldof.RowSize(group);
680 opd.nb = 1;
681 opd.buf = const_cast<T*>(buf);
682
683 switch (layout)
684 {
685 case 1:
686 {
687 MFEM_ABORT("layout 1 is not supported");
688 T *dest = ldata + group_ldof.GetI()[group];
689 for (int j = 0; j < opd.nldofs; j++)
690 {
691 dest[j] += buf[j];
692 }
693 break;
694 }
695 case 2:
696 {
697 opd.ldofs = const_cast<int*>(group_ltdof.GetRow(group));
698 Op(opd);
699 break;
700 }
701 default:
702 {
703 opd.ldofs = const_cast<int*>(group_ldof.GetRow(group));
704 Op(opd);
705 break;
706 }
707 }
708 return buf + opd.nldofs;
709}
710
711template <class T>
712void GroupCommunicator::BcastBegin(T *ldata, int layout) const
713{
714 MFEM_VERIFY(comm_lock == 0, "object is already in use");
715
716 if (group_buf_size == 0) { return; }
717
718 int request_counter = 0;
719 switch (mode)
720 {
721 case byGroup: // ***** Communication by groups *****
722 {
723 T *buf;
724 if (layout != 1)
725 {
727 buf = (T *)group_buf.GetData();
728 MFEM_VERIFY(layout != 2 || group_ltdof.Size() == group_ldof.Size(),
729 "'group_ltdof' is not set, use SetLTDofTable()");
730 }
731 else
732 {
733 buf = ldata;
734 }
735
736 for (int gr = 1; gr < group_ldof.Size(); gr++)
737 {
738 const int nldofs = group_ldof.RowSize(gr);
739
740 // ignore groups without dofs
741 if (nldofs == 0) { continue; }
742
743 if (!gtopo.IAmMaster(gr)) // we are not the master
744 {
745 MPI_Irecv(buf,
746 nldofs,
749 40822 + gtopo.GetGroupMasterGroup(gr),
750 gtopo.GetComm(),
751 &requests[request_counter]);
752 request_marker[request_counter] = gr;
753 request_counter++;
754 }
755 else // we are the master
756 {
757 if (layout != 1)
758 {
759 CopyGroupToBuffer(ldata, buf, gr, layout);
760 }
761 const int gs = gtopo.GetGroupSize(gr);
762 const int *nbs = gtopo.GetGroup(gr);
763 for (int i = 0; i < gs; i++)
764 {
765 if (nbs[i] != 0)
766 {
767 MPI_Isend(buf,
768 nldofs,
770 gtopo.GetNeighborRank(nbs[i]),
771 40822 + gtopo.GetGroupMasterGroup(gr),
772 gtopo.GetComm(),
773 &requests[request_counter]);
774 request_marker[request_counter] = -1; // mark as send req.
775 request_counter++;
776 }
777 }
778 }
779 buf += nldofs;
780 }
781 break;
782 }
783
784 case byNeighbor: // ***** Communication by neighbors *****
785 {
787 T *buf = (T *)group_buf.GetData();
788 for (int nbr = 1; nbr < nbr_send_groups.Size(); nbr++)
789 {
790 const int num_send_groups = nbr_send_groups.RowSize(nbr);
791 if (num_send_groups > 0)
792 {
793 // Possible optimization:
794 // if (num_send_groups == 1) and (layout == 1) then we do not
795 // need to copy the data in order to send it.
796 T *buf_start = buf;
797 const int *grp_list = nbr_send_groups.GetRow(nbr);
798 for (int i = 0; i < num_send_groups; i++)
799 {
800 buf = CopyGroupToBuffer(ldata, buf, grp_list[i], layout);
801 }
802 MPI_Isend(buf_start,
803 buf - buf_start,
806 40822,
807 gtopo.GetComm(),
808 &requests[request_counter]);
809 request_marker[request_counter] = -1; // mark as send request
810 request_counter++;
811 }
812
813 const int num_recv_groups = nbr_recv_groups.RowSize(nbr);
814 if (num_recv_groups > 0)
815 {
816 // Possible optimization (requires interface change):
817 // if (num_recv_groups == 1) and the (output layout == 1) then
818 // we can receive directly in the output buffer; however, at
819 // this point we do not have that information.
820 const int *grp_list = nbr_recv_groups.GetRow(nbr);
821 int recv_size = 0;
822 for (int i = 0; i < num_recv_groups; i++)
823 {
824 recv_size += group_ldof.RowSize(grp_list[i]);
825 }
826 MPI_Irecv(buf,
827 recv_size,
830 40822,
831 gtopo.GetComm(),
832 &requests[request_counter]);
833 request_marker[request_counter] = nbr;
834 request_counter++;
835 buf_offsets[nbr] = buf - (T*)group_buf.GetData();
836 buf += recv_size;
837 }
838 }
839 MFEM_ASSERT(buf - (T*)group_buf.GetData() == group_buf_size, "");
840 break;
841 }
842 }
843
844 comm_lock = 1; // 1 - locked for Bcast
845 num_requests = request_counter;
846}
847
848template <class T>
849void GroupCommunicator::BcastEnd(T *ldata, int layout) const
850{
851 if (comm_lock == 0) { return; }
852 // The above also handles the case (group_buf_size == 0).
853 MFEM_VERIFY(comm_lock == 1, "object is NOT locked for Bcast");
854
855 switch (mode)
856 {
857 case byGroup: // ***** Communication by groups *****
858 {
859 if (layout == 1)
860 {
861 MPI_Waitall(num_requests, requests, MPI_STATUSES_IGNORE);
862 }
863 else if (layout == 0)
864 {
865 // copy the received data from the buffer to ldata, as it arrives
866 int idx;
867 while (MPI_Waitany(num_requests, requests, &idx, MPI_STATUS_IGNORE),
868 idx != MPI_UNDEFINED)
869 {
870 int gr = request_marker[idx];
871 if (gr == -1) { continue; } // skip send requests
872
873 // groups without dofs are skipped, so here nldofs > 0.
874 T *buf = (T *)group_buf.GetData() + group_ldof.GetI()[gr];
875 CopyGroupFromBuffer(buf, ldata, gr, layout);
876 }
877 }
878 break;
879 }
880
881 case byNeighbor: // ***** Communication by neighbors *****
882 {
883 // copy the received data from the buffer to ldata, as it arrives
884 int idx;
885 while (MPI_Waitany(num_requests, requests, &idx, MPI_STATUS_IGNORE),
886 idx != MPI_UNDEFINED)
887 {
888 int nbr = request_marker[idx];
889 if (nbr == -1) { continue; } // skip send requests
890
891 const int num_recv_groups = nbr_recv_groups.RowSize(nbr);
892 if (num_recv_groups > 0)
893 {
894 const int *grp_list = nbr_recv_groups.GetRow(nbr);
895 const T *buf = (T*)group_buf.GetData() + buf_offsets[nbr];
896 for (int i = 0; i < num_recv_groups; i++)
897 {
898 buf = CopyGroupFromBuffer(buf, ldata, grp_list[i], layout);
899 }
900 }
901 }
902 break;
903 }
904 }
905
906 comm_lock = 0; // 0 - no lock
907 num_requests = 0;
908}
909
910template <class T>
911void GroupCommunicator::ReduceBegin(const T *ldata) const
912{
913 MFEM_VERIFY(comm_lock == 0, "object is already in use");
914
915 if (group_buf_size == 0) { return; }
916
917 int request_counter = 0;
919 T *buf = (T *)group_buf.GetData();
920 switch (mode)
921 {
922 case byGroup: // ***** Communication by groups *****
923 {
924 for (int gr = 1; gr < group_ldof.Size(); gr++)
925 {
926 const int nldofs = group_ldof.RowSize(gr);
927 // ignore groups without dofs
928 if (nldofs == 0) { continue; }
929
930 if (!gtopo.IAmMaster(gr)) // we are not the master
931 {
932 const int layout = 0;
933 CopyGroupToBuffer(ldata, buf, gr, layout);
934 MPI_Isend(buf,
935 nldofs,
938 43822 + gtopo.GetGroupMasterGroup(gr),
939 gtopo.GetComm(),
940 &requests[request_counter]);
941 request_marker[request_counter] = -1; // mark as send request
942 request_counter++;
943 buf += nldofs;
944 }
945 else // we are the master
946 {
947 const int gs = gtopo.GetGroupSize(gr);
948 const int *nbs = gtopo.GetGroup(gr);
949 buf_offsets[gr] = buf - (T *)group_buf.GetData();
950 for (int i = 0; i < gs; i++)
951 {
952 if (nbs[i] != 0)
953 {
954 MPI_Irecv(buf,
955 nldofs,
957 gtopo.GetNeighborRank(nbs[i]),
958 43822 + gtopo.GetGroupMasterGroup(gr),
959 gtopo.GetComm(),
960 &requests[request_counter]);
961 request_marker[request_counter] = gr;
962 request_counter++;
963 buf += nldofs;
964 }
965 }
966 }
967 }
968 break;
969 }
970
971 case byNeighbor: // ***** Communication by neighbors *****
972 {
973 for (int nbr = 1; nbr < nbr_send_groups.Size(); nbr++)
974 {
975 // In Reduce operation: send_groups <--> recv_groups
976 const int num_send_groups = nbr_recv_groups.RowSize(nbr);
977 if (num_send_groups > 0)
978 {
979 T *buf_start = buf;
980 const int *grp_list = nbr_recv_groups.GetRow(nbr);
981 for (int i = 0; i < num_send_groups; i++)
982 {
983 const int layout = 0; // ldata is an array on all ldofs
984 buf = CopyGroupToBuffer(ldata, buf, grp_list[i], layout);
985 }
986 MPI_Isend(buf_start,
987 buf - buf_start,
990 43822,
991 gtopo.GetComm(),
992 &requests[request_counter]);
993 request_marker[request_counter] = -1; // mark as send request
994 request_counter++;
995 }
996
997 // In Reduce operation: send_groups <--> recv_groups
998 const int num_recv_groups = nbr_send_groups.RowSize(nbr);
999 if (num_recv_groups > 0)
1000 {
1001 const int *grp_list = nbr_send_groups.GetRow(nbr);
1002 int recv_size = 0;
1003 for (int i = 0; i < num_recv_groups; i++)
1004 {
1005 recv_size += group_ldof.RowSize(grp_list[i]);
1006 }
1007 MPI_Irecv(buf,
1008 recv_size,
1011 43822,
1012 gtopo.GetComm(),
1013 &requests[request_counter]);
1014 request_marker[request_counter] = nbr;
1015 request_counter++;
1016 buf_offsets[nbr] = buf - (T*)group_buf.GetData();
1017 buf += recv_size;
1018 }
1019 }
1020 MFEM_ASSERT(buf - (T*)group_buf.GetData() == group_buf_size, "");
1021 break;
1022 }
1023 }
1024
1025 comm_lock = 2;
1026 num_requests = request_counter;
1027}
1028
1029template <class T>
1030void GroupCommunicator::ReduceEnd(T *ldata, int layout,
1031 void (*Op)(OpData<T>)) const
1032{
1033 if (comm_lock == 0) { return; }
1034 // The above also handles the case (group_buf_size == 0).
1035 MFEM_VERIFY(comm_lock == 2, "object is NOT locked for Reduce");
1036
1037 switch (mode)
1038 {
1039 case byGroup: // ***** Communication by groups *****
1040 {
1041 OpData<T> opd;
1042 opd.ldata = ldata;
1043 Array<int> group_num_req(group_ldof.Size());
1044 for (int gr = 1; gr < group_ldof.Size(); gr++)
1045 {
1046 group_num_req[gr] =
1047 gtopo.IAmMaster(gr) ? gtopo.GetGroupSize(gr)-1 : 0;
1048 }
1049 int idx;
1050 while (MPI_Waitany(num_requests, requests, &idx, MPI_STATUS_IGNORE),
1051 idx != MPI_UNDEFINED)
1052 {
1053 int gr = request_marker[idx];
1054 if (gr == -1) { continue; } // skip send requests
1055
1056 // Delay the processing of a group until all receive requests, for
1057 // that group, are done:
1058 if ((--group_num_req[gr]) != 0) { continue; }
1059
1060 opd.nldofs = group_ldof.RowSize(gr);
1061 // groups without dofs are skipped, so here nldofs > 0.
1062
1063 opd.buf = (T *)group_buf.GetData() + buf_offsets[gr];
1064 opd.ldofs = (layout == 0) ?
1066 opd.nb = gtopo.GetGroupSize(gr)-1;
1067 Op(opd);
1068 }
1069 break;
1070 }
1071
1072 case byNeighbor: // ***** Communication by neighbors *****
1073 {
1074 MPI_Waitall(num_requests, requests, MPI_STATUSES_IGNORE);
1075
1076 for (int nbr = 1; nbr < nbr_send_groups.Size(); nbr++)
1077 {
1078 // In Reduce operation: send_groups <--> recv_groups
1079 const int num_recv_groups = nbr_send_groups.RowSize(nbr);
1080 if (num_recv_groups > 0)
1081 {
1082 const int *grp_list = nbr_send_groups.GetRow(nbr);
1083 const T *buf = (T*)group_buf.GetData() + buf_offsets[nbr];
1084 for (int i = 0; i < num_recv_groups; i++)
1085 {
1086 buf = ReduceGroupFromBuffer(buf, ldata, grp_list[i],
1087 layout, Op);
1088 }
1089 }
1090 }
1091 break;
1092 }
1093 }
1094
1095 comm_lock = 0; // 0 - no lock
1096 num_requests = 0;
1097}
1098
1099template <class T>
1101{
1102 if (opd.nb == 1)
1103 {
1104 for (int i = 0; i < opd.nldofs; i++)
1105 {
1106 opd.ldata[opd.ldofs[i]] += opd.buf[i];
1107 }
1108 }
1109 else
1110 {
1111 for (int i = 0; i < opd.nldofs; i++)
1112 {
1113 T data = opd.ldata[opd.ldofs[i]];
1114 for (int j = 0; j < opd.nb; j++)
1115 {
1116 data += opd.buf[j*opd.nldofs+i];
1117 }
1118 opd.ldata[opd.ldofs[i]] = data;
1119 }
1120 }
1121}
1122
1123template <class T>
1125{
1126 for (int i = 0; i < opd.nldofs; i++)
1127 {
1128 T data = opd.ldata[opd.ldofs[i]];
1129 for (int j = 0; j < opd.nb; j++)
1130 {
1131 T b = opd.buf[j*opd.nldofs+i];
1132 if (data > b)
1133 {
1134 data = b;
1135 }
1136 }
1137 opd.ldata[opd.ldofs[i]] = data;
1138 }
1139}
1140
1141template <class T>
1143{
1144 for (int i = 0; i < opd.nldofs; i++)
1145 {
1146 T data = opd.ldata[opd.ldofs[i]];
1147 for (int j = 0; j < opd.nb; j++)
1148 {
1149 T b = opd.buf[j*opd.nldofs+i];
1150 if (data < b)
1151 {
1152 data = b;
1153 }
1154 }
1155 opd.ldata[opd.ldofs[i]] = data;
1156 }
1157}
1158
1159template <class T>
1161{
1162 for (int i = 0; i < opd.nldofs; i++)
1163 {
1164 T data = opd.ldata[opd.ldofs[i]];
1165 for (int j = 0; j < opd.nb; j++)
1166 {
1167 data |= opd.buf[j*opd.nldofs+i];
1168 }
1169 opd.ldata[opd.ldofs[i]] = data;
1170 }
1171}
1172
1173void GroupCommunicator::PrintInfo(std::ostream &os) const
1174{
1175 char c = '\0';
1176 const int tag = 46800;
1177 const int myid = gtopo.MyRank();
1178
1179 int num_sends = 0, num_recvs = 0;
1180 size_t mem_sends = 0, mem_recvs = 0;
1181 int num_master_groups = 0, num_empty_groups = 0;
1182 int num_active_neighbors = 0; // for mode == byNeighbor
1183 switch (mode)
1184 {
1185 case byGroup:
1186 for (int gr = 1; gr < group_ldof.Size(); gr++)
1187 {
1188 const int nldofs = group_ldof.RowSize(gr);
1189 if (nldofs == 0)
1190 {
1191 num_empty_groups++;
1192 continue;
1193 }
1194 if (gtopo.IAmMaster(gr))
1195 {
1196 num_sends += (gtopo.GetGroupSize(gr)-1);
1197 mem_sends += sizeof(double)*nldofs*(gtopo.GetGroupSize(gr)-1);
1198 num_master_groups++;
1199 }
1200 else
1201 {
1202 num_recvs++;
1203 mem_recvs += sizeof(double)*nldofs;
1204 }
1205 }
1206 break;
1207
1208 case byNeighbor:
1209 for (int gr = 1; gr < group_ldof.Size(); gr++)
1210 {
1211 const int nldofs = group_ldof.RowSize(gr);
1212 if (nldofs == 0)
1213 {
1214 num_empty_groups++;
1215 continue;
1216 }
1217 if (gtopo.IAmMaster(gr))
1218 {
1219 num_master_groups++;
1220 }
1221 }
1222 for (int nbr = 1; nbr < nbr_send_groups.Size(); nbr++)
1223 {
1224 const int num_send_groups = nbr_send_groups.RowSize(nbr);
1225 if (num_send_groups > 0)
1226 {
1227 const int *grp_list = nbr_send_groups.GetRow(nbr);
1228 for (int i = 0; i < num_send_groups; i++)
1229 {
1230 mem_sends += sizeof(double)*group_ldof.RowSize(grp_list[i]);
1231 }
1232 num_sends++;
1233 }
1234
1235 const int num_recv_groups = nbr_recv_groups.RowSize(nbr);
1236 if (num_recv_groups > 0)
1237 {
1238 const int *grp_list = nbr_recv_groups.GetRow(nbr);
1239 for (int i = 0; i < num_recv_groups; i++)
1240 {
1241 mem_recvs += sizeof(double)*group_ldof.RowSize(grp_list[i]);
1242 }
1243 num_recvs++;
1244 }
1245 if (num_send_groups > 0 || num_recv_groups > 0)
1246 {
1247 num_active_neighbors++;
1248 }
1249 }
1250 break;
1251 }
1252 if (myid != 0)
1253 {
1254 MPI_Recv(&c, 1, MPI_CHAR, myid-1, tag, gtopo.GetComm(),
1255 MPI_STATUS_IGNORE);
1256 }
1257 else
1258 {
1259 os << "\nGroupCommunicator:\n";
1260 }
1261 os << "Rank " << myid << ":\n"
1262 " mode = " <<
1263 (mode == byGroup ? "byGroup" : "byNeighbor") << "\n"
1264 " number of sends = " << num_sends <<
1265 " (" << mem_sends << " bytes)\n"
1266 " number of recvs = " << num_recvs <<
1267 " (" << mem_recvs << " bytes)\n";
1268 os <<
1269 " num groups = " << group_ldof.Size() << " = " <<
1270 num_master_groups << " + " <<
1271 group_ldof.Size()-num_master_groups-num_empty_groups << " + " <<
1272 num_empty_groups << " (master + slave + empty)\n";
1273 if (mode == byNeighbor)
1274 {
1275 os <<
1276 " num neighbors = " << nbr_send_groups.Size() << " = " <<
1277 num_active_neighbors << " + " <<
1278 nbr_send_groups.Size()-num_active_neighbors <<
1279 " (active + inactive)\n";
1280 }
1281 if (myid != gtopo.NRanks()-1)
1282 {
1283 os << std::flush;
1284 MPI_Send(&c, 1, MPI_CHAR, myid+1, tag, gtopo.GetComm());
1285 }
1286 else
1287 {
1288 os << std::endl;
1289 }
1290 MPI_Barrier(gtopo.GetComm());
1291}
1292
1294{
1295 delete [] buf_offsets;
1296 delete [] request_marker;
1297 // delete [] statuses;
1298 delete [] requests;
1299}
1300
1301// @cond DOXYGEN_SKIP
1302
1303// instantiate GroupCommunicator::Bcast and Reduce for int and double
1304template void GroupCommunicator::BcastBegin<int>(int *, int) const;
1305template void GroupCommunicator::BcastEnd<int>(int *, int) const;
1306template void GroupCommunicator::ReduceBegin<int>(const int *) const;
1308 int *, int, void (*)(OpData<int>)) const;
1309
1310template void GroupCommunicator::BcastBegin<double>(double *, int) const;
1311template void GroupCommunicator::BcastEnd<double>(double *, int) const;
1312template void GroupCommunicator::ReduceBegin<double>(const double *) const;
1314 double *, int, void (*)(OpData<double>)) const;
1315
1316template void GroupCommunicator::BcastBegin<float>(float *, int) const;
1317template void GroupCommunicator::BcastEnd<float>(float *, int) const;
1318template void GroupCommunicator::ReduceBegin<float>(const float *) const;
1320 float *, int, void (*)(OpData<float>)) const;
1321
1322// @endcond
1323
1324// instantiate reduce operators for int and double
1325template void GroupCommunicator::Sum<int>(OpData<int>);
1326template void GroupCommunicator::Min<int>(OpData<int>);
1327template void GroupCommunicator::Max<int>(OpData<int>);
1328template void GroupCommunicator::BitOR<int>(OpData<int>);
1329
1330template void GroupCommunicator::Sum<double>(OpData<double>);
1331template void GroupCommunicator::Min<double>(OpData<double>);
1332template void GroupCommunicator::Max<double>(OpData<double>);
1333
1334template void GroupCommunicator::Sum<float>(OpData<float>);
1335template void GroupCommunicator::Min<float>(OpData<float>);
1336template void GroupCommunicator::Max<float>(OpData<float>);
1337
1338
1339#ifdef __bgq__
1340static void DebugRankCoords(int** coords, int dim, int size)
1341{
1342 for (int i = 0; i < size; i++)
1343 {
1344 mfem::out << "Rank " << i << " coords: ";
1345 for (int j = 0; j < dim; j++)
1346 {
1347 mfem::out << coords[i][j] << " ";
1348 }
1349 mfem::out << endl;
1350 }
1351}
1352
1353struct CompareCoords
1354{
1355 CompareCoords(int coord) : coord(coord) {}
1356 int coord;
1357
1358 bool operator()(int* const &a, int* const &b) const
1359 { return a[coord] < b[coord]; }
1360};
1361
1362void KdTreeSort(int** coords, int d, int dim, int size)
1363{
1364 if (size > 1)
1365 {
1366 bool all_same = true;
1367 for (int i = 1; i < size && all_same; i++)
1368 {
1369 for (int j = 0; j < dim; j++)
1370 {
1371 if (coords[i][j] != coords[0][j]) { all_same = false; break; }
1372 }
1373 }
1374 if (all_same) { return; }
1375
1376 // sort by coordinate 'd'
1377 std::sort(coords, coords + size, CompareCoords(d));
1378 int next = (d + 1) % dim;
1379
1380 if (coords[0][d] < coords[size-1][d])
1381 {
1382 KdTreeSort(coords, next, dim, size/2);
1383 KdTreeSort(coords + size/2, next, dim, size - size/2);
1384 }
1385 else
1386 {
1387 // skip constant dimension
1388 KdTreeSort(coords, next, dim, size);
1389 }
1390 }
1391}
1392
1393MPI_Comm ReorderRanksZCurve(MPI_Comm comm)
1394{
1395 MPI_Status status;
1396
1397 int rank, size;
1398 MPI_Comm_rank(comm, &rank);
1399 MPI_Comm_size(comm, &size);
1400
1401 int dim;
1402 MPIX_Torus_ndims(&dim);
1403
1404 int* mycoords = new int[dim + 1];
1405 MPIX_Rank2torus(rank, mycoords);
1406
1407 MPI_Send(mycoords, dim, MPI_INT, 0, 111, comm);
1408 delete [] mycoords;
1409
1410 if (rank == 0)
1411 {
1412 int** coords = new int*[size];
1413 for (int i = 0; i < size; i++)
1414 {
1415 coords[i] = new int[dim + 1];
1416 coords[i][dim] = i;
1417 MPI_Recv(coords[i], dim, MPI_INT, i, 111, comm, &status);
1418 }
1419
1420 KdTreeSort(coords, 0, dim, size);
1421
1422 // DebugRankCoords(coords, dim, size);
1423
1424 for (int i = 0; i < size; i++)
1425 {
1426 MPI_Send(&coords[i][dim], 1, MPI_INT, i, 112, comm);
1427 delete [] coords[i];
1428 }
1429 delete [] coords;
1430 }
1431
1432 int new_rank;
1433 MPI_Recv(&new_rank, 1, MPI_INT, 0, 112, comm, &status);
1434
1435 MPI_Comm new_comm;
1436 MPI_Comm_split(comm, 0, new_rank, &new_comm);
1437 return new_comm;
1438}
1439
1440#else // __bgq__
1441
1442MPI_Comm ReorderRanksZCurve(MPI_Comm comm)
1443{
1444 // pass
1445 return comm;
1446}
1447#endif // __bgq__
1448
1449} // namespace mfem
1450
1451#endif
void Sort()
Sorts the array in ascending order. This requires operator< to be defined for T.
Definition array.hpp:261
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition array.hpp:697
int Size() const
Return the logical size of the array.
Definition array.hpp:144
int Append(const T &el)
Append element 'el' to array, resize if necessary.
Definition array.hpp:769
T * GetData()
Returns the data.
Definition array.hpp:118
void Copy(Array &copy) const
Create a copy of the internal array to the provided copy.
Definition array.hpp:874
void Unique()
Removes duplicities from a sorted array. This requires operator== to be defined for T.
Definition array.hpp:269
GroupCommunicator(const GroupTopology &gt, Mode m=byNeighbor)
Construct a GroupCommunicator object.
void GetNeighborLTDofTable(Table &nbr_ltdof) const
Dofs to be sent to communication neighbors.
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().
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'.
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.
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'.
A set of integers.
Definition sets.hpp:24
void Recreate(const int n, const int *p)
Create an integer set from C-array 'p' of 'n' integers. Overwrites any existing set data.
Definition sets.cpp:68
List of integer sets.
Definition sets.hpp:63
int Insert(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...
Definition sets.cpp:91
void AsTable(Table &t)
Write the list of sets into table 't'.
Definition sets.cpp:116
int Lookup(IntegerSet &s)
Definition sets.cpp:104
int PickElementInSet(int i)
Return the value of the first element of the ith set.
Definition sets.hpp:73
static MFEM_EXPORT int default_thread_required
Default level of thread support for MPI_Init_thread.
int * GetJ()
Definition table.hpp:114
void AddConnections(int r, const int *c, int nc)
Definition table.cpp:104
int RowSize(int i) const
Definition table.hpp:108
void ShiftUpI()
Definition table.cpp:115
void GetRow(int i, Array< int > &row) const
Return row i in array row (the Table must be finalized)
Definition table.cpp:187
void AddConnection(int r, int c)
Definition table.hpp:80
void MakeI(int nrows)
Next 7 methods are used together with the default constructor.
Definition table.cpp:81
int Size() const
Returns the number of TYPE I elements.
Definition table.hpp:92
int Size_of_connections() const
Definition table.hpp:98
void AddColumnsInRow(int r, int ncol)
Definition table.hpp:78
void MakeFromList(int nrows, const Array< Connection > &list)
Definition table.cpp:280
void Copy(Table &copy) const
Definition table.cpp:390
void MakeJ()
Definition table.cpp:91
int * GetI()
Definition table.hpp:113
void AddAColumnInRow(int r)
Definition table.hpp:77
void SetDims(int rows, int nnz)
Definition table.cpp:140
int dim
Definition ex24.cpp:53
int index(int i, int j, int nx, int ny)
Definition life.cpp:235
real_t b
Definition lissajous.cpp:42
real_t a
Definition lissajous.cpp:41
void Swap(Array< T > &, Array< T > &)
Definition array.hpp:648
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
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.
Definition text.hpp:31
real_t p(const Vector &x, real_t t)
RefCoord s[3]
Helper struct for defining a connectivity table, see Table::MakeFromList.
Definition table.hpp:28
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.