Dune Core Modules (2.7.1)

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/simd/simd.hh>
14 #include <dune/common/unused.hh>
16 
17 #include <dune/istl/solverfactory.hh>
18 #include "preconditioner.hh"
19 #include "solver.hh"
20 #include "solvercategory.hh"
21 #include "istlexception.hh"
22 #include "matrixutils.hh"
23 #include "gsetc.hh"
24 #include "ildl.hh"
25 #include "ilu.hh"
26 
27 
28 namespace Dune {
71  template<class O, int c = -1>
73  public Preconditioner<typename O::domain_type, typename O::range_type>
74  {
75  public:
77  typedef typename O::domain_type domain_type;
79  typedef typename O::range_type range_type;
81  typedef typename range_type::field_type field_type;
85  typedef O InverseOperator;
86 
92  : inverse_operator_(inverse_operator)
93  {
94  if(c != -1 && SolverCategory::category(inverse_operator_) != c)
95  DUNE_THROW(InvalidStateException, "User supplied solver category does not match that of the supplied iverser operator");
96  }
97 
98  virtual void pre(domain_type&,range_type&)
99  {}
100 
101  virtual void apply(domain_type& v, const range_type& d)
102  {
104  range_type copy(d);
105  inverse_operator_.apply(v, copy, res);
106  }
107 
108  virtual void post(domain_type&)
109  {}
110 
113  {
114  return SolverCategory::category(inverse_operator_);
115  }
116 
117  private:
118  InverseOperator& inverse_operator_;
119  };
120 
121  //=====================================================================
122  // Implementation of this interface for sequential ISTL-preconditioners
123  //=====================================================================
124 
125 
137  template<class M, class X, class Y, int l=1>
138  class SeqSSOR : public Preconditioner<X,Y> {
139  public:
141  typedef M matrix_type;
143  typedef X domain_type;
145  typedef Y range_type;
147  typedef typename X::field_type field_type;
150 
158  SeqSSOR (const M& A, int n, scalar_field_type w)
159  : _A_(A), _n(n), _w(w)
160  {
162  }
163 
177  SeqSSOR (const M& A, const ParameterTree& configuration)
178  : _A_(A)
179  {
180  _n = configuration.get<int>("iterations",1);
181  _w = configuration.get<scalar_field_type>("relaxation",1.0);
183  }
184 
190  virtual void pre (X& x, Y& b)
191  {
194 
195  }
196 
202  virtual void apply (X& v, const Y& d)
203  {
204  for (int i=0; i<_n; i++) {
205  bsorf(_A_,v,d,_w,BL<l>());
206  bsorb(_A_,v,d,_w,BL<l>());
207  }
208  }
209 
215  virtual void post (X& x)
216  {
218  }
219 
222  {
224  }
225 
226  private:
228  const M& _A_;
230  int _n;
233  };
234  DUNE_REGISTER_PRECONDITIONER("ssor", defaultPreconditionerBlockLevelCreator<Dune::SeqSSOR>());
235 
236 
248  template<class M, class X, class Y, int l=1>
249  class SeqSOR : public Preconditioner<X,Y> {
250  public:
252  typedef M matrix_type;
254  typedef X domain_type;
256  typedef Y range_type;
258  typedef typename X::field_type field_type;
261 
269  SeqSOR (const M& A, int n, scalar_field_type w)
270  : _A_(A), _n(n), _w(w)
271  {
273  }
274 
288  SeqSOR (const M& A, const ParameterTree& configuration)
289  : _A_(A)
290  {
291  _n = configuration.get<int>("iterations",1);
292  _w = configuration.get<scalar_field_type>("relaxation",1.0);
294  }
295 
301  virtual void pre (X& x, Y& b)
302  {
305  }
306 
312  virtual void apply (X& v, const Y& d)
313  {
314  this->template apply<true>(v,d);
315  }
316 
325  template<bool forward>
326  void apply(X& v, const Y& d)
327  {
328  if(forward)
329  for (int i=0; i<_n; i++) {
330  bsorf(_A_,v,d,_w,BL<l>());
331  }
332  else
333  for (int i=0; i<_n; i++) {
334  bsorb(_A_,v,d,_w,BL<l>());
335  }
336  }
337 
343  virtual void post (X& x)
344  {
346  }
347 
350  {
352  }
353 
354  private:
356  const M& _A_;
358  int _n;
361  };
362  DUNE_REGISTER_PRECONDITIONER("sor", defaultPreconditionerBlockLevelCreator<Dune::SeqSOR>());
363 
364 
375  template<class M, class X, class Y, int l=1>
377  DUNE_REGISTER_PRECONDITIONER("gs", defaultPreconditionerBlockLevelCreator<Dune::SeqGS>());
378 
389  template<class M, class X, class Y, int l=1>
390  class SeqJac : public Preconditioner<X,Y> {
391  public:
393  typedef M matrix_type;
395  typedef X domain_type;
397  typedef Y range_type;
399  typedef typename X::field_type field_type;
402 
410  SeqJac (const M& A, int n, scalar_field_type w)
411  : _A_(A), _n(n), _w(w)
412  {
414  }
415 
429  SeqJac (const M& A, const ParameterTree& configuration)
430  : _A_(A)
431  {
432  _n = configuration.get<int>("iterations",1);
433  _w = configuration.get<scalar_field_type>("relaxation",1.0);
435  }
436 
442  virtual void pre (X& x, Y& b)
443  {
446  }
447 
453  virtual void apply (X& v, const Y& d)
454  {
455  for (int i=0; i<_n; i++) {
456  dbjac(_A_,v,d,_w,BL<l>());
457  }
458  }
459 
465  virtual void post (X& x)
466  {
468  }
469 
472  {
474  }
475 
476  private:
478  const M& _A_;
480  int _n;
482  scalar_field_type _w;
483  };
484  DUNE_REGISTER_PRECONDITIONER("jac", defaultPreconditionerBlockLevelCreator<Dune::SeqJac>());
485 
486 
487 
499  template<class M, class X, class Y, int l=1>
500  class SeqILU : public Preconditioner<X,Y> {
501  public:
503  typedef typename std::remove_const<M>::type matrix_type;
505  typedef typename matrix_type :: block_type block_type;
507  typedef X domain_type;
509  typedef Y range_type;
510 
512  typedef typename X::field_type field_type;
513 
516 
518  typedef typename ILU::CRS< block_type , typename M::allocator_type> CRS;
519 
527  SeqILU (const M& A, scalar_field_type w, const bool resort = false )
528  : SeqILU( A, 0, w, resort ) // construct ILU(0)
529  {
530  }
531 
546  SeqILU(const M& A, const ParameterTree& config)
547  : SeqILU(A, config.get("n", 0),
548  config.get<scalar_field_type>("relaxation", 1.0),
549  config.get("resort", false))
550  {}
551 
560  SeqILU (const M& A, int n, scalar_field_type w, const bool resort = false )
561  : ILU_(),
562  lower_(),
563  upper_(),
564  inv_(),
565  w_(w),
566  wNotIdentity_([w]{using std::abs; return abs(w - scalar_field_type(1)) > 1e-15;}() )
567  {
568  if( n == 0 )
569  {
570  // copy A
571  ILU_.reset( new matrix_type( A ) );
572  // create ILU(0) decomposition
574  }
575  else
576  {
577  // create matrix in build mode
578  ILU_.reset( new matrix_type( A.N(), A.M(), matrix_type::row_wise) );
579  // create ILU(n) decomposition
580  bilu_decomposition( A, n, *ILU_ );
581  }
582 
583  if( resort )
584  {
585  // store ILU in simple CRS format
586  ILU::convertToCRS( *ILU_, lower_, upper_, inv_ );
587  ILU_.reset();
588  }
589  }
590 
596  virtual void pre (X& x, Y& b)
597  {
600  }
601 
607  virtual void apply (X& v, const Y& d)
608  {
609  if( ILU_ )
610  {
611  bilu_backsolve( *ILU_, v, d);
612  }
613  else
614  {
615  ILU::bilu_backsolve(lower_, upper_, inv_, v, d);
616  }
617 
618  if( wNotIdentity_ )
619  {
620  v *= w_;
621  }
622  }
623 
629  virtual void post (X& x)
630  {
632  }
633 
636  {
638  }
639 
640  protected:
642  std::unique_ptr< matrix_type > ILU_;
643 
646  CRS upper_;
647  std::vector< block_type, typename matrix_type::allocator_type > inv_;
648 
652  const bool wNotIdentity_;
653  };
654  DUNE_REGISTER_PRECONDITIONER("ilu", defaultPreconditionerBlockLevelCreator<Dune::SeqILU>());
655 
669  template<class M, class X, class Y, int l=1>
670  class DUNE_DEPRECATED_MSG("Use SeqILU instead of SeqILU0!") SeqILU0: public Preconditioner<X,Y> {
671  public:
673  typedef typename std::remove_const<M>::type matrix_type;
675  typedef X domain_type;
677  typedef Y range_type;
679  typedef typename X::field_type field_type;
682 
689  SeqILU0 (const M& A, scalar_field_type w)
690  : _w(w),
691  ILU(A) // copy A
692 
693  {
694  bilu0_decomposition(ILU);
695  }
696 
709  SeqILU0 (const M& A, const ParameterTree& configuration)
710  : ILU(A) // copy A
711  {
712  _w = configuration.get<scalar_field_type>("relaxation");
713  bilu0_decomposition(ILU);
714  }
715 
721  virtual void pre (X& x, Y& b)
722  {
725  }
726 
732  virtual void apply (X& v, const Y& d)
733  {
734  bilu_backsolve(ILU,v,d);
735  v *= _w;
736  }
737 
743  virtual void post (X& x)
744  {
746  }
747 
750  {
752  }
753 
754  private:
756  scalar_field_type _w;
758  matrix_type ILU;
759  };
760 
761 
777  template<class M, class X, class Y, int l=1>
778  class DUNE_DEPRECATED_MSG("Use SeqILU instead of SeqILUn!") SeqILUn : public Preconditioner<X,Y> {
779  public:
781  typedef typename std::remove_const<M>::type matrix_type;
783  typedef X domain_type;
785  typedef Y range_type;
787  typedef typename X::field_type field_type;
790 
798  SeqILUn (const M& A, int n, scalar_field_type w)
799  : ILU(A.N(),A.M(),M::row_wise),
800  _n(n),
801  _w(w)
802  {
803  bilu_decomposition(A,n,ILU);
804  }
805 
809  SeqILUn (const M& A, const ParameterTree& configuration)
810  : ILU(A.N(),A.M(),M::row_wise)
811  {
812  _n = configuration.get<int>("iterations");
813  _w = configuration.get<scalar_field_type>("relaxation");
814  bilu_decomposition(A,_n,ILU);
815  }
816 
822  virtual void pre (X& x, Y& b)
823  {
826  }
827 
833  virtual void apply (X& v, const Y& d)
834  {
835  bilu_backsolve(ILU,v,d);
836  v *= _w;
837  }
838 
844  virtual void post (X& x)
845  {
847  }
848 
851  {
853  }
854 
855  private:
857  matrix_type ILU;
859  int _n;
861  scalar_field_type _w;
862  };
863 
864 
865 
874  template<class X, class Y>
875  class Richardson : public Preconditioner<X,Y> {
876  public:
878  typedef X domain_type;
880  typedef Y range_type;
882  typedef typename X::field_type field_type;
885 
892  _w(w)
893  {}
894 
906  Richardson (const ParameterTree& configuration)
907  {
908  _w = configuration.get<scalar_field_type>("relaxation");
909  }
910 
916  virtual void pre (X& x, Y& b)
917  {
920  }
921 
927  virtual void apply (X& v, const Y& d)
928  {
929  v = d;
930  v *= _w;
931  }
932 
938  virtual void post (X& x)
939  {
941  }
942 
945  {
947  }
948 
949  private:
952  };
953  DUNE_REGISTER_PRECONDITIONER("richardson", [](auto tl, const auto& mat, const ParameterTree& config){
954  using D = typename Dune::TypeListElement<1, decltype(tl)>::type;
955  using R = typename Dune::TypeListElement<2, decltype(tl)>::type;
956  return std::make_shared<Richardson<D,R>>(config);
957  });
958 
959 
970  template< class M, class X, class Y >
971  class SeqILDL
972  : public Preconditioner< X, Y >
973  {
974  typedef SeqILDL< M, X, Y > This;
976 
977  public:
979  typedef std::remove_const_t< M > matrix_type;
981  typedef X domain_type;
983  typedef Y range_type;
985  typedef typename X::field_type field_type;
988 
1001  SeqILDL(const matrix_type& A, const ParameterTree& config)
1002  : SeqILDL(A, config.template get<scalar_field_type>("relaxation", 1.0))
1003  {}
1004 
1013  explicit SeqILDL ( const matrix_type &A, scalar_field_type relax = scalar_field_type( 1 ) )
1014  : decomposition_( A.N(), A.M(), matrix_type::random ),
1015  relax_( relax )
1016  {
1017  // setup row sizes for lower triangular matrix
1018  for( auto i = A.begin(), iend = A.end(); i != iend; ++i )
1019  {
1020  const auto &A_i = *i;
1021  const auto ij = A_i.find( i.index() );
1022  if( ij != A_i.end() )
1023  decomposition_.setrowsize( i.index(), ij.offset()+1 );
1024  else
1025  DUNE_THROW( ISTLError, "diagonal entry missing" );
1026  }
1027  decomposition_.endrowsizes();
1028 
1029  // setup row indices for lower triangular matrix
1030  for( auto i = A.begin(), iend = A.end(); i != iend; ++i )
1031  {
1032  const auto &A_i = *i;
1033  for( auto ij = A_i.begin(); ij.index() < i.index() ; ++ij )
1034  decomposition_.addindex( i.index(), ij.index() );
1035  decomposition_.addindex( i.index(), i.index() );
1036  }
1037  decomposition_.endindices();
1038 
1039  // copy values of lower triangular matrix
1040  auto i = A.begin();
1041  for( auto row = decomposition_.begin(), rowend = decomposition_.end(); row != rowend; ++row, ++i )
1042  {
1043  auto ij = i->begin();
1044  for( auto col = row->begin(), colend = row->end(); col != colend; ++col, ++ij )
1045  *col = *ij;
1046  }
1047 
1048  // perform ILDL decomposition
1049  bildl_decompose( decomposition_ );
1050  }
1051 
1053  void pre ( X &x, Y &b ) override
1054  {
1055  DUNE_UNUSED_PARAMETER( x );
1056  DUNE_UNUSED_PARAMETER( b );
1057  }
1058 
1060  void apply ( X &v, const Y &d ) override
1061  {
1062  bildl_backsolve( decomposition_, v, d, true );
1063  v *= relax_;
1064  }
1065 
1067  void post ( X &x ) override
1068  {
1069  DUNE_UNUSED_PARAMETER( x );
1070  }
1071 
1074 
1075  private:
1076  matrix_type decomposition_;
1077  scalar_field_type relax_;
1078  };
1079  DUNE_REGISTER_PRECONDITIONER("ildl", defaultPreconditionerCreator<Dune::SeqILDL>());
1080 
1083 } // end namespace
1084 
1085 
1086 #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:74
O::range_type range_type
The range type of the preconditioner.
Definition: preconditioners.hh:79
O::domain_type domain_type
The domain type of the preconditioner.
Definition: preconditioners.hh:77
virtual void post(domain_type &)
Clean up.
Definition: preconditioners.hh:108
range_type::field_type field_type
The field type of the preconditioner.
Definition: preconditioners.hh:81
Simd::Scalar< field_type > scalar_field_type
scalar type underlying the field_type
Definition: preconditioners.hh:83
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: preconditioners.hh:112
virtual void pre(domain_type &, range_type &)
Prepare the preconditioner.
Definition: preconditioners.hh:98
InverseOperator2Preconditioner(InverseOperator &inverse_operator)
Construct the preconditioner from the solver.
Definition: preconditioners.hh:91
O InverseOperator
type of the wrapped inverse operator
Definition: preconditioners.hh:85
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:101
Abstract base class for all solvers.
Definition: solver.hh:97
Hierarchical structure of string parameters.
Definition: parametertree.hh:35
std::string get(const std::string &key, const std::string &defaultValue) const
get value as string
Definition: parametertree.cc:183
Base class for matrix free definition of preconditioners.
Definition: preconditioner.hh:30
Richardson preconditioner.
Definition: preconditioners.hh:875
X::field_type field_type
The field type of the preconditioner.
Definition: preconditioners.hh:882
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: preconditioners.hh:944
Y range_type
The range type of the preconditioner.
Definition: preconditioners.hh:880
virtual void pre(X &x, Y &b)
Prepare the preconditioner.
Definition: preconditioners.hh:916
virtual void post(X &x)
Clean up.
Definition: preconditioners.hh:938
Richardson(scalar_field_type w=1.0)
Constructor.
Definition: preconditioners.hh:891
Simd::Scalar< field_type > scalar_field_type
scalar type underlying the field_type
Definition: preconditioners.hh:884
Richardson(const ParameterTree &configuration)
Constructor.
Definition: preconditioners.hh:906
X domain_type
The domain type of the preconditioner.
Definition: preconditioners.hh:878
virtual void apply(X &v, const Y &d)
Apply the precondioner.
Definition: preconditioners.hh:927
sequential ILDL preconditioner
Definition: preconditioners.hh:973
SeqILDL(const matrix_type &A, const ParameterTree &config)
Constructor.
Definition: preconditioners.hh:1001
SeqILDL(const matrix_type &A, scalar_field_type relax=scalar_field_type(1))
constructor
Definition: preconditioners.hh:1013
X domain_type
domain type of the preconditioner
Definition: preconditioners.hh:981
void post(X &x) override
Clean up.
Definition: preconditioners.hh:1067
Y range_type
range type of the preconditioner
Definition: preconditioners.hh:983
std::remove_const_t< M > matrix_type
type of matrix the preconditioner is for
Definition: preconditioners.hh:979
void apply(X &v, const Y &d) override
Apply one step of the preconditioner to the system A(v)=d.
Definition: preconditioners.hh:1060
void pre(X &x, Y &b) override
Prepare the preconditioner.
Definition: preconditioners.hh:1053
Simd::Scalar< field_type > scalar_field_type
scalar type underlying the field_type
Definition: preconditioners.hh:987
X::field_type field_type
field type of the preconditioner
Definition: preconditioners.hh:985
SolverCategory::Category category() const override
Category of the preconditioner (see SolverCategory::Category)
Definition: preconditioners.hh:1073
Sequential ILU0 preconditioner.
Definition: preconditioners.hh:670
virtual void post(X &x)
Clean up.
Definition: preconditioners.hh:743
Simd::Scalar< field_type > scalar_field_type
scalar type underlying the field_type
Definition: preconditioners.hh:681
virtual void pre(X &x, Y &b)
Prepare the preconditioner.
Definition: preconditioners.hh:721
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: preconditioners.hh:749
virtual void apply(X &v, const Y &d)
Apply the preconditioner.
Definition: preconditioners.hh:732
Y range_type
The range type of the preconditioner.
Definition: preconditioners.hh:677
SeqILU0(const M &A, scalar_field_type w)
Constructor.
Definition: preconditioners.hh:689
SeqILU0(const M &A, const ParameterTree &configuration)
Constructor.
Definition: preconditioners.hh:709
X domain_type
The domain type of the preconditioner.
Definition: preconditioners.hh:675
X::field_type field_type
The field type of the preconditioner.
Definition: preconditioners.hh:679
std::remove_const< M >::type matrix_type
The matrix type the preconditioner is for.
Definition: preconditioners.hh:673
Sequential ILU preconditioner.
Definition: preconditioners.hh:500
virtual void post(X &x)
Clean up.
Definition: preconditioners.hh:629
virtual void pre(X &x, Y &b)
Prepare the preconditioner.
Definition: preconditioners.hh:596
matrix_type ::block_type block_type
block type of matrix
Definition: preconditioners.hh:505
virtual void apply(X &v, const Y &d)
Apply the preconditioner.
Definition: preconditioners.hh:607
ILU::CRS< block_type, typename M::allocator_type > CRS
type of ILU storage
Definition: preconditioners.hh:518
Y range_type
The range type of the preconditioner.
Definition: preconditioners.hh:509
const scalar_field_type w_
The relaxation factor to use.
Definition: preconditioners.hh:650
CRS lower_
The ILU(n) decomposition of the matrix. As storage a CRS structure is used.
Definition: preconditioners.hh:645
const bool wNotIdentity_
true if w != 1.0
Definition: preconditioners.hh:652
SeqILU(const M &A, const ParameterTree &config)
Constructor.
Definition: preconditioners.hh:546
SeqILU(const M &A, int n, scalar_field_type w, const bool resort=false)
Constructor.
Definition: preconditioners.hh:560
std::remove_const< M >::type matrix_type
The matrix type the preconditioner is for.
Definition: preconditioners.hh:503
SeqILU(const M &A, scalar_field_type w, const bool resort=false)
Constructor.
Definition: preconditioners.hh:527
X::field_type field_type
The field type of the preconditioner.
Definition: preconditioners.hh:512
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: preconditioners.hh:635
X domain_type
The domain type of the preconditioner.
Definition: preconditioners.hh:507
Simd::Scalar< field_type > scalar_field_type
scalar type underlying the field_type
Definition: preconditioners.hh:515
std::unique_ptr< matrix_type > ILU_
The ILU(n) decomposition of the matrix. As storage a BCRSMatrix is used.
Definition: preconditioners.hh:642
Sequential ILU(n) preconditioner.
Definition: preconditioners.hh:778
Simd::Scalar< field_type > scalar_field_type
scalar type underlying the field_type
Definition: preconditioners.hh:789
SeqILUn(const M &A, int n, scalar_field_type w)
Constructor.
Definition: preconditioners.hh:798
virtual void pre(X &x, Y &b)
Prepare the preconditioner.
Definition: preconditioners.hh:822
std::remove_const< M >::type matrix_type
The matrix type the preconditioner is for.
Definition: preconditioners.hh:781
virtual void post(X &x)
Clean up.
Definition: preconditioners.hh:844
X domain_type
The domain type of the preconditioner.
Definition: preconditioners.hh:783
SeqILUn(const M &A, const ParameterTree &configuration)
Constructor.
Definition: preconditioners.hh:809
Y range_type
The range type of the preconditioner.
Definition: preconditioners.hh:785
virtual void apply(X &v, const Y &d)
Apply the precondioner.
Definition: preconditioners.hh:833
X::field_type field_type
The field type of the preconditioner.
Definition: preconditioners.hh:787
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: preconditioners.hh:850
The sequential jacobian preconditioner.
Definition: preconditioners.hh:390
virtual void post(X &x)
Clean up.
Definition: preconditioners.hh:465
SeqJac(const M &A, const ParameterTree &configuration)
Constructor.
Definition: preconditioners.hh:429
virtual void apply(X &v, const Y &d)
Apply the preconditioner.
Definition: preconditioners.hh:453
M matrix_type
The matrix type the preconditioner is for.
Definition: preconditioners.hh:393
Simd::Scalar< field_type > scalar_field_type
scalar type underlying the field_type
Definition: preconditioners.hh:401
X::field_type field_type
The field type of the preconditioner.
Definition: preconditioners.hh:399
SeqJac(const M &A, int n, scalar_field_type w)
Constructor.
Definition: preconditioners.hh:410
virtual void pre(X &x, Y &b)
Prepare the preconditioner.
Definition: preconditioners.hh:442
X domain_type
The domain type of the preconditioner.
Definition: preconditioners.hh:395
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: preconditioners.hh:471
Y range_type
The range type of the preconditioner.
Definition: preconditioners.hh:397
Sequential SOR preconditioner.
Definition: preconditioners.hh:249
M matrix_type
The matrix type the preconditioner is for.
Definition: preconditioners.hh:252
void apply(X &v, const Y &d)
Apply the preconditioner in a special direction.
Definition: preconditioners.hh:326
X domain_type
The domain type of the preconditioner.
Definition: preconditioners.hh:254
virtual void post(X &x)
Clean up.
Definition: preconditioners.hh:343
SeqSOR(const M &A, int n, scalar_field_type w)
Constructor.
Definition: preconditioners.hh:269
virtual void pre(X &x, Y &b)
Prepare the preconditioner.
Definition: preconditioners.hh:301
Simd::Scalar< field_type > scalar_field_type
scalar type underlying the field_type
Definition: preconditioners.hh:260
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: preconditioners.hh:349
virtual void apply(X &v, const Y &d)
Apply the preconditioner.
Definition: preconditioners.hh:312
Y range_type
The range type of the preconditioner.
Definition: preconditioners.hh:256
SeqSOR(const M &A, const ParameterTree &configuration)
Constructor.
Definition: preconditioners.hh:288
X::field_type field_type
The field type of the preconditioner.
Definition: preconditioners.hh:258
Sequential SSOR preconditioner.
Definition: preconditioners.hh:138
virtual void post(X &x)
Clean up.
Definition: preconditioners.hh:215
SeqSSOR(const M &A, const ParameterTree &configuration)
Constructor.
Definition: preconditioners.hh:177
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: preconditioners.hh:221
X::field_type field_type
The field type of the preconditioner.
Definition: preconditioners.hh:147
Simd::Scalar< field_type > scalar_field_type
scalar type underlying the field_type
Definition: preconditioners.hh:149
SeqSSOR(const M &A, int n, scalar_field_type w)
Constructor.
Definition: preconditioners.hh:158
X domain_type
The domain type of the preconditioner.
Definition: preconditioners.hh:143
M matrix_type
The matrix type the preconditioner is for.
Definition: preconditioners.hh:141
virtual void apply(X &v, const Y &d)
Apply the preconditioner.
Definition: preconditioners.hh:202
virtual void pre(X &x, Y &b)
Prepare the preconditioner.
Definition: preconditioners.hh:190
Y range_type
The range type of the preconditioner.
Definition: preconditioners.hh:145
#define DUNE_DEPRECATED_MSG(text)
Mark some entity as deprecated.
Definition: deprecated.hh:169
#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:640
void dbjac(const M &A, X &x, const Y &b, const K &w)
Jacobi step.
Definition: gsetc.hh:652
void bilu_backsolve(const M &A, X &v, const Y &d)
LU backsolve with stored inverse.
Definition: ilu.hh:95
void bilu_decomposition(const M &A, int n, M &ILU)
Definition: ilu.hh:169
void bilu0_decomposition(M &A)
compute ILU decomposition of A. A is overwritten by its decomposition
Definition: ilu.hh:34
void bsorf(const M &A, X &x, const Y &b, const K &w)
SOR step.
Definition: gsetc.hh:628
typename Overloads::ScalarType< std::decay_t< V > >::type Scalar
Element type of some SIMD type.
Definition: interface.hh:233
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:14
void bildl_decompose(Matrix &A)
compute ILDL decomposition of a symmetric matrix A
Definition: ildl.hh:86
A hierarchical structure of string parameters.
Include file for users of the SIMD abstraction layer.
Define general, extensible interface for inverse operators.
compile-time parameter for block recursion depth
Definition: gsetc.hh:43
static void check(const Matrix &mat)
Check whether the a matrix has diagonal values on blocklevel recursion levels.
Definition: matrixutils.hh:52
Statistics about the application of an inverse operator.
Definition: solver.hh:46
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 8, 22:30, 2024)