3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_LINEAR_PARALLELHELPERS_HH
25#define DUMUX_LINEAR_PARALLELHELPERS_HH
26
27#include <dune/common/exceptions.hh>
28#include <dune/geometry/dimension.hh>
29#include <dune/grid/common/datahandleif.hh>
30#include <dune/grid/common/partitionset.hh>
31#include <dune/istl/owneroverlapcopy.hh>
32#include <dune/istl/paamg/pinfo.hh>
35
36namespace Dumux {
37
43// operator that resets result to zero at constrained DOFS
44template<class LinearSolverTraits>
46{
47 using GridView = typename LinearSolverTraits::GridView;
48 using DofMapper = typename LinearSolverTraits::DofMapper;
49 static constexpr int dofCodim = LinearSolverTraits::dofCodim;
50
51 // TODO: this is some large number (replace by limits?)
52 static constexpr std::size_t ghostMarker_ = 1<<24;
53
54 class BaseGatherScatter
55 {
56 public:
57 BaseGatherScatter(const DofMapper& mapper)
58 : mapper_(mapper) {}
59
60 template<class EntityType>
61 int index(const EntityType& e) const
62 { return mapper_.index(e); }
63
64 bool contains(int dim, int codim) const
65 { return dofCodim == codim; }
66
68 bool fixedSize(int dim, int codim) const
69 { return true; }
70
71 template<class EntityType>
72 std::size_t size(EntityType& e) const
73 { return 1; }
74
75 template<class EntityType>
76 bool isNeitherInteriorNorBorderEntity(EntityType& e) const
77 { return e.partitionType() != Dune::InteriorEntity && e.partitionType() != Dune::BorderEntity; }
78
79 private:
80 const DofMapper& mapper_;
81 };
82
89 class GhostGatherScatter
90 : public BaseGatherScatter
91 , public Dune::CommDataHandleIF<GhostGatherScatter, std::size_t>
92 {
93 public:
94 using DataType = std::size_t;
95 using BaseGatherScatter::contains;
96 using BaseGatherScatter::fixedSize;
97 using BaseGatherScatter::size;
98
99 GhostGatherScatter(std::vector<std::size_t>& ranks, const DofMapper& mapper)
100 : BaseGatherScatter(mapper)
101 , ranks_(ranks)
102 {}
103
104 template<class MessageBuffer, class EntityType>
105 void gather(MessageBuffer& buff, const EntityType& e) const
106 {
107 auto& data = ranks_[this->index(e)];
108 if (this->isNeitherInteriorNorBorderEntity(e))
109 data = ghostMarker_;
110 buff.write(data);
111 }
112
113 template<class MessageBuffer, class EntityType>
114 void scatter(MessageBuffer& buff, const EntityType& e, std::size_t n)
115 {
116 std::size_t x;
117 buff.read(x);
118 auto& data = ranks_[this->index(e)];
119 if (this->isNeitherInteriorNorBorderEntity(e))
120 data = ghostMarker_;
121 }
122 private:
123 std::vector<std::size_t>& ranks_;
124 };
125
132 class InteriorBorderGatherScatter
133 : public BaseGatherScatter
134 , public Dune::CommDataHandleIF<InteriorBorderGatherScatter, std::size_t>
135 {
136 public:
137 using DataType = std::size_t;
138 using BaseGatherScatter::contains;
139 using BaseGatherScatter::fixedSize;
140 using BaseGatherScatter::size;
141
142 InteriorBorderGatherScatter(std::vector<std::size_t>& ranks, const DofMapper& mapper)
143 : BaseGatherScatter(mapper)
144 , ranks_(ranks)
145 {}
146
150 template<class MessageBuffer, class EntityType>
151 void gather(MessageBuffer& buff, const EntityType& e) const
152 {
153 auto& data = ranks_[this->index(e)];
154 if (this->isNeitherInteriorNorBorderEntity(e))
155 data = ghostMarker_;
156 buff.write(data);
157 }
158
162 template<class MessageBuffer, class EntityType>
163 void scatter(MessageBuffer& buff, const EntityType& e, std::size_t n)
164 {
165 std::size_t x;
166 buff.read(x);
167 auto& data = ranks_[this->index(e)];
168
169 // we leave ghost unchanged
170 // for other dofs, the process with the lowest rank
171 // is assigned to be the (unique) owner
172 using std::min;
173 data = this->isNeitherInteriorNorBorderEntity(e) ? x : min(data, x);
174 }
175 private:
176 std::vector<std::size_t>& ranks_;
177 };
178
183 struct NeighbourGatherScatter
184 : public BaseGatherScatter
185 , public Dune::CommDataHandleIF<NeighbourGatherScatter, int>
186 {
187 using DataType = int;
188 using BaseGatherScatter::contains;
189 using BaseGatherScatter::fixedSize;
190 using BaseGatherScatter::size;
191
192 NeighbourGatherScatter(const DofMapper& mapper, int rank, std::set<int>& neighbours)
193 : BaseGatherScatter(mapper)
194 , rank_(rank)
195 , neighbours_(neighbours)
196 {}
197
198 template<class MessageBuffer, class EntityType>
199 void gather(MessageBuffer& buff, const EntityType& e) const
200 {
201 buff.write(rank_);
202 }
203
204 template<class MessageBuffer, class EntityType>
205 void scatter(MessageBuffer& buff, const EntityType& e, std::size_t n)
206 {
207 int x;
208 buff.read(x);
209 neighbours_.insert(x);
210 }
211 private:
212 int rank_;
213 std::set<int>& neighbours_;
214 };
215
216
221 struct SharedGatherScatter
222 : public BaseGatherScatter
223 , public Dune::CommDataHandleIF<SharedGatherScatter, int>
224 {
225 using DataType = int;
226 using BaseGatherScatter::contains;
227 using BaseGatherScatter::fixedSize;
228 using BaseGatherScatter::size;
229
230 SharedGatherScatter(std::vector<int>& shared, const DofMapper& mapper)
231 : BaseGatherScatter(mapper)
232 , shared_(shared)
233 {}
234
235 template<class MessageBuffer, class EntityType>
236 void gather(MessageBuffer& buff, EntityType& e) const
237 {
238 int data = true;
239 buff.write(data);
240 }
241
242 template<class MessageBuffer, class EntityType>
243 void scatter(MessageBuffer& buff, const EntityType &e, std::size_t n)
244 {
245 int x;
246 buff.read(x);
247 auto& data = shared_[this->index(e)];
248 data = data || x;
249 }
250 private:
251 std::vector<int>& shared_;
252 };
253
258 template<typename GlobalIndex>
259 struct GlobalIndexGatherScatter
260 : public BaseGatherScatter
261 , public Dune::CommDataHandleIF<GlobalIndexGatherScatter<GlobalIndex>, GlobalIndex>
262 {
263 using DataType = GlobalIndex;
264 using BaseGatherScatter::contains;
265 using BaseGatherScatter::fixedSize;
266 using BaseGatherScatter::size;
267
268 GlobalIndexGatherScatter(std::vector<GlobalIndex>& globalIndices, const DofMapper& mapper)
269 : BaseGatherScatter(mapper)
270 , globalIndices_(globalIndices)
271 {}
272
273 template<class MessageBuffer, class EntityType>
274 void gather(MessageBuffer& buff, const EntityType& e) const
275 {
276 buff.write(globalIndices_[this->index(e)]);
277 }
278
279 template<class MessageBuffer, class EntityType>
280 void scatter(MessageBuffer& buff, const EntityType& e, std::size_t n)
281 {
282 DataType x;
283 buff.read(x);
284 using std::min;
285 auto& data = globalIndices_[this->index(e)];
286 data = min(data, x);
287 }
288 private:
289 std::vector<GlobalIndex>& globalIndices_;
290 };
291
292public:
293
294 ParallelISTLHelper(const GridView& gridView, const DofMapper& mapper)
295 : gridView_(gridView), mapper_(mapper)
296 {
297 if constexpr (Detail::canCommunicate<typename GridView::Traits::Grid, dofCodim>)
298 initGhostsAndOwners_();
299 else
300 DUNE_THROW(Dune::InvalidStateException,
301 "Cannot initialize parallel helper for a grid that cannot communicate codim-" << dofCodim << "-entities."
302 );
303 }
304
305 bool isGhost(std::size_t i) const
306 { return isGhost_[i] == ghostMarker_; }
307
308 bool isOwned(std::size_t i) const
309 { return isOwned_[i] == 1; }
310
322 template<class Comm>
323 void createParallelIndexSet(Comm& comm) const
324 {
325 if constexpr (Detail::canCommunicate<typename GridView::Traits::Grid, dofCodim>)
326 {
327 if (gridView_.comm().size() <= 1)
328 {
329 comm.remoteIndices().template rebuild<false>();
330 return;
331 }
332
333 // First find out which dofs we share with other processors
334 std::vector<int> isShared(mapper_.size(), false);
335 SharedGatherScatter sgs(isShared, mapper_);
336 gridView_.communicate(sgs, Dune::All_All_Interface, Dune::ForwardCommunication);
337
338 // Count shared dofs that we own
339 using GlobalIndex = typename Comm::ParallelIndexSet::GlobalIndex;
340 GlobalIndex count = 0;
341 for (std::size_t i = 0; i < isShared.size(); ++i)
342 if (isShared[i] && isOwned_[i] == 1)
343 ++count;
344
345 std::vector<GlobalIndex> counts(gridView_.comm().size());
346 gridView_.comm().allgather(&count, 1, counts.data());
347
348 // Compute start index start_p = \sum_{i=0}^{i<p} counts_i
349 const int rank = gridView_.comm().rank();
350 auto start = std::accumulate(counts.begin(), counts.begin() + rank, GlobalIndex(0));
351
352 // starting from start (offset), number global indices owned by this process consecutively
353 std::vector<GlobalIndex> globalIndices(mapper_.size(), std::numeric_limits<GlobalIndex>::max());
354 for (std::size_t i = 0; i < globalIndices.size(); ++i)
355 {
356 if (isOwned_[i] == 1 && isShared[i])
357 {
358 globalIndices[i] = start;
359 ++start;
360 }
361 }
362
363 // publish global indices for the shared DOFS to other processors.
364 GlobalIndexGatherScatter<GlobalIndex> gigs(globalIndices, mapper_);
365 gridView_.communicate(gigs, Dune::All_All_Interface, Dune::ForwardCommunication);
366
367 // put the information into the parallel index set
368 resizeAndFillIndexSet_(comm, globalIndices);
369
370 // Compute neighbours using communication
371 std::set<int> neighbours;
372 NeighbourGatherScatter ngs(mapper_, gridView_.comm().rank(), neighbours);
373 gridView_.communicate(ngs, Dune::All_All_Interface, Dune::ForwardCommunication);
374 comm.remoteIndices().setNeighbours(neighbours);
375
376 comm.remoteIndices().template rebuild<false>();
377 }
378 else
379 DUNE_THROW(Dune::InvalidStateException,
380 "Cannot build parallel index set for a grid that cannot communicate codim-" << dofCodim << "-entities."
381 );
382 }
383
385 const DofMapper& dofMapper() const
386 { return mapper_; }
387
389 const GridView& gridView() const
390 { return gridView_; }
391
392private:
401 void initGhostsAndOwners_()
402 {
403 const auto rank = gridView_.comm().rank();
404 isOwned_.assign(mapper_.size(), rank);
405
406 // find out about ghosts
407 GhostGatherScatter ggs(isOwned_, mapper_);
408
409 if (gridView_.comm().size() > 1)
410 gridView_.communicate(ggs, Dune::InteriorBorder_All_Interface, Dune::ForwardCommunication);
411
412 // isGhost_ contains rank when the dof is not on a ghost entity (interior or border)
413 // and ghostMarker_ when the dof is on a ghost
414 isGhost_ = isOwned_;
415
416 // partition interior/border uniquely
417 // after the communication each dof is assigned to one unique owner process rank
418 InteriorBorderGatherScatter dh(isOwned_, mapper_);
419
420 if (gridView_.comm().size() > 1)
421 gridView_.communicate(dh, Dune::InteriorBorder_InteriorBorder_Interface, Dune::ForwardCommunication);
422
423 // convert vector into mask vector
424 for (auto& v : isOwned_)
425 v = (v == rank) ? 1 : 0;
426 }
427
428 template<class Comm, class GlobalIndices>
429 void resizeAndFillIndexSet_(Comm& comm, const GlobalIndices& globalIndices) const
430 {
431 comm.indexSet().beginResize();
432
433 // add triplets characterizing each dof in the parallel index set
434 // (globalIndex, (localIndex, attribute)) where attribute is owner, overlap or copy
435 for (std::size_t localIndex = 0; localIndex < globalIndices.size(); ++localIndex)
436 {
437 const auto globalIndex = globalIndices[localIndex];
438 if (globalIndex != std::numeric_limits<typename GlobalIndices::value_type>::max())
439 {
440 const bool isOwned = isOwned_[localIndex] > 0;
441 const auto attr = getAttribute_(comm, isOwned, isGhost(localIndex));
442 using LocalIndex = typename Comm::ParallelIndexSet::LocalIndex;
443 comm.indexSet().add(globalIndex, LocalIndex{localIndex, attr});
444 }
445 }
446
447 comm.indexSet().endResize();
448 }
449
450 template<class Comm>
451 Dune::OwnerOverlapCopyAttributeSet::AttributeSet
452 getAttribute_(const Comm& comm, const bool isOwned, const bool isGhost) const
453 {
454 if (isOwned) // this process owns the dof (uniquely)
455 return Dune::OwnerOverlapCopyAttributeSet::owner;
456 else if (isGhost && (comm.category() == Dune::SolverCategory::nonoverlapping) )
457 return Dune::OwnerOverlapCopyAttributeSet::overlap;
458 else
459 return Dune::OwnerOverlapCopyAttributeSet::copy;
460 }
461
462 const GridView gridView_;
463 const DofMapper& mapper_;
465 std::vector<std::size_t> isOwned_;
467 std::vector<std::size_t> isGhost_;
468
469}; // class ParallelISTLHelper
470
471template<class GridView, class DofMapper, int dofCodim>
473{
474public:
475 ParallelVectorHelper(const GridView& gridView, const DofMapper& mapper)
476 : gridView_(gridView), mapper_(mapper)
477 {}
478
479 // \brief Make a vector of the box model consistent.
480 template<class Block, class Alloc>
481 void makeNonOverlappingConsistent(Dune::BlockVector<Block, Alloc>& v) const
482 {
483 if constexpr (Detail::canCommunicate<typename GridView::Traits::Grid, dofCodim>)
484 {
486 if (gridView_.comm().size() > 1)
487 gridView_.communicate(gs, Dune::InteriorBorder_InteriorBorder_Interface,
488 Dune::ForwardCommunication);
489 }
490 else
491 DUNE_THROW(Dune::InvalidStateException, "Cannot call makeNonOverlappingConsistent for a grid that cannot communicate codim-" << dofCodim << "-entities.");
492 }
493
494private:
495 const GridView gridView_;
496 const DofMapper& mapper_;
497};
498
507template<class Matrix, class GridView,
508 class RowDofMapper, int rowDofCodim>
510{
511 using IdType = typename GridView::Traits::Grid::Traits::GlobalIdSet::IdType;
512
534 template<class ColIsGhostFunc>
535 struct MatrixPatternExchange
536 : public Dune::CommDataHandleIF<MatrixPatternExchange<ColIsGhostFunc>, IdType>
537 {
539 using DataType = IdType;
540
551 MatrixPatternExchange(const RowDofMapper& rowEntityMapper,
552 const std::map<IdType,int>& globalToLocal,
553 const std::map<int,IdType>& localToGlobal,
554 Matrix& A,
555 std::vector<std::set<int>>& sparsityPattern,
556 const ColIsGhostFunc& isGhostColumDof)
557 : rowEntityMapper_(rowEntityMapper), idToIndex_(globalToLocal), indexToID_(localToGlobal)
558 , sparsityPattern_(sparsityPattern), A_(A), isGhostColumDof_(isGhostColumDof)
559 {
560 sparsityPattern_.resize(A.N());
561 }
562
566 bool contains (int dim, int codim) const
567 { return (codim == rowDofCodim); }
568
570 bool fixedSize(int dim, int codim) const
571 { return false; }
572
576 template<class EntityType>
577 std::size_t size(EntityType& e) const
578 {
579 const auto rowIdx = rowEntityMapper_.index(e);
580
581 // all column entity indices of this row that are in the index set
582 std::size_t n = 0;
583 for (auto colIt = A_[rowIdx].begin(); colIt != A_[rowIdx].end(); ++colIt)
584 if (indexToID_.count(colIt.index()))
585 n++;
586
587 return n;
588 }
589
593 template<class MessageBuffer, class EntityType>
594 void gather(MessageBuffer& buff, const EntityType& e) const
595 {
596 const auto rowIdx = rowEntityMapper_.index(e);
597
598 // send all column entity ids of this row that are in the index set
599 for (auto colIt = A_[rowIdx].begin(); colIt != A_[rowIdx].end(); ++colIt)
600 if (auto it = indexToID_.find(colIt.index()); it != indexToID_.end())
601 buff.write(it->second);
602 }
603
607 template<class MessageBuffer, class EntityType>
608 void scatter(MessageBuffer& buff, const EntityType& e, std::size_t n)
609 {
610 const auto rowIdx = rowEntityMapper_.index(e);
611 for (std::size_t k = 0; k < n; k++)
612 {
613 // receive n column entity IDs
614 IdType id;
615 buff.read(id);
616
617 // only add entries that are contained in the index set
618 if (const auto it = idToIndex_.find(id); it != idToIndex_.end())
619 {
620 const auto colIdx = it->second; // get local column index
621 // add this entry (if it doesn't exist yet)
622 // (and only if the column dof is not associated with a ghost entity)
623 if (!isGhostColumDof_(colIdx))
624 // std::set takes care that each index only is inserted once
625 sparsityPattern_[rowIdx].insert(colIdx);
626 }
627 }
628 }
629
630 private:
631 const RowDofMapper& rowEntityMapper_;
632 const std::map<IdType, int>& idToIndex_;
633 const std::map<int, IdType>& indexToID_;
634 std::vector<std::set<int>>& sparsityPattern_;
635 Matrix& A_;
636 const ColIsGhostFunc& isGhostColumDof_;
637
638 }; // class MatrixPatternExchange
639
641 struct MatrixEntry
642 {
643 IdType first;
644 typename Matrix::block_type second;
645 };
646
648 struct MatrixEntryExchange
649 : public Dune::CommDataHandleIF<MatrixEntryExchange, MatrixEntry>
650 {
652 using DataType = MatrixEntry;
653
654 MatrixEntryExchange(const RowDofMapper& rowEntityMapper,
655 const std::map<IdType, int>& globalToLocal,
656 const std::map<int, IdType>& localToGlobal,
657 Matrix& A)
658 : rowEntityMapper_(rowEntityMapper), idToIndex_(globalToLocal), indexToID_(localToGlobal), A_(A)
659 {}
660
664 bool contains(int dim, int codim) const
665 { return (codim == rowDofCodim); }
666
668 bool fixedSize(int dim, int codim) const
669 { return false; }
670
674 template<class EntityType>
675 std::size_t size(EntityType& e) const
676 {
677 const auto rowIdx = rowEntityMapper_.index(e);
678 std::size_t n = 0;
679 for (auto colIt = A_[rowIdx].begin(); colIt != A_[rowIdx].end(); ++colIt)
680 if (indexToID_.count(colIt.index()))
681 n++;
682
683 return n;
684 }
685
689 template<class MessageBuffer, class EntityType>
690 void gather(MessageBuffer& buff, const EntityType& e) const
691 {
692 const auto rowIdx = rowEntityMapper_.index(e);
693 // send all matrix entries for which the column index is in the index set
694 for (auto colIt = A_[rowIdx].begin(); colIt != A_[rowIdx].end(); ++colIt)
695 if (auto it = indexToID_.find(colIt.index()); it != indexToID_.end())
696 buff.write(MatrixEntry{ it->second, *colIt });
697 }
698
702 template<class MessageBuffer, class EntityType>
703 void scatter(MessageBuffer& buff, const EntityType& e, std::size_t n)
704 {
705 const auto rowIdx = rowEntityMapper_.index(e);
706 for (std::size_t k = 0; k < n; k++)
707 {
708 MatrixEntry m;
709 buff.read(m);
710 const auto& [colDofID, matrixBlock] = m; // unpack
711
712 // only add entries in the index set and for which A has allocated memory
713 if (auto it = idToIndex_.find(colDofID); it != idToIndex_.end())
714 if (A_[rowIdx].find(it->second) != A_[rowIdx].end())
715 A_[rowIdx][it->second] += matrixBlock;
716 }
717 }
718
719 private:
720 const RowDofMapper& rowEntityMapper_;
721 const std::map<IdType, int>& idToIndex_;
722 const std::map<int, IdType>& indexToID_;
723 Matrix& A_;
724
725 }; // class MatrixEntryExchange
726
727public:
728
729 ParallelMatrixHelper(const GridView& gridView, const RowDofMapper& mapper)
730 : gridView_(gridView), mapper_(mapper)
731 {
732 idToIndex_.clear();
733 indexToID_.clear();
734
735 for (const auto& entity : entities(gridView_, Dune::Codim<rowDofCodim>()))
736 {
737 if (entity.partitionType() == Dune::BorderEntity)
738 {
739 const int localIdx = mapper_.index(entity);
740 IdType dofIdxGlobal = gridView_.grid().globalIdSet().id(entity);
741
742 idToIndex_.emplace(dofIdxGlobal, localIdx);
743 indexToID_.emplace(localIdx, dofIdxGlobal);
744 }
745 }
746 }
747
753 template<class IsGhostFunc>
754 void extendMatrix(Matrix& A, const IsGhostFunc& isGhost)
755 {
756 if constexpr (Detail::canCommunicate<typename GridView::Grid, rowDofCodim>)
757 {
758 if (gridView_.comm().size() <= 1)
759 return;
760
761 // make a copy of the matrix as originally assembled
762 Matrix matrixAsAssembled(A);
763
764 // first get entries in sparsity pattern from other processes
765 std::size_t numNonZeroEntries = 0;
766 std::vector<std::set<int>> sparsityPattern; // column indices for every row
767 MatrixPatternExchange<IsGhostFunc> dataHandle(mapper_, idToIndex_, indexToID_, A, sparsityPattern, isGhost);
768 gridView_.communicate(
769 dataHandle, Dune::InteriorBorder_InteriorBorder_Interface, Dune::ForwardCommunication
770 );
771
772 // add own entries to the sparsity pattern and count number of non-zeros
773 for (auto rowIt = A.begin(); rowIt != A.end(); ++rowIt)
774 {
775 const auto colEndIt = A[rowIt.index()].end();
776 for (auto colIt = rowIt->begin(); colIt != colEndIt; ++colIt)
777 if (!sparsityPattern[rowIt.index()].count(colIt.index()))
778 sparsityPattern[rowIt.index()].insert(colIt.index());
779
780 numNonZeroEntries += sparsityPattern[rowIt.index()].size();
781 }
782
783 // insert any additional entries into the matrix
784 A.setSize(matrixAsAssembled.N(), matrixAsAssembled.M(), numNonZeroEntries);
785 A.setBuildMode(Matrix::row_wise);
786 auto citer = A.createbegin();
787 for (auto i = sparsityPattern.begin(), end = sparsityPattern.end(); i!=end; ++i, ++citer)
788 for (auto si = i->begin(), send = i->end(); si!=send; ++si)
789 citer.insert(*si);
790
791 // reset matrix to contain original values
792 A = 0;
793 const auto rowEndIt = matrixAsAssembled.end();
794 for (auto rowIt = matrixAsAssembled.begin(); rowIt != rowEndIt; ++rowIt)
795 for (auto colIt = matrixAsAssembled[rowIt.index()].begin(); colIt != matrixAsAssembled[rowIt.index()].end(); ++colIt)
796 A[rowIt.index()][colIt.index()] = *colIt;
797
798 // The matrix A has now a possible extended sparsity pattern but any additional entries are
799 // initialized to zero and have to be filled by the sumEntries function
800 }
801 else
802 DUNE_THROW(Dune::InvalidStateException, "Cannot call extendMatrix for a grid that cannot communicate codim-" << rowDofCodim << "-entities.");
803 }
804
818 void sumEntries(Matrix& A)
819 {
820 if constexpr (Detail::canCommunicate<typename GridView::Grid, rowDofCodim>)
821 {
822 if (gridView_.comm().size() <= 1)
823 return;
824
825 MatrixEntryExchange dataHandle(mapper_, idToIndex_, indexToID_, A);
826 gridView_.communicate(
827 dataHandle, Dune::InteriorBorder_InteriorBorder_Interface, Dune::ForwardCommunication
828 );
829 }
830 else
831 DUNE_THROW(Dune::InvalidStateException,
832 "Cannot call sumEntries for a grid that cannot communicate codim-" << rowDofCodim << "-entities."
833 );
834 }
835
836private:
837 const GridView gridView_;
838 const RowDofMapper& mapper_;
839 std::map<IdType, int> idToIndex_;
840 std::map<int, IdType> indexToID_;
841
842};
843
847template<class LinearSolverTraits, class ParallelTraits,
848 class Matrix, class Vector, class ParallelHelper>
849void prepareLinearAlgebraParallel(Matrix& A, Vector& b, ParallelHelper& pHelper)
850{
851 if constexpr (ParallelTraits::isNonOverlapping)
852 {
853 // extend the matrix pattern such that it is usable for a parallel solver
854 // and make right-hand side consistent
855 using GridView = typename LinearSolverTraits::GridView;
856 using DofMapper = typename LinearSolverTraits::DofMapper;
857 static constexpr int dofCodim = LinearSolverTraits::dofCodim;
858 ParallelMatrixHelper<Matrix, GridView, DofMapper, dofCodim> matrixHelper(pHelper.gridView(), pHelper.dofMapper());
859 matrixHelper.extendMatrix(A, [&pHelper](auto idx){ return pHelper.isGhost(idx); });
860 matrixHelper.sumEntries(A);
861
862 ParallelVectorHelper<GridView, DofMapper, dofCodim> vectorHelper(pHelper.gridView(), pHelper.dofMapper());
863 vectorHelper.makeNonOverlappingConsistent(b);
864 }
865}
866
870template<class LinearSolverTraits, class ParallelTraits,
871 class Matrix, class Vector, class ParallelHelper>
872void prepareLinearAlgebraParallel(Matrix& A, Vector& b,
873 std::shared_ptr<typename ParallelTraits::Comm>& comm,
874 std::shared_ptr<typename ParallelTraits::LinearOperator>& fop,
875 std::shared_ptr<typename ParallelTraits::ScalarProduct>& sp,
876 ParallelHelper& pHelper)
877{
878 prepareLinearAlgebraParallel<LinearSolverTraits, ParallelTraits>(A, b, pHelper);
879 const auto category = ParallelTraits::isNonOverlapping ?
880 Dune::SolverCategory::nonoverlapping : Dune::SolverCategory::overlapping;
881
882 comm = std::make_shared<typename ParallelTraits::Comm>(pHelper.gridView().comm(), category);
883 pHelper.createParallelIndexSet(*comm);
884 fop = std::make_shared<typename ParallelTraits::LinearOperator>(A, *comm);
885 sp = std::make_shared<typename ParallelTraits::ScalarProduct>(*comm);
886}
887
888} // end namespace Dumux
889
890#endif
dune-grid capabilities compatibility layer
Contains a class to exchange entries of a vector.
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
void prepareLinearAlgebraParallel(Matrix &A, Vector &b, ParallelHelper &pHelper)
Prepare linear algebra variables for parallel solvers.
Definition: parallelhelpers.hh:849
LinearSolverTraitsImpl< GridGeometry, typename GridGeometry::DiscretizationMethod > LinearSolverTraits
The type traits required for using the IstlFactoryBackend.
Definition: linearsolvertraits.hh:49
A parallel helper class providing a nonoverlapping decomposition of all degrees of freedom.
Definition: parallelhelpers.hh:46
bool isGhost(std::size_t i) const
Definition: parallelhelpers.hh:305
void createParallelIndexSet(Comm &comm) const
Creates a parallel index set for a communicator.
Definition: parallelhelpers.hh:323
const DofMapper & dofMapper() const
Return the dofMapper.
Definition: parallelhelpers.hh:385
bool isOwned(std::size_t i) const
Definition: parallelhelpers.hh:308
const GridView & gridView() const
Return the gridView.
Definition: parallelhelpers.hh:389
ParallelISTLHelper(const GridView &gridView, const DofMapper &mapper)
Definition: parallelhelpers.hh:294
Definition: parallelhelpers.hh:473
ParallelVectorHelper(const GridView &gridView, const DofMapper &mapper)
Definition: parallelhelpers.hh:475
void makeNonOverlappingConsistent(Dune::BlockVector< Block, Alloc > &v) const
Definition: parallelhelpers.hh:481
Helper class for adding up matrix entries for border entities.
Definition: parallelhelpers.hh:510
ParallelMatrixHelper(const GridView &gridView, const RowDofMapper &mapper)
Definition: parallelhelpers.hh:729
void sumEntries(Matrix &A)
Sums up the entries corresponding to border entities (usually vertices or faces)
Definition: parallelhelpers.hh:818
void extendMatrix(Matrix &A, const IsGhostFunc &isGhost)
communicates values for the sparsity pattern of the new matrix.
Definition: parallelhelpers.hh:754
A data handle class to exchange entries of a vector.
Definition: vectorcommdatahandle.hh:78