DUNE PDELab (2.7)

matrixredistribute.hh
Go to the documentation of this file.
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3#ifndef DUNE_ISTL_MATRIXREDISTRIBUTE_HH
4#define DUNE_ISTL_MATRIXREDISTRIBUTE_HH
5#include <memory>
6#include "repartition.hh"
11#include <dune/istl/paamg/pinfo.hh>
17namespace Dune
18{
19 template<typename T>
20 struct RedistributeInformation
21 {
22 bool isSetup() const
23 {
24 return false;
25 }
26 template<class D>
27 void redistribute(const D& from, D& to) const
28 {
31 }
32
33 template<class D>
34 void redistributeBackward(D& from, const D& to) const
35 {
38 }
39
40 void resetSetup()
41 {}
42
43 void setNoRows(std::size_t size)
44 {
46 }
47
48 void setNoCopyRows(std::size_t size)
49 {
51 }
52
53 void setNoBackwardsCopyRows(std::size_t size)
54 {
56 }
57
58 std::size_t getRowSize(std::size_t index) const
59 {
61 return -1;
62 }
63
64 std::size_t getCopyRowSize(std::size_t index) const
65 {
67 return -1;
68 }
69
70 std::size_t getBackwardsCopyRowSize(std::size_t index) const
71 {
73 return -1;
74 }
75
76 };
77
78#if HAVE_MPI
79 template<typename T, typename T1>
80 class RedistributeInformation<OwnerOverlapCopyCommunication<T,T1> >
81 {
82 public:
83 typedef OwnerOverlapCopyCommunication<T,T1> Comm;
84
85 RedistributeInformation()
86 : interface(), setup_(false)
87 {}
88
89 RedistributeInterface& getInterface()
90 {
91 return interface;
92 }
93 template<typename IS>
94 void checkInterface(const IS& source,
95 const IS& target, MPI_Comm comm)
96 {
97 auto ri = std::make_unique<RemoteIndices<IS> >(source, target, comm);
98 ri->template rebuild<true>();
99 Interface inf;
100 typename OwnerOverlapCopyCommunication<int>::OwnerSet flags;
101 int rank;
102 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
103 inf.free();
104 inf.build(*ri, flags, flags);
105
106
107#ifdef DEBUG_REPART
108 if(inf!=interface) {
109
110 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
111 if(rank==0)
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;
115
116 throw "autsch!";
117 }
118#endif
119 }
120 void setSetup()
121 {
122 setup_=true;
123 interface.strip();
124 }
125
126 void resetSetup()
127 {
128 setup_=false;
129 }
130
131 template<class GatherScatter, class D>
132 void redistribute(const D& from, D& to) const
133 {
134 BufferedCommunicator communicator;
135 communicator.template build<D>(from,to, interface);
136 communicator.template forward<GatherScatter>(from, to);
137 communicator.free();
138 }
139 template<class GatherScatter, class D>
140 void redistributeBackward(D& from, const D& to) const
141 {
142
143 BufferedCommunicator communicator;
144 communicator.template build<D>(from,to, interface);
145 communicator.template backward<GatherScatter>(from, to);
146 communicator.free();
147 }
148
149 template<class D>
150 void redistribute(const D& from, D& to) const
151 {
152 redistribute<CopyGatherScatter<D> >(from,to);
153 }
154 template<class D>
155 void redistributeBackward(D& from, const D& to) const
156 {
157 redistributeBackward<CopyGatherScatter<D> >(from,to);
158 }
159 bool isSetup() const
160 {
161 return setup_;
162 }
163
164 void reserve(std::size_t size)
165 {}
166
167 std::size_t& getRowSize(std::size_t index)
168 {
169 return rowSize[index];
170 }
171
172 std::size_t getRowSize(std::size_t index) const
173 {
174 return rowSize[index];
175 }
176
177 std::size_t& getCopyRowSize(std::size_t index)
178 {
179 return copyrowSize[index];
180 }
181
182 std::size_t getCopyRowSize(std::size_t index) const
183 {
184 return copyrowSize[index];
185 }
186
187 std::size_t& getBackwardsCopyRowSize(std::size_t index)
188 {
189 return backwardscopyrowSize[index];
190 }
191
192 std::size_t getBackwardsCopyRowSize(std::size_t index) const
193 {
194 return backwardscopyrowSize[index];
195 }
196
197 void setNoRows(std::size_t rows)
198 {
199 rowSize.resize(rows, 0);
200 }
201
202 void setNoCopyRows(std::size_t rows)
203 {
204 copyrowSize.resize(rows, 0);
205 }
206
207 void setNoBackwardsCopyRows(std::size_t rows)
208 {
209 backwardscopyrowSize.resize(rows, 0);
210 }
211
212 private:
213 std::vector<std::size_t> rowSize;
214 std::vector<std::size_t> copyrowSize;
215 std::vector<std::size_t> backwardscopyrowSize;
216 RedistributeInterface interface;
217 bool setup_;
218 };
219
228 template<class M, class RI>
230 {
231 // Make the default communication policy work.
232 typedef typename M::size_type value_type;
233 typedef typename M::size_type size_type;
234
240 CommMatrixRowSize(const M& m_, RI& rowsize_)
241 : matrix(m_), rowsize(rowsize_)
242 {}
243 const M& matrix;
244 RI& rowsize;
245
246 };
247
248
257 template<class M, class I>
259 {
260 typedef typename M::size_type size_type;
261
268 CommMatrixSparsityPattern(const M& m_, const Dune::GlobalLookupIndexSet<I>& idxset_, const I& aggidxset_)
269 : matrix(m_), idxset(idxset_), aggidxset(aggidxset_), rowsize()
270 {}
271
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_)
282 {}
283
291 {
292 // insert diagonal to overlap rows
295 std::size_t nnz=0;
296#ifdef DEBUG_REPART
297 int rank;
298
299 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
300#endif
301 for(IIter i= aggidxset.begin(), end=aggidxset.end(); i!=end; ++i) {
302 if(!OwnerSet::contains(i->local().attribute())) {
303#ifdef DEBUG_REPART
304 std::cout<<rank<<" Inserting diagonal for"<<i->local()<<std::endl;
305#endif
306 sparsity[i->local()].insert(i->local());
307 }
308
309 nnz+=sparsity[i->local()].size();
310 }
311 assert( aggidxset.size()==sparsity.size());
312
313 if(nnz>0) {
314 m.setSize(aggidxset.size(), aggidxset.size(), nnz);
315 m.setBuildMode(M::row_wise);
316 typename M::CreateIterator citer=m.createbegin();
317#ifdef DEBUG_REPART
318 std::size_t idx=0;
319 bool correct=true;
320 Dune::GlobalLookupIndexSet<I> global(aggidxset);
321#endif
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)
324 {
325 typedef typename std::set<size_type>::const_iterator SIter;
326 for(SIter si=i->begin(), send=i->end(); si!=send; ++si)
327 citer.insert(*si);
328#ifdef DEBUG_REPART
329 if(i->find(idx)==i->end()) {
330 const typename I::IndexPair* gi=global.pair(idx);
331 assert(gi);
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;
335 correct=false;
336 }
337 ++idx;
338#endif
339 }
340#ifdef DEBUG_REPART
341 if(!correct)
342 throw "bla";
343#endif
344 }
345 }
346
354 void completeSparsityPattern(std::vector<std::set<size_type> > add_sparsity)
355 {
356 for (unsigned int i = 0; i != sparsity.size(); ++i) {
357 if (add_sparsity[i].size() != 0) {
358 typedef std::set<size_type> Set;
359 Set tmp_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);
364 }
365 }
366 }
367
368 const M& matrix;
369 typedef Dune::GlobalLookupIndexSet<I> LookupIndexSet;
370 const Dune::GlobalLookupIndexSet<I>& idxset;
371 const I& aggidxset;
372 std::vector<std::set<size_type> > sparsity;
373 const std::vector<size_type>* rowsize;
374 };
375
376 template<class M, class I>
377 struct CommPolicy<CommMatrixSparsityPattern<M,I> >
378 {
379 typedef CommMatrixSparsityPattern<M,I> Type;
380
385 typedef typename I::GlobalIndex IndexedType;
386
388 typedef VariableSize IndexedTypeFlag;
389
390 static typename M::size_type getSize(const Type& t, std::size_t i)
391 {
392 if(!t.rowsize)
393 return t.matrix[i].size();
394 else
395 {
396 assert((*t.rowsize)[i]>0);
397 return (*t.rowsize)[i];
398 }
399 }
400 };
401
408 template<class M, class I>
410 {
419 CommMatrixRow(M& m_, const Dune::GlobalLookupIndexSet<I>& idxset_, const I& aggidxset_)
420 : matrix(m_), idxset(idxset_), aggidxset(aggidxset_), rowsize()
421 {}
422
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_)
429 {}
436 {
439
440 for(Iter i= aggidxset.begin(), end=aggidxset.end(); i!=end; ++i)
441 if(!OwnerSet::contains(i->local().attribute())) {
442 // Set to Dirchlet
443 typedef typename M::ColIterator CIter;
444 for(CIter c=matrix[i->local()].begin(), cend= matrix[i->local()].end();
445 c!= cend; ++c)
446 {
447 *c=0;
448 if(c.index()==i->local()) {
449 auto setDiagonal = [](auto&& scalarOrMatrix, const auto& value) {
450 auto&& matrix = Dune::Impl::asMatrix(scalarOrMatrix);
451 for (auto rowIt = matrix.begin(); rowIt != matrix.end(); ++rowIt)
452 (*rowIt)[rowIt.index()] = value;
453 };
454 setDiagonal(*c, 1);
455 }
456 }
457 }
458 }
464 const I& aggidxset;
466 std::vector<size_t>* rowsize; // row sizes differ from sender side in overlap!
467 };
468
469 template<class M, class I>
470 struct CommPolicy<CommMatrixRow<M,I> >
471 {
472 typedef CommMatrixRow<M,I> Type;
473
478 typedef std::pair<typename I::GlobalIndex,typename M::block_type> IndexedType;
479
482
483 static std::size_t getSize(const Type& t, std::size_t i)
484 {
485 if(!t.rowsize)
486 return t.matrix[i].size();
487 else
488 {
489 assert((*t.rowsize)[i]>0);
490 return (*t.rowsize)[i];
491 }
492 }
493 };
494
495 template<class M, class I, class RI>
496 struct MatrixRowSizeGatherScatter
497 {
498 typedef CommMatrixRowSize<M,RI> Container;
499
500 static const typename M::size_type gather(const Container& cont, std::size_t i)
501 {
502 return cont.matrix[i].size();
503 }
504 static void scatter(Container& cont, const typename M::size_type& rowsize, std::size_t i)
505 {
506 assert(rowsize);
507 cont.rowsize.getRowSize(i)=rowsize;
508 }
509
510 };
511
512 template<class M, class I, class RI>
513 struct MatrixCopyRowSizeGatherScatter
514 {
515 typedef CommMatrixRowSize<M,RI> Container;
516
517 static const typename M::size_type gather(const Container& cont, std::size_t i)
518 {
519 return cont.matrix[i].size();
520 }
521 static void scatter(Container& cont, const typename M::size_type& rowsize, std::size_t i)
522 {
523 assert(rowsize);
524 if (rowsize > cont.rowsize.getCopyRowSize(i))
525 cont.rowsize.getCopyRowSize(i)=rowsize;
526 }
527
528 };
529
530 template<class M, class I>
531 struct MatrixSparsityPatternGatherScatter
532 {
533 typedef typename I::GlobalIndex GlobalIndex;
534 typedef CommMatrixSparsityPattern<M,I> Container;
535 typedef typename M::ConstColIterator ColIter;
536
537 static ColIter col;
538 static GlobalIndex numlimits;
539
540 static const GlobalIndex& gather(const Container& cont, std::size_t i, std::size_t j)
541 {
542 if(j==0)
543 col=cont.matrix[i].begin();
544 else if (col!=cont.matrix[i].end())
545 ++col;
546
547 //copy communication: different row sizes for copy rows with the same global index
548 //are possible. If all values of current matrix row are sent, send
549 //std::numeric_limits<GlobalIndex>::max()
550 //and receiver will ignore it.
551 if (col==cont.matrix[i].end()) {
553 return numlimits;
554 }
555 else {
556 const typename I::IndexPair* index=cont.idxset.pair(col.index());
557 assert(index);
558 // Only send index if col is no ghost
559 if ( index->local().attribute() != 2)
560 return index->global();
561 else {
563 return numlimits;
564 }
565 }
566 }
567 static void scatter(Container& cont, const GlobalIndex& gi, std::size_t i, std::size_t j)
568 {
570 try{
572 const typename I::IndexPair& ip=cont.aggidxset.at(gi);
573 assert(ip.global()==gi);
574 std::size_t column = ip.local();
575 cont.sparsity[i].insert(column);
576
578 if(!OwnerSet::contains(ip.local().attribute()))
579 // preserve symmetry for overlap
580 cont.sparsity[column].insert(i);
581 }
582 }
583 catch(const Dune::RangeError&) {
584 // Entry not present in the new index set. Ignore!
585#ifdef DEBUG_REPART
586 typedef typename Container::LookupIndexSet GlobalLookup;
587 typedef typename GlobalLookup::IndexPair IndexPair;
589
590 GlobalLookup lookup(cont.aggidxset);
591 const IndexPair* pi=lookup.pair(i);
592 assert(pi);
593 if(OwnerSet::contains(pi->local().attribute())) {
594 int rank;
595 MPI_Comm_rank(MPI_COMM_WORLD,&rank);
596 std::cout<<rank<<cont.aggidxset<<std::endl;
597 std::cout<<rank<<": row "<<i<<" (global="<<gi <<") not in index set for owner index "<<pi->global()<<std::endl;
598 throw;
599 }
600#endif
601 }
602 }
603
604 };
605 template<class M, class I>
606 typename MatrixSparsityPatternGatherScatter<M,I>::ColIter MatrixSparsityPatternGatherScatter<M,I>::col;
607
608 template<class M, class I>
609 typename MatrixSparsityPatternGatherScatter<M,I>::GlobalIndex MatrixSparsityPatternGatherScatter<M,I>::numlimits;
610
611
612 template<class M, class I>
613 struct MatrixRowGatherScatter
614 {
615 typedef typename I::GlobalIndex GlobalIndex;
616 typedef CommMatrixRow<M,I> Container;
617 typedef typename M::ConstColIterator ColIter;
618 typedef typename std::pair<GlobalIndex,typename M::block_type> Data;
619 static ColIter col;
620 static Data datastore;
621 static GlobalIndex numlimits;
622
623 static const Data& gather(const Container& cont, std::size_t i, std::size_t j)
624 {
625 if(j==0)
626 col=cont.matrix[i].begin();
627 else if (col!=cont.matrix[i].end())
628 ++col;
629 // copy communication: different row sizes for copy rows with the same global index
630 // are possible. If all values of current matrix row are sent, send
631 // std::numeric_limits<GlobalIndex>::max()
632 // and receiver will ignore it.
633 if (col==cont.matrix[i].end()) {
635 datastore = std::make_pair(numlimits,*col);
636 return datastore;
637 }
638 else {
639 // convert local column index to global index
640 const typename I::IndexPair* index=cont.idxset.pair(col.index());
641 assert(index);
642 // Store the data to prevent reference to temporary
643 // Only send index if col is no ghost
644 if ( index->local().attribute() != 2)
645 datastore = std::make_pair(index->global(),*col);
646 else {
648 datastore = std::make_pair(numlimits,*col);
649 }
650 return datastore;
651 }
652 }
653 static void scatter(Container& cont, const Data& data, std::size_t i, std::size_t j)
654 {
656 try{
657 if (data.first != std::numeric_limits<GlobalIndex>::max()) {
658 typename M::size_type column=cont.aggidxset.at(data.first).local();
659 cont.matrix[i][column]=data.second;
660 }
661 }
662 catch(const Dune::RangeError&) {
663 // This an overlap row and might therefore lack some entries!
664 }
665
666 }
667 };
668
669 template<class M, class I>
670 typename MatrixRowGatherScatter<M,I>::ColIter MatrixRowGatherScatter<M,I>::col;
671
672 template<class M, class I>
673 typename MatrixRowGatherScatter<M,I>::Data MatrixRowGatherScatter<M,I>::datastore;
674
675 template<class M, class I>
676 typename MatrixRowGatherScatter<M,I>::GlobalIndex MatrixRowGatherScatter<M,I>::numlimits;
677
678 template<typename M, typename C>
679 void redistributeSparsityPattern(M& origMatrix, M& newMatrix, C& origComm, C& newComm,
680 RedistributeInformation<C>& ri)
681 {
682 typename C::CopySet copyflags;
683 typename C::OwnerSet ownerflags;
684 typedef typename C::ParallelIndexSet IndexSet;
685 typedef RedistributeInformation<C> RI;
686 std::vector<typename M::size_type> rowsize(newComm.indexSet().size(), 0);
687 std::vector<typename M::size_type> copyrowsize(newComm.indexSet().size(), 0);
688 std::vector<typename M::size_type> backwardscopyrowsize(origComm.indexSet().size(), 0);
689
690 // get owner rowsizes
691 CommMatrixRowSize<M,RI> commRowSize(origMatrix, ri);
692 ri.template redistribute<MatrixRowSizeGatherScatter<M,IndexSet,RI> >(commRowSize,commRowSize);
693
694 origComm.buildGlobalLookup();
695
696 for (std::size_t i=0; i < newComm.indexSet().size(); i++) {
697 rowsize[i] = ri.getRowSize(i);
698 }
699 // get sparsity pattern from owner rows
700 CommMatrixSparsityPattern<M,IndexSet>
701 origsp(origMatrix, origComm.globalLookup(), newComm.indexSet());
702 CommMatrixSparsityPattern<M,IndexSet>
703 newsp(origMatrix, origComm.globalLookup(), newComm.indexSet(), rowsize);
704
705 ri.template redistribute<MatrixSparsityPatternGatherScatter<M,IndexSet> >(origsp,newsp);
706
707 // build copy to owner interface to get missing matrix values for novlp case
709 RemoteIndices<IndexSet> *ris = new RemoteIndices<IndexSet>(origComm.indexSet(),
710 newComm.indexSet(),
711 origComm.communicator());
712 ris->template rebuild<true>();
713
714 ri.getInterface().free();
715 ri.getInterface().build(*ris,copyflags,ownerflags);
716
717 // get copy rowsizes
718 CommMatrixRowSize<M,RI> commRowSize_copy(origMatrix, ri);
719 ri.template redistribute<MatrixCopyRowSizeGatherScatter<M,IndexSet,RI> >(commRowSize_copy,
720 commRowSize_copy);
721
722 for (std::size_t i=0; i < newComm.indexSet().size(); i++) {
723 copyrowsize[i] = ri.getCopyRowSize(i);
724 }
725 //get copy rowsizes for sender
726 ri.redistributeBackward(backwardscopyrowsize,copyrowsize);
727 for (std::size_t i=0; i < origComm.indexSet().size(); i++) {
728 ri.getBackwardsCopyRowSize(i) = backwardscopyrowsize[i];
729 }
730
731 // get sparsity pattern from copy rows
732 CommMatrixSparsityPattern<M,IndexSet> origsp_copy(origMatrix,
733 origComm.globalLookup(),
734 newComm.indexSet(),
735 backwardscopyrowsize);
736 CommMatrixSparsityPattern<M,IndexSet> newsp_copy(origMatrix, origComm.globalLookup(),
737 newComm.indexSet(), copyrowsize);
738 ri.template redistribute<MatrixSparsityPatternGatherScatter<M,IndexSet> >(origsp_copy,
739 newsp_copy);
740
741 newsp.completeSparsityPattern(newsp_copy.sparsity);
742 newsp.storeSparsityPattern(newMatrix);
743 }
744 else
745 newsp.storeSparsityPattern(newMatrix);
746
747#ifdef DUNE_ISTL_WITH_CHECKING
748 // Check for symmetry
749 int ret=0;
750 typedef typename M::ConstRowIterator RIter;
751 for(RIter row=newMatrix.begin(), rend=newMatrix.end(); row != rend; ++row) {
752 typedef typename M::ConstColIterator CIter;
753 for(CIter col=row->begin(), cend=row->end(); col!=cend; ++col)
754 {
755 try{
756 newMatrix[col.index()][row.index()];
757 }catch(const Dune::ISTLError&) {
758 std::cerr<<newComm.communicator().rank()<<": entry ("
759 <<col.index()<<","<<row.index()<<") missing! for symmetry!"<<std::endl;
760 ret=1;
761
762 }
763
764 }
765 }
766
767 if(ret)
768 DUNE_THROW(ISTLError, "Matrix not symmetric!");
769#endif
770 }
771
772 template<typename M, typename C>
773 void redistributeMatrixEntries(M& origMatrix, M& newMatrix, C& origComm, C& newComm,
774 RedistributeInformation<C>& ri)
775 {
776 typedef typename C::ParallelIndexSet IndexSet;
777 typename C::OwnerSet ownerflags;
778 std::vector<typename M::size_type> rowsize(newComm.indexSet().size(), 0);
779 std::vector<typename M::size_type> copyrowsize(newComm.indexSet().size(), 0);
780 std::vector<typename M::size_type> backwardscopyrowsize(origComm.indexSet().size(), 0);
781
782 for (std::size_t i=0; i < newComm.indexSet().size(); i++) {
783 rowsize[i] = ri.getRowSize(i);
785 copyrowsize[i] = ri.getCopyRowSize(i);
786 }
787 }
788
789 for (std::size_t i=0; i < origComm.indexSet().size(); i++)
791 backwardscopyrowsize[i] = ri.getBackwardsCopyRowSize(i);
792
793
795 // fill sparsity pattern from copy rows
796 CommMatrixRow<M,IndexSet> origrow_copy(origMatrix, origComm.globalLookup(),
797 newComm.indexSet(), backwardscopyrowsize);
798 CommMatrixRow<M,IndexSet> newrow_copy(newMatrix, origComm.globalLookup(),
799 newComm.indexSet(),copyrowsize);
800 ri.template redistribute<MatrixRowGatherScatter<M,IndexSet> >(origrow_copy,
801 newrow_copy);
802 ri.getInterface().free();
803 RemoteIndices<IndexSet> *ris = new RemoteIndices<IndexSet>(origComm.indexSet(),
804 newComm.indexSet(),
805 origComm.communicator());
806 ris->template rebuild<true>();
807 ri.getInterface().build(*ris,ownerflags,ownerflags);
808 }
809
810 CommMatrixRow<M,IndexSet>
811 origrow(origMatrix, origComm.globalLookup(), newComm.indexSet());
812 CommMatrixRow<M,IndexSet>
813 newrow(newMatrix, origComm.globalLookup(), newComm.indexSet(),rowsize);
814 ri.template redistribute<MatrixRowGatherScatter<M,IndexSet> >(origrow,newrow);
815 if (SolverCategory::category(origComm) != static_cast<int>(SolverCategory::nonoverlapping))
816 newrow.setOverlapRowsToDirichlet();
817 }
818
835 template<typename M, typename C>
836 void redistributeMatrix(M& origMatrix, M& newMatrix, C& origComm, C& newComm,
837 RedistributeInformation<C>& ri)
838 {
839 ri.setNoRows(newComm.indexSet().size());
840 ri.setNoCopyRows(newComm.indexSet().size());
841 ri.setNoBackwardsCopyRows(origComm.indexSet().size());
842 redistributeSparsityPattern(origMatrix, newMatrix, origComm, newComm, ri);
843 redistributeMatrixEntries(origMatrix, newMatrix, origComm, newComm, ri);
844 }
845#endif
846
847template<typename M>
848 void redistributeMatrixEntries(M& origMatrix, M& newMatrix,
849 Dune::Amg::SequentialInformation& origComm,
850 Dune::Amg::SequentialInformation& newComm,
851 RedistributeInformation<Dune::Amg::SequentialInformation>& ri)
852 {
853 DUNE_THROW(InvalidStateException, "Trying to redistribute in sequential program!");
854 }
855 template<typename M>
856 void redistributeMatrix(M& origMatrix, M& newMatrix,
857 Dune::Amg::SequentialInformation& origComm,
858 Dune::Amg::SequentialInformation& newComm,
859 RedistributeInformation<Dune::Amg::SequentialInformation>& ri)
860 {
861 DUNE_THROW(InvalidStateException, "Trying to redistribute in sequential program!");
862 }
863}
864#endif
A constant random access iterator for the Dune::ArrayList class.
Definition: arraylist.hh:377
A set consisting only of one item.
Definition: enumset.hh:60
Decorates an index set with the possibility to find a global index that is mapped to a specific local...
Definition: indexset.hh:508
derive error class from the base class in common
Definition: istlexception.hh:16
Default exception class for range errors.
Definition: exceptions.hh:252
A few common exception classes.
const IndexPair * pair(const std::size_t &local) const
Get the index pair corresponding to a local index.
#define DUNE_UNUSED_PARAMETER(parm)
A macro to mark intentionally unused function parameters with.
Definition: unused.hh:25
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
auto max(ADLTag< 0 >, const V &v1, const V &v2)
implements binary Simd::max()
Definition: defaults.hh:79
Provides a map between global and local indices.
Dune namespace.
Definition: alignedallocator.hh:14
void redistributeMatrix(M &origMatrix, M &newMatrix, C &origComm, C &newComm, RedistributeInformation< C > &ri)
Redistribute a matrix according to given domain decompositions.
Definition: matrixredistribute.hh:836
Classes providing communication interfaces for overlapping Schwarz methods.
Functionality for redistributing a parallel index set using graph partitioning.
Utility class to communicate and set the row sizes of a redistributed matrix.
Definition: matrixredistribute.hh:230
CommMatrixRowSize(const M &m_, RI &rowsize_)
Constructor.
Definition: matrixredistribute.hh:240
Utility class for comunicating the matrix entries.
Definition: matrixredistribute.hh:410
std::vector< size_t > * rowsize
row size information for the receiving side.
Definition: matrixredistribute.hh:466
M & matrix
The matrix to communicate the values of.
Definition: matrixredistribute.hh:460
CommMatrixRow(M &m_, const Dune::GlobalLookupIndexSet< I > &idxset_, const I &aggidxset_, std::vector< size_t > &rowsize_)
Constructor.
Definition: matrixredistribute.hh:426
const Dune::GlobalLookupIndexSet< I > & idxset
Index set for the original matrix.
Definition: matrixredistribute.hh:462
void setOverlapRowsToDirichlet()
Sets the non-owner rows correctly as Dirichlet boundaries.
Definition: matrixredistribute.hh:435
const I & aggidxset
Index set for the redistributed matrix.
Definition: matrixredistribute.hh:464
CommMatrixRow(M &m_, const Dune::GlobalLookupIndexSet< I > &idxset_, const I &aggidxset_)
Constructor.
Definition: matrixredistribute.hh:419
Utility class to communicate and build the sparsity pattern of a redistributed matrix.
Definition: matrixredistribute.hh:259
void storeSparsityPattern(M &m)
Creates and stores the sparsity pattern of the redistributed matrix.
Definition: matrixredistribute.hh:290
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:354
CommMatrixSparsityPattern(const M &m_, const Dune::GlobalLookupIndexSet< I > &idxset_, const I &aggidxset_)
Constructor for the original side.
Definition: matrixredistribute.hh:268
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:279
Default policy used for communicating an indexed type.
Definition: communicator.hh:127
V::value_type IndexedType
The type we get at each index with operator[].
Definition: communicator.hh:146
static int getSize(const V &, int index)
Get the number of primitve elements at that index.
SizeOne IndexedTypeFlag
Whether the indexed type has variable size or there is always one value at each index.
Definition: communicator.hh:152
V Type
The type the policy is for.
Definition: communicator.hh:139
@ nonoverlapping
Category for non-overlapping solvers.
Definition: solvercategory.hh:25
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:32
Flag for marking indexed data structures where the data at each index may be a variable multiple of a...
Definition: communicator.hh:117
Definition of the DUNE_UNUSED macro for the case that config.h is not available.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)