3 #ifndef DUNE_ISTL_MATRIXREDISTRIBUTE_HH
4 #define DUNE_ISTL_MATRIXREDISTRIBUTE_HH
6 #include <dune/common/exceptions.hh>
7 #include <dune/common/parallel/indexset.hh>
8 #include <dune/common/unused.hh>
27 DUNE_UNUSED_PARAMETER(from);
28 DUNE_UNUSED_PARAMETER(to);
34 DUNE_UNUSED_PARAMETER(from);
35 DUNE_UNUSED_PARAMETER(to);
43 DUNE_UNUSED_PARAMETER(size);
48 DUNE_UNUSED_PARAMETER(size);
53 DUNE_UNUSED_PARAMETER(size);
58 DUNE_UNUSED_PARAMETER(index);
64 DUNE_UNUSED_PARAMETER(index);
70 DUNE_UNUSED_PARAMETER(index);
77 template<
typename T,
typename T1>
84 : interface(), setup_(false)
93 const IS& target, MPI_Comm comm)
95 RemoteIndices<IS> *ri=
new RemoteIndices<IS>(source, target, comm);
96 ri->template rebuild<true>();
100 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
102 inf.build(*ri, flags, flags);
108 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
110 std::cout<<
"Interfaces do not match!"<<std::endl;
111 std::cout<<rank<<
": redist interface new :"<<inf<<std::endl;
112 std::cout<<rank<<
": redist interface :"<<interface<<std::endl;
132 template<
class GatherScatter,
class D>
135 BufferedCommunicator communicator;
136 communicator.template build<D>(from,to, interface);
137 communicator.template forward<GatherScatter>(from, to);
140 template<
class GatherScatter,
class D>
144 BufferedCommunicator communicator;
145 communicator.template build<D>(from,to, interface);
146 communicator.template backward<GatherScatter>(from, to);
153 redistribute<CopyGatherScatter<D> >(from,to);
158 redistributeBackward<CopyGatherScatter<D> >(from,to);
170 return rowSize[index];
175 return rowSize[index];
180 return copyrowSize[index];
185 return copyrowSize[index];
190 return backwardscopyrowSize[index];
195 return backwardscopyrowSize[index];
200 rowSize.resize(rows, 0);
205 copyrowSize.resize(rows, 0);
210 backwardscopyrowSize.resize(rows, 0);
214 std::vector<std::size_t> rowSize;
215 std::vector<std::size_t> copyrowSize;
216 std::vector<std::size_t> backwardscopyrowSize;
229 template<
class M,
class RI>
258 template<
class M,
class I>
281 const std::vector<typename M::size_type>& rowsize_)
294 typedef typename Dune::GlobalLookupIndexSet<I>::const_iterator IIter;
300 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
303 if(!OwnerSet::contains(i->local().attribute())) {
305 std::cout<<rank<<
" Inserting diagonal for"<<i->local()<<std::endl;
307 sparsity[i->local()].insert(i->local());
316 m.setBuildMode(M::row_wise);
317 typename M::CreateIterator citer=m.createbegin();
321 Dune::GlobalLookupIndexSet<I> global(
aggidxset);
323 typedef typename std::vector<std::set<size_type> >::const_iterator Iter;
326 typedef typename std::set<size_type>::const_iterator SIter;
327 for(SIter si=i->begin(), send=i->end(); si!=send; ++si)
330 if(i->find(idx)==i->end()) {
331 const typename I::IndexPair* gi=global.pair(idx);
333 std::cout<<rank<<
": row "<<idx<<
" is missing a diagonal entry! global="<<gi->global()<<
" attr="<<gi->local().attribute()<<
" "<<
334 OwnerSet::contains(gi->local().attribute())<<
335 " row size="<<i->size()<<std::endl;
357 for (
unsigned int i = 0; i !=
sparsity.size(); ++i) {
358 if (add_sparsity[i].size() != 0) {
359 typedef std::set<size_type> Set;
361 std::insert_iterator<Set> tmp_insert (tmp_set, tmp_set.begin());
362 std::set_union(add_sparsity[i].begin(), add_sparsity[i].end(),
371 const Dune::GlobalLookupIndexSet<I>&
idxset;
377 template<
class M,
class I>
391 static typename M::size_type
getSize(
const Type& t, std::size_t i)
394 return t.
matrix[i].size();
409 template<
class M,
class I>
420 CommMatrixRow(M& m_,
const Dune::GlobalLookupIndexSet<I>& idxset_,
const I& aggidxset_)
427 CommMatrixRow(M& m_,
const Dune::GlobalLookupIndexSet<I>& idxset_,
const I& aggidxset_,
428 std::vector<size_t>& rowsize_)
438 typedef typename Dune::GlobalLookupIndexSet<I>::const_iterator Iter;
442 if(!OwnerSet::contains(i->local().attribute())) {
444 typedef typename M::ColIterator CIter;
445 for(CIter c=
matrix[i->local()].begin(), cend=
matrix[i->local()].end();
449 if(c.index()==i->local()) {
450 typedef typename M::block_type::RowIterator RIter;
451 for(RIter r=c->begin(), rend=c->end();
461 const Dune::GlobalLookupIndexSet<I>&
idxset;
468 template<
class M,
class I>
477 typedef std::pair<typename I::GlobalIndex,typename M::block_type>
IndexedType;
485 return t.
matrix[i].size();
494 template<
class M,
class I,
class RI>
501 return cont.
matrix[i].size();
506 cont.
rowsize.getRowSize(i)=rowsize;
511 template<
class M,
class I,
class RI>
518 return cont.
matrix[i].size();
523 if (rowsize > cont.
rowsize.getCopyRowSize(i))
524 cont.
rowsize.getCopyRowSize(i)=rowsize;
529 template<
class M,
class I>
551 numlimits = std::numeric_limits<GlobalIndex>::max();
555 const typename I::IndexPair* index=cont.
idxset.pair(
col.index());
558 if ( index->local().attribute() != 2)
559 return index->global();
561 numlimits = std::numeric_limits<GlobalIndex>::max();
568 DUNE_UNUSED_PARAMETER(j);
570 if (gi != std::numeric_limits<GlobalIndex>::max()) {
571 const typename I::IndexPair& ip=cont.
aggidxset.at(gi);
572 assert(ip.global()==gi);
573 std::size_t column = ip.local();
577 if(!OwnerSet::contains(ip.local().attribute()))
582 catch(Dune::RangeError er) {
586 typedef typename GlobalLookup::IndexPair IndexPair;
590 const IndexPair* pi=lookup.pair(i);
592 if(OwnerSet::contains(pi->local().attribute())) {
594 MPI_Comm_rank(MPI_COMM_WORLD,&rank);
595 std::cout<<rank<<cont.
aggidxset<<std::endl;
596 std::cout<<rank<<
": row "<<i<<
" (global="<<gi <<
") not in index set for owner index "<<pi->global()<<std::endl;
604 template<
class M,
class I>
607 template<
class M,
class I>
608 typename MatrixSparsityPatternGatherScatter<M,I>::GlobalIndex MatrixSparsityPatternGatherScatter<M,I>::numlimits;
611 template<
class M,
class I>
617 typedef typename std::pair<GlobalIndex,typename M::block_type>
Data;
633 numlimits = std::numeric_limits<GlobalIndex>::max();
639 const typename I::IndexPair* index=cont.
idxset.pair(
col.index());
643 if ( index->local().attribute() != 2)
646 numlimits = std::numeric_limits<GlobalIndex>::max();
654 DUNE_UNUSED_PARAMETER(j);
656 if (data.first != std::numeric_limits<GlobalIndex>::max()) {
657 typename M::size_type column=cont.
aggidxset.at(data.first).local();
658 cont.
matrix[i][column]=data.second;
661 catch(Dune::RangeError er) {
668 template<
class M,
class I>
671 template<
class M,
class I>
672 typename MatrixRowGatherScatter<M,I>::Data MatrixRowGatherScatter<M,I>::datastore;
674 template<
class M,
class I>
675 typename MatrixRowGatherScatter<M,I>::GlobalIndex MatrixRowGatherScatter<M,I>::numlimits;
677 template<
typename M,
typename C>
681 typename C::CopySet copyflags;
682 typename C::OwnerSet ownerflags;
683 typedef typename C::ParallelIndexSet IndexSet;
685 std::vector<typename M::size_type> rowsize(newComm.indexSet().size(), 0);
686 std::vector<typename M::size_type> copyrowsize(newComm.indexSet().size(), 0);
687 std::vector<typename M::size_type> backwardscopyrowsize(origComm.indexSet().size(), 0);
691 ri.template redistribute<MatrixRowSizeGatherScatter<M,IndexSet,RI> >(commRowSize,commRowSize);
693 origComm.buildGlobalLookup();
695 for (std::size_t i=0; i < newComm.indexSet().size(); i++) {
700 origsp(origMatrix, origComm.globalLookup(), newComm.indexSet());
702 newsp(origMatrix, origComm.globalLookup(), newComm.indexSet(), rowsize);
704 ri.template redistribute<MatrixSparsityPatternGatherScatter<M,IndexSet> >(origsp,newsp);
708 RemoteIndices<IndexSet> *ris =
new RemoteIndices<IndexSet>(origComm.indexSet(),
710 origComm.communicator());
711 ris->template rebuild<true>();
713 ri.getInterface().free();
714 ri.getInterface().build(*ris,copyflags,ownerflags);
718 ri.template redistribute<MatrixCopyRowSizeGatherScatter<M,IndexSet,RI> >(commRowSize_copy,
721 for (std::size_t i=0; i < newComm.indexSet().size(); i++) {
726 for (std::size_t i=0; i < origComm.indexSet().size(); i++) {
732 origComm.globalLookup(),
734 backwardscopyrowsize);
736 newComm.indexSet(), copyrowsize);
737 ri.template redistribute<MatrixSparsityPatternGatherScatter<M,IndexSet> >(origsp_copy,
740 newsp.completeSparsityPattern(newsp_copy.sparsity);
741 newsp.storeSparsityPattern(newMatrix);
744 newsp.storeSparsityPattern(newMatrix);
746 #ifdef DUNE_ISTL_WITH_CHECKING
749 typedef typename M::ConstRowIterator RIter;
750 for(RIter
row=newMatrix.begin(), rend=newMatrix.end();
row != rend; ++
row) {
751 typedef typename M::ConstColIterator CIter;
755 newMatrix[
col.index()][
row.index()];
757 std::cerr<<newComm.communicator().rank()<<
": entry ("
758 <<
col.index()<<
","<<
row.index()<<
") missing! for symmetry!"<<std::endl;
767 DUNE_THROW(
ISTLError,
"Matrix not symmetric!");
771 template<
typename M,
typename C>
775 typedef typename C::ParallelIndexSet IndexSet;
776 typename C::OwnerSet ownerflags;
777 std::vector<typename M::size_type> rowsize(newComm.indexSet().size(), 0);
778 std::vector<typename M::size_type> copyrowsize(newComm.indexSet().size(), 0);
779 std::vector<typename M::size_type> backwardscopyrowsize(origComm.indexSet().size(), 0);
781 for (std::size_t i=0; i < newComm.indexSet().size(); i++) {
788 for (std::size_t i=0; i < origComm.indexSet().size(); i++)
796 newComm.indexSet(), backwardscopyrowsize);
798 newComm.indexSet(),copyrowsize);
799 ri.template redistribute<MatrixRowGatherScatter<M,IndexSet> >(origrow_copy,
801 ri.getInterface().free();
802 RemoteIndices<IndexSet> *ris =
new RemoteIndices<IndexSet>(origComm.indexSet(),
804 origComm.communicator());
805 ris->template rebuild<true>();
806 ri.getInterface().build(*ris,ownerflags,ownerflags);
810 origrow(origMatrix, origComm.globalLookup(), newComm.indexSet());
812 newrow(newMatrix, origComm.globalLookup(), newComm.indexSet(),rowsize);
813 ri.template redistribute<MatrixRowGatherScatter<M,IndexSet> >(origrow,newrow);
815 newrow.setOverlapRowsToDirichlet();
816 if(newMatrix.N()>0&&newMatrix.N()<20)
817 printmatrix(std::cout, newMatrix,
"redist",
"row");
836 template<
typename M,
typename C>
RI & rowsize
Definition: matrixredistribute.hh:245
std::size_t & getCopyRowSize(std::size_t index)
Definition: matrixredistribute.hh:178
void redistribute(const D &from, D &to) const
Definition: matrixredistribute.hh:151
static M::size_type getSize(const Type &t, std::size_t i)
Definition: matrixredistribute.hh:391
CommMatrixSparsityPattern(const M &m_, const Dune::GlobalLookupIndexSet< I > &idxset_, const I &aggidxset_)
Constructor for the original side.
Definition: matrixredistribute.hh:269
void checkInterface(const IS &source, const IS &target, MPI_Comm comm)
Definition: matrixredistribute.hh:92
Row row
Definition: matrixmatrix.hh:345
M::ConstColIterator ColIter
Definition: matrixredistribute.hh:616
const M & matrix
Definition: matrixredistribute.hh:244
void resetSetup()
Definition: matrixredistribute.hh:127
Classes providing communication interfaces for overlapping Schwarz methods.
Definition: matrixredistribute.hh:530
static void scatter(Container &cont, const Data &data, std::size_t i, std::size_t j)
Definition: matrixredistribute.hh:652
Definition: matrixredistribute.hh:495
M::ConstColIterator ColIter
Definition: matrixredistribute.hh:534
void reserve(std::size_t size)
Definition: matrixredistribute.hh:165
RedistributeInterface & getInterface()
Definition: matrixredistribute.hh:87
Functionality for redistributing a parallel index set using graph partitioning.
void setOverlapRowsToDirichlet()
Sets the non-owner rows correctly as Dirichlet boundaries.
Definition: matrixredistribute.hh:436
derive error class from the base class in common
Definition: istlexception.hh:16
Col col
Definition: matrixmatrix.hh:347
static const Data & gather(const Container &cont, std::size_t i, std::size_t j)
Definition: matrixredistribute.hh:622
bool isSetup() const
Definition: matrixredistribute.hh:160
Definition: matrixredistribute.hh:512
static Data datastore
Definition: matrixredistribute.hh:619
I::GlobalIndex GlobalIndex
Definition: matrixredistribute.hh:614
CommMatrixRowSize< M, RI > Container
Definition: matrixredistribute.hh:497
std::size_t getRowSize(std::size_t index) const
Definition: matrixredistribute.hh:173
OwnerOverlapCopyCommunication< T, T1 > Comm
Definition: matrixredistribute.hh:81
const I & aggidxset
Index set for the redistributed matrix.
Definition: matrixredistribute.hh:463
Definition: repartition.hh:253
static const GlobalIndex & gather(const Container &cont, std::size_t i, std::size_t j)
Definition: matrixredistribute.hh:539
Definition: matrixredistribute.hh:18
static std::size_t getSize(const Type &t, std::size_t i)
Definition: matrixredistribute.hh:482
M::size_type size_type
Definition: matrixredistribute.hh:234
static GlobalIndex numlimits
Definition: matrixredistribute.hh:537
void setNoBackwardsCopyRows(std::size_t rows)
Definition: matrixredistribute.hh:208
void setNoCopyRows(std::size_t rows)
Definition: matrixredistribute.hh:203
static ColIter col
Definition: matrixredistribute.hh:536
static const M::size_type gather(const Container &cont, std::size_t i)
Definition: matrixredistribute.hh:499
Utility class to communicate and build the sparsity pattern of a redistributed matrix.
Definition: matrixredistribute.hh:259
I::GlobalIndex IndexedType
The indexed type we send. This is the global index indentitfying the column.
Definition: matrixredistribute.hh:386
CommMatrixSparsityPattern< M, I > Container
Definition: matrixredistribute.hh:533
std::size_t getRowSize(std::size_t index) const
Definition: matrixredistribute.hh:56
RedistributeInformation()
Definition: matrixredistribute.hh:83
void redistributeMatrix(M &origMatrix, M &newMatrix, C &origComm, C &newComm, RedistributeInformation< C > &ri)
Redistribute a matrix according to given domain decompositions.
Definition: matrixredistribute.hh:837
static GlobalIndex numlimits
Definition: matrixredistribute.hh:620
void setSetup()
Definition: matrixredistribute.hh:121
I::GlobalIndex GlobalIndex
Definition: matrixredistribute.hh:532
VariableSize IndexedTypeFlag
Each row varies in size.
Definition: matrixredistribute.hh:480
const I & aggidxset
Definition: matrixredistribute.hh:372
void storeSparsityPattern(M &m)
Creates and stores the sparsity pattern of the redistributed matrix.
Definition: matrixredistribute.hh:291
M::size_type value_type
Definition: matrixredistribute.hh:233
CommMatrixRowSize(const M &m_, RI &rowsize_)
Constructor.
Definition: matrixredistribute.hh:241
static ColIter col
Definition: matrixredistribute.hh:618
void redistributeMatrixEntries(M &origMatrix, M &newMatrix, C &origComm, C &newComm, RedistributeInformation< C > &ri)
Definition: matrixredistribute.hh:772
std::vector< size_t > * rowsize
row size information for the receiving side.
Definition: matrixredistribute.hh:465
static void scatter(Container &cont, const typename M::size_type &rowsize, std::size_t i)
Definition: matrixredistribute.hh:520
VariableSize IndexedTypeFlag
Each row varies in size.
Definition: matrixredistribute.hh:389
CommMatrixRow< M, I > Container
Definition: matrixredistribute.hh:615
void resetSetup()
Definition: matrixredistribute.hh:38
std::pair< GlobalIndex, typename M::block_type > Data
Definition: matrixredistribute.hh:617
const Dune::GlobalLookupIndexSet< I > & idxset
Definition: matrixredistribute.hh:371
bool isSetup() const
Definition: matrixredistribute.hh:20
EnumItem< AttributeSet, OwnerOverlapCopyAttributeSet::owner > OwnerSet
Definition: owneroverlapcopy.hh:193
void redistribute(const D &from, D &to) const
Definition: matrixredistribute.hh:133
CommMatrixRow(M &m_, const Dune::GlobalLookupIndexSet< I > &idxset_, const I &aggidxset_, std::vector< size_t > &rowsize_)
Constructor.
Definition: matrixredistribute.hh:427
CommMatrixRow< M, I > Type
Definition: matrixredistribute.hh:471
void redistributeBackward(D &from, const D &to) const
Definition: matrixredistribute.hh:32
M & matrix
The matrix to communicate the values of.
Definition: matrixredistribute.hh:459
std::size_t getCopyRowSize(std::size_t index) const
Definition: matrixredistribute.hh:183
Utility class to communicate and set the row sizes of a redistributed matrix.
Definition: matrixredistribute.hh:230
std::size_t getBackwardsCopyRowSize(std::size_t index) const
Definition: matrixredistribute.hh:193
void printmatrix(std::ostream &s, const M &A, std::string title, std::string rowtext, int width=10, int precision=2)
Print a generic block matrix.
Definition: io.hh:259
Dune::GlobalLookupIndexSet< I > LookupIndexSet
Definition: matrixredistribute.hh:370
static void scatter(Container &cont, const typename M::size_type &rowsize, std::size_t i)
Definition: matrixredistribute.hh:503
std::size_t & getBackwardsCopyRowSize(std::size_t index)
Definition: matrixredistribute.hh:188
Category for on overlapping solvers.
Definition: solvercategory.hh:23
A class setting up standard communication for a two-valued attribute set with owner/overlap/copy sema...
Definition: owneroverlapcopy.hh:172
CommMatrixRowSize< M, RI > Container
Definition: matrixredistribute.hh:514
void redistributeBackward(D &from, const D &to) const
Definition: matrixredistribute.hh:156
void redistributeBackward(D &from, const D &to) const
Definition: matrixredistribute.hh:141
void setNoBackwardsCopyRows(std::size_t size)
Definition: matrixredistribute.hh:51
CommMatrixRow(M &m_, const Dune::GlobalLookupIndexSet< I > &idxset_, const I &aggidxset_)
Constructor.
Definition: matrixredistribute.hh:420
static void scatter(Container &cont, const GlobalIndex &gi, std::size_t i, std::size_t j)
Definition: matrixredistribute.hh:566
Definition: matrixredistribute.hh:612
const Dune::GlobalLookupIndexSet< I > & idxset
Index set for the original matrix.
Definition: matrixredistribute.hh:461
const std::vector< size_type > * rowsize
Definition: matrixredistribute.hh:374
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:280
static const M::size_type gather(const Container &cont, std::size_t i)
Definition: matrixredistribute.hh:516
Utility class for comunicating the matrix entries.
Definition: matrixredistribute.hh:410
void setNoCopyRows(std::size_t size)
Definition: matrixredistribute.hh:46
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:477
void setNoRows(std::size_t size)
Definition: matrixredistribute.hh:41
void redistribute(const D &from, D &to) const
Definition: matrixredistribute.hh:25
std::size_t & getRowSize(std::size_t index)
Definition: matrixredistribute.hh:168
std::vector< std::set< size_type > > sparsity
Definition: matrixredistribute.hh:373
std::size_t getCopyRowSize(std::size_t index) const
Definition: matrixredistribute.hh:62
void setNoRows(std::size_t rows)
Definition: matrixredistribute.hh:198
void redistributeSparsityPattern(M &origMatrix, M &newMatrix, C &origComm, C &newComm, RedistributeInformation< C > &ri)
Definition: matrixredistribute.hh:678
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:355
std::size_t getBackwardsCopyRowSize(std::size_t index) const
Definition: matrixredistribute.hh:68
const M & matrix
Definition: matrixredistribute.hh:369
M::size_type size_type
Definition: matrixredistribute.hh:261
CommMatrixSparsityPattern< M, I > Type
Definition: matrixredistribute.hh:380