version 3.11-dev
Loading...
Searching...
No Matches
parallelhelpers.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_LINEAR_PARALLELHELPERS_HH
13#define DUMUX_LINEAR_PARALLELHELPERS_HH
14
15#include <algorithm>
16#include <dune/common/exceptions.hh>
17#include <dune/geometry/dimension.hh>
18#include <dune/grid/common/datahandleif.hh>
19#include <dune/grid/common/partitionset.hh>
20#include <dune/istl/owneroverlapcopy.hh>
21#include <dune/istl/paamg/pinfo.hh>
22#include <dune/istl/bvector.hh>
23#include <dune/istl/multitypeblockvector.hh>
28
29#include <string>
30#include <utility>
31
32namespace Dumux::Detail {
33
34template<class LinearSolverTraits, bool canCommunicate = false>
36
37template<class LinearSolverTraits>
39{
40 using GridView = typename LinearSolverTraits::GridView;
41 using DofMapper = typename LinearSolverTraits::DofMapper;
42
43 // TODO: this is some large number (replace by limits?)
44 static constexpr std::size_t ghostMarker_ = 1<<24;
45
46 class BaseGatherScatter
47 {
48 static constexpr int numCodims_ = GridView::dimension + 1;
49
50 static constexpr auto getActiveCodims_()
51 {
52 if constexpr (requires { LinearSolverTraits::dofCodims; })
53 return LinearSolverTraits::dofCodims;
54 else
55 return std::bitset<numCodims_>{ 1UL << LinearSolverTraits::dofCodim };
56 }
57
58 public:
59 BaseGatherScatter(const DofMapper& mapper)
60 : mapper_(mapper) {}
61
62 template<class EntityType>
63 int index(const EntityType& e) const
64 { return mapper_.index(e); }
65
66 template<class EntityType>
67 auto indices(const EntityType& e) const
68 { return asMultiMapper(mapper_).indices(e); }
69
70 bool contains(int dim, int codim) const
71 { return codim >= 0 && codim < numCodims_ && activeCodims_.test(codim); }
72
75 bool fixedSize(int dim, int codim) const
76 {
77 if (codim < 0 || codim >= numCodims_ || !activeCodims_.test(codim)) return true;
78 const auto& gtypes = mapper_.types(codim);
79 if (gtypes.empty()) return true;
80 const auto nDofs = mapper_.size(gtypes[0]);
81 for (const auto& gt : gtypes)
82 if (mapper_.size(gt) != nDofs) return false;
83 return true;
84 }
85
87 template<class EntityType>
88 std::size_t size(EntityType& e) const
89 { return asMultiMapper(mapper_).indices(e).size(); }
90
91 template<class EntityType>
92 bool isNeitherInteriorNorBorderEntity(EntityType& e) const
93 { return e.partitionType() != Dune::InteriorEntity && e.partitionType() != Dune::BorderEntity; }
94
95 protected:
96 const DofMapper& mapper_;
97 private:
98 static constexpr auto activeCodims_ = getActiveCodims_();
99 };
100
107 class GhostGatherScatter
108 : public BaseGatherScatter
109 , public Dune::CommDataHandleIF<GhostGatherScatter, std::size_t>
110 {
111 public:
112 using DataType = std::size_t;
113 using BaseGatherScatter::contains;
114 using BaseGatherScatter::fixedSize;
115 using BaseGatherScatter::size;
116
117 GhostGatherScatter(std::vector<std::size_t>& ranks, const DofMapper& mapper)
118 : BaseGatherScatter(mapper)
119 , ranks_(ranks)
120 {}
121
122 template<class MessageBuffer, class EntityType>
123 void gather(MessageBuffer& buff, const EntityType& e) const
124 {
125 const bool isGhost = this->isNeitherInteriorNorBorderEntity(e);
126 for (auto i : this->indices(e)) {
127 auto& data = ranks_[i];
128 if (isGhost) data = ghostMarker_;
129 buff.write(data);
130 }
131 }
132
133 template<class MessageBuffer, class EntityType>
134 void scatter(MessageBuffer& buff, const EntityType& e, std::size_t n)
135 {
136 const bool isGhost = this->isNeitherInteriorNorBorderEntity(e);
137 for (auto i : this->indices(e)) {
138 std::size_t x;
139 buff.read(x);
140 if (isGhost) ranks_[i] = ghostMarker_;
141 }
142 }
143 private:
144 std::vector<std::size_t>& ranks_;
145 };
146
153 class InteriorBorderGatherScatter
154 : public BaseGatherScatter
155 , public Dune::CommDataHandleIF<InteriorBorderGatherScatter, std::size_t>
156 {
157 public:
158 using DataType = std::size_t;
159 using BaseGatherScatter::contains;
160 using BaseGatherScatter::fixedSize;
161 using BaseGatherScatter::size;
162
163 InteriorBorderGatherScatter(std::vector<std::size_t>& ranks, const DofMapper& mapper)
164 : BaseGatherScatter(mapper)
165 , ranks_(ranks)
166 {}
167
171 template<class MessageBuffer, class EntityType>
172 void gather(MessageBuffer& buff, const EntityType& e) const
173 {
174 const bool isGhost = this->isNeitherInteriorNorBorderEntity(e);
175 for (auto i : this->indices(e)) {
176 auto& data = ranks_[i];
177 if (isGhost) data = ghostMarker_;
178 buff.write(data);
179 }
180 }
181
185 template<class MessageBuffer, class EntityType>
186 void scatter(MessageBuffer& buff, const EntityType& e, std::size_t n)
187 {
188 const bool isGhost = this->isNeitherInteriorNorBorderEntity(e);
189 for (auto i : this->indices(e)) {
190 std::size_t x;
191 buff.read(x);
192 auto& data = ranks_[i];
193 // we leave ghost unchanged
194 // for other dofs, the process with the lowest rank
195 // is assigned to be the (unique) owner
196 using std::min;
197 data = isGhost ? x : min(data, x);
198 }
199 }
200 private:
201 std::vector<std::size_t>& ranks_;
202 };
203
208 struct NeighbourGatherScatter
209 : public BaseGatherScatter
210 , public Dune::CommDataHandleIF<NeighbourGatherScatter, int>
211 {
212 using DataType = int;
213 using BaseGatherScatter::contains;
214 using BaseGatherScatter::fixedSize;
215 using BaseGatherScatter::size;
216
217 NeighbourGatherScatter(const DofMapper& mapper, int rank, std::set<int>& neighbours)
218 : BaseGatherScatter(mapper)
219 , rank_(rank)
220 , neighbours_(neighbours)
221 {}
222
223 // NeighbourGatherScatter sends exactly 1 rank value per entity (not per DOF).
224 template<class EntityType>
225 std::size_t size(EntityType& e) const { return 1; }
226 bool fixedSize(int dim, int codim) const { return true; }
227
228 template<class MessageBuffer, class EntityType>
229 void gather(MessageBuffer& buff, const EntityType& e) const
230 {
231 buff.write(rank_);
232 }
233
234 template<class MessageBuffer, class EntityType>
235 void scatter(MessageBuffer& buff, const EntityType& e, std::size_t n)
236 {
237 int x;
238 buff.read(x);
239 neighbours_.insert(x);
240 }
241 private:
242 int rank_;
243 std::set<int>& neighbours_;
244 };
245
246
251 struct SharedGatherScatter
252 : public BaseGatherScatter
253 , public Dune::CommDataHandleIF<SharedGatherScatter, int>
254 {
255 using DataType = int;
256 using BaseGatherScatter::contains;
257 using BaseGatherScatter::fixedSize;
258 using BaseGatherScatter::size;
259
260 SharedGatherScatter(std::vector<int>& shared, const DofMapper& mapper)
261 : BaseGatherScatter(mapper)
262 , shared_(shared)
263 {}
264
265 template<class MessageBuffer, class EntityType>
266 void gather(MessageBuffer& buff, EntityType& e) const
267 {
268 for (std::size_t k = 0; k < this->indices(e).size(); ++k)
269 buff.write(int(true));
270 }
271
272 template<class MessageBuffer, class EntityType>
273 void scatter(MessageBuffer& buff, const EntityType &e, std::size_t n)
274 {
275 for (auto i : this->indices(e)) {
276 int x;
277 buff.read(x);
278 shared_[i] = shared_[i] || x;
279 }
280 }
281 private:
282 std::vector<int>& shared_;
283 };
284
289 template<typename GlobalIndex>
290 struct GlobalIndexGatherScatter
291 : public BaseGatherScatter
292 , public Dune::CommDataHandleIF<GlobalIndexGatherScatter<GlobalIndex>, GlobalIndex>
293 {
294 using DataType = GlobalIndex;
295 using BaseGatherScatter::contains;
296 using BaseGatherScatter::fixedSize;
297 using BaseGatherScatter::size;
298
299 GlobalIndexGatherScatter(std::vector<GlobalIndex>& globalIndices, const DofMapper& mapper)
300 : BaseGatherScatter(mapper)
301 , globalIndices_(globalIndices)
302 {}
303
304 template<class MessageBuffer, class EntityType>
305 void gather(MessageBuffer& buff, const EntityType& e) const
306 {
307 for (auto i : this->indices(e))
308 buff.write(globalIndices_[i]);
309 }
310
311 template<class MessageBuffer, class EntityType>
312 void scatter(MessageBuffer& buff, const EntityType& e, std::size_t n)
313 {
314 for (auto i : this->indices(e)) {
315 DataType x;
316 buff.read(x);
317 using std::min;
318 globalIndices_[i] = min(globalIndices_[i], x);
319 }
320 }
321 private:
322 std::vector<GlobalIndex>& globalIndices_;
323 };
324
325public:
326
327 ParallelISTLHelperImpl(const GridView& gridView, const DofMapper& mapper)
328 : gridView_(gridView), mapper_(mapper)
329 {
330 if constexpr (LinearSolverTraits::canCommunicate)
331 initGhostsAndOwners_();
332 else
333 DUNE_THROW(Dune::InvalidStateException,
334 "Cannot initialize parallel helper for a grid that cannot communicate."
335 );
336 }
337
338 bool isGhost(std::size_t i) const
339 { return isGhost_[i] == ghostMarker_; }
340
341 bool isOwned(std::size_t i) const
342 { return isOwned_[i] == 1; }
343
355 template<class Comm>
356 void createParallelIndexSet(Comm& comm) const
357 {
358 if constexpr (LinearSolverTraits::canCommunicate)
359 {
360 if (gridView_.comm().size() <= 1)
361 {
362 comm.remoteIndices().template rebuild<false>();
363 return;
364 }
365
366 // First find out which dofs we share with other processors
367 std::vector<int> isShared(mapper_.size(), false);
368 SharedGatherScatter sgs(isShared, mapper_);
369 gridView_.communicate(sgs, Dune::All_All_Interface, Dune::ForwardCommunication);
370
371 // Count shared dofs that we own
372 using GlobalIndex = typename Comm::ParallelIndexSet::GlobalIndex;
373 GlobalIndex count = 0;
374 for (std::size_t i = 0; i < isShared.size(); ++i)
375 if (isShared[i] && isOwned_[i] == 1)
376 ++count;
377
378 std::vector<GlobalIndex> counts(gridView_.comm().size());
379 gridView_.comm().allgather(&count, 1, counts.data());
380
381 // Compute start index start_p = \sum_{i=0}^{i<p} counts_i
382 const int rank = gridView_.comm().rank();
383 auto start = std::accumulate(counts.begin(), counts.begin() + rank, GlobalIndex(0));
384
385 // starting from start (offset), number global indices owned by this process consecutively
386 std::vector<GlobalIndex> globalIndices(mapper_.size(), std::numeric_limits<GlobalIndex>::max());
387 for (std::size_t i = 0; i < globalIndices.size(); ++i)
388 {
389 if (isOwned_[i] == 1 && isShared[i])
390 {
391 globalIndices[i] = start;
392 ++start;
393 }
394 }
395
396 // publish global indices for the shared DOFS to other processors.
397 GlobalIndexGatherScatter<GlobalIndex> gigs(globalIndices, mapper_);
398 gridView_.communicate(gigs, Dune::All_All_Interface, Dune::ForwardCommunication);
399
400 // put the information into the parallel index set
401 resizeAndFillIndexSet_(comm, globalIndices);
402
403 // Compute neighbours using communication
404 std::set<int> neighbours;
405 NeighbourGatherScatter ngs(mapper_, gridView_.comm().rank(), neighbours);
406 gridView_.communicate(ngs, Dune::All_All_Interface, Dune::ForwardCommunication);
407 comm.remoteIndices().setNeighbours(neighbours);
408
409 comm.remoteIndices().template rebuild<false>();
410 }
411 else
412 DUNE_THROW(Dune::InvalidStateException,
413 "Cannot build parallel index set for a grid that cannot communicate."
414 );
415 }
416
418 const DofMapper& dofMapper() const
419 { return mapper_; }
420
422 const GridView& gridView() const
423 { return gridView_; }
424
425private:
434 void initGhostsAndOwners_()
435 {
436 const auto rank = gridView_.comm().rank();
437 isOwned_.assign(mapper_.size(), rank);
438
439 // find out about ghosts
440 GhostGatherScatter ggs(isOwned_, mapper_);
441
442 if (gridView_.comm().size() > 1)
443 gridView_.communicate(ggs, Dune::InteriorBorder_All_Interface, Dune::ForwardCommunication);
444
445 // isGhost_ contains rank when the dof is not on a ghost entity (interior or border)
446 // and ghostMarker_ when the dof is on a ghost
447 isGhost_ = isOwned_;
448
449 // partition interior/border uniquely
450 // after the communication each dof is assigned to one unique owner process rank
451 InteriorBorderGatherScatter dh(isOwned_, mapper_);
452
453 if (gridView_.comm().size() > 1)
454 gridView_.communicate(dh, Dune::InteriorBorder_InteriorBorder_Interface, Dune::ForwardCommunication);
455
456 // convert vector into mask vector
457 for (auto& v : isOwned_)
458 v = (v == rank) ? 1 : 0;
459 }
460
461 template<class Comm, class GlobalIndices>
462 void resizeAndFillIndexSet_(Comm& comm, const GlobalIndices& globalIndices) const
463 {
464 comm.indexSet().beginResize();
465
466 // add triplets characterizing each dof in the parallel index set
467 // (globalIndex, (localIndex, attribute)) where attribute is owner, overlap or copy
468 for (std::size_t localIndex = 0; localIndex < globalIndices.size(); ++localIndex)
469 {
470 const auto globalIndex = globalIndices[localIndex];
471 if (globalIndex != std::numeric_limits<typename GlobalIndices::value_type>::max())
472 {
473 const bool isOwned = isOwned_[localIndex] > 0;
474 const auto attr = getAttribute_(comm, isOwned, isGhost(localIndex));
475 using LocalIndex = typename Comm::ParallelIndexSet::LocalIndex;
476 comm.indexSet().add(globalIndex, LocalIndex{localIndex, attr});
477 }
478 }
479
480 comm.indexSet().endResize();
481 }
482
483 template<class Comm>
484 Dune::OwnerOverlapCopyAttributeSet::AttributeSet
485 getAttribute_(const Comm& comm, const bool isOwned, const bool isGhost) const
486 {
487 if (isOwned) // this process owns the dof (uniquely)
488 return Dune::OwnerOverlapCopyAttributeSet::owner;
489 else if (isGhost && (comm.category() == Dune::SolverCategory::nonoverlapping) )
490 return Dune::OwnerOverlapCopyAttributeSet::overlap;
491 else
492 return Dune::OwnerOverlapCopyAttributeSet::copy;
493 }
494
495 const GridView gridView_;
496 const DofMapper& mapper_;
498 std::vector<std::size_t> isOwned_;
500 std::vector<std::size_t> isGhost_;
501
502};
503
504} // end namespace Dumux::Detail
505
506namespace Dumux {
507
513template<class LinearSolverTraits>
516 LinearSolverTraits, LinearSolverTraits::canCommunicate
517 >;
518
519
520template<class GridView, class DofMapper, int dofCodim>
522{
523public:
524 ParallelVectorHelper(const GridView& gridView, const DofMapper& mapper)
525 : gridView_(gridView), mapper_(mapper)
526 {}
527
529 template<class Block, class Alloc>
530 void makeNonOverlappingConsistent(Dune::BlockVector<Block, Alloc>& v) const
531 {
532 if constexpr (Dune::Capabilities::canCommunicate<typename GridView::Traits::Grid, dofCodim>::v)
533 {
535 if (gridView_.comm().size() > 1)
536 gridView_.communicate(gs, Dune::InteriorBorder_InteriorBorder_Interface,
537 Dune::ForwardCommunication);
538 }
539 else
540 DUNE_THROW(Dune::InvalidStateException, "Cannot call makeNonOverlappingConsistent for a grid that cannot communicate codim-" << dofCodim << "-entities.");
541 }
542
544 template<class Block, class Alloc>
545 void makeOverlappingConsistent(Dune::BlockVector<Block, Alloc>& v) const
546 {
547 if constexpr (Dune::Capabilities::canCommunicate<typename GridView::Traits::Grid, dofCodim>::v)
548 {
549 // Iterate over the gridView_ to make overlap entries 0. This is necessary to make
550 // `VectorCommDataHandleSum` work properly.
551 for (const auto& element : elements(gridView_))
552 if (element.partitionType() == Dune::OverlapEntity)
553 v[gridView_.indexSet().index(element)] = 0;
554
555 if (gridView_.comm().size() > 1)
556 {
558 gridView_.communicate(gs, Dune::Overlap_All_Interface,
559 Dune::ForwardCommunication);
560 }
561 }
562 else
563 DUNE_THROW(Dune::InvalidStateException, "Cannot call makeNonOverlappingConsistent for a grid that cannot communicate codim-" << dofCodim << "-entities.");
564 }
565
567 template<class... Blocks>
569 {
570 DUNE_THROW(Dune::NotImplemented, "makeNonOverlappingConsistent for Dune::MultiTypeBlockVector");
571 }
572
573private:
574 const GridView gridView_;
575 const DofMapper& mapper_;
576};
577
578template<class GridView, class DofMapper>
580{
581public:
582 MultiCodimParallelVectorHelper(const GridView& gridView, const DofMapper& mapper)
583 : gridView_(gridView), mapper_(mapper)
584 {}
585
587 template<class Block, class Alloc, std::size_t numCodims>
589 Dune::BlockVector<Block, Alloc>& v,
590 const std::bitset<numCodims>& activeCodims
591 ) const
592 {
593 static_assert(numCodims == GridView::dimension + 1,
594 "Number of codims must match GridView::dimension+1");
595
596 if (gridView_.comm().size() > 1)
597 {
599 mapper_, v, activeCodims
600 );
601 gridView_.communicate(gs, Dune::InteriorBorder_InteriorBorder_Interface, Dune::ForwardCommunication);
602 }
603 }
604
606 template<class Block, class Alloc, std::size_t numCodims>
608 Dune::BlockVector<Block, Alloc>& v,
609 const std::bitset<numCodims>& activeCodims
610 ) const
611 {
612 static_assert(numCodims == GridView::dimension + 1,
613 "Number of codims must match GridView::dimension+1");
614
615 zeroOverlapEntries_<0>(v, activeCodims);
616
617 if (gridView_.comm().size() > 1)
618 {
620 mapper_, v, activeCodims
621 );
622 gridView_.communicate(gs, Dune::Overlap_All_Interface, Dune::ForwardCommunication);
623 }
624 }
625
626private:
627 template<int codim, class Block, class Alloc, std::size_t numCodims>
628 void zeroOverlapEntries_(Dune::BlockVector<Block, Alloc>& v,
629 const std::bitset<numCodims>& activeCodims) const
630 {
631 if constexpr (codim <= GridView::dimension)
632 {
633 if (activeCodims.test(codim))
634 {
635 for (const auto& element : elements(gridView_))
636 {
637 for (int i = 0; i < element.subEntities(codim); ++i)
638 {
639 const auto entity = element.template subEntity<codim>(i);
640 if (entity.partitionType() == Dune::OverlapEntity)
641 for (auto i : asMultiMapper(mapper_).indices(entity))
642 v[i] = 0;
643 }
644 }
645 }
646
647 zeroOverlapEntries_<codim + 1>(v, activeCodims);
648 }
649 }
650
651 const GridView gridView_;
652 const DofMapper& mapper_;
653};
654
663template<class Matrix, class GridView,
664 class RowDofMapper, int rowDofCodim>
666{
667 using IdType = typename GridView::Traits::Grid::Traits::GlobalIdSet::IdType;
668
669 struct DofIdKey
670 {
671 int codim;
672 IdType id;
673
674 friend bool operator<(const DofIdKey& a, const DofIdKey& b)
675 {
676 if (a.codim != b.codim)
677 return a.codim < b.codim;
678 return a.id < b.id;
679 }
680 };
681
703 template<class ColIsGhostFunc>
704 struct MatrixPatternExchange
705 : public Dune::CommDataHandleIF<MatrixPatternExchange<ColIsGhostFunc>, DofIdKey>
706 {
708 using DataType = DofIdKey;
709
720 MatrixPatternExchange(const RowDofMapper& rowEntityMapper,
721 const std::map<DofIdKey,int>& globalToLocal,
722 const std::map<int,DofIdKey>& localToGlobal,
723 Matrix& A,
724 std::vector<std::set<int>>& sparsityPattern,
725 const ColIsGhostFunc& isGhostColumDof)
726 : rowEntityMapper_(rowEntityMapper)
727 , idToIndex_(globalToLocal), indexToID_(localToGlobal)
728 , sparsityPattern_(sparsityPattern), A_(A), isGhostColumDof_(isGhostColumDof)
729 {
730 sparsityPattern_.resize(A.N());
731 }
732
736 bool contains (int dim, int codim) const
737 { return (codim == rowDofCodim); }
738
740 bool fixedSize(int dim, int codim) const
741 { return false; }
742
746 template<class EntityType>
747 std::size_t size(EntityType& e) const
748 {
749 const auto rowIdx = rowEntityMapper_.index(e);
750
751 // all column entity indices of this row that are in the index set
752 std::size_t n = 0;
753 for (auto colIt = A_[rowIdx].begin(); colIt != A_[rowIdx].end(); ++colIt)
754 {
755 if (indexToID_.count(colIt.index()))
756 ++n;
757 }
758
759 return n;
760 }
761
765 template<class MessageBuffer, class EntityType>
766 void gather(MessageBuffer& buff, const EntityType& e) const
767 {
768 const auto rowIdx = rowEntityMapper_.index(e);
769
770 // send all column entity ids of this row that are in the index set
771 for (auto colIt = A_[rowIdx].begin(); colIt != A_[rowIdx].end(); ++colIt)
772 {
773 if (auto it = indexToID_.find(colIt.index()); it != indexToID_.end())
774 buff.write(it->second);
775 }
776 }
777
781 template<class MessageBuffer, class EntityType>
782 void scatter(MessageBuffer& buff, const EntityType& e, std::size_t n)
783 {
784 const auto rowIdx = rowEntityMapper_.index(e);
785 for (std::size_t k = 0; k < n; k++)
786 {
787 // receive n column entity IDs
788 DofIdKey key;
789 buff.read(key);
790
791 // only add entries that are contained in the index set
792 if (const auto it = idToIndex_.find(key); it != idToIndex_.end())
793 {
794 const auto colIdx = it->second;
795 if (!isGhostColumDof_(colIdx))
796 sparsityPattern_[rowIdx].insert(colIdx);
797 }
798 }
799 }
800
801 private:
802 const RowDofMapper& rowEntityMapper_;
803 const std::map<DofIdKey, int>& idToIndex_;
804 const std::map<int, DofIdKey>& indexToID_;
805 std::vector<std::set<int>>& sparsityPattern_;
806 Matrix& A_;
807 const ColIsGhostFunc& isGhostColumDof_;
808
809 }; // class MatrixPatternExchange
810
812 struct MatrixEntry
813 {
814 DofIdKey first;
815 typename Matrix::block_type second;
816 };
817
819 struct MatrixEntryExchange
820 : public Dune::CommDataHandleIF<MatrixEntryExchange, MatrixEntry>
821 {
823 using DataType = MatrixEntry;
824
825 MatrixEntryExchange(const RowDofMapper& rowEntityMapper,
826 const std::map<DofIdKey, int>& globalToLocal,
827 const std::map<int, DofIdKey>& localToGlobal,
828 Matrix& A)
829 : rowEntityMapper_(rowEntityMapper)
830 , idToIndex_(globalToLocal), indexToID_(localToGlobal)
831 , A_(A)
832 {}
833
837 bool contains(int dim, int codim) const
838 { return (codim == rowDofCodim); }
839
841 bool fixedSize(int dim, int codim) const
842 { return false; }
843
847 template<class EntityType>
848 std::size_t size(EntityType& e) const
849 {
850 const auto rowIdx = rowEntityMapper_.index(e);
851 std::size_t n = 0;
852 for (auto colIt = A_[rowIdx].begin(); colIt != A_[rowIdx].end(); ++colIt)
853 {
854 if (indexToID_.count(colIt.index()))
855 ++n;
856 }
857
858 return n;
859 }
860
864 template<class MessageBuffer, class EntityType>
865 void gather(MessageBuffer& buff, const EntityType& e) const
866 {
867 const auto rowIdx = rowEntityMapper_.index(e);
868 // send all matrix entries for which the column index is in the index set
869 for (auto colIt = A_[rowIdx].begin(); colIt != A_[rowIdx].end(); ++colIt)
870 {
871 if (auto it = indexToID_.find(colIt.index()); it != indexToID_.end())
872 buff.write(MatrixEntry{ it->second, *colIt });
873 }
874 }
875
879 template<class MessageBuffer, class EntityType>
880 void scatter(MessageBuffer& buff, const EntityType& e, std::size_t n)
881 {
882 const auto rowIdx = rowEntityMapper_.index(e);
883 for (std::size_t k = 0; k < n; k++)
884 {
885 MatrixEntry m;
886 buff.read(m);
887 const auto& [colDofID, matrixBlock] = m; // unpack
888
889 // only add entries in the index set and for which A has allocated memory
890 if (auto it = idToIndex_.find(colDofID); it != idToIndex_.end())
891 {
892 if (A_[rowIdx].find(it->second) != A_[rowIdx].end())
893 A_[rowIdx][it->second] += matrixBlock;
894 }
895 }
896 }
897
898 private:
899 const RowDofMapper& rowEntityMapper_;
900 const std::map<DofIdKey, int>& idToIndex_;
901 const std::map<int, DofIdKey>& indexToID_;
902 Matrix& A_;
903
904 }; // class MatrixEntryExchange
905
906public:
907
908 ParallelMatrixHelper(const GridView& gridView, const RowDofMapper& mapper)
909 : gridView_(gridView), mapper_(mapper)
910 {
911 idToIndex_.clear();
912 indexToID_.clear();
913
915 }
916
917 template<int codim>
919 {
920 if constexpr (codim <= GridView::dimension)
921 {
922 for (const auto& element : elements(gridView_))
923 {
924 for (int i = 0; i < element.subEntities(codim); ++i)
925 {
926 const auto entity = element.template subEntity<codim>(i);
927 if (entity.partitionType() == Dune::BorderEntity)
928 {
929 const auto localIdx = mapper_.index(entity);
930 const DofIdKey key{codim, gridView_.grid().globalIdSet().id(entity)};
931 idToIndex_.emplace(key, localIdx);
932 indexToID_.emplace(localIdx, key);
933 }
934 }
935 }
936
938 }
939 }
940
946 template<class IsGhostFunc>
947 void extendMatrix(Matrix& A, const IsGhostFunc& isGhost)
948 {
949 if constexpr (Dune::Capabilities::canCommunicate<typename GridView::Grid, rowDofCodim>::v)
950 {
951 if (gridView_.comm().size() <= 1)
952 return;
953
954 // make a copy of the matrix as originally assembled
955 Matrix matrixAsAssembled(A);
956
957 // first get entries in sparsity pattern from other processes
958 std::size_t numNonZeroEntries = 0;
959 std::vector<std::set<int>> sparsityPattern; // column indices for every row
960 MatrixPatternExchange<IsGhostFunc> dataHandle(mapper_, idToIndex_, indexToID_, A, sparsityPattern, isGhost);
961 gridView_.communicate(
962 dataHandle, Dune::InteriorBorder_InteriorBorder_Interface, Dune::ForwardCommunication
963 );
964
965 // add own entries to the sparsity pattern and count number of non-zeros
966 for (auto rowIt = A.begin(); rowIt != A.end(); ++rowIt)
967 {
968 const auto colEndIt = A[rowIt.index()].end();
969 for (auto colIt = rowIt->begin(); colIt != colEndIt; ++colIt)
970 if (!sparsityPattern[rowIt.index()].count(colIt.index()))
971 sparsityPattern[rowIt.index()].insert(colIt.index());
972
973 numNonZeroEntries += sparsityPattern[rowIt.index()].size();
974 }
975
976 // insert any additional entries into the matrix
977 A.setSize(matrixAsAssembled.N(), matrixAsAssembled.M(), numNonZeroEntries);
978 A.setBuildMode(Matrix::row_wise);
979 auto citer = A.createbegin();
980 for (auto i = sparsityPattern.begin(), end = sparsityPattern.end(); i!=end; ++i, ++citer)
981 for (auto si = i->begin(), send = i->end(); si!=send; ++si)
982 citer.insert(*si);
983
984 // reset matrix to contain original values
985 A = 0;
986 const auto rowEndIt = matrixAsAssembled.end();
987 for (auto rowIt = matrixAsAssembled.begin(); rowIt != rowEndIt; ++rowIt)
988 for (auto colIt = matrixAsAssembled[rowIt.index()].begin(); colIt != matrixAsAssembled[rowIt.index()].end(); ++colIt)
989 A[rowIt.index()][colIt.index()] = *colIt;
990
991 // The matrix A has now a possible extended sparsity pattern but any additional entries are
992 // initialized to zero and have to be filled by the sumEntries function
993 }
994 else
995 DUNE_THROW(Dune::InvalidStateException, "Cannot call extendMatrix for a grid that cannot communicate codim-" << rowDofCodim << "-entities.");
996 }
997
1011 void sumEntries(Matrix& A)
1012 {
1013 if constexpr (Dune::Capabilities::canCommunicate<typename GridView::Grid, rowDofCodim>::v)
1014 {
1015 if (gridView_.comm().size() <= 1)
1016 return;
1017
1018 MatrixEntryExchange dataHandle(mapper_, idToIndex_, indexToID_, A);
1019 gridView_.communicate(
1020 dataHandle, Dune::InteriorBorder_InteriorBorder_Interface, Dune::ForwardCommunication
1021 );
1022 }
1023 else
1024 DUNE_THROW(Dune::InvalidStateException,
1025 "Cannot call sumEntries for a grid that cannot communicate codim-" << rowDofCodim << "-entities."
1026 );
1027 }
1028
1029private:
1030 const GridView gridView_;
1031 const RowDofMapper& mapper_;
1032 std::map<DofIdKey, int> idToIndex_;
1033 std::map<int, DofIdKey> indexToID_;
1034
1035};
1036
1041template<class Matrix, class GridView, class DofMapper, std::size_t numCodims>
1043{
1044 using IdType = typename GridView::Traits::Grid::Traits::GlobalIdSet::IdType;
1045
1046 struct DofIdKey
1047 {
1048 int codim;
1049 IdType id;
1050 int subIndex;
1051
1052 friend bool operator<(const DofIdKey& a, const DofIdKey& b)
1053 {
1054 if (a.codim != b.codim)
1055 return a.codim < b.codim;
1056 if (a.id != b.id)
1057 return a.id < b.id;
1058 return a.subIndex < b.subIndex;
1059 }
1060 };
1061
1062 // Each sent item encodes which row DOF within the entity and which column DOF.
1063 // This handles multi-DOF entities (e.g. PQ3 edges with 2 DOFs) correctly.
1064 struct PatternEntry
1065 {
1066 int rowSubIndex;
1067 DofIdKey colKey;
1068 };
1069
1070 template<class ColIsGhostFunc>
1071 struct MatrixPatternExchange
1072 : public Dune::CommDataHandleIF<MatrixPatternExchange<ColIsGhostFunc>, PatternEntry>
1073 {
1074 using DataType = PatternEntry;
1075
1076 MatrixPatternExchange(const DofMapper& dofMapper,
1077 const std::map<DofIdKey, int>& globalToLocal,
1078 const std::map<int, DofIdKey>& localToGlobal,
1079 Matrix& A,
1080 std::vector<std::set<int>>& sparsityPattern,
1081 const ColIsGhostFunc& isGhostColumDof,
1082 const std::bitset<numCodims>& activeCodims)
1083 : dofMapper_(dofMapper)
1084 , idToIndex_(globalToLocal)
1085 , indexToID_(localToGlobal)
1086 , sparsityPattern_(sparsityPattern)
1087 , A_(A)
1088 , isGhostColumDof_(isGhostColumDof)
1089 , activeCodims_(activeCodims)
1090 {
1091 sparsityPattern_.resize(A.N());
1092 }
1093
1094 bool contains(int dim, int codim) const
1095 {
1096 if (codim < 0 || codim >= static_cast<int>(numCodims))
1097 return false;
1098 return activeCodims_.test(codim);
1099 }
1100
1101 bool fixedSize(int dim, int codim) const
1102 { return false; }
1103
1105 template<class EntityType>
1106 std::size_t size(EntityType& e) const
1107 {
1108 std::size_t n = 0;
1109 for (auto rowIdx : asMultiMapper(dofMapper_).indices(e))
1110 for (auto colIt = A_[rowIdx].begin(); colIt != A_[rowIdx].end(); ++colIt)
1111 if (indexToID_.count(colIt.index()))
1112 ++n;
1113 return n;
1114 }
1115
1117 template<class MessageBuffer, class EntityType>
1118 void gather(MessageBuffer& buff, const EntityType& e) const
1119 {
1120 int rowSub = 0;
1121 for (auto rowIdx : asMultiMapper(dofMapper_).indices(e))
1122 {
1123 for (auto colIt = A_[rowIdx].begin(); colIt != A_[rowIdx].end(); ++colIt)
1124 if (auto it = indexToID_.find(colIt.index()); it != indexToID_.end())
1125 buff.write(PatternEntry{rowSub, it->second});
1126 ++rowSub;
1127 }
1128 }
1129
1131 template<class MessageBuffer, class EntityType>
1132 void scatter(MessageBuffer& buff, const EntityType& e, std::size_t n)
1133 {
1134 const auto baseIdx = dofMapper_.index(e);
1135 for (std::size_t k = 0; k < n; k++)
1136 {
1137 PatternEntry entry;
1138 buff.read(entry);
1139 const auto rowIdx = baseIdx + entry.rowSubIndex;
1140 if (const auto it = idToIndex_.find(entry.colKey); it != idToIndex_.end())
1141 {
1142 const auto colIdx = it->second;
1143 if (!isGhostColumDof_(colIdx))
1144 sparsityPattern_[rowIdx].insert(colIdx);
1145 }
1146 }
1147 }
1148
1149 private:
1150 const DofMapper& dofMapper_;
1151 const std::map<DofIdKey, int>& idToIndex_;
1152 const std::map<int, DofIdKey>& indexToID_;
1153 std::vector<std::set<int>>& sparsityPattern_;
1154 Matrix& A_;
1155 const ColIsGhostFunc& isGhostColumDof_;
1156 const std::bitset<numCodims>& activeCodims_;
1157 };
1158
1159 struct MatrixEntry
1160 {
1161 int rowSubIndex;
1162 DofIdKey colKey;
1163 typename Matrix::block_type block;
1164 };
1165
1166 struct MatrixEntryExchange
1167 : public Dune::CommDataHandleIF<MatrixEntryExchange, MatrixEntry>
1168 {
1169 using DataType = MatrixEntry;
1170
1171 MatrixEntryExchange(const DofMapper& dofMapper,
1172 const std::map<DofIdKey, int>& globalToLocal,
1173 const std::map<int, DofIdKey>& localToGlobal,
1174 Matrix& A,
1175 const std::bitset<numCodims>& activeCodims)
1176 : dofMapper_(dofMapper)
1177 , idToIndex_(globalToLocal)
1178 , indexToID_(localToGlobal)
1179 , A_(A)
1180 , activeCodims_(activeCodims)
1181 {}
1182
1183 bool contains(int dim, int codim) const
1184 {
1185 if (codim < 0 || codim >= static_cast<int>(numCodims))
1186 return false;
1187 return activeCodims_.test(codim);
1188 }
1189
1190 bool fixedSize(int dim, int codim) const
1191 { return false; }
1192
1194 template<class EntityType>
1195 std::size_t size(EntityType& e) const
1196 {
1197 std::size_t n = 0;
1198 for (auto rowIdx : asMultiMapper(dofMapper_).indices(e))
1199 for (auto colIt = A_[rowIdx].begin(); colIt != A_[rowIdx].end(); ++colIt)
1200 if (indexToID_.count(colIt.index()))
1201 ++n;
1202 return n;
1203 }
1204
1206 template<class MessageBuffer, class EntityType>
1207 void gather(MessageBuffer& buff, const EntityType& e) const
1208 {
1209 int rowSub = 0;
1210 for (auto rowIdx : asMultiMapper(dofMapper_).indices(e))
1211 {
1212 for (auto colIt = A_[rowIdx].begin(); colIt != A_[rowIdx].end(); ++colIt)
1213 if (auto it = indexToID_.find(colIt.index()); it != indexToID_.end())
1214 buff.write(MatrixEntry{rowSub, it->second, *colIt});
1215 ++rowSub;
1216 }
1217 }
1218
1220 template<class MessageBuffer, class EntityType>
1221 void scatter(MessageBuffer& buff, const EntityType& e, std::size_t n)
1222 {
1223 const auto baseIdx = dofMapper_.index(e);
1224 for (std::size_t k = 0; k < n; k++)
1225 {
1226 MatrixEntry m;
1227 buff.read(m);
1228 const auto rowIdx = baseIdx + m.rowSubIndex;
1229 if (auto it = idToIndex_.find(m.colKey); it != idToIndex_.end())
1230 {
1231 if (A_[rowIdx].find(it->second) != A_[rowIdx].end())
1232 A_[rowIdx][it->second] += m.block;
1233 }
1234 }
1235 }
1236
1237 private:
1238 const DofMapper& dofMapper_;
1239 const std::map<DofIdKey, int>& idToIndex_;
1240 const std::map<int, DofIdKey>& indexToID_;
1241 Matrix& A_;
1242 const std::bitset<numCodims>& activeCodims_;
1243 };
1244
1245public:
1246 // Standard constructor: only communicates border-entity DOF couplings.
1247 // Suitable for iterative solvers (ISTL AMG) where ghost DOF equations
1248 // are handled separately by the parallel operator mechanism.
1249 MultiCodimParallelMatrixHelper(const GridView& gridView,
1250 const DofMapper& mapper,
1251 const std::bitset<numCodims>& activeCodims)
1252 : gridView_(gridView)
1253 , mapper_(mapper)
1254 , activeCodims_(activeCodims)
1255 , includeGhostAndAdjacent_(false)
1256 , commInterface_(Dune::InteriorBorder_InteriorBorder_Interface)
1257 {
1258 initMapsForCodim_<0>();
1259 }
1260
1261 // Extended constructor for direct solvers (e.g. Trilinos/MUMPS):
1262 // also includes ghost entities and boundary-adjacent interior entities
1263 // so that element-interior DOF column entries in border DOF matrix rows
1264 // are correctly communicated across ranks.
1265 // Uses All_All_Interface so that ghost entities (GhostEntity on non-owner
1266 // rank) can participate in the scatter step and receive data from owners.
1267 MultiCodimParallelMatrixHelper(const GridView& gridView,
1268 const DofMapper& mapper,
1269 const std::bitset<numCodims>& activeCodims,
1270 bool includeGhostAndAdjacent)
1271 : gridView_(gridView)
1272 , mapper_(mapper)
1273 , activeCodims_(activeCodims)
1274 , includeGhostAndAdjacent_(includeGhostAndAdjacent)
1275 // InteriorBorder_All_Interface: owner (Interior/Border) gathers and
1276 // sends; ghost (All) scatters and receives. This propagates the
1277 // owner's physical stiffness to ghost DOFs without sending ghost-
1278 // treatment (identity) entries from the non-owner back to the owner.
1279 , commInterface_(includeGhostAndAdjacent ? Dune::InteriorBorder_All_Interface
1280 : Dune::InteriorBorder_InteriorBorder_Interface)
1281 {
1282 initMapsForCodim_<0>();
1283 }
1284
1285 template<class IsGhostFunc>
1286 void extendMatrix(Matrix& A, const IsGhostFunc& isGhost)
1287 {
1288 if (gridView_.comm().size() <= 1)
1289 return;
1290
1291 Matrix matrixAsAssembled(A);
1292 std::size_t numNonZeroEntries = 0;
1293 std::vector<std::set<int>> sparsityPattern;
1294 MatrixPatternExchange<IsGhostFunc> dataHandle(mapper_, idToIndex_, indexToID_, A, sparsityPattern, isGhost, activeCodims_);
1295 gridView_.communicate(dataHandle, commInterface_, Dune::ForwardCommunication);
1296
1297 for (auto rowIt = A.begin(); rowIt != A.end(); ++rowIt)
1298 {
1299 const auto colEndIt = A[rowIt.index()].end();
1300 for (auto colIt = rowIt->begin(); colIt != colEndIt; ++colIt)
1301 {
1302 if (!sparsityPattern[rowIt.index()].count(colIt.index()))
1303 sparsityPattern[rowIt.index()].insert(colIt.index());
1304 }
1305
1306 numNonZeroEntries += sparsityPattern[rowIt.index()].size();
1307 }
1308
1309 A.setSize(matrixAsAssembled.N(), matrixAsAssembled.M(), numNonZeroEntries);
1310 A.setBuildMode(Matrix::row_wise);
1311 auto citer = A.createbegin();
1312 for (auto i = sparsityPattern.begin(), end = sparsityPattern.end(); i != end; ++i, ++citer)
1313 for (auto si = i->begin(), send = i->end(); si != send; ++si)
1314 citer.insert(*si);
1315
1316 A = 0;
1317 const auto rowEndIt = matrixAsAssembled.end();
1318 for (auto rowIt = matrixAsAssembled.begin(); rowIt != rowEndIt; ++rowIt)
1319 for (auto colIt = matrixAsAssembled[rowIt.index()].begin(); colIt != matrixAsAssembled[rowIt.index()].end(); ++colIt)
1320 A[rowIt.index()][colIt.index()] = *colIt;
1321 }
1322
1323 void sumEntries(Matrix& A)
1324 {
1325 if (gridView_.comm().size() <= 1)
1326 return;
1327
1328 MatrixEntryExchange dataHandle(mapper_, idToIndex_, indexToID_, A, activeCodims_);
1329 gridView_.communicate(dataHandle, commInterface_, Dune::ForwardCommunication);
1330 }
1331
1332private:
1333 template<std::size_t codim>
1334 void initMapsForCodim_()
1335 {
1336 if constexpr (codim <= GridView::dimension)
1337 {
1338 if (activeCodims_.test(codim))
1339 {
1340 for (const auto& element : elements(gridView_))
1341 {
1342 for (int i = 0; i < element.subEntities(codim); ++i)
1343 {
1344 const auto entity = element.template subEntity<codim>(i);
1345 const auto pt = entity.partitionType();
1346
1347 // Standard: always include BorderEntity DOFs.
1348 bool shouldProcess = (pt == Dune::BorderEntity);
1349
1350 // Extended mode (for direct solvers): also include ghost entities
1351 // and interior entities adjacent to the partition boundary so that
1352 // element-interior DOF column entries in border rows are communicated.
1353 if (!shouldProcess && includeGhostAndAdjacent_) {
1354 if constexpr (codim == 0) {
1355 // Elements are never BorderEntity; include GhostEntity elements
1356 // AND InteriorEntity elements adjacent to a border sub-entity.
1357 shouldProcess = (pt == Dune::GhostEntity)
1358 || (pt == Dune::InteriorEntity
1359 && elementHasBorderSubEntity_(entity));
1360 } else {
1361 shouldProcess = (pt == Dune::GhostEntity);
1362 }
1363 }
1364
1365 if (shouldProcess)
1366 {
1367 // Use asMultiMapper(mapper_).indices(entity) to handle multi-DOF entities
1368 // (e.g. PQ3 with 2 DOFs per edge). Each DOF gets a unique
1369 // DofIdKey that encodes both the entity ID and the DOF sub-index.
1370 const auto globalId = gridView_.grid().globalIdSet().id(entity);
1371 int subIdx = 0;
1372 for (auto dofIdx : asMultiMapper(mapper_).indices(entity))
1373 {
1374 const DofIdKey key{static_cast<int>(codim), globalId, subIdx};
1375 idToIndex_.emplace(key, dofIdx);
1376 indexToID_.emplace(dofIdx, key);
1377 ++subIdx;
1378 }
1379 }
1380 }
1381 }
1382 }
1383
1384 initMapsForCodim_<codim + 1>();
1385 }
1386 }
1387
1388 // Helper: returns true if element has at least one BorderEntity sub-entity
1389 template<class Element>
1390 static bool elementHasBorderSubEntity_(const Element& element)
1391 {
1392 for (int i = 0; i < element.subEntities(1); ++i)
1393 if (element.template subEntity<1>(i).partitionType() == Dune::BorderEntity)
1394 return true;
1395 if constexpr (Element::mydimension >= 2)
1396 for (int i = 0; i < element.subEntities(2); ++i)
1397 if (element.template subEntity<2>(i).partitionType() == Dune::BorderEntity)
1398 return true;
1399 if constexpr (Element::mydimension >= 3)
1400 for (int i = 0; i < element.subEntities(3); ++i)
1401 if (element.template subEntity<3>(i).partitionType() == Dune::BorderEntity)
1402 return true;
1403 return false;
1404 }
1405
1406 const GridView gridView_;
1407 const DofMapper& mapper_;
1408 std::bitset<numCodims> activeCodims_;
1409 bool includeGhostAndAdjacent_;
1410 Dune::InterfaceType commInterface_;
1411 std::map<DofIdKey, int> idToIndex_;
1412 std::map<int, DofIdKey> indexToID_;
1413};
1414
1418template<class LinearSolverTraits, class ParallelTraits,
1419 class Matrix, class ParallelHelper>
1420void prepareMatrixParallel(Matrix& A, ParallelHelper& pHelper)
1421{
1422 if constexpr (ParallelTraits::isNonOverlapping)
1423 {
1424 using GridView = typename LinearSolverTraits::GridView;
1425 using DofMapper = typename LinearSolverTraits::DofMapper;
1426
1427 if constexpr (requires { LinearSolverTraits::dofCodims; })
1428 {
1429 constexpr auto dofCodims = LinearSolverTraits::dofCodims;
1431 pHelper.gridView(), pHelper.dofMapper(), dofCodims
1432 );
1433 matrixHelper.extendMatrix(A, [&pHelper](auto idx){ return pHelper.isGhost(idx); });
1434 matrixHelper.sumEntries(A);
1435 }
1436 else
1437 {
1438 static constexpr int dofCodim = LinearSolverTraits::dofCodim;
1439 ParallelMatrixHelper<Matrix, GridView, DofMapper, dofCodim> matrixHelper(pHelper.gridView(), pHelper.dofMapper());
1440 matrixHelper.extendMatrix(A, [&pHelper](auto idx){ return pHelper.isGhost(idx); });
1441 matrixHelper.sumEntries(A);
1442 }
1443 }
1444}
1445
1449template<class LinearSolverTraits, class ParallelTraits,
1450 class Vector, class ParallelHelper>
1451void prepareVectorParallel(Vector& b, ParallelHelper& pHelper)
1452{
1453 if constexpr (ParallelTraits::isNonOverlapping)
1454 {
1455 using GridView = typename LinearSolverTraits::GridView;
1456 using DofMapper = typename LinearSolverTraits::DofMapper;
1457
1458 if constexpr (requires { LinearSolverTraits::dofCodims; })
1459 {
1460 MultiCodimParallelVectorHelper<GridView, DofMapper> vectorHelper(pHelper.gridView(), pHelper.dofMapper());
1461 vectorHelper.makeNonOverlappingConsistent(b, LinearSolverTraits::dofCodims);
1462 }
1463 else
1464 {
1465 static constexpr int dofCodim = LinearSolverTraits::dofCodim;
1466 ParallelVectorHelper<GridView, DofMapper, dofCodim> vectorHelper(pHelper.gridView(), pHelper.dofMapper());
1467 vectorHelper.makeNonOverlappingConsistent(b);
1468 }
1469 }
1470}
1471
1475template<class LinearSolverTraits, class ParallelTraits,
1476 class Matrix, class Vector, class ParallelHelper>
1477void prepareLinearAlgebraParallel(Matrix& A, Vector& b, ParallelHelper& pHelper)
1478{
1481}
1482
1486template<class LinearSolverTraits, class ParallelTraits,
1487 class Matrix, class Vector, class ParallelHelper>
1488void prepareLinearAlgebraParallel(Matrix& A, Vector& b,
1489 std::shared_ptr<typename ParallelTraits::Comm>& comm,
1490 std::shared_ptr<typename ParallelTraits::LinearOperator>& fop,
1491 std::shared_ptr<typename ParallelTraits::ScalarProduct>& sp,
1492 ParallelHelper& pHelper)
1493{
1495 const auto category = ParallelTraits::isNonOverlapping ?
1496 Dune::SolverCategory::nonoverlapping : Dune::SolverCategory::overlapping;
1497
1498 comm = std::make_shared<typename ParallelTraits::Comm>(pHelper.gridView().comm(), category);
1499 pHelper.createParallelIndexSet(*comm);
1500 fop = std::make_shared<typename ParallelTraits::LinearOperator>(A, *comm);
1501 sp = std::make_shared<typename ParallelTraits::ScalarProduct>(*comm);
1502}
1503
1504} // end namespace Dumux
1505
1506#endif
const DofMapper & dofMapper() const
Return the dofMapper.
Definition parallelhelpers.hh:418
ParallelISTLHelperImpl(const GridView &gridView, const DofMapper &mapper)
Definition parallelhelpers.hh:327
bool isOwned(std::size_t i) const
Definition parallelhelpers.hh:341
const GridView & gridView() const
Return the gridView.
Definition parallelhelpers.hh:422
bool isGhost(std::size_t i) const
Definition parallelhelpers.hh:338
void createParallelIndexSet(Comm &comm) const
Creates a parallel index set for a communicator.
Definition parallelhelpers.hh:356
Definition parallelhelpers.hh:35
Dedicated helper that exchanges matrix pattern and entries for multiple active codims in one pass.
Definition parallelhelpers.hh:1043
void sumEntries(Matrix &A)
Definition parallelhelpers.hh:1323
MultiCodimParallelMatrixHelper(const GridView &gridView, const DofMapper &mapper, const std::bitset< numCodims > &activeCodims)
Definition parallelhelpers.hh:1249
MultiCodimParallelMatrixHelper(const GridView &gridView, const DofMapper &mapper, const std::bitset< numCodims > &activeCodims, bool includeGhostAndAdjacent)
Definition parallelhelpers.hh:1267
void extendMatrix(Matrix &A, const IsGhostFunc &isGhost)
Definition parallelhelpers.hh:1286
Definition parallelhelpers.hh:580
MultiCodimParallelVectorHelper(const GridView &gridView, const DofMapper &mapper)
Definition parallelhelpers.hh:582
void makeOverlappingConsistent(Dune::BlockVector< Block, Alloc > &v, const std::bitset< numCodims > &activeCodims) const
Make a vector consistent for overlapping domain decomposition methods.
Definition parallelhelpers.hh:607
void makeNonOverlappingConsistent(Dune::BlockVector< Block, Alloc > &v, const std::bitset< numCodims > &activeCodims) const
Make a vector consistent for non-overlapping domain decomposition methods.
Definition parallelhelpers.hh:588
Helper class for adding up matrix entries for border entities.
Definition parallelhelpers.hh:666
ParallelMatrixHelper(const GridView &gridView, const RowDofMapper &mapper)
Definition parallelhelpers.hh:908
void sumEntries(Matrix &A)
Sums up the entries corresponding to border entities (usually vertices or faces).
Definition parallelhelpers.hh:1011
void initMapsForCodim_()
Definition parallelhelpers.hh:918
void extendMatrix(Matrix &A, const IsGhostFunc &isGhost)
communicates values for the sparsity pattern of the new matrix.
Definition parallelhelpers.hh:947
Definition parallelhelpers.hh:522
void makeNonOverlappingConsistent(Dune::MultiTypeBlockVector< Blocks... > &v) const
Make a vector consistent for non-overlapping domain decomposition methods.
Definition parallelhelpers.hh:568
ParallelVectorHelper(const GridView &gridView, const DofMapper &mapper)
Definition parallelhelpers.hh:524
void makeOverlappingConsistent(Dune::BlockVector< Block, Alloc > &v) const
Make a vector consistent for overlapping domain decomposition methods.
Definition parallelhelpers.hh:545
void makeNonOverlappingConsistent(Dune::BlockVector< Block, Alloc > &v) const
Make a vector consistent for non-overlapping domain decomposition methods.
Definition parallelhelpers.hh:530
Definition variablesbackend.hh:34
dune-grid capabilities compatibility layer
@ element
Definition fieldtype.hh:23
Detail::ParallelISTLHelperImpl< LinearSolverTraits, LinearSolverTraits::canCommunicate > ParallelISTLHelper
A parallel helper class providing a parallel decomposition of all degrees of freedom.
Definition parallelhelpers.hh:514
Adapter to expose a multi-DOF mapper interface for single- and multi-DOF mappers.
Definition cvfelocalresidual.hh:25
Definition adapt.hh:17
VectorCommDataHandle< Mapper, Vector, codim, Detail::Sum, DataType > VectorCommDataHandleSum
Definition vectorcommdatahandle.hh:214
void prepareMatrixParallel(Matrix &A, ParallelHelper &pHelper)
Prepare a matrix for parallel solvers.
Definition parallelhelpers.hh:1420
LinearSolverTraitsImpl< GridGeometry, typename GridGeometry::DiscretizationMethod > LinearSolverTraits
The type traits required for using the IstlFactoryBackend.
Definition linearsolvertraits.hh:39
void prepareLinearAlgebraParallel(Matrix &A, Vector &b, ParallelHelper &pHelper)
Prepare linear algebra variables for parallel solvers.
Definition parallelhelpers.hh:1477
constexpr auto asMultiMapper(const Mapper &mapper)
Definition multimapperview.hh:48
void prepareVectorParallel(Vector &b, ParallelHelper &pHelper)
Prepare a vector for parallel solvers.
Definition parallelhelpers.hh:1451
MultiCodimVectorCommDataHandle< Mapper, Vector, dim, Detail::Sum, DataType > MultiCodimVectorCommDataHandleSum
Definition vectorcommdatahandle.hh:226
Definition common/pdesolver.hh:24
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Contains a class to exchange entries of a vector.