5#ifndef DUNE_ISTL_MATRIXREDISTRIBUTE_HH
6#define DUNE_ISTL_MATRIXREDISTRIBUTE_HH
9#include <dune/common/exceptions.hh>
10#include <dune/common/parallel/indexset.hh>
28 void redistribute([[maybe_unused]]
const D& from, [[maybe_unused]] D& to)
const
47 std::size_t
getRowSize([[maybe_unused]] std::size_t index)
const
65 template<
typename T,
typename T1>
72 : interface(), setup_(false)
81 const IS& target, MPI_Comm comm)
83 auto ri = std::make_unique<RemoteIndices<IS> >(source, target, comm);
84 ri->template rebuild<true>();
88 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
90 inf.build(*ri, flags, flags);
96 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
98 std::cout<<
"Interfaces do not match!"<<std::endl;
99 std::cout<<rank<<
": redist interface new :"<<inf<<std::endl;
100 std::cout<<rank<<
": redist interface :"<<interface<<std::endl;
117 template<
class GatherScatter,
class D>
120 BufferedCommunicator communicator;
121 communicator.template build<D>(from,to, interface);
122 communicator.template forward<GatherScatter>(from, to);
125 template<
class GatherScatter,
class D>
129 BufferedCommunicator communicator;
130 communicator.template build<D>(from,to, interface);
131 communicator.template backward<GatherScatter>(from, to);
155 return rowSize[index];
160 return rowSize[index];
165 return copyrowSize[index];
170 return copyrowSize[index];
175 return backwardscopyrowSize[index];
180 return backwardscopyrowSize[index];
185 rowSize.resize(rows, 0);
190 copyrowSize.resize(rows, 0);
195 backwardscopyrowSize.resize(rows, 0);
199 std::vector<std::size_t> rowSize;
200 std::vector<std::size_t> copyrowSize;
201 std::vector<std::size_t> backwardscopyrowSize;
214 template<
class M,
class RI>
243 template<
class M,
class I>
266 const std::vector<typename M::size_type>& rowsize_)
279 typedef typename Dune::GlobalLookupIndexSet<I>::const_iterator IIter;
285 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
288 if(!OwnerSet::contains(i->local().attribute())) {
290 std::cout<<rank<<
" Inserting diagonal for"<<i->local()<<std::endl;
292 sparsity[i->local()].insert(i->local());
301 m.setBuildMode(M::row_wise);
302 typename M::CreateIterator citer=m.createbegin();
306 Dune::GlobalLookupIndexSet<I> global(
aggidxset);
308 typedef typename std::vector<std::set<size_type> >::const_iterator Iter;
311 typedef typename std::set<size_type>::const_iterator SIter;
312 for(SIter si=i->begin(), send=i->end(); si!=send; ++si)
315 if(i->find(idx)==i->end()) {
316 const typename I::IndexPair* gi=global.pair(idx);
318 std::cout<<rank<<
": row "<<idx<<
" is missing a diagonal entry! global="<<gi->global()<<
" attr="<<gi->local().attribute()<<
" "<<
319 OwnerSet::contains(gi->local().attribute())<<
320 " row size="<<i->size()<<std::endl;
342 for (
unsigned int i = 0; i !=
sparsity.size(); ++i) {
343 if (add_sparsity[i].size() != 0) {
344 typedef std::set<size_type> Set;
346 std::insert_iterator<Set> tmp_insert (tmp_set, tmp_set.begin());
347 std::set_union(add_sparsity[i].begin(), add_sparsity[i].end(),
356 const Dune::GlobalLookupIndexSet<I>&
idxset;
362 template<
class M,
class I>
376 static typename M::size_type
getSize(
const Type& t, std::size_t i)
379 return t.
matrix[i].size();
394 template<
class M,
class I>
405 CommMatrixRow(M& m_,
const Dune::GlobalLookupIndexSet<I>& idxset_,
const I& aggidxset_)
412 CommMatrixRow(M& m_,
const Dune::GlobalLookupIndexSet<I>& idxset_,
const I& aggidxset_,
413 std::vector<size_t>& rowsize_)
423 typedef typename Dune::GlobalLookupIndexSet<I>::const_iterator Iter;
427 if(!OwnerSet::contains(i->local().attribute())) {
429 typedef typename M::ColIterator CIter;
430 for(CIter c=
matrix[i->local()].begin(), cend=
matrix[i->local()].end();
434 if(c.index()==i->local()) {
435 auto setDiagonal = [](
auto&& scalarOrMatrix,
const auto& value) {
436 auto&& matrixView = Dune::Impl::asMatrix(scalarOrMatrix);
437 for (
auto rowIt = matrixView.begin(); rowIt != matrixView.end(); ++rowIt)
438 (*rowIt)[rowIt.index()] = value;
448 const Dune::GlobalLookupIndexSet<I>&
idxset;
455 template<
class M,
class I>
464 typedef std::pair<typename I::GlobalIndex,typename M::block_type>
IndexedType;
472 return t.
matrix[i].size();
481 template<
class M,
class I,
class RI>
488 return cont.
matrix[i].size();
493 cont.
rowsize.getRowSize(i)=rowsize;
498 template<
class M,
class I,
class RI>
505 return cont.
matrix[i].size();
510 if (rowsize > cont.
rowsize.getCopyRowSize(i))
511 cont.
rowsize.getCopyRowSize(i)=rowsize;
516 template<
class M,
class I>
538 numlimits = std::numeric_limits<GlobalIndex>::max();
542 const typename I::IndexPair* index=cont.
idxset.pair(
col.index());
545 if ( index->local().attribute() != 2)
546 return index->global();
548 numlimits = std::numeric_limits<GlobalIndex>::max();
556 if (gi != std::numeric_limits<GlobalIndex>::max()) {
557 const typename I::IndexPair& ip=cont.
aggidxset.at(gi);
558 assert(ip.global()==gi);
559 std::size_t column = ip.local();
563 if(!OwnerSet::contains(ip.local().attribute()))
568 catch(
const Dune::RangeError&) {
572 typedef typename GlobalLookup::IndexPair IndexPair;
576 const IndexPair* pi=lookup.pair(i);
578 if(OwnerSet::contains(pi->local().attribute())) {
580 MPI_Comm_rank(MPI_COMM_WORLD,&rank);
581 std::cout<<rank<<cont.
aggidxset<<std::endl;
582 std::cout<<rank<<
": row "<<i<<
" (global="<<gi <<
") not in index set for owner index "<<pi->global()<<std::endl;
590 template<
class M,
class I>
593 template<
class M,
class I>
597 template<
class M,
class I>
603 typedef typename std::pair<GlobalIndex,typename M::block_type>
Data;
619 numlimits = std::numeric_limits<GlobalIndex>::max();
625 const typename I::IndexPair* index=cont.
idxset.pair(
col.index());
629 if ( index->local().attribute() != 2)
632 numlimits = std::numeric_limits<GlobalIndex>::max();
641 if (data.first != std::numeric_limits<GlobalIndex>::max()) {
642 typename M::size_type column=cont.
aggidxset.at(data.first).local();
643 cont.
matrix[i][column]=data.second;
646 catch(
const Dune::RangeError&) {
653 template<
class M,
class I>
656 template<
class M,
class I>
659 template<
class M,
class I>
662 template<
typename M,
typename C>
666 typename C::CopySet copyflags;
667 typename C::OwnerSet ownerflags;
668 typedef typename C::ParallelIndexSet IndexSet;
670 std::vector<typename M::size_type> rowsize(newComm.indexSet().size(), 0);
671 std::vector<typename M::size_type> copyrowsize(newComm.indexSet().size(), 0);
672 std::vector<typename M::size_type> backwardscopyrowsize(origComm.indexSet().size(), 0);
676 ri.template redistribute<MatrixRowSizeGatherScatter<M,IndexSet,RI> >(commRowSize,commRowSize);
678 origComm.buildGlobalLookup();
680 for (std::size_t i=0; i < newComm.indexSet().size(); i++) {
685 origsp(origMatrix, origComm.globalLookup(), newComm.indexSet());
687 newsp(origMatrix, origComm.globalLookup(), newComm.indexSet(), rowsize);
689 ri.template redistribute<MatrixSparsityPatternGatherScatter<M,IndexSet> >(origsp,newsp);
693 RemoteIndices<IndexSet> *ris =
new RemoteIndices<IndexSet>(origComm.indexSet(),
695 origComm.communicator());
696 ris->template rebuild<true>();
698 ri.getInterface().free();
699 ri.getInterface().build(*ris,copyflags,ownerflags);
703 ri.template redistribute<MatrixCopyRowSizeGatherScatter<M,IndexSet,RI> >(commRowSize_copy,
706 for (std::size_t i=0; i < newComm.indexSet().size(); i++) {
711 for (std::size_t i=0; i < origComm.indexSet().size(); i++) {
717 origComm.globalLookup(),
719 backwardscopyrowsize);
721 newComm.indexSet(), copyrowsize);
722 ri.template redistribute<MatrixSparsityPatternGatherScatter<M,IndexSet> >(origsp_copy,
731#ifdef DUNE_ISTL_WITH_CHECKING
734 typedef typename M::ConstRowIterator RIter;
735 for(RIter row=newMatrix.begin(), rend=newMatrix.end(); row != rend; ++row) {
736 typedef typename M::ConstColIterator CIter;
737 for(CIter
col=row->begin(), cend=row->end();
col!=cend; ++
col)
740 newMatrix[
col.index()][row.index()];
742 std::cerr<<newComm.communicator().rank()<<
": entry ("
743 <<
col.index()<<
","<<row.index()<<
") missing! for symmetry!"<<std::endl;
752 DUNE_THROW(
ISTLError,
"Matrix not symmetric!");
756 template<
typename M,
typename C>
760 typedef typename C::ParallelIndexSet IndexSet;
761 typename C::OwnerSet ownerflags;
762 std::vector<typename M::size_type> rowsize(newComm.indexSet().size(), 0);
763 std::vector<typename M::size_type> copyrowsize(newComm.indexSet().size(), 0);
764 std::vector<typename M::size_type> backwardscopyrowsize(origComm.indexSet().size(), 0);
766 for (std::size_t i=0; i < newComm.indexSet().size(); i++) {
773 for (std::size_t i=0; i < origComm.indexSet().size(); i++)
781 newComm.indexSet(), backwardscopyrowsize);
783 newComm.indexSet(),copyrowsize);
784 ri.template redistribute<MatrixRowGatherScatter<M,IndexSet> >(origrow_copy,
786 ri.getInterface().free();
787 RemoteIndices<IndexSet> *ris =
new RemoteIndices<IndexSet>(origComm.indexSet(),
789 origComm.communicator());
790 ris->template rebuild<true>();
791 ri.getInterface().build(*ris,ownerflags,ownerflags);
795 origrow(origMatrix, origComm.globalLookup(), newComm.indexSet());
797 newrow(newMatrix, origComm.globalLookup(), newComm.indexSet(),rowsize);
798 ri.template redistribute<MatrixRowGatherScatter<M,IndexSet> >(origrow,newrow);
819 template<
typename M,
typename C>
837 DUNE_THROW(InvalidStateException,
"Trying to redistribute in sequential program!");
845 DUNE_THROW(InvalidStateException,
"Trying to redistribute in sequential program!");
Classes providing communication interfaces for overlapping Schwarz methods.
Functionality for redistributing a parallel index set using graph partitioning.
Col col
Definition matrixmatrix.hh:351
Definition allocator.hh:11
void redistributeMatrixEntries(M &origMatrix, M &newMatrix, C &origComm, C &newComm, RedistributeInformation< C > &ri)
Definition matrixredistribute.hh:757
void redistributeSparsityPattern(M &origMatrix, M &newMatrix, C &origComm, C &newComm, RedistributeInformation< C > &ri)
Definition matrixredistribute.hh:663
void redistributeMatrix(M &origMatrix, M &newMatrix, C &origComm, C &newComm, RedistributeInformation< C > &ri)
Redistribute a matrix according to given domain decompositions.
Definition matrixredistribute.hh:820
derive error class from the base class in common
Definition istlexception.hh:19
Definition matrixredistribute.hh:22
void setNoBackwardsCopyRows(std::size_t size)
Definition matrixredistribute.hh:44
void redistribute(const D &from, D &to) const
Definition matrixredistribute.hh:28
void resetSetup()
Definition matrixredistribute.hh:35
void setNoCopyRows(std::size_t size)
Definition matrixredistribute.hh:41
bool isSetup() const
Definition matrixredistribute.hh:23
void setNoRows(std::size_t size)
Definition matrixredistribute.hh:38
void redistributeBackward(D &from, const D &to) const
Definition matrixredistribute.hh:32
std::size_t getBackwardsCopyRowSize(std::size_t index) const
Definition matrixredistribute.hh:57
std::size_t getRowSize(std::size_t index) const
Definition matrixredistribute.hh:47
std::size_t getCopyRowSize(std::size_t index) const
Definition matrixredistribute.hh:52
std::size_t getRowSize(std::size_t index) const
Definition matrixredistribute.hh:158
std::size_t & getBackwardsCopyRowSize(std::size_t index)
Definition matrixredistribute.hh:173
RedistributeInterface & getInterface()
Definition matrixredistribute.hh:75
void redistribute(const D &from, D &to) const
Definition matrixredistribute.hh:136
void setNoBackwardsCopyRows(std::size_t rows)
Definition matrixredistribute.hh:193
std::size_t & getCopyRowSize(std::size_t index)
Definition matrixredistribute.hh:163
RedistributeInformation()
Definition matrixredistribute.hh:71
std::size_t getCopyRowSize(std::size_t index) const
Definition matrixredistribute.hh:168
void setNoRows(std::size_t rows)
Definition matrixredistribute.hh:183
void reserve(std::size_t size)
Definition matrixredistribute.hh:150
void setNoCopyRows(std::size_t rows)
Definition matrixredistribute.hh:188
void redistributeBackward(D &from, const D &to) const
Definition matrixredistribute.hh:141
void setSetup()
Definition matrixredistribute.hh:106
OwnerOverlapCopyCommunication< T, T1 > Comm
Definition matrixredistribute.hh:69
std::size_t getBackwardsCopyRowSize(std::size_t index) const
Definition matrixredistribute.hh:178
void redistribute(const D &from, D &to) const
Definition matrixredistribute.hh:118
void resetSetup()
Definition matrixredistribute.hh:112
void redistributeBackward(D &from, const D &to) const
Definition matrixredistribute.hh:126
void checkInterface(const IS &source, const IS &target, MPI_Comm comm)
Definition matrixredistribute.hh:80
bool isSetup() const
Definition matrixredistribute.hh:145
std::size_t & getRowSize(std::size_t index)
Definition matrixredistribute.hh:153
Utility class to communicate and set the row sizes of a redistributed matrix.
Definition matrixredistribute.hh:216
M::size_type value_type
Definition matrixredistribute.hh:218
RI & rowsize
Definition matrixredistribute.hh:230
const M & matrix
Definition matrixredistribute.hh:229
M::size_type size_type
Definition matrixredistribute.hh:219
CommMatrixRowSize(const M &m_, RI &rowsize_)
Constructor.
Definition matrixredistribute.hh:226
Utility class to communicate and build the sparsity pattern of a redistributed matrix.
Definition matrixredistribute.hh:245
const Dune::GlobalLookupIndexSet< I > & idxset
Definition matrixredistribute.hh:356
M::size_type size_type
Definition matrixredistribute.hh:246
void storeSparsityPattern(M &m)
Creates and stores the sparsity pattern of the redistributed matrix.
Definition matrixredistribute.hh:276
const I & aggidxset
Definition matrixredistribute.hh:357
const std::vector< size_type > * rowsize
Definition matrixredistribute.hh:359
void completeSparsityPattern(std::vector< std::set< size_type > > add_sparsity)
Completes the sparsity pattern of the redistributed matrix with data from copy rows for the novlp cas...
Definition matrixredistribute.hh:340
CommMatrixSparsityPattern(const M &m_, const Dune::GlobalLookupIndexSet< I > &idxset_, const I &aggidxset_)
Constructor for the original side.
Definition matrixredistribute.hh:254
Dune::GlobalLookupIndexSet< I > LookupIndexSet
Definition matrixredistribute.hh:355
const M & matrix
Definition matrixredistribute.hh:354
CommMatrixSparsityPattern(const M &m_, const Dune::GlobalLookupIndexSet< I > &idxset_, const I &aggidxset_, const std::vector< typename M::size_type > &rowsize_)
Constructor for the redistruted side.
Definition matrixredistribute.hh:265
std::vector< std::set< size_type > > sparsity
Definition matrixredistribute.hh:358
VariableSize IndexedTypeFlag
Each row varies in size.
Definition matrixredistribute.hh:374
I::GlobalIndex IndexedType
The indexed type we send. This is the global index indentitfying the column.
Definition matrixredistribute.hh:371
CommMatrixSparsityPattern< M, I > Type
Definition matrixredistribute.hh:365
static M::size_type getSize(const Type &t, std::size_t i)
Definition matrixredistribute.hh:376
Utility class for comunicating the matrix entries.
Definition matrixredistribute.hh:396
std::vector< size_t > * rowsize
row size information for the receiving side.
Definition matrixredistribute.hh:452
M & matrix
The matrix to communicate the values of.
Definition matrixredistribute.hh:446
CommMatrixRow(M &m_, const Dune::GlobalLookupIndexSet< I > &idxset_, const I &aggidxset_, std::vector< size_t > &rowsize_)
Constructor.
Definition matrixredistribute.hh:412
const Dune::GlobalLookupIndexSet< I > & idxset
Index set for the original matrix.
Definition matrixredistribute.hh:448
void setOverlapRowsToDirichlet()
Sets the non-owner rows correctly as Dirichlet boundaries.
Definition matrixredistribute.hh:421
const I & aggidxset
Index set for the redistributed matrix.
Definition matrixredistribute.hh:450
CommMatrixRow(M &m_, const Dune::GlobalLookupIndexSet< I > &idxset_, const I &aggidxset_)
Constructor.
Definition matrixredistribute.hh:405
static std::size_t getSize(const Type &t, std::size_t i)
Definition matrixredistribute.hh:469
VariableSize IndexedTypeFlag
Each row varies in size.
Definition matrixredistribute.hh:467
CommMatrixRow< M, I > Type
Definition matrixredistribute.hh:458
std::pair< typename I::GlobalIndex, typename M::block_type > IndexedType
The indexed type we send. This is the pair of global index indentitfying the column and the value its...
Definition matrixredistribute.hh:464
Definition matrixredistribute.hh:483
static void scatter(Container &cont, const typename M::size_type &rowsize, std::size_t i)
Definition matrixredistribute.hh:490
static const M::size_type gather(const Container &cont, std::size_t i)
Definition matrixredistribute.hh:486
CommMatrixRowSize< M, RI > Container
Definition matrixredistribute.hh:484
Definition matrixredistribute.hh:500
static const M::size_type gather(const Container &cont, std::size_t i)
Definition matrixredistribute.hh:503
static void scatter(Container &cont, const typename M::size_type &rowsize, std::size_t i)
Definition matrixredistribute.hh:507
CommMatrixRowSize< M, RI > Container
Definition matrixredistribute.hh:501
Definition matrixredistribute.hh:518
static void scatter(Container &cont, const GlobalIndex &gi, std::size_t i, std::size_t j)
Definition matrixredistribute.hh:553
M::ConstColIterator ColIter
Definition matrixredistribute.hh:521
static GlobalIndex numlimits
Definition matrixredistribute.hh:524
static ColIter col
Definition matrixredistribute.hh:523
I::GlobalIndex GlobalIndex
Definition matrixredistribute.hh:519
CommMatrixSparsityPattern< M, I > Container
Definition matrixredistribute.hh:520
static const GlobalIndex & gather(const Container &cont, std::size_t i, std::size_t j)
Definition matrixredistribute.hh:526
Definition matrixredistribute.hh:599
static Data datastore
Definition matrixredistribute.hh:605
static GlobalIndex numlimits
Definition matrixredistribute.hh:606
M::ConstColIterator ColIter
Definition matrixredistribute.hh:602
static const Data & gather(const Container &cont, std::size_t i, std::size_t j)
Definition matrixredistribute.hh:608
I::GlobalIndex GlobalIndex
Definition matrixredistribute.hh:600
static void scatter(Container &cont, const Data &data, std::size_t i, std::size_t j)
Definition matrixredistribute.hh:638
std::pair< GlobalIndex, typename M::block_type > Data
Definition matrixredistribute.hh:603
static ColIter col
Definition matrixredistribute.hh:604
CommMatrixRow< M, I > Container
Definition matrixredistribute.hh:601
A class setting up standard communication for a two-valued attribute set with owner/overlap/copy sema...
Definition globalaggregates.hh:142
EnumItem< AttributeSet, OwnerOverlapCopyAttributeSet::owner > OwnerSet
Definition owneroverlapcopy.hh:194
Definition repartition.hh:260
@ nonoverlapping
Category for non-overlapping solvers.
Definition solvercategory.hh:27
static Category category(const OP &op, decltype(op.category()) *=nullptr)
Helperfunction to extract the solver category either from an enum, or from the newly introduced virtu...
Definition solvercategory.hh:34