3#ifndef DUNE_ISTL_ILU_HH
4#define DUNE_ISTL_ILU_HH
16#include "istlexception.hh"
35 typedef typename M::RowIterator rowiterator;
36 typedef typename M::ColIterator coliterator;
37 typedef typename M::block_type block;
40 rowiterator endi=A.end();
41 for (rowiterator i=A.begin(); i!=endi; ++i)
44 coliterator endij=(*i).end();
48 for (ij=(*i).begin(); ij.index()<i.index(); ++ij)
51 coliterator jj = A[ij.index()].find(ij.index());
54 Impl::asMatrix(*ij).rightmultiply(Impl::asMatrix(*jj));
57 coliterator endjk=A[ij.index()].end();
58 coliterator jk=jj; ++jk;
59 coliterator ik=ij; ++ik;
60 while (ik!=endij && jk!=endjk)
61 if (ik.index()==jk.index())
64 Impl::asMatrix(B).leftmultiply(Impl::asMatrix(*ij));
70 if (ik.index()<jk.index())
78 if (ij.index()!=i.index())
81 Impl::asMatrix(*ij).invert();
85 << i.index() <<
"][" << ij.index() <<
"]" << e.what();
86 th__ex.r=i.index(); th__ex.c=ij.index(););
92 template<
class M,
class X,
class Y>
96 typedef typename M::ConstRowIterator rowiterator;
97 typedef typename M::ConstColIterator coliterator;
98 typedef typename Y::block_type dblock;
99 typedef typename X::block_type vblock;
102 rowiterator endi=A.end();
103 for (rowiterator i=A.begin(); i!=endi; ++i)
111 dblock rhsValue(d[i.index()]);
112 auto&& rhs = Impl::asVector(rhsValue);
113 for (coliterator j=(*i).begin(); j.index()<i.index(); ++j)
114 Impl::asMatrix(*j).mmv(Impl::asVector(v[j.index()]),rhs);
115 Impl::asVector(v[i.index()]) = rhs;
119 rowiterator rendi=A.beforeBegin();
120 for (rowiterator i=A.beforeEnd(); i!=rendi; --i)
128 vblock rhsValue(v[i.index()]);
129 auto&& rhs = Impl::asVector(rhsValue);
131 for (j=(*i).beforeEnd(); j.index()>i.index(); --j)
132 Impl::asMatrix(*j).mmv(Impl::asVector(v[j.index()]),rhs);
133 auto&& vi = Impl::asVector(v[i.index()]);
134 Impl::asMatrix(*j).mv(rhs,vi);
140 typename M::field_type& firstMatrixElement (M& A,
143 return firstMatrixElement(*(A.begin()->begin()));
147 K& firstMatrixElement (K& A,
153 template<
class K,
int n,
int m>
154 K& firstMatrixElement (FieldMatrix<K,n,m>& A)
169 typedef typename M::ColIterator coliterator;
170 typedef typename M::ConstRowIterator crowiterator;
171 typedef typename M::ConstColIterator ccoliterator;
172 typedef typename M::CreateIterator createiterator;
173 typedef typename M::field_type K;
174 typedef std::map<size_t, int> map;
175 typedef typename map::iterator mapiterator;
178 crowiterator endi=A.end();
179 createiterator ci=ILU.createbegin();
180 for (crowiterator i=A.begin(); i!=endi; ++i)
185 for (ccoliterator j=(*i).begin(); j!=(*i).end(); ++j)
186 rowpattern[j.index()] = 0;
189 for (mapiterator ik=rowpattern.begin(); (*ik).first<i.index(); ++ik)
193 coliterator endk = ILU[(*ik).first].end();
194 coliterator kj = ILU[(*ik).first].find((*ik).first);
195 for (++kj; kj!=endk; ++kj)
200 int generation = (int)
Simd::lane(0, abs( firstMatrixElement(*kj) ));
203 mapiterator ij = rowpattern.find(kj.index());
204 if (ij==rowpattern.end())
206 rowpattern[kj.index()] = generation+1;
214 for (mapiterator ik=rowpattern.begin(); ik!=rowpattern.end(); ++ik)
215 ci.insert((*ik).first);
219 coliterator endILUij = ILU[i.index()].end();;
220 for (coliterator ILUij=ILU[i.index()].begin(); ILUij!=endILUij; ++ILUij)
225 for (crowiterator i=A.begin(); i!=endi; ++i)
228 coliterator endILUij = ILU[i.index()].end();;
229 for (ILUij=ILU[i.index()].begin(); ILUij!=endILUij; ++ILUij)
231 ccoliterator Aij = (*i).begin();
232 ccoliterator endAij = (*i).end();
233 ILUij = ILU[i.index()].begin();
234 while (Aij!=endAij && ILUij!=endILUij)
236 if (Aij.index()==ILUij.index())
243 if (Aij.index()<ILUij.index())
256 template <
class B,
class Alloc = std::allocator<B>>
259 typedef B block_type;
260 typedef size_t size_type;
262 CRS() : nRows_( 0 ) {}
264 size_type rows()
const {
return nRows_; }
266 size_type nonZeros()
const
268 assert( rows_[ rows() ] != size_type(-1) );
269 return rows_[ rows() ];
272 void resize(
const size_type nRows )
274 if( nRows_ != nRows )
277 rows_.resize( nRows_+1, size_type(-1) );
281 void reserveAdditional(
const size_type nonZeros )
283 const size_type needed = values_.size() + nonZeros ;
284 if( values_.capacity() < needed )
286 const size_type estimate = needed * 1.1;
287 values_.reserve( estimate );
288 cols_.reserve( estimate );
292 void push_back(
const block_type& value,
const size_type index )
294 values_.push_back( value );
295 cols_.push_back( index );
298 std::vector< size_type > rows_;
299 std::vector< block_type, Alloc> values_;
300 std::vector< size_type > cols_;
305 template<
class M,
class CRS,
class InvVector>
308 typedef typename M :: size_type size_type;
310 lower.resize( A.N() );
311 upper.resize( A.N() );
315 const size_t memEstimate = (A.nonzeroes() - A.N())/2;
317 assert( A.nonzeroes() != 0 );
318 lower.reserveAdditional( memEstimate );
319 upper.reserveAdditional( memEstimate );
321 const auto endi = A.end();
323 size_type colcount = 0;
324 lower.rows_[ 0 ] = colcount;
325 for (
auto i=A.begin(); i!=endi; ++i, ++row)
327 const size_type iIndex = i.index();
330 for (
auto j=(*i).begin(); j.index() < iIndex; ++j )
332 lower.push_back( (*j), j.index() );
335 lower.rows_[ iIndex+1 ] = colcount;
338 const auto rendi = A.beforeBegin();
341 upper.rows_[ 0 ] = colcount ;
345 for (
auto i=A.beforeEnd(); i!=rendi; --i, ++ row )
347 const auto endij=(*i).beforeBegin();
349 const size_type iIndex = i.index();
352 for (
auto j=(*i).beforeEnd(); j != endij; --j )
354 const size_type jIndex = j.index();
355 if( j.index() == iIndex )
360 else if ( j.index() >= i.index() )
362 upper.push_back( (*j), jIndex );
366 upper.rows_[ row+1 ] = colcount;
370 template<
class CRS,
class InvVector,
class X,
class Y>
376 template<
class CRS,
class InvVector,
class X,
class Y>
379 const InvVector& inv,
383 typedef typename Y :: block_type dblock;
384 typedef typename X :: block_type vblock;
385 typedef typename X :: size_type size_type ;
387 const size_type iEnd = lower.rows();
388 const size_type lastRow = iEnd - 1;
389 if( iEnd != upper.rows() )
395 for( size_type i=0; i<iEnd; ++ i )
397 dblock rhsValue( d[ i ] );
398 auto&& rhs = Impl::asVector(rhsValue);
399 const size_type rowI = lower.rows_[ i ];
400 const size_type rowINext = lower.rows_[ i+1 ];
402 for( size_type col = rowI; col < rowINext; ++ col )
403 Impl::asMatrix(lower.values_[ col ]).mmv( Impl::asVector(v[ lower.cols_[ col ] ] ), rhs );
405 Impl::asVector(v[ i ]) = rhs;
409 for( size_type i=0; i<iEnd; ++ i )
411 auto&& vBlock = Impl::asVector(v[ lastRow - i ]);
412 vblock rhsValue ( v[ lastRow - i ] );
413 auto&& rhs = Impl::asVector(rhsValue);
414 const size_type rowI = upper.rows_[ i ];
415 const size_type rowINext = upper.rows_[ i+1 ];
417 for( size_type col = rowI; col < rowINext; ++ col )
418 Impl::asMatrix(upper.values_[ col ]).mmv( Impl::asVector(v[ upper.cols_[ col ] ]), rhs );
421 Impl::asMatrix(inv[ i ]).mv(rhs, vBlock);
433 template<
class M,
class X,
class Y>
435 void
bilu_backsolve (const M& A, X& v, const Y& d) { ILU::blockILUBacksolve(A, v, d); }
443 void
bilu_decomposition (const M& A,
int n, M& ilu) {
return ILU::blockILUDecomposition(A, n, ilu); }
Error thrown if operations of a FieldMatrix fail.
Definition: densematrix.hh:151
derive error class from the base class in common
Definition: istlexception.hh:17
Error when performing an operation on a matrix block.
Definition: istlexception.hh:59
Definition of the DUNE_DEPRECATED macro for the case that config.h is not available.
Implements a matrix constructed from a given type representing a field and compile-time given number ...
#define DUNE_DEPRECATED_MSG(text)
Mark some entity as deprecated.
Definition: deprecated.hh:32
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
decltype(auto) lane(std::size_t l, V &&v)
Extract an element of a SIMD type.
Definition: interface.hh:322
typename Overloads::ScalarType< std::decay_t< V > >::type Scalar
Element type of some SIMD type.
Definition: interface.hh:233
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:306
void blockILUBacksolve(const M &A, X &v, const Y &d)
LU backsolve with stored inverse.
Definition: ilu.hh:93
void blockILU0Decomposition(M &A)
compute ILU decomposition of A. A is overwritten by its decomposition
Definition: ilu.hh:32
void blockILUDecomposition(const M &A, int n, M &ILU)
Definition: ilu.hh:166
Dune namespace.
Definition: alignedallocator.hh:11
decltype(auto) firstmatrixelement(M &A)
Definition: ilu.hh:439
void bilu_backsolve(const M &A, X &v, const Y &d)
Definition: ilu.hh:435
void bilu_decomposition(const M &A, int n, M &ilu)
Definition: ilu.hh:443
void bilu0_decomposition(M &A)
Definition: ilu.hh:431
Implements a scalar matrix view wrapper around an existing scalar.
Implements a scalar vector view wrapper around an existing scalar.
a simple compressed row storage matrix class
Definition: ilu.hh:258
Whether this type acts as a scalar in the context of (hierarchically blocked) containers.
Definition: typetraits.hh:194