5#ifndef DUNE_ISTL_ILU_HH
6#define DUNE_ISTL_ILU_HH
15#include <dune/common/fmatrix.hh>
16#include <dune/common/scalarvectorview.hh>
17#include <dune/common/scalarmatrixview.hh>
38 typedef typename M::RowIterator rowiterator;
39 typedef typename M::ColIterator coliterator;
40 typedef typename M::block_type block;
43 rowiterator endi=A.end();
44 for (rowiterator i=A.begin(); i!=endi; ++i)
47 coliterator endij=(*i).end();
51 for (ij=(*i).begin(); ij.index()<i.index(); ++ij)
54 coliterator jj = A[ij.index()].find(ij.index());
57 Impl::asMatrix(*ij).rightmultiply(Impl::asMatrix(*jj));
60 coliterator endjk=A[ij.index()].end();
61 coliterator jk=jj; ++jk;
62 coliterator ik=ij; ++ik;
63 while (ik!=endij && jk!=endjk)
64 if (ik.index()==jk.index())
67 Impl::asMatrix(B).leftmultiply(Impl::asMatrix(*ij));
73 if (ik.index()<jk.index())
81 if (ij.index()!=i.index())
82 DUNE_THROW(
ISTLError,
"diagonal entry missing");
84 Impl::asMatrix(*ij).invert();
86 catch (Dune::FMatrixError & e) {
87 std::ostringstream sstream;
90 <<
"ILU failed to invert matrix block A["
91 << i.index() <<
"][" << ij.index() <<
"]" << e.what();
93 ex.message(sstream.str());
102 template<
class M,
class X,
class Y>
106 typedef typename M::ConstRowIterator rowiterator;
107 typedef typename M::ConstColIterator coliterator;
108 typedef typename Y::block_type dblock;
109 typedef typename X::block_type vblock;
112 rowiterator endi=A.end();
113 for (rowiterator i=A.begin(); i!=endi; ++i)
121 dblock rhsValue(d[i.index()]);
122 auto&& rhs = Impl::asVector(rhsValue);
123 for (coliterator j=(*i).begin(); j.index()<i.index(); ++j)
124 Impl::asMatrix(*j).mmv(Impl::asVector(v[j.index()]),rhs);
125 Impl::asVector(v[i.index()]) = rhs;
129 rowiterator rendi=A.beforeBegin();
130 for (rowiterator i=A.beforeEnd(); i!=rendi; --i)
138 vblock rhsValue(v[i.index()]);
139 auto&& rhs = Impl::asVector(rhsValue);
141 for (j=(*i).beforeEnd(); j.index()>i.index(); --j)
142 Impl::asMatrix(*j).mmv(Impl::asVector(v[j.index()]),rhs);
143 auto&& vi = Impl::asVector(v[i.index()]);
144 Impl::asMatrix(*j).mv(rhs,vi);
151 [[maybe_unused]]
typename std::enable_if_t<!Dune::IsNumber<M>::value>* sfinae =
nullptr)
158 [[maybe_unused]]
typename std::enable_if_t<Dune::IsNumber<K>::value>* sfinae =
nullptr)
163 template<
class K,
int n,
int m>
179 typedef typename M::ColIterator coliterator;
180 typedef typename M::ConstRowIterator crowiterator;
181 typedef typename M::ConstColIterator ccoliterator;
182 typedef typename M::CreateIterator createiterator;
183 typedef typename M::field_type K;
184 typedef std::map<size_t, int> map;
185 typedef typename map::iterator mapiterator;
188 crowiterator endi=A.end();
189 createiterator ci=
ILU.createbegin();
190 for (crowiterator i=A.begin(); i!=endi; ++i)
195 for (ccoliterator j=(*i).begin(); j!=(*i).end(); ++j)
196 rowpattern[j.index()] = 0;
199 for (mapiterator ik=rowpattern.begin(); (*ik).first<i.index(); ++ik)
203 coliterator endk =
ILU[(*ik).first].end();
204 coliterator kj =
ILU[(*ik).first].find((*ik).first);
205 for (++kj; kj!=endk; ++kj)
213 mapiterator ij = rowpattern.find(kj.index());
214 if (ij==rowpattern.end())
216 rowpattern[kj.index()] = generation+1;
224 for (mapiterator ik=rowpattern.begin(); ik!=rowpattern.end(); ++ik)
225 ci.insert((*ik).first);
229 coliterator endILUij =
ILU[i.index()].end();;
230 for (coliterator ILUij=
ILU[i.index()].begin(); ILUij!=endILUij; ++ILUij)
231 Simd::lane(0,
firstMatrixElement(*ILUij)) = (Simd::Scalar<K>) rowpattern[ILUij.index()];
235 for (crowiterator i=A.begin(); i!=endi; ++i)
238 coliterator endILUij =
ILU[i.index()].end();;
239 for (ILUij=
ILU[i.index()].begin(); ILUij!=endILUij; ++ILUij)
241 ccoliterator Aij = (*i).begin();
242 ccoliterator endAij = (*i).end();
243 ILUij =
ILU[i.index()].begin();
244 while (Aij!=endAij && ILUij!=endILUij)
246 if (Aij.index()==ILUij.index())
253 if (Aij.index()<ILUij.index())
266 template <
class B,
class Alloc = std::allocator<B>>
294 if(
values_.capacity() < needed )
298 cols_.reserve( estimate );
305 cols_.push_back( index );
315 template<
class M,
class CRS,
class InvVector>
318 typedef typename M :: size_type
size_type;
325 const size_t memEstimate = (A.nonzeroes() - A.N())/2;
327 assert( A.nonzeroes() != 0 );
331 const auto endi = A.end();
334 lower.
rows_[ 0 ] = colcount;
335 for (
auto i=A.begin(); i!=endi; ++i, ++row)
340 for (
auto j=(*i).begin(); j.index() < iIndex; ++j )
345 lower.
rows_[ iIndex+1 ] = colcount;
351 upper.
rows_[ 0 ] = colcount ;
355 for (
auto i=A.beforeEnd(); i!=rendi; --i, ++ row )
357 const auto endij=(*i).beforeBegin();
362 for (
auto j=(*i).beforeEnd(); j != endij; --j )
365 if( j.index() == iIndex )
370 else if ( j.index() >= i.index() )
376 upper.
rows_[ row+1 ] = colcount;
381 template<
class CRS,
class InvVector,
class X,
class Y>
384 const InvVector& inv,
388 typedef typename Y :: block_type dblock;
389 typedef typename X :: block_type vblock;
390 typedef typename X :: size_type
size_type ;
394 if( iEnd != upper.
rows() )
396 DUNE_THROW(
ISTLError,
"ILU::blockILUBacksolve: lower and upper rows must be the same");
402 dblock rhsValue( d[ i ] );
403 auto&& rhs = Impl::asVector(rhsValue);
408 Impl::asMatrix(lower.
values_[
col ]).mmv( Impl::asVector(v[ lower.
cols_[
col ] ] ), rhs );
410 Impl::asVector(v[ i ]) = rhs;
416 auto&& vBlock = Impl::asVector(v[ lastRow - i ]);
417 vblock rhsValue ( v[ lastRow - i ] );
418 auto&& rhs = Impl::asVector(rhsValue);
423 Impl::asMatrix(upper.
values_[
col ]).mmv( Impl::asVector(v[ upper.
cols_[
col ] ]), rhs );
426 Impl::asMatrix(inv[ i ]).mv(rhs, vBlock);
Definition allocator.hh:11
MatrixSparsityPatternGatherScatter< M, I >::ColIter MatrixSparsityPatternGatherScatter< M, I >::col
Definition matrixredistribute.hh:591
void convertToCRS(const M &A, CRS &lower, CRS &upper, InvVector &inv)
convert ILU decomposition into CRS format for lower and upper triangular and inverse.
Definition ilu.hh:316
void blockILUBacksolve(const M &A, X &v, const Y &d)
LU backsolve with stored inverse.
Definition ilu.hh:103
M::field_type & firstMatrixElement(M &A, typename std::enable_if_t<!Dune::IsNumber< M >::value > *sfinae=nullptr)
Definition ilu.hh:150
void blockILU0Decomposition(M &A)
compute ILU decomposition of A. A is overwritten by its decomposition
Definition ilu.hh:35
void blockILUDecomposition(const M &A, int n, M &ILU)
Definition ilu.hh:176
Iterator beforeBegin()
Definition bcrsmatrix.hh:692
a simple compressed row storage matrix class
Definition ilu.hh:268
std::vector< size_type > cols_
Definition ilu.hh:310
size_type nonZeros() const
Definition ilu.hh:276
void resize(const size_type nRows)
Definition ilu.hh:282
size_type rows() const
Definition ilu.hh:274
CRS()
Definition ilu.hh:272
void reserveAdditional(const size_type nonZeros)
Definition ilu.hh:291
B block_type
Definition ilu.hh:269
std::vector< block_type, typename M::allocator_type > values_
Definition ilu.hh:309
size_type nRows_
Definition ilu.hh:311
size_t size_type
Definition ilu.hh:270
std::vector< size_type > rows_
Definition ilu.hh:308
void push_back(const block_type &value, const size_type index)
Definition ilu.hh:302
derive error class from the base class in common
Definition istlexception.hh:19
Error when performing an operation on a matrix block.
Definition istlexception.hh:52
int c
Definition istlexception.hh:54
int r
Definition istlexception.hh:54
Definition matrixutils.hh:27