5#ifndef DUNE_ISTL_OVERLAPPINGSCHWARZ_HH
6#define DUNE_ISTL_OVERLAPPINGSCHWARZ_HH
13#include <dune/common/dynmatrix.hh>
14#include <dune/common/sllist.hh>
39 template<
class M,
class X,
class TM,
class TD,
class TA>
45 template<
class I,
class S,
class D>
54 typedef typename AtomInitializer::Matrix
Matrix;
55 typedef typename Matrix::const_iterator
Iter;
56 typedef typename Matrix::row_type::const_iterator
CIter;
82 typedef std::map<size_type,size_type> Map;
83 typedef typename Map::iterator iterator;
84 typedef typename Map::const_iterator const_iterator;
96 const_iterator
begin()
const;
100 const_iterator
end()
const;
103 std::map<size_type,size_type> map_;
108 typedef typename InitializerList::iterator InitIterator;
109 typedef typename IndexSet::const_iterator IndexIteratur;
112 mutable std::vector<IndexMap> indexMaps;
139 template<
class M,
class X,
class Y>
143 template<
class K,
class Al,
class X,
class Y>
156 static constexpr size_t n = std::decay_t<decltype(Impl::asMatrix(std::declval<K>()))>::rows;
162 void apply (DynamicVector<field_type>& v, DynamicVector<field_type>& d)
164 assert(v.size() > 0);
165 assert(v.size() == d.size());
166 assert(A.rows() <= v.size());
167 assert(A.cols() <= v.size());
168 size_t sz = A.rows();
172 v.resize(v.capacity());
173 d.resize(d.capacity());
186 size_t sz = rowset.size();
188 typedef typename S::const_iterator SIter;
190 for(SIter rowIdx = rowset.begin(), rowEnd=rowset.end();
191 rowIdx!= rowEnd; ++rowIdx, r++)
194 for(SIter colIdx = rowset.begin(), colEnd=rowset.end();
195 colIdx!= colEnd; ++colIdx, c++)
197 if (BCRS[*rowIdx].find(*colIdx) == BCRS[*rowIdx].end())
199 for (
size_t i=0; i<
n; i++)
201 for (
size_t j=0; j<
n; j++)
203 A[r*
n+i][c*
n+j] = Impl::asMatrix(BCRS[*rowIdx][*colIdx])[i][j];
213 template<
typename T,
bool tag>
221 template<
class K,
class Al,
class X,
class Y>
230 static constexpr size_t n = std::decay_t<decltype(Impl::asMatrix(std::declval<K>()))>::rows;
257 DynamicVector<field_type> &
lhs();
264 DynamicVector<field_type> &
rhs();
296 DynamicVector<field_type> * rhs_;
299 DynamicVector<field_type> * lhs_;
307 std::size_t maxlength_;
310#if HAVE_SUPERLU || HAVE_SUITESPARSE_UMFPACK
311 template<
template<
class>
class S,
typename T,
typename A>
321 static constexpr size_t n = std::decay_t<decltype(Impl::asMatrix(std::declval<T>()))>::rows;
322 static constexpr size_t m = std::decay_t<decltype(Impl::asMatrix(std::declval<T>()))>::cols;
393 std::size_t maxlength_;
398 template<
class M,
class X,
class Y>
481 template<
class M,
class X,
class Y>
500 template<
class M,
class X,
class Y>
518 template<
typename S,
typename T>
522 template<
typename S,
typename T,
typename A>
531 static constexpr size_t n = std::decay_t<decltype(Impl::asVector(std::declval<T>()))>::dimension;
540 template<
typename S,
typename T>
544 template<
typename S,
typename T,
typename A>
553 static constexpr size_t n = std::decay_t<decltype(Impl::asVector(std::declval<T>()))>::dimension;
570 template<
typename T,
class X,
class S>
574 template<
class X,
class S>
580 template<
class X,
class S>
586 template<
class X,
class S>
603 template<
typename T1,
typename T2,
bool forward>
631 template<
typename T1,
typename T2>
679 template<
class M,
class X,
class TD,
class TA>
692 template<
class T,
bool tag>
699 template<
class K,
class Al,
class X,
class Y>
703 static constexpr size_t n = std::decay_t<decltype(Impl::asMatrix(std::declval<K>()))>::rows;
704 template<
class RowToDomain,
class Solvers,
class SubDomains>
706 Solvers& solvers,
const SubDomains& domains,
710 template<
template<
class>
class S,
typename T,
typename A>
714 static constexpr size_t n = std::decay_t<decltype(Impl::asMatrix(std::declval<T>()))>::rows;
715 template<
class RowToDomain,
class Solvers,
class SubDomains>
717 Solvers& solvers,
const SubDomains& domains,
721 template<
class M,
class X,
class Y>
725 template<
class RowToDomain,
class Solvers,
class SubDomains>
727 Solvers& solvers,
const SubDomains& domains,
731 template<
class M,
class X,
class Y>
736 template<
class M,
class X,
class Y>
792 typedef std::set<size_type, std::less<size_type>,
793 typename std::allocator_traits<TA>::template rebind_alloc<size_type> >
797 typedef std::vector<subdomain_type, typename std::allocator_traits<TA>::template rebind_alloc<subdomain_type> >
subdomain_vector;
800 typedef SLList<size_type, typename std::allocator_traits<TA>::template rebind_alloc<size_type> >
subdomain_list;
803 typedef std::vector<subdomain_list, typename std::allocator_traits<TA>::template rebind_alloc<subdomain_list> >
rowtodomain_vector;
809 typedef std::vector<slu, typename std::allocator_traits<TA>::template rebind_alloc<slu> >
slu_vector;
825 field_type relaxationFactor=1,
bool onTheFly_=
true);
839 field_type relaxationFactor=1,
bool onTheFly_=
true);
846 virtual void pre ([[maybe_unused]] X& x, [[maybe_unused]] X& b)
854 virtual void apply (X& v,
const X& d);
861 virtual void post ([[maybe_unused]] X& x)
864 template<
bool forward>
879 typename M::size_type maxlength;
886 template<
class I,
class S,
class D>
890 : initializers(&il), indices(&idx), indexMaps(il.size()), domains(domains_)
894 template<
class I,
class S,
class D>
897 typedef typename IndexSet::value_type::const_iterator iterator;
898 for(iterator domain=(*indices)[row.index()].begin(); domain != (*indices)[row.index()].end(); ++domain) {
899 (*initializers)[*domain].addRowNnz(row, domains[*domain]);
900 indexMaps[*domain].insert(row.index());
904 template<
class I,
class S,
class D>
907 for(
auto&& i: *initializers)
908 i.allocateMatrixStorage();
909 for(
auto&& i: *initializers)
913 template<
class I,
class S,
class D>
916 typedef typename IndexSet::value_type::const_iterator iterator;
917 for(iterator domain=(*indices)[row.index()].begin(); domain != (*indices)[row.index()].end(); ++domain) {
918 typename std::map<size_type,size_type>::const_iterator v = indexMaps[*domain].find(
col.index());
919 if(v!= indexMaps[*domain].end()) {
920 (*initializers)[*domain].countEntries(indexMaps[*domain].find(
col.index())->second);
925 template<
class I,
class S,
class D>
928 for(
auto&& i : *initializers)
932 template<
class I,
class S,
class D>
935 typedef typename IndexSet::value_type::const_iterator iterator;
936 for(iterator domain=(*indices)[row.index()].begin(); domain!= (*indices)[row.index()].end(); ++domain) {
937 typename std::map<size_type,size_type>::const_iterator v = indexMaps[*domain].find(
col.index());
938 if(v!= indexMaps[*domain].end()) {
939 assert(indexMaps[*domain].end()!=indexMaps[*domain].find(row.index()));
940 (*initializers)[*domain].copyValue(
col, indexMaps[*domain].find(row.index())->second,
946 template<
class I,
class S,
class D>
949 std::vector<IndexMap>().swap(indexMaps);
950 for(
auto&& i: *initializers)
954 template<
class I,
class S,
class D>
959 template<
class I,
class S,
class D>
962 assert(map_.find(grow)==map_.end());
963 map_.insert(std::make_pair(grow, row++));
966 template<
class I,
class S,
class D>
967 typename OverlappingSchwarzInitializer<I,S,D>::IndexMap::const_iterator
970 return map_.find(grow);
973 template<
class I,
class S,
class D>
974 typename OverlappingSchwarzInitializer<I,S,D>::IndexMap::iterator
977 return map_.find(grow);
980 template<
class I,
class S,
class D>
981 typename OverlappingSchwarzInitializer<I,S,D>::IndexMap::const_iterator
987 template<
class I,
class S,
class D>
988 typename OverlappingSchwarzInitializer<I,S,D>::IndexMap::iterator
994 template<
class I,
class S,
class D>
995 typename OverlappingSchwarzInitializer<I,S,D>::IndexMap::const_iterator
1001 template<
class I,
class S,
class D>
1002 typename OverlappingSchwarzInitializer<I,S,D>::IndexMap::iterator
1005 return map_.begin();
1008 template<
class M,
class X,
class TM,
class TD,
class TA>
1011 : mat(mat_), relax(relaxationFactor), onTheFly(fly)
1013 typedef typename rowtodomain_vector::const_iterator RowDomainIterator;
1014 typedef typename subdomain_list::const_iterator DomainIterator;
1015#ifdef DUNE_ISTL_WITH_CHECKING
1016 assert(rowToDomain.size()==mat.N());
1017 assert(rowToDomain.size()==mat.M());
1019 for(RowDomainIterator iter=rowToDomain.begin(); iter != rowToDomain.end(); ++iter)
1020 assert(iter->size()>0);
1025 for(RowDomainIterator iter=rowToDomain.begin(); iter != rowToDomain.end(); ++iter)
1026 for(DomainIterator d=iter->begin(); d != iter->end(); ++d)
1027 domains=std::max(domains, *d);
1030 solvers.resize(domains);
1031 subDomains.resize(domains);
1035 for(RowDomainIterator iter=rowToDomain.begin(); iter != rowToDomain.end(); ++iter, ++row)
1036 for(DomainIterator d=iter->begin(); d != iter->end(); ++d)
1037 subDomains[*d].insert(row);
1039#ifdef DUNE_ISTL_WITH_CHECKING
1041 typedef typename subdomain_vector::const_iterator iterator;
1042 for(iterator iter=subDomains.begin(); iter != subDomains.end(); ++iter) {
1043 typedef typename subdomain_type::const_iterator entry_iterator;
1044 Dune::dvverb<<
"domain "<<i++<<
":";
1045 for(entry_iterator entry = iter->begin(); entry != iter->end(); ++entry) {
1046 Dune::dvverb<<
" "<<*entry;
1048 Dune::dvverb<<std::endl;
1055 template<
class M,
class X,
class TM,
class TD,
class TA>
1060 : mat(mat_), solvers(sd.size()), subDomains(sd), relax(relaxationFactor),
1063 typedef typename subdomain_vector::const_iterator DomainIterator;
1065#ifdef DUNE_ISTL_WITH_CHECKING
1068 for(DomainIterator d=sd.begin(); d != sd.end(); ++d,++i) {
1070 assert(d->size()>0);
1071 typedef typename DomainIterator::value_type::const_iterator entry_iterator;
1072 Dune::dvverb<<
"domain "<<i<<
":";
1073 for(entry_iterator entry = d->begin(); entry != d->end(); ++entry) {
1074 Dune::dvverb<<
" "<<*entry;
1076 Dune::dvverb<<std::endl;
1082 rowtodomain_vector rowToDomain(mat.N());
1086 for(DomainIterator domain=sd.begin(); domain != sd.end(); ++domain, ++domainId) {
1087 typedef typename subdomain_type::const_iterator iterator;
1088 for(iterator row=domain->begin(); row != domain->end(); ++row)
1089 rowToDomain[*row].push_back(domainId);
1092 maxlength = SeqOverlappingSchwarzAssembler<slu>
1093 ::assembleLocalProblems(rowToDomain, mat, solvers, subDomains, onTheFly);
1105 template<
typename T,
typename A>
1108 static constexpr size_t n = std::decay_t<decltype(Impl::asMatrix(std::declval<T>()))>::rows;
1109 static constexpr size_t m = std::decay_t<decltype(Impl::asMatrix(std::declval<T>()))>::cols;
1110 template<
class Domain>
1118 template<
class K,
class Al,
class X,
class Y>
1119 template<
class RowToDomain,
class Solvers,
class SubDomains>
1122 assembleLocalProblems([[maybe_unused]]
const RowToDomain& rowToDomain,
1124 [[maybe_unused]] Solvers& solvers,
1125 const SubDomains& subDomains,
1126 [[maybe_unused]]
bool onTheFly)
1128 typedef typename SubDomains::const_iterator DomainIterator;
1129 std::size_t maxlength = 0;
1133 for(DomainIterator domain=subDomains.begin(); domain!=subDomains.end(); ++domain)
1134 maxlength=std::max(maxlength, domain->size());
1140#if HAVE_SUPERLU || HAVE_SUITESPARSE_UMFPACK
1141 template<
template<
class>
class S,
typename T,
typename A>
1142 template<
class RowToDomain,
class Solvers,
class SubDomains>
1146 const SubDomains& subDomains,
1149 typedef typename S<BCRSMatrix<T,A>>::MatrixInitializer MatrixInitializer;
1150 typedef typename std::vector<MatrixInitializer>::iterator InitializerIterator;
1151 typedef typename SubDomains::const_iterator DomainIterator;
1152 typedef typename Solvers::iterator SolverIterator;
1153 std::size_t maxlength = 0;
1156 for(DomainIterator domain=subDomains.begin(); domain!=subDomains.end(); ++domain)
1157 maxlength=std::max(maxlength, domain->size());
1158 maxlength*=Impl::asMatrix(*mat[0].begin()).N();
1161 DomainIterator domain=subDomains.begin();
1164 std::vector<MatrixInitializer> initializers(subDomains.size());
1166 SolverIterator solver=solvers.begin();
1167 for(InitializerIterator initializer=initializers.begin(); initializer!=initializers.end();
1168 ++initializer, ++solver, ++domain) {
1172 *initializer=MatrixInitializer(solver->getInternalMatrix());
1177 RowToDomain, SubDomains> Initializer;
1179 Initializer initializer(initializers, rowToDomain, subDomains);
1180 copyToBCCSMatrix(initializer, mat);
1183 for(
auto&& s: solvers)
1185 for (SolverIterator solverIt = solvers.begin(); solverIt != solvers.end(); ++solverIt)
1187 assert(solverIt->getInternalMatrix().N() == solverIt->getInternalMatrix().M());
1188 maxlength = std::max(maxlength, solverIt->getInternalMatrix().N());
1196 template<
class M,
class X,
class Y>
1197 template<
class RowToDomain,
class Solvers,
class SubDomains>
1201 const SubDomains& subDomains,
1204 typedef typename SubDomains::const_iterator DomainIterator;
1205 typedef typename Solvers::iterator SolverIterator;
1206 std::size_t maxlength = 0;
1209 for(DomainIterator domain=subDomains.begin(); domain!=subDomains.end(); ++domain)
1210 maxlength=std::max(maxlength, domain->size());
1213 SolverIterator solver=solvers.begin();
1214 for(DomainIterator domain=subDomains.begin(); domain!=subDomains.end();
1215 ++domain, ++solver) {
1216 solver->setSubMatrix(mat, *domain);
1217 maxlength=std::max(maxlength, domain->size());
1226 template<
class M,
class X,
class TM,
class TD,
class TA>
1232 template<
class M,
class X,
class TM,
class TD,
class TA>
1233 template<
bool forward>
1249 Adder adder(v, x, assigner, relax);
1253 std::for_each(domain->begin(), domain->end(), assigner);
1254 assigner.resetIndexForNextDomain();
1258 sdsolver.setSubMatrix(mat, *domain);
1260 sdsolver.apply(assigner.lhs(), assigner.rhs());
1262 solver->apply(assigner.lhs(), assigner.rhs());
1267 std::for_each(domain->begin(), domain->end(), adder);
1268 assigner.resetIndexForNextDomain();
1273 assigner.deallocate();
1276 template<
class K,
class Al,
class X,
class Y>
1279 const X& b_, Y& x_) :
1281 rhs_( new DynamicVector<
field_type>(maxlength, 42) ),
1282 lhs_( new DynamicVector<
field_type>(maxlength, -42) ),
1286 maxlength_(maxlength)
1289 template<
class K,
class Al,
class X,
class Y>
1298 template<
class K,
class Al,
class X,
class Y>
1306 template<
class K,
class Al,
class X,
class Y>
1307 DynamicVector<typename X::field_type> &
1314 template<
class K,
class Al,
class X,
class Y>
1315 DynamicVector<typename X::field_type> &
1322 template<
class K,
class Al,
class X,
class Y>
1330 template<
class K,
class Al,
class X,
class Y>
1339 assert(i<maxlength_);
1340 rhs()[i]=(*b)[domainIndex][j];
1344 typedef typename matrix_type::ConstColIterator col_iterator;
1347 for(col_iterator
col=(*mat)[domainIndex].begin();
col!=(*mat)[domainIndex].end(); ++
col) {
1349 (*col).mv((*x)[
col.index()], tmp);
1352 assert(i<maxlength_);
1359 assert(i<maxlength_);
1360 rhs()[i]=Impl::asVector((*b)[domainIndex])[j];
1363 typedef typename matrix_type::ConstColIterator col_iterator;
1366 for(col_iterator
col=(*mat)[domainIndex].begin();
col!=(*mat)[domainIndex].end(); ++
col) {
1368 rhs()[i]-=Impl::asMatrix(*
col)[j][k] * Impl::asVector((*x)[
col.index()])[k];
1375 template<
class K,
class Al,
class X,
class Y>
1382 assert(i<maxlength_);
1383 Impl::asVector(res)[j]+=
lhs()[i];
1387#if HAVE_SUPERLU || HAVE_SUITESPARSE_UMFPACK
1389 template<
template<
class>
class S,
typename T,
typename A>
1397 x(&x_), i(0), maxlength_(maxlength)
1404 template<
template<
class>
class S,
typename T,
typename A>
1411 template<
template<
class>
class S,
typename T,
typename A>
1418 assert(i<maxlength_);
1419 rhs_[i]=Impl::asVector((*b)[domainIndex])[j];
1424 typedef typename matrix_type::ConstColIterator col_iterator;
1427 for(col_iterator
col=(*mat)[domainIndex].begin();
col!=(*mat)[domainIndex].end(); ++
col) {
1429 Impl::asMatrix(*col).mv((*x)[
col.index()], tmp);
1432 assert(i<maxlength_);
1433 rhs_[i]-=Impl::asVector(tmp)[j];
1440 template<
template<
class>
class S,
typename T,
typename A>
1444 assert(i<maxlength_);
1450 template<
template<
class>
class S,
typename T,
typename A>
1455 assert(i<maxlength_);
1456 Impl::asVector(res)[j]+=lhs_[i];
1460 template<
template<
class>
class S,
typename T,
typename A>
1466 template<
template<
class>
class S,
typename T,
typename A>
1473 template<
template<
class>
class S,
typename T,
typename A>
1482 template<
class M,
class X,
class Y>
1491 rhs_=
new Y(maxlength);
1492 lhs_ =
new X(maxlength);
1495 template<
class M,
class X,
class Y>
1502 template<
class M,
class X,
class Y>
1505 (*rhs_)[i]=(*b)[domainIndex];
1508 typedef typename matrix_type::ConstColIterator col_iterator;
1511 for(col_iterator
col=(*mat)[domainIndex].begin();
col!=(*mat)[domainIndex].end(); ++
col) {
1512 Impl::asMatrix(*col).mmv((*x)[
col.index()], (*rhs_)[i]);
1518 template<
class M,
class X,
class Y>
1524 template<
class M,
class X,
class Y>
1530 template<
class M,
class X,
class Y>
1536 template<
class M,
class X,
class Y>
1542 template<
class M,
class X,
class Y>
1548 template<
typename S,
typename T,
typename A>
1553 : v(&v_), x(&x_), assigner(&assigner_), relax(relax_)
1556 template<
typename S,
typename T,
typename A>
1560 assigner->assignResult((*v)[domainIndex]);
1564 template<
typename S,
typename T,
typename A>
1572 template<
typename S,
typename T,
typename A>
1577 : x(&x_), assigner(&assigner_), relax(relax_)
1581 template<
typename S,
typename T,
typename A>
1585 assigner->relaxResult(relax);
1586 assigner->assignResult((*x)[domainIndex]);
1590 template<
typename S,
typename T,
typename A>
Classes for using SuperLU with ISTL matrices.
Templates characterizing the type of a solver.
This file implements a vector space as a tensor product of a given vector space. The number of compon...
Implementation of the BCRSMatrix class.
Classes for using UMFPack with ISTL matrices.
Various local subdomain solvers based on ILU for SeqOverlappingSchwarz.
Define general preconditioner interface.
*The type for the index access and the size typedef A::size_type size_type
Definition bcrsmatrix.hh:497
typename Imp::BlockTraits< B >::field_type field_type
Definition bcrsmatrix.hh:488
void addRowNnz(const Iter &row)
Definition overlappingschwarz.hh:895
void apply(X &v, const X &d)
Apply one step of the preconditioner to the system A(v)=d.
Definition overlappingschwarz.hh:1234
X & lhs()
Get the local left hand side.
Definition overlappingschwarz.hh:1531
void calcColstart() const
Definition overlappingschwarz.hh:926
MultiplicativeAdder(BlockVector< T, A > &v, BlockVector< T, A > &x, OverlappingAssigner< S > &assigner_, const field_type &relax_)
Definition overlappingschwarz.hh:1574
Y & rhs()
Get the local right hand side.
Definition overlappingschwarz.hh:1537
iterator end()
Definition overlappingschwarz.hh:989
void operator()(const size_type &domain)
Definition overlappingschwarz.hh:1582
OverlappingAssignerHelper(std::size_t maxlength, const BCRSMatrix< K, Al > &mat_, const X &b_, Y &x_)
Constructor.
Definition overlappingschwarz.hh:1278
void resetIndexForNextDomain()
Resets the local index to zero.
Definition overlappingschwarz.hh:1543
AdditiveAdder(BlockVector< T, A > &v, BlockVector< T, A > &x, OverlappingAssigner< S > &assigner, const field_type &relax_)
Definition overlappingschwarz.hh:1549
void copyValue(const Iter &row, const CIter &col) const
Definition overlappingschwarz.hh:933
void operator()(const size_type &domainIndex)
calculate one entry of the local defect.
Definition overlappingschwarz.hh:1333
void createMatrix() const
Definition overlappingschwarz.hh:947
field_type * lhs()
Get the local left hand side.
Definition overlappingschwarz.hh:1468
iterator begin()
Definition overlappingschwarz.hh:1003
OverlappingSchwarzInitializer(InitializerList &il, const IndexSet &indices, const subdomain_vector &domains)
Definition overlappingschwarz.hh:887
IndexMap()
Definition overlappingschwarz.hh:955
virtual void apply(X &v, const X &d)
Apply the preconditioner.
Definition overlappingschwarz.hh:1227
field_type * rhs()
Get the local right hand side.
Definition overlappingschwarz.hh:1475
void operator()(const size_type &domain)
Definition overlappingschwarz.hh:1557
OverlappingAssignerILUBase(std::size_t maxlength, const M &mat, const Y &b, X &x)
Constructor.
Definition overlappingschwarz.hh:1483
static std::size_t assembleLocalProblems(const RowToDomain &rowToDomain, const matrix_type &mat, Solvers &solvers, const SubDomains &domains, bool onTheFly)
Definition overlappingschwarz.hh:1122
const_iterator find(size_type grow) const
Definition overlappingschwarz.hh:968
void relaxResult(field_type relax)
relax the result.
Definition overlappingschwarz.hh:1441
void operator()(const size_type &domain)
calculate one entry of the local defect.
Definition overlappingschwarz.hh:1503
SeqOverlappingSchwarz(const matrix_type &mat, const subdomain_vector &subDomains, field_type relaxationFactor=1, bool onTheFly_=true)
Construct the overlapping Schwarz method.
Definition overlappingschwarz.hh:1056
void allocate()
Definition overlappingschwarz.hh:905
void deallocate()
Deallocates memory of the local vector.
Definition overlappingschwarz.hh:1496
void axpy()
Definition overlappingschwarz.hh:1565
DynamicVector< field_type > & lhs()
Get the local left hand side.
Definition overlappingschwarz.hh:1309
void axpy()
Definition overlappingschwarz.hh:1591
void deallocate()
Deallocates memory of the local vector.
Definition overlappingschwarz.hh:1405
SeqOverlappingSchwarzAssemblerHelper< T, Dune::StoresColumnCompressed< T >::value > SeqOverlappingSchwarzAssembler
Definition overlappingschwarz.hh:697
void resetIndexForNextDomain()
Resets the local index to zero.
Definition overlappingschwarz.hh:1301
void assignResult(block_type &res)
Assigns the block to the current local index. At the same time the local defect is calculated for the...
Definition overlappingschwarz.hh:1378
void insert(size_type grow)
Definition overlappingschwarz.hh:960
void countEntries(const Iter &row, const CIter &col) const
Definition overlappingschwarz.hh:914
void operator()(const size_type &domain)
calculate one entry of the local defect.
Definition overlappingschwarz.hh:1412
void assignResult(block_type &res)
Assigns the block to the current local index. At the same time the local defect is calculated for the...
Definition overlappingschwarz.hh:1451
static std::size_t assembleLocalProblems(const RowToDomain &rowToDomain, const matrix_type &mat, Solvers &solvers, const SubDomains &domains, bool onTheFly)
Definition overlappingschwarz.hh:1143
static std::size_t assembleLocalProblems(const RowToDomain &rowToDomain, const matrix_type &mat, Solvers &solvers, const SubDomains &domains, bool onTheFly)
Definition overlappingschwarz.hh:1198
void relaxResult(field_type relax)
relax the result.
Definition overlappingschwarz.hh:1519
OverlappingAssignerHelper< T, Dune::StoresColumnCompressed< T >::value > OverlappingAssigner
Definition overlappingschwarz.hh:218
void resetIndexForNextDomain()
Definition overlappingschwarz.hh:1461
OverlappingAssignerHelper(std::size_t maxlength, const matrix_type &mat, const range_type &b, range_type &x)
Constructor.
Definition overlappingschwarz.hh:1391
void relaxResult(field_type relax)
relax the result.
Definition overlappingschwarz.hh:1325
void deallocate()
Deallocates memory of the local vector.
Definition overlappingschwarz.hh:1292
void assignResult(block_type &res)
Assigns the block to the current local index. At the same time the local defect is calculated for the...
Definition overlappingschwarz.hh:1525
SeqOverlappingSchwarz(const matrix_type &mat, const rowtodomain_vector &rowToDomain, field_type relaxationFactor=1, bool onTheFly_=true)
Definition overlappingschwarz.hh:1009
DynamicVector< field_type > & rhs()
Get the local right hand side.
Definition overlappingschwarz.hh:1317
Definition allocator.hh:11
MatrixSparsityPatternGatherScatter< M, I >::ColIter MatrixSparsityPatternGatherScatter< M, I >::col
Definition matrixredistribute.hh:591
Initializer for SuperLU Matrices representing the subdomains.
Definition overlappingschwarz.hh:47
Matrix::row_type::const_iterator CIter
Definition overlappingschwarz.hh:56
S IndexSet
Definition overlappingschwarz.hh:58
Matrix::const_iterator Iter
Definition overlappingschwarz.hh:55
IndexSet::size_type size_type
Definition overlappingschwarz.hh:59
I InitializerList
Definition overlappingschwarz.hh:52
AtomInitializer::Matrix Matrix
Definition overlappingschwarz.hh:54
InitializerList::value_type AtomInitializer
Definition overlappingschwarz.hh:53
D subdomain_vector
The vector type containing the subdomain to row index mapping.
Definition overlappingschwarz.hh:50
A vector of blocks with memory management.
Definition bvector.hh:392
Exact subdomain solver using ILU(p) with appropriate p.
Definition ilusubdomainsolver.hh:78
Definition ilusubdomainsolver.hh:111
Sequential overlapping Schwarz preconditioner.
Definition overlappingschwarz.hh:755
X::field_type field_type
Definition overlappingschwarz.hh:783
SLList< size_type, typename std::allocator_traits< TA >::template rebind_alloc< size_type > > subdomain_list
Definition overlappingschwarz.hh:800
std::vector< slu, typename std::allocator_traits< TA >::template rebind_alloc< slu > > slu_vector
Definition overlappingschwarz.hh:809
M matrix_type
Definition overlappingschwarz.hh:760
SymmetricMultiplicativeSchwarzMode Mode
Definition overlappingschwarz.hh:778
X range_type
Definition overlappingschwarz.hh:770
std::set< size_type, std::less< size_type >, typename std::allocator_traits< TA >::template rebind_alloc< size_type > > subdomain_type
Definition overlappingschwarz.hh:794
X domain_type
Definition overlappingschwarz.hh:765
TD slu
Definition overlappingschwarz.hh:806
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category).
Definition overlappingschwarz.hh:868
virtual void pre(X &x, X &b)
Prepare the preconditioner.
Definition overlappingschwarz.hh:846
virtual void post(X &x)
Postprocess the preconditioner.
Definition overlappingschwarz.hh:861
TA allocator
Definition overlappingschwarz.hh:789
std::vector< subdomain_type, typename std::allocator_traits< TA >::template rebind_alloc< subdomain_type > > subdomain_vector
Definition overlappingschwarz.hh:797
std::vector< subdomain_list, typename std::allocator_traits< TA >::template rebind_alloc< subdomain_list > > rowtodomain_vector
Definition overlappingschwarz.hh:803
matrix_type::size_type size_type
Definition overlappingschwarz.hh:786
Definition overlappingschwarz.hh:694
Definition matrixutils.hh:24
Tag that the tells the Schwarz method to be additive.
Definition overlappingschwarz.hh:120
Tag that tells the Schwarz method to be multiplicative.
Definition overlappingschwarz.hh:126
Tag that tells the Schwarz method to be multiplicative and symmetric.
Definition overlappingschwarz.hh:133
Exact subdomain solver using Dune::DynamicMatrix<T>::solve.
Definition overlappingschwarz.hh:140
std::remove_const< M >::type matrix_type
The matrix type the preconditioner is for.
Definition overlappingschwarz.hh:149
static constexpr size_t n
Definition overlappingschwarz.hh:156
X::field_type field_type
Definition overlappingschwarz.hh:150
X domain_type
The domain type of the preconditioner.
Definition overlappingschwarz.hh:153
Y range_type
The range type of the preconditioner.
Definition overlappingschwarz.hh:155
void setSubMatrix(const M &BCRS, S &rowset)
Set the data of the local problem.
Definition overlappingschwarz.hh:184
void apply(DynamicVector< field_type > &v, DynamicVector< field_type > &d)
Apply the subdomain solver.
Definition overlappingschwarz.hh:162
std::remove_const< M >::type rilu_type
Definition overlappingschwarz.hh:151
Definition overlappingschwarz.hh:215
matrix_type::size_type size_type
Definition overlappingschwarz.hh:229
X::field_type field_type
Definition overlappingschwarz.hh:226
BCRSMatrix< K, Al > matrix_type
Definition overlappingschwarz.hh:225
Y range_type
Definition overlappingschwarz.hh:227
static constexpr size_t n
Definition overlappingschwarz.hh:230
range_type::block_type block_type
Definition overlappingschwarz.hh:228
S< BCRSMatrix< T, A > >::range_type range_type
Definition overlappingschwarz.hh:315
range_type::block_type block_type
Definition overlappingschwarz.hh:317
static constexpr size_t m
Definition overlappingschwarz.hh:322
range_type::field_type field_type
Definition overlappingschwarz.hh:316
static constexpr size_t n
Definition overlappingschwarz.hh:321
matrix_type::size_type size_type
Definition overlappingschwarz.hh:319
BCRSMatrix< T, A > matrix_type
Definition overlappingschwarz.hh:314
matrix_type::size_type size_type
Definition overlappingschwarz.hh:408
Y::field_type field_type
Definition overlappingschwarz.hh:404
Y::block_type block_type
Definition overlappingschwarz.hh:406
M matrix_type
Definition overlappingschwarz.hh:402
OverlappingAssignerHelper(std::size_t maxlength, const M &mat, const Y &b, X &x)
Constructor.
Definition overlappingschwarz.hh:493
OverlappingAssignerHelper(std::size_t maxlength, const M &mat, const Y &b, X &x)
Constructor.
Definition overlappingschwarz.hh:512
Definition overlappingschwarz.hh:520
std::decay_t< decltype(Impl::asVector(std::declval< T >()))>::field_type field_type
Definition overlappingschwarz.hh:526
static constexpr size_t n
Definition overlappingschwarz.hh:531
A::size_type size_type
Definition overlappingschwarz.hh:525
Definition overlappingschwarz.hh:542
static constexpr size_t n
Definition overlappingschwarz.hh:553
A::size_type size_type
Definition overlappingschwarz.hh:547
std::decay_t< decltype(Impl::asVector(std::declval< T >()))>::field_type field_type
Definition overlappingschwarz.hh:548
template meta program for choosing how to add the correction.
Definition overlappingschwarz.hh:572
AdditiveAdder< S, X > Adder
Definition overlappingschwarz.hh:577
MultiplicativeAdder< S, X > Adder
Definition overlappingschwarz.hh:583
MultiplicativeAdder< S, X > Adder
Definition overlappingschwarz.hh:589
Helper template meta program for application of overlapping Schwarz.
Definition overlappingschwarz.hh:605
static solver_iterator begin(solver_vector &sv)
Definition overlappingschwarz.hh:611
solver_vector::iterator solver_iterator
Definition overlappingschwarz.hh:607
static domain_iterator end(const subdomain_vector &sv)
Definition overlappingschwarz.hh:625
subdomain_vector::const_iterator domain_iterator
Definition overlappingschwarz.hh:609
T1 solver_vector
Definition overlappingschwarz.hh:606
static domain_iterator begin(const subdomain_vector &sv)
Definition overlappingschwarz.hh:620
T2 subdomain_vector
Definition overlappingschwarz.hh:608
static solver_iterator end(solver_vector &sv)
Definition overlappingschwarz.hh:616
T2 subdomain_vector
Definition overlappingschwarz.hh:636
static solver_iterator end(solver_vector &sv)
Definition overlappingschwarz.hh:644
solver_vector::reverse_iterator solver_iterator
Definition overlappingschwarz.hh:635
subdomain_vector::const_reverse_iterator domain_iterator
Definition overlappingschwarz.hh:637
static solver_iterator begin(solver_vector &sv)
Definition overlappingschwarz.hh:639
T1 solver_vector
Definition overlappingschwarz.hh:634
static domain_iterator begin(const subdomain_vector &sv)
Definition overlappingschwarz.hh:648
static domain_iterator end(const subdomain_vector &sv)
Definition overlappingschwarz.hh:653
Helper template meta program for application of overlapping Schwarz.
Definition overlappingschwarz.hh:669
smoother::range_type range_type
Definition overlappingschwarz.hh:671
T smoother
Definition overlappingschwarz.hh:670
static void apply(smoother &sm, range_type &v, const range_type &b)
Definition overlappingschwarz.hh:673
static void apply(smoother &sm, range_type &v, const range_type &b)
Definition overlappingschwarz.hh:685
SeqOverlappingSchwarz< M, X, SymmetricMultiplicativeSchwarzMode, TD, TA > smoother
Definition overlappingschwarz.hh:682
smoother::range_type range_type
Definition overlappingschwarz.hh:683
BCRSMatrix< K, Al > matrix_type
Definition overlappingschwarz.hh:702
static constexpr size_t n
Definition overlappingschwarz.hh:703
BCRSMatrix< T, A > matrix_type
Definition overlappingschwarz.hh:713
static constexpr size_t n
Definition overlappingschwarz.hh:714
Definition overlappingschwarz.hh:723
M matrix_type
Definition overlappingschwarz.hh:724
Definition overlappingschwarz.hh:1103
static int size(const Domain &d)
Definition overlappingschwarz.hh:1111
static constexpr size_t m
Definition overlappingschwarz.hh:1109
static constexpr size_t n
Definition overlappingschwarz.hh:1108
Base class for matrix free definition of preconditioners.
Definition preconditioner.hh:33
Category
Definition solvercategory.hh:23
@ sequential
Category for sequential solvers.
Definition solvercategory.hh:25