5#ifndef DUNE_ISTL_ILU_HH
6#define DUNE_ISTL_ILU_HH
13#include <dune/common/fmatrix.hh>
14#include <dune/common/scalarvectorview.hh>
15#include <dune/common/scalarmatrixview.hh>
36 typedef typename M::RowIterator rowiterator;
37 typedef typename M::ColIterator coliterator;
38 typedef typename M::block_type block;
41 rowiterator endi=A.end();
42 for (rowiterator i=A.begin(); i!=endi; ++i)
45 coliterator endij=(*i).end();
49 for (ij=(*i).begin(); ij.index()<i.index(); ++ij)
52 coliterator jj = A[ij.index()].find(ij.index());
55 Impl::asMatrix(*ij).rightmultiply(Impl::asMatrix(*jj));
58 coliterator endjk=A[ij.index()].end();
59 coliterator jk=jj; ++jk;
60 coliterator ik=ij; ++ik;
61 while (ik!=endij && jk!=endjk)
62 if (ik.index()==jk.index())
65 Impl::asMatrix(B).leftmultiply(Impl::asMatrix(*ij));
71 if (ik.index()<jk.index())
79 if (ij.index()!=i.index())
80 DUNE_THROW(
ISTLError,
"diagonal entry missing");
82 Impl::asMatrix(*ij).invert();
84 catch (Dune::FMatrixError & e) {
86 << i.index() <<
"][" << ij.index() <<
"]" << e.what();
87 th__ex.
r=i.index(); th__ex.c=ij.index(););
93 template<
class M,
class X,
class Y>
97 typedef typename M::ConstRowIterator rowiterator;
98 typedef typename M::ConstColIterator coliterator;
99 typedef typename Y::block_type dblock;
100 typedef typename X::block_type vblock;
103 rowiterator endi=A.end();
104 for (rowiterator i=A.begin(); i!=endi; ++i)
112 dblock rhsValue(d[i.index()]);
113 auto&& rhs = Impl::asVector(rhsValue);
114 for (coliterator j=(*i).begin(); j.index()<i.index(); ++j)
115 Impl::asMatrix(*j).mmv(Impl::asVector(v[j.index()]),rhs);
116 Impl::asVector(v[i.index()]) = rhs;
120 rowiterator rendi=A.beforeBegin();
121 for (rowiterator i=A.beforeEnd(); i!=rendi; --i)
129 vblock rhsValue(v[i.index()]);
130 auto&& rhs = Impl::asVector(rhsValue);
132 for (j=(*i).beforeEnd(); j.index()>i.index(); --j)
133 Impl::asMatrix(*j).mmv(Impl::asVector(v[j.index()]),rhs);
134 auto&& vi = Impl::asVector(v[i.index()]);
135 Impl::asMatrix(*j).mv(rhs,vi);
142 [[maybe_unused]]
typename std::enable_if_t<!Dune::IsNumber<M>::value>* sfinae =
nullptr)
149 [[maybe_unused]]
typename std::enable_if_t<Dune::IsNumber<K>::value>* sfinae =
nullptr)
154 template<
class K,
int n,
int m>
170 typedef typename M::ColIterator coliterator;
171 typedef typename M::ConstRowIterator crowiterator;
172 typedef typename M::ConstColIterator ccoliterator;
173 typedef typename M::CreateIterator createiterator;
174 typedef typename M::field_type K;
175 typedef std::map<size_t, int> map;
176 typedef typename map::iterator mapiterator;
179 crowiterator endi=A.end();
180 createiterator ci=ILU.createbegin();
181 for (crowiterator i=A.begin(); i!=endi; ++i)
186 for (ccoliterator j=(*i).begin(); j!=(*i).end(); ++j)
187 rowpattern[j.index()] = 0;
190 for (mapiterator ik=rowpattern.begin(); (*ik).first<i.index(); ++ik)
194 coliterator endk = ILU[(*ik).first].end();
195 coliterator kj = ILU[(*ik).first].find((*ik).first);
196 for (++kj; kj!=endk; ++kj)
204 mapiterator ij = rowpattern.find(kj.index());
205 if (ij==rowpattern.end())
207 rowpattern[kj.index()] = generation+1;
215 for (mapiterator ik=rowpattern.begin(); ik!=rowpattern.end(); ++ik)
216 ci.insert((*ik).first);
220 coliterator endILUij = ILU[i.index()].end();;
221 for (coliterator ILUij=ILU[i.index()].begin(); ILUij!=endILUij; ++ILUij)
222 Simd::lane(0,
firstMatrixElement(*ILUij)) = (Simd::Scalar<K>) rowpattern[ILUij.index()];
226 for (crowiterator i=A.begin(); i!=endi; ++i)
229 coliterator endILUij = ILU[i.index()].end();;
230 for (ILUij=ILU[i.index()].begin(); ILUij!=endILUij; ++ILUij)
232 ccoliterator Aij = (*i).begin();
233 ccoliterator endAij = (*i).end();
234 ILUij = ILU[i.index()].begin();
235 while (Aij!=endAij && ILUij!=endILUij)
237 if (Aij.index()==ILUij.index())
244 if (Aij.index()<ILUij.index())
257 template <
class B,
class Alloc = std::allocator<B>>
285 if(
values_.capacity() < needed )
289 cols_.reserve( estimate );
296 cols_.push_back( index );
306 template<
class M,
class CRS,
class InvVector>
309 typedef typename M :: size_type size_type;
316 const size_t memEstimate = (A.nonzeroes() - A.N())/2;
318 assert( A.nonzeroes() != 0 );
322 const auto endi = A.end();
324 size_type colcount = 0;
325 lower.
rows_[ 0 ] = colcount;
326 for (
auto i=A.begin(); i!=endi; ++i, ++row)
328 const size_type iIndex = i.index();
331 for (
auto j=(*i).begin(); j.index() < iIndex; ++j )
336 lower.
rows_[ iIndex+1 ] = colcount;
339 const auto rendi = A.beforeBegin();
342 upper.
rows_[ 0 ] = colcount ;
346 for (
auto i=A.beforeEnd(); i!=rendi; --i, ++ row )
348 const auto endij=(*i).beforeBegin();
350 const size_type iIndex = i.index();
353 for (
auto j=(*i).beforeEnd(); j != endij; --j )
355 const size_type jIndex = j.index();
356 if( j.index() == iIndex )
361 else if ( j.index() >= i.index() )
367 upper.
rows_[ row+1 ] = colcount;
372 template<
class CRS,
class InvVector,
class X,
class Y>
375 const InvVector& inv,
379 typedef typename Y :: block_type dblock;
380 typedef typename X :: block_type vblock;
381 typedef typename X :: size_type size_type ;
383 const size_type iEnd = lower.
rows();
384 const size_type lastRow = iEnd - 1;
385 if( iEnd != upper.
rows() )
387 DUNE_THROW(
ISTLError,
"ILU::blockILUBacksolve: lower and upper rows must be the same");
391 for( size_type i=0; i<iEnd; ++ i )
393 dblock rhsValue( d[ i ] );
394 auto&& rhs = Impl::asVector(rhsValue);
395 const size_type rowI = lower.
rows_[ i ];
396 const size_type rowINext = lower.
rows_[ i+1 ];
398 for( size_type
col = rowI;
col < rowINext; ++
col )
399 Impl::asMatrix(lower.
values_[
col ]).mmv( Impl::asVector(v[ lower.
cols_[
col ] ] ), rhs );
401 Impl::asVector(v[ i ]) = rhs;
405 for( size_type i=0; i<iEnd; ++ i )
407 auto&& vBlock = Impl::asVector(v[ lastRow - i ]);
408 vblock rhsValue ( v[ lastRow - i ] );
409 auto&& rhs = Impl::asVector(rhsValue);
410 const size_type rowI = upper.
rows_[ i ];
411 const size_type rowINext = upper.
rows_[ i+1 ];
413 for( size_type
col = rowI;
col < rowINext; ++
col )
414 Impl::asMatrix(upper.
values_[
col ]).mmv( Impl::asVector(v[ upper.
cols_[
col ] ]), rhs );
417 Impl::asMatrix(inv[ i ]).mv(rhs, vBlock);
Col col
Definition: matrixmatrix.hh:351
Definition: allocator.hh:11
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:307
void blockILUBacksolve(const M &A, X &v, const Y &d)
LU backsolve with stored inverse.
Definition: ilu.hh:94
M::field_type & firstMatrixElement(M &A, typename std::enable_if_t<!Dune::IsNumber< M >::value > *sfinae=nullptr)
Definition: ilu.hh:141
void blockILU0Decomposition(M &A)
compute ILU decomposition of A. A is overwritten by its decomposition
Definition: ilu.hh:33
void blockILUDecomposition(const M &A, int n, M &ILU)
Definition: ilu.hh:167
a simple compressed row storage matrix class
Definition: ilu.hh:259
std::vector< size_type > cols_
Definition: ilu.hh:301
size_type nonZeros() const
Definition: ilu.hh:267
void resize(const size_type nRows)
Definition: ilu.hh:273
size_type rows() const
Definition: ilu.hh:265
CRS()
Definition: ilu.hh:263
void reserveAdditional(const size_type nonZeros)
Definition: ilu.hh:282
B block_type
Definition: ilu.hh:260
std::vector< block_type, Alloc > values_
Definition: ilu.hh:300
size_type nRows_
Definition: ilu.hh:302
size_t size_type
Definition: ilu.hh:261
std::vector< size_type > rows_
Definition: ilu.hh:299
void push_back(const block_type &value, const size_type index)
Definition: ilu.hh:293
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 r
Definition: istlexception.hh:54
Definition: matrixutils.hh:27