version 3.11-dev
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 <dune/common/exceptions.hh>
16#include <dune/geometry/dimension.hh>
17#include <dune/grid/common/datahandleif.hh>
18#include <dune/grid/common/partitionset.hh>
19#include <dune/istl/owneroverlapcopy.hh>
20#include <dune/istl/paamg/pinfo.hh>
21#include <dune/istl/bvector.hh>
22#include <dune/istl/multitypeblockvector.hh>
25
26namespace Dumux::Detail {
27
28template<class LinearSolverTraits, bool canCommunicate = false>
30
31template<class LinearSolverTraits>
33{
34 using GridView = typename LinearSolverTraits::GridView;
35 using DofMapper = typename LinearSolverTraits::DofMapper;
36 static constexpr int dofCodim = LinearSolverTraits::dofCodim;
37
38 // TODO: this is some large number (replace by limits?)
39 static constexpr std::size_t ghostMarker_ = 1<<24;
40
41 class BaseGatherScatter
42 {
43 public:
44 BaseGatherScatter(const DofMapper& mapper)
45 : mapper_(mapper) {}
46
47 template<class EntityType>
48 int index(const EntityType& e) const
49 { return mapper_.index(e); }
50
51 bool contains(int dim, int codim) const
52 { return dofCodim == codim; }
53
55 bool fixedSize(int dim, int codim) const
56 { return true; }
57
58 template<class EntityType>
59 std::size_t size(EntityType& e) const
60 { return 1; }
61
62 template<class EntityType>
63 bool isNeitherInteriorNorBorderEntity(EntityType& e) const
64 { return e.partitionType() != Dune::InteriorEntity && e.partitionType() != Dune::BorderEntity; }
65
66 private:
67 const DofMapper& mapper_;
68 };
69
76 class GhostGatherScatter
77 : public BaseGatherScatter
78 , public Dune::CommDataHandleIF<GhostGatherScatter, std::size_t>
79 {
80 public:
81 using DataType = std::size_t;
82 using BaseGatherScatter::contains;
83 using BaseGatherScatter::fixedSize;
84 using BaseGatherScatter::size;
85
86 GhostGatherScatter(std::vector<std::size_t>& ranks, const DofMapper& mapper)
87 : BaseGatherScatter(mapper)
88 , ranks_(ranks)
89 {}
90
91 template<class MessageBuffer, class EntityType>
92 void gather(MessageBuffer& buff, const EntityType& e) const
93 {
94 auto& data = ranks_[this->index(e)];
95 if (this->isNeitherInteriorNorBorderEntity(e))
96 data = ghostMarker_;
97 buff.write(data);
98 }
99
100 template<class MessageBuffer, class EntityType>
101 void scatter(MessageBuffer& buff, const EntityType& e, std::size_t n)
102 {
103 std::size_t x;
104 buff.read(x);
105 auto& data = ranks_[this->index(e)];
106 if (this->isNeitherInteriorNorBorderEntity(e))
107 data = ghostMarker_;
108 }
109 private:
110 std::vector<std::size_t>& ranks_;
111 };
112
119 class InteriorBorderGatherScatter
120 : public BaseGatherScatter
121 , public Dune::CommDataHandleIF<InteriorBorderGatherScatter, std::size_t>
122 {
123 public:
124 using DataType = std::size_t;
125 using BaseGatherScatter::contains;
126 using BaseGatherScatter::fixedSize;
127 using BaseGatherScatter::size;
128
129 InteriorBorderGatherScatter(std::vector<std::size_t>& ranks, const DofMapper& mapper)
130 : BaseGatherScatter(mapper)
131 , ranks_(ranks)
132 {}
133
137 template<class MessageBuffer, class EntityType>
138 void gather(MessageBuffer& buff, const EntityType& e) const
139 {
140 auto& data = ranks_[this->index(e)];
141 if (this->isNeitherInteriorNorBorderEntity(e))
142 data = ghostMarker_;
143 buff.write(data);
144 }
145
149 template<class MessageBuffer, class EntityType>
150 void scatter(MessageBuffer& buff, const EntityType& e, std::size_t n)
151 {
152 std::size_t x;
153 buff.read(x);
154 auto& data = ranks_[this->index(e)];
155
156 // we leave ghost unchanged
157 // for other dofs, the process with the lowest rank
158 // is assigned to be the (unique) owner
159 using std::min;
160 data = this->isNeitherInteriorNorBorderEntity(e) ? x : min(data, x);
161 }
162 private:
163 std::vector<std::size_t>& ranks_;
164 };
165
170 struct NeighbourGatherScatter
171 : public BaseGatherScatter
172 , public Dune::CommDataHandleIF<NeighbourGatherScatter, int>
173 {
174 using DataType = int;
175 using BaseGatherScatter::contains;
176 using BaseGatherScatter::fixedSize;
177 using BaseGatherScatter::size;
178
179 NeighbourGatherScatter(const DofMapper& mapper, int rank, std::set<int>& neighbours)
180 : BaseGatherScatter(mapper)
181 , rank_(rank)
182 , neighbours_(neighbours)
183 {}
184
185 template<class MessageBuffer, class EntityType>
186 void gather(MessageBuffer& buff, const EntityType& e) const
187 {
188 buff.write(rank_);
189 }
190
191 template<class MessageBuffer, class EntityType>
192 void scatter(MessageBuffer& buff, const EntityType& e, std::size_t n)
193 {
194 int x;
195 buff.read(x);
196 neighbours_.insert(x);
197 }
198 private:
199 int rank_;
200 std::set<int>& neighbours_;
201 };
202
203
208 struct SharedGatherScatter
209 : public BaseGatherScatter
210 , public Dune::CommDataHandleIF<SharedGatherScatter, int>
211 {
212 using DataType = int;
213 using BaseGatherScatter::contains;
214 using BaseGatherScatter::fixedSize;
215 using BaseGatherScatter::size;
216
217 SharedGatherScatter(std::vector<int>& shared, const DofMapper& mapper)
218 : BaseGatherScatter(mapper)
219 , shared_(shared)
220 {}
221
222 template<class MessageBuffer, class EntityType>
223 void gather(MessageBuffer& buff, EntityType& e) const
224 {
225 int data = true;
226 buff.write(data);
227 }
228
229 template<class MessageBuffer, class EntityType>
230 void scatter(MessageBuffer& buff, const EntityType &e, std::size_t n)
231 {
232 int x;
233 buff.read(x);
234 auto& data = shared_[this->index(e)];
235 data = data || x;
236 }
237 private:
238 std::vector<int>& shared_;
239 };
240
245 template<typename GlobalIndex>
246 struct GlobalIndexGatherScatter
247 : public BaseGatherScatter
248 , public Dune::CommDataHandleIF<GlobalIndexGatherScatter<GlobalIndex>, GlobalIndex>
249 {
250 using DataType = GlobalIndex;
251 using BaseGatherScatter::contains;
252 using BaseGatherScatter::fixedSize;
253 using BaseGatherScatter::size;
254
255 GlobalIndexGatherScatter(std::vector<GlobalIndex>& globalIndices, const DofMapper& mapper)
256 : BaseGatherScatter(mapper)
257 , globalIndices_(globalIndices)
258 {}
259
260 template<class MessageBuffer, class EntityType>
261 void gather(MessageBuffer& buff, const EntityType& e) const
262 {
263 buff.write(globalIndices_[this->index(e)]);
264 }
265
266 template<class MessageBuffer, class EntityType>
267 void scatter(MessageBuffer& buff, const EntityType& e, std::size_t n)
268 {
269 DataType x;
270 buff.read(x);
271 using std::min;
272 DataType& data = globalIndices_[this->index(e)];
273 data = min(data, x);
274 }
275 private:
276 std::vector<GlobalIndex>& globalIndices_;
277 };
278
279public:
280
281 ParallelISTLHelperImpl(const GridView& gridView, const DofMapper& mapper)
282 : gridView_(gridView), mapper_(mapper)
283 {
284 if constexpr (Dune::Capabilities::canCommunicate<typename GridView::Traits::Grid, dofCodim>::v)
285 initGhostsAndOwners_();
286 else
287 DUNE_THROW(Dune::InvalidStateException,
288 "Cannot initialize parallel helper for a grid that cannot communicate codim-" << dofCodim << "-entities."
289 );
290 }
291
292 bool isGhost(std::size_t i) const
293 { return isGhost_[i] == ghostMarker_; }
294
295 bool isOwned(std::size_t i) const
296 { return isOwned_[i] == 1; }
297
309 template<class Comm>
310 void createParallelIndexSet(Comm& comm) const
311 {
312 if constexpr (Dune::Capabilities::canCommunicate<typename GridView::Traits::Grid, dofCodim>::v)
313 {
314 if (gridView_.comm().size() <= 1)
315 {
316 comm.remoteIndices().template rebuild<false>();
317 return;
318 }
319
320 // First find out which dofs we share with other processors
321 std::vector<int> isShared(mapper_.size(), false);
322 SharedGatherScatter sgs(isShared, mapper_);
323 gridView_.communicate(sgs, Dune::All_All_Interface, Dune::ForwardCommunication);
324
325 // Count shared dofs that we own
326 using GlobalIndex = typename Comm::ParallelIndexSet::GlobalIndex;
327 GlobalIndex count = 0;
328 for (std::size_t i = 0; i < isShared.size(); ++i)
329 if (isShared[i] && isOwned_[i] == 1)
330 ++count;
331
332 std::vector<GlobalIndex> counts(gridView_.comm().size());
333 gridView_.comm().allgather(&count, 1, counts.data());
334
335 // Compute start index start_p = \sum_{i=0}^{i<p} counts_i
336 const int rank = gridView_.comm().rank();
337 auto start = std::accumulate(counts.begin(), counts.begin() + rank, GlobalIndex(0));
338
339 // starting from start (offset), number global indices owned by this process consecutively
340 std::vector<GlobalIndex> globalIndices(mapper_.size(), std::numeric_limits<GlobalIndex>::max());
341 for (std::size_t i = 0; i < globalIndices.size(); ++i)
342 {
343 if (isOwned_[i] == 1 && isShared[i])
344 {
345 globalIndices[i] = start;
346 ++start;
347 }
348 }
349
350 // publish global indices for the shared DOFS to other processors.
351 GlobalIndexGatherScatter<GlobalIndex> gigs(globalIndices, mapper_);
352 gridView_.communicate(gigs, Dune::All_All_Interface, Dune::ForwardCommunication);
353
354 // put the information into the parallel index set
355 resizeAndFillIndexSet_(comm, globalIndices);
356
357 // Compute neighbours using communication
358 std::set<int> neighbours;
359 NeighbourGatherScatter ngs(mapper_, gridView_.comm().rank(), neighbours);
360 gridView_.communicate(ngs, Dune::All_All_Interface, Dune::ForwardCommunication);
361 comm.remoteIndices().setNeighbours(neighbours);
362
363 comm.remoteIndices().template rebuild<false>();
364 }
365 else
366 DUNE_THROW(Dune::InvalidStateException,
367 "Cannot build parallel index set for a grid that cannot communicate codim-" << dofCodim << "-entities."
368 );
369 }
370
372 const DofMapper& dofMapper() const
373 { return mapper_; }
374
376 const GridView& gridView() const
377 { return gridView_; }
378
379private:
388 void initGhostsAndOwners_()
389 {
390 const auto rank = gridView_.comm().rank();
391 isOwned_.assign(mapper_.size(), rank);
392
393 // find out about ghosts
394 GhostGatherScatter ggs(isOwned_, mapper_);
395
396 if (gridView_.comm().size() > 1)
397 gridView_.communicate(ggs, Dune::InteriorBorder_All_Interface, Dune::ForwardCommunication);
398
399 // isGhost_ contains rank when the dof is not on a ghost entity (interior or border)
400 // and ghostMarker_ when the dof is on a ghost
401 isGhost_ = isOwned_;
402
403 // partition interior/border uniquely
404 // after the communication each dof is assigned to one unique owner process rank
405 InteriorBorderGatherScatter dh(isOwned_, mapper_);
406
407 if (gridView_.comm().size() > 1)
408 gridView_.communicate(dh, Dune::InteriorBorder_InteriorBorder_Interface, Dune::ForwardCommunication);
409
410 // convert vector into mask vector
411 for (auto& v : isOwned_)
412 v = (v == rank) ? 1 : 0;
413 }
414
415 template<class Comm, class GlobalIndices>
416 void resizeAndFillIndexSet_(Comm& comm, const GlobalIndices& globalIndices) const
417 {
418 comm.indexSet().beginResize();
419
420 // add triplets characterizing each dof in the parallel index set
421 // (globalIndex, (localIndex, attribute)) where attribute is owner, overlap or copy
422 for (std::size_t localIndex = 0; localIndex < globalIndices.size(); ++localIndex)
423 {
424 const auto globalIndex = globalIndices[localIndex];
425 if (globalIndex != std::numeric_limits<typename GlobalIndices::value_type>::max())
426 {
427 const bool isOwned = isOwned_[localIndex] > 0;
428 const auto attr = getAttribute_(comm, isOwned, isGhost(localIndex));
429 using LocalIndex = typename Comm::ParallelIndexSet::LocalIndex;
430 comm.indexSet().add(globalIndex, LocalIndex{localIndex, attr});
431 }
432 }
433
434 comm.indexSet().endResize();
435 }
436
437 template<class Comm>
438 Dune::OwnerOverlapCopyAttributeSet::AttributeSet
439 getAttribute_(const Comm& comm, const bool isOwned, const bool isGhost) const
440 {
441 if (isOwned) // this process owns the dof (uniquely)
442 return Dune::OwnerOverlapCopyAttributeSet::owner;
443 else if (isGhost && (comm.category() == Dune::SolverCategory::nonoverlapping) )
444 return Dune::OwnerOverlapCopyAttributeSet::overlap;
445 else
446 return Dune::OwnerOverlapCopyAttributeSet::copy;
447 }
448
449 const GridView gridView_;
450 const DofMapper& mapper_;
452 std::vector<std::size_t> isOwned_;
454 std::vector<std::size_t> isGhost_;
455
456};
457
458} // end namespace Dumux::Detail
459
460namespace Dumux {
461
467template<class LinearSolverTraits>
470 LinearSolverTraits, LinearSolverTraits::canCommunicate
471 >;
472
473
474template<class GridView, class DofMapper, int dofCodim>
476{
477public:
478 ParallelVectorHelper(const GridView& gridView, const DofMapper& mapper)
479 : gridView_(gridView), mapper_(mapper)
480 {}
481
483 template<class Block, class Alloc>
484 void makeNonOverlappingConsistent(Dune::BlockVector<Block, Alloc>& v) const
485 {
486 if constexpr (Dune::Capabilities::canCommunicate<typename GridView::Traits::Grid, dofCodim>::v)
487 {
489 if (gridView_.comm().size() > 1)
490 gridView_.communicate(gs, Dune::InteriorBorder_InteriorBorder_Interface,
491 Dune::ForwardCommunication);
492 }
493 else
494 DUNE_THROW(Dune::InvalidStateException, "Cannot call makeNonOverlappingConsistent for a grid that cannot communicate codim-" << dofCodim << "-entities.");
495 }
496
498 template<class Block, class Alloc>
499 void makeOverlappingConsistent(Dune::BlockVector<Block, Alloc>& v) const
500 {
501 if constexpr (Dune::Capabilities::canCommunicate<typename GridView::Traits::Grid, dofCodim>::v)
502 {
503 // Iterate over the gridView_ to make overlap entries 0. This is necessary to make
504 // `VectorCommDataHandleSum` work properly.
505 for (const auto& element : elements(gridView_))
506 if (element.partitionType() == Dune::OverlapEntity)
507 v[gridView_.indexSet().index(element)] = 0;
508
509 if (gridView_.comm().size() > 1)
510 {
512 gridView_.communicate(gs, Dune::Overlap_All_Interface,
513 Dune::ForwardCommunication);
514 }
515 }
516 else
517 DUNE_THROW(Dune::InvalidStateException, "Cannot call makeNonOverlappingConsistent for a grid that cannot communicate codim-" << dofCodim << "-entities.");
518 }
519
521 template<class... Blocks>
523 {
524 DUNE_THROW(Dune::NotImplemented, "makeNonOverlappingConsistent for Dune::MultiTypeBlockVector");
525 }
526
527private:
528 const GridView gridView_;
529 const DofMapper& mapper_;
530};
531
540template<class Matrix, class GridView,
541 class RowDofMapper, int rowDofCodim>
543{
544 using IdType = typename GridView::Traits::Grid::Traits::GlobalIdSet::IdType;
545
567 template<class ColIsGhostFunc>
568 struct MatrixPatternExchange
569 : public Dune::CommDataHandleIF<MatrixPatternExchange<ColIsGhostFunc>, IdType>
570 {
572 using DataType = IdType;
573
584 MatrixPatternExchange(const RowDofMapper& rowEntityMapper,
585 const std::map<IdType,int>& globalToLocal,
586 const std::map<int,IdType>& localToGlobal,
587 Matrix& A,
588 std::vector<std::set<int>>& sparsityPattern,
589 const ColIsGhostFunc& isGhostColumDof)
590 : rowEntityMapper_(rowEntityMapper), idToIndex_(globalToLocal), indexToID_(localToGlobal)
591 , sparsityPattern_(sparsityPattern), A_(A), isGhostColumDof_(isGhostColumDof)
592 {
593 sparsityPattern_.resize(A.N());
594 }
595
599 bool contains (int dim, int codim) const
600 { return (codim == rowDofCodim); }
601
603 bool fixedSize(int dim, int codim) const
604 { return false; }
605
609 template<class EntityType>
610 std::size_t size(EntityType& e) const
611 {
612 const auto rowIdx = rowEntityMapper_.index(e);
613
614 // all column entity indices of this row that are in the index set
615 std::size_t n = 0;
616 for (auto colIt = A_[rowIdx].begin(); colIt != A_[rowIdx].end(); ++colIt)
617 if (indexToID_.count(colIt.index()))
618 n++;
619
620 return n;
621 }
622
626 template<class MessageBuffer, class EntityType>
627 void gather(MessageBuffer& buff, const EntityType& e) const
628 {
629 const auto rowIdx = rowEntityMapper_.index(e);
630
631 // send all column entity ids of this row that are in the index set
632 for (auto colIt = A_[rowIdx].begin(); colIt != A_[rowIdx].end(); ++colIt)
633 if (auto it = indexToID_.find(colIt.index()); it != indexToID_.end())
634 buff.write(it->second);
635 }
636
640 template<class MessageBuffer, class EntityType>
641 void scatter(MessageBuffer& buff, const EntityType& e, std::size_t n)
642 {
643 const auto rowIdx = rowEntityMapper_.index(e);
644 for (std::size_t k = 0; k < n; k++)
645 {
646 // receive n column entity IDs
647 IdType id;
648 buff.read(id);
649
650 // only add entries that are contained in the index set
651 if (const auto it = idToIndex_.find(id); it != idToIndex_.end())
652 {
653 const auto colIdx = it->second; // get local column index
654 // add this entry (if it doesn't exist yet)
655 // (and only if the column dof is not associated with a ghost entity)
656 if (!isGhostColumDof_(colIdx))
657 // std::set takes care that each index only is inserted once
658 sparsityPattern_[rowIdx].insert(colIdx);
659 }
660 }
661 }
662
663 private:
664 const RowDofMapper& rowEntityMapper_;
665 const std::map<IdType, int>& idToIndex_;
666 const std::map<int, IdType>& indexToID_;
667 std::vector<std::set<int>>& sparsityPattern_;
668 Matrix& A_;
669 const ColIsGhostFunc& isGhostColumDof_;
670
671 }; // class MatrixPatternExchange
672
674 struct MatrixEntry
675 {
676 IdType first;
677 typename Matrix::block_type second;
678 };
679
681 struct MatrixEntryExchange
682 : public Dune::CommDataHandleIF<MatrixEntryExchange, MatrixEntry>
683 {
685 using DataType = MatrixEntry;
686
687 MatrixEntryExchange(const RowDofMapper& rowEntityMapper,
688 const std::map<IdType, int>& globalToLocal,
689 const std::map<int, IdType>& localToGlobal,
690 Matrix& A)
691 : rowEntityMapper_(rowEntityMapper), idToIndex_(globalToLocal), indexToID_(localToGlobal), A_(A)
692 {}
693
697 bool contains(int dim, int codim) const
698 { return (codim == rowDofCodim); }
699
701 bool fixedSize(int dim, int codim) const
702 { return false; }
703
707 template<class EntityType>
708 std::size_t size(EntityType& e) const
709 {
710 const auto rowIdx = rowEntityMapper_.index(e);
711 std::size_t n = 0;
712 for (auto colIt = A_[rowIdx].begin(); colIt != A_[rowIdx].end(); ++colIt)
713 if (indexToID_.count(colIt.index()))
714 n++;
715
716 return n;
717 }
718
722 template<class MessageBuffer, class EntityType>
723 void gather(MessageBuffer& buff, const EntityType& e) const
724 {
725 const auto rowIdx = rowEntityMapper_.index(e);
726 // send all matrix entries for which the column index is in the index set
727 for (auto colIt = A_[rowIdx].begin(); colIt != A_[rowIdx].end(); ++colIt)
728 if (auto it = indexToID_.find(colIt.index()); it != indexToID_.end())
729 buff.write(MatrixEntry{ it->second, *colIt });
730 }
731
735 template<class MessageBuffer, class EntityType>
736 void scatter(MessageBuffer& buff, const EntityType& e, std::size_t n)
737 {
738 const auto rowIdx = rowEntityMapper_.index(e);
739 for (std::size_t k = 0; k < n; k++)
740 {
741 MatrixEntry m;
742 buff.read(m);
743 const auto& [colDofID, matrixBlock] = m; // unpack
744
745 // only add entries in the index set and for which A has allocated memory
746 if (auto it = idToIndex_.find(colDofID); it != idToIndex_.end())
747 if (A_[rowIdx].find(it->second) != A_[rowIdx].end())
748 A_[rowIdx][it->second] += matrixBlock;
749 }
750 }
751
752 private:
753 const RowDofMapper& rowEntityMapper_;
754 const std::map<IdType, int>& idToIndex_;
755 const std::map<int, IdType>& indexToID_;
756 Matrix& A_;
757
758 }; // class MatrixEntryExchange
759
760public:
761
762 ParallelMatrixHelper(const GridView& gridView, const RowDofMapper& mapper)
763 : gridView_(gridView), mapper_(mapper)
764 {
765 idToIndex_.clear();
766 indexToID_.clear();
767
768 std::vector<bool> handledDof(gridView_.size(rowDofCodim), false);
769 for (const auto& element : elements(gridView_))
770 {
771 for (int i = 0; i < element.subEntities(rowDofCodim); ++i)
772 {
773 const auto entity = element.template subEntity<rowDofCodim>(i);
774 if (entity.partitionType() == Dune::BorderEntity)
775 {
776 const auto localRowIdx = mapper_.index(entity);
777 if (!handledDof[localRowIdx])
778 {
779 IdType dofIdxGlobal = gridView_.grid().globalIdSet().id(entity);
780 idToIndex_.emplace(dofIdxGlobal, localRowIdx);
781 indexToID_.emplace(localRowIdx, dofIdxGlobal);
782 }
783 }
784 }
785 }
786 }
787
793 template<class IsGhostFunc>
794 void extendMatrix(Matrix& A, const IsGhostFunc& isGhost)
795 {
796 if constexpr (Dune::Capabilities::canCommunicate<typename GridView::Grid, rowDofCodim>::v)
797 {
798 if (gridView_.comm().size() <= 1)
799 return;
800
801 // make a copy of the matrix as originally assembled
802 Matrix matrixAsAssembled(A);
803
804 // first get entries in sparsity pattern from other processes
805 std::size_t numNonZeroEntries = 0;
806 std::vector<std::set<int>> sparsityPattern; // column indices for every row
807 MatrixPatternExchange<IsGhostFunc> dataHandle(mapper_, idToIndex_, indexToID_, A, sparsityPattern, isGhost);
808 gridView_.communicate(
809 dataHandle, Dune::InteriorBorder_InteriorBorder_Interface, Dune::ForwardCommunication
810 );
811
812 // add own entries to the sparsity pattern and count number of non-zeros
813 for (auto rowIt = A.begin(); rowIt != A.end(); ++rowIt)
814 {
815 const auto colEndIt = A[rowIt.index()].end();
816 for (auto colIt = rowIt->begin(); colIt != colEndIt; ++colIt)
817 if (!sparsityPattern[rowIt.index()].count(colIt.index()))
818 sparsityPattern[rowIt.index()].insert(colIt.index());
819
820 numNonZeroEntries += sparsityPattern[rowIt.index()].size();
821 }
822
823 // insert any additional entries into the matrix
824 A.setSize(matrixAsAssembled.N(), matrixAsAssembled.M(), numNonZeroEntries);
825 A.setBuildMode(Matrix::row_wise);
826 auto citer = A.createbegin();
827 for (auto i = sparsityPattern.begin(), end = sparsityPattern.end(); i!=end; ++i, ++citer)
828 for (auto si = i->begin(), send = i->end(); si!=send; ++si)
829 citer.insert(*si);
830
831 // reset matrix to contain original values
832 A = 0;
833 const auto rowEndIt = matrixAsAssembled.end();
834 for (auto rowIt = matrixAsAssembled.begin(); rowIt != rowEndIt; ++rowIt)
835 for (auto colIt = matrixAsAssembled[rowIt.index()].begin(); colIt != matrixAsAssembled[rowIt.index()].end(); ++colIt)
836 A[rowIt.index()][colIt.index()] = *colIt;
837
838 // The matrix A has now a possible extended sparsity pattern but any additional entries are
839 // initialized to zero and have to be filled by the sumEntries function
840 }
841 else
842 DUNE_THROW(Dune::InvalidStateException, "Cannot call extendMatrix for a grid that cannot communicate codim-" << rowDofCodim << "-entities.");
843 }
844
858 void sumEntries(Matrix& A)
859 {
860 if constexpr (Dune::Capabilities::canCommunicate<typename GridView::Grid, rowDofCodim>::v)
861 {
862 if (gridView_.comm().size() <= 1)
863 return;
864
865 MatrixEntryExchange dataHandle(mapper_, idToIndex_, indexToID_, A);
866 gridView_.communicate(
867 dataHandle, Dune::InteriorBorder_InteriorBorder_Interface, Dune::ForwardCommunication
868 );
869 }
870 else
871 DUNE_THROW(Dune::InvalidStateException,
872 "Cannot call sumEntries for a grid that cannot communicate codim-" << rowDofCodim << "-entities."
873 );
874 }
875
876private:
877 const GridView gridView_;
878 const RowDofMapper& mapper_;
879 std::map<IdType, int> idToIndex_;
880 std::map<int, IdType> indexToID_;
881
882};
883
887template<class LinearSolverTraits, class ParallelTraits,
888 class Matrix, class ParallelHelper>
889void prepareMatrixParallel(Matrix& A, ParallelHelper& pHelper)
890{
891 if constexpr (ParallelTraits::isNonOverlapping)
892 {
893 // extend the matrix pattern such that it is usable for a parallel solver
894 // and make right-hand side consistent
895 using GridView = typename LinearSolverTraits::GridView;
896 using DofMapper = typename LinearSolverTraits::DofMapper;
897 static constexpr int dofCodim = LinearSolverTraits::dofCodim;
898 ParallelMatrixHelper<Matrix, GridView, DofMapper, dofCodim> matrixHelper(pHelper.gridView(), pHelper.dofMapper());
899 matrixHelper.extendMatrix(A, [&pHelper](auto idx){ return pHelper.isGhost(idx); });
900 matrixHelper.sumEntries(A);
901 }
902}
903
907template<class LinearSolverTraits, class ParallelTraits,
908 class Vector, class ParallelHelper>
909void prepareVectorParallel(Vector& b, ParallelHelper& pHelper)
910{
911 if constexpr (ParallelTraits::isNonOverlapping)
912 {
913 // extend the matrix pattern such that it is usable for a parallel solver
914 // and make right-hand side consistent
915 using GridView = typename LinearSolverTraits::GridView;
916 using DofMapper = typename LinearSolverTraits::DofMapper;
917 static constexpr int dofCodim = LinearSolverTraits::dofCodim;
918 ParallelVectorHelper<GridView, DofMapper, dofCodim> vectorHelper(pHelper.gridView(), pHelper.dofMapper());
919 vectorHelper.makeNonOverlappingConsistent(b);
920 }
921}
922
926template<class LinearSolverTraits, class ParallelTraits,
927 class Matrix, class Vector, class ParallelHelper>
928void prepareLinearAlgebraParallel(Matrix& A, Vector& b, ParallelHelper& pHelper)
929{
930 prepareMatrixParallel<LinearSolverTraits, ParallelTraits>(A, pHelper);
931 prepareVectorParallel<LinearSolverTraits, ParallelTraits>(b, pHelper);
932}
933
937template<class LinearSolverTraits, class ParallelTraits,
938 class Matrix, class Vector, class ParallelHelper>
939void prepareLinearAlgebraParallel(Matrix& A, Vector& b,
940 std::shared_ptr<typename ParallelTraits::Comm>& comm,
941 std::shared_ptr<typename ParallelTraits::LinearOperator>& fop,
942 std::shared_ptr<typename ParallelTraits::ScalarProduct>& sp,
943 ParallelHelper& pHelper)
944{
945 prepareLinearAlgebraParallel<LinearSolverTraits, ParallelTraits>(A, b, pHelper);
946 const auto category = ParallelTraits::isNonOverlapping ?
947 Dune::SolverCategory::nonoverlapping : Dune::SolverCategory::overlapping;
948
949 comm = std::make_shared<typename ParallelTraits::Comm>(pHelper.gridView().comm(), category);
950 pHelper.createParallelIndexSet(*comm);
951 fop = std::make_shared<typename ParallelTraits::LinearOperator>(A, *comm);
952 sp = std::make_shared<typename ParallelTraits::ScalarProduct>(*comm);
953}
954
955} // end namespace Dumux
956
957#endif
const DofMapper & dofMapper() const
Return the dofMapper.
Definition: parallelhelpers.hh:372
ParallelISTLHelperImpl(const GridView &gridView, const DofMapper &mapper)
Definition: parallelhelpers.hh:281
bool isOwned(std::size_t i) const
Definition: parallelhelpers.hh:295
const GridView & gridView() const
Return the gridView.
Definition: parallelhelpers.hh:376
bool isGhost(std::size_t i) const
Definition: parallelhelpers.hh:292
void createParallelIndexSet(Comm &comm) const
Creates a parallel index set for a communicator.
Definition: parallelhelpers.hh:310
Definition: parallelhelpers.hh:29
Helper class for adding up matrix entries for border entities.
Definition: parallelhelpers.hh:543
ParallelMatrixHelper(const GridView &gridView, const RowDofMapper &mapper)
Definition: parallelhelpers.hh:762
void sumEntries(Matrix &A)
Sums up the entries corresponding to border entities (usually vertices or faces)
Definition: parallelhelpers.hh:858
void extendMatrix(Matrix &A, const IsGhostFunc &isGhost)
communicates values for the sparsity pattern of the new matrix.
Definition: parallelhelpers.hh:794
Definition: parallelhelpers.hh:476
void makeNonOverlappingConsistent(Dune::MultiTypeBlockVector< Blocks... > &v) const
Make a vector consistent for non-overlapping domain decomposition methods.
Definition: parallelhelpers.hh:522
ParallelVectorHelper(const GridView &gridView, const DofMapper &mapper)
Definition: parallelhelpers.hh:478
void makeOverlappingConsistent(Dune::BlockVector< Block, Alloc > &v) const
Make a vector consistent for overlapping domain decomposition methods.
Definition: parallelhelpers.hh:499
void makeNonOverlappingConsistent(Dune::BlockVector< Block, Alloc > &v) const
Make a vector consistent for non-overlapping domain decomposition methods.
Definition: parallelhelpers.hh:484
A data handle class to exchange entries of a vector.
Definition: vectorcommdatahandle.hh:66
Definition: variablesbackend.hh:34
dune-grid capabilities compatibility layer
Definition: cvfelocalresidual.hh:28
Definition: adapt.hh:17
void prepareMatrixParallel(Matrix &A, ParallelHelper &pHelper)
Prepare a matrix for parallel solvers.
Definition: parallelhelpers.hh:889
void prepareLinearAlgebraParallel(Matrix &A, Vector &b, ParallelHelper &pHelper)
Prepare linear algebra variables for parallel solvers.
Definition: parallelhelpers.hh:928
LinearSolverTraitsImpl< GridGeometry, typename GridGeometry::DiscretizationMethod > LinearSolverTraits
The type traits required for using the IstlFactoryBackend.
Definition: linearsolvertraits.hh:37
void prepareVectorParallel(Vector &b, ParallelHelper &pHelper)
Prepare a vector for parallel solvers.
Definition: parallelhelpers.hh:909
The implementation is specialized for the different discretizations.
Definition: linearsolvertraits.hh:31
Contains a class to exchange entries of a vector.