3 #ifndef DUNE_ISTL_MATRIXREDISTRIBUTE_HH
4 #define DUNE_ISTL_MATRIXREDISTRIBUTE_HH
7 #include <dune/common/exceptions.hh>
8 #include <dune/common/parallel/indexset.hh>
9 #include <dune/common/unused.hh>
29 DUNE_UNUSED_PARAMETER(from);
30 DUNE_UNUSED_PARAMETER(to);
36 DUNE_UNUSED_PARAMETER(from);
37 DUNE_UNUSED_PARAMETER(to);
45 DUNE_UNUSED_PARAMETER(size);
50 DUNE_UNUSED_PARAMETER(size);
55 DUNE_UNUSED_PARAMETER(size);
60 DUNE_UNUSED_PARAMETER(index);
66 DUNE_UNUSED_PARAMETER(index);
72 DUNE_UNUSED_PARAMETER(index);
79 template<
typename T,
typename T1>
80 class RedistributeInformation<OwnerOverlapCopyCommunication<T,T1> >
83 typedef OwnerOverlapCopyCommunication<T,T1> Comm;
85 RedistributeInformation()
86 : interface(), setup_(false)
89 RedistributeInterface& getInterface()
94 void checkInterface(
const IS& source,
95 const IS& target, MPI_Comm comm)
97 auto ri = std::make_unique<RemoteIndices<IS> >(source, target, comm);
98 ri->template rebuild<true>();
100 typename OwnerOverlapCopyCommunication<int>::OwnerSet flags;
102 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
104 inf.build(*ri, flags, flags);
110 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
112 std::cout<<
"Interfaces do not match!"<<std::endl;
113 std::cout<<rank<<
": redist interface new :"<<inf<<std::endl;
114 std::cout<<rank<<
": redist interface :"<<interface<<std::endl;
131 template<
class GatherScatter,
class D>
134 BufferedCommunicator communicator;
135 communicator.template build<D>(from,to, interface);
136 communicator.template forward<GatherScatter>(from, to);
139 template<
class GatherScatter,
class D>
143 BufferedCommunicator communicator;
144 communicator.template build<D>(from,to, interface);
145 communicator.template backward<GatherScatter>(from, to);
152 redistribute<CopyGatherScatter<D> >(from,to);
157 redistributeBackward<CopyGatherScatter<D> >(from,to);
164 void reserve(std::size_t size)
169 return rowSize[index];
172 std::size_t
getRowSize(std::size_t index)
const
174 return rowSize[index];
179 return copyrowSize[index];
184 return copyrowSize[index];
189 return backwardscopyrowSize[index];
194 return backwardscopyrowSize[index];
199 rowSize.resize(rows, 0);
204 copyrowSize.resize(rows, 0);
209 backwardscopyrowSize.resize(rows, 0);
213 std::vector<std::size_t> rowSize;
214 std::vector<std::size_t> copyrowSize;
215 std::vector<std::size_t> backwardscopyrowSize;
216 RedistributeInterface interface;
228 template<
class M,
class RI>
229 struct CommMatrixRowSize
232 typedef typename M::size_type value_type;
233 typedef typename M::size_type size_type;
240 CommMatrixRowSize(
const M& m_, RI& rowsize_)
241 : matrix(m_), rowsize(rowsize_)
257 template<
class M,
class I>
258 struct CommMatrixSparsityPattern
260 typedef typename M::size_type size_type;
268 CommMatrixSparsityPattern(
const M& m_,
const Dune::GlobalLookupIndexSet<I>& idxset_,
const I& aggidxset_)
269 : matrix(m_), idxset(idxset_), aggidxset(aggidxset_), rowsize()
279 CommMatrixSparsityPattern(
const M& m_,
const Dune::GlobalLookupIndexSet<I>& idxset_,
const I& aggidxset_,
280 const std::vector<typename M::size_type>& rowsize_)
281 : matrix(m_), idxset(idxset_), aggidxset(aggidxset_), sparsity(aggidxset_.size()), rowsize(&rowsize_)
290 void storeSparsityPattern(M& m)
293 typedef typename Dune::GlobalLookupIndexSet<I>::const_iterator IIter;
294 typedef typename Dune::OwnerOverlapCopyCommunication<int>::OwnerSet OwnerSet;
299 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
301 for(IIter i= aggidxset.begin(), end=aggidxset.end(); i!=end; ++i) {
302 if(!OwnerSet::contains(i->local().attribute())) {
304 std::cout<<rank<<
" Inserting diagonal for"<<i->local()<<std::endl;
306 sparsity[i->local()].insert(i->local());
309 nnz+=sparsity[i->local()].size();
311 assert( aggidxset.size()==sparsity.size());
314 m.setSize(aggidxset.size(), aggidxset.size(), nnz);
315 m.setBuildMode(M::row_wise);
316 typename M::CreateIterator citer=m.createbegin();
320 Dune::GlobalLookupIndexSet<I> global(aggidxset);
322 typedef typename std::vector<std::set<size_type> >::const_iterator Iter;
323 for(Iter i=sparsity.begin(), end=sparsity.end(); i!=end; ++i, ++citer)
325 typedef typename std::set<size_type>::const_iterator SIter;
326 for(SIter si=i->begin(), send=i->end(); si!=send; ++si)
329 if(i->find(idx)==i->end()) {
330 const typename I::IndexPair* gi=global.pair(idx);
332 std::cout<<rank<<
": row "<<idx<<
" is missing a diagonal entry! global="<<gi->global()<<
" attr="<<gi->local().attribute()<<
" "<<
333 OwnerSet::contains(gi->local().attribute())<<
334 " row size="<<i->size()<<std::endl;
354 void completeSparsityPattern(std::vector<std::set<size_type> > add_sparsity)
356 for (
unsigned int i = 0; i != sparsity.size(); ++i) {
357 if (add_sparsity[i].size() != 0) {
358 typedef std::set<size_type> Set;
360 std::insert_iterator<Set> tmp_insert (tmp_set, tmp_set.begin());
361 std::set_union(add_sparsity[i].begin(), add_sparsity[i].end(),
362 sparsity[i].begin(), sparsity[i].end(), tmp_insert);
363 sparsity[i].swap(tmp_set);
369 typedef Dune::GlobalLookupIndexSet<I> LookupIndexSet;
370 const Dune::GlobalLookupIndexSet<I>& idxset;
372 std::vector<std::set<size_type> > sparsity;
373 const std::vector<size_type>* rowsize;
376 template<
class M,
class I>
377 struct CommPolicy<CommMatrixSparsityPattern<M,I> >
379 typedef CommMatrixSparsityPattern<M,I> Type;
385 typedef typename I::GlobalIndex IndexedType;
388 typedef VariableSize IndexedTypeFlag;
390 static typename M::size_type getSize(
const Type& t, std::size_t i)
393 return t.matrix[i].size();
396 assert((*t.rowsize)[i]>0);
397 return (*t.rowsize)[i];
408 template<
class M,
class I>
419 CommMatrixRow(M& m_,
const Dune::GlobalLookupIndexSet<I>& idxset_,
const I& aggidxset_)
420 : matrix(m_), idxset(idxset_), aggidxset(aggidxset_), rowsize()
426 CommMatrixRow(M& m_,
const Dune::GlobalLookupIndexSet<I>& idxset_,
const I& aggidxset_,
427 std::vector<size_t>& rowsize_)
428 : matrix(m_), idxset(idxset_), aggidxset(aggidxset_), rowsize(&rowsize_)
435 void setOverlapRowsToDirichlet()
437 typedef typename Dune::GlobalLookupIndexSet<I>::const_iterator Iter;
438 typedef typename Dune::OwnerOverlapCopyCommunication<int>::OwnerSet OwnerSet;
440 for(Iter i= aggidxset.begin(), end=aggidxset.end(); i!=end; ++i)
441 if(!OwnerSet::contains(i->local().attribute())) {
443 typedef typename M::ColIterator CIter;
444 for(CIter c=matrix[i->local()].begin(), cend= matrix[i->local()].end();
448 if(c.index()==i->local()) {
449 typedef typename M::block_type::RowIterator RIter;
450 for(RIter r=c->begin(), rend=c->end();
460 const Dune::GlobalLookupIndexSet<I>& idxset;
464 std::vector<size_t>* rowsize;
467 template<
class M,
class I>
468 struct CommPolicy<CommMatrixRow<M,I> >
470 typedef CommMatrixRow<M,I> Type;
476 typedef std::pair<typename I::GlobalIndex,typename M::block_type> IndexedType;
479 typedef VariableSize IndexedTypeFlag;
481 static std::size_t getSize(
const Type& t, std::size_t i)
484 return t.matrix[i].size();
487 assert((*t.rowsize)[i]>0);
488 return (*t.rowsize)[i];
493 template<
class M,
class I,
class RI>
494 struct MatrixRowSizeGatherScatter
496 typedef CommMatrixRowSize<M,RI> Container;
498 static const typename M::size_type gather(
const Container& cont, std::size_t i)
500 return cont.matrix[i].size();
502 static void scatter(Container& cont,
const typename M::size_type& rowsize, std::size_t i)
505 cont.rowsize.getRowSize(i)=rowsize;
510 template<
class M,
class I,
class RI>
511 struct MatrixCopyRowSizeGatherScatter
513 typedef CommMatrixRowSize<M,RI> Container;
515 static const typename M::size_type gather(
const Container& cont, std::size_t i)
517 return cont.matrix[i].size();
519 static void scatter(Container& cont,
const typename M::size_type& rowsize, std::size_t i)
522 if (rowsize > cont.rowsize.getCopyRowSize(i))
523 cont.rowsize.getCopyRowSize(i)=rowsize;
528 template<
class M,
class I>
529 struct MatrixSparsityPatternGatherScatter
531 typedef typename I::GlobalIndex GlobalIndex;
532 typedef CommMatrixSparsityPattern<M,I> Container;
533 typedef typename M::ConstColIterator ColIter;
536 static GlobalIndex numlimits;
538 static const GlobalIndex& gather(
const Container& cont, std::size_t i, std::size_t j)
541 col=cont.matrix[i].begin();
542 else if (
col!=cont.matrix[i].end())
549 if (
col==cont.matrix[i].end()) {
550 numlimits = std::numeric_limits<GlobalIndex>::max();
554 const typename I::IndexPair* index=cont.idxset.pair(
col.index());
557 if ( index->local().attribute() != 2)
558 return index->global();
560 numlimits = std::numeric_limits<GlobalIndex>::max();
565 static void scatter(Container& cont,
const GlobalIndex& gi, std::size_t i, std::size_t j)
567 DUNE_UNUSED_PARAMETER(j);
569 if (gi != std::numeric_limits<GlobalIndex>::max()) {
570 const typename I::IndexPair& ip=cont.aggidxset.at(gi);
571 assert(ip.global()==gi);
572 std::size_t column = ip.local();
573 cont.sparsity[i].insert(column);
575 typedef typename Dune::OwnerOverlapCopyCommunication<int>::OwnerSet OwnerSet;
576 if(!OwnerSet::contains(ip.local().attribute()))
578 cont.sparsity[column].insert(i);
581 catch(
const Dune::RangeError&) {
584 typedef typename Container::LookupIndexSet GlobalLookup;
585 typedef typename GlobalLookup::IndexPair IndexPair;
586 typedef typename Dune::OwnerOverlapCopyCommunication<int>::OwnerSet OwnerSet;
588 GlobalLookup lookup(cont.aggidxset);
589 const IndexPair* pi=lookup.pair(i);
591 if(OwnerSet::contains(pi->local().attribute())) {
593 MPI_Comm_rank(MPI_COMM_WORLD,&rank);
594 std::cout<<rank<<cont.aggidxset<<std::endl;
595 std::cout<<rank<<
": row "<<i<<
" (global="<<gi <<
") not in index set for owner index "<<pi->global()<<std::endl;
603 template<
class M,
class I>
606 template<
class M,
class I>
607 typename MatrixSparsityPatternGatherScatter<M,I>::GlobalIndex MatrixSparsityPatternGatherScatter<M,I>::numlimits;
610 template<
class M,
class I>
611 struct MatrixRowGatherScatter
613 typedef typename I::GlobalIndex GlobalIndex;
614 typedef CommMatrixRow<M,I> Container;
615 typedef typename M::ConstColIterator ColIter;
616 typedef typename std::pair<GlobalIndex,typename M::block_type> Data;
618 static Data datastore;
619 static GlobalIndex numlimits;
621 static const Data& gather(
const Container& cont, std::size_t i, std::size_t j)
624 col=cont.matrix[i].begin();
625 else if (
col!=cont.matrix[i].end())
631 if (
col==cont.matrix[i].end()) {
632 numlimits = std::numeric_limits<GlobalIndex>::max();
633 datastore = std::make_pair(numlimits,*
col);
638 const typename I::IndexPair* index=cont.idxset.pair(
col.index());
642 if ( index->local().attribute() != 2)
643 datastore = std::make_pair(index->global(),*
col);
645 numlimits = std::numeric_limits<GlobalIndex>::max();
646 datastore = std::make_pair(numlimits,*
col);
651 static void scatter(Container& cont,
const Data& data, std::size_t i, std::size_t j)
653 DUNE_UNUSED_PARAMETER(j);
655 if (data.first != std::numeric_limits<GlobalIndex>::max()) {
656 typename M::size_type column=cont.aggidxset.at(data.first).local();
657 cont.matrix[i][column]=data.second;
660 catch(
const Dune::RangeError&) {
667 template<
class M,
class I>
670 template<
class M,
class I>
671 typename MatrixRowGatherScatter<M,I>::Data MatrixRowGatherScatter<M,I>::datastore;
673 template<
class M,
class I>
674 typename MatrixRowGatherScatter<M,I>::GlobalIndex MatrixRowGatherScatter<M,I>::numlimits;
676 template<
typename M,
typename C>
677 void redistributeSparsityPattern(M& origMatrix, M& newMatrix, C& origComm, C& newComm,
678 RedistributeInformation<C>& ri)
680 typename C::CopySet copyflags;
681 typename C::OwnerSet ownerflags;
682 typedef typename C::ParallelIndexSet IndexSet;
683 typedef RedistributeInformation<C> RI;
684 std::vector<typename M::size_type> rowsize(newComm.indexSet().size(), 0);
685 std::vector<typename M::size_type> copyrowsize(newComm.indexSet().size(), 0);
686 std::vector<typename M::size_type> backwardscopyrowsize(origComm.indexSet().size(), 0);
689 CommMatrixRowSize<M,RI> commRowSize(origMatrix, ri);
690 ri.template redistribute<MatrixRowSizeGatherScatter<M,IndexSet,RI> >(commRowSize,commRowSize);
692 origComm.buildGlobalLookup();
694 for (std::size_t i=0; i < newComm.indexSet().size(); i++) {
695 rowsize[i] = ri.getRowSize(i);
698 CommMatrixSparsityPattern<M,IndexSet>
699 origsp(origMatrix, origComm.globalLookup(), newComm.indexSet());
700 CommMatrixSparsityPattern<M,IndexSet>
701 newsp(origMatrix, origComm.globalLookup(), newComm.indexSet(), rowsize);
703 ri.template redistribute<MatrixSparsityPatternGatherScatter<M,IndexSet> >(origsp,newsp);
707 RemoteIndices<IndexSet> *ris =
new RemoteIndices<IndexSet>(origComm.indexSet(),
709 origComm.communicator());
710 ris->template rebuild<true>();
712 ri.getInterface().free();
713 ri.getInterface().build(*ris,copyflags,ownerflags);
716 CommMatrixRowSize<M,RI> commRowSize_copy(origMatrix, ri);
717 ri.template redistribute<MatrixCopyRowSizeGatherScatter<M,IndexSet,RI> >(commRowSize_copy,
720 for (std::size_t i=0; i < newComm.indexSet().size(); i++) {
721 copyrowsize[i] = ri.getCopyRowSize(i);
724 ri.redistributeBackward(backwardscopyrowsize,copyrowsize);
725 for (std::size_t i=0; i < origComm.indexSet().size(); i++) {
726 ri.getBackwardsCopyRowSize(i) = backwardscopyrowsize[i];
730 CommMatrixSparsityPattern<M,IndexSet> origsp_copy(origMatrix,
731 origComm.globalLookup(),
733 backwardscopyrowsize);
734 CommMatrixSparsityPattern<M,IndexSet> newsp_copy(origMatrix, origComm.globalLookup(),
735 newComm.indexSet(), copyrowsize);
736 ri.template redistribute<MatrixSparsityPatternGatherScatter<M,IndexSet> >(origsp_copy,
739 newsp.completeSparsityPattern(newsp_copy.sparsity);
740 newsp.storeSparsityPattern(newMatrix);
743 newsp.storeSparsityPattern(newMatrix);
745 #ifdef DUNE_ISTL_WITH_CHECKING
748 typedef typename M::ConstRowIterator RIter;
749 for(RIter row=newMatrix.begin(), rend=newMatrix.end(); row != rend; ++row) {
750 typedef typename M::ConstColIterator CIter;
751 for(CIter
col=row->begin(), cend=row->end();
col!=cend; ++
col)
754 newMatrix[
col.index()][row.index()];
756 std::cerr<<newComm.communicator().rank()<<
": entry ("
757 <<
col.index()<<
","<<row.index()<<
") missing! for symmetry!"<<std::endl;
766 DUNE_THROW(ISTLError,
"Matrix not symmetric!");
770 template<
typename M,
typename C>
772 RedistributeInformation<C>& ri)
774 typedef typename C::ParallelIndexSet IndexSet;
775 typename C::OwnerSet ownerflags;
776 std::vector<typename M::size_type> rowsize(newComm.indexSet().size(), 0);
777 std::vector<typename M::size_type> copyrowsize(newComm.indexSet().size(), 0);
778 std::vector<typename M::size_type> backwardscopyrowsize(origComm.indexSet().size(), 0);
780 for (std::size_t i=0; i < newComm.indexSet().size(); i++) {
781 rowsize[i] = ri.getRowSize(i);
783 copyrowsize[i] = ri.getCopyRowSize(i);
787 for (std::size_t i=0; i < origComm.indexSet().size(); i++)
789 backwardscopyrowsize[i] = ri.getBackwardsCopyRowSize(i);
794 CommMatrixRow<M,IndexSet> origrow_copy(origMatrix, origComm.globalLookup(),
795 newComm.indexSet(), backwardscopyrowsize);
796 CommMatrixRow<M,IndexSet> newrow_copy(newMatrix, origComm.globalLookup(),
797 newComm.indexSet(),copyrowsize);
798 ri.template redistribute<MatrixRowGatherScatter<M,IndexSet> >(origrow_copy,
800 ri.getInterface().free();
801 RemoteIndices<IndexSet> *ris =
new RemoteIndices<IndexSet>(origComm.indexSet(),
803 origComm.communicator());
804 ris->template rebuild<true>();
805 ri.getInterface().build(*ris,ownerflags,ownerflags);
808 CommMatrixRow<M,IndexSet>
809 origrow(origMatrix, origComm.globalLookup(), newComm.indexSet());
810 CommMatrixRow<M,IndexSet>
811 newrow(newMatrix, origComm.globalLookup(), newComm.indexSet(),rowsize);
812 ri.template redistribute<MatrixRowGatherScatter<M,IndexSet> >(origrow,newrow);
814 newrow.setOverlapRowsToDirichlet();
833 template<
typename M,
typename C>
835 RedistributeInformation<C>& ri)
837 ri.setNoRows(newComm.indexSet().size());
838 ri.setNoCopyRows(newComm.indexSet().size());
839 ri.setNoBackwardsCopyRows(origComm.indexSet().size());
840 redistributeSparsityPattern(origMatrix, newMatrix, origComm, newComm, ri);
851 DUNE_THROW(InvalidStateException,
"Trying to redistribute in sequential program!");
859 DUNE_THROW(InvalidStateException,
"Trying to redistribute in sequential program!");