dune-istl  2.2.1
matrixredistribute.hh
Go to the documentation of this file.
1 // -*- tab-width: 8; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2 // vi: set ts=8 sw=2 et sts=2:
3 #ifndef DUNE_MATRIXREDIST_HH
4 #define DUNE_MATRIXREDIST_HH
5 #include"repartition.hh"
6 #include<dune/common/exceptions.hh>
7 #include<dune/common/parallel/indexset.hh>
14 namespace Dune
15 {
16  template<typename T>
18  {
19  bool isSetup() const
20  {
21  return false;
22  }
23  template<class D>
24  void redistribute(const D& from, D& to) const
25  {}
26 
27  template<class D>
28  void redistributeBackward(D& from, const D& to) const
29  {}
30 
31  void resetSetup()
32  {}
33 
34  void setNoRows(std::size_t size)
35  {}
36 
37  void setNoCopyRows(std::size_t size)
38  {}
39 
40  void setNoBackwardsCopyRows(std::size_t size)
41  {}
42 
43  std::size_t getRowSize(std::size_t index) const
44  {
45  return -1;
46  }
47 
48  std::size_t getCopyRowSize(std::size_t index) const
49  {
50  return -1;
51  }
52 
53  std::size_t getBackwardsCopyRowSize(std::size_t index) const
54  {
55  return -1;
56  }
57 
58  };
59 
60 #if HAVE_MPI
61  template<typename T, typename T1>
63  {
64  public:
66 
68  : interface(), setup_(false)
69  {}
70 
71  RedistributeInterface& getInterface()
72  {
73  return interface;
74  }
75  template<typename IS>
76  void checkInterface(const IS& source,
77  const IS& target, MPI_Comm comm)
78  {
79  RemoteIndices<IS> *ri=new RemoteIndices<IS>(source, target, comm);
80  ri->template rebuild<true>();
81  Interface inf;
83  int rank;
84  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
85  inf.free();
86  inf.build(*ri, flags, flags);
87 
88 
89 #ifdef DEBUG_REPART
90  if(inf!=interface){
91 
92  int rank;
93  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
94  if(rank==0)
95  std::cout<<"Interfaces do not match!"<<std::endl;
96  std::cout<<rank<<": redist interface new :"<<inf<<std::endl;
97  std::cout<<rank<<": redist interface :"<<interface<<std::endl;
98 
99  throw "autsch!";
100  delete ri;
101  }else
102 
103 #endif
104  delete ri;
105  }
106  void setSetup()
107  {
108  setup_=true;
109  interface.strip();
110  }
111 
112  void resetSetup()
113  {
114  setup_=false;
115  }
116 
117  template<class GatherScatter, class D>
118  void redistribute(const D& from, D& to) const
119  {
120  BufferedCommunicator communicator;
121  communicator.template build<D>(from,to, interface);
122  communicator.template forward<GatherScatter>(from, to);
123  communicator.free();
124  }
125  template<class GatherScatter, class D>
126  void redistributeBackward(D& from, const D& to) const
127  {
128 
129  BufferedCommunicator communicator;
130  communicator.template build<D>(from,to, interface);
131  communicator.template backward<GatherScatter>(from, to);
132  communicator.free();
133  }
134 
135  template<class D>
136  void redistribute(const D& from, D& to) const
137  {
138  redistribute<CopyGatherScatter<D> >(from,to);
139  }
140  template<class D>
141  void redistributeBackward(D& from, const D& to) const
142  {
143  redistributeBackward<CopyGatherScatter<D> >(from,to);
144  }
145  bool isSetup() const
146  {
147  return setup_;
148  }
149 
150  void reserve(std::size_t size)
151  {}
152 
153  std::size_t& getRowSize(std::size_t index)
154  {
155  return rowSize[index];
156  }
157 
158  std::size_t getRowSize(std::size_t index) const
159  {
160  return rowSize[index];
161  }
162 
163  std::size_t& getCopyRowSize(std::size_t index)
164  {
165  return copyrowSize[index];
166  }
167 
168  std::size_t getCopyRowSize(std::size_t index) const
169  {
170  return copyrowSize[index];
171  }
172 
173  std::size_t& getBackwardsCopyRowSize(std::size_t index)
174  {
175  return backwardscopyrowSize[index];
176  }
177 
178  std::size_t getBackwardsCopyRowSize(std::size_t index) const
179  {
180  return backwardscopyrowSize[index];
181  }
182 
183  void setNoRows(std::size_t rows)
184  {
185  rowSize.resize(rows, 0);
186  }
187 
188  void setNoCopyRows(std::size_t rows)
189  {
190  copyrowSize.resize(rows, 0);
191  }
192 
193  void setNoBackwardsCopyRows(std::size_t rows)
194  {
195  backwardscopyrowSize.resize(rows, 0);
196  }
197 
198  private:
199  std::vector<std::size_t> rowSize;
200  std::vector<std::size_t> copyrowSize;
201  std::vector<std::size_t> backwardscopyrowSize;
202  RedistributeInterface interface;
203  bool setup_;
204  };
205 
214  template<class M, class RI>
216  {
217  // Make the default communication policy work.
218  typedef typename M::size_type value_type;
219  typedef typename M::size_type size_type;
220 
226  CommMatrixRowSize(const M& m_, RI& rowsize_)
227  : matrix(m_), rowsize(rowsize_)
228  {}
229  const M& matrix;
230  RI& rowsize;
231 
232  };
233 
234 
243  template<class M, class I>
245  {
246  typedef typename M::size_type size_type;
247 
254  CommMatrixSparsityPattern(const M& m_, const Dune::GlobalLookupIndexSet<I>& idxset_, const I& aggidxset_)
255  : matrix(m_), idxset(idxset_), aggidxset(aggidxset_), rowsize()
256  {}
257 
265  CommMatrixSparsityPattern(const M& m_, const Dune::GlobalLookupIndexSet<I>& idxset_, const I& aggidxset_,
266  const std::vector<typename M::size_type>& rowsize_)
267  : matrix(m_), idxset(idxset_), aggidxset(aggidxset_), sparsity(aggidxset_.size()), rowsize(&rowsize_)
268  {}
269 
277  {
278  // insert diagonal to overlap rows
279  typedef typename Dune::GlobalLookupIndexSet<I>::const_iterator IIter;
280  typedef typename Dune::OwnerOverlapCopyCommunication<int>::OwnerSet OwnerSet;
281  std::size_t nnz=0;
282 #ifdef DEBUG_REPART
283  int rank;
284 
285  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
286 #endif
287  for(IIter i= aggidxset.begin(), end=aggidxset.end(); i!=end; ++i){
288  if(!OwnerSet::contains(i->local().attribute())){
289 #ifdef DEBUG_REPART
290  std::cout<<rank<<" Inserting diagonal for"<<i->local()<<std::endl;
291 #endif
292  sparsity[i->local()].insert(i->local());
293  }
294 
295  nnz+=sparsity[i->local()].size();
296  }
297  assert( aggidxset.size()==sparsity.size());
298 
299  if(nnz>0){
300  m.setSize(aggidxset.size(), aggidxset.size(), nnz);
301  m.setBuildMode(M::row_wise);
302  typename M::CreateIterator citer=m.createbegin();
303 #ifdef DEBUG_REPART
304  std::size_t idx=0;
305  bool correct=true;
306  Dune::GlobalLookupIndexSet<I> global(aggidxset);
307 #endif
308  typedef typename std::vector<std::set<size_type> >::const_iterator Iter;
309  for(Iter i=sparsity.begin(), end=sparsity.end(); i!=end; ++i, ++citer)
310  {
311  typedef typename std::set<size_type>::const_iterator SIter;
312  for(SIter si=i->begin(), send=i->end(); si!=send; ++si)
313  citer.insert(*si);
314 #ifdef DEBUG_REPART
315  if(i->find(idx)==i->end()){
316  const typename I::IndexPair* gi=global.pair(idx);
317  assert(gi);
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;
321  correct=false;
322  }
323  ++idx;
324 #endif
325  }
326 #ifdef DEBUG_REPART
327  if(!correct)
328  throw "bla";
329 #endif
330  }
331  }
332 
340  void completeSparsityPattern(std::vector<std::set<size_type> > add_sparsity)
341  {
342  for (unsigned int i = 0; i != sparsity.size(); ++i) {
343  if (add_sparsity[i].size() != 0){
344  typedef std::set<size_type> Set;
345  Set tmp_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(),
348  sparsity[i].begin(), sparsity[i].end(), tmp_insert);
349  sparsity[i].swap(tmp_set);
350  }
351  }
352  }
353 
354  const M& matrix;
355  typedef Dune::GlobalLookupIndexSet<I> LookupIndexSet;
356  const Dune::GlobalLookupIndexSet<I>& idxset;
357  const I& aggidxset;
358  std::vector<std::set<size_type> > sparsity;
359  const std::vector<size_type>* rowsize;
360  };
361 
362  template<class M, class I>
363  struct CommPolicy<CommMatrixSparsityPattern<M,I> >
364  {
366 
371  typedef typename I::GlobalIndex IndexedType;
372 
374  typedef VariableSize IndexedTypeFlag;
375 
376  static typename M::size_type getSize(const Type& t, std::size_t i)
377  {
378  if(!t.rowsize)
379  return t.matrix[i].size();
380  else
381  {
382  assert((*t.rowsize)[i]>0);
383  return (*t.rowsize)[i];
384  }
385  }
386  };
387 
394  template<class M, class I>
396  {
405  CommMatrixRow(M& m_, const Dune::GlobalLookupIndexSet<I>& idxset_, const I& aggidxset_)
406  : matrix(m_), idxset(idxset_), aggidxset(aggidxset_), rowsize()
407  {}
408 
412  CommMatrixRow(M& m_, const Dune::GlobalLookupIndexSet<I>& idxset_, const I& aggidxset_,
413  std::vector<size_t>& rowsize_)
414  : matrix(m_), idxset(idxset_), aggidxset(aggidxset_), rowsize(&rowsize_)
415  {}
422  {
423  typedef typename Dune::GlobalLookupIndexSet<I>::const_iterator Iter;
424  typedef typename Dune::OwnerOverlapCopyCommunication<int>::OwnerSet OwnerSet;
425 
426  for(Iter i= aggidxset.begin(), end=aggidxset.end(); i!=end; ++i)
427  if(!OwnerSet::contains(i->local().attribute())){
428  // Set to Dirchlet
429  typedef typename M::ColIterator CIter;
430  for(CIter c=matrix[i->local()].begin(), cend= matrix[i->local()].end();
431  c!= cend; ++c)
432  {
433  *c=0;
434  if(c.index()==i->local()){
435  typedef typename M::block_type::RowIterator RIter;
436  for(RIter r=c->begin(), rend=c->end();
437  r != rend; ++r)
438  (*r)[r.index()]=1;
439  }
440  }
441  }
442  }
444  M& matrix;
446  const Dune::GlobalLookupIndexSet<I>& idxset;
448  const I& aggidxset;
450  std::vector<size_t>* rowsize; // row sizes differ from sender side in overlap!
451  };
452 
453  template<class M, class I>
454  struct CommPolicy<CommMatrixRow<M,I> >
455  {
457 
462  typedef std::pair<typename I::GlobalIndex,typename M::block_type> IndexedType;
463 
465  typedef VariableSize IndexedTypeFlag;
466 
467  static std::size_t getSize(const Type& t, std::size_t i)
468  {
469  if(!t.rowsize)
470  return t.matrix[i].size();
471  else
472  {
473  assert((*t.rowsize)[i]>0);
474  return (*t.rowsize)[i];
475  }
476  }
477  };
478 
479  template<class M, class I, class RI>
481  {
483 
484  static const typename M::size_type gather(const Container& cont, std::size_t i)
485  {
486  return cont.matrix[i].size();
487  }
488  static void scatter(Container& cont, const typename M::size_type& rowsize, std::size_t i)
489  {
490  assert(rowsize);
491  cont.rowsize.getRowSize(i)=rowsize;
492  }
493 
494  };
495 
496  template<class M, class I, class RI>
498  {
500 
501  static const typename M::size_type gather(const Container& cont, std::size_t i)
502  {
503  return cont.matrix[i].size();
504  }
505  static void scatter(Container& cont, const typename M::size_type& rowsize, std::size_t i)
506  {
507  assert(rowsize);
508  if (rowsize > cont.rowsize.getCopyRowSize(i))
509  cont.rowsize.getCopyRowSize(i)=rowsize;
510  }
511 
512  };
513 
514  template<class M, class I>
516  {
517  typedef typename I::GlobalIndex GlobalIndex;
519  typedef typename M::ConstColIterator ColIter;
520 
521  static ColIter col;
523 
524  static const GlobalIndex& gather(const Container& cont, std::size_t i, std::size_t j)
525  {
526  if(j==0)
527  col=cont.matrix[i].begin();
528  else if (col!=cont.matrix[i].end())
529  ++col;
530 
531  //copy communication: different row sizes for copy rows with the same global index
532  //are possible. If all values of current matrix row are sent, send
533  //std::numeric_limits<GlobalIndex>::max()
534  //and receiver will ignore it.
535  if (col==cont.matrix[i].end()){
536  numlimits = std::numeric_limits<GlobalIndex>::max();
537  return numlimits;
538  }
539  else {
540  const typename I::IndexPair* index=cont.idxset.pair(col.index());
541  assert(index);
542  // Only send index if col is no ghost
543  if ( index->local().attribute() != 2)
544  return index->global();
545  else {
546  numlimits = std::numeric_limits<GlobalIndex>::max();
547  return numlimits;
548  }
549  }
550  }
551  static void scatter(Container& cont, const GlobalIndex& gi, std::size_t i, std::size_t j)
552  {
553  try{
554  if (gi != std::numeric_limits<GlobalIndex>::max()){
555  const typename I::IndexPair& ip=cont.aggidxset.at(gi);
556  assert(ip.global()==gi);
557  std::size_t col = ip.local();
558  cont.sparsity[i].insert(col);
559 
560  typedef typename Dune::OwnerOverlapCopyCommunication<int>::OwnerSet OwnerSet;
561  if(!OwnerSet::contains(ip.local().attribute()))
562  // preserve symmetry for overlap
563  cont.sparsity[col].insert(i);
564  }
565  }
566  catch(Dune::RangeError er){
567  // Entry not present in the new index set. Ignore!
568 #ifdef DEBUG_REPART
569  typedef typename Container::LookupIndexSet GlobalLookup;
570  typedef typename GlobalLookup::IndexPair IndexPair;
571  typedef typename Dune::OwnerOverlapCopyCommunication<int>::OwnerSet OwnerSet;
572 
573  GlobalLookup lookup(cont.aggidxset);
574  const IndexPair* pi=lookup.pair(i);
575  assert(pi);
576  if(OwnerSet::contains(pi->local().attribute())){
577  int rank;
578  MPI_Comm_rank(MPI_COMM_WORLD,&rank);
579  std::cout<<rank<<cont.aggidxset<<std::endl;
580  std::cout<<rank<<": row "<<i<<" (global="<<gi <<") not in index set for owner index "<<pi->global()<<std::endl;
581  throw er;
582  }
583 #endif
584  }
585  }
586 
587  };
588  template<class M, class I>
589  typename MatrixSparsityPatternGatherScatter<M,I>::ColIter MatrixSparsityPatternGatherScatter<M,I>::col;
590 
591  template<class M, class I>
592  typename MatrixSparsityPatternGatherScatter<M,I>::GlobalIndex MatrixSparsityPatternGatherScatter<M,I>::numlimits;
593 
594 
595  template<class M, class I>
597  {
598  typedef typename I::GlobalIndex GlobalIndex;
600  typedef typename M::ConstColIterator ColIter;
601  typedef typename std::pair<GlobalIndex,typename M::block_type> Data;
602  static ColIter col;
603  static Data datastore;
605 
606  static const Data& gather(const Container& cont, std::size_t i, std::size_t j)
607  {
608  if(j==0)
609  col=cont.matrix[i].begin();
610  else if (col!=cont.matrix[i].end())
611  ++col;
612  // copy communication: different row sizes for copy rows with the same global index
613  // are possible. If all values of current matrix row are sent, send
614  // std::numeric_limits<GlobalIndex>::max()
615  // and receiver will ignore it.
616  if (col==cont.matrix[i].end()){
617  numlimits = std::numeric_limits<GlobalIndex>::max();
618  datastore = std::make_pair(numlimits,*col);
619  return datastore;
620  }
621  else {
622  // convert local column index to global index
623  const typename I::IndexPair* index=cont.idxset.pair(col.index());
624  assert(index);
625  // Store the data to prevent reference to temporary
626  // Only send index if col is no ghost
627  if ( index->local().attribute() != 2)
628  datastore = std::make_pair(index->global(),*col);
629  else {
630  numlimits = std::numeric_limits<GlobalIndex>::max();
631  datastore = std::make_pair(numlimits,*col);
632  }
633  return datastore;
634  }
635  }
636  static void scatter(Container& cont, const Data& data, std::size_t i, std::size_t j)
637  {
638  try{
639  if (data.first != std::numeric_limits<GlobalIndex>::max()){
640  typename M::size_type column=cont.aggidxset.at(data.first).local();
641  cont.matrix[i][column]=data.second;
642  }
643  }
644  catch(Dune::RangeError er){
645  // This an overlap row and might therefore lack some entries!
646  }
647 
648  }
649  };
650 
651  template<class M, class I>
652  typename MatrixRowGatherScatter<M,I>::ColIter MatrixRowGatherScatter<M,I>::col;
653 
654  template<class M, class I>
655  typename MatrixRowGatherScatter<M,I>::Data MatrixRowGatherScatter<M,I>::datastore;
656 
657  template<class M, class I>
658  typename MatrixRowGatherScatter<M,I>::GlobalIndex MatrixRowGatherScatter<M,I>::numlimits;
659 
660  template<typename M, typename C>
661  void redistributeSparsityPattern(M& origMatrix, M& newMatrix, C& origComm, C& newComm,
663  {
664  typename C::CopySet copyflags;
665  typename C::OwnerSet ownerflags;
666  typedef typename C::ParallelIndexSet IndexSet;
667  typedef RedistributeInformation<C> RI;
668  std::vector<typename M::size_type> rowsize(newComm.indexSet().size(), 0);
669  std::vector<typename M::size_type> copyrowsize(newComm.indexSet().size(), 0);
670  std::vector<typename M::size_type> backwardscopyrowsize(origComm.indexSet().size(), 0);
671 
672  // get owner rowsizes
673  CommMatrixRowSize<M,RI> commRowSize(origMatrix, ri);
674  ri.template redistribute<MatrixRowSizeGatherScatter<M,IndexSet,RI> >(commRowSize,commRowSize);
675 
676  origComm.buildGlobalLookup();
677 
678  for (std::size_t i=0; i < newComm.indexSet().size(); i++){
679  rowsize[i] = ri.getRowSize(i);
680  }
681  // get sparsity pattern from owner rows
683  origsp(origMatrix, origComm.globalLookup(), newComm.indexSet());
685  newsp(origMatrix, origComm.globalLookup(), newComm.indexSet(), rowsize);
686 
687  ri.template redistribute<MatrixSparsityPatternGatherScatter<M,IndexSet> >(origsp,newsp);
688 
689  // build copy to owner interface to get missing matrix values for novlp case
690  if (origComm.getSolverCategory() == SolverCategory::nonoverlapping){
691  RemoteIndices<IndexSet> *ris = new RemoteIndices<IndexSet>(origComm.indexSet(),
692  newComm.indexSet(),
693  origComm.communicator());
694  ris->template rebuild<true>();
695 
696  ri.getInterface().free();
697  ri.getInterface().build(*ris,copyflags,ownerflags);
698 
699  // get copy rowsizes
700  CommMatrixRowSize<M,RI> commRowSize_copy(origMatrix, ri);
701  ri.template redistribute<MatrixCopyRowSizeGatherScatter<M,IndexSet,RI> >(commRowSize_copy,
702  commRowSize_copy);
703 
704  for (std::size_t i=0; i < newComm.indexSet().size(); i++){
705  copyrowsize[i] = ri.getCopyRowSize(i);
706  }
707  //get copy rowsizes for sender
708  ri.redistributeBackward(backwardscopyrowsize,copyrowsize);
709  for (std::size_t i=0; i < origComm.indexSet().size(); i++){
710  ri.getBackwardsCopyRowSize(i) = backwardscopyrowsize[i];
711  }
712 
713  // get sparsity pattern from copy rows
714  CommMatrixSparsityPattern<M,IndexSet> origsp_copy(origMatrix,
715  origComm.globalLookup(),
716  newComm.indexSet(),
717  backwardscopyrowsize);
718  CommMatrixSparsityPattern<M,IndexSet> newsp_copy(origMatrix, origComm.globalLookup(),
719  newComm.indexSet(), copyrowsize);
720  ri.template redistribute<MatrixSparsityPatternGatherScatter<M,IndexSet> >(origsp_copy,
721  newsp_copy);
722 
723  newsp.completeSparsityPattern(newsp_copy.sparsity);
724  newsp.storeSparsityPattern(newMatrix);
725  }
726  else
727  newsp.storeSparsityPattern(newMatrix);
728 
729 #ifdef DUNE_ISTL_WITH_CHECKING
730  // Check for symmetry
731  int ret=0;
732  typedef typename M::ConstRowIterator RIter;
733  for(RIter row=newMatrix.begin(), rend=newMatrix.end(); row != rend; ++row){
734  typedef typename M::ConstColIterator CIter;
735  for(CIter col=row->begin(), cend=row->end(); col!=cend; ++col)
736  {
737  try{
738  newMatrix[col.index()][row.index()];
739  }catch(Dune::ISTLError e){
740  std::cerr<<newComm.communicator().rank()<<": entry ("
741  <<col.index()<<","<<row.index()<<") missing! for symmetry!"<<std::endl;
742  ret=1;
743 
744  }
745 
746  }
747  }
748 
749  if(ret)
750  DUNE_THROW(ISTLError, "Matrix not symmetric!");
751 #endif
752  }
753 
754  template<typename M, typename C>
755  void redistributeMatrixEntries(M& origMatrix, M& newMatrix, C& origComm, C& newComm,
757  {
758  typedef typename C::ParallelIndexSet IndexSet;
759  typename C::OwnerSet ownerflags;
760  std::vector<typename M::size_type> rowsize(newComm.indexSet().size(), 0);
761  std::vector<typename M::size_type> copyrowsize(newComm.indexSet().size(), 0);
762  std::vector<typename M::size_type> backwardscopyrowsize(origComm.indexSet().size(), 0);
763 
764  for (std::size_t i=0; i < newComm.indexSet().size(); i++){
765  rowsize[i] = ri.getRowSize(i);
766  if (origComm.getSolverCategory() == SolverCategory::nonoverlapping){
767  copyrowsize[i] = ri.getCopyRowSize(i);
768  }
769  }
770 
771  for (std::size_t i=0; i < origComm.indexSet().size(); i++)
772  if (origComm.getSolverCategory() == SolverCategory::nonoverlapping)
773  backwardscopyrowsize[i] = ri.getBackwardsCopyRowSize(i);
774 
775 
776  if (origComm.getSolverCategory() == SolverCategory::nonoverlapping){
777  // fill sparsity pattern from copy rows
778  CommMatrixRow<M,IndexSet> origrow_copy(origMatrix, origComm.globalLookup(),
779  newComm.indexSet(), backwardscopyrowsize);
780  CommMatrixRow<M,IndexSet> newrow_copy(newMatrix, origComm.globalLookup(),
781  newComm.indexSet(),copyrowsize);
782  ri.template redistribute<MatrixRowGatherScatter<M,IndexSet> >(origrow_copy,
783  newrow_copy);
784  ri.getInterface().free();
785  RemoteIndices<IndexSet> *ris = new RemoteIndices<IndexSet>(origComm.indexSet(),
786  newComm.indexSet(),
787  origComm.communicator());
788  ris->template rebuild<true>();
789  ri.getInterface().build(*ris,ownerflags,ownerflags);
790  }
791 
793  origrow(origMatrix, origComm.globalLookup(), newComm.indexSet());
795  newrow(newMatrix, origComm.globalLookup(), newComm.indexSet(),rowsize);
796  ri.template redistribute<MatrixRowGatherScatter<M,IndexSet> >(origrow,newrow);
797  if (origComm.getSolverCategory() != static_cast<int>(SolverCategory::nonoverlapping))
798  newrow.setOverlapRowsToDirichlet();
799  if(newMatrix.N()>0&&newMatrix.N()<20)
800  printmatrix(std::cout, newMatrix, "redist", "row");
801  }
802 
819  template<typename M, typename C>
820  void redistributeMatrix(M& origMatrix, M& newMatrix, C& origComm, C& newComm,
822  {
823  ri.setNoRows(newComm.indexSet().size());
824  ri.setNoCopyRows(newComm.indexSet().size());
825  ri.setNoBackwardsCopyRows(origComm.indexSet().size());
826  redistributeSparsityPattern(origMatrix, newMatrix, origComm, newComm, ri);
827  redistributeMatrixEntries(origMatrix, newMatrix, origComm, newComm, ri);
828  }
829 #endif
830 }
831 #endif