Dune Core Modules (2.6.0)

preconditioners.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_PRECONDITIONERS_HH
4 #define DUNE_ISTL_PRECONDITIONERS_HH
5 
6 #include <cmath>
7 #include <complex>
8 #include <iostream>
9 #include <iomanip>
10 #include <memory>
11 #include <string>
12 
13 #include <dune/common/unused.hh>
14 
15 #include "preconditioner.hh"
16 #include "solver.hh"
17 #include "solvercategory.hh"
18 #include "istlexception.hh"
19 #include "matrixutils.hh"
20 #include "gsetc.hh"
21 #include "ildl.hh"
22 #include "ilu.hh"
23 
24 
25 namespace Dune {
68  template<class O, int c = -1>
70  public Preconditioner<typename O::domain_type, typename O::range_type>
71  {
72  public:
74  typedef typename O::domain_type domain_type;
76  typedef typename O::range_type range_type;
78  typedef typename range_type::field_type field_type;
80  typedef SimdScalar<field_type> scalar_field_type;
82  typedef O InverseOperator;
83 
89  : inverse_operator_(inverse_operator)
90  {
91  if(c != -1 && SolverCategory::category(inverse_operator_) != c)
92  DUNE_THROW(InvalidStateException, "User supplied solver category does not match that of the supplied iverser operator");
93  }
94 
95  virtual void pre(domain_type&,range_type&)
96  {}
97 
98  virtual void apply(domain_type& v, const range_type& d)
99  {
101  range_type copy(d);
102  inverse_operator_.apply(v, copy, res);
103  }
104 
105  virtual void post(domain_type&)
106  {}
107 
110  {
111  return SolverCategory::category(inverse_operator_);
112  }
113 
114  private:
115  InverseOperator& inverse_operator_;
116  };
117 
118  //=====================================================================
119  // Implementation of this interface for sequential ISTL-preconditioners
120  //=====================================================================
121 
122 
134  template<class M, class X, class Y, int l=1>
135  class SeqSSOR : public Preconditioner<X,Y> {
136  public:
138  typedef M matrix_type;
140  typedef X domain_type;
142  typedef Y range_type;
144  typedef typename X::field_type field_type;
146  typedef SimdScalar<field_type> scalar_field_type;
147 
155  SeqSSOR (const M& A, int n, scalar_field_type w)
156  : _A_(A), _n(n), _w(w)
157  {
159  }
160 
166  virtual void pre (X& x, Y& b)
167  {
170 
171  }
172 
178  virtual void apply (X& v, const Y& d)
179  {
180  for (int i=0; i<_n; i++) {
181  bsorf(_A_,v,d,_w,BL<l>());
182  bsorb(_A_,v,d,_w,BL<l>());
183  }
184  }
185 
191  virtual void post (X& x)
192  {
194  }
195 
198  {
200  }
201 
202  private:
204  const M& _A_;
206  int _n;
209  };
210 
211 
212 
224  template<class M, class X, class Y, int l=1>
225  class SeqSOR : public Preconditioner<X,Y> {
226  public:
228  typedef M matrix_type;
230  typedef X domain_type;
232  typedef Y range_type;
234  typedef typename X::field_type field_type;
236  typedef SimdScalar<field_type> scalar_field_type;
237 
245  SeqSOR (const M& A, int n, scalar_field_type w)
246  : _A_(A), _n(n), _w(w)
247  {
249  }
250 
256  virtual void pre (X& x, Y& b)
257  {
260  }
261 
267  virtual void apply (X& v, const Y& d)
268  {
269  this->template apply<true>(v,d);
270  }
271 
280  template<bool forward>
281  void apply(X& v, const Y& d)
282  {
283  if(forward)
284  for (int i=0; i<_n; i++) {
285  bsorf(_A_,v,d,_w,BL<l>());
286  }
287  else
288  for (int i=0; i<_n; i++) {
289  bsorb(_A_,v,d,_w,BL<l>());
290  }
291  }
292 
298  virtual void post (X& x)
299  {
301  }
302 
305  {
307  }
308 
309  private:
311  const M& _A_;
313  int _n;
316  };
317 
318 
329  template<class M, class X, class Y, int l=1>
330  class SeqGS : public Preconditioner<X,Y> {
331  public:
333  typedef M matrix_type;
335  typedef X domain_type;
337  typedef Y range_type;
339  typedef typename X::field_type field_type;
341  typedef SimdScalar<field_type> scalar_field_type;
342 
350  SeqGS (const M& A, int n, scalar_field_type w)
351  : _A_(A), _n(n), _w(w)
352  {
354  }
355 
361  virtual void pre (X& x, Y& b)
362  {
365  }
366 
372  virtual void apply (X& v, const Y& d)
373  {
374  for (int i=0; i<_n; i++) {
375  dbgs(_A_,v,d,_w,BL<l>());
376  }
377  }
378 
384  virtual void post (X& x)
385  {
387  }
388 
391  {
393  }
394 
395  private:
397  const M& _A_;
399  int _n;
401  scalar_field_type _w;
402  };
403 
404 
415  template<class M, class X, class Y, int l=1>
416  class SeqJac : public Preconditioner<X,Y> {
417  public:
419  typedef M matrix_type;
421  typedef X domain_type;
423  typedef Y range_type;
425  typedef typename X::field_type field_type;
427  typedef SimdScalar<field_type> scalar_field_type;
428 
436  SeqJac (const M& A, int n, scalar_field_type w)
437  : _A_(A), _n(n), _w(w)
438  {
440  }
441 
447  virtual void pre (X& x, Y& b)
448  {
451  }
452 
458  virtual void apply (X& v, const Y& d)
459  {
460  for (int i=0; i<_n; i++) {
461  dbjac(_A_,v,d,_w,BL<l>());
462  }
463  }
464 
470  virtual void post (X& x)
471  {
473  }
474 
477  {
479  }
480 
481  private:
483  const M& _A_;
485  int _n;
487  scalar_field_type _w;
488  };
489 
490 
491 
503  template<class M, class X, class Y, int l=1>
504  class SeqILU : public Preconditioner<X,Y> {
505  public:
507  typedef typename std::remove_const<M>::type matrix_type;
509  typedef typename matrix_type :: block_type block_type;
511  typedef X domain_type;
513  typedef Y range_type;
514 
516  typedef typename X::field_type field_type;
517 
519  typedef SimdScalar<field_type> scalar_field_type;
520 
522  typedef typename ILU::CRS< block_type > CRS;
523 
531  SeqILU (const M& A, scalar_field_type w, const bool resort = false )
532  : SeqILU( A, 0, w, resort ) // construct ILU(0)
533  {
534  }
535 
544  SeqILU (const M& A, int n, scalar_field_type w, const bool resort = false )
545  : ILU_(),
546  lower_(),
547  upper_(),
548  inv_(),
549  w_(w),
550  wNotIdentity_( std::abs( w_ - scalar_field_type(1) ) > 1e-15 )
551  {
552  if( n == 0 )
553  {
554  // copy A
555  ILU_.reset( new matrix_type( A ) );
556  // create ILU(0) decomposition
557  bilu0_decomposition( *ILU_ );
558  }
559  else
560  {
561  // create matrix in build mode
562  ILU_.reset( new matrix_type( A.N(), A.M(), matrix_type::row_wise) );
563  // create ILU(n) decomposition
564  bilu_decomposition( A, n, *ILU_ );
565  }
566 
567  if( resort )
568  {
569  // store ILU in simple CRS format
570  ILU::convertToCRS( *ILU_, lower_, upper_, inv_ );
571  ILU_.reset();
572  }
573  }
574 
580  virtual void pre (X& x, Y& b)
581  {
584  }
585 
591  virtual void apply (X& v, const Y& d)
592  {
593  if( ILU_ )
594  {
595  bilu_backsolve( *ILU_, v, d);
596  }
597  else
598  {
599  ILU::bilu_backsolve(lower_, upper_, inv_, v, d);
600  }
601 
602  if( wNotIdentity_ )
603  {
604  v *= w_;
605  }
606  }
607 
613  virtual void post (X& x)
614  {
616  }
617 
620  {
622  }
623 
624  protected:
626  std::unique_ptr< matrix_type > ILU_;
627 
630  CRS upper_;
631  std::vector< block_type > inv_;
632 
636  const bool wNotIdentity_;
637  };
638 
639 
651  template<class M, class X, class Y, int l=1>
652  class SeqILU0 : public Preconditioner<X,Y> {
653  public:
655  typedef typename std::remove_const<M>::type matrix_type;
657  typedef X domain_type;
659  typedef Y range_type;
661  typedef typename X::field_type field_type;
663  typedef SimdScalar<field_type> scalar_field_type;
664 
671  SeqILU0 (const M& A, scalar_field_type w)
672  : _w(w),
673  ILU(A) // copy A
674 
675  {
676  bilu0_decomposition(ILU);
677  }
678 
684  virtual void pre (X& x, Y& b)
685  {
688  }
689 
695  virtual void apply (X& v, const Y& d)
696  {
697  bilu_backsolve(ILU,v,d);
698  v *= _w;
699  }
700 
706  virtual void post (X& x)
707  {
709  }
710 
713  {
715  }
716 
717  private:
719  scalar_field_type _w;
721  matrix_type ILU;
722  };
723 
724 
738  template<class M, class X, class Y, int l=1>
739  class SeqILUn : public Preconditioner<X,Y> {
740  public:
742  typedef typename std::remove_const<M>::type matrix_type;
744  typedef X domain_type;
746  typedef Y range_type;
748  typedef typename X::field_type field_type;
750  typedef SimdScalar<field_type> scalar_field_type;
751 
759  SeqILUn (const M& A, int n, scalar_field_type w)
760  : ILU(A.N(),A.M(),M::row_wise),
761  _n(n),
762  _w(w)
763  {
764  bilu_decomposition(A,n,ILU);
765  }
766 
772  virtual void pre (X& x, Y& b)
773  {
776  }
777 
783  virtual void apply (X& v, const Y& d)
784  {
785  bilu_backsolve(ILU,v,d);
786  v *= _w;
787  }
788 
794  virtual void post (X& x)
795  {
797  }
798 
801  {
803  }
804 
805  private:
807  matrix_type ILU;
809  int _n;
811  scalar_field_type _w;
812  };
813 
814 
815 
824  template<class X, class Y>
826  public:
828  typedef X domain_type;
830  typedef Y range_type;
832  typedef typename X::field_type field_type;
834  typedef SimdScalar<field_type> scalar_field_type;
835 
842  _w(w)
843  {}
844 
850  virtual void pre (X& x, Y& b)
851  {
854  }
855 
861  virtual void apply (X& v, const Y& d)
862  {
863  v = d;
864  v *= _w;
865  }
866 
872  virtual void post (X& x)
873  {
875  }
876 
879  {
881  }
882 
883  private:
885  scalar_field_type _w;
886  };
887 
888 
899  template< class M, class X, class Y >
900  class SeqILDL
901  : public Preconditioner< X, Y >
902  {
903  typedef SeqILDL< M, X, Y > This;
905 
906  public:
908  typedef std::remove_const_t< M > matrix_type;
910  typedef X domain_type;
912  typedef Y range_type;
914  typedef typename X::field_type field_type;
916  typedef SimdScalar<field_type> scalar_field_type;
917 
926  explicit SeqILDL ( const matrix_type &A, scalar_field_type relax = scalar_field_type( 1 ) )
927  : decomposition_( A.N(), A.M(), matrix_type::random ),
928  relax_( relax )
929  {
930  // setup row sizes for lower triangular matrix
931  for( auto i = A.begin(), iend = A.end(); i != iend; ++i )
932  {
933  const auto &A_i = *i;
934  const auto ij = A_i.find( i.index() );
935  if( ij != A_i.end() )
936  decomposition_.setrowsize( i.index(), ij.offset()+1 );
937  else
938  DUNE_THROW( ISTLError, "diagonal entry missing" );
939  }
940  decomposition_.endrowsizes();
941 
942  // setup row indices for lower triangular matrix
943  for( auto i = A.begin(), iend = A.end(); i != iend; ++i )
944  {
945  const auto &A_i = *i;
946  for( auto ij = A_i.begin(); ij.index() < i.index() ; ++ij )
947  decomposition_.addindex( i.index(), ij.index() );
948  decomposition_.addindex( i.index(), i.index() );
949  }
950  decomposition_.endindices();
951 
952  // copy values of lower triangular matrix
953  auto i = A.begin();
954  for( auto row = decomposition_.begin(), rowend = decomposition_.end(); row != rowend; ++row, ++i )
955  {
956  auto ij = i->begin();
957  for( auto col = row->begin(), colend = row->end(); col != colend; ++col, ++ij )
958  *col = *ij;
959  }
960 
961  // perform ILDL decomposition
962  bildl_decompose( decomposition_ );
963  }
964 
966  void pre ( X &x, Y &b ) override
967  {
970  }
971 
973  void apply ( X &v, const Y &d ) override
974  {
975  bildl_backsolve( decomposition_, v, d, true );
976  v *= relax_;
977  }
978 
980  void post ( X &x ) override
981  {
983  }
984 
987 
988  private:
989  matrix_type decomposition_;
990  scalar_field_type relax_;
991  };
992 
995 } // end namespace
996 
997 #endif
derive error class from the base class in common
Definition: istlexception.hh:16
Default exception if a function was called while the object is not in a valid state for that function...
Definition: exceptions.hh:279
Turns an InverseOperator into a Preconditioner.
Definition: preconditioners.hh:71
O::range_type range_type
The range type of the preconditioner.
Definition: preconditioners.hh:76
O::domain_type domain_type
The domain type of the preconditioner.
Definition: preconditioners.hh:74
virtual void post(domain_type &)
Clean up.
Definition: preconditioners.hh:105
range_type::field_type field_type
The field type of the preconditioner.
Definition: preconditioners.hh:78
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: preconditioners.hh:109
virtual void pre(domain_type &, range_type &)
Prepare the preconditioner.
Definition: preconditioners.hh:95
SimdScalar< field_type > scalar_field_type
scalar type underlying the field_type
Definition: preconditioners.hh:80
InverseOperator2Preconditioner(InverseOperator &inverse_operator)
Construct the preconditioner from the solver.
Definition: preconditioners.hh:88
O InverseOperator
type of the wrapped inverse operator
Definition: preconditioners.hh:82
virtual void apply(domain_type &v, const range_type &d)
Apply one step of the preconditioner to the system A(v)=d.
Definition: preconditioners.hh:98
Abstract base class for all solvers.
Definition: solver.hh:91
RowIterator end()
Get iterator to one beyond last row.
Definition: matrix.hh:615
RowIterator begin()
Get iterator to first row.
Definition: matrix.hh:609
size_type M() const
Return the number of columns.
Definition: matrix.hh:695
size_type N() const
Return the number of rows.
Definition: matrix.hh:690
Base class for matrix free definition of preconditioners.
Definition: preconditioner.hh:30
Richardson preconditioner.
Definition: preconditioners.hh:825
X::field_type field_type
The field type of the preconditioner.
Definition: preconditioners.hh:832
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: preconditioners.hh:878
Y range_type
The range type of the preconditioner.
Definition: preconditioners.hh:830
virtual void pre(X &x, Y &b)
Prepare the preconditioner.
Definition: preconditioners.hh:850
SimdScalar< field_type > scalar_field_type
scalar type underlying the field_type
Definition: preconditioners.hh:834
virtual void post(X &x)
Clean up.
Definition: preconditioners.hh:872
Richardson(scalar_field_type w=1.0)
Constructor.
Definition: preconditioners.hh:841
X domain_type
The domain type of the preconditioner.
Definition: preconditioners.hh:828
virtual void apply(X &v, const Y &d)
Apply the precondioner.
Definition: preconditioners.hh:861
Sequential Gauss Seidel preconditioner.
Definition: preconditioners.hh:330
virtual void post(X &x)
Clean up.
Definition: preconditioners.hh:384
virtual void pre(X &x, Y &b)
Prepare the preconditioner.
Definition: preconditioners.hh:361
SeqGS(const M &A, int n, scalar_field_type w)
Constructor.
Definition: preconditioners.hh:350
Y range_type
The range type of the preconditioner.
Definition: preconditioners.hh:337
SimdScalar< field_type > scalar_field_type
scalar type underlying the field_type
Definition: preconditioners.hh:341
X::field_type field_type
The field type of the preconditioner.
Definition: preconditioners.hh:339
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: preconditioners.hh:390
virtual void apply(X &v, const Y &d)
Apply the preconditioner.
Definition: preconditioners.hh:372
M matrix_type
The matrix type the preconditioner is for.
Definition: preconditioners.hh:333
X domain_type
The domain type of the preconditioner.
Definition: preconditioners.hh:335
sequential ILDL preconditioner
Definition: preconditioners.hh:902
SeqILDL(const matrix_type &A, scalar_field_type relax=scalar_field_type(1))
constructor
Definition: preconditioners.hh:926
X domain_type
domain type of the preconditioner
Definition: preconditioners.hh:910
void post(X &x) override
Clean up.
Definition: preconditioners.hh:980
Y range_type
range type of the preconditioner
Definition: preconditioners.hh:912
std::remove_const_t< M > matrix_type
type of matrix the preconditioner is for
Definition: preconditioners.hh:908
SimdScalar< field_type > scalar_field_type
scalar type underlying the field_type
Definition: preconditioners.hh:916
void apply(X &v, const Y &d) override
Apply one step of the preconditioner to the system A(v)=d.
Definition: preconditioners.hh:973
void pre(X &x, Y &b) override
Prepare the preconditioner.
Definition: preconditioners.hh:966
X::field_type field_type
field type of the preconditioner
Definition: preconditioners.hh:914
SolverCategory::Category category() const override
Category of the preconditioner (see SolverCategory::Category)
Definition: preconditioners.hh:986
Sequential ILU0 preconditioner.
Definition: preconditioners.hh:652
virtual void post(X &x)
Clean up.
Definition: preconditioners.hh:706
virtual void pre(X &x, Y &b)
Prepare the preconditioner.
Definition: preconditioners.hh:684
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: preconditioners.hh:712
virtual void apply(X &v, const Y &d)
Apply the preconditoner.
Definition: preconditioners.hh:695
Y range_type
The range type of the preconditioner.
Definition: preconditioners.hh:659
SeqILU0(const M &A, scalar_field_type w)
Constructor.
Definition: preconditioners.hh:671
X domain_type
The domain type of the preconditioner.
Definition: preconditioners.hh:657
X::field_type field_type
The field type of the preconditioner.
Definition: preconditioners.hh:661
std::remove_const< M >::type matrix_type
The matrix type the preconditioner is for.
Definition: preconditioners.hh:655
SimdScalar< field_type > scalar_field_type
scalar type underlying the field_type
Definition: preconditioners.hh:663
Sequential ILU preconditioner.
Definition: preconditioners.hh:504
virtual void post(X &x)
Clean up.
Definition: preconditioners.hh:613
virtual void pre(X &x, Y &b)
Prepare the preconditioner.
Definition: preconditioners.hh:580
matrix_type ::block_type block_type
block type of matrix
Definition: preconditioners.hh:509
virtual void apply(X &v, const Y &d)
Apply the preconditoner.
Definition: preconditioners.hh:591
Y range_type
The range type of the preconditioner.
Definition: preconditioners.hh:513
const scalar_field_type w_
The relaxation factor to use.
Definition: preconditioners.hh:634
CRS lower_
The ILU(n) decomposition of the matrix. As storage a CRS structure is used.
Definition: preconditioners.hh:629
const bool wNotIdentity_
true if w != 1.0
Definition: preconditioners.hh:636
SeqILU(const M &A, int n, scalar_field_type w, const bool resort=false)
Constructor.
Definition: preconditioners.hh:544
std::remove_const< M >::type matrix_type
The matrix type the preconditioner is for.
Definition: preconditioners.hh:507
SeqILU(const M &A, scalar_field_type w, const bool resort=false)
Constructor.
Definition: preconditioners.hh:531
SimdScalar< field_type > scalar_field_type
scalar type underlying the field_type
Definition: preconditioners.hh:519
X::field_type field_type
The field type of the preconditioner.
Definition: preconditioners.hh:516
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: preconditioners.hh:619
X domain_type
The domain type of the preconditioner.
Definition: preconditioners.hh:511
ILU::CRS< block_type > CRS
type of ILU storage
Definition: preconditioners.hh:522
std::unique_ptr< matrix_type > ILU_
The ILU(n) decomposition of the matrix. As storage a BCRSMatrix is used.
Definition: preconditioners.hh:626
Sequential ILU(n) preconditioner.
Definition: preconditioners.hh:739
SeqILUn(const M &A, int n, scalar_field_type w)
Constructor.
Definition: preconditioners.hh:759
virtual void pre(X &x, Y &b)
Prepare the preconditioner.
Definition: preconditioners.hh:772
std::remove_const< M >::type matrix_type
The matrix type the preconditioner is for.
Definition: preconditioners.hh:742
virtual void post(X &x)
Clean up.
Definition: preconditioners.hh:794
SimdScalar< field_type > scalar_field_type
scalar type underlying the field_type
Definition: preconditioners.hh:750
X domain_type
The domain type of the preconditioner.
Definition: preconditioners.hh:744
Y range_type
The range type of the preconditioner.
Definition: preconditioners.hh:746
virtual void apply(X &v, const Y &d)
Apply the precondioner.
Definition: preconditioners.hh:783
X::field_type field_type
The field type of the preconditioner.
Definition: preconditioners.hh:748
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: preconditioners.hh:800
The sequential jacobian preconditioner.
Definition: preconditioners.hh:416
virtual void post(X &x)
Clean up.
Definition: preconditioners.hh:470
virtual void apply(X &v, const Y &d)
Apply the preconditioner.
Definition: preconditioners.hh:458
M matrix_type
The matrix type the preconditioner is for.
Definition: preconditioners.hh:419
X::field_type field_type
The field type of the preconditioner.
Definition: preconditioners.hh:425
SeqJac(const M &A, int n, scalar_field_type w)
Constructor.
Definition: preconditioners.hh:436
virtual void pre(X &x, Y &b)
Prepare the preconditioner.
Definition: preconditioners.hh:447
X domain_type
The domain type of the preconditioner.
Definition: preconditioners.hh:421
SimdScalar< field_type > scalar_field_type
scalar type underlying the field_type
Definition: preconditioners.hh:427
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: preconditioners.hh:476
Y range_type
The range type of the preconditioner.
Definition: preconditioners.hh:423
Sequential SOR preconditioner.
Definition: preconditioners.hh:225
M matrix_type
The matrix type the preconditioner is for.
Definition: preconditioners.hh:228
SimdScalar< field_type > scalar_field_type
scalar type underlying the field_type
Definition: preconditioners.hh:236
void apply(X &v, const Y &d)
Apply the preconditioner in a special direction.
Definition: preconditioners.hh:281
X domain_type
The domain type of the preconditioner.
Definition: preconditioners.hh:230
virtual void post(X &x)
Clean up.
Definition: preconditioners.hh:298
SeqSOR(const M &A, int n, scalar_field_type w)
Constructor.
Definition: preconditioners.hh:245
virtual void pre(X &x, Y &b)
Prepare the preconditioner.
Definition: preconditioners.hh:256
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: preconditioners.hh:304
virtual void apply(X &v, const Y &d)
Apply the preconditioner.
Definition: preconditioners.hh:267
Y range_type
The range type of the preconditioner.
Definition: preconditioners.hh:232
X::field_type field_type
The field type of the preconditioner.
Definition: preconditioners.hh:234
Sequential SSOR preconditioner.
Definition: preconditioners.hh:135
virtual void post(X &x)
Clean up.
Definition: preconditioners.hh:191
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: preconditioners.hh:197
SimdScalar< field_type > scalar_field_type
scalar type underlying the field_type
Definition: preconditioners.hh:146
X::field_type field_type
The field type of the preconditioner.
Definition: preconditioners.hh:144
SeqSSOR(const M &A, int n, scalar_field_type w)
Constructor.
Definition: preconditioners.hh:155
X domain_type
The domain type of the preconditioner.
Definition: preconditioners.hh:140
M matrix_type
The matrix type the preconditioner is for.
Definition: preconditioners.hh:138
virtual void apply(X &v, const Y &d)
Apply the preconditioner.
Definition: preconditioners.hh:178
virtual void pre(X &x, Y &b)
Prepare the preconditioner.
Definition: preconditioners.hh:166
Y range_type
The range type of the preconditioner.
Definition: preconditioners.hh:142
#define DUNE_UNUSED_PARAMETER(parm)
A macro to mark intentionally unused function parameters with.
Definition: unused.hh:25
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
void bsorb(const M &A, X &x, const Y &b, const K &w)
SSOR step.
Definition: gsetc.hh:591
void dbjac(const M &A, X &x, const Y &b, const K &w)
Jacobi step.
Definition: gsetc.hh:603
void dbgs(const M &A, X &x, const Y &b, const K &w)
GS step.
Definition: gsetc.hh:567
void bilu_backsolve(const M &A, X &v, const Y &d)
LU backsolve with stored inverse.
Definition: ilu.hh:97
void bilu_decomposition(const M &A, int n, M &ILU)
Definition: ilu.hh:157
void bilu0_decomposition(M &A)
compute ILU decomposition of A. A is overwritten by its decomposition
Definition: ilu.hh:36
void bsorf(const M &A, X &x, const Y &b, const K &w)
SOR step.
Definition: gsetc.hh:579
Simple iterative methods like Jacobi, Gauss-Seidel, SOR, SSOR, etc. in a generic way.
Incomplete LDL decomposition.
Some handy generic functions for ISTL matrices.
Dune namespace.
Definition: alignedallocator.hh:10
void bildl_decompose(Matrix &A)
compute ILDL decomposition of a symmetric matrix A
Definition: ildl.hh:76
Define general, extensible interface for inverse operators.
compile-time parameter for block recursion depth
Definition: gsetc.hh:40
static void check(const Matrix &mat)
Check whether the a matrix has diagonal values on blocklevel recursion levels.
Definition: matrixutils.hh:88
Statistics about the application of an inverse operator.
Definition: solver.hh:41
Category
Definition: solvercategory.hh:21
@ sequential
Category for sequential solvers.
Definition: solvercategory.hh:23
static Category category(const OP &op, decltype(op.category()) *=nullptr)
Helperfunction to extract the solver category either from an enum, or from the newly introduced virtu...
Definition: solvercategory.hh:32
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.80.0 (May 2, 22:35, 2024)