DUNE PDELab (2.7)

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 <memory>
9#include <vector>
10#include <set>
12#include <dune/common/sllist.hh>
13#include <dune/common/unused.hh>
14#include "preconditioners.hh"
15#include "superlu.hh"
16#include "umfpack.hh"
17#include "bvector.hh"
18#include "bcrsmatrix.hh"
19#include "ilusubdomainsolver.hh"
21
22namespace Dune
23{
24
36 template<class M, class X, class TM, class TD, class TA>
37 class SeqOverlappingSchwarz;
38
42 template<class I, class S, class D>
44 {
45 public:
48
49 typedef I InitializerList;
50 typedef typename InitializerList::value_type AtomInitializer;
51 typedef typename AtomInitializer::Matrix Matrix;
52 typedef typename Matrix::const_iterator Iter;
53 typedef typename Matrix::row_type::const_iterator CIter;
54
55 typedef S IndexSet;
56 typedef typename IndexSet::size_type size_type;
57
58 OverlappingSchwarzInitializer(InitializerList& il,
59 const IndexSet& indices,
60 const subdomain_vector& domains);
61
62
63 void addRowNnz(const Iter& row);
64
65 void allocate();
66
67 void countEntries(const Iter& row, const CIter& col) const;
68
69 void calcColstart() const;
70
71 void copyValue(const Iter& row, const CIter& col) const;
72
73 void createMatrix() const;
74 private:
75 class IndexMap
76 {
77 public:
78 typedef typename S::size_type size_type;
79 typedef std::map<size_type,size_type> Map;
80 typedef typename Map::iterator iterator;
81 typedef typename Map::const_iterator const_iterator;
82
83 IndexMap();
84
85 void insert(size_type grow);
86
87 const_iterator find(size_type grow) const;
88
89 iterator find(size_type grow);
90
91 iterator begin();
92
93 const_iterator begin() const;
94
95 iterator end();
96
97 const_iterator end() const;
98
99 private:
100 std::map<size_type,size_type> map_;
101 size_type row;
102 };
103
104
105 typedef typename InitializerList::iterator InitIterator;
106 typedef typename IndexSet::const_iterator IndexIteratur;
107 InitializerList* initializers;
108 const IndexSet *indices;
109 mutable std::vector<IndexMap> indexMaps;
110 const subdomain_vector& domains;
111 };
112
117 {};
118
123 {};
124
130 {};
131
136 template<class M, class X, class Y>
138
139 // Specialization for BCRSMatrix
140 template<class K, class Al, class X, class Y>
141 class DynamicMatrixSubdomainSolver< BCRSMatrix< K, Al>, X, Y >
142 {
143 typedef BCRSMatrix< K, Al> M;
144 public:
146 typedef typename std::remove_const<M>::type matrix_type;
147 typedef typename X::field_type field_type;
148 typedef typename std::remove_const<M>::type rilu_type;
150 typedef X domain_type;
152 typedef Y range_type;
153 static constexpr size_t n = std::decay_t<decltype(Impl::asMatrix(std::declval<K>()))>::rows;
154
160 {
161 assert(v.size() > 0);
162 assert(v.size() == d.size());
163 assert(A.rows() <= v.size());
164 assert(A.cols() <= v.size());
165 size_t sz = A.rows();
166 v.resize(sz);
167 d.resize(sz);
168 A.solve(v,d);
169 v.resize(v.capacity());
170 d.resize(d.capacity());
171 }
172
180 template<class S>
181 void setSubMatrix(const M& BCRS, S& rowset)
182 {
183 size_t sz = rowset.size();
184 A.resize(sz*n,sz*n);
185 typedef typename S::const_iterator SIter;
186 size_t r = 0, c = 0;
187 for(SIter rowIdx = rowset.begin(), rowEnd=rowset.end();
188 rowIdx!= rowEnd; ++rowIdx, r++)
189 {
190 c = 0;
191 for(SIter colIdx = rowset.begin(), colEnd=rowset.end();
192 colIdx!= colEnd; ++colIdx, c++)
193 {
194 if (BCRS[*rowIdx].find(*colIdx) == BCRS[*rowIdx].end())
195 continue;
196 for (size_t i=0; i<n; i++)
197 {
198 for (size_t j=0; j<n; j++)
199 {
200 A[r*n+i][c*n+j] = Impl::asMatrix(BCRS[*rowIdx][*colIdx])[i][j];
201 }
202 }
203 }
204 }
205 }
206 private:
207 DynamicMatrix<K> A;
208 };
209
210 template<typename T, bool tag>
211 class OverlappingAssignerHelper
212 {};
213
214 template<typename T>
215 using OverlappingAssigner = OverlappingAssignerHelper<T, Dune::StoresColumnCompressed<T>::value>;
216
217 // specialization for DynamicMatrix
218 template<class K, class Al, class X, class Y>
219 class OverlappingAssignerHelper< DynamicMatrixSubdomainSolver< BCRSMatrix<K, Al>, X, Y >,false>
220 {
221 public:
222 typedef BCRSMatrix< K, Al> matrix_type;
223 typedef typename X::field_type field_type;
224 typedef Y range_type;
225 typedef typename range_type::block_type block_type;
226 typedef typename matrix_type::size_type size_type;
227 static constexpr size_t n = std::decay_t<decltype(Impl::asMatrix(std::declval<K>()))>::rows;
235 OverlappingAssignerHelper(std::size_t maxlength, const BCRSMatrix<K, Al>& mat_, const X& b_, Y& x_);
236
240 inline
241 void deallocate();
242
246 inline
247 void resetIndexForNextDomain();
248
253 inline
254 DynamicVector<field_type> & lhs();
255
260 inline
261 DynamicVector<field_type> & rhs();
262
267 inline
268 void relaxResult(field_type relax);
269
274 void operator()(const size_type& domainIndex);
275
283 inline
284 void assignResult(block_type& res);
285
286 private:
290 const matrix_type* mat;
292 // we need a pointer, because we have to avoid deep copies
293 DynamicVector<field_type> * rhs_;
295 // we need a pointer, because we have to avoid deep copies
296 DynamicVector<field_type> * lhs_;
298 const range_type* b;
300 range_type* x;
302 std::size_t i;
304 std::size_t maxlength_;
305 };
306
307#if HAVE_SUPERLU || HAVE_SUITESPARSE_UMFPACK
308 template<template<class> class S, typename T, typename A>
309 struct OverlappingAssignerHelper<S<BCRSMatrix<T, A>>, true>
310 {
311 typedef BCRSMatrix<T, A> matrix_type;
312 typedef typename S<BCRSMatrix<T, A>>::range_type range_type;
313 typedef typename range_type::field_type field_type;
314 typedef typename range_type::block_type block_type;
315
316 typedef typename matrix_type::size_type size_type;
317
318 static constexpr size_t n = std::decay_t<decltype(Impl::asMatrix(std::declval<T>()))>::rows;
319 static constexpr size_t m = std::decay_t<decltype(Impl::asMatrix(std::declval<T>()))>::cols;
327 OverlappingAssignerHelper(std::size_t maxlength, const matrix_type& mat,
328 const range_type& b, range_type& x);
334 void deallocate();
335
336 /*
337 * @brief Resets the local index to zero.
338 */
339 void resetIndexForNextDomain();
340
345 field_type* lhs();
346
351 field_type* rhs();
352
357 void relaxResult(field_type relax);
358
363 void operator()(const size_type& domain);
364
372 void assignResult(block_type& res);
373
374 private:
378 const matrix_type* mat;
380 field_type* rhs_;
382 field_type* lhs_;
384 const range_type* b;
386 range_type* x;
388 std::size_t i;
390 std::size_t maxlength_;
391 };
392
393#endif // HAVE_SUPERLU || HAVE_SUITESPARSE_UMFPACK
394
395 template<class M, class X, class Y>
396 class OverlappingAssignerILUBase
397 {
398 public:
399 typedef M matrix_type;
400
401 typedef typename Y::field_type field_type;
402
403 typedef typename Y::block_type block_type;
404
405 typedef typename matrix_type::size_type size_type;
413 OverlappingAssignerILUBase(std::size_t maxlength, const M& mat,
414 const Y& b, X& x);
420 void deallocate();
421
425 void resetIndexForNextDomain();
426
431 X& lhs();
432
437 Y& rhs();
438
443 void relaxResult(field_type relax);
444
449 void operator()(const size_type& domain);
450
458 void assignResult(block_type& res);
459
460 private:
464 const M* mat;
466 X* lhs_;
468 Y* rhs_;
470 const Y* b;
472 X* x;
474 size_type i;
475 };
476
477 // specialization for ILU0
478 template<class M, class X, class Y>
479 class OverlappingAssignerHelper<ILU0SubdomainSolver<M,X,Y>, false>
480 : public OverlappingAssignerILUBase<M,X,Y>
481 {
482 public:
490 OverlappingAssignerHelper(std::size_t maxlength, const M& mat,
491 const Y& b, X& x)
492 : OverlappingAssignerILUBase<M,X,Y>(maxlength, mat,b,x)
493 {}
494 };
495
496 // specialization for ILUN
497 template<class M, class X, class Y>
498 class OverlappingAssignerHelper<ILUNSubdomainSolver<M,X,Y>,false>
499 : public OverlappingAssignerILUBase<M,X,Y>
500 {
501 public:
509 OverlappingAssignerHelper(std::size_t maxlength, const M& mat,
510 const Y& b, X& x)
511 : OverlappingAssignerILUBase<M,X,Y>(maxlength, mat,b,x)
512 {}
513 };
514
515 template<typename S, typename T>
516 struct AdditiveAdder
517 {};
518
519 template<typename S, typename T, typename A>
520 struct AdditiveAdder<S, BlockVector<T,A> >
521 {
522 typedef typename A::size_type size_type;
523 typedef typename std::decay_t<decltype(Impl::asVector(std::declval<T>()))>::field_type field_type;
524 AdditiveAdder(BlockVector<T,A>& v, BlockVector<T,A>& x,
525 OverlappingAssigner<S>& assigner, const field_type& relax_);
526 void operator()(const size_type& domain);
527 void axpy();
528 static constexpr size_t n = std::decay_t<decltype(Impl::asVector(std::declval<T>()))>::dimension;
529
530 private:
531 BlockVector<T,A>* v;
532 BlockVector<T,A>* x;
533 OverlappingAssigner<S>* assigner;
534 field_type relax;
535 };
536
537 template<typename S,typename T>
538 struct MultiplicativeAdder
539 {};
540
541 template<typename S, typename T, typename A>
542 struct MultiplicativeAdder<S, BlockVector<T,A> >
543 {
544 typedef typename A::size_type size_type;
545 typedef typename std::decay_t<decltype(Impl::asVector(std::declval<T>()))>::field_type field_type;
546 MultiplicativeAdder(BlockVector<T,A>& v, BlockVector<T,A>& x,
547 OverlappingAssigner<S>& assigner_, const field_type& relax_);
548 void operator()(const size_type& domain);
549 void axpy();
550 static constexpr size_t n = std::decay_t<decltype(Impl::asVector(std::declval<T>()))>::dimension;
551
552 private:
553 BlockVector<T,A>* x;
554 OverlappingAssigner<S>* assigner;
555 field_type relax;
556 };
557
567 template<typename T, class X, class S>
569 {};
570
571 template<class X, class S>
573 {
574 typedef AdditiveAdder<S,X> Adder;
575 };
576
577 template<class X, class S>
578 struct AdderSelector<MultiplicativeSchwarzMode,X,S>
579 {
580 typedef MultiplicativeAdder<S,X> Adder;
581 };
582
583 template<class X, class S>
584 struct AdderSelector<SymmetricMultiplicativeSchwarzMode,X,S>
585 {
586 typedef MultiplicativeAdder<S,X> Adder;
587 };
588
600 template<typename T1, typename T2, bool forward>
602 {
603 typedef T1 solver_vector;
604 typedef typename solver_vector::iterator solver_iterator;
605 typedef T2 subdomain_vector;
606 typedef typename subdomain_vector::const_iterator domain_iterator;
607
608 static solver_iterator begin(solver_vector& sv)
609 {
610 return sv.begin();
611 }
612
613 static solver_iterator end(solver_vector& sv)
614 {
615 return sv.end();
616 }
617 static domain_iterator begin(const subdomain_vector& sv)
618 {
619 return sv.begin();
620 }
621
622 static domain_iterator end(const subdomain_vector& sv)
623 {
624 return sv.end();
625 }
626 };
627
628 template<typename T1, typename T2>
629 struct IteratorDirectionSelector<T1,T2,false>
630 {
631 typedef T1 solver_vector;
632 typedef typename solver_vector::reverse_iterator solver_iterator;
633 typedef T2 subdomain_vector;
634 typedef typename subdomain_vector::const_reverse_iterator domain_iterator;
635
636 static solver_iterator begin(solver_vector& sv)
637 {
638 return sv.rbegin();
639 }
640
641 static solver_iterator end(solver_vector& sv)
642 {
643 return sv.rend();
644 }
645 static domain_iterator begin(const subdomain_vector& sv)
646 {
647 return sv.rbegin();
648 }
649
650 static domain_iterator end(const subdomain_vector& sv)
651 {
652 return sv.rend();
653 }
654 };
655
664 template<class T>
666 {
667 typedef T smoother;
668 typedef typename smoother::range_type range_type;
669
670 static void apply(smoother& sm, range_type& v, const range_type& b)
671 {
672 sm.template apply<true>(v, b);
673 }
674 };
675
676 template<class M, class X, class TD, class TA>
678 {
680 typedef typename smoother::range_type range_type;
681
682 static void apply(smoother& sm, range_type& v, const range_type& b)
683 {
684 sm.template apply<true>(v, b);
685 sm.template apply<false>(v, b);
686 }
687 };
688
689 template<class T, bool tag>
690 struct SeqOverlappingSchwarzAssemblerHelper
691 {};
692
693 template<class T>
694 using SeqOverlappingSchwarzAssembler = SeqOverlappingSchwarzAssemblerHelper<T,Dune::StoresColumnCompressed<T>::value>;
695
696 template<class K, class Al, class X, class Y>
697 struct SeqOverlappingSchwarzAssemblerHelper< DynamicMatrixSubdomainSolver< BCRSMatrix< K, Al>, X, Y >,false>
698 {
699 typedef BCRSMatrix< K, Al> matrix_type;
700 static constexpr size_t n = std::decay_t<decltype(Impl::asMatrix(std::declval<K>()))>::rows;
701 template<class RowToDomain, class Solvers, class SubDomains>
702 static std::size_t assembleLocalProblems(const RowToDomain& rowToDomain, const matrix_type& mat,
703 Solvers& solvers, const SubDomains& domains,
704 bool onTheFly);
705 };
706
707 template<template<class> class S, typename T, typename A>
708 struct SeqOverlappingSchwarzAssemblerHelper<S<BCRSMatrix<T,A>>,true>
709 {
710 typedef BCRSMatrix<T,A> matrix_type;
711 static constexpr size_t n = std::decay_t<decltype(Impl::asMatrix(std::declval<T>()))>::rows;
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 SeqOverlappingSchwarzAssemblerILUBase
720 {
721 typedef M matrix_type;
722 template<class RowToDomain, class Solvers, class SubDomains>
723 static std::size_t assembleLocalProblems(const RowToDomain& rowToDomain, const matrix_type& mat,
724 Solvers& solvers, const SubDomains& domains,
725 bool onTheFly);
726 };
727
728 template<class M,class X, class Y>
729 struct SeqOverlappingSchwarzAssemblerHelper<ILU0SubdomainSolver<M,X,Y>,false>
730 : public SeqOverlappingSchwarzAssemblerILUBase<M,X,Y>
731 {};
732
733 template<class M,class X, class Y>
734 struct SeqOverlappingSchwarzAssemblerHelper<ILUNSubdomainSolver<M,X,Y>,false>
735 : public SeqOverlappingSchwarzAssemblerILUBase<M,X,Y>
736 {};
737
748 template<class M, class X, class TM=AdditiveSchwarzMode,
749 class TD=ILU0SubdomainSolver<M,X,X>, class TA=std::allocator<X> >
751 : public Preconditioner<X,X>
752 {
753 public:
757 typedef M matrix_type;
758
762 typedef X domain_type;
763
767 typedef X range_type;
768
775 typedef TM Mode;
776
780 typedef typename X::field_type field_type;
781
783 typedef typename matrix_type::size_type size_type;
784
786 typedef TA allocator;
787
789 typedef std::set<size_type, std::less<size_type>,
790 typename std::allocator_traits<TA>::template rebind_alloc<size_type> >
792
794 typedef std::vector<subdomain_type, typename std::allocator_traits<TA>::template rebind_alloc<subdomain_type> > subdomain_vector;
795
798
800 typedef std::vector<subdomain_list, typename std::allocator_traits<TA>::template rebind_alloc<subdomain_list> > rowtodomain_vector;
801
803 typedef TD slu;
804
806 typedef std::vector<slu, typename std::allocator_traits<TA>::template rebind_alloc<slu> > slu_vector;
807
822 field_type relaxationFactor=1, bool onTheFly_=true);
823
836 field_type relaxationFactor=1, bool onTheFly_=true);
837
843 virtual void pre (X& x, X& b)
844 {
847 }
848
854 virtual void apply (X& v, const X& d);
855
861 virtual void post (X& x)
862 {
864 }
865
866 template<bool forward>
867 void apply(X& v, const X& d);
868
871 {
873 }
874
875 private:
876 const M& mat;
877 slu_vector solvers;
878 subdomain_vector subDomains;
879 field_type relax;
880
881 typename M::size_type maxlength;
882
883 bool onTheFly;
884 };
885
886
887
888 template<class I, class S, class D>
889 OverlappingSchwarzInitializer<I,S,D>::OverlappingSchwarzInitializer(InitializerList& il,
890 const IndexSet& idx,
891 const subdomain_vector& domains_)
892 : initializers(&il), indices(&idx), indexMaps(il.size()), domains(domains_)
893 {}
894
895
896 template<class I, class S, class D>
897 void OverlappingSchwarzInitializer<I,S,D>::addRowNnz(const Iter& row)
898 {
899 typedef typename IndexSet::value_type::const_iterator iterator;
900 for(iterator domain=(*indices)[row.index()].begin(); domain != (*indices)[row.index()].end(); ++domain) {
901 (*initializers)[*domain].addRowNnz(row, domains[*domain]);
902 indexMaps[*domain].insert(row.index());
903 }
904 }
905
906 template<class I, class S, class D>
907 void OverlappingSchwarzInitializer<I,S,D>::allocate()
908 {
909 for(auto&& i: *initializers)
910 i.allocateMatrixStorage();
911 for(auto&& i: *initializers)
912 i.allocateMarker();
913 }
914
915 template<class I, class S, class D>
916 void OverlappingSchwarzInitializer<I,S,D>::countEntries(const Iter& row, const CIter& col) const
917 {
918 typedef typename IndexSet::value_type::const_iterator iterator;
919 for(iterator domain=(*indices)[row.index()].begin(); domain != (*indices)[row.index()].end(); ++domain) {
920 typename std::map<size_type,size_type>::const_iterator v = indexMaps[*domain].find(col.index());
921 if(v!= indexMaps[*domain].end()) {
922 (*initializers)[*domain].countEntries(indexMaps[*domain].find(col.index())->second);
923 }
924 }
925 }
926
927 template<class I, class S, class D>
928 void OverlappingSchwarzInitializer<I,S,D>::calcColstart() const
929 {
930 for(auto&& i : *initializers)
931 i.calcColstart();
932 }
933
934 template<class I, class S, class D>
935 void OverlappingSchwarzInitializer<I,S,D>::copyValue(const Iter& row, const CIter& col) const
936 {
937 typedef typename IndexSet::value_type::const_iterator iterator;
938 for(iterator domain=(*indices)[row.index()].begin(); domain!= (*indices)[row.index()].end(); ++domain) {
939 typename std::map<size_type,size_type>::const_iterator v = indexMaps[*domain].find(col.index());
940 if(v!= indexMaps[*domain].end()) {
941 assert(indexMaps[*domain].end()!=indexMaps[*domain].find(row.index()));
942 (*initializers)[*domain].copyValue(col, indexMaps[*domain].find(row.index())->second,
943 v->second);
944 }
945 }
946 }
947
948 template<class I, class S, class D>
949 void OverlappingSchwarzInitializer<I,S,D>::createMatrix() const
950 {
951 std::vector<IndexMap>().swap(indexMaps);
952 for(auto&& i: *initializers)
953 i.createMatrix();
954 }
955
956 template<class I, class S, class D>
957 OverlappingSchwarzInitializer<I,S,D>::IndexMap::IndexMap()
958 : row(0)
959 {}
960
961 template<class I, class S, class D>
962 void OverlappingSchwarzInitializer<I,S,D>::IndexMap::insert(size_type grow)
963 {
964 assert(map_.find(grow)==map_.end());
965 map_.insert(std::make_pair(grow, row++));
966 }
967
968 template<class I, class S, class D>
969 typename OverlappingSchwarzInitializer<I,S,D>::IndexMap::const_iterator
970 OverlappingSchwarzInitializer<I,S,D>::IndexMap::find(size_type grow) const
971 {
972 return map_.find(grow);
973 }
974
975 template<class I, class S, class D>
976 typename OverlappingSchwarzInitializer<I,S,D>::IndexMap::iterator
977 OverlappingSchwarzInitializer<I,S,D>::IndexMap::find(size_type grow)
978 {
979 return map_.find(grow);
980 }
981
982 template<class I, class S, class D>
983 typename OverlappingSchwarzInitializer<I,S,D>::IndexMap::const_iterator
984 OverlappingSchwarzInitializer<I,S,D>::IndexMap::end() const
985 {
986 return map_.end();
987 }
988
989 template<class I, class S, class D>
990 typename OverlappingSchwarzInitializer<I,S,D>::IndexMap::iterator
991 OverlappingSchwarzInitializer<I,S,D>::IndexMap::end()
992 {
993 return map_.end();
994 }
995
996 template<class I, class S, class D>
997 typename OverlappingSchwarzInitializer<I,S,D>::IndexMap::const_iterator
998 OverlappingSchwarzInitializer<I,S,D>::IndexMap::begin() const
999 {
1000 return map_.begin();
1001 }
1002
1003 template<class I, class S, class D>
1004 typename OverlappingSchwarzInitializer<I,S,D>::IndexMap::iterator
1005 OverlappingSchwarzInitializer<I,S,D>::IndexMap::begin()
1006 {
1007 return map_.begin();
1008 }
1009
1010 template<class M, class X, class TM, class TD, class TA>
1012 field_type relaxationFactor, bool fly)
1013 : mat(mat_), relax(relaxationFactor), onTheFly(fly)
1014 {
1015 typedef typename rowtodomain_vector::const_iterator RowDomainIterator;
1016 typedef typename subdomain_list::const_iterator DomainIterator;
1017#ifdef DUNE_ISTL_WITH_CHECKING
1018 assert(rowToDomain.size()==mat.N());
1019 assert(rowToDomain.size()==mat.M());
1020
1021 for(RowDomainIterator iter=rowToDomain.begin(); iter != rowToDomain.end(); ++iter)
1022 assert(iter->size()>0);
1023
1024#endif
1025 // calculate the number of domains
1026 size_type domains=0;
1027 for(RowDomainIterator iter=rowToDomain.begin(); iter != rowToDomain.end(); ++iter)
1028 for(DomainIterator d=iter->begin(); d != iter->end(); ++d)
1029 domains=std::max(domains, *d);
1030 ++domains;
1031
1032 solvers.resize(domains);
1033 subDomains.resize(domains);
1034
1035 // initialize subdomains to row mapping from row to subdomain mapping
1036 size_type row=0;
1037 for(RowDomainIterator iter=rowToDomain.begin(); iter != rowToDomain.end(); ++iter, ++row)
1038 for(DomainIterator d=iter->begin(); d != iter->end(); ++d)
1039 subDomains[*d].insert(row);
1040
1041#ifdef DUNE_ISTL_WITH_CHECKING
1042 size_type i=0;
1043 typedef typename subdomain_vector::const_iterator iterator;
1044 for(iterator iter=subDomains.begin(); iter != subDomains.end(); ++iter) {
1045 typedef typename subdomain_type::const_iterator entry_iterator;
1046 Dune::dvverb<<"domain "<<i++<<":";
1047 for(entry_iterator entry = iter->begin(); entry != iter->end(); ++entry) {
1048 Dune::dvverb<<" "<<*entry;
1049 }
1050 Dune::dvverb<<std::endl;
1051 }
1052#endif
1053 maxlength = SeqOverlappingSchwarzAssembler<slu>
1054 ::assembleLocalProblems(rowToDomain, mat, solvers, subDomains, onTheFly);
1055 }
1056
1057 template<class M, class X, class TM, class TD, class TA>
1059 const subdomain_vector& sd,
1060 field_type relaxationFactor,
1061 bool fly)
1062 : mat(mat_), solvers(sd.size()), subDomains(sd), relax(relaxationFactor),
1063 onTheFly(fly)
1064 {
1065 typedef typename subdomain_vector::const_iterator DomainIterator;
1066
1067#ifdef DUNE_ISTL_WITH_CHECKING
1068 size_type i=0;
1069
1070 for(DomainIterator d=sd.begin(); d != sd.end(); ++d,++i) {
1071 //std::cout<<i<<": "<<d->size()<<std::endl;
1072 assert(d->size()>0);
1073 typedef typename DomainIterator::value_type::const_iterator entry_iterator;
1074 Dune::dvverb<<"domain "<<i<<":";
1075 for(entry_iterator entry = d->begin(); entry != d->end(); ++entry) {
1076 Dune::dvverb<<" "<<*entry;
1077 }
1078 Dune::dvverb<<std::endl;
1079 }
1080
1081#endif
1082
1083 // Create a row to subdomain mapping
1084 rowtodomain_vector rowToDomain(mat.N());
1085
1086 size_type domainId=0;
1087
1088 for(DomainIterator domain=sd.begin(); domain != sd.end(); ++domain, ++domainId) {
1089 typedef typename subdomain_type::const_iterator iterator;
1090 for(iterator row=domain->begin(); row != domain->end(); ++row)
1091 rowToDomain[*row].push_back(domainId);
1092 }
1093
1094 maxlength = SeqOverlappingSchwarzAssembler<slu>
1095 ::assembleLocalProblems(rowToDomain, mat, solvers, subDomains, onTheFly);
1096 }
1097
1104 template<class M>
1106
1107 template<typename T, typename A>
1109 {
1110 static constexpr size_t n = std::decay_t<decltype(Impl::asMatrix(std::declval<T>()))>::rows;
1111 static constexpr size_t m = std::decay_t<decltype(Impl::asMatrix(std::declval<T>()))>::cols;
1112 template<class Domain>
1113 static int size(const Domain & d)
1114 {
1115 assert(n==m);
1116 return m*d.size();
1117 }
1118 };
1119
1120 template<class K, class Al, class X, class Y>
1121 template<class RowToDomain, class Solvers, class SubDomains>
1122 std::size_t
1123 SeqOverlappingSchwarzAssemblerHelper< DynamicMatrixSubdomainSolver< BCRSMatrix< K, Al>, X, Y >,false>::
1124 assembleLocalProblems(const RowToDomain& rowToDomain,
1125 const matrix_type& mat,
1126 Solvers& solvers,
1127 const SubDomains& subDomains,
1128 bool onTheFly)
1129 {
1130 DUNE_UNUSED_PARAMETER(onTheFly);
1131 DUNE_UNUSED_PARAMETER(rowToDomain);
1133 DUNE_UNUSED_PARAMETER(solvers);
1134 typedef typename SubDomains::const_iterator DomainIterator;
1135 std::size_t maxlength = 0;
1136
1137 assert(onTheFly);
1138
1139 for(DomainIterator domain=subDomains.begin(); domain!=subDomains.end(); ++domain)
1140 maxlength=std::max(maxlength, domain->size());
1141 maxlength*=n;
1142
1143 return maxlength;
1144 }
1145
1146#if HAVE_SUPERLU || HAVE_SUITESPARSE_UMFPACK
1147 template<template<class> class S, typename T, typename A>
1148 template<class RowToDomain, class Solvers, class SubDomains>
1149 std::size_t SeqOverlappingSchwarzAssemblerHelper<S<BCRSMatrix<T,A>>,true>::assembleLocalProblems(const RowToDomain& rowToDomain,
1150 const matrix_type& mat,
1151 Solvers& solvers,
1152 const SubDomains& subDomains,
1153 bool onTheFly)
1154 {
1155 typedef typename S<BCRSMatrix<T,A>>::MatrixInitializer MatrixInitializer;
1156 typedef typename std::vector<MatrixInitializer>::iterator InitializerIterator;
1157 typedef typename SubDomains::const_iterator DomainIterator;
1158 typedef typename Solvers::iterator SolverIterator;
1159 std::size_t maxlength = 0;
1160
1161 if(onTheFly) {
1162 for(DomainIterator domain=subDomains.begin(); domain!=subDomains.end(); ++domain)
1163 maxlength=std::max(maxlength, domain->size());
1164 maxlength*=Impl::asMatrix(*mat[0].begin()).N();
1165 }else{
1166 // initialize the initializers
1167 DomainIterator domain=subDomains.begin();
1168
1169 // Create the initializers list.
1170 std::vector<MatrixInitializer> initializers(subDomains.size());
1171
1172 SolverIterator solver=solvers.begin();
1173 for(InitializerIterator initializer=initializers.begin(); initializer!=initializers.end();
1174 ++initializer, ++solver, ++domain) {
1175 solver->getInternalMatrix().N_=SeqOverlappingSchwarzDomainSize<matrix_type>::size(*domain);
1176 solver->getInternalMatrix().M_=SeqOverlappingSchwarzDomainSize<matrix_type>::size(*domain);
1177 //solver->setVerbosity(true);
1178 *initializer=MatrixInitializer(solver->getInternalMatrix());
1179 }
1180
1181 // Set up the supermatrices according to the subdomains
1182 typedef OverlappingSchwarzInitializer<std::vector<MatrixInitializer>,
1183 RowToDomain, SubDomains> Initializer;
1184
1185 Initializer initializer(initializers, rowToDomain, subDomains);
1186 copyToColCompMatrix(initializer, mat);
1187
1188 // Calculate the LU decompositions
1189 for(auto&& s: solvers)
1190 s.decompose();
1191 for (SolverIterator solverIt = solvers.begin(); solverIt != solvers.end(); ++solverIt)
1192 {
1193 assert(solverIt->getInternalMatrix().N() == solverIt->getInternalMatrix().M());
1194 maxlength = std::max(maxlength, solverIt->getInternalMatrix().N());
1195 }
1196 }
1197 return maxlength;
1198 }
1199
1200#endif // HAVE_SUPERLU || HAVE_SUITESPARSE_UMFPACK
1201
1202 template<class M,class X,class Y>
1203 template<class RowToDomain, class Solvers, class SubDomains>
1204 std::size_t SeqOverlappingSchwarzAssemblerILUBase<M,X,Y>::assembleLocalProblems(const RowToDomain& rowToDomain,
1205 const matrix_type& mat,
1206 Solvers& solvers,
1207 const SubDomains& subDomains,
1208 bool onTheFly)
1209 {
1210 DUNE_UNUSED_PARAMETER(rowToDomain);
1211 typedef typename SubDomains::const_iterator DomainIterator;
1212 typedef typename Solvers::iterator SolverIterator;
1213 std::size_t maxlength = 0;
1214
1215 if(onTheFly) {
1216 for(DomainIterator domain=subDomains.begin(); domain!=subDomains.end(); ++domain)
1217 maxlength=std::max(maxlength, domain->size());
1218 }else{
1219 // initialize the solvers of the local prolems.
1220 SolverIterator solver=solvers.begin();
1221 for(DomainIterator domain=subDomains.begin(); domain!=subDomains.end();
1222 ++domain, ++solver) {
1223 solver->setSubMatrix(mat, *domain);
1224 maxlength=std::max(maxlength, domain->size());
1225 }
1226 }
1227
1228 return maxlength;
1229
1230 }
1231
1232
1233 template<class M, class X, class TM, class TD, class TA>
1235 {
1237 }
1238
1239 template<class M, class X, class TM, class TD, class TA>
1240 template<bool forward>
1241 void SeqOverlappingSchwarz<M,X,TM,TD,TA>::apply(X& x, const X& b)
1242 {
1243 typedef slu_vector solver_vector;
1244 typedef typename IteratorDirectionSelector<solver_vector,subdomain_vector,forward>::solver_iterator iterator;
1245 typedef typename IteratorDirectionSelector<solver_vector,subdomain_vector,forward>::domain_iterator
1246 domain_iterator;
1247
1248 OverlappingAssigner<TD> assigner(maxlength, mat, b, x);
1249
1252 X v(x); // temporary for the update
1253 v=0;
1254
1255 typedef typename AdderSelector<TM,X,TD >::Adder Adder;
1256 Adder adder(v, x, assigner, relax);
1257
1258 for(; domain != IteratorDirectionSelector<solver_vector,subdomain_vector,forward>::end(subDomains); ++domain) {
1259 //Copy rhs to C-array for SuperLU
1260 std::for_each(domain->begin(), domain->end(), assigner);
1261 assigner.resetIndexForNextDomain();
1262 if(onTheFly) {
1263 // Create the subdomain solver
1264 slu sdsolver;
1265 sdsolver.setSubMatrix(mat, *domain);
1266 // Apply
1267 sdsolver.apply(assigner.lhs(), assigner.rhs());
1268 }else{
1269 solver->apply(assigner.lhs(), assigner.rhs());
1270 ++solver;
1271 }
1272
1273 //Add relaxed correction to from SuperLU to v
1274 std::for_each(domain->begin(), domain->end(), adder);
1275 assigner.resetIndexForNextDomain();
1276
1277 }
1278
1279 adder.axpy();
1280 assigner.deallocate();
1281 }
1282
1283 template<class K, class Al, class X, class Y>
1284 OverlappingAssignerHelper< DynamicMatrixSubdomainSolver< BCRSMatrix< K, Al>, X, Y >,false>
1285 ::OverlappingAssignerHelper(std::size_t maxlength, const BCRSMatrix<K, Al>& mat_,
1286 const X& b_, Y& x_) :
1287 mat(&mat_),
1288 rhs_( new DynamicVector<field_type>(maxlength, 42) ),
1289 lhs_( new DynamicVector<field_type>(maxlength, -42) ),
1290 b(&b_),
1291 x(&x_),
1292 i(0),
1293 maxlength_(maxlength)
1294 {}
1295
1296 template<class K, class Al, class X, class Y>
1297 void
1298 OverlappingAssignerHelper< DynamicMatrixSubdomainSolver< BCRSMatrix< K, Al>, X, Y >,false>
1299 ::deallocate()
1300 {
1301 delete rhs_;
1302 delete lhs_;
1303 }
1304
1305 template<class K, class Al, class X, class Y>
1306 void
1307 OverlappingAssignerHelper< DynamicMatrixSubdomainSolver< BCRSMatrix< K, Al>, X, Y >,false>
1308 ::resetIndexForNextDomain()
1309 {
1310 i=0;
1311 }
1312
1313 template<class K, class Al, class X, class Y>
1315 OverlappingAssignerHelper< DynamicMatrixSubdomainSolver< BCRSMatrix< K, Al>, X, Y >,false>
1316 ::lhs()
1317 {
1318 return *lhs_;
1319 }
1320
1321 template<class K, class Al, class X, class Y>
1323 OverlappingAssignerHelper< DynamicMatrixSubdomainSolver< BCRSMatrix< K, Al>, X, Y >,false>
1324 ::rhs()
1325 {
1326 return *rhs_;
1327 }
1328
1329 template<class K, class Al, class X, class Y>
1330 void
1331 OverlappingAssignerHelper< DynamicMatrixSubdomainSolver< BCRSMatrix< K, Al>, X, Y >,false>
1332 ::relaxResult(field_type relax)
1333 {
1334 lhs() *= relax;
1335 }
1336
1337 template<class K, class Al, class X, class Y>
1338 void
1339 OverlappingAssignerHelper< DynamicMatrixSubdomainSolver< BCRSMatrix< K, Al>, X, Y >,false>
1340 ::operator()(const size_type& domainIndex)
1341 {
1342 lhs() = 0.0;
1343#if 0
1344 //assign right hand side of current domainindex block
1345 for(size_type j=0; j<n; ++j, ++i) {
1346 assert(i<maxlength_);
1347 rhs()[i]=(*b)[domainIndex][j];
1348 }
1349
1350 // loop over all Matrix row entries and calculate defect.
1351 typedef typename matrix_type::ConstColIterator col_iterator;
1352
1353 // calculate defect for current row index block
1354 for(col_iterator col=(*mat)[domainIndex].begin(); col!=(*mat)[domainIndex].end(); ++col) {
1355 block_type tmp(0.0);
1356 (*col).mv((*x)[col.index()], tmp);
1357 i-=n;
1358 for(size_type j=0; j<n; ++j, ++i) {
1359 assert(i<maxlength_);
1360 rhs()[i]-=tmp[j];
1361 }
1362 }
1363#else
1364 //assign right hand side of current domainindex block
1365 for(size_type j=0; j<n; ++j, ++i) {
1366 assert(i<maxlength_);
1367 rhs()[i]=Impl::asVector((*b)[domainIndex])[j];
1368
1369 // loop over all Matrix row entries and calculate defect.
1370 typedef typename matrix_type::ConstColIterator col_iterator;
1371
1372 // calculate defect for current row index block
1373 for(col_iterator col=(*mat)[domainIndex].begin(); col!=(*mat)[domainIndex].end(); ++col) {
1374 for(size_type k=0; k<n; ++k) {
1375 rhs()[i]-=Impl::asMatrix(*col)[j][k] * Impl::asVector((*x)[col.index()])[k];
1376 }
1377 }
1378 }
1379#endif
1380 }
1381
1382 template<class K, class Al, class X, class Y>
1383 void
1384 OverlappingAssignerHelper< DynamicMatrixSubdomainSolver< BCRSMatrix< K, Al>, X, Y >,false>
1385 ::assignResult(block_type& res)
1386 {
1387 // assign the result of the local solve to the global vector
1388 for(size_type j=0; j<n; ++j, ++i) {
1389 assert(i<maxlength_);
1390 Impl::asVector(res)[j]+=lhs()[i];
1391 }
1392 }
1393
1394#if HAVE_SUPERLU || HAVE_SUITESPARSE_UMFPACK
1395
1396 template<template<class> class S, typename T, typename A>
1397 OverlappingAssignerHelper<S<BCRSMatrix<T,A>>,true>
1398 ::OverlappingAssignerHelper(std::size_t maxlength,
1399 const BCRSMatrix<T,A>& mat_,
1400 const range_type& b_,
1401 range_type& x_)
1402 : mat(&mat_),
1403 b(&b_),
1404 x(&x_), i(0), maxlength_(maxlength)
1405 {
1406 rhs_ = new field_type[maxlength];
1407 lhs_ = new field_type[maxlength];
1408
1409 }
1410
1411 template<template<class> class S, typename T, typename A>
1412 void OverlappingAssignerHelper<S<BCRSMatrix<T,A> >,true>::deallocate()
1413 {
1414 delete[] rhs_;
1415 delete[] lhs_;
1416 }
1417
1418 template<template<class> class S, typename T, typename A>
1419 void OverlappingAssignerHelper<S<BCRSMatrix<T,A>>,true>::operator()(const size_type& domainIndex)
1420 {
1421 //assign right hand side of current domainindex block
1422 // rhs is an array of doubles!
1423 // rhs[starti] = b[domainindex]
1424 for(size_type j=0; j<n; ++j, ++i) {
1425 assert(i<maxlength_);
1426 rhs_[i]=Impl::asVector((*b)[domainIndex])[j];
1427 }
1428
1429
1430 // loop over all Matrix row entries and calculate defect.
1431 typedef typename matrix_type::ConstColIterator col_iterator;
1432
1433 // calculate defect for current row index block
1434 for(col_iterator col=(*mat)[domainIndex].begin(); col!=(*mat)[domainIndex].end(); ++col) {
1435 block_type tmp;
1436 Impl::asMatrix(*col).mv((*x)[col.index()], tmp);
1437 i-=n;
1438 for(size_type j=0; j<n; ++j, ++i) {
1439 assert(i<maxlength_);
1440 rhs_[i]-=Impl::asVector(tmp)[j];
1441 }
1442
1443 }
1444
1445 }
1446
1447 template<template<class> class S, typename T, typename A>
1448 void OverlappingAssignerHelper<S<BCRSMatrix<T,A>>,true>::relaxResult(field_type relax)
1449 {
1450 for(size_type j=i+n; i<j; ++i) {
1451 assert(i<maxlength_);
1452 lhs_[i]*=relax;
1453 }
1454 i-=n;
1455 }
1456
1457 template<template<class> class S, typename T, typename A>
1458 void OverlappingAssignerHelper<S<BCRSMatrix<T,A>>,true>::assignResult(block_type& res)
1459 {
1460 // assign the result of the local solve to the global vector
1461 for(size_type j=0; j<n; ++j, ++i) {
1462 assert(i<maxlength_);
1463 Impl::asVector(res)[j]+=lhs_[i];
1464 }
1465 }
1466
1467 template<template<class> class S, typename T, typename A>
1468 void OverlappingAssignerHelper<S<BCRSMatrix<T,A>>,true>::resetIndexForNextDomain()
1469 {
1470 i=0;
1471 }
1472
1473 template<template<class> class S, typename T, typename A>
1474 typename OverlappingAssignerHelper<S<BCRSMatrix<T,A>>,true>::field_type*
1475 OverlappingAssignerHelper<S<BCRSMatrix<T,A>>,true>::lhs()
1476 {
1477 return lhs_;
1478 }
1479
1480 template<template<class> class S, typename T, typename A>
1481 typename OverlappingAssignerHelper<S<BCRSMatrix<T,A>>,true>::field_type*
1482 OverlappingAssignerHelper<S<BCRSMatrix<T,A>>,true>::rhs()
1483 {
1484 return rhs_;
1485 }
1486
1487#endif // HAVE_SUPERLU || HAVE_SUITESPARSE_UMFPACK
1488
1489 template<class M, class X, class Y>
1490 OverlappingAssignerILUBase<M,X,Y>::OverlappingAssignerILUBase(std::size_t maxlength,
1491 const M& mat_,
1492 const Y& b_,
1493 X& x_)
1494 : mat(&mat_),
1495 b(&b_),
1496 x(&x_), i(0)
1497 {
1498 rhs_= new Y(maxlength);
1499 lhs_ = new X(maxlength);
1500 }
1501
1502 template<class M, class X, class Y>
1503 void OverlappingAssignerILUBase<M,X,Y>::deallocate()
1504 {
1505 delete rhs_;
1506 delete lhs_;
1507 }
1508
1509 template<class M, class X, class Y>
1510 void OverlappingAssignerILUBase<M,X,Y>::operator()(const size_type& domainIndex)
1511 {
1512 (*rhs_)[i]=(*b)[domainIndex];
1513
1514 // loop over all Matrix row entries and calculate defect.
1515 typedef typename matrix_type::ConstColIterator col_iterator;
1516
1517 // calculate defect for current row index block
1518 for(col_iterator col=(*mat)[domainIndex].begin(); col!=(*mat)[domainIndex].end(); ++col) {
1519 Impl::asMatrix(*col).mmv((*x)[col.index()], (*rhs_)[i]);
1520 }
1521 // Goto next local index
1522 ++i;
1523 }
1524
1525 template<class M, class X, class Y>
1526 void OverlappingAssignerILUBase<M,X,Y>::relaxResult(field_type relax)
1527 {
1528 (*lhs_)[i]*=relax;
1529 }
1530
1531 template<class M, class X, class Y>
1532 void OverlappingAssignerILUBase<M,X,Y>::assignResult(block_type& res)
1533 {
1534 res+=(*lhs_)[i++];
1535 }
1536
1537 template<class M, class X, class Y>
1538 X& OverlappingAssignerILUBase<M,X,Y>::lhs()
1539 {
1540 return *lhs_;
1541 }
1542
1543 template<class M, class X, class Y>
1544 Y& OverlappingAssignerILUBase<M,X,Y>::rhs()
1545 {
1546 return *rhs_;
1547 }
1548
1549 template<class M, class X, class Y>
1550 void OverlappingAssignerILUBase<M,X,Y>::resetIndexForNextDomain()
1551 {
1552 i=0;
1553 }
1554
1555 template<typename S, typename T, typename A>
1556 AdditiveAdder<S,BlockVector<T,A> >::AdditiveAdder(BlockVector<T,A>& v_,
1557 BlockVector<T,A>& x_,
1558 OverlappingAssigner<S>& assigner_,
1559 const field_type& relax_)
1560 : v(&v_), x(&x_), assigner(&assigner_), relax(relax_)
1561 {}
1562
1563 template<typename S, typename T, typename A>
1564 void AdditiveAdder<S,BlockVector<T,A> >::operator()(const size_type& domainIndex)
1565 {
1566 // add the result of the local solve to the current update
1567 assigner->assignResult((*v)[domainIndex]);
1568 }
1569
1570
1571 template<typename S, typename T, typename A>
1572 void AdditiveAdder<S,BlockVector<T,A> >::axpy()
1573 {
1574 // relax the update and add it to the current guess.
1575 x->axpy(relax,*v);
1576 }
1577
1578
1579 template<typename S, typename T, typename A>
1580 MultiplicativeAdder<S,BlockVector<T,A> >
1581 ::MultiplicativeAdder(BlockVector<T,A>& v_,
1582 BlockVector<T,A>& x_,
1583 OverlappingAssigner<S>& assigner_, const field_type& relax_)
1584 : x(&x_), assigner(&assigner_), relax(relax_)
1585 {
1587 }
1588
1589
1590 template<typename S,typename T, typename A>
1591 void MultiplicativeAdder<S,BlockVector<T,A> >::operator()(const size_type& domainIndex)
1592 {
1593 // add the result of the local solve to the current guess
1594 assigner->relaxResult(relax);
1595 assigner->assignResult((*x)[domainIndex]);
1596 }
1597
1598
1599 template<typename S,typename T, typename A>
1600 void MultiplicativeAdder<S,BlockVector<T,A> >::axpy()
1601 {
1602 // nothing to do, as the corrections already relaxed and added in operator()
1603 }
1604
1605
1607}
1608
1609#endif
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:426
row_type::ConstIterator ConstColIterator
Const iterator to the entries of a row.
Definition: bcrsmatrix.hh:700
Iterator begin()
begin iterator
Definition: densevector.hh:348
Iterator end()
end iterator
Definition: densevector.hh:354
Iterator find(size_type i)
return iterator to given element or end()
Definition: densevector.hh:374
Exact subdomain solver using Dune::DynamicMatrix<T>::solve.
Definition: overlappingschwarz.hh:137
size_type capacity() const
Number of elements for which memory has been allocated.
Definition: dynvector.hh:135
Index Set Interface base class.
Definition: indexidset.hh:76
Initializer for SuperLU Matrices representing the subdomains.
Definition: overlappingschwarz.hh:44
D subdomain_vector
The vector type containing the subdomain to row index mapping.
Definition: overlappingschwarz.hh:47
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:752
X::field_type field_type
The field type of the preconditioner.
Definition: overlappingschwarz.hh:780
void apply(X &v, const X &d)
Apply one step of the preconditioner to the system A(v)=d.
SLList< size_type, typename std::allocator_traits< TA >::template rebind_alloc< size_type > > subdomain_list
The type for the row to subdomain mapping.
Definition: overlappingschwarz.hh:797
std::vector< slu, typename std::allocator_traits< TA >::template rebind_alloc< slu > > slu_vector
The vector type containing subdomain solvers.
Definition: overlappingschwarz.hh:806
M matrix_type
The type of the matrix to precondition.
Definition: overlappingschwarz.hh:757
TM Mode
The mode (additive or multiplicative) of the Schwarz method.
Definition: overlappingschwarz.hh:775
X range_type
The range type of the preconditioner.
Definition: overlappingschwarz.hh:767
std::set< size_type, std::less< size_type >, typename std::allocator_traits< TA >::template rebind_alloc< size_type > > subdomain_type
The type for the subdomain to row index mapping.
Definition: overlappingschwarz.hh:791
X domain_type
The domain type of the preconditioner.
Definition: overlappingschwarz.hh:762
TD slu
The type for the subdomain solver in use.
Definition: overlappingschwarz.hh:803
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: overlappingschwarz.hh:870
virtual void pre(X &x, X &b)
Prepare the preconditioner.
Definition: overlappingschwarz.hh:843
virtual void post(X &x)
Postprocess the preconditioner.
Definition: overlappingschwarz.hh:861
TA allocator
The allocator to use.
Definition: overlappingschwarz.hh:786
std::vector< subdomain_type, typename std::allocator_traits< TA >::template rebind_alloc< subdomain_type > > subdomain_vector
The vector type containing the subdomain to row index mapping.
Definition: overlappingschwarz.hh:794
std::vector< subdomain_list, typename std::allocator_traits< TA >::template rebind_alloc< subdomain_list > > rowtodomain_vector
The vector type containing the row index to subdomain mapping.
Definition: overlappingschwarz.hh:800
matrix_type::size_type size_type
The return type of the size method.
Definition: overlappingschwarz.hh:783
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:46
#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:1234
SeqOverlappingSchwarz(const matrix_type &mat, const subdomain_vector &subDomains, field_type relaxationFactor=1, bool onTheFly_=true)
Construct the overlapping Schwarz method.
Definition: overlappingschwarz.hh:1058
SeqOverlappingSchwarz(const matrix_type &mat, const rowtodomain_vector &rowToDomain, field_type relaxationFactor=1, bool onTheFly_=true)
Definition: overlappingschwarz.hh:1011
auto max(ADLTag< 0 >, const V &v1, const V &v2)
implements binary Simd::max()
Definition: defaults.hh:79
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:14
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:569
Tag that the tells the Schwarz method to be additive.
Definition: overlappingschwarz.hh:117
Helper template meta program for application of overlapping Schwarz.
Definition: overlappingschwarz.hh:602
Tag that tells the Schwarz method to be multiplicative.
Definition: overlappingschwarz.hh:123
Helper template meta program for application of overlapping Schwarz.
Definition: overlappingschwarz.hh:666
Definition: overlappingschwarz.hh:1105
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:130
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.111.3 (Jul 15, 22:36, 2024)