Dune Core Modules (2.6.0)

overlappingschwarz.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_OVERLAPPINGSCHWARZ_HH
4 #define DUNE_ISTL_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 <dune/common/unused.hh>
13 #include "preconditioners.hh"
14 #include "superlu.hh"
15 #include "umfpack.hh"
16 #include "bvector.hh"
17 #include "bcrsmatrix.hh"
18 #include "ilusubdomainsolver.hh"
19 #include <dune/istl/solvertype.hh>
20 
21 namespace Dune
22 {
23 
35  template<class M, class X, class TM, class TD, class TA>
36  class SeqOverlappingSchwarz;
37 
41  template<class I, class S, class D>
43  {
44  public:
46  typedef D subdomain_vector;
47 
48  typedef I InitializerList;
49  typedef typename InitializerList::value_type AtomInitializer;
50  typedef typename AtomInitializer::Matrix Matrix;
51  typedef typename Matrix::const_iterator Iter;
52  typedef typename Matrix::row_type::const_iterator CIter;
53 
54  typedef S IndexSet;
55  typedef typename IndexSet::size_type size_type;
56 
57  OverlappingSchwarzInitializer(InitializerList& il,
58  const IndexSet& indices,
59  const subdomain_vector& domains);
60 
61 
62  void addRowNnz(const Iter& row);
63 
64  void allocate();
65 
66  void countEntries(const Iter& row, const CIter& col) const;
67 
68  void calcColstart() const;
69 
70  void copyValue(const Iter& row, const CIter& col) const;
71 
72  void createMatrix() const;
73  private:
74  class IndexMap
75  {
76  public:
77  typedef typename S::size_type size_type;
78  typedef std::map<size_type,size_type> Map;
79  typedef typename Map::iterator iterator;
80  typedef typename Map::const_iterator const_iterator;
81 
82  IndexMap();
83 
84  void insert(size_type grow);
85 
86  const_iterator find(size_type grow) const;
87 
88  iterator find(size_type grow);
89 
90  iterator begin();
91 
92  const_iterator begin() const;
93 
94  iterator end();
95 
96  const_iterator end() const;
97 
98  private:
99  std::map<size_type,size_type> map_;
100  size_type row;
101  };
102 
103 
104  typedef typename InitializerList::iterator InitIterator;
105  typedef typename IndexSet::const_iterator IndexIteratur;
106  InitializerList* initializers;
107  const IndexSet *indices;
108  mutable std::vector<IndexMap> indexMaps;
109  const subdomain_vector& domains;
110  };
111 
116  {};
117 
122  {};
123 
129  {};
130 
135  template<class M, class X, class Y>
137 
138  // Specialization for BCRSMatrix
139  template<class K, int n, class Al, class X, class Y>
140  class DynamicMatrixSubdomainSolver< BCRSMatrix< FieldMatrix<K,n,n>, Al>, X, Y >
141  {
142  typedef BCRSMatrix< FieldMatrix<K,n,n>, Al> M;
143  public:
145  typedef typename std::remove_const<M>::type matrix_type;
146  typedef K field_type;
147  typedef typename std::remove_const<M>::type rilu_type;
149  typedef X domain_type;
151  typedef Y range_type;
152 
158  {
159  assert(v.size() > 0);
160  assert(v.size() == d.size());
161  assert(A.rows() <= v.size());
162  assert(A.cols() <= v.size());
163  size_t sz = A.rows();
164  v.resize(sz);
165  d.resize(sz);
166  A.solve(v,d);
167  v.resize(v.capacity());
168  d.resize(d.capacity());
169  }
170 
178  template<class S>
179  void setSubMatrix(const M& BCRS, S& rowset)
180  {
181  size_t sz = rowset.size();
182  A.resize(sz*n,sz*n);
183  typedef typename S::const_iterator SIter;
184  size_t r = 0, c = 0;
185  for(SIter rowIdx = rowset.begin(), rowEnd=rowset.end();
186  rowIdx!= rowEnd; ++rowIdx, r++)
187  {
188  c = 0;
189  for(SIter colIdx = rowset.begin(), colEnd=rowset.end();
190  colIdx!= colEnd; ++colIdx, c++)
191  {
192  if (BCRS[*rowIdx].find(*colIdx) == BCRS[*rowIdx].end())
193  continue;
194  for (size_t i=0; i<n; i++)
195  {
196  for (size_t j=0; j<n; j++)
197  {
198  A[r*n+i][c*n+j] = BCRS[*rowIdx][*colIdx][i][j];
199  }
200  }
201  }
202  }
203  }
204  private:
205  DynamicMatrix<K> A;
206  };
207 
208  template<typename T, bool tag>
209  class OverlappingAssignerHelper
210  {};
211 
212  template<typename T>
213  using OverlappingAssigner = OverlappingAssignerHelper<T, Dune::StoresColumnCompressed<T>::value>;
214 
215  // specialization for DynamicMatrix
216  template<class K, int n, class Al, class X, class Y>
217  class OverlappingAssignerHelper< DynamicMatrixSubdomainSolver< BCRSMatrix< FieldMatrix<K,n,n>, Al>, X, Y >,false>
218  {
219  public:
220  typedef BCRSMatrix< FieldMatrix<K,n,n>, Al> matrix_type;
221  typedef K field_type;
222  typedef Y range_type;
223  typedef typename range_type::block_type block_type;
224  typedef typename matrix_type::size_type size_type;
225 
233  OverlappingAssignerHelper(std::size_t maxlength, const BCRSMatrix<FieldMatrix<K,n,n>, Al>& mat_, const X& b_, Y& x_);
234 
238  inline
239  void deallocate();
240 
244  inline
245  void resetIndexForNextDomain();
246 
251  inline
252  DynamicVector<K> & lhs();
253 
258  inline
259  DynamicVector<K> & rhs();
260 
265  inline
266  void relaxResult(field_type relax);
267 
272  void operator()(const size_type& domainIndex);
273 
281  inline
282  void assignResult(block_type& res);
283 
284  private:
288  const matrix_type* mat;
290  // we need a pointer, because we have to avoid deep copies
291  DynamicVector<field_type> * rhs_;
293  // we need a pointer, because we have to avoid deep copies
294  DynamicVector<field_type> * lhs_;
296  const range_type* b;
298  range_type* x;
300  std::size_t i;
302  std::size_t maxlength_;
303  };
304 
305 #if HAVE_SUPERLU || HAVE_SUITESPARSE_UMFPACK
306  template<template<class> class S, int n, int m, typename T, typename A>
307  struct OverlappingAssignerHelper<S<BCRSMatrix<FieldMatrix<T,n,m>, A> >, true>
308  {
309  typedef BCRSMatrix<FieldMatrix<T,n,m>, A> matrix_type;
310  typedef typename S<BCRSMatrix<FieldMatrix<T,n,m>, A> >::range_type range_type;
311  typedef typename range_type::field_type field_type;
312  typedef typename range_type::block_type block_type;
313 
314  typedef typename matrix_type::size_type size_type;
315 
323  OverlappingAssignerHelper(std::size_t maxlength, const matrix_type& mat,
324  const range_type& b, range_type& x);
330  void deallocate();
331 
332  /*
333  * @brief Resets the local index to zero.
334  */
335  void resetIndexForNextDomain();
336 
341  field_type* lhs();
342 
347  field_type* rhs();
348 
353  void relaxResult(field_type relax);
354 
359  void operator()(const size_type& domain);
360 
368  void assignResult(block_type& res);
369 
370  private:
374  const matrix_type* mat;
376  field_type* rhs_;
378  field_type* lhs_;
380  const range_type* b;
382  range_type* x;
384  std::size_t i;
386  std::size_t maxlength_;
387  };
388 
389 #endif // HAVE_SUPERLU || HAVE_SUITESPARSE_UMFPACK
390 
391  template<class M, class X, class Y>
392  class OverlappingAssignerILUBase
393  {
394  public:
395  typedef M matrix_type;
396 
397  typedef typename M::field_type field_type;
398 
399  typedef typename Y::block_type block_type;
400 
401  typedef typename matrix_type::size_type size_type;
409  OverlappingAssignerILUBase(std::size_t maxlength, const M& mat,
410  const Y& b, X& x);
416  void deallocate();
417 
421  void resetIndexForNextDomain();
422 
427  X& lhs();
428 
433  Y& rhs();
434 
439  void relaxResult(field_type relax);
440 
445  void operator()(const size_type& domain);
446 
454  void assignResult(block_type& res);
455 
456  private:
460  const M* mat;
462  X* lhs_;
464  Y* rhs_;
466  const Y* b;
468  X* x;
470  size_type i;
471  };
472 
473  // specialization for ILU0
474  template<class M, class X, class Y>
475  class OverlappingAssignerHelper<ILU0SubdomainSolver<M,X,Y>, false>
476  : public OverlappingAssignerILUBase<M,X,Y>
477  {
478  public:
486  OverlappingAssignerHelper(std::size_t maxlength, const M& mat,
487  const Y& b, X& x)
488  : OverlappingAssignerILUBase<M,X,Y>(maxlength, mat,b,x)
489  {}
490  };
491 
492  // specialization for ILUN
493  template<class M, class X, class Y>
494  class OverlappingAssignerHelper<ILUNSubdomainSolver<M,X,Y>,false>
495  : public OverlappingAssignerILUBase<M,X,Y>
496  {
497  public:
505  OverlappingAssignerHelper(std::size_t maxlength, const M& mat,
506  const Y& b, X& x)
507  : OverlappingAssignerILUBase<M,X,Y>(maxlength, mat,b,x)
508  {}
509  };
510 
511  template<typename S, typename T>
512  struct AdditiveAdder
513  {};
514 
515  template<typename S, typename T, typename A, int n>
516  struct AdditiveAdder<S, BlockVector<FieldVector<T,n>,A> >
517  {
518  typedef typename A::size_type size_type;
519  AdditiveAdder(BlockVector<FieldVector<T,n>,A>& v, BlockVector<FieldVector<T,n>,A>& x,
520  OverlappingAssigner<S>& assigner, const T& relax_);
521  void operator()(const size_type& domain);
522  void axpy();
523 
524  private:
525  BlockVector<FieldVector<T,n>,A>* v;
526  BlockVector<FieldVector<T,n>,A>* x;
527  OverlappingAssigner<S>* assigner;
528  T relax;
529  };
530 
531  template<typename S,typename T>
532  struct MultiplicativeAdder
533  {};
534 
535  template<typename S, typename T, typename A, int n>
536  struct MultiplicativeAdder<S, BlockVector<FieldVector<T,n>,A> >
537  {
538  typedef typename A::size_type size_type;
539  MultiplicativeAdder(BlockVector<FieldVector<T,n>,A>& v, BlockVector<FieldVector<T,n>,A>& x,
540  OverlappingAssigner<S>& assigner_, const T& relax_);
541  void operator()(const size_type& domain);
542  void axpy();
543 
544  private:
545  BlockVector<FieldVector<T,n>,A>* x;
546  OverlappingAssigner<S>* assigner;
547  T relax;
548  };
549 
559  template<typename T, class X, class S>
561  {};
562 
563  template<class X, class S>
565  {
566  typedef AdditiveAdder<S,X> Adder;
567  };
568 
569  template<class X, class S>
570  struct AdderSelector<MultiplicativeSchwarzMode,X,S>
571  {
572  typedef MultiplicativeAdder<S,X> Adder;
573  };
574 
575  template<class X, class S>
576  struct AdderSelector<SymmetricMultiplicativeSchwarzMode,X,S>
577  {
578  typedef MultiplicativeAdder<S,X> Adder;
579  };
580 
592  template<typename T1, typename T2, bool forward>
594  {
595  typedef T1 solver_vector;
596  typedef typename solver_vector::iterator solver_iterator;
597  typedef T2 subdomain_vector;
598  typedef typename subdomain_vector::const_iterator domain_iterator;
599 
600  static solver_iterator begin(solver_vector& sv)
601  {
602  return sv.begin();
603  }
604 
605  static solver_iterator end(solver_vector& sv)
606  {
607  return sv.end();
608  }
609  static domain_iterator begin(const subdomain_vector& sv)
610  {
611  return sv.begin();
612  }
613 
614  static domain_iterator end(const subdomain_vector& sv)
615  {
616  return sv.end();
617  }
618  };
619 
620  template<typename T1, typename T2>
621  struct IteratorDirectionSelector<T1,T2,false>
622  {
623  typedef T1 solver_vector;
624  typedef typename solver_vector::reverse_iterator solver_iterator;
625  typedef T2 subdomain_vector;
626  typedef typename subdomain_vector::const_reverse_iterator domain_iterator;
627 
628  static solver_iterator begin(solver_vector& sv)
629  {
630  return sv.rbegin();
631  }
632 
633  static solver_iterator end(solver_vector& sv)
634  {
635  return sv.rend();
636  }
637  static domain_iterator begin(const subdomain_vector& sv)
638  {
639  return sv.rbegin();
640  }
641 
642  static domain_iterator end(const subdomain_vector& sv)
643  {
644  return sv.rend();
645  }
646  };
647 
656  template<class T>
658  {
659  typedef T smoother;
660  typedef typename smoother::range_type range_type;
661 
662  static void apply(smoother& sm, range_type& v, const range_type& b)
663  {
664  sm.template apply<true>(v, b);
665  }
666  };
667 
668  template<class M, class X, class TD, class TA>
670  {
672  typedef typename smoother::range_type range_type;
673 
674  static void apply(smoother& sm, range_type& v, const range_type& b)
675  {
676  sm.template apply<true>(v, b);
677  sm.template apply<false>(v, b);
678  }
679  };
680 
681  template<class T, bool tag>
682  struct SeqOverlappingSchwarzAssemblerHelper
683  {};
684 
685  template<class T>
686  using SeqOverlappingSchwarzAssembler = SeqOverlappingSchwarzAssemblerHelper<T,Dune::StoresColumnCompressed<T>::value>;
687 
688  template<class K, int n, class Al, class X, class Y>
689  struct SeqOverlappingSchwarzAssemblerHelper< DynamicMatrixSubdomainSolver< BCRSMatrix< FieldMatrix<K,n,n>, Al>, X, Y >,false>
690  {
691  typedef BCRSMatrix< FieldMatrix<K,n,n>, Al> matrix_type;
692  template<class RowToDomain, class Solvers, class SubDomains>
693  static std::size_t assembleLocalProblems(const RowToDomain& rowToDomain, const matrix_type& mat,
694  Solvers& solvers, const SubDomains& domains,
695  bool onTheFly);
696  };
697 
698  template<template<class> class S, typename T, typename A, int m, int n>
699  struct SeqOverlappingSchwarzAssemblerHelper<S<BCRSMatrix<FieldMatrix<T,m,n>,A> >,true>
700  {
701  typedef BCRSMatrix<FieldMatrix<T,m,n>,A> matrix_type;
702  template<class RowToDomain, class Solvers, class SubDomains>
703  static std::size_t assembleLocalProblems(const RowToDomain& rowToDomain, const matrix_type& mat,
704  Solvers& solvers, const SubDomains& domains,
705  bool onTheFly);
706  };
707 
708  template<class M,class X, class Y>
709  struct SeqOverlappingSchwarzAssemblerILUBase
710  {
711  typedef M matrix_type;
712  template<class RowToDomain, class Solvers, class SubDomains>
713  static std::size_t assembleLocalProblems(const RowToDomain& rowToDomain, const matrix_type& mat,
714  Solvers& solvers, const SubDomains& domains,
715  bool onTheFly);
716  };
717 
718  template<class M,class X, class Y>
719  struct SeqOverlappingSchwarzAssemblerHelper<ILU0SubdomainSolver<M,X,Y>,false>
720  : public SeqOverlappingSchwarzAssemblerILUBase<M,X,Y>
721  {};
722 
723  template<class M,class X, class Y>
724  struct SeqOverlappingSchwarzAssemblerHelper<ILUNSubdomainSolver<M,X,Y>,false>
725  : public SeqOverlappingSchwarzAssemblerILUBase<M,X,Y>
726  {};
727 
738  template<class M, class X, class TM=AdditiveSchwarzMode,
739  class TD=ILU0SubdomainSolver<M,X,X>, class TA=std::allocator<X> >
741  : public Preconditioner<X,X>
742  {
743  public:
747  typedef M matrix_type;
748 
752  typedef X domain_type;
753 
757  typedef X range_type;
758 
765  typedef TM Mode;
766 
770  typedef typename X::field_type field_type;
771 
773  typedef typename matrix_type::size_type size_type;
774 
776  typedef TA allocator;
777 
779  typedef std::set<size_type, std::less<size_type>,
780  typename TA::template rebind<size_type>::other>
782 
784  typedef std::vector<subdomain_type, typename TA::template rebind<subdomain_type>::other> subdomain_vector;
785 
788 
790  typedef std::vector<subdomain_list, typename TA::template rebind<subdomain_list>::other > rowtodomain_vector;
791 
793  typedef TD slu;
794 
796  typedef std::vector<slu, typename TA::template rebind<slu>::other> slu_vector;
797 
811  SeqOverlappingSchwarz(const matrix_type& mat, const subdomain_vector& subDomains,
812  field_type relaxationFactor=1, bool onTheFly_=true);
813 
825  SeqOverlappingSchwarz(const matrix_type& mat, const rowtodomain_vector& rowToDomain,
826  field_type relaxationFactor=1, bool onTheFly_=true);
827 
833  virtual void pre (X& x, X& b)
834  {
837  }
838 
844  virtual void apply (X& v, const X& d);
845 
851  virtual void post (X& x)
852  {
854  }
855 
856  template<bool forward>
857  void apply(X& v, const X& d);
858 
861  {
863  }
864 
865  private:
866  const M& mat;
867  slu_vector solvers;
868  subdomain_vector subDomains;
869  field_type relax;
870 
871  typename M::size_type maxlength;
872 
873  bool onTheFly;
874  };
875 
876 
877 
878  template<class I, class S, class D>
879  OverlappingSchwarzInitializer<I,S,D>::OverlappingSchwarzInitializer(InitializerList& il,
880  const IndexSet& idx,
881  const subdomain_vector& domains_)
882  : initializers(&il), indices(&idx), indexMaps(il.size()), domains(domains_)
883  {}
884 
885 
886  template<class I, class S, class D>
887  void OverlappingSchwarzInitializer<I,S,D>::addRowNnz(const Iter& row)
888  {
889  typedef typename IndexSet::value_type::const_iterator iterator;
890  for(iterator domain=(*indices)[row.index()].begin(); domain != (*indices)[row.index()].end(); ++domain) {
891  (*initializers)[*domain].addRowNnz(row, domains[*domain]);
892  indexMaps[*domain].insert(row.index());
893  }
894  }
895 
896  template<class I, class S, class D>
897  void OverlappingSchwarzInitializer<I,S,D>::allocate()
898  {
899  for(auto&& i: *initializers)
900  i.allocateMatrixStorage();
901  for(auto&& i: *initializers)
902  i.allocateMarker();
903  }
904 
905  template<class I, class S, class D>
906  void OverlappingSchwarzInitializer<I,S,D>::countEntries(const Iter& row, const CIter& col) const
907  {
908  typedef typename IndexSet::value_type::const_iterator iterator;
909  for(iterator domain=(*indices)[row.index()].begin(); domain != (*indices)[row.index()].end(); ++domain) {
910  typename std::map<size_type,size_type>::const_iterator v = indexMaps[*domain].find(col.index());
911  if(v!= indexMaps[*domain].end()) {
912  (*initializers)[*domain].countEntries(indexMaps[*domain].find(col.index())->second);
913  }
914  }
915  }
916 
917  template<class I, class S, class D>
918  void OverlappingSchwarzInitializer<I,S,D>::calcColstart() const
919  {
920  for(auto&& i : *initializers)
921  i.calcColstart();
922  }
923 
924  template<class I, class S, class D>
925  void OverlappingSchwarzInitializer<I,S,D>::copyValue(const Iter& row, const CIter& col) const
926  {
927  typedef typename IndexSet::value_type::const_iterator iterator;
928  for(iterator domain=(*indices)[row.index()].begin(); domain!= (*indices)[row.index()].end(); ++domain) {
929  typename std::map<size_type,size_type>::const_iterator v = indexMaps[*domain].find(col.index());
930  if(v!= indexMaps[*domain].end()) {
931  assert(indexMaps[*domain].end()!=indexMaps[*domain].find(row.index()));
932  (*initializers)[*domain].copyValue(col, indexMaps[*domain].find(row.index())->second,
933  v->second);
934  }
935  }
936  }
937 
938  template<class I, class S, class D>
939  void OverlappingSchwarzInitializer<I,S,D>::createMatrix() const
940  {
941  std::vector<IndexMap>().swap(indexMaps);
942  for(auto&& i: *initializers)
943  i.createMatrix();
944  }
945 
946  template<class I, class S, class D>
947  OverlappingSchwarzInitializer<I,S,D>::IndexMap::IndexMap()
948  : row(0)
949  {}
950 
951  template<class I, class S, class D>
952  void OverlappingSchwarzInitializer<I,S,D>::IndexMap::insert(size_type grow)
953  {
954  assert(map_.find(grow)==map_.end());
955  map_.insert(std::make_pair(grow, row++));
956  }
957 
958  template<class I, class S, class D>
959  typename OverlappingSchwarzInitializer<I,S,D>::IndexMap::const_iterator
960  OverlappingSchwarzInitializer<I,S,D>::IndexMap::find(size_type grow) const
961  {
962  return map_.find(grow);
963  }
964 
965  template<class I, class S, class D>
966  typename OverlappingSchwarzInitializer<I,S,D>::IndexMap::iterator
967  OverlappingSchwarzInitializer<I,S,D>::IndexMap::find(size_type grow)
968  {
969  return map_.find(grow);
970  }
971 
972  template<class I, class S, class D>
973  typename OverlappingSchwarzInitializer<I,S,D>::IndexMap::const_iterator
974  OverlappingSchwarzInitializer<I,S,D>::IndexMap::end() const
975  {
976  return map_.end();
977  }
978 
979  template<class I, class S, class D>
980  typename OverlappingSchwarzInitializer<I,S,D>::IndexMap::iterator
981  OverlappingSchwarzInitializer<I,S,D>::IndexMap::end()
982  {
983  return map_.end();
984  }
985 
986  template<class I, class S, class D>
987  typename OverlappingSchwarzInitializer<I,S,D>::IndexMap::const_iterator
988  OverlappingSchwarzInitializer<I,S,D>::IndexMap::begin() const
989  {
990  return map_.begin();
991  }
992 
993  template<class I, class S, class D>
994  typename OverlappingSchwarzInitializer<I,S,D>::IndexMap::iterator
995  OverlappingSchwarzInitializer<I,S,D>::IndexMap::begin()
996  {
997  return map_.begin();
998  }
999 
1000  template<class M, class X, class TM, class TD, class TA>
1002  field_type relaxationFactor, bool fly)
1003  : mat(mat_), relax(relaxationFactor), onTheFly(fly)
1004  {
1005  typedef typename rowtodomain_vector::const_iterator RowDomainIterator;
1006  typedef typename subdomain_list::const_iterator DomainIterator;
1007 #ifdef DUNE_ISTL_WITH_CHECKING
1008  assert(rowToDomain.size()==mat.N());
1009  assert(rowToDomain.size()==mat.M());
1010 
1011  for(RowDomainIterator iter=rowToDomain.begin(); iter != rowToDomain.end(); ++iter)
1012  assert(iter->size()>0);
1013 
1014 #endif
1015  // calculate the number of domains
1016  size_type domains=0;
1017  for(RowDomainIterator iter=rowToDomain.begin(); iter != rowToDomain.end(); ++iter)
1018  for(DomainIterator d=iter->begin(); d != iter->end(); ++d)
1019  domains=std::max(domains, *d);
1020  ++domains;
1021 
1022  solvers.resize(domains);
1023  subDomains.resize(domains);
1024 
1025  // initialize subdomains to row mapping from row to subdomain mapping
1026  size_type row=0;
1027  for(RowDomainIterator iter=rowToDomain.begin(); iter != rowToDomain.end(); ++iter, ++row)
1028  for(DomainIterator d=iter->begin(); d != iter->end(); ++d)
1029  subDomains[*d].insert(row);
1030 
1031 #ifdef DUNE_ISTL_WITH_CHECKING
1032  size_type i=0;
1033  typedef typename subdomain_vector::const_iterator iterator;
1034  for(iterator iter=subDomains.begin(); iter != subDomains.end(); ++iter) {
1035  typedef typename subdomain_type::const_iterator entry_iterator;
1036  Dune::dvverb<<"domain "<<i++<<":";
1037  for(entry_iterator entry = iter->begin(); entry != iter->end(); ++entry) {
1038  Dune::dvverb<<" "<<*entry;
1039  }
1040  Dune::dvverb<<std::endl;
1041  }
1042 #endif
1043  maxlength = SeqOverlappingSchwarzAssembler<slu>
1044  ::assembleLocalProblems(rowToDomain, mat, solvers, subDomains, onTheFly);
1045  }
1046 
1047  template<class M, class X, class TM, class TD, class TA>
1049  const subdomain_vector& sd,
1050  field_type relaxationFactor,
1051  bool fly)
1052  : mat(mat_), solvers(sd.size()), subDomains(sd), relax(relaxationFactor),
1053  onTheFly(fly)
1054  {
1055  typedef typename subdomain_vector::const_iterator DomainIterator;
1056 
1057 #ifdef DUNE_ISTL_WITH_CHECKING
1058  size_type i=0;
1059 
1060  for(DomainIterator d=sd.begin(); d != sd.end(); ++d,++i) {
1061  //std::cout<<i<<": "<<d->size()<<std::endl;
1062  assert(d->size()>0);
1063  typedef typename DomainIterator::value_type::const_iterator entry_iterator;
1064  Dune::dvverb<<"domain "<<i<<":";
1065  for(entry_iterator entry = d->begin(); entry != d->end(); ++entry) {
1066  Dune::dvverb<<" "<<*entry;
1067  }
1068  Dune::dvverb<<std::endl;
1069  }
1070 
1071 #endif
1072 
1073  // Create a row to subdomain mapping
1074  rowtodomain_vector rowToDomain(mat.N());
1075 
1076  size_type domainId=0;
1077 
1078  for(DomainIterator domain=sd.begin(); domain != sd.end(); ++domain, ++domainId) {
1079  typedef typename subdomain_type::const_iterator iterator;
1080  for(iterator row=domain->begin(); row != domain->end(); ++row)
1081  rowToDomain[*row].push_back(domainId);
1082  }
1083 
1084  maxlength = SeqOverlappingSchwarzAssembler<slu>
1085  ::assembleLocalProblems(rowToDomain, mat, solvers, subDomains, onTheFly);
1086  }
1087 
1094  template<class M>
1096 
1097  template<typename T, typename A, int n, int m>
1099  {
1100  template<class Domain>
1101  static int size(const Domain & d)
1102  {
1103  assert(n==m);
1104  return m*d.size();
1105  }
1106  };
1107 
1108  template<class K, int n, class Al, class X, class Y>
1109  template<class RowToDomain, class Solvers, class SubDomains>
1110  std::size_t
1111  SeqOverlappingSchwarzAssemblerHelper< DynamicMatrixSubdomainSolver< BCRSMatrix< FieldMatrix<K,n,n>, Al>, X, Y >,false>::
1112  assembleLocalProblems(const RowToDomain& rowToDomain,
1113  const matrix_type& mat,
1114  Solvers& solvers,
1115  const SubDomains& subDomains,
1116  bool onTheFly)
1117  {
1118  DUNE_UNUSED_PARAMETER(onTheFly);
1119  DUNE_UNUSED_PARAMETER(rowToDomain);
1120  DUNE_UNUSED_PARAMETER(mat);
1121  DUNE_UNUSED_PARAMETER(solvers);
1122  typedef typename SubDomains::const_iterator DomainIterator;
1123  std::size_t maxlength = 0;
1124 
1125  assert(onTheFly);
1126 
1127  for(DomainIterator domain=subDomains.begin(); domain!=subDomains.end(); ++domain)
1128  maxlength=std::max(maxlength, domain->size());
1129  maxlength*=n;
1130 
1131  return maxlength;
1132  }
1133 
1134 #if HAVE_SUPERLU || HAVE_SUITESPARSE_UMFPACK
1135  template<template<class> class S, typename T, typename A, int m, int n>
1136  template<class RowToDomain, class Solvers, class SubDomains>
1137  std::size_t SeqOverlappingSchwarzAssemblerHelper<S<BCRSMatrix<FieldMatrix<T,m,n>,A> >,true>::assembleLocalProblems(const RowToDomain& rowToDomain,
1138  const matrix_type& mat,
1139  Solvers& solvers,
1140  const SubDomains& subDomains,
1141  bool onTheFly)
1142  {
1143  typedef typename S<BCRSMatrix<FieldMatrix<T,m,n>,A> >::MatrixInitializer MatrixInitializer;
1144  typedef typename std::vector<MatrixInitializer>::iterator InitializerIterator;
1145  typedef typename SubDomains::const_iterator DomainIterator;
1146  typedef typename Solvers::iterator SolverIterator;
1147  std::size_t maxlength = 0;
1148 
1149  if(onTheFly) {
1150  for(DomainIterator domain=subDomains.begin(); domain!=subDomains.end(); ++domain)
1151  maxlength=std::max(maxlength, domain->size());
1152  maxlength*=mat[0].begin()->N();
1153  }else{
1154  // initialize the initializers
1155  DomainIterator domain=subDomains.begin();
1156 
1157  // Create the initializers list.
1158  std::vector<MatrixInitializer> initializers(subDomains.size());
1159 
1160  SolverIterator solver=solvers.begin();
1161  for(InitializerIterator initializer=initializers.begin(); initializer!=initializers.end();
1162  ++initializer, ++solver, ++domain) {
1163  solver->getInternalMatrix().N_=SeqOverlappingSchwarzDomainSize<matrix_type>::size(*domain);
1164  solver->getInternalMatrix().M_=SeqOverlappingSchwarzDomainSize<matrix_type>::size(*domain);
1165  //solver->setVerbosity(true);
1166  *initializer=MatrixInitializer(solver->getInternalMatrix());
1167  }
1168 
1169  // Set up the supermatrices according to the subdomains
1170  typedef OverlappingSchwarzInitializer<std::vector<MatrixInitializer>,
1171  RowToDomain, SubDomains> Initializer;
1172 
1173  Initializer initializer(initializers, rowToDomain, subDomains);
1174  copyToColCompMatrix(initializer, mat);
1175 
1176  // Calculate the LU decompositions
1177  for(auto&& s: solvers)
1178  s.decompose();
1179  for (SolverIterator solverIt = solvers.begin(); solverIt != solvers.end(); ++solverIt)
1180  {
1181  assert(solverIt->getInternalMatrix().N() == solverIt->getInternalMatrix().M());
1182  maxlength = std::max(maxlength, solverIt->getInternalMatrix().N());
1183  }
1184  }
1185  return maxlength;
1186  }
1187 
1188 #endif // HAVE_SUPERLU || HAVE_SUITESPARSE_UMFPACK
1189 
1190  template<class M,class X,class Y>
1191  template<class RowToDomain, class Solvers, class SubDomains>
1192  std::size_t SeqOverlappingSchwarzAssemblerILUBase<M,X,Y>::assembleLocalProblems(const RowToDomain& rowToDomain,
1193  const matrix_type& mat,
1194  Solvers& solvers,
1195  const SubDomains& subDomains,
1196  bool onTheFly)
1197  {
1198  DUNE_UNUSED_PARAMETER(rowToDomain);
1199  typedef typename SubDomains::const_iterator DomainIterator;
1200  typedef typename Solvers::iterator SolverIterator;
1201  std::size_t maxlength = 0;
1202 
1203  if(onTheFly) {
1204  for(DomainIterator domain=subDomains.begin(); domain!=subDomains.end(); ++domain)
1205  maxlength=std::max(maxlength, domain->size());
1206  }else{
1207  // initialize the solvers of the local prolems.
1208  SolverIterator solver=solvers.begin();
1209  for(DomainIterator domain=subDomains.begin(); domain!=subDomains.end();
1210  ++domain, ++solver) {
1211  solver->setSubMatrix(mat, *domain);
1212  maxlength=std::max(maxlength, domain->size());
1213  }
1214  }
1215 
1216  return maxlength;
1217 
1218  }
1219 
1220 
1221  template<class M, class X, class TM, class TD, class TA>
1223  {
1225  }
1226 
1227  template<class M, class X, class TM, class TD, class TA>
1228  template<bool forward>
1229  void SeqOverlappingSchwarz<M,X,TM,TD,TA>::apply(X& x, const X& b)
1230  {
1231  typedef slu_vector solver_vector;
1232  typedef typename IteratorDirectionSelector<solver_vector,subdomain_vector,forward>::solver_iterator iterator;
1233  typedef typename IteratorDirectionSelector<solver_vector,subdomain_vector,forward>::domain_iterator
1234  domain_iterator;
1235 
1236  OverlappingAssigner<TD> assigner(maxlength, mat, b, x);
1237 
1240  X v(x); // temporary for the update
1241  v=0;
1242 
1243  typedef typename AdderSelector<TM,X,TD >::Adder Adder;
1244  Adder adder(v, x, assigner, relax);
1245 
1246  for(; domain != IteratorDirectionSelector<solver_vector,subdomain_vector,forward>::end(subDomains); ++domain) {
1247  //Copy rhs to C-array for SuperLU
1248  std::for_each(domain->begin(), domain->end(), assigner);
1249  assigner.resetIndexForNextDomain();
1250  if(onTheFly) {
1251  // Create the subdomain solver
1252  slu sdsolver;
1253  sdsolver.setSubMatrix(mat, *domain);
1254  // Apply
1255  sdsolver.apply(assigner.lhs(), assigner.rhs());
1256  }else{
1257  solver->apply(assigner.lhs(), assigner.rhs());
1258  ++solver;
1259  }
1260 
1261  //Add relaxed correction to from SuperLU to v
1262  std::for_each(domain->begin(), domain->end(), adder);
1263  assigner.resetIndexForNextDomain();
1264 
1265  }
1266 
1267  adder.axpy();
1268  assigner.deallocate();
1269  }
1270 
1271  template<class K, int n, class Al, class X, class Y>
1272  OverlappingAssignerHelper< DynamicMatrixSubdomainSolver< BCRSMatrix< FieldMatrix<K,n,n>, Al>, X, Y >,false>
1273  ::OverlappingAssignerHelper(std::size_t maxlength, const BCRSMatrix<FieldMatrix<K,n,n>, Al>& mat_,
1274  const X& b_, Y& x_) :
1275  mat(&mat_),
1276  rhs_( new DynamicVector<field_type>(maxlength, 42) ),
1277  lhs_( new DynamicVector<field_type>(maxlength, -42) ),
1278  b(&b_),
1279  x(&x_),
1280  i(0),
1281  maxlength_(maxlength)
1282  {}
1283 
1284  template<class K, int n, class Al, class X, class Y>
1285  void
1286  OverlappingAssignerHelper< DynamicMatrixSubdomainSolver< BCRSMatrix< FieldMatrix<K,n,n>, Al>, X, Y >,false>
1287  ::deallocate()
1288  {
1289  delete rhs_;
1290  delete lhs_;
1291  }
1292 
1293  template<class K, int n, class Al, class X, class Y>
1294  void
1295  OverlappingAssignerHelper< DynamicMatrixSubdomainSolver< BCRSMatrix< FieldMatrix<K,n,n>, Al>, X, Y >,false>
1296  ::resetIndexForNextDomain()
1297  {
1298  i=0;
1299  }
1300 
1301  template<class K, int n, class Al, class X, class Y>
1303  OverlappingAssignerHelper< DynamicMatrixSubdomainSolver< BCRSMatrix< FieldMatrix<K,n,n>, Al>, X, Y >,false>
1304  ::lhs()
1305  {
1306  return *lhs_;
1307  }
1308 
1309  template<class K, int n, class Al, class X, class Y>
1311  OverlappingAssignerHelper< DynamicMatrixSubdomainSolver< BCRSMatrix< FieldMatrix<K,n,n>, Al>, X, Y >,false>
1312  ::rhs()
1313  {
1314  return *rhs_;
1315  }
1316 
1317  template<class K, int n, class Al, class X, class Y>
1318  void
1319  OverlappingAssignerHelper< DynamicMatrixSubdomainSolver< BCRSMatrix< FieldMatrix<K,n,n>, Al>, X, Y >,false>
1320  ::relaxResult(field_type relax)
1321  {
1322  lhs() *= relax;
1323  }
1324 
1325  template<class K, int n, class Al, class X, class Y>
1326  void
1327  OverlappingAssignerHelper< DynamicMatrixSubdomainSolver< BCRSMatrix< FieldMatrix<K,n,n>, Al>, X, Y >,false>
1328  ::operator()(const size_type& domainIndex)
1329  {
1330  lhs() = 0.0;
1331 #if 0
1332  //assign right hand side of current domainindex block
1333  for(size_type j=0; j<n; ++j, ++i) {
1334  assert(i<maxlength_);
1335  rhs()[i]=(*b)[domainIndex][j];
1336  }
1337 
1338  // loop over all Matrix row entries and calculate defect.
1339  typedef typename matrix_type::ConstColIterator col_iterator;
1340 
1341  // calculate defect for current row index block
1342  for(col_iterator col=(*mat)[domainIndex].begin(); col!=(*mat)[domainIndex].end(); ++col) {
1343  block_type tmp(0.0);
1344  (*col).mv((*x)[col.index()], tmp);
1345  i-=n;
1346  for(size_type j=0; j<n; ++j, ++i) {
1347  assert(i<maxlength_);
1348  rhs()[i]-=tmp[j];
1349  }
1350  }
1351 #else
1352  //assign right hand side of current domainindex block
1353  for(size_type j=0; j<n; ++j, ++i) {
1354  assert(i<maxlength_);
1355  rhs()[i]=(*b)[domainIndex][j];
1356 
1357  // loop over all Matrix row entries and calculate defect.
1358  typedef typename matrix_type::ConstColIterator col_iterator;
1359 
1360  // calculate defect for current row index block
1361  for(col_iterator col=(*mat)[domainIndex].begin(); col!=(*mat)[domainIndex].end(); ++col) {
1362  for(size_type k=0; k<n; ++k) {
1363  rhs()[i]-=(*col)[j][k] * (*x)[col.index()][k];
1364  }
1365  }
1366  }
1367 #endif
1368  }
1369 
1370  template<class K, int n, class Al, class X, class Y>
1371  void
1372  OverlappingAssignerHelper< DynamicMatrixSubdomainSolver< BCRSMatrix< FieldMatrix<K,n,n>, Al>, X, Y >,false>
1373  ::assignResult(block_type& res)
1374  {
1375  // assign the result of the local solve to the global vector
1376  for(size_type j=0; j<n; ++j, ++i) {
1377  assert(i<maxlength_);
1378  res[j]+=lhs()[i];
1379  }
1380  }
1381 
1382 #if HAVE_SUPERLU || HAVE_SUITESPARSE_UMFPACK
1383 
1384  template<template<class> class S, int n, int m, typename T, typename A>
1385  OverlappingAssignerHelper<S<BCRSMatrix<FieldMatrix<T,n,m>,A> >,true>
1386  ::OverlappingAssignerHelper(std::size_t maxlength,
1387  const BCRSMatrix<FieldMatrix<T,n,m>,A>& mat_,
1388  const range_type& b_,
1389  range_type& x_)
1390  : mat(&mat_),
1391  b(&b_),
1392  x(&x_), i(0), maxlength_(maxlength)
1393  {
1394  rhs_ = new field_type[maxlength];
1395  lhs_ = new field_type[maxlength];
1396 
1397  }
1398 
1399  template<template<class> class S, int n, int m, typename T, typename A>
1400  void OverlappingAssignerHelper<S<BCRSMatrix<FieldMatrix<T,n,m>,A> >,true>::deallocate()
1401  {
1402  delete[] rhs_;
1403  delete[] lhs_;
1404  }
1405 
1406  template<template<class> class S, int n, int m, typename T, typename A>
1407  void OverlappingAssignerHelper<S<BCRSMatrix<FieldMatrix<T,n,m>,A> >,true>::operator()(const size_type& domainIndex)
1408  {
1409  //assign right hand side of current domainindex block
1410  // rhs is an array of doubles!
1411  // rhs[starti] = b[domainindex]
1412  for(size_type j=0; j<n; ++j, ++i) {
1413  assert(i<maxlength_);
1414  rhs_[i]=(*b)[domainIndex][j];
1415  }
1416 
1417 
1418  // loop over all Matrix row entries and calculate defect.
1419  typedef typename matrix_type::ConstColIterator col_iterator;
1420 
1421  // calculate defect for current row index block
1422  for(col_iterator col=(*mat)[domainIndex].begin(); col!=(*mat)[domainIndex].end(); ++col) {
1423  block_type tmp;
1424  (*col).mv((*x)[col.index()], tmp);
1425  i-=n;
1426  for(size_type j=0; j<n; ++j, ++i) {
1427  assert(i<maxlength_);
1428  rhs_[i]-=tmp[j];
1429  }
1430 
1431  }
1432 
1433  }
1434 
1435  template<template<class> class S, int n, int m, typename T, typename A>
1436  void OverlappingAssignerHelper<S<BCRSMatrix<FieldMatrix<T,n,m>,A> >,true>::relaxResult(field_type relax)
1437  {
1438  for(size_type j=i+n; i<j; ++i) {
1439  assert(i<maxlength_);
1440  lhs_[i]*=relax;
1441  }
1442  i-=n;
1443  }
1444 
1445  template<template<class> class S, int n, int m, typename T, typename A>
1446  void OverlappingAssignerHelper<S<BCRSMatrix<FieldMatrix<T,n,m>,A> >,true>::assignResult(block_type& res)
1447  {
1448  // assign the result of the local solve to the global vector
1449  for(size_type j=0; j<n; ++j, ++i) {
1450  assert(i<maxlength_);
1451  res[j]+=lhs_[i];
1452  }
1453  }
1454 
1455  template<template<class> class S, int n, int m, typename T, typename A>
1456  void OverlappingAssignerHelper<S<BCRSMatrix<FieldMatrix<T,n,m>,A> >,true>::resetIndexForNextDomain()
1457  {
1458  i=0;
1459  }
1460 
1461  template<template<class> class S, int n, int m, typename T, typename A>
1462  typename OverlappingAssignerHelper<S<BCRSMatrix<FieldMatrix<T,n,m>,A> >,true>::field_type*
1463  OverlappingAssignerHelper<S<BCRSMatrix<FieldMatrix<T,n,m>,A> >,true>::lhs()
1464  {
1465  return lhs_;
1466  }
1467 
1468  template<template<class> class S, int n, int m, typename T, typename A>
1469  typename OverlappingAssignerHelper<S<BCRSMatrix<FieldMatrix<T,n,m>,A> >,true>::field_type*
1470  OverlappingAssignerHelper<S<BCRSMatrix<FieldMatrix<T,n,m>,A> >,true>::rhs()
1471  {
1472  return rhs_;
1473  }
1474 
1475 #endif // HAVE_SUPERLU || HAVE_SUITESPARSE_UMFPACK
1476 
1477  template<class M, class X, class Y>
1478  OverlappingAssignerILUBase<M,X,Y>::OverlappingAssignerILUBase(std::size_t maxlength,
1479  const M& mat_,
1480  const Y& b_,
1481  X& x_)
1482  : mat(&mat_),
1483  b(&b_),
1484  x(&x_), i(0)
1485  {
1486  rhs_= new Y(maxlength);
1487  lhs_ = new X(maxlength);
1488  }
1489 
1490  template<class M, class X, class Y>
1491  void OverlappingAssignerILUBase<M,X,Y>::deallocate()
1492  {
1493  delete rhs_;
1494  delete lhs_;
1495  }
1496 
1497  template<class M, class X, class Y>
1498  void OverlappingAssignerILUBase<M,X,Y>::operator()(const size_type& domainIndex)
1499  {
1500  (*rhs_)[i]=(*b)[domainIndex];
1501 
1502  // loop over all Matrix row entries and calculate defect.
1503  typedef typename matrix_type::ConstColIterator col_iterator;
1504 
1505  // calculate defect for current row index block
1506  for(col_iterator col=(*mat)[domainIndex].begin(); col!=(*mat)[domainIndex].end(); ++col) {
1507  (*col).mmv((*x)[col.index()], (*rhs_)[i]);
1508  }
1509  // Goto next local index
1510  ++i;
1511  }
1512 
1513  template<class M, class X, class Y>
1514  void OverlappingAssignerILUBase<M,X,Y>::relaxResult(field_type relax)
1515  {
1516  (*lhs_)[i]*=relax;
1517  }
1518 
1519  template<class M, class X, class Y>
1520  void OverlappingAssignerILUBase<M,X,Y>::assignResult(block_type& res)
1521  {
1522  res+=(*lhs_)[i++];
1523  }
1524 
1525  template<class M, class X, class Y>
1526  X& OverlappingAssignerILUBase<M,X,Y>::lhs()
1527  {
1528  return *lhs_;
1529  }
1530 
1531  template<class M, class X, class Y>
1532  Y& OverlappingAssignerILUBase<M,X,Y>::rhs()
1533  {
1534  return *rhs_;
1535  }
1536 
1537  template<class M, class X, class Y>
1538  void OverlappingAssignerILUBase<M,X,Y>::resetIndexForNextDomain()
1539  {
1540  i=0;
1541  }
1542 
1543  template<typename S, typename T, typename A, int n>
1544  AdditiveAdder<S,BlockVector<FieldVector<T,n>,A> >::AdditiveAdder(BlockVector<FieldVector<T,n>,A>& v_,
1546  OverlappingAssigner<S>& assigner_,
1547  const T& relax_)
1548  : v(&v_), x(&x_), assigner(&assigner_), relax(relax_)
1549  {}
1550 
1551  template<typename S, typename T, typename A, int n>
1552  void AdditiveAdder<S,BlockVector<FieldVector<T,n>,A> >::operator()(const size_type& domainIndex)
1553  {
1554  // add the result of the local solve to the current update
1555  assigner->assignResult((*v)[domainIndex]);
1556  }
1557 
1558 
1559  template<typename S, typename T, typename A, int n>
1560  void AdditiveAdder<S,BlockVector<FieldVector<T,n>,A> >::axpy()
1561  {
1562  // relax the update and add it to the current guess.
1563  x->axpy(relax,*v);
1564  }
1565 
1566 
1567  template<typename S, typename T, typename A, int n>
1568  MultiplicativeAdder<S,BlockVector<FieldVector<T,n>,A> >
1569  ::MultiplicativeAdder(BlockVector<FieldVector<T,n>,A>& v_,
1570  BlockVector<FieldVector<T,n>,A>& x_,
1571  OverlappingAssigner<S>& assigner_, const T& relax_)
1572  : x(&x_), assigner(&assigner_), relax(relax_)
1573  {
1575  }
1576 
1577 
1578  template<typename S,typename T, typename A, int n>
1579  void MultiplicativeAdder<S,BlockVector<FieldVector<T,n>,A> >::operator()(const size_type& domainIndex)
1580  {
1581  // add the result of the local solve to the current guess
1582  assigner->relaxResult(relax);
1583  assigner->assignResult((*x)[domainIndex]);
1584  }
1585 
1586 
1587  template<typename S,typename T, typename A, int n>
1588  void MultiplicativeAdder<S,BlockVector<FieldVector<T,n>,A> >::axpy()
1589  {
1590  // nothing to do, as the corrections already relaxed and added in operator()
1591  }
1592 
1593 
1595 }
1596 
1597 #endif
Implementation of the BCRSMatrix class.
This file implements a vector space as a tensor product of a given vector space. The number of compon...
A sparse block matrix with compressed row storage.
Definition: bcrsmatrix.hh:423
row_type::ConstIterator ConstColIterator
Const iterator to the entries of a row.
Definition: bcrsmatrix.hh:700
A vector of blocks with memory management.
Definition: bvector.hh:317
Iterator begin()
begin iterator
Definition: densevector.hh:308
Iterator end()
end iterator
Definition: densevector.hh:314
Iterator find(size_type i)
return iterator to given element or end()
Definition: densevector.hh:334
Exact subdomain solver using Dune::DynamicMatrix<T>::solve.
Definition: overlappingschwarz.hh:136
size_type capacity() const
Number of elements for which memory has been allocated.
Definition: dynvector.hh:135
A dense n x m matrix.
Definition: fmatrix.hh:68
vector space out of a tensor product of fields.
Definition: fvector.hh:93
Index Set Interface base class.
Definition: indexidset.hh:76
RowIterator begin()
Get iterator to first row.
Definition: matrix.hh:609
Initializer for SuperLU Matrices representing the subdomains.
Definition: overlappingschwarz.hh:43
D subdomain_vector
The vector type containing the subdomain to row index mapping.
Definition: overlappingschwarz.hh:46
Base class for matrix free definition of preconditioners.
Definition: preconditioner.hh:30
A constant iterator for the SLList.
Definition: sllist.hh:369
A single linked list.
Definition: sllist.hh:42
Sequential overlapping Schwarz preconditioner.
Definition: overlappingschwarz.hh:742
X::field_type field_type
The field type of the preconditioner.
Definition: overlappingschwarz.hh:770
void apply(X &v, const X &d)
Apply one step of the preconditioner to the system A(v)=d.
M matrix_type
The type of the matrix to precondition.
Definition: overlappingschwarz.hh:747
TM Mode
The mode (additive or multiplicative) of the Schwarz method.
Definition: overlappingschwarz.hh:765
X range_type
The range type of the preconditioner.
Definition: overlappingschwarz.hh:757
X domain_type
The domain type of the preconditioner.
Definition: overlappingschwarz.hh:752
TD slu
The type for the subdomain solver in use.
Definition: overlappingschwarz.hh:793
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: overlappingschwarz.hh:860
virtual void pre(X &x, X &b)
Prepare the preconditioner.
Definition: overlappingschwarz.hh:833
virtual void post(X &x)
Postprocess the preconditioner.
Definition: overlappingschwarz.hh:851
std::vector< subdomain_type, typename TA::template rebind< subdomain_type >::other > subdomain_vector
The vector type containing the subdomain to row index mapping.
Definition: overlappingschwarz.hh:784
TA allocator
The allocator to use.
Definition: overlappingschwarz.hh:776
SLList< size_type, typename TA::template rebind< size_type >::other > subdomain_list
The type for the row to subdomain mapping.
Definition: overlappingschwarz.hh:787
std::vector< subdomain_list, typename TA::template rebind< subdomain_list >::other > rowtodomain_vector
The vector type containing the row index to subdomain mapping.
Definition: overlappingschwarz.hh:790
std::vector< slu, typename TA::template rebind< slu >::other > slu_vector
The vector type containing subdomain solvers.
Definition: overlappingschwarz.hh:796
matrix_type::size_type size_type
The return type of the size method.
Definition: overlappingschwarz.hh:773
std::set< size_type, std::less< size_type >, typename TA::template rebind< size_type >::other > subdomain_type
The type for the subdomain to row index mapping.
Definition: overlappingschwarz.hh:781
This file implements a dense matrix with dynamic numbers of rows and columns.
decltype(auto) apply(F &&f, ArgTuple &&args)
Apply function with arguments given as tuple.
Definition: apply.hh:58
#define DUNE_UNUSED_PARAMETER(parm)
A macro to mark intentionally unused function parameters with.
Definition: unused.hh:25
virtual void apply(X &v, const X &d)
Apply the precondtioner.
Definition: overlappingschwarz.hh:1222
SeqOverlappingSchwarz(const matrix_type &mat, const subdomain_vector &subDomains, field_type relaxationFactor=1, bool onTheFly_=true)
Construct the overlapping Schwarz method.
Definition: overlappingschwarz.hh:1048
SeqOverlappingSchwarz(const matrix_type &mat, const rowtodomain_vector &rowToDomain, field_type relaxationFactor=1, bool onTheFly_=true)
Definition: overlappingschwarz.hh:1001
DVVerbType dvverb(std::cout)
stream for very verbose output.
Definition: stdstreams.hh:93
Various local subdomain solvers based on ILU for SeqOverlappingSchwarz.
Dune namespace.
Definition: alignedallocator.hh:10
Define general preconditioner interface.
Implements a singly linked list together with the necessary iterators.
Templates characterizing the type of a solver.
template meta program for choosing how to add the correction.
Definition: overlappingschwarz.hh:561
Tag that the tells the schwarz method to be additive.
Definition: overlappingschwarz.hh:116
Helper template meta program for application of overlapping schwarz.
Definition: overlappingschwarz.hh:594
Tag that tells the Schwarz method to be multiplicative.
Definition: overlappingschwarz.hh:122
Helper template meta program for application of overlapping schwarz.
Definition: overlappingschwarz.hh:658
Definition: overlappingschwarz.hh:1095
Category
Definition: solvercategory.hh:21
@ sequential
Category for sequential solvers.
Definition: solvercategory.hh:23
Tag that tells the Schwarz method to be multiplicative and symmetric.
Definition: overlappingschwarz.hh:129
Classes for using SuperLU with ISTL matrices.
Classes for using UMFPack with ISTL matrices.
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.80.0 (Apr 26, 22:29, 2024)