Dune Core Modules (2.9.0)

preconditioners.hh
Go to the documentation of this file.
1 // SPDX-FileCopyrightText: Copyright (C) DUNE Project contributors, see file LICENSE.md in module root
2 // SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
3 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
4 // vi: set et ts=4 sw=2 sts=2:
5 #ifndef DUNE_ISTL_PRECONDITIONERS_HH
6 #define DUNE_ISTL_PRECONDITIONERS_HH
7 
8 #include <cmath>
9 #include <complex>
10 #include <iostream>
11 #include <iomanip>
12 #include <memory>
13 #include <string>
14 
15 #include <dune/common/simd/simd.hh>
17 
18 #include <dune/istl/solverregistry.hh>
19 #include "preconditioner.hh"
20 #include "solver.hh"
21 #include "solvercategory.hh"
22 #include "istlexception.hh"
23 #include "matrixutils.hh"
24 #include "gsetc.hh"
25 #include "ildl.hh"
26 #include "ilu.hh"
27 
28 
29 namespace Dune {
72  template<class O, int c = -1>
74  public Preconditioner<typename O::domain_type, typename O::range_type>
75  {
76  public:
78  typedef typename O::domain_type domain_type;
80  typedef typename O::range_type range_type;
82  typedef typename range_type::field_type field_type;
86  typedef typename FieldTraits<scalar_field_type>::real_type real_field_type;
88  typedef O InverseOperator;
89 
95  : inverse_operator_(inverse_operator)
96  {
97  if(c != -1 && SolverCategory::category(inverse_operator_) != c)
98  DUNE_THROW(InvalidStateException, "User-supplied solver category does not match that of the given inverse operator");
99  }
100 
101  virtual void pre(domain_type&,range_type&)
102  {}
103 
104  virtual void apply(domain_type& v, const range_type& d)
105  {
107  range_type copy(d);
108  inverse_operator_.apply(v, copy, res);
109  }
110 
111  virtual void post(domain_type&)
112  {}
113 
116  {
117  return SolverCategory::category(inverse_operator_);
118  }
119 
120  private:
121  InverseOperator& inverse_operator_;
122  };
123 
124  //=====================================================================
125  // Implementation of this interface for sequential ISTL-preconditioners
126  //=====================================================================
127 
128 
140  template<class M, class X, class Y, int l=1>
141  class SeqSSOR : public Preconditioner<X,Y> {
142  public:
144  typedef M matrix_type;
146  typedef X domain_type;
148  typedef Y range_type;
150  typedef typename X::field_type field_type;
154  typedef typename FieldTraits<scalar_field_type>::real_type real_field_type;
155 
163  SeqSSOR (const M& A, int n, real_field_type w)
164  : _A_(A), _n(n), _w(w)
165  {
167  }
168 
182  SeqSSOR (const std::shared_ptr<const AssembledLinearOperator<M,X,Y>>& A, const ParameterTree& configuration)
183  : SeqSSOR(A->getmat(), configuration)
184  {}
185 
199  SeqSSOR (const M& A, const ParameterTree& configuration)
200  : SeqSSOR(A, configuration.get<int>("iterations",1), configuration.get<real_field_type>("relaxation",1.0))
201  {}
202 
208  virtual void pre ([[maybe_unused]] X& x, [[maybe_unused]] Y& b)
209  {}
210 
216  virtual void apply (X& v, const Y& d)
217  {
218  for (int i=0; i<_n; i++) {
219  bsorf(_A_,v,d,_w,BL<l>());
220  bsorb(_A_,v,d,_w,BL<l>());
221  }
222  }
223 
229  virtual void post ([[maybe_unused]] X& x)
230  {}
231 
234  {
236  }
237 
238  private:
240  const M& _A_;
242  int _n;
244  real_field_type _w;
245  };
246  DUNE_REGISTER_PRECONDITIONER("ssor", defaultPreconditionerBlockLevelCreator<Dune::SeqSSOR>());
247 
248 
260  template<class M, class X, class Y, int l=1>
261  class SeqSOR : public Preconditioner<X,Y> {
262  public:
264  typedef M matrix_type;
266  typedef X domain_type;
268  typedef Y range_type;
270  typedef typename X::field_type field_type;
274  typedef typename FieldTraits<scalar_field_type>::real_type real_field_type;
275 
283  SeqSOR (const M& A, int n, real_field_type w)
284  : _A_(A), _n(n), _w(w)
285  {
287  }
288 
302  SeqSOR (const std::shared_ptr<const AssembledLinearOperator<M,X,Y>>& A, const ParameterTree& configuration)
303  : SeqSOR(A->getmat(), configuration)
304  {}
305 
319  SeqSOR (const M& A, const ParameterTree& configuration)
320  : SeqSOR(A, configuration.get<int>("iterations",1), configuration.get<real_field_type>("relaxation",1.0))
321  {}
322 
328  virtual void pre ([[maybe_unused]] X& x, [[maybe_unused]] Y& b)
329  {}
330 
336  virtual void apply (X& v, const Y& d)
337  {
338  this->template apply<true>(v,d);
339  }
340 
349  template<bool forward>
350  void apply(X& v, const Y& d)
351  {
352  if(forward)
353  for (int i=0; i<_n; i++) {
354  bsorf(_A_,v,d,_w,BL<l>());
355  }
356  else
357  for (int i=0; i<_n; i++) {
358  bsorb(_A_,v,d,_w,BL<l>());
359  }
360  }
361 
367  virtual void post ([[maybe_unused]] X& x)
368  {}
369 
372  {
374  }
375 
376  private:
378  const M& _A_;
380  int _n;
382  real_field_type _w;
383  };
384  DUNE_REGISTER_PRECONDITIONER("sor", defaultPreconditionerBlockLevelCreator<Dune::SeqSOR>());
385 
386 
397  template<class M, class X, class Y, int l=1>
399  DUNE_REGISTER_PRECONDITIONER("gs", defaultPreconditionerBlockLevelCreator<Dune::SeqGS>());
400 
411  template<class M, class X, class Y, int l=1>
412  class SeqJac : public Preconditioner<X,Y> {
413  public:
415  typedef M matrix_type;
417  typedef X domain_type;
419  typedef Y range_type;
421  typedef typename X::field_type field_type;
425  typedef typename FieldTraits<scalar_field_type>::real_type real_field_type;
426 
434  SeqJac (const M& A, int n, real_field_type w)
435  : _A_(A), _n(n), _w(w)
436  {
438  }
439 
453  SeqJac (const std::shared_ptr<const AssembledLinearOperator<M,X,Y>>& A, const ParameterTree& configuration)
454  : SeqJac(A->getmat(), configuration)
455  {}
456 
470  SeqJac (const M& A, const ParameterTree& configuration)
471  : SeqJac(A, configuration.get<int>("iterations",1), configuration.get<real_field_type>("relaxation",1.0))
472  {}
473 
479  virtual void pre ([[maybe_unused]] X& x, [[maybe_unused]] Y& b)
480  {}
481 
487  virtual void apply (X& v, const Y& d)
488  {
489  for (int i=0; i<_n; i++) {
490  dbjac(_A_,v,d,_w,BL<l>());
491  }
492  }
493 
499  virtual void post ([[maybe_unused]] X& x)
500  {}
501 
504  {
506  }
507 
508  private:
510  const M& _A_;
512  int _n;
514  real_field_type _w;
515  };
516  DUNE_REGISTER_PRECONDITIONER("jac", defaultPreconditionerBlockLevelCreator<Dune::SeqJac>());
517 
518 
519 
531  template<class M, class X, class Y, int l=1>
532  class SeqILU : public Preconditioner<X,Y> {
533  public:
535  typedef typename std::remove_const<M>::type matrix_type;
537  typedef typename matrix_type :: block_type block_type;
539  typedef X domain_type;
541  typedef Y range_type;
542 
544  typedef typename X::field_type field_type;
545 
549  typedef typename FieldTraits<scalar_field_type>::real_type real_field_type;
550 
553 
561  SeqILU (const M& A, real_field_type w, const bool resort = false )
562  : SeqILU( A, 0, w, resort ) // construct ILU(0)
563  {
564  }
565 
580  SeqILU (const std::shared_ptr<const AssembledLinearOperator<M,X,Y>>& A, const ParameterTree& configuration)
581  : SeqILU(A->getmat(), configuration)
582  {}
583 
598  SeqILU(const M& A, const ParameterTree& config)
599  : SeqILU(A, config.get("n", 0),
600  config.get<real_field_type>("relaxation", 1.0),
601  config.get("resort", false))
602  {}
603 
612  SeqILU (const M& A, int n, real_field_type w, const bool resort = false )
613  : ILU_(),
614  lower_(),
615  upper_(),
616  inv_(),
617  w_(w),
618  wNotIdentity_([w]{using std::abs; return abs(w - real_field_type(1)) > 1e-15;}() )
619  {
620  if( n == 0 )
621  {
622  // copy A
623  ILU_.reset( new matrix_type( A ) );
624  // create ILU(0) decomposition
625  ILU::blockILU0Decomposition( *ILU_ );
626  }
627  else
628  {
629  // create matrix in build mode
630  ILU_.reset( new matrix_type( A.N(), A.M(), matrix_type::row_wise) );
631  // create ILU(n) decomposition
632  ILU::blockILUDecomposition( A, n, *ILU_ );
633  }
634 
635  if( resort )
636  {
637  // store ILU in simple CRS format
638  ILU::convertToCRS( *ILU_, lower_, upper_, inv_ );
639  ILU_.reset();
640  }
641  }
642 
648  virtual void pre ([[maybe_unused]] X& x, [[maybe_unused]] Y& b)
649  {}
650 
656  virtual void apply (X& v, const Y& d)
657  {
658  if( ILU_ )
659  {
660  ILU::blockILUBacksolve( *ILU_, v, d);
661  }
662  else
663  {
664  ILU::blockILUBacksolve(lower_, upper_, inv_, v, d);
665  }
666 
667  if( wNotIdentity_ )
668  {
669  v *= w_;
670  }
671  }
672 
678  virtual void post ([[maybe_unused]] X& x)
679  {}
680 
683  {
685  }
686 
687  protected:
689  std::unique_ptr< matrix_type > ILU_;
690 
693  CRS upper_;
694  std::vector< block_type, typename matrix_type::allocator_type > inv_;
695 
699  const bool wNotIdentity_;
700  };
701  DUNE_REGISTER_PRECONDITIONER("ilu", defaultPreconditionerBlockLevelCreator<Dune::SeqILU>());
702 
703 
712  template<class X, class Y>
713  class Richardson : public Preconditioner<X,Y> {
714  public:
716  typedef X domain_type;
718  typedef Y range_type;
720  typedef typename X::field_type field_type;
724  typedef typename FieldTraits<scalar_field_type>::real_type real_field_type;
725 
732  _w(w)
733  {}
734 
746  Richardson (const ParameterTree& configuration)
747  : Richardson(configuration.get<real_field_type>("relaxation", 1.0))
748  {}
749 
755  virtual void pre ([[maybe_unused]] X& x, [[maybe_unused]] Y& b)
756  {}
757 
763  virtual void apply (X& v, const Y& d)
764  {
765  v = d;
766  v *= _w;
767  }
768 
774  virtual void post ([[maybe_unused]] X& x)
775  {}
776 
779  {
781  }
782 
783  private:
785  real_field_type _w;
786  };
787  DUNE_REGISTER_PRECONDITIONER("richardson", [](auto tl, const auto& /* mat */, const ParameterTree& config){
788  using D = typename Dune::TypeListElement<1, decltype(tl)>::type;
789  using R = typename Dune::TypeListElement<2, decltype(tl)>::type;
790  return std::make_shared<Richardson<D,R>>(config);
791  });
792 
793 
804  template< class M, class X, class Y >
805  class SeqILDL
806  : public Preconditioner< X, Y >
807  {
808  typedef SeqILDL< M, X, Y > This;
810 
811  public:
813  typedef std::remove_const_t< M > matrix_type;
815  typedef X domain_type;
817  typedef Y range_type;
819  typedef typename X::field_type field_type;
823  typedef typename FieldTraits<scalar_field_type>::real_type real_field_type;
824 
837  SeqILDL (const std::shared_ptr<const AssembledLinearOperator<M,X,Y>>& A, const ParameterTree& configuration)
838  : SeqILDL(A->getmat(), configuration)
839  {}
840 
853  SeqILDL(const matrix_type& A, const ParameterTree& config)
854  : SeqILDL(A, config.get<real_field_type>("relaxation", 1.0))
855  {}
856 
865  explicit SeqILDL ( const matrix_type &A, real_field_type relax = real_field_type( 1 ) )
866  : decomposition_( A.N(), A.M(), matrix_type::random ),
867  relax_( relax )
868  {
869  // setup row sizes for lower triangular matrix
870  for( auto i = A.begin(), iend = A.end(); i != iend; ++i )
871  {
872  const auto &A_i = *i;
873  const auto ij = A_i.find( i.index() );
874  if( ij != A_i.end() )
875  decomposition_.setrowsize( i.index(), ij.offset()+1 );
876  else
877  DUNE_THROW( ISTLError, "diagonal entry missing" );
878  }
879  decomposition_.endrowsizes();
880 
881  // setup row indices for lower triangular matrix
882  for( auto i = A.begin(), iend = A.end(); i != iend; ++i )
883  {
884  const auto &A_i = *i;
885  for( auto ij = A_i.begin(); ij.index() < i.index() ; ++ij )
886  decomposition_.addindex( i.index(), ij.index() );
887  decomposition_.addindex( i.index(), i.index() );
888  }
889  decomposition_.endindices();
890 
891  // copy values of lower triangular matrix
892  auto i = A.begin();
893  for( auto row = decomposition_.begin(), rowend = decomposition_.end(); row != rowend; ++row, ++i )
894  {
895  auto ij = i->begin();
896  for( auto col = row->begin(), colend = row->end(); col != colend; ++col, ++ij )
897  *col = *ij;
898  }
899 
900  // perform ILDL decomposition
901  bildl_decompose( decomposition_ );
902  }
903 
905  void pre ([[maybe_unused]] X &x, [[maybe_unused]] Y &b) override
906  {}
907 
909  void apply ( X &v, const Y &d ) override
910  {
911  bildl_backsolve( decomposition_, v, d, true );
912  v *= relax_;
913  }
914 
916  void post ([[maybe_unused]] X &x) override
917  {}
918 
921 
922  private:
923  matrix_type decomposition_;
924  real_field_type relax_;
925  };
926  DUNE_REGISTER_PRECONDITIONER("ildl", defaultPreconditionerCreator<Dune::SeqILDL>());
927 
930 } // end namespace
931 
932 
933 #endif
A linear operator exporting itself in matrix form.
Definition: operators.hh:109
derive error class from the base class in common
Definition: istlexception.hh:19
Default exception if a function was called while the object is not in a valid state for that function...
Definition: exceptions.hh:281
Turns an InverseOperator into a Preconditioner.
Definition: preconditioners.hh:75
O::range_type range_type
The range type of the preconditioner.
Definition: preconditioners.hh:80
O::domain_type domain_type
The domain type of the preconditioner.
Definition: preconditioners.hh:78
virtual void post(domain_type &)
Clean up.
Definition: preconditioners.hh:111
range_type::field_type field_type
The field type of the preconditioner.
Definition: preconditioners.hh:82
FieldTraits< scalar_field_type >::real_type real_field_type
real scalar type underlying the field_type
Definition: preconditioners.hh:86
Simd::Scalar< field_type > scalar_field_type
scalar type underlying the field_type
Definition: preconditioners.hh:84
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: preconditioners.hh:115
virtual void pre(domain_type &, range_type &)
Prepare the preconditioner.
Definition: preconditioners.hh:101
InverseOperator2Preconditioner(InverseOperator &inverse_operator)
Construct the preconditioner from the solver.
Definition: preconditioners.hh:94
O InverseOperator
type of the wrapped inverse operator
Definition: preconditioners.hh:88
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:104
Abstract base class for all solvers.
Definition: solver.hh:99
size_type M() const
Return the number of columns.
Definition: matrix.hh:700
size_type N() const
Return the number of rows.
Definition: matrix.hh:695
Hierarchical structure of string parameters.
Definition: parametertree.hh:37
Base class for matrix free definition of preconditioners.
Definition: preconditioner.hh:32
Richardson preconditioner.
Definition: preconditioners.hh:713
virtual void post([[maybe_unused]] X &x)
Clean up.
Definition: preconditioners.hh:774
X::field_type field_type
The field type of the preconditioner.
Definition: preconditioners.hh:720
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: preconditioners.hh:778
Y range_type
The range type of the preconditioner.
Definition: preconditioners.hh:718
Richardson(real_field_type w=1.0)
Constructor.
Definition: preconditioners.hh:731
FieldTraits< scalar_field_type >::real_type real_field_type
real scalar type underlying the field_type
Definition: preconditioners.hh:724
Simd::Scalar< field_type > scalar_field_type
scalar type underlying the field_type
Definition: preconditioners.hh:722
Richardson(const ParameterTree &configuration)
Constructor.
Definition: preconditioners.hh:746
X domain_type
The domain type of the preconditioner.
Definition: preconditioners.hh:716
virtual void pre([[maybe_unused]] X &x, [[maybe_unused]] Y &b)
Prepare the preconditioner.
Definition: preconditioners.hh:755
virtual void apply(X &v, const Y &d)
Apply the precondioner.
Definition: preconditioners.hh:763
sequential ILDL preconditioner
Definition: preconditioners.hh:807
SeqILDL(const matrix_type &A, const ParameterTree &config)
Constructor.
Definition: preconditioners.hh:853
SeqILDL(const matrix_type &A, real_field_type relax=real_field_type(1))
constructor
Definition: preconditioners.hh:865
void pre([[maybe_unused]] X &x, [[maybe_unused]] Y &b) override
Prepare the preconditioner.
Definition: preconditioners.hh:905
void post([[maybe_unused]] X &x) override
Clean up.
Definition: preconditioners.hh:916
X domain_type
domain type of the preconditioner
Definition: preconditioners.hh:815
Y range_type
range type of the preconditioner
Definition: preconditioners.hh:817
std::remove_const_t< M > matrix_type
type of matrix the preconditioner is for
Definition: preconditioners.hh:813
void apply(X &v, const Y &d) override
Apply one step of the preconditioner to the system A(v)=d.
Definition: preconditioners.hh:909
FieldTraits< scalar_field_type >::real_type real_field_type
real scalar type underlying the field_type
Definition: preconditioners.hh:823
SeqILDL(const std::shared_ptr< const AssembledLinearOperator< M, X, Y >> &A, const ParameterTree &configuration)
Constructor.
Definition: preconditioners.hh:837
Simd::Scalar< field_type > scalar_field_type
scalar type underlying the field_type
Definition: preconditioners.hh:821
X::field_type field_type
field type of the preconditioner
Definition: preconditioners.hh:819
SolverCategory::Category category() const override
Category of the preconditioner (see SolverCategory::Category)
Definition: preconditioners.hh:920
Sequential ILU preconditioner.
Definition: preconditioners.hh:532
SeqILU(const M &A, int n, real_field_type w, const bool resort=false)
Constructor.
Definition: preconditioners.hh:612
matrix_type ::block_type block_type
block type of matrix
Definition: preconditioners.hh:537
virtual void apply(X &v, const Y &d)
Apply the preconditioner.
Definition: preconditioners.hh:656
ILU::CRS< block_type, typename M::allocator_type > CRS
type of ILU storage
Definition: preconditioners.hh:552
Y range_type
The range type of the preconditioner.
Definition: preconditioners.hh:541
CRS lower_
The ILU(n) decomposition of the matrix. As storage a CRS structure is used.
Definition: preconditioners.hh:692
const bool wNotIdentity_
true if w != 1.0
Definition: preconditioners.hh:699
SeqILU(const M &A, const ParameterTree &config)
Constructor.
Definition: preconditioners.hh:598
virtual void post([[maybe_unused]] X &x)
Clean up.
Definition: preconditioners.hh:678
std::remove_const< M >::type matrix_type
The matrix type the preconditioner is for.
Definition: preconditioners.hh:535
FieldTraits< scalar_field_type >::real_type real_field_type
real scalar type underlying the field_type
Definition: preconditioners.hh:549
virtual void pre([[maybe_unused]] X &x, [[maybe_unused]] Y &b)
Prepare the preconditioner.
Definition: preconditioners.hh:648
X::field_type field_type
The field type of the preconditioner.
Definition: preconditioners.hh:544
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: preconditioners.hh:682
SeqILU(const M &A, real_field_type w, const bool resort=false)
Constructor.
Definition: preconditioners.hh:561
const real_field_type w_
The relaxation factor to use.
Definition: preconditioners.hh:697
X domain_type
The domain type of the preconditioner.
Definition: preconditioners.hh:539
Simd::Scalar< field_type > scalar_field_type
scalar type underlying the field_type
Definition: preconditioners.hh:547
SeqILU(const std::shared_ptr< const AssembledLinearOperator< M, X, Y >> &A, const ParameterTree &configuration)
Constructor.
Definition: preconditioners.hh:580
std::unique_ptr< matrix_type > ILU_
The ILU(n) decomposition of the matrix. As storage a BCRSMatrix is used.
Definition: preconditioners.hh:689
The sequential jacobian preconditioner.
Definition: preconditioners.hh:412
virtual void post([[maybe_unused]] X &x)
Clean up.
Definition: preconditioners.hh:499
SeqJac(const M &A, const ParameterTree &configuration)
Constructor.
Definition: preconditioners.hh:470
virtual void apply(X &v, const Y &d)
Apply the preconditioner.
Definition: preconditioners.hh:487
M matrix_type
The matrix type the preconditioner is for.
Definition: preconditioners.hh:415
Simd::Scalar< field_type > scalar_field_type
scalar type underlying the field_type
Definition: preconditioners.hh:423
SeqJac(const std::shared_ptr< const AssembledLinearOperator< M, X, Y >> &A, const ParameterTree &configuration)
Constructor.
Definition: preconditioners.hh:453
X::field_type field_type
The field type of the preconditioner.
Definition: preconditioners.hh:421
virtual void pre([[maybe_unused]] X &x, [[maybe_unused]] Y &b)
Prepare the preconditioner.
Definition: preconditioners.hh:479
X domain_type
The domain type of the preconditioner.
Definition: preconditioners.hh:417
FieldTraits< scalar_field_type >::real_type real_field_type
real scalar type underlying the field_type
Definition: preconditioners.hh:425
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: preconditioners.hh:503
SeqJac(const M &A, int n, real_field_type w)
Constructor.
Definition: preconditioners.hh:434
Y range_type
The range type of the preconditioner.
Definition: preconditioners.hh:419
Sequential SOR preconditioner.
Definition: preconditioners.hh:261
M matrix_type
The matrix type the preconditioner is for.
Definition: preconditioners.hh:264
FieldTraits< scalar_field_type >::real_type real_field_type
real scalar type underlying the field_type
Definition: preconditioners.hh:274
void apply(X &v, const Y &d)
Apply the preconditioner in a special direction.
Definition: preconditioners.hh:350
X domain_type
The domain type of the preconditioner.
Definition: preconditioners.hh:266
Simd::Scalar< field_type > scalar_field_type
scalar type underlying the field_type
Definition: preconditioners.hh:272
SeqSOR(const std::shared_ptr< const AssembledLinearOperator< M, X, Y >> &A, const ParameterTree &configuration)
Constructor.
Definition: preconditioners.hh:302
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: preconditioners.hh:371
virtual void pre([[maybe_unused]] X &x, [[maybe_unused]] Y &b)
Prepare the preconditioner.
Definition: preconditioners.hh:328
virtual void apply(X &v, const Y &d)
Apply the preconditioner.
Definition: preconditioners.hh:336
Y range_type
The range type of the preconditioner.
Definition: preconditioners.hh:268
SeqSOR(const M &A, const ParameterTree &configuration)
Constructor.
Definition: preconditioners.hh:319
X::field_type field_type
The field type of the preconditioner.
Definition: preconditioners.hh:270
SeqSOR(const M &A, int n, real_field_type w)
Constructor.
Definition: preconditioners.hh:283
virtual void post([[maybe_unused]] X &x)
Clean up.
Definition: preconditioners.hh:367
Sequential SSOR preconditioner.
Definition: preconditioners.hh:141
SeqSSOR(const std::shared_ptr< const AssembledLinearOperator< M, X, Y >> &A, const ParameterTree &configuration)
Constructor.
Definition: preconditioners.hh:182
SeqSSOR(const M &A, const ParameterTree &configuration)
Constructor.
Definition: preconditioners.hh:199
virtual void pre([[maybe_unused]] X &x, [[maybe_unused]] Y &b)
Prepare the preconditioner.
Definition: preconditioners.hh:208
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: preconditioners.hh:233
X::field_type field_type
The field type of the preconditioner.
Definition: preconditioners.hh:150
Simd::Scalar< field_type > scalar_field_type
scalar type underlying the field_type
Definition: preconditioners.hh:152
X domain_type
The domain type of the preconditioner.
Definition: preconditioners.hh:146
M matrix_type
The matrix type the preconditioner is for.
Definition: preconditioners.hh:144
virtual void post([[maybe_unused]] X &x)
Clean up.
Definition: preconditioners.hh:229
virtual void apply(X &v, const Y &d)
Apply the preconditioner.
Definition: preconditioners.hh:216
Y range_type
The range type of the preconditioner.
Definition: preconditioners.hh:148
FieldTraits< scalar_field_type >::real_type real_field_type
real scalar type underlying the field_type
Definition: preconditioners.hh:154
SeqSSOR(const M &A, int n, real_field_type w)
Constructor.
Definition: preconditioners.hh:163
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218
void bsorb(const M &A, X &x, const Y &b, const K &w)
SSOR step.
Definition: gsetc.hh:646
void dbjac(const M &A, X &x, const Y &b, const K &w)
Jacobi step.
Definition: gsetc.hh:658
void bsorf(const M &A, X &x, const Y &b, const K &w)
SOR step.
Definition: gsetc.hh:634
typename Overloads::ScalarType< std::decay_t< V > >::type Scalar
Element type of some SIMD type.
Definition: interface.hh:235
Simple iterative methods like Jacobi, Gauss-Seidel, SOR, SSOR, etc. in a generic way.
Incomplete LDL decomposition.
The incomplete LU factorization kernels.
Some handy generic functions for ISTL matrices.
Dune namespace.
Definition: alignedallocator.hh:13
void bildl_decompose(Matrix &A)
compute ILDL decomposition of a symmetric matrix A
Definition: ildl.hh:88
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:45
static void check([[maybe_unused]] const Matrix &mat)
Check whether the a matrix has diagonal values on blocklevel recursion levels.
Definition: matrixutils.hh:53
Statistics about the application of an inverse operator.
Definition: solver.hh:48
Category
Definition: solvercategory.hh:23
@ sequential
Category for sequential solvers.
Definition: solvercategory.hh:25
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:34
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (Apr 26, 22:29, 2024)