25#ifndef DUMUX_LINEAR_PARALLELHELPERS_HH
26#define DUMUX_LINEAR_PARALLELHELPERS_HH
28#include <dune/common/exceptions.hh>
29#include <dune/geometry/dimension.hh>
30#include <dune/grid/common/datahandleif.hh>
31#include <dune/grid/common/partitionset.hh>
32#include <dune/istl/owneroverlapcopy.hh>
33#include <dune/istl/paamg/pinfo.hh>
45template<
class LinearSolverTraits>
48 using GridView =
typename LinearSolverTraits::GridView;
49 using DofMapper =
typename LinearSolverTraits::DofMapper;
50 enum { dofCodim = LinearSolverTraits::dofCodim };
53 static constexpr std::size_t ghostMarker_ = 1<<24;
55 class BaseGatherScatter
58 BaseGatherScatter(
const DofMapper& mapper)
61 template<
class EntityType>
62 int index(
const EntityType& e)
const
63 {
return mapper_.index(e); }
65 bool contains(
int dim,
int codim)
const
66 {
return dofCodim == codim; }
69 bool fixedSize(
int dim,
int codim)
const
72 template<
class EntityType>
73 std::size_t size(EntityType& e)
const
76 template<
class EntityType>
77 bool isNeitherInteriorNorBorderEntity(EntityType& e)
const
78 {
return e.partitionType() != Dune::InteriorEntity && e.partitionType() != Dune::BorderEntity; }
81 const DofMapper& mapper_;
90 class GhostGatherScatter
91 :
public BaseGatherScatter,
92 public Dune::CommDataHandleIF<GhostGatherScatter,std::size_t>
95 using DataType = std::size_t;
96 using BaseGatherScatter::contains;
97 using BaseGatherScatter::fixedSize;
98 using BaseGatherScatter::size;
100 GhostGatherScatter(std::vector<std::size_t>& ranks,
const DofMapper& mapper)
101 : BaseGatherScatter(mapper), ranks_(ranks)
104 template<
class MessageBuffer,
class EntityType>
105 void gather(MessageBuffer& buff,
const EntityType& e)
const
107 auto& data = ranks_[this->index(e)];
108 if (this->isNeitherInteriorNorBorderEntity(e))
113 template<
class MessageBuffer,
class EntityType>
114 void scatter(MessageBuffer& buff,
const EntityType& e, std::size_t n)
118 auto& data = ranks_[this->index(e)];
119 if (this->isNeitherInteriorNorBorderEntity(e))
123 std::vector<std::size_t>& ranks_;
132 class InteriorBorderGatherScatter
133 :
public BaseGatherScatter,
134 public Dune::CommDataHandleIF<InteriorBorderGatherScatter,std::size_t>
137 using DataType = std::size_t;
138 using BaseGatherScatter::contains;
139 using BaseGatherScatter::fixedSize;
140 using BaseGatherScatter::size;
142 InteriorBorderGatherScatter(std::vector<std::size_t>& ranks,
const DofMapper& mapper)
143 : BaseGatherScatter(mapper), ranks_(ranks)
146 template<
class MessageBuffer,
class EntityType>
147 void gather(MessageBuffer& buff,
const EntityType& e)
const
149 auto& data = ranks_[this->index(e)];
150 if (this->isNeitherInteriorNorBorderEntity(e))
155 template<
class MessageBuffer,
class EntityType>
156 void scatter(MessageBuffer& buff,
const EntityType& e, std::size_t n)
160 auto& data = ranks_[this->index(e)];
162 data = this->isNeitherInteriorNorBorderEntity(e) ? x : min(data,x);
165 std::vector<std::size_t>& ranks_;
172 struct NeighbourGatherScatter
173 :
public BaseGatherScatter,
174 public Dune::CommDataHandleIF<NeighbourGatherScatter,int>
176 using DataType = int;
177 using BaseGatherScatter::contains;
178 using BaseGatherScatter::fixedSize;
179 using BaseGatherScatter::size;
181 NeighbourGatherScatter(
const DofMapper& mapper,
int rank, std::set<int>& neighbours)
182 : BaseGatherScatter(mapper), rank_(rank), neighbours_(neighbours)
185 template<
class MessageBuffer,
class EntityType>
186 void gather(MessageBuffer& buff,
const EntityType& e)
const
191 template<
class MessageBuffer,
class EntityType>
192 void scatter(MessageBuffer& buff,
const EntityType& e, std::size_t n)
196 neighbours_.insert(x);
200 std::set<int>& neighbours_;
208 struct SharedGatherScatter
209 :
public BaseGatherScatter,
210 public Dune::CommDataHandleIF<SharedGatherScatter,int>
212 using DataType = int;
213 using BaseGatherScatter::contains;
214 using BaseGatherScatter::fixedSize;
215 using BaseGatherScatter::size;
217 SharedGatherScatter(std::vector<int>& shared,
const DofMapper& mapper)
218 : BaseGatherScatter(mapper), shared_(shared)
221 template<
class MessageBuffer,
class EntityType>
222 void gather(MessageBuffer& buff, EntityType& e)
const
228 template<
class MessageBuffer,
class EntityType>
229 void scatter(MessageBuffer& buff,
const EntityType &e, std::size_t n)
233 auto& data = shared_[this->index(e)];
237 std::vector<int>& shared_;
244 template<
typename GlobalIndex>
245 struct GlobalIndexGatherScatter
246 :
public BaseGatherScatter,
247 public Dune::CommDataHandleIF<GlobalIndexGatherScatter<GlobalIndex>, GlobalIndex>
249 using DataType = GlobalIndex;
250 using BaseGatherScatter::contains;
251 using BaseGatherScatter::fixedSize;
252 using BaseGatherScatter::size;
254 GlobalIndexGatherScatter(std::vector<GlobalIndex>& globalIndices,
const DofMapper& mapper)
255 : BaseGatherScatter(mapper), globalIndices_(globalIndices)
258 template<
class MessageBuffer,
class EntityType>
259 void gather(MessageBuffer& buff,
const EntityType& e)
const
261 buff.write(globalIndices_[this->index(e)]);
264 template<
class MessageBuffer,
class EntityType>
265 void scatter(MessageBuffer& buff,
const EntityType& e, std::size_t n)
270 auto& data = globalIndices_[this->index(e)];
274 std::vector<GlobalIndex>& globalIndices_;
280 : gridView_(
gridView), mapper_(mapper)
282 if constexpr (Detail::canCommunicate<typename GridView::Traits::Grid, dofCodim>)
283 initGhostsAndOwners_();
285 DUNE_THROW(Dune::InvalidStateException,
"Cannot initialize parallel helper for a grid that cannot communicate codim-" << dofCodim <<
"-entities.");
289 {
return isGhost_[i] == ghostMarker_; }
300 if constexpr (Detail::canCommunicate<typename GridView::Traits::Grid, dofCodim>)
302 if (gridView_.comm().size() <= 1)
304 comm.remoteIndices().template rebuild<false>();
309 std::vector<int> isShared(mapper_.size(),
false);
311 SharedGatherScatter sgs(isShared, mapper_);
312 gridView_.communicate(sgs, Dune::All_All_Interface, Dune::ForwardCommunication);
315 using GlobalIndex =
typename Comm::ParallelIndexSet::GlobalIndex;
316 GlobalIndex count = 0;
317 for (std::size_t i = 0; i < isShared.size(); ++i)
318 if (isShared[i] && isOwned_[i] == 1)
321 std::vector<GlobalIndex> counts(gridView_.comm().size());
322 gridView_.comm().allgather(&count, 1, counts.data());
325 const int rank = gridView_.comm().rank();
326 auto start = std::accumulate(counts.begin(), counts.begin() + rank, GlobalIndex(0));
328 std::vector<GlobalIndex> globalIndices(mapper_.size(), std::numeric_limits<GlobalIndex>::max());
330 for (std::size_t i = 0; i < globalIndices.size(); ++i)
332 if (isOwned_[i] == 1 && isShared[i])
334 globalIndices[i] = start;
340 GlobalIndexGatherScatter<GlobalIndex> gigs(globalIndices, mapper_);
341 gridView_.communicate(gigs, Dune::All_All_Interface, Dune::ForwardCommunication);
343 resizeIndexSet_(comm, globalIndices);
346 std::set<int> neighbours;
347 NeighbourGatherScatter ngs(mapper_, gridView_.comm().rank(), neighbours);
348 gridView_.communicate(ngs, Dune::All_All_Interface, Dune::ForwardCommunication);
349 comm.remoteIndices().setNeighbours(neighbours);
351 comm.remoteIndices().template rebuild<false>();
354 DUNE_THROW(Dune::InvalidStateException,
"Cannot build parallel index set for a grid that cannot communicate codim-" << dofCodim <<
"-entities.");
363 {
return gridView_; }
366 void initGhostsAndOwners_()
368 const auto rank = gridView_.comm().rank();
369 isOwned_.assign(mapper_.size(), rank);
371 GhostGatherScatter ggs(isOwned_, mapper_);
373 if (gridView_.comm().size() > 1)
374 gridView_.communicate(ggs, Dune::InteriorBorder_All_Interface, Dune::ForwardCommunication);
379 InteriorBorderGatherScatter dh(isOwned_, mapper_);
381 if (gridView_.comm().size() > 1)
382 gridView_.communicate(dh, Dune::InteriorBorder_InteriorBorder_Interface, Dune::ForwardCommunication);
385 for (
auto& v : isOwned_)
386 v = (v == rank) ? 1 : 0;
389 template<
class Comm,
class GlobalIndices>
390 void resizeIndexSet_(Comm& comm,
const GlobalIndices& globalIndices)
const
392 comm.indexSet().beginResize();
394 for (std::size_t localIndex = 0; localIndex < globalIndices.size(); ++localIndex)
396 const auto globalIndex = globalIndices[localIndex];
397 if (globalIndex != std::numeric_limits<typename GlobalIndices::value_type>::max())
399 const bool isOwned = isOwned_[localIndex] > 0;
400 const auto attr = getAttribute_(comm, isOwned,
isGhost(localIndex));
401 using LocalIndex =
typename Comm::ParallelIndexSet::LocalIndex;
402 comm.indexSet().add(globalIndex, LocalIndex{localIndex, attr});
406 comm.indexSet().endResize();
410 Dune::OwnerOverlapCopyAttributeSet::AttributeSet
411 getAttribute_(
const Comm& comm,
const bool isOwned,
const bool isGhost)
const
414 return Dune::OwnerOverlapCopyAttributeSet::owner;
415 else if (
isGhost && (comm.category() ==
static_cast<int>(Dune::SolverCategory::nonoverlapping)) )
416 return Dune::OwnerOverlapCopyAttributeSet::overlap;
418 return Dune::OwnerOverlapCopyAttributeSet::copy;
421 const GridView gridView_;
422 const DofMapper& mapper_;
423 std::vector<std::size_t> isOwned_;
424 std::vector<std::size_t> isGhost_;
428template<
class Gr
idView,
class DofMapper,
int dofCodim>
433 : gridView_(gridView), mapper_(mapper)
437 template<
class Block,
class Alloc>
440 if constexpr (Detail::canCommunicate<typename GridView::Traits::Grid, dofCodim>)
443 if (gridView_.comm().size() > 1)
444 gridView_.communicate(gs, Dune::InteriorBorder_InteriorBorder_Interface,
445 Dune::ForwardCommunication);
448 DUNE_THROW(Dune::InvalidStateException,
"Cannot call makeNonOverlappingConsistent for a grid that cannot communicate codim-" << dofCodim <<
"-entities.");
452 const GridView gridView_;
453 const DofMapper& mapper_;
460template<
class Matrix,
class Gr
idView,
class DofMapper,
int dofCodim>
463 static constexpr int dim = GridView::dimension;
464 using Grid =
typename GridView::Traits::Grid;
465 using BlockType =
typename Matrix::block_type;
466 using IDS =
typename Grid::Traits::GlobalIdSet;
467 using IdType =
typename IDS::IdType;
490 template<
class IsGhostFunc>
491 struct MatrixPatternExchange
492 :
public Dune::CommDataHandleIF<MatrixPatternExchange<IsGhostFunc>, IdType>
495 using DataType = IdType;
497 MatrixPatternExchange(
const DofMapper& mapper,
498 const std::map<IdType,int>& globalToLocal,
499 const std::map<int,IdType>& localToGlobal,
501 std::vector<std::set<int>>& sparsity,
502 const IsGhostFunc& isGhost)
503 : mapper_(mapper), idToIndex_(globalToLocal), indexToID_(localToGlobal)
504 , sparsity_(sparsity), A_(A), isGhost_(isGhost)
506 sparsity_.resize(A.N());
512 bool contains (
int dim,
int codim)
const
513 {
return (codim == dofCodim); }
516 bool fixedSize(
int dim,
int codim)
const
522 template<
class EntityType>
523 std::size_t size(EntityType& e)
const
525 const auto rowIdx = mapper_.index(e);
527 for (
auto colIt = A_[rowIdx].begin(); colIt != A_[rowIdx].end(); ++colIt)
528 if (indexToID_.count(colIt.index()))
537 template<
class MessageBuffer,
class EntityType>
538 void gather(MessageBuffer& buff,
const EntityType& e)
const
540 const auto rowIdx = mapper_.index(e);
541 for (
auto colIt = A_[rowIdx].begin(); colIt != A_[rowIdx].end(); ++colIt)
543 auto it = indexToID_.find(colIt.index());
544 if (it != indexToID_.end())
545 buff.write(it->second);
552 template<
class MessageBuffer,
class EntityType>
553 void scatter(MessageBuffer& buff,
const EntityType& e, std::size_t n)
555 const auto rowIdx = mapper_.index(e);
556 for (std::size_t k = 0; k < n; k++)
561 const auto it = idToIndex_.find(
id);
562 if (it != idToIndex_.end())
564 const auto colIdx = it->second;
565 if (!sparsity_[rowIdx].count(colIdx) && !isGhost_(colIdx))
566 sparsity_[rowIdx].insert(colIdx);
572 const DofMapper& mapper_;
573 const std::map<IdType,int>& idToIndex_;
574 const std::map<int,IdType>& indexToID_;
575 std::vector<std::set<int>>& sparsity_;
577 const IsGhostFunc& isGhost_;
586 MatrixEntry (
const IdType& f,
const BlockType& s) : first(f), second(s) {}
591 struct MatrixEntryExchange
592 :
public Dune::CommDataHandleIF<MatrixEntryExchange, MatrixEntry>
595 using DataType = MatrixEntry;
597 MatrixEntryExchange(
const DofMapper& mapper,
598 const std::map<IdType,int>& globalToLocal,
599 const std::map<int,IdType>& localToGlobal,
601 : mapper_(mapper), idToIndex_(globalToLocal), indexToID_(localToGlobal), A_(A)
607 bool contains(
int dim,
int codim)
const
608 {
return (codim == dofCodim); }
611 bool fixedSize(
int dim,
int codim)
const
617 template<
class EntityType>
618 std::size_t size(EntityType& e)
const
620 const auto rowIdx = mapper_.index(e);
622 for (
auto colIt = A_[rowIdx].begin(); colIt != A_[rowIdx].end(); ++colIt)
623 if (indexToID_.count(colIt.index()))
632 template<
class MessageBuffer,
class EntityType>
633 void gather(MessageBuffer& buff,
const EntityType& e)
const
635 const auto rowIdx = mapper_.index(e);
636 for (
auto colIt = A_[rowIdx].begin(); colIt != A_[rowIdx].end(); ++colIt)
638 auto it = indexToID_.find(colIt.index());
639 if (it != indexToID_.end())
640 buff.write(MatrixEntry(it->second,*colIt));
647 template<
class MessageBuffer,
class EntityType>
648 void scatter(MessageBuffer& buff,
const EntityType& e, std::size_t n)
650 const auto rowIdx = mapper_.index(e);
651 for (std::size_t k = 0; k < n; k++)
656 auto it = idToIndex_.find(m.first);
657 if (it != idToIndex_.end())
658 if (A_[rowIdx].find(it->second) != A_[rowIdx].end())
659 A_[rowIdx][it->second] += m.second;
664 const DofMapper& mapper_;
665 const std::map<IdType,int>& idToIndex_;
666 const std::map<int,IdType>& indexToID_;
674 : gridView_(gridView), mapper_(mapper)
679 for (
const auto& entity : entities(gridView_, Dune::Codim<dofCodim>()))
681 if (entity.partitionType() == Dune::BorderEntity)
683 const int localIdx = mapper_.index(entity);
684 IdType dofIdxGlobal = gridView_.grid().globalIdSet().id(entity);
686 idToIndex_.emplace(dofIdxGlobal, localIdx);
687 indexToID_.emplace(localIdx, dofIdxGlobal);
697 template<
class IsGhostFunc>
700 if constexpr (Detail::canCommunicate<Grid, dofCodim>)
702 if (gridView_.comm().size() <= 1)
708 std::vector<std::set<int>> sparsity;
709 MatrixPatternExchange<IsGhostFunc> datahandle(mapper_, idToIndex_, indexToID_, A, sparsity, isGhost);
710 gridView_.communicate(datahandle, Dune::InteriorBorder_InteriorBorder_Interface,
711 Dune::ForwardCommunication);
713 for (
auto rowIt = A.begin(); rowIt != A.end(); ++rowIt)
715 for (
auto colIt = A[rowIt.index()].begin(); colIt != A[rowIt.index()].end(); ++colIt)
716 if (!sparsity[rowIt.index()].count(colIt.index()))
717 sparsity[rowIt.index()].insert(colIt.index());
719 nnz += sparsity[rowIt.index()].size();
722 A.setSize(tmp.N(), tmp.N(), nnz);
723 A.setBuildMode(Matrix::row_wise);
724 auto citer = A.createbegin();
725 for (
auto i = sparsity.begin(), end = sparsity.end(); i!=end; ++i, ++citer)
727 for (
auto si = i->begin(), send = i->end(); si!=send; ++si)
733 for (
auto rowIt = tmp.begin(); rowIt != tmp.end(); ++rowIt)
734 for (
auto colIt = tmp[rowIt.index()].begin(); colIt != tmp[rowIt.index()].end(); ++colIt)
735 A[rowIt.index()][colIt.index()] = tmp[rowIt.index()][colIt.index()];
738 DUNE_THROW(Dune::InvalidStateException,
"Cannot call extendMatrix for a grid that cannot communicate codim-" << dofCodim <<
"-entities.");
747 if constexpr (Detail::canCommunicate<Grid, dofCodim>)
749 if (gridView_.comm().size() <= 1)
752 MatrixEntryExchange datahandle(mapper_, idToIndex_, indexToID_, A);
753 gridView_.communicate(datahandle, Dune::InteriorBorder_InteriorBorder_Interface,
754 Dune::ForwardCommunication);
757 DUNE_THROW(Dune::InvalidStateException,
"Cannot call sumEntries for a grid that cannot communicate codim-" << dofCodim <<
"-entities.");
761 const GridView gridView_;
762 const DofMapper& mapper_;
763 std::map<IdType, int> idToIndex_;
764 std::map<int, IdType> indexToID_;
772 class Matrix,
class Vector,
class ParallelHelper>
775 if constexpr (ParallelTraits::isNonOverlapping)
779 using GridView =
typename LinearSolverTraits::GridView;
780 using DofMapper =
typename LinearSolverTraits::DofMapper;
781 static constexpr int dofCodim = LinearSolverTraits::dofCodim;
783 matrixHelper.
extendMatrix(A, [&pHelper](
auto idx){
return pHelper.isGhost(idx); });
795 class Matrix,
class Vector,
class ParallelHelper>
797 std::shared_ptr<typename ParallelTraits::Comm>& comm,
798 std::shared_ptr<typename ParallelTraits::LinearOperator>& fop,
799 std::shared_ptr<typename ParallelTraits::ScalarProduct>& sp,
800 ParallelHelper& pHelper)
802 prepareLinearAlgebraParallel<LinearSolverTraits, ParallelTraits>(A, b, pHelper);
803 const auto category = ParallelTraits::isNonOverlapping ?
804 Dune::SolverCategory::nonoverlapping : Dune::SolverCategory::overlapping;
806 comm = std::make_shared<typename ParallelTraits::Comm>(pHelper.gridView().comm(), category);
807 pHelper.createParallelIndexSet(*comm);
808 fop = std::make_shared<typename ParallelTraits::LinearOperator>(A, *comm);
809 sp = std::make_shared<typename ParallelTraits::ScalarProduct>(*comm);
dune-grid capabilities compatibility layer
Contains a class to exchange entries of a vector.
void prepareLinearAlgebraParallel(Matrix &A, Vector &b, ParallelHelper &pHelper)
Prepare linear algebra variables for parallel solvers.
Definition: parallelhelpers.hh:773
LinearSolverTraitsImpl< GridGeometry, typename GridGeometry::DiscretizationMethod > LinearSolverTraits
Definition: linearsolvertraits.hh:48
A parallel helper class providing a nonoverlapping decomposition of all degrees of freedom.
Definition: parallelhelpers.hh:47
bool isGhost(std::size_t i) const
Definition: parallelhelpers.hh:288
void createParallelIndexSet(Comm &comm) const
Creates a parallel index set.
Definition: parallelhelpers.hh:298
const DofMapper & dofMapper() const
Return the dofMapper.
Definition: parallelhelpers.hh:358
const GridView & gridView() const
Return the gridView.
Definition: parallelhelpers.hh:362
ParallelISTLHelper(const GridView &gridView, const DofMapper &mapper)
Definition: parallelhelpers.hh:279
Definition: parallelhelpers.hh:430
ParallelVectorHelper(const GridView &gridView, const DofMapper &mapper)
Definition: parallelhelpers.hh:432
void makeNonOverlappingConsistent(Dune::BlockVector< Block, Alloc > &v) const
Definition: parallelhelpers.hh:438
Helper class for adding up matrix entries on border.
Definition: parallelhelpers.hh:462
ParallelMatrixHelper(const GridView &gridView, const DofMapper &mapper)
Definition: parallelhelpers.hh:673
void extendMatrix(Matrix &A, const IsGhostFunc &isGhost)
communicates values for the sparsity pattern of the new matrix.
Definition: parallelhelpers.hh:698
void sumEntries(Matrix &A)
Sums up the entries corresponding to border vertices.
Definition: parallelhelpers.hh:745
A data handle class to exchange entries of a vector.
Definition: vectorcommdatahandle.hh:78