MFEM v4.9.0
Finite element discretization library
Loading...
Searching...
No Matches
communication.cpp
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#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/// \cond DO_NOT_DOCUMENT
351// Initialize the static mpi_type for the specializations of MPITypeMap:
352const MPI_Datatype MPITypeMap<bool>::mpi_type = MFEM_MPI_CXX_BOOL;
353const MPI_Datatype MPITypeMap<char>::mpi_type = MPI_CHAR;
354const MPI_Datatype MPITypeMap<unsigned char>::mpi_type = MPI_UNSIGNED_CHAR;
355const MPI_Datatype MPITypeMap<short>::mpi_type = MPI_SHORT;
356const MPI_Datatype MPITypeMap<unsigned short>::mpi_type = MPI_UNSIGNED_SHORT;
357const MPI_Datatype MPITypeMap<int>::mpi_type = MPI_INT;
358const MPI_Datatype MPITypeMap<unsigned int>::mpi_type = MPI_UNSIGNED;
359const MPI_Datatype MPITypeMap<long>::mpi_type = MPI_LONG;
360const MPI_Datatype MPITypeMap<unsigned long>::mpi_type = MPI_UNSIGNED_LONG;
361const MPI_Datatype MPITypeMap<long long>::mpi_type = MPI_LONG_LONG;
363 MPI_UNSIGNED_LONG_LONG;
364const MPI_Datatype MPITypeMap<float>::mpi_type = MPI_FLOAT;
365const MPI_Datatype MPITypeMap<double>::mpi_type = MPI_DOUBLE;
366/// \endcond DO_NOT_DOCUMENT
367
369 : gtopo(gt), mode(m)
370{
371 group_buf_size = 0;
372 requests = NULL;
373 // statuses = NULL;
374 comm_lock = 0;
375 num_requests = 0;
376 request_marker = NULL;
377 buf_offsets = NULL;
378}
379
381{
383 for (int i = 0; i < ldof_group.Size(); i++)
384 {
385 int group = ldof_group[i];
386 if (group != 0)
387 {
389 }
390 }
392
393 for (int i = 0; i < ldof_group.Size(); i++)
394 {
395 int group = ldof_group[i];
396 if (group != 0)
397 {
398 group_ldof.AddConnection(group, i);
399 }
400 }
402
403 Finalize();
404}
405
407{
408 int request_counter = 0;
409
410 // size buf_offsets = max(number of groups, number of neighbors)
411 buf_offsets = new int[max(group_ldof.Size(), gtopo.GetNumNeighbors())];
412 buf_offsets[0] = 0;
413 for (int gr = 1; gr < group_ldof.Size(); gr++)
414 {
415 if (group_ldof.RowSize(gr) != 0)
416 {
417 int gr_requests;
418 if (!gtopo.IAmMaster(gr)) // we are not the master
419 {
420 gr_requests = 1;
421 }
422 else
423 {
424 gr_requests = gtopo.GetGroupSize(gr)-1;
425 }
426
427 request_counter += gr_requests;
428 group_buf_size += gr_requests * group_ldof.RowSize(gr);
429 }
430 }
431
432 requests = new MPI_Request[request_counter];
433 // statuses = new MPI_Status[request_counter];
434 request_marker = new int[request_counter];
435
436 // Construct nbr_send_groups and nbr_recv_groups: (nbr 0 = me)
439 for (int gr = 1; gr < group_ldof.Size(); gr++)
440 {
441 const int nldofs = group_ldof.RowSize(gr);
442 if (nldofs == 0) { continue; }
443
444 if (!gtopo.IAmMaster(gr)) // we are not the master
445 {
447 }
448 else // we are the master
449 {
450 const int grp_size = gtopo.GetGroupSize(gr);
451 const int *grp_nbr_list = gtopo.GetGroup(gr);
452 for (int i = 0; i < grp_size; i++)
453 {
454 if (grp_nbr_list[i] != 0)
455 {
456 nbr_send_groups.AddAColumnInRow(grp_nbr_list[i]);
457 }
458 }
459 }
460 }
463 for (int gr = 1; gr < group_ldof.Size(); gr++)
464 {
465 const int nldofs = group_ldof.RowSize(gr);
466 if (nldofs == 0) { continue; }
467
468 if (!gtopo.IAmMaster(gr)) // we are not the master
469 {
471 }
472 else // we are the master
473 {
474 const int grp_size = gtopo.GetGroupSize(gr);
475 const int *grp_nbr_list = gtopo.GetGroup(gr);
476 for (int i = 0; i < grp_size; i++)
477 {
478 if (grp_nbr_list[i] != 0)
479 {
480 nbr_send_groups.AddConnection(grp_nbr_list[i], gr);
481 }
482 }
483 }
484 }
487 // The above construction creates the Tables with the column indices
488 // sorted, i.e. the group lists are sorted. To coordinate this order between
489 // processors, we will sort the group lists in the nbr_recv_groups Table
490 // according to their indices in the master. This does not require any
491 // communication because we have access to the group indices in the master
492 // by calling: master_group_id = gtopo.GetGroupMasterGroup(my_group_id).
493 Array<Pair<int,int> > group_ids;
494 for (int nbr = 1; nbr < nbr_recv_groups.Size(); nbr++)
495 {
496 const int num_recv_groups = nbr_recv_groups.RowSize(nbr);
497 if (num_recv_groups > 0)
498 {
499 int *grp_list = nbr_recv_groups.GetRow(nbr);
500 group_ids.SetSize(num_recv_groups);
501 for (int i = 0; i < num_recv_groups; i++)
502 {
503 group_ids[i].one = gtopo.GetGroupMasterGroup(grp_list[i]);
504 group_ids[i].two = grp_list[i]; // my_group_id
505 }
506 group_ids.Sort();
507 for (int i = 0; i < num_recv_groups; i++)
508 {
509 grp_list[i] = group_ids[i].two;
510 }
511 }
512 }
513}
514
516{
517 if (group_ltdof.Size() == group_ldof.Size()) { return; }
518
520 for (int gr = 1; gr < group_ldof.Size(); gr++)
521 {
522 if (gtopo.IAmMaster(gr))
523 {
525 }
526 }
528 for (int gr = 1; gr < group_ldof.Size(); gr++)
529 {
530 if (gtopo.IAmMaster(gr))
531 {
532 const int *ldofs = group_ldof.GetRow(gr);
533 const int nldofs = group_ldof.RowSize(gr);
534 for (int i = 0; i < nldofs; i++)
535 {
536 group_ltdof.AddConnection(gr, ldof_ltdof[ldofs[i]]);
537 }
538 }
539 }
541}
542
544{
545 nbr_ltdof.MakeI(nbr_send_groups.Size());
546 for (int nbr = 1; nbr < nbr_send_groups.Size(); nbr++)
547 {
548 const int num_send_groups = nbr_send_groups.RowSize(nbr);
549 if (num_send_groups > 0)
550 {
551 const int *grp_list = nbr_send_groups.GetRow(nbr);
552 for (int i = 0; i < num_send_groups; i++)
553 {
554 const int group = grp_list[i];
555 const int nltdofs = group_ltdof.RowSize(group);
556 nbr_ltdof.AddColumnsInRow(nbr, nltdofs);
557 }
558 }
559 }
560 nbr_ltdof.MakeJ();
561 for (int nbr = 1; nbr < nbr_send_groups.Size(); nbr++)
562 {
563 const int num_send_groups = nbr_send_groups.RowSize(nbr);
564 if (num_send_groups > 0)
565 {
566 const int *grp_list = nbr_send_groups.GetRow(nbr);
567 for (int i = 0; i < num_send_groups; i++)
568 {
569 const int group = grp_list[i];
570 const int nltdofs = group_ltdof.RowSize(group);
571 const int *ltdofs = group_ltdof.GetRow(group);
572 nbr_ltdof.AddConnections(nbr, ltdofs, nltdofs);
573 }
574 }
575 }
576 nbr_ltdof.ShiftUpI();
577}
578
580{
581 nbr_ldof.MakeI(nbr_recv_groups.Size());
582 for (int nbr = 1; nbr < nbr_recv_groups.Size(); nbr++)
583 {
584 const int num_recv_groups = nbr_recv_groups.RowSize(nbr);
585 if (num_recv_groups > 0)
586 {
587 const int *grp_list = nbr_recv_groups.GetRow(nbr);
588 for (int i = 0; i < num_recv_groups; i++)
589 {
590 const int group = grp_list[i];
591 const int nldofs = group_ldof.RowSize(group);
592 nbr_ldof.AddColumnsInRow(nbr, nldofs);
593 }
594 }
595 }
596 nbr_ldof.MakeJ();
597 for (int nbr = 1; nbr < nbr_recv_groups.Size(); nbr++)
598 {
599 const int num_recv_groups = nbr_recv_groups.RowSize(nbr);
600 if (num_recv_groups > 0)
601 {
602 const int *grp_list = nbr_recv_groups.GetRow(nbr);
603 for (int i = 0; i < num_recv_groups; i++)
604 {
605 const int group = grp_list[i];
606 const int nldofs = group_ldof.RowSize(group);
607 const int *ldofs = group_ldof.GetRow(group);
608 nbr_ldof.AddConnections(nbr, ldofs, nldofs);
609 }
610 }
611 }
612 nbr_ldof.ShiftUpI();
613}
614
615template <class T>
616T *GroupCommunicator::CopyGroupToBuffer(const T *ldata, T *buf, int group,
617 int layout) const
618{
619 switch (layout)
620 {
621 case 1:
622 {
623 return std::copy(ldata + group_ldof.GetI()[group],
624 ldata + group_ldof.GetI()[group+1],
625 buf);
626 }
627 case 2:
628 {
629 const int nltdofs = group_ltdof.RowSize(group);
630 const int *ltdofs = group_ltdof.GetRow(group);
631 for (int j = 0; j < nltdofs; j++)
632 {
633 buf[j] = ldata[ltdofs[j]];
634 }
635 return buf + nltdofs;
636 }
637 default:
638 {
639 const int nldofs = group_ldof.RowSize(group);
640 const int *ldofs = group_ldof.GetRow(group);
641 for (int j = 0; j < nldofs; j++)
642 {
643 buf[j] = ldata[ldofs[j]];
644 }
645 return buf + nldofs;
646 }
647 }
648}
649
650template <class T>
651const T *GroupCommunicator::CopyGroupFromBuffer(const T *buf, T *ldata,
652 int group, int layout) const
653{
654 const int nldofs = group_ldof.RowSize(group);
655 switch (layout)
656 {
657 case 1:
658 {
659 std::copy(buf, buf + nldofs, ldata + group_ldof.GetI()[group]);
660 break;
661 }
662 case 2:
663 {
664 const int *ltdofs = group_ltdof.GetRow(group);
665 for (int j = 0; j < nldofs; j++)
666 {
667 ldata[ltdofs[j]] = buf[j];
668 }
669 break;
670 }
671 default:
672 {
673 const int *ldofs = group_ldof.GetRow(group);
674 for (int j = 0; j < nldofs; j++)
675 {
676 ldata[ldofs[j]] = buf[j];
677 }
678 break;
679 }
680 }
681 return buf + nldofs;
682}
683
684template <class T>
685const T *GroupCommunicator::ReduceGroupFromBuffer(const T *buf, T *ldata,
686 int group, int layout,
687 void (*Op)(OpData<T>)) const
688{
689 OpData<T> opd;
690 opd.ldata = ldata;
691 opd.nldofs = group_ldof.RowSize(group);
692 opd.nb = 1;
693 opd.buf = const_cast<T*>(buf);
694
695 switch (layout)
696 {
697 case 1:
698 {
699 MFEM_ABORT("layout 1 is not supported");
700 T *dest = ldata + group_ldof.GetI()[group];
701 for (int j = 0; j < opd.nldofs; j++)
702 {
703 dest[j] += buf[j];
704 }
705 break;
706 }
707 case 2:
708 {
709 opd.ldofs = const_cast<int*>(group_ltdof.GetRow(group));
710 Op(opd);
711 break;
712 }
713 default:
714 {
715 opd.ldofs = const_cast<int*>(group_ldof.GetRow(group));
716 Op(opd);
717 break;
718 }
719 }
720 return buf + opd.nldofs;
721}
722
723template <class T>
724void GroupCommunicator::BcastBegin(T *ldata, int layout) const
725{
726 MFEM_VERIFY(comm_lock == 0, "object is already in use");
727
728 if (group_buf_size == 0) { return; }
729
730 int request_counter = 0;
731 switch (mode)
732 {
733 case byGroup: // ***** Communication by groups *****
734 {
735 T *buf;
736 if (layout != 1)
737 {
739 buf = (T *)group_buf.GetData();
740 MFEM_VERIFY(layout != 2 || group_ltdof.Size() == group_ldof.Size(),
741 "'group_ltdof' is not set, use SetLTDofTable()");
742 }
743 else
744 {
745 buf = ldata;
746 }
747
748 for (int gr = 1; gr < group_ldof.Size(); gr++)
749 {
750 const int nldofs = group_ldof.RowSize(gr);
751
752 // ignore groups without dofs
753 if (nldofs == 0) { continue; }
754
755 if (!gtopo.IAmMaster(gr)) // we are not the master
756 {
757 MPI_Irecv(buf,
758 nldofs,
761 40822 + gtopo.GetGroupMasterGroup(gr),
762 gtopo.GetComm(),
763 &requests[request_counter]);
764 request_marker[request_counter] = gr;
765 request_counter++;
766 }
767 else // we are the master
768 {
769 if (layout != 1)
770 {
771 CopyGroupToBuffer(ldata, buf, gr, layout);
772 }
773 const int gs = gtopo.GetGroupSize(gr);
774 const int *nbs = gtopo.GetGroup(gr);
775 for (int i = 0; i < gs; i++)
776 {
777 if (nbs[i] != 0)
778 {
779 MPI_Isend(buf,
780 nldofs,
782 gtopo.GetNeighborRank(nbs[i]),
783 40822 + gtopo.GetGroupMasterGroup(gr),
784 gtopo.GetComm(),
785 &requests[request_counter]);
786 request_marker[request_counter] = -1; // mark as send req.
787 request_counter++;
788 }
789 }
790 }
791 buf += nldofs;
792 }
793 break;
794 }
795
796 case byNeighbor: // ***** Communication by neighbors *****
797 {
799 T *buf = (T *)group_buf.GetData();
800 for (int nbr = 1; nbr < nbr_send_groups.Size(); nbr++)
801 {
802 const int num_send_groups = nbr_send_groups.RowSize(nbr);
803 if (num_send_groups > 0)
804 {
805 // Possible optimization:
806 // if (num_send_groups == 1) and (layout == 1) then we do not
807 // need to copy the data in order to send it.
808 T *buf_start = buf;
809 const int *grp_list = nbr_send_groups.GetRow(nbr);
810 for (int i = 0; i < num_send_groups; i++)
811 {
812 buf = CopyGroupToBuffer(ldata, buf, grp_list[i], layout);
813 }
814 MPI_Isend(buf_start,
815 buf - buf_start,
818 40822,
819 gtopo.GetComm(),
820 &requests[request_counter]);
821 request_marker[request_counter] = -1; // mark as send request
822 request_counter++;
823 }
824
825 const int num_recv_groups = nbr_recv_groups.RowSize(nbr);
826 if (num_recv_groups > 0)
827 {
828 // Possible optimization (requires interface change):
829 // if (num_recv_groups == 1) and the (output layout == 1) then
830 // we can receive directly in the output buffer; however, at
831 // this point we do not have that information.
832 const int *grp_list = nbr_recv_groups.GetRow(nbr);
833 int recv_size = 0;
834 for (int i = 0; i < num_recv_groups; i++)
835 {
836 recv_size += group_ldof.RowSize(grp_list[i]);
837 }
838 MPI_Irecv(buf,
839 recv_size,
842 40822,
843 gtopo.GetComm(),
844 &requests[request_counter]);
845 request_marker[request_counter] = nbr;
846 request_counter++;
847 buf_offsets[nbr] = buf - (T*)group_buf.GetData();
848 buf += recv_size;
849 }
850 }
851 MFEM_ASSERT(buf - (T*)group_buf.GetData() == group_buf_size, "");
852 break;
853 }
854 }
855
856 comm_lock = 1; // 1 - locked for Bcast
857 num_requests = request_counter;
858}
859
860template <class T>
861void GroupCommunicator::BcastEnd(T *ldata, int layout) const
862{
863 if (comm_lock == 0) { return; }
864 // The above also handles the case (group_buf_size == 0).
865 MFEM_VERIFY(comm_lock == 1, "object is NOT locked for Bcast");
866
867 switch (mode)
868 {
869 case byGroup: // ***** Communication by groups *****
870 {
871 if (layout == 1)
872 {
873 MPI_Waitall(num_requests, requests, MPI_STATUSES_IGNORE);
874 }
875 else if (layout == 0)
876 {
877 // copy the received data from the buffer to ldata, as it arrives
878 int idx;
879 while (MPI_Waitany(num_requests, requests, &idx, MPI_STATUS_IGNORE),
880 idx != MPI_UNDEFINED)
881 {
882 int gr = request_marker[idx];
883 if (gr == -1) { continue; } // skip send requests
884
885 // groups without dofs are skipped, so here nldofs > 0.
886 T *buf = (T *)group_buf.GetData() + group_ldof.GetI()[gr];
887 CopyGroupFromBuffer(buf, ldata, gr, layout);
888 }
889 }
890 break;
891 }
892
893 case byNeighbor: // ***** Communication by neighbors *****
894 {
895 // copy the received data from the buffer to ldata, as it arrives
896 int idx;
897 while (MPI_Waitany(num_requests, requests, &idx, MPI_STATUS_IGNORE),
898 idx != MPI_UNDEFINED)
899 {
900 int nbr = request_marker[idx];
901 if (nbr == -1) { continue; } // skip send requests
902
903 const int num_recv_groups = nbr_recv_groups.RowSize(nbr);
904 if (num_recv_groups > 0)
905 {
906 const int *grp_list = nbr_recv_groups.GetRow(nbr);
907 const T *buf = (T*)group_buf.GetData() + buf_offsets[nbr];
908 for (int i = 0; i < num_recv_groups; i++)
909 {
910 buf = CopyGroupFromBuffer(buf, ldata, grp_list[i], layout);
911 }
912 }
913 }
914 break;
915 }
916 }
917
918 comm_lock = 0; // 0 - no lock
919 num_requests = 0;
920}
921
922template <class T>
923void GroupCommunicator::ReduceBegin(const T *ldata) const
924{
925 MFEM_VERIFY(comm_lock == 0, "object is already in use");
926
927 if (group_buf_size == 0) { return; }
928
929 int request_counter = 0;
931 T *buf = (T *)group_buf.GetData();
932 switch (mode)
933 {
934 case byGroup: // ***** Communication by groups *****
935 {
936 for (int gr = 1; gr < group_ldof.Size(); gr++)
937 {
938 const int nldofs = group_ldof.RowSize(gr);
939 // ignore groups without dofs
940 if (nldofs == 0) { continue; }
941
942 if (!gtopo.IAmMaster(gr)) // we are not the master
943 {
944 const int layout = 0;
945 CopyGroupToBuffer(ldata, buf, gr, layout);
946 MPI_Isend(buf,
947 nldofs,
950 43822 + gtopo.GetGroupMasterGroup(gr),
951 gtopo.GetComm(),
952 &requests[request_counter]);
953 request_marker[request_counter] = -1; // mark as send request
954 request_counter++;
955 buf += nldofs;
956 }
957 else // we are the master
958 {
959 const int gs = gtopo.GetGroupSize(gr);
960 const int *nbs = gtopo.GetGroup(gr);
961 buf_offsets[gr] = buf - (T *)group_buf.GetData();
962 for (int i = 0; i < gs; i++)
963 {
964 if (nbs[i] != 0)
965 {
966 MPI_Irecv(buf,
967 nldofs,
969 gtopo.GetNeighborRank(nbs[i]),
970 43822 + gtopo.GetGroupMasterGroup(gr),
971 gtopo.GetComm(),
972 &requests[request_counter]);
973 request_marker[request_counter] = gr;
974 request_counter++;
975 buf += nldofs;
976 }
977 }
978 }
979 }
980 break;
981 }
982
983 case byNeighbor: // ***** Communication by neighbors *****
984 {
985 for (int nbr = 1; nbr < nbr_send_groups.Size(); nbr++)
986 {
987 // In Reduce operation: send_groups <--> recv_groups
988 const int num_send_groups = nbr_recv_groups.RowSize(nbr);
989 if (num_send_groups > 0)
990 {
991 T *buf_start = buf;
992 const int *grp_list = nbr_recv_groups.GetRow(nbr);
993 for (int i = 0; i < num_send_groups; i++)
994 {
995 const int layout = 0; // ldata is an array on all ldofs
996 buf = CopyGroupToBuffer(ldata, buf, grp_list[i], layout);
997 }
998 MPI_Isend(buf_start,
999 buf - buf_start,
1002 43822,
1003 gtopo.GetComm(),
1004 &requests[request_counter]);
1005 request_marker[request_counter] = -1; // mark as send request
1006 request_counter++;
1007 }
1008
1009 // In Reduce operation: send_groups <--> recv_groups
1010 const int num_recv_groups = nbr_send_groups.RowSize(nbr);
1011 if (num_recv_groups > 0)
1012 {
1013 const int *grp_list = nbr_send_groups.GetRow(nbr);
1014 int recv_size = 0;
1015 for (int i = 0; i < num_recv_groups; i++)
1016 {
1017 recv_size += group_ldof.RowSize(grp_list[i]);
1018 }
1019 MPI_Irecv(buf,
1020 recv_size,
1023 43822,
1024 gtopo.GetComm(),
1025 &requests[request_counter]);
1026 request_marker[request_counter] = nbr;
1027 request_counter++;
1028 buf_offsets[nbr] = buf - (T*)group_buf.GetData();
1029 buf += recv_size;
1030 }
1031 }
1032 MFEM_ASSERT(buf - (T*)group_buf.GetData() == group_buf_size, "");
1033 break;
1034 }
1035 }
1036
1037 comm_lock = 2;
1038 num_requests = request_counter;
1039}
1040
1041template <class T>
1042void GroupCommunicator::ReduceEnd(T *ldata, int layout,
1043 void (*Op)(OpData<T>)) const
1044{
1045 if (comm_lock == 0) { return; }
1046 // The above also handles the case (group_buf_size == 0).
1047 MFEM_VERIFY(comm_lock == 2, "object is NOT locked for Reduce");
1048
1049 switch (mode)
1050 {
1051 case byGroup: // ***** Communication by groups *****
1052 {
1053 OpData<T> opd;
1054 opd.ldata = ldata;
1055 Array<int> group_num_req(group_ldof.Size());
1056 for (int gr = 1; gr < group_ldof.Size(); gr++)
1057 {
1058 group_num_req[gr] =
1059 gtopo.IAmMaster(gr) ? gtopo.GetGroupSize(gr)-1 : 0;
1060 }
1061 int idx;
1062 while (MPI_Waitany(num_requests, requests, &idx, MPI_STATUS_IGNORE),
1063 idx != MPI_UNDEFINED)
1064 {
1065 int gr = request_marker[idx];
1066 if (gr == -1) { continue; } // skip send requests
1067
1068 // Delay the processing of a group until all receive requests, for
1069 // that group, are done:
1070 if ((--group_num_req[gr]) != 0) { continue; }
1071
1072 opd.nldofs = group_ldof.RowSize(gr);
1073 // groups without dofs are skipped, so here nldofs > 0.
1074
1075 opd.buf = (T *)group_buf.GetData() + buf_offsets[gr];
1076 opd.ldofs = (layout == 0) ?
1078 opd.nb = gtopo.GetGroupSize(gr)-1;
1079 Op(opd);
1080 }
1081 break;
1082 }
1083
1084 case byNeighbor: // ***** Communication by neighbors *****
1085 {
1086 MPI_Waitall(num_requests, requests, MPI_STATUSES_IGNORE);
1087
1088 for (int nbr = 1; nbr < nbr_send_groups.Size(); nbr++)
1089 {
1090 // In Reduce operation: send_groups <--> recv_groups
1091 const int num_recv_groups = nbr_send_groups.RowSize(nbr);
1092 if (num_recv_groups > 0)
1093 {
1094 const int *grp_list = nbr_send_groups.GetRow(nbr);
1095 const T *buf = (T*)group_buf.GetData() + buf_offsets[nbr];
1096 for (int i = 0; i < num_recv_groups; i++)
1097 {
1098 buf = ReduceGroupFromBuffer(buf, ldata, grp_list[i],
1099 layout, Op);
1100 }
1101 }
1102 }
1103 break;
1104 }
1105 }
1106
1107 comm_lock = 0; // 0 - no lock
1108 num_requests = 0;
1109}
1110
1111template <class T>
1113{
1114 if (opd.nb == 1)
1115 {
1116 for (int i = 0; i < opd.nldofs; i++)
1117 {
1118 opd.ldata[opd.ldofs[i]] += opd.buf[i];
1119 }
1120 }
1121 else
1122 {
1123 for (int i = 0; i < opd.nldofs; i++)
1124 {
1125 T data = opd.ldata[opd.ldofs[i]];
1126 for (int j = 0; j < opd.nb; j++)
1127 {
1128 data += opd.buf[j*opd.nldofs+i];
1129 }
1130 opd.ldata[opd.ldofs[i]] = data;
1131 }
1132 }
1133}
1134
1135template <class T>
1137{
1138 for (int i = 0; i < opd.nldofs; i++)
1139 {
1140 T data = opd.ldata[opd.ldofs[i]];
1141 for (int j = 0; j < opd.nb; j++)
1142 {
1143 T b = opd.buf[j*opd.nldofs+i];
1144 if (data > b)
1145 {
1146 data = b;
1147 }
1148 }
1149 opd.ldata[opd.ldofs[i]] = data;
1150 }
1151}
1152
1153template <class T>
1155{
1156 for (int i = 0; i < opd.nldofs; i++)
1157 {
1158 T data = opd.ldata[opd.ldofs[i]];
1159 for (int j = 0; j < opd.nb; j++)
1160 {
1161 T b = opd.buf[j*opd.nldofs+i];
1162 if (data < b)
1163 {
1164 data = b;
1165 }
1166 }
1167 opd.ldata[opd.ldofs[i]] = data;
1168 }
1169}
1170
1171template <class T>
1173{
1174 for (int i = 0; i < opd.nldofs; i++)
1175 {
1176 T data = opd.ldata[opd.ldofs[i]];
1177 for (int j = 0; j < opd.nb; j++)
1178 {
1179 data |= opd.buf[j*opd.nldofs+i];
1180 }
1181 opd.ldata[opd.ldofs[i]] = data;
1182 }
1183}
1184
1185void GroupCommunicator::PrintInfo(std::ostream &os) const
1186{
1187 char c = '\0';
1188 const int tag = 46800;
1189 const int myid = gtopo.MyRank();
1190
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; // for mode == byNeighbor
1195 switch (mode)
1196 {
1197 case byGroup:
1198 for (int gr = 1; gr < group_ldof.Size(); gr++)
1199 {
1200 const int nldofs = group_ldof.RowSize(gr);
1201 if (nldofs == 0)
1202 {
1203 num_empty_groups++;
1204 continue;
1205 }
1206 if (gtopo.IAmMaster(gr))
1207 {
1208 num_sends += (gtopo.GetGroupSize(gr)-1);
1209 mem_sends += sizeof(double)*nldofs*(gtopo.GetGroupSize(gr)-1);
1210 num_master_groups++;
1211 }
1212 else
1213 {
1214 num_recvs++;
1215 mem_recvs += sizeof(double)*nldofs;
1216 }
1217 }
1218 break;
1219
1220 case byNeighbor:
1221 for (int gr = 1; gr < group_ldof.Size(); gr++)
1222 {
1223 const int nldofs = group_ldof.RowSize(gr);
1224 if (nldofs == 0)
1225 {
1226 num_empty_groups++;
1227 continue;
1228 }
1229 if (gtopo.IAmMaster(gr))
1230 {
1231 num_master_groups++;
1232 }
1233 }
1234 for (int nbr = 1; nbr < nbr_send_groups.Size(); nbr++)
1235 {
1236 const int num_send_groups = nbr_send_groups.RowSize(nbr);
1237 if (num_send_groups > 0)
1238 {
1239 const int *grp_list = nbr_send_groups.GetRow(nbr);
1240 for (int i = 0; i < num_send_groups; i++)
1241 {
1242 mem_sends += sizeof(double)*group_ldof.RowSize(grp_list[i]);
1243 }
1244 num_sends++;
1245 }
1246
1247 const int num_recv_groups = nbr_recv_groups.RowSize(nbr);
1248 if (num_recv_groups > 0)
1249 {
1250 const int *grp_list = nbr_recv_groups.GetRow(nbr);
1251 for (int i = 0; i < num_recv_groups; i++)
1252 {
1253 mem_recvs += sizeof(double)*group_ldof.RowSize(grp_list[i]);
1254 }
1255 num_recvs++;
1256 }
1257 if (num_send_groups > 0 || num_recv_groups > 0)
1258 {
1259 num_active_neighbors++;
1260 }
1261 }
1262 break;
1263 }
1264 if (myid != 0)
1265 {
1266 MPI_Recv(&c, 1, MPI_CHAR, myid-1, tag, gtopo.GetComm(),
1267 MPI_STATUS_IGNORE);
1268 }
1269 else
1270 {
1271 os << "\nGroupCommunicator:\n";
1272 }
1273 os << "Rank " << myid << ":\n"
1274 " mode = " <<
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";
1280 os <<
1281 " num groups = " << group_ldof.Size() << " = " <<
1282 num_master_groups << " + " <<
1283 group_ldof.Size()-num_master_groups-num_empty_groups << " + " <<
1284 num_empty_groups << " (master + slave + empty)\n";
1285 if (mode == byNeighbor)
1286 {
1287 os <<
1288 " num neighbors = " << nbr_send_groups.Size() << " = " <<
1289 num_active_neighbors << " + " <<
1290 nbr_send_groups.Size()-num_active_neighbors <<
1291 " (active + inactive)\n";
1292 }
1293 if (myid != gtopo.NRanks()-1)
1294 {
1295 os << std::flush;
1296 MPI_Send(&c, 1, MPI_CHAR, myid+1, tag, gtopo.GetComm());
1297 }
1298 else
1299 {
1300 os << std::endl;
1301 }
1302 MPI_Barrier(gtopo.GetComm());
1303}
1304
1306{
1307 delete [] buf_offsets;
1308 delete [] request_marker;
1309 // delete [] statuses;
1310 delete [] requests;
1311}
1312
1313// @cond DOXYGEN_SKIP
1314
1315// instantiate GroupCommunicator::Bcast and Reduce for int and double
1316template void GroupCommunicator::BcastBegin<int>(int *, int) const;
1317template void GroupCommunicator::BcastEnd<int>(int *, int) const;
1318template void GroupCommunicator::ReduceBegin<int>(const int *) const;
1320 int *, int, void (*)(OpData<int>)) const;
1321
1322template void GroupCommunicator::BcastBegin<double>(double *, int) const;
1323template void GroupCommunicator::BcastEnd<double>(double *, int) const;
1324template void GroupCommunicator::ReduceBegin<double>(const double *) const;
1326 double *, int, void (*)(OpData<double>)) const;
1327
1328template void GroupCommunicator::BcastBegin<float>(float *, int) const;
1329template void GroupCommunicator::BcastEnd<float>(float *, int) const;
1330template void GroupCommunicator::ReduceBegin<float>(const float *) const;
1332 float *, int, void (*)(OpData<float>)) const;
1333
1334// @endcond
1335
1336// instantiate reduce operators for int and double
1337template void GroupCommunicator::Sum<int>(OpData<int>);
1338template void GroupCommunicator::Min<int>(OpData<int>);
1339template void GroupCommunicator::Max<int>(OpData<int>);
1340template void GroupCommunicator::BitOR<int>(OpData<int>);
1341
1342template void GroupCommunicator::Sum<double>(OpData<double>);
1343template void GroupCommunicator::Min<double>(OpData<double>);
1344template void GroupCommunicator::Max<double>(OpData<double>);
1345
1346template void GroupCommunicator::Sum<float>(OpData<float>);
1347template void GroupCommunicator::Min<float>(OpData<float>);
1348template void GroupCommunicator::Max<float>(OpData<float>);
1349
1350
1351#ifdef __bgq__
1352static void DebugRankCoords(int** coords, int dim, int size)
1353{
1354 for (int i = 0; i < size; i++)
1355 {
1356 mfem::out << "Rank " << i << " coords: ";
1357 for (int j = 0; j < dim; j++)
1358 {
1359 mfem::out << coords[i][j] << " ";
1360 }
1361 mfem::out << endl;
1362 }
1363}
1364
1365struct CompareCoords
1366{
1367 CompareCoords(int coord) : coord(coord) {}
1368 int coord;
1369
1370 bool operator()(int* const &a, int* const &b) const
1371 { return a[coord] < b[coord]; }
1372};
1373
1374void KdTreeSort(int** coords, int d, int dim, int size)
1375{
1376 if (size > 1)
1377 {
1378 bool all_same = true;
1379 for (int i = 1; i < size && all_same; i++)
1380 {
1381 for (int j = 0; j < dim; j++)
1382 {
1383 if (coords[i][j] != coords[0][j]) { all_same = false; break; }
1384 }
1385 }
1386 if (all_same) { return; }
1387
1388 // sort by coordinate 'd'
1389 std::sort(coords, coords + size, CompareCoords(d));
1390 int next = (d + 1) % dim;
1391
1392 if (coords[0][d] < coords[size-1][d])
1393 {
1394 KdTreeSort(coords, next, dim, size/2);
1395 KdTreeSort(coords + size/2, next, dim, size - size/2);
1396 }
1397 else
1398 {
1399 // skip constant dimension
1400 KdTreeSort(coords, next, dim, size);
1401 }
1402 }
1403}
1404
1405MPI_Comm ReorderRanksZCurve(MPI_Comm comm)
1406{
1407 MPI_Status status;
1408
1409 int rank, size;
1410 MPI_Comm_rank(comm, &rank);
1411 MPI_Comm_size(comm, &size);
1412
1413 int dim;
1414 MPIX_Torus_ndims(&dim);
1415
1416 int* mycoords = new int[dim + 1];
1417 MPIX_Rank2torus(rank, mycoords);
1418
1419 MPI_Send(mycoords, dim, MPI_INT, 0, 111, comm);
1420 delete [] mycoords;
1421
1422 if (rank == 0)
1423 {
1424 int** coords = new int*[size];
1425 for (int i = 0; i < size; i++)
1426 {
1427 coords[i] = new int[dim + 1];
1428 coords[i][dim] = i;
1429 MPI_Recv(coords[i], dim, MPI_INT, i, 111, comm, &status);
1430 }
1431
1432 KdTreeSort(coords, 0, dim, size);
1433
1434 // DebugRankCoords(coords, dim, size);
1435
1436 for (int i = 0; i < size; i++)
1437 {
1438 MPI_Send(&coords[i][dim], 1, MPI_INT, i, 112, comm);
1439 delete [] coords[i];
1440 }
1441 delete [] coords;
1442 }
1443
1444 int new_rank;
1445 MPI_Recv(&new_rank, 1, MPI_INT, 0, 112, comm, &status);
1446
1447 MPI_Comm new_comm;
1448 MPI_Comm_split(comm, 0, new_rank, &new_comm);
1449 return new_comm;
1450}
1451
1452#else // __bgq__
1453
1454MPI_Comm ReorderRanksZCurve(MPI_Comm comm)
1455{
1456 // pass
1457 return comm;
1458}
1459#endif // __bgq__
1460
1461} // namespace mfem
1462
1463#endif
void Sort()
Sorts the array in ascending order. This requires operator< to be defined for T.
Definition array.hpp:312
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition array.hpp:840
int Size() const
Return the logical size of the array.
Definition array.hpp:166
int Append(const T &el)
Append element 'el' to array, resize if necessary.
Definition array.hpp:912
T * GetData()
Returns the data.
Definition array.hpp:140
void Copy(Array &copy) const
Create a copy of the internal array to the provided copy.
Definition array.hpp:1042
void Unique()
Removes duplicities from a sorted array. This requires operator== to be defined for T.
Definition array.hpp:320
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'.
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:33
List of integer sets.
Definition sets.hpp:51
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...
Definition sets.cpp:56
void AsTable(Table &t) const
Write the list of sets into table 't'.
Definition sets.cpp:81
int PickElementInSet(int i) const
Return the value of the first element of the ith set.
Definition sets.hpp:61
int Lookup(const IntegerSet &s) const
Definition sets.cpp:69
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,...
Definition table.hpp:43
int * GetJ()
Definition table.hpp:128
void AddConnections(int r, const int *c, int nc)
Definition table.cpp:152
int RowSize(int i) const
Definition table.hpp:122
void ShiftUpI()
Definition table.cpp:163
void GetRow(int i, Array< int > &row) const
Return row i in array row (the Table must be finalized)
Definition table.cpp:233
void AddConnection(int r, int c)
Definition table.hpp:89
void MakeI(int nrows)
Definition table.cpp:130
int Size() const
Returns the number of TYPE I elements.
Definition table.hpp:103
int Size_of_connections() const
Returns the number of connections in the table.
Definition table.hpp:110
void AddColumnsInRow(int r, int ncol)
Definition table.hpp:87
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 ...
Definition table.cpp:322
void Copy(Table &copy) const
Definition table.cpp:427
void MakeJ()
Definition table.cpp:140
int * GetI()
Definition table.hpp:127
void AddAColumnInRow(int r)
Definition table.hpp:86
void SetDims(int rows, int nnz)
Set the rows and the number of all connections for the table.
Definition table.cpp:188
int dim
Definition ex24.cpp:53
int index(int i, int j, int nx, int ny)
Definition life.cpp:236
real_t b
Definition lissajous.cpp:42
real_t a
Definition lissajous.cpp:41
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 Swap(T &a, T &b)
Swap objects of type T. The operation is performed using the most specialized swap function from the ...
Definition array.hpp:738
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)
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.