dune-istl  2.2.0
overlappingschwarz.hh
Go to the documentation of this file.
1 // -*- tab-width: 8; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2 // vi: set et ts=8 sw=2 sts=2:
3 #ifndef DUNE_OVERLAPPINGSCHWARZ_HH
4 #define DUNE_OVERLAPPINGSCHWARZ_HH
5 #include<cassert>
6 #include<algorithm>
7 #include<functional>
8 #include<vector>
9 #include<set>
10 #include<dune/common/dynmatrix.hh>
11 #include<dune/common/sllist.hh>
12 #include"preconditioners.hh"
13 #include"superlu.hh"
14 #include"bvector.hh"
15 #include"bcrsmatrix.hh"
16 #include"ilusubdomainsolver.hh"
17 
18 namespace Dune
19 {
20 
32  template<class M, class X, class TM, class TD, class TA>
33  class SeqOverlappingSchwarz;
34 
38  template<class I, class S, class D>
40  {
41  public:
43  typedef D subdomain_vector;
44 
45  typedef I InitializerList;
46  typedef typename InitializerList::value_type AtomInitializer;
47  typedef typename AtomInitializer::Matrix Matrix;
48  typedef typename Matrix::const_iterator Iter;
49  typedef typename Matrix::row_type::const_iterator CIter;
50 
51  typedef S IndexSet;
52  typedef typename IndexSet::size_type size_type;
53 
55  const IndexSet& indices,
56  const subdomain_vector& domains);
57 
58 
59  void addRowNnz(const Iter& row);
60 
61  void allocate();
62 
63  void countEntries(const Iter& row, const CIter& col) const;
64 
65  void calcColstart() const;
66 
67  void copyValue(const Iter& row, const CIter& col) const;
68 
69  void createMatrix() const;
70  private:
71  class IndexMap
72  {
73  public:
74  typedef typename S::size_type size_type;
75  typedef std::map<size_type,size_type> Map;
76  typedef typename Map::iterator iterator;
77  typedef typename Map::const_iterator const_iterator;
78 
79  IndexMap();
80 
81  void insert(size_type grow);
82 
83  const_iterator find(size_type grow) const;
84 
85  iterator find(size_type grow);
86 
87  iterator begin();
88 
89  const_iterator begin()const;
90 
91  iterator end();
92 
93  const_iterator end() const;
94 
95  private:
96  std::map<size_type,size_type> map_;
97  size_type row;
98  };
99 
100 
101  typedef typename InitializerList::iterator InitIterator;
102  typedef typename IndexSet::const_iterator IndexIteratur;
103  InitializerList* initializers;
104  const IndexSet *indices;
105  mutable std::vector<IndexMap> indexMaps;
106  const subdomain_vector& domains;
107  };
108 
113  {};
114 
119  {};
120 
126  {};
127 
132  template<class M, class X, class Y>
133  class DynamicMatrixSubdomainSolver;
134 
135  // Specialization for BCRSMatrix
136  template<class K, int n, class Al, class X, class Y>
137  class DynamicMatrixSubdomainSolver< BCRSMatrix< FieldMatrix<K,n,n>, Al>, X, Y >
138  {
139  typedef BCRSMatrix< FieldMatrix<K,n,n>, Al> M;
140  public:
143  typedef K field_type;
146  typedef X domain_type;
148  typedef Y range_type;
149 
154  void apply (DynamicVector<field_type>& v, DynamicVector<field_type>& d)
155  {
156  assert(v.size() > 0);
157  assert(v.size() == d.size());
158  assert(A.rows() <= v.size());
159  assert(A.cols() <= v.size());
160  size_t sz = A.rows();
161  v.resize(sz);
162  d.resize(sz);
163  A.solve(v,d);
164  v.resize(v.capacity());
165  d.resize(d.capacity());
166  }
167 
175  template<class S>
176  void setSubMatrix(const M& BCRS, S& rowset)
177  {
178  size_t sz = rowset.size();
179  A.resize(sz*n,sz*n);
180  typename DynamicMatrix<K>::RowIterator rIt = A.begin();
181  typedef typename S::const_iterator SIter;
182  size_t r = 0, c = 0;
183  for(SIter rowIdx = rowset.begin(), rowEnd=rowset.end();
184  rowIdx!= rowEnd; ++rowIdx, r++)
185  {
186  c = 0;
187  for(SIter colIdx = rowset.begin(), colEnd=rowset.end();
188  colIdx!= colEnd; ++colIdx, c++)
189  {
190  if (BCRS[*rowIdx].find(*colIdx) == BCRS[*rowIdx].end())
191  continue;
192  for (size_t i=0; i<n; i++)
193  {
194  for (size_t j=0; j<n; j++)
195  {
196  A[r*n+i][c*n+j] = BCRS[*rowIdx][*colIdx][i][j];
197  }
198  }
199  }
200  }
201  }
202  private:
203  DynamicMatrix<K> A;
204  };
205 
206  template<typename T>
208  {
209  };
210 
211  // specialization for DynamicMatrix
212  template<class K, int n, class Al, class X, class Y>
213  class OverlappingAssigner< DynamicMatrixSubdomainSolver< BCRSMatrix< FieldMatrix<K,n,n>, Al>, X, Y > >
214  {
215  public:
217  typedef K field_type;
218  typedef Y range_type;
219  typedef typename range_type::block_type block_type;
221 
229  OverlappingAssigner(std::size_t maxlength, const BCRSMatrix<FieldMatrix<K,n,n>, Al>& mat_, const X& b_, Y& x_);
230 
234  inline
235  void deallocate();
236 
240  inline
241  void resetIndexForNextDomain();
242 
247  inline
248  DynamicVector<K> & lhs();
249 
254  inline
255  DynamicVector<K> & rhs();
256 
261  inline
262  void relaxResult(field_type relax);
263 
268  void operator()(const size_type& domainIndex);
269 
277  inline
278  void assignResult(block_type& res);
279 
280  private:
284  const matrix_type* mat;
286  // we need a pointer, because we have to avoid deep copies
287  DynamicVector<field_type> * rhs_;
289  // we need a pointer, because we have to avoid deep copies
290  DynamicVector<field_type> * lhs_;
292  const range_type* b;
294  range_type* x;
296  std::size_t i;
298  std::size_t maxlength_;
299  };
300 
301 #if HAVE_SUPERLU
302  template<typename T, typename A, int n, int m>
303  struct OverlappingAssigner<SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> > >
304  {
307  typedef typename range_type::field_type field_type;
308  typedef typename range_type::block_type block_type;
309 
311 
319  OverlappingAssigner(std::size_t maxlength, const BCRSMatrix<FieldMatrix<T,n,m>,A>& mat,
320  const range_type& b, range_type& x);
326  void deallocate();
327 
328  /*
329  * @brief Resets the local index to zero.
330  */
331  void resetIndexForNextDomain();
332 
337  field_type* lhs();
338 
343  field_type* rhs();
344 
349  void relaxResult(field_type relax);
350 
355  void operator()(const size_type& domain);
356 
364  void assignResult(block_type& res);
365 
366  private:
370  const matrix_type* mat;
372  field_type* rhs_;
374  field_type* lhs_;
376  const range_type* b;
378  range_type* x;
380  std::size_t i;
382  std::size_t maxlength_;
383  };
384 
385 #endif
386 
387  template<class M, class X, class Y>
389  {
390  public:
391  typedef M matrix_type;
392 
393  typedef typename M::field_type field_type;
394 
395  typedef typename Y::block_type block_type;
396 
397  typedef typename matrix_type::size_type size_type;
405  OverlappingAssignerILUBase(std::size_t maxlength, const M& mat,
406  const Y& b, X& x);
412  void deallocate();
413 
418 
423  X& lhs();
424 
429  Y& rhs();
430 
435  void relaxResult(field_type relax);
436 
441  void operator()(const size_type& domain);
442 
450  void assignResult(block_type& res);
451 
452  private:
456  const M* mat;
458  X* lhs_;
460  Y* rhs_;
462  const Y* b;
464  X* x;
466  size_type i;
467  };
468 
469  // specialization for ILU0
470  template<class M, class X, class Y>
472  : public OverlappingAssignerILUBase<M,X,Y>
473  {
474  public:
482  OverlappingAssigner(std::size_t maxlength, const M& mat,
483  const Y& b, X& x)
484  : OverlappingAssignerILUBase<M,X,Y>(maxlength, mat,b,x)
485  {}
486  };
487 
488  // specialization for ILUN
489  template<class M, class X, class Y>
491  : public OverlappingAssignerILUBase<M,X,Y>
492  {
493  public:
501  OverlappingAssigner(std::size_t maxlength, const M& mat,
502  const Y& b, X& x)
503  : OverlappingAssignerILUBase<M,X,Y>(maxlength, mat,b,x)
504  {}
505  };
506 
507  template<typename S, typename T>
509  {
510  };
511 
512  template<typename S, typename T, typename A, int n>
513  struct AdditiveAdder<S, BlockVector<FieldVector<T,n>,A> >
514  {
515  typedef typename A::size_type size_type;
516  AdditiveAdder(BlockVector<FieldVector<T,n>,A>& v, BlockVector<FieldVector<T,n>,A>& x,
517  OverlappingAssigner<S>& assigner, const T& relax_);
518  void operator()(const size_type& domain);
519  void axpy();
520 
521  private:
524  OverlappingAssigner<S>* assigner;
525  T relax;
526  };
527 
528  template<typename S,typename T>
530  {
531  };
532 
533  template<typename S, typename T, typename A, int n>
534  struct MultiplicativeAdder<S, BlockVector<FieldVector<T,n>,A> >
535  {
536  typedef typename A::size_type size_type;
537  MultiplicativeAdder(BlockVector<FieldVector<T,n>,A>& v, BlockVector<FieldVector<T,n>,A>& x,
538  OverlappingAssigner<S>& assigner_, const T& relax_);
539  void operator()(const size_type& domain);
540  void axpy();
541 
542  private:
544  OverlappingAssigner<S>* assigner;
545  T relax;
546  };
547 
557  template<typename T, class X, class S>
559  {};
560 
561  template<class X, class S>
563  {
565  };
566 
567  template<class X, class S>
569  {
571  };
572 
573  template<class X, class S>
575  {
577  };
578 
590  template<typename T1, typename T2, bool forward>
592  {
593  typedef T1 solver_vector;
594  typedef typename solver_vector::iterator solver_iterator;
595  typedef T2 subdomain_vector;
596  typedef typename subdomain_vector::const_iterator domain_iterator;
597 
599  {
600  return sv.begin();
601  }
602 
604  {
605  return sv.end();
606  }
608  {
609  return sv.begin();
610  }
611 
613  {
614  return sv.end();
615  }
616  };
617 
618  template<typename T1, typename T2>
619  struct IteratorDirectionSelector<T1,T2,false>
620  {
621  typedef T1 solver_vector;
622  typedef typename solver_vector::reverse_iterator solver_iterator;
623  typedef T2 subdomain_vector;
624  typedef typename subdomain_vector::const_reverse_iterator domain_iterator;
625 
627  {
628  return sv.rbegin();
629  }
630 
632  {
633  return sv.rend();
634  }
636  {
637  return sv.rbegin();
638  }
639 
641  {
642  return sv.rend();
643  }
644  };
645 
654  template<class T>
656  {
657  typedef T smoother;
658  typedef typename smoother::range_type range_type;
659 
660  static void apply(smoother& sm, range_type& v, const range_type& b)
661  {
662  sm.template apply<true>(v, b);
663  }
664  };
665 
666  template<class M, class X, class TD, class TA>
668  {
671 
672  static void apply(smoother& sm, range_type& v, const range_type& b)
673  {
674  sm.template apply<true>(v, b);
675  sm.template apply<false>(v, b);
676  }
677  };
678 
679  template<class T>
681  {};
682 
683 
684  template<class K, int n, class Al, class X, class Y>
685  struct SeqOverlappingSchwarzAssembler< DynamicMatrixSubdomainSolver< BCRSMatrix< FieldMatrix<K,n,n>, Al>, X, Y > >
686  {
688  template<class RowToDomain, class Solvers, class SubDomains>
689  static std::size_t assembleLocalProblems(const RowToDomain& rowToDomain, const matrix_type& mat,
690  Solvers& solvers, const SubDomains& domains,
691  bool onTheFly);
692  };
693 
694 #if HAVE_SUPERLU
695  template<class T>
697  {
698  typedef T matrix_type;
699  template<class RowToDomain, class Solvers, class SubDomains>
700  static std::size_t assembleLocalProblems(const RowToDomain& rowToDomain, const matrix_type& mat,
701  Solvers& solvers, const SubDomains& domains,
702  bool onTheFly);
703  };
704 #endif
705 
706  template<class M,class X, class Y>
708  {
709  typedef M matrix_type;
710  template<class RowToDomain, class Solvers, class SubDomains>
711  static std::size_t assembleLocalProblems(const RowToDomain& rowToDomain, const matrix_type& mat,
712  Solvers& solvers, const SubDomains& domains,
713  bool onTheFly);
714  };
715 
716  template<class M,class X, class Y>
719  {};
720 
721  template<class M,class X, class Y>
724  {};
725 
736  template<class M, class X, class TM=AdditiveSchwarzMode,
737  class TD=ILU0SubdomainSolver<M,X,X>, class TA=std::allocator<X> >
739  : public Preconditioner<X,X>
740  {
741  public:
745  typedef M matrix_type;
746 
750  typedef X domain_type;
751 
755  typedef X range_type;
756 
763  typedef TM Mode;
764 
768  typedef typename X::field_type field_type;
769 
771  typedef typename matrix_type::size_type size_type;
772 
774  typedef TA allocator;
775 
777  typedef std::set<size_type, std::less<size_type>,
778  typename TA::template rebind<size_type>::other>
780 
782  typedef std::vector<subdomain_type, typename TA::template rebind<subdomain_type>::other> subdomain_vector;
783 
785  typedef SLList<size_type, typename TA::template rebind<size_type>::other> subdomain_list;
786 
788  typedef std::vector<subdomain_list, typename TA::template rebind<subdomain_list>::other > rowtodomain_vector;
789 
791  typedef TD slu;
792 
794  typedef std::vector<slu, typename TA::template rebind<slu>::other> slu_vector;
795 
796  enum{
799  };
800 
814  SeqOverlappingSchwarz(const matrix_type& mat, const subdomain_vector& subDomains,
815  field_type relaxationFactor=1, bool onTheFly_=true);
816 
828  SeqOverlappingSchwarz(const matrix_type& mat, const rowtodomain_vector& rowToDomain,
829  field_type relaxationFactor=1, bool onTheFly_=true);
830 
836  virtual void pre (X& x, X& b) {}
837 
843  virtual void apply (X& v, const X& d);
844 
845  template<bool forward>
846  void apply(X& v, const X& d);
847 
853  virtual void post (X& x) {
854  Dune::dverb<<" avg nnz over subdomain is "<<nnz<<std::endl;
855  }
856 
857  private:
858  const M& mat;
859  slu_vector solvers;
860  subdomain_vector subDomains;
861  field_type relax;
862 
863  typename M::size_type maxlength;
864  std::size_t nnz;
865 
866  bool onTheFly;
867  };
868 
869 
870 
871  template<class I, class S, class D>
873  const IndexSet& idx,
874  const subdomain_vector& domains_)
875  : initializers(&il), indices(&idx), indexMaps(il.size()), domains(domains_)
876  {
877  }
878 
879 
880  template<class I, class S, class D>
882  {
883  typedef typename IndexSet::value_type::const_iterator iterator;
884  for(iterator domain=(*indices)[row.index()].begin(); domain != (*indices)[row.index()].end(); ++domain){
885  (*initializers)[*domain].addRowNnz(row, domains[*domain]);
886  indexMaps[*domain].insert(row.index());
887  }
888  }
889 
890  template<class I, class S, class D>
892  {
893  std::for_each(initializers->begin(), initializers->end(),
894  std::mem_fun_ref(&AtomInitializer::allocateMatrixStorage));
895  std::for_each(initializers->begin(), initializers->end(),
896  std::mem_fun_ref(&AtomInitializer::allocateMarker));
897  }
898 
899  template<class I, class S, class D>
901  {
902  typedef typename IndexSet::value_type::const_iterator iterator;
903  for(iterator domain=(*indices)[row.index()].begin(); domain != (*indices)[row.index()].end(); ++domain){
904  typename std::map<size_type,size_type>::const_iterator v = indexMaps[*domain].find(col.index());
905  if(v!= indexMaps[*domain].end()){
906  (*initializers)[*domain].countEntries(indexMaps[*domain].find(col.index())->second);
907  }
908  }
909  }
910 
911  template<class I, class S, class D>
913  {
914  std::for_each(initializers->begin(), initializers->end(),
915  std::mem_fun_ref(&AtomInitializer::calcColstart));
916  }
917 
918  template<class I, class S, class D>
920  {
921  typedef typename IndexSet::value_type::const_iterator iterator;
922  for(iterator domain=(*indices)[row.index()].begin(); domain!= (*indices)[row.index()].end(); ++domain){
923  typename std::map<size_type,size_type>::const_iterator v = indexMaps[*domain].find(col.index());
924  if(v!= indexMaps[*domain].end()){
925  assert(indexMaps[*domain].end()!=indexMaps[*domain].find(row.index()));
926  (*initializers)[*domain].copyValue(col, indexMaps[*domain].find(row.index())->second,
927  v->second);
928  }
929  }
930  }
931 
932  template<class I, class S, class D>
934  {
935  indexMaps.clear();
936  indexMaps.swap(std::vector<IndexMap>(indexMaps));
937  std::for_each(initializers->begin(), initializers->end(),
938  std::mem_fun_ref(&AtomInitializer::createMatrix));
939  }
940 
941  template<class I, class S, class D>
943  : row(0)
944  {}
945 
946  template<class I, class S, class D>
948  {
949  assert(map_.find(grow)==map_.end());
950  map_.insert(std::make_pair(grow, row++));
951  }
952 
953  template<class I, class S, class D>
954  typename OverlappingSchwarzInitializer<I,S,D>::IndexMap::const_iterator
956  {
957  return map_.find(grow);
958  }
959 
960  template<class I, class S, class D>
961  typename OverlappingSchwarzInitializer<I,S,D>::IndexMap::iterator
963  {
964  return map_.find(grow);
965  }
966 
967  template<class I, class S, class D>
968  typename OverlappingSchwarzInitializer<I,S,D>::IndexMap::const_iterator
970  {
971  return map_.end();
972  }
973 
974  template<class I, class S, class D>
975  typename OverlappingSchwarzInitializer<I,S,D>::IndexMap::iterator
977  {
978  return map_.end();
979  }
980 
981  template<class I, class S, class D>
982  typename OverlappingSchwarzInitializer<I,S,D>::IndexMap::const_iterator
984  {
985  return map_.begin();
986  }
987 
988  template<class I, class S, class D>
989  typename OverlappingSchwarzInitializer<I,S,D>::IndexMap::iterator
991  {
992  return map_.begin();
993  }
994 
995  template<class M, class X, class TM, class TD, class TA>
997  field_type relaxationFactor, bool fly)
998  : mat(mat_), relax(relaxationFactor), onTheFly(fly)
999  {
1000  typedef typename rowtodomain_vector::const_iterator RowDomainIterator;
1001  typedef typename subdomain_list::const_iterator DomainIterator;
1002 #ifdef DUNE_ISTL_WITH_CHECKING
1003  assert(rowToDomain.size()==mat.N());
1004  assert(rowToDomain.size()==mat.M());
1005 
1006  for(RowDomainIterator iter=rowToDomain.begin(); iter != rowToDomain.end(); ++iter)
1007  assert(iter->size()>0);
1008 
1009 #endif
1010  // calculate the number of domains
1011  size_type domains=0;
1012  for(RowDomainIterator iter=rowToDomain.begin(); iter != rowToDomain.end(); ++iter)
1013  for(DomainIterator d=iter->begin(); d != iter->end(); ++d)
1014  domains=std::max(domains, *d);
1015  ++domains;
1016 
1017  solvers.resize(domains);
1018  subDomains.resize(domains);
1019 
1020  // initialize subdomains to row mapping from row to subdomain mapping
1021  size_type row=0;
1022  for(RowDomainIterator iter=rowToDomain.begin(); iter != rowToDomain.end(); ++iter, ++row)
1023  for(DomainIterator d=iter->begin(); d != iter->end(); ++d)
1024  subDomains[*d].insert(row);
1025 
1026 #ifdef DUNE_ISTL_WITH_CHECKING
1027  size_type i=0;
1028  typedef typename subdomain_vector::const_iterator iterator;
1029  for(iterator iter=subDomains.begin(); iter != subDomains.end(); ++iter){
1030  typedef typename subdomain_type::const_iterator entry_iterator;
1031  Dune::dvverb<<"domain "<<i++<<":";
1032  for(entry_iterator entry = iter->begin(); entry != iter->end(); ++entry){
1033  Dune::dvverb<<" "<<*entry;
1034  }
1035  Dune::dvverb<<std::endl;
1036  }
1037 #endif
1039  ::assembleLocalProblems(rowToDomain, mat, solvers, subDomains, onTheFly);
1040  }
1041 
1042  template<class M, class X, class TM, class TD, class TA>
1044  const subdomain_vector& sd,
1045  field_type relaxationFactor,
1046  bool fly)
1047  : mat(mat_), solvers(sd.size()), subDomains(sd), relax(relaxationFactor),
1048  onTheFly(fly)
1049  {
1050  typedef typename subdomain_vector::const_iterator DomainIterator;
1051 
1052 #ifdef DUNE_ISTL_WITH_CHECKING
1053  size_type i=0;
1054 
1055  for(DomainIterator d=sd.begin(); d != sd.end(); ++d,++i){
1056  //std::cout<<i<<": "<<d->size()<<std::endl;
1057  assert(d->size()>0);
1058  typedef typename DomainIterator::value_type::const_iterator entry_iterator;
1059  Dune::dvverb<<"domain "<<i<<":";
1060  for(entry_iterator entry = d->begin(); entry != d->end(); ++entry){
1061  Dune::dvverb<<" "<<*entry;
1062  }
1063  Dune::dvverb<<std::endl;
1064  }
1065 
1066 #endif
1067 
1068  // Create a row to subdomain mapping
1069  rowtodomain_vector rowToDomain(mat.N());
1070 
1071  size_type domainId=0;
1072 
1073  for(DomainIterator domain=sd.begin(); domain != sd.end(); ++domain, ++domainId){
1074  typedef typename subdomain_type::const_iterator iterator;
1075  for(iterator row=domain->begin(); row != domain->end(); ++row)
1076  rowToDomain[*row].push_back(domainId);
1077  }
1078 
1080  ::assembleLocalProblems(rowToDomain, mat, solvers, subDomains, onTheFly);
1081  }
1082 
1089  template<class M>
1091 
1092  template<typename T, typename A, int n, int m>
1093  struct SeqOverlappingSchwarzDomainSize<BCRSMatrix<FieldMatrix<T,n,m>,A > >
1094  {
1095  template<class Domain>
1096  static int size(const Domain & d)
1097  {
1098  assert(n==m);
1099  return m*d.size();
1100  }
1101  };
1102 
1103  template<class K, int n, class Al, class X, class Y>
1104  template<class RowToDomain, class Solvers, class SubDomains>
1105  std::size_t
1106  SeqOverlappingSchwarzAssembler< DynamicMatrixSubdomainSolver< BCRSMatrix< FieldMatrix<K,n,n>, Al>, X, Y > >::
1107  assembleLocalProblems(const RowToDomain& rowToDomain,
1108  const matrix_type& mat,
1109  Solvers& solvers,
1110  const SubDomains& subDomains,
1111  bool onTheFly)
1112  {
1113  typedef typename SubDomains::const_iterator DomainIterator;
1114  std::size_t maxlength = 0;
1115 
1116  assert(onTheFly);
1117 
1118  for(DomainIterator domain=subDomains.begin();domain!=subDomains.end();++domain)
1119  maxlength=std::max(maxlength, domain->size());
1120  maxlength*=n;
1121 
1122  return maxlength;
1123  }
1124 
1125 #if HAVE_SUPERLU
1126  template<class T>
1127  template<class RowToDomain, class Solvers, class SubDomains>
1128  std::size_t SeqOverlappingSchwarzAssembler<SuperLU<T> >::assembleLocalProblems(const RowToDomain& rowToDomain,
1129  const matrix_type& mat,
1130  Solvers& solvers,
1131  const SubDomains& subDomains,
1132  bool onTheFly)
1133  {
1134  typedef typename std::vector<SuperMatrixInitializer<matrix_type> >::iterator InitializerIterator;
1135  typedef typename SubDomains::const_iterator DomainIterator;
1136  typedef typename Solvers::iterator SolverIterator;
1137  std::size_t maxlength = 0;
1138 
1139  if(onTheFly){
1140  for(DomainIterator domain=subDomains.begin();domain!=subDomains.end();++domain)
1141  maxlength=std::max(maxlength, domain->size());
1142  maxlength*=mat[0].begin()->N();
1143  }else{
1144  // initialize the initializers
1145  DomainIterator domain=subDomains.begin();
1146 
1147  // Create the initializers list.
1148  std::vector<SuperMatrixInitializer<matrix_type> > initializers(subDomains.size());
1149 
1150  SolverIterator solver=solvers.begin();
1151  for(InitializerIterator initializer=initializers.begin(); initializer!=initializers.end();
1152  ++initializer, ++solver, ++domain){
1155  //solver->setVerbosity(true);
1156  *initializer=SuperMatrixInitializer<matrix_type>(solver->mat);
1157  }
1158 
1159  // Set up the supermatrices according to the subdomains
1161  RowToDomain, SubDomains> Initializer;
1162 
1163  Initializer initializer(initializers, rowToDomain, subDomains);
1164  copyToSuperMatrix(initializer, mat);
1165  if(solvers.size()==1)
1166  assert(solvers[0].mat==mat);
1167 
1168  /* for(SolverIterator solver=solvers.begin(); solver!=solvers.end(); ++solver)
1169  dPrint_CompCol_Matrix("superlu", &static_cast<SuperMatrix&>(solver->mat)); */
1170 
1171  // Calculate the LU decompositions
1172  std::for_each(solvers.begin(), solvers.end(), std::mem_fun_ref(&SuperLU<matrix_type>::decompose));
1173  for(SolverIterator solver=solvers.begin(); solver!=solvers.end(); ++solver){
1174  assert(solver->mat.N()==solver->mat.M());
1175  maxlength=std::max(maxlength, solver->mat.N());
1176  //writeCompColMatrixToMatlab(static_cast<SuperLUMatrix<M>&>(solver->mat), std::cout);
1177  }
1178  }
1179  return maxlength;
1180  }
1181 
1182 #endif
1183 
1184  template<class M,class X,class Y>
1185  template<class RowToDomain, class Solvers, class SubDomains>
1186  std::size_t SeqOverlappingSchwarzAssemblerILUBase<M,X,Y>::assembleLocalProblems(const RowToDomain& rowToDomain,
1187  const matrix_type& mat,
1188  Solvers& solvers,
1189  const SubDomains& subDomains,
1190  bool onTheFly)
1191  {
1192  typedef typename SubDomains::const_iterator DomainIterator;
1193  typedef typename Solvers::iterator SolverIterator;
1194  std::size_t maxlength = 0;
1195 
1196  if(onTheFly){
1197  for(DomainIterator domain=subDomains.begin();domain!=subDomains.end();++domain)
1198  maxlength=std::max(maxlength, domain->size());
1199  }else{
1200  // initialize the solvers of the local prolems.
1201  SolverIterator solver=solvers.begin();
1202  for(DomainIterator domain=subDomains.begin(); domain!=subDomains.end();
1203  ++domain, ++solver){
1204  solver->setSubMatrix(mat, *domain);
1205  maxlength=std::max(maxlength, domain->size());
1206  }
1207  }
1208 
1209  return maxlength;
1210 
1211  }
1212 
1213 
1214  template<class M, class X, class TM, class TD, class TA>
1216  {
1218  }
1219 
1220  template<class M, class X, class TM, class TD, class TA>
1221  template<bool forward>
1222  void SeqOverlappingSchwarz<M,X,TM,TD,TA>::apply(X& x, const X& b)
1223  {
1224  typedef typename X::block_type block;
1225  typedef slu_vector solver_vector;
1228  domain_iterator;
1229 
1230  OverlappingAssigner<TD> assigner(maxlength, mat, b, x);
1231 
1234  X v(x); // temporary for the update
1235  v=0;
1236 
1237  typedef typename AdderSelector<TM,X,TD >::Adder Adder;
1238  Adder adder(v, x, assigner, relax);
1239 
1240  nnz=0;
1241  std::size_t no=0;
1242  for(;domain != IteratorDirectionSelector<solver_vector,subdomain_vector,forward>::end(subDomains); ++domain){
1243  //Copy rhs to C-array for SuperLU
1244  std::for_each(domain->begin(), domain->end(), assigner);
1245  assigner.resetIndexForNextDomain();
1246  if(onTheFly){
1247  // Create the subdomain solver
1248  slu sdsolver;
1249  sdsolver.setSubMatrix(mat, *domain);
1250  // Apply
1251  sdsolver.apply(assigner.lhs(), assigner.rhs());
1252  //nnz+=sdsolver.nnz();
1253  }else{
1254  solver->apply(assigner.lhs(), assigner.rhs());
1255  //nnz+=solver->nnz();
1256  ++solver;
1257  }
1258  ++no;
1259  //Add relaxed correction to from SuperLU to v
1260  std::for_each(domain->begin(), domain->end(), adder);
1261  assigner.resetIndexForNextDomain();
1262 
1263  }
1264  nnz/=no;
1265 
1266  adder.axpy();
1267  assigner.deallocate();
1268  }
1269 
1270  template<class K, int n, class Al, class X, class Y>
1271  OverlappingAssigner< DynamicMatrixSubdomainSolver< BCRSMatrix< FieldMatrix<K,n,n>, Al>, X, Y > >
1272  ::OverlappingAssigner(std::size_t maxlength, const BCRSMatrix<FieldMatrix<K,n,n>, Al>& mat_,
1273  const X& b_, Y& x_) :
1274  mat(&mat_),
1275  rhs_( new DynamicVector<field_type>(maxlength, 42) ),
1276  lhs_( new DynamicVector<field_type>(maxlength, -42) ),
1277  b(&b_),
1278  x(&x_),
1279  i(0),
1280  maxlength_(maxlength)
1281  {}
1282 
1283  template<class K, int n, class Al, class X, class Y>
1284  void
1286  ::deallocate()
1287  {
1288  delete rhs_;
1289  delete lhs_;
1290  }
1291 
1292  template<class K, int n, class Al, class X, class Y>
1293  void
1295  ::resetIndexForNextDomain()
1296  {
1297  i=0;
1298  }
1299 
1300  template<class K, int n, class Al, class X, class Y>
1301  DynamicVector<K> &
1303  ::lhs()
1304  {
1305  return *lhs_;
1306  }
1307 
1308  template<class K, int n, class Al, class X, class Y>
1309  DynamicVector<K> &
1311  ::rhs()
1312  {
1313  return *rhs_;
1314  }
1315 
1316  template<class K, int n, class Al, class X, class Y>
1317  void
1319  ::relaxResult(field_type relax)
1320  {
1321  lhs() *= relax;
1322  }
1323 
1324  template<class K, int n, class Al, class X, class Y>
1325  void
1327  ::operator()(const size_type& domainIndex)
1328  {
1329  lhs() = 0.0;
1330 #if 0
1331  //assign right hand side of current domainindex block
1332  for(size_type j=0; j<n; ++j, ++i){
1333  assert(i<maxlength_);
1334  rhs()[i]=(*b)[domainIndex][j];
1335  }
1336 
1337  // loop over all Matrix row entries and calculate defect.
1338  typedef typename matrix_type::ConstColIterator col_iterator;
1339 
1340  // calculate defect for current row index block
1341  for(col_iterator col=(*mat)[domainIndex].begin(); col!=(*mat)[domainIndex].end(); ++col){
1342  block_type tmp(0.0);
1343  (*col).mv((*x)[col.index()], tmp);
1344  i-=n;
1345  for(size_type j=0; j<n; ++j, ++i){
1346  assert(i<maxlength_);
1347  rhs()[i]-=tmp[j];
1348  }
1349  }
1350 #else
1351  //assign right hand side of current domainindex block
1352  for(size_type j=0; j<n; ++j, ++i){
1353  assert(i<maxlength_);
1354  rhs()[i]=(*b)[domainIndex][j];
1355 
1356  // loop over all Matrix row entries and calculate defect.
1357  typedef typename matrix_type::ConstColIterator col_iterator;
1358 
1359  // calculate defect for current row index block
1360  for(col_iterator col=(*mat)[domainIndex].begin(); col!=(*mat)[domainIndex].end(); ++col){
1361  for(size_type k=0; k<n; ++k){
1362  rhs()[i]-=(*col)[j][k] * (*x)[col.index()][k];
1363  }
1364  }
1365  }
1366 #endif
1367  }
1368 
1369  template<class K, int n, class Al, class X, class Y>
1370  void
1372  ::assignResult(block_type& res)
1373  {
1374  // assign the result of the local solve to the global vector
1375  for(size_type j=0; j<n; ++j, ++i){
1376  assert(i<maxlength_);
1377  res[j]+=lhs()[i];
1378  }
1379  }
1380 
1381 #if HAVE_SUPERLU
1382 
1383  template<typename T, typename A, int n, int m>
1385  ::OverlappingAssigner(std::size_t maxlength,
1386  const BCRSMatrix<FieldMatrix<T,n,m>,A>& mat_,
1387  const range_type& b_,
1388  range_type& x_)
1389  : mat(&mat_),
1390  b(&b_),
1391  x(&x_), i(0), maxlength_(maxlength)
1392  {
1393  rhs_ = new field_type[maxlength];
1394  lhs_ = new field_type[maxlength];
1395 
1396  }
1397 
1398  template<typename T, typename A, int n, int m>
1400  {
1401  delete[] rhs_;
1402  delete[] lhs_;
1403  }
1404 
1405  template<typename T, typename A, int n, int m>
1407  {
1408  //assign right hand side of current domainindex block
1409  // rhs is an array of doubles!
1410  // rhs[starti] = b[domainindex]
1411  for(size_type j=0; j<n; ++j, ++i){
1412  assert(i<maxlength_);
1413  rhs_[i]=(*b)[domainIndex][j];
1414  }
1415 
1416 
1417  // loop over all Matrix row entries and calculate defect.
1418  typedef typename matrix_type::ConstColIterator col_iterator;
1419 
1420  // calculate defect for current row index block
1421  for(col_iterator col=(*mat)[domainIndex].begin(); col!=(*mat)[domainIndex].end(); ++col){
1422  block_type tmp;
1423  (*col).mv((*x)[col.index()], tmp);
1424  i-=n;
1425  for(size_type j=0; j<n; ++j, ++i){
1426  assert(i<maxlength_);
1427  rhs_[i]-=tmp[j];
1428  }
1429 
1430  }
1431 
1432  }
1433  template<typename T, typename A, int n, int m>
1435  {
1436  for(size_type j=i+n; i<j; ++i){
1437  assert(i<maxlength_);
1438  lhs_[i]*=relax;
1439  }
1440  i-=n;
1441  }
1442 
1443  template<typename T, typename A, int n, int m>
1445  {
1446  // assign the result of the local solve to the global vector
1447  for(size_type j=0; j<n; ++j, ++i){
1448  assert(i<maxlength_);
1449  res[j]+=lhs_[i];
1450  }
1451  }
1452 
1453  template<typename T, typename A, int n, int m>
1455  {
1456  i=0;
1457  }
1458 
1459  template<typename T, typename A, int n, int m>
1462  {
1463  return lhs_;
1464  }
1465 
1466  template<typename T, typename A, int n, int m>
1469  {
1470  return rhs_;
1471  }
1472 
1473 #endif
1474 
1475  template<class M, class X, class Y>
1477  const M& mat_,
1478  const Y& b_,
1479  X& x_)
1480  : mat(&mat_),
1481  b(&b_),
1482  x(&x_), i(0)
1483  {
1484  rhs_= new Y(maxlength);
1485  lhs_ = new X(maxlength);
1486  }
1487 
1488  template<class M, class X, class Y>
1490  {
1491  delete rhs_;
1492  delete lhs_;
1493  }
1494 
1495  template<class M, class X, class Y>
1497  {
1498  (*rhs_)[i]=(*b)[domainIndex];
1499 
1500  // loop over all Matrix row entries and calculate defect.
1501  typedef typename matrix_type::ConstColIterator col_iterator;
1502 
1503  // calculate defect for current row index block
1504  for(col_iterator col=(*mat)[domainIndex].begin(); col!=(*mat)[domainIndex].end(); ++col){
1505  (*col).mmv((*x)[col.index()], (*rhs_)[i]);
1506  }
1507  // Goto next local index
1508  ++i;
1509  }
1510 
1511  template<class M, class X, class Y>
1513  {
1514  (*lhs_)[i]*=relax;
1515  }
1516 
1517  template<class M, class X, class Y>
1519  {
1520  res+=(*lhs_)[i++];
1521  }
1522 
1523  template<class M, class X, class Y>
1525  {
1526  return *lhs_;
1527  }
1528 
1529  template<class M, class X, class Y>
1531  {
1532  return *rhs_;
1533  }
1534 
1535  template<class M, class X, class Y>
1537  {
1538  i=0;
1539  }
1540 
1541  template<typename S, typename T, typename A, int n>
1543  BlockVector<FieldVector<T,n>,A>& x_,
1544  OverlappingAssigner<S>& assigner_,
1545  const T& relax_)
1546  : v(&v_), x(&x_), assigner(&assigner_), relax(relax_)
1547  {}
1548 
1549  template<typename S, typename T, typename A, int n>
1550  void AdditiveAdder<S,BlockVector<FieldVector<T,n>,A> >::operator()(const size_type& domainIndex)
1551  {
1552  // add the result of the local solve to the current update
1553  assigner->assignResult((*v)[domainIndex]);
1554  }
1555 
1556 
1557  template<typename S, typename T, typename A, int n>
1559  {
1560  // relax the update and add it to the current guess.
1561  x->axpy(relax,*v);
1562  }
1563 
1564 
1565  template<typename S, typename T, typename A, int n>
1567  ::MultiplicativeAdder(BlockVector<FieldVector<T,n>,A>& v_,
1568  BlockVector<FieldVector<T,n>,A>& x_,
1569  OverlappingAssigner<S>& assigner_, const T& relax_)
1570  : x(&x_), assigner(&assigner_), relax(relax_)
1571  {}
1572 
1573 
1574  template<typename S,typename T, typename A, int n>
1575  void MultiplicativeAdder<S,BlockVector<FieldVector<T,n>,A> >::operator()(const size_type& domainIndex)
1576  {
1577  // add the result of the local solve to the current guess
1578  assigner->relaxResult(relax);
1579  assigner->assignResult((*x)[domainIndex]);
1580  }
1581 
1582 
1583  template<typename S,typename T, typename A, int n>
1585  {
1586  // nothing to do, as the corrections already relaxed and added in operator()
1587  }
1588 
1589 
1591 }
1592 
1593 #endif