Dune Core Modules (2.8.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
15
16#include <dune/istl/solverregistry.hh>
17#include "preconditioner.hh"
18#include "solver.hh"
19#include "solvercategory.hh"
20#include "istlexception.hh"
21#include "matrixutils.hh"
22#include "gsetc.hh"
23#include "ildl.hh"
24#include "ilu.hh"
25
26
27namespace Dune {
70 template<class O, int c = -1>
72 public Preconditioner<typename O::domain_type, typename O::range_type>
73 {
74 public:
76 typedef typename O::domain_type domain_type;
78 typedef typename O::range_type range_type;
80 typedef typename range_type::field_type field_type;
84 typedef typename FieldTraits<scalar_field_type>::real_type real_field_type;
86 typedef O InverseOperator;
87
93 : inverse_operator_(inverse_operator)
94 {
95 if(c != -1 && SolverCategory::category(inverse_operator_) != c)
96 DUNE_THROW(InvalidStateException, "User-supplied solver category does not match that of the given inverse operator");
97 }
98
99 virtual void pre(domain_type&,range_type&)
100 {}
101
102 virtual void apply(domain_type& v, const range_type& d)
103 {
105 range_type copy(d);
106 inverse_operator_.apply(v, copy, res);
107 }
108
109 virtual void post(domain_type&)
110 {}
111
114 {
115 return SolverCategory::category(inverse_operator_);
116 }
117
118 private:
119 InverseOperator& inverse_operator_;
120 };
121
122 //=====================================================================
123 // Implementation of this interface for sequential ISTL-preconditioners
124 //=====================================================================
125
126
138 template<class M, class X, class Y, int l=1>
139 class SeqSSOR : public Preconditioner<X,Y> {
140 public:
142 typedef M matrix_type;
144 typedef X domain_type;
146 typedef Y range_type;
148 typedef typename X::field_type field_type;
152 typedef typename FieldTraits<scalar_field_type>::real_type real_field_type;
153
161 SeqSSOR (const M& A, int n, real_field_type w)
162 : _A_(A), _n(n), _w(w)
163 {
165 }
166
180 SeqSSOR (const std::shared_ptr<const AssembledLinearOperator<M,X,Y>>& A, const ParameterTree& configuration)
181 : SeqSSOR(A->getmat(), configuration)
182 {}
183
197 SeqSSOR (const M& A, const ParameterTree& configuration)
198 : SeqSSOR(A, configuration.get<int>("iterations",1), configuration.get<real_field_type>("relaxation",1.0))
199 {}
200
206 virtual void pre ([[maybe_unused]] X& x, [[maybe_unused]] Y& b)
207 {}
208
214 virtual void apply (X& v, const Y& d)
215 {
216 for (int i=0; i<_n; i++) {
217 bsorf(_A_,v,d,_w,BL<l>());
218 bsorb(_A_,v,d,_w,BL<l>());
219 }
220 }
221
227 virtual void post ([[maybe_unused]] X& x)
228 {}
229
232 {
234 }
235
236 private:
238 const M& _A_;
240 int _n;
243 };
244 DUNE_REGISTER_PRECONDITIONER("ssor", defaultPreconditionerBlockLevelCreator<Dune::SeqSSOR>());
245
246
258 template<class M, class X, class Y, int l=1>
259 class SeqSOR : public Preconditioner<X,Y> {
260 public:
262 typedef M matrix_type;
264 typedef X domain_type;
266 typedef Y range_type;
268 typedef typename X::field_type field_type;
272 typedef typename FieldTraits<scalar_field_type>::real_type real_field_type;
273
281 SeqSOR (const M& A, int n, real_field_type w)
282 : _A_(A), _n(n), _w(w)
283 {
285 }
286
300 SeqSOR (const std::shared_ptr<const AssembledLinearOperator<M,X,Y>>& A, const ParameterTree& configuration)
301 : SeqSOR(A->getmat(), configuration)
302 {}
303
317 SeqSOR (const M& A, const ParameterTree& configuration)
318 : SeqSOR(A, configuration.get<int>("iterations",1), configuration.get<real_field_type>("relaxation",1.0))
319 {}
320
326 virtual void pre ([[maybe_unused]] X& x, [[maybe_unused]] Y& b)
327 {}
328
334 virtual void apply (X& v, const Y& d)
335 {
336 this->template apply<true>(v,d);
337 }
338
347 template<bool forward>
348 void apply(X& v, const Y& d)
349 {
350 if(forward)
351 for (int i=0; i<_n; i++) {
352 bsorf(_A_,v,d,_w,BL<l>());
353 }
354 else
355 for (int i=0; i<_n; i++) {
356 bsorb(_A_,v,d,_w,BL<l>());
357 }
358 }
359
365 virtual void post ([[maybe_unused]] X& x)
366 {}
367
370 {
372 }
373
374 private:
376 const M& _A_;
378 int _n;
381 };
382 DUNE_REGISTER_PRECONDITIONER("sor", defaultPreconditionerBlockLevelCreator<Dune::SeqSOR>());
383
384
395 template<class M, class X, class Y, int l=1>
397 DUNE_REGISTER_PRECONDITIONER("gs", defaultPreconditionerBlockLevelCreator<Dune::SeqGS>());
398
409 template<class M, class X, class Y, int l=1>
411 public:
413 typedef M matrix_type;
415 typedef X domain_type;
417 typedef Y range_type;
419 typedef typename X::field_type field_type;
423 typedef typename FieldTraits<scalar_field_type>::real_type real_field_type;
424
432 SeqJac (const M& A, int n, real_field_type w)
433 : _A_(A), _n(n), _w(w)
434 {
436 }
437
451 SeqJac (const std::shared_ptr<const AssembledLinearOperator<M,X,Y>>& A, const ParameterTree& configuration)
452 : SeqJac(A->getmat(), configuration)
453 {}
454
468 SeqJac (const M& A, const ParameterTree& configuration)
469 : SeqJac(A, configuration.get<int>("iterations",1), configuration.get<real_field_type>("relaxation",1.0))
470 {}
471
477 virtual void pre ([[maybe_unused]] X& x, [[maybe_unused]] Y& b)
478 {}
479
485 virtual void apply (X& v, const Y& d)
486 {
487 for (int i=0; i<_n; i++) {
488 dbjac(_A_,v,d,_w,BL<l>());
489 }
490 }
491
497 virtual void post ([[maybe_unused]] X& x)
498 {}
499
502 {
504 }
505
506 private:
508 const M& _A_;
510 int _n;
512 real_field_type _w;
513 };
514 DUNE_REGISTER_PRECONDITIONER("jac", defaultPreconditionerBlockLevelCreator<Dune::SeqJac>());
515
516
517
529 template<class M, class X, class Y, int l=1>
531 public:
533 typedef typename std::remove_const<M>::type matrix_type;
535 typedef typename matrix_type :: block_type block_type;
537 typedef X domain_type;
539 typedef Y range_type;
540
542 typedef typename X::field_type field_type;
543
547 typedef typename FieldTraits<scalar_field_type>::real_type real_field_type;
548
551
559 SeqILU (const M& A, real_field_type w, const bool resort = false )
560 : SeqILU( A, 0, w, resort ) // construct ILU(0)
561 {
562 }
563
578 SeqILU (const std::shared_ptr<const AssembledLinearOperator<M,X,Y>>& A, const ParameterTree& configuration)
579 : SeqILU(A->getmat(), configuration)
580 {}
581
596 SeqILU(const M& A, const ParameterTree& config)
597 : SeqILU(A, config.get("n", 0),
598 config.get<real_field_type>("relaxation", 1.0),
599 config.get("resort", false))
600 {}
601
610 SeqILU (const M& A, int n, real_field_type w, const bool resort = false )
611 : ILU_(),
612 lower_(),
613 upper_(),
614 inv_(),
615 w_(w),
616 wNotIdentity_([w]{using std::abs; return abs(w - real_field_type(1)) > 1e-15;}() )
617 {
618 if( n == 0 )
619 {
620 // copy A
621 ILU_.reset( new matrix_type( A ) );
622 // create ILU(0) decomposition
623 ILU::blockILU0Decomposition( *ILU_ );
624 }
625 else
626 {
627 // create matrix in build mode
628 ILU_.reset( new matrix_type( A.N(), A.M(), matrix_type::row_wise) );
629 // create ILU(n) decomposition
630 ILU::blockILUDecomposition( A, n, *ILU_ );
631 }
632
633 if( resort )
634 {
635 // store ILU in simple CRS format
636 ILU::convertToCRS( *ILU_, lower_, upper_, inv_ );
637 ILU_.reset();
638 }
639 }
640
646 virtual void pre ([[maybe_unused]] X& x, [[maybe_unused]] Y& b)
647 {}
648
654 virtual void apply (X& v, const Y& d)
655 {
656 if( ILU_ )
657 {
658 ILU::blockILUBacksolve( *ILU_, v, d);
659 }
660 else
661 {
662 ILU::blockILUBacksolve(lower_, upper_, inv_, v, d);
663 }
664
665 if( wNotIdentity_ )
666 {
667 v *= w_;
668 }
669 }
670
676 virtual void post ([[maybe_unused]] X& x)
677 {}
678
681 {
683 }
684
685 protected:
687 std::unique_ptr< matrix_type > ILU_;
688
691 CRS upper_;
692 std::vector< block_type, typename matrix_type::allocator_type > inv_;
693
697 const bool wNotIdentity_;
698 };
699 DUNE_REGISTER_PRECONDITIONER("ilu", defaultPreconditionerBlockLevelCreator<Dune::SeqILU>());
700
701
710 template<class X, class Y>
711 class Richardson : public Preconditioner<X,Y> {
712 public:
714 typedef X domain_type;
716 typedef Y range_type;
718 typedef typename X::field_type field_type;
722 typedef typename FieldTraits<scalar_field_type>::real_type real_field_type;
723
730 _w(w)
731 {}
732
744 Richardson (const ParameterTree& configuration)
745 : Richardson(configuration.get<real_field_type>("relaxation", 1.0))
746 {}
747
753 virtual void pre ([[maybe_unused]] X& x, [[maybe_unused]] Y& b)
754 {}
755
761 virtual void apply (X& v, const Y& d)
762 {
763 v = d;
764 v *= _w;
765 }
766
772 virtual void post ([[maybe_unused]] X& x)
773 {}
774
777 {
779 }
780
781 private:
784 };
785 DUNE_REGISTER_PRECONDITIONER("richardson", [](auto tl, const auto& /* mat */, const ParameterTree& config){
786 using D = typename Dune::TypeListElement<1, decltype(tl)>::type;
787 using R = typename Dune::TypeListElement<2, decltype(tl)>::type;
788 return std::make_shared<Richardson<D,R>>(config);
789 });
790
791
802 template< class M, class X, class Y >
804 : public Preconditioner< X, Y >
805 {
806 typedef SeqILDL< M, X, Y > This;
808
809 public:
811 typedef std::remove_const_t< M > matrix_type;
813 typedef X domain_type;
815 typedef Y range_type;
817 typedef typename X::field_type field_type;
821 typedef typename FieldTraits<scalar_field_type>::real_type real_field_type;
822
835 SeqILDL (const std::shared_ptr<const AssembledLinearOperator<M,X,Y>>& A, const ParameterTree& configuration)
836 : SeqILDL(A->getmat(), configuration)
837 {}
838
851 SeqILDL(const matrix_type& A, const ParameterTree& config)
852 : SeqILDL(A, config.get<real_field_type>("relaxation", 1.0))
853 {}
854
863 explicit SeqILDL ( const matrix_type &A, real_field_type relax = real_field_type( 1 ) )
864 : decomposition_( A.N(), A.M(), matrix_type::random ),
865 relax_( relax )
866 {
867 // setup row sizes for lower triangular matrix
868 for( auto i = A.begin(), iend = A.end(); i != iend; ++i )
869 {
870 const auto &A_i = *i;
871 const auto ij = A_i.find( i.index() );
872 if( ij != A_i.end() )
873 decomposition_.setrowsize( i.index(), ij.offset()+1 );
874 else
875 DUNE_THROW( ISTLError, "diagonal entry missing" );
876 }
877 decomposition_.endrowsizes();
878
879 // setup row indices for lower triangular matrix
880 for( auto i = A.begin(), iend = A.end(); i != iend; ++i )
881 {
882 const auto &A_i = *i;
883 for( auto ij = A_i.begin(); ij.index() < i.index() ; ++ij )
884 decomposition_.addindex( i.index(), ij.index() );
885 decomposition_.addindex( i.index(), i.index() );
886 }
887 decomposition_.endindices();
888
889 // copy values of lower triangular matrix
890 auto i = A.begin();
891 for( auto row = decomposition_.begin(), rowend = decomposition_.end(); row != rowend; ++row, ++i )
892 {
893 auto ij = i->begin();
894 for( auto col = row->begin(), colend = row->end(); col != colend; ++col, ++ij )
895 *col = *ij;
896 }
897
898 // perform ILDL decomposition
899 bildl_decompose( decomposition_ );
900 }
901
903 void pre ([[maybe_unused]] X &x, [[maybe_unused]] Y &b) override
904 {}
905
907 void apply ( X &v, const Y &d ) override
908 {
909 bildl_backsolve( decomposition_, v, d, true );
910 v *= relax_;
911 }
912
914 void post ([[maybe_unused]] X &x) override
915 {}
916
919
920 private:
921 matrix_type decomposition_;
922 real_field_type relax_;
923 };
924 DUNE_REGISTER_PRECONDITIONER("ildl", defaultPreconditionerCreator<Dune::SeqILDL>());
925
928} // end namespace
929
930
931#endif
A linear operator exporting itself in matrix form.
Definition: operators.hh:107
derive error class from the base class in common
Definition: istlexception.hh:17
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:73
O::range_type range_type
The range type of the preconditioner.
Definition: preconditioners.hh:78
O::domain_type domain_type
The domain type of the preconditioner.
Definition: preconditioners.hh:76
virtual void post(domain_type &)
Clean up.
Definition: preconditioners.hh:109
range_type::field_type field_type
The field type of the preconditioner.
Definition: preconditioners.hh:80
FieldTraits< scalar_field_type >::real_type real_field_type
real scalar type underlying the field_type
Definition: preconditioners.hh:84
Simd::Scalar< field_type > scalar_field_type
scalar type underlying the field_type
Definition: preconditioners.hh:82
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: preconditioners.hh:113
virtual void pre(domain_type &, range_type &)
Prepare the preconditioner.
Definition: preconditioners.hh:99
InverseOperator2Preconditioner(InverseOperator &inverse_operator)
Construct the preconditioner from the solver.
Definition: preconditioners.hh:92
O InverseOperator
type of the wrapped inverse operator
Definition: preconditioners.hh:86
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:102
Abstract base class for all solvers.
Definition: solver.hh:97
size_type M() const
Return the number of columns.
Definition: matrix.hh:698
size_type N() const
Return the number of rows.
Definition: matrix.hh:693
Hierarchical structure of string parameters.
Definition: parametertree.hh:35
Base class for matrix free definition of preconditioners.
Definition: preconditioner.hh:30
Richardson preconditioner.
Definition: preconditioners.hh:711
X::field_type field_type
The field type of the preconditioner.
Definition: preconditioners.hh:718
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: preconditioners.hh:776
Y range_type
The range type of the preconditioner.
Definition: preconditioners.hh:716
virtual void pre(X &x, Y &b)
Prepare the preconditioner.
Definition: preconditioners.hh:753
Richardson(real_field_type w=1.0)
Constructor.
Definition: preconditioners.hh:729
virtual void post(X &x)
Clean up.
Definition: preconditioners.hh:772
FieldTraits< scalar_field_type >::real_type real_field_type
real scalar type underlying the field_type
Definition: preconditioners.hh:722
Simd::Scalar< field_type > scalar_field_type
scalar type underlying the field_type
Definition: preconditioners.hh:720
Richardson(const ParameterTree &configuration)
Constructor.
Definition: preconditioners.hh:744
X domain_type
The domain type of the preconditioner.
Definition: preconditioners.hh:714
virtual void apply(X &v, const Y &d)
Apply the precondioner.
Definition: preconditioners.hh:761
sequential ILDL preconditioner
Definition: preconditioners.hh:805
SeqILDL(const matrix_type &A, const ParameterTree &config)
Constructor.
Definition: preconditioners.hh:851
SeqILDL(const matrix_type &A, real_field_type relax=real_field_type(1))
constructor
Definition: preconditioners.hh:863
X domain_type
domain type of the preconditioner
Definition: preconditioners.hh:813
void post(X &x) override
Clean up.
Definition: preconditioners.hh:914
Y range_type
range type of the preconditioner
Definition: preconditioners.hh:815
std::remove_const_t< M > matrix_type
type of matrix the preconditioner is for
Definition: preconditioners.hh:811
void apply(X &v, const Y &d) override
Apply one step of the preconditioner to the system A(v)=d.
Definition: preconditioners.hh:907
FieldTraits< scalar_field_type >::real_type real_field_type
real scalar type underlying the field_type
Definition: preconditioners.hh:821
SeqILDL(const std::shared_ptr< const AssembledLinearOperator< M, X, Y > > &A, const ParameterTree &configuration)
Constructor.
Definition: preconditioners.hh:835
void pre(X &x, Y &b) override
Prepare the preconditioner.
Definition: preconditioners.hh:903
Simd::Scalar< field_type > scalar_field_type
scalar type underlying the field_type
Definition: preconditioners.hh:819
X::field_type field_type
field type of the preconditioner
Definition: preconditioners.hh:817
SolverCategory::Category category() const override
Category of the preconditioner (see SolverCategory::Category)
Definition: preconditioners.hh:918
Sequential ILU preconditioner.
Definition: preconditioners.hh:530
virtual void post(X &x)
Clean up.
Definition: preconditioners.hh:676
SeqILU(const M &A, int n, real_field_type w, const bool resort=false)
Constructor.
Definition: preconditioners.hh:610
virtual void pre(X &x, Y &b)
Prepare the preconditioner.
Definition: preconditioners.hh:646
virtual void apply(X &v, const Y &d)
Apply the preconditioner.
Definition: preconditioners.hh:654
ILU::CRS< block_type, typename M::allocator_type > CRS
type of ILU storage
Definition: preconditioners.hh:550
Y range_type
The range type of the preconditioner.
Definition: preconditioners.hh:539
CRS lower_
The ILU(n) decomposition of the matrix. As storage a CRS structure is used.
Definition: preconditioners.hh:690
const bool wNotIdentity_
true if w != 1.0
Definition: preconditioners.hh:697
SeqILU(const M &A, const ParameterTree &config)
Constructor.
Definition: preconditioners.hh:596
std::remove_const< M >::type matrix_type
The matrix type the preconditioner is for.
Definition: preconditioners.hh:533
matrix_type::block_type block_type
block type of matrix
Definition: preconditioners.hh:535
FieldTraits< scalar_field_type >::real_type real_field_type
real scalar type underlying the field_type
Definition: preconditioners.hh:547
X::field_type field_type
The field type of the preconditioner.
Definition: preconditioners.hh:542
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: preconditioners.hh:680
SeqILU(const M &A, real_field_type w, const bool resort=false)
Constructor.
Definition: preconditioners.hh:559
const real_field_type w_
The relaxation factor to use.
Definition: preconditioners.hh:695
SeqILU(const std::shared_ptr< const AssembledLinearOperator< M, X, Y > > &A, const ParameterTree &configuration)
Constructor.
Definition: preconditioners.hh:578
X domain_type
The domain type of the preconditioner.
Definition: preconditioners.hh:537
Simd::Scalar< field_type > scalar_field_type
scalar type underlying the field_type
Definition: preconditioners.hh:545
std::unique_ptr< matrix_type > ILU_
The ILU(n) decomposition of the matrix. As storage a BCRSMatrix is used.
Definition: preconditioners.hh:687
The sequential jacobian preconditioner.
Definition: preconditioners.hh:410
virtual void post(X &x)
Clean up.
Definition: preconditioners.hh:497
SeqJac(const M &A, const ParameterTree &configuration)
Constructor.
Definition: preconditioners.hh:468
virtual void apply(X &v, const Y &d)
Apply the preconditioner.
Definition: preconditioners.hh:485
M matrix_type
The matrix type the preconditioner is for.
Definition: preconditioners.hh:413
Simd::Scalar< field_type > scalar_field_type
scalar type underlying the field_type
Definition: preconditioners.hh:421
SeqJac(const std::shared_ptr< const AssembledLinearOperator< M, X, Y > > &A, const ParameterTree &configuration)
Constructor.
Definition: preconditioners.hh:451
X::field_type field_type
The field type of the preconditioner.
Definition: preconditioners.hh:419
virtual void pre(X &x, Y &b)
Prepare the preconditioner.
Definition: preconditioners.hh:477
X domain_type
The domain type of the preconditioner.
Definition: preconditioners.hh:415
FieldTraits< scalar_field_type >::real_type real_field_type
real scalar type underlying the field_type
Definition: preconditioners.hh:423
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: preconditioners.hh:501
SeqJac(const M &A, int n, real_field_type w)
Constructor.
Definition: preconditioners.hh:432
Y range_type
The range type of the preconditioner.
Definition: preconditioners.hh:417
Sequential SOR preconditioner.
Definition: preconditioners.hh:259
SeqSOR(const std::shared_ptr< const AssembledLinearOperator< M, X, Y > > &A, const ParameterTree &configuration)
Constructor.
Definition: preconditioners.hh:300
M matrix_type
The matrix type the preconditioner is for.
Definition: preconditioners.hh:262
FieldTraits< scalar_field_type >::real_type real_field_type
real scalar type underlying the field_type
Definition: preconditioners.hh:272
void apply(X &v, const Y &d)
Apply the preconditioner in a special direction.
Definition: preconditioners.hh:348
X domain_type
The domain type of the preconditioner.
Definition: preconditioners.hh:264
virtual void post(X &x)
Clean up.
Definition: preconditioners.hh:365
virtual void pre(X &x, Y &b)
Prepare the preconditioner.
Definition: preconditioners.hh:326
Simd::Scalar< field_type > scalar_field_type
scalar type underlying the field_type
Definition: preconditioners.hh:270
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: preconditioners.hh:369
virtual void apply(X &v, const Y &d)
Apply the preconditioner.
Definition: preconditioners.hh:334
Y range_type
The range type of the preconditioner.
Definition: preconditioners.hh:266
SeqSOR(const M &A, const ParameterTree &configuration)
Constructor.
Definition: preconditioners.hh:317
X::field_type field_type
The field type of the preconditioner.
Definition: preconditioners.hh:268
SeqSOR(const M &A, int n, real_field_type w)
Constructor.
Definition: preconditioners.hh:281
Sequential SSOR preconditioner.
Definition: preconditioners.hh:139
virtual void post(X &x)
Clean up.
Definition: preconditioners.hh:227
SeqSSOR(const std::shared_ptr< const AssembledLinearOperator< M, X, Y > > &A, const ParameterTree &configuration)
Constructor.
Definition: preconditioners.hh:180
SeqSSOR(const M &A, const ParameterTree &configuration)
Constructor.
Definition: preconditioners.hh:197
virtual SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: preconditioners.hh:231
X::field_type field_type
The field type of the preconditioner.
Definition: preconditioners.hh:148
Simd::Scalar< field_type > scalar_field_type
scalar type underlying the field_type
Definition: preconditioners.hh:150
X domain_type
The domain type of the preconditioner.
Definition: preconditioners.hh:144
M matrix_type
The matrix type the preconditioner is for.
Definition: preconditioners.hh:142
virtual void apply(X &v, const Y &d)
Apply the preconditioner.
Definition: preconditioners.hh:214
virtual void pre(X &x, Y &b)
Prepare the preconditioner.
Definition: preconditioners.hh:206
Y range_type
The range type of the preconditioner.
Definition: preconditioners.hh:146
FieldTraits< scalar_field_type >::real_type real_field_type
real scalar type underlying the field_type
Definition: preconditioners.hh:152
SeqSSOR(const M &A, int n, real_field_type w)
Constructor.
Definition: preconditioners.hh:161
#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:644
void dbjac(const M &A, X &x, const Y &b, const K &w)
Jacobi step.
Definition: gsetc.hh:656
void bsorf(const M &A, X &x, const Y &b, const K &w)
SOR step.
Definition: gsetc.hh:632
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.
The incomplete LU factorization kernels.
Some handy generic functions for ISTL matrices.
Dune namespace.
Definition: alignedallocator.hh:11
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:51
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
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 21, 23:30, 2024)