DUNE-FEM (unstable)

krylovinverseoperators.hh
1#ifndef DUNE_FEM_SOLVER_INVERSEOPERATORS_HH
2#define DUNE_FEM_SOLVER_INVERSEOPERATORS_HH
3
4#include <dune/fem/function/adaptivefunction.hh>
5#include <dune/fem/operator/common/operator.hh>
6#include <dune/fem/io/parameter.hh>
7#include <dune/fem/solver/parameter.hh>
8
9#include <dune/fem/solver/cginverseoperator.hh>
10#include <dune/fem/solver/fempreconditioning.hh>
11
12#include <dune/fem/solver/linear/gmres.hh>
13#include <dune/fem/solver/linear/bicgstab.hh>
14#include <dune/fem/solver/linear/cg.hh>
15#include <dune/fem/solver/inverseoperatorinterface.hh>
16
17#include <dune/fem/misc/mpimanager.hh>
18
19namespace Dune
20{
21
22 namespace Fem
23 {
24
25 // KrylovInverseOperator
26 // ---------------------
27
28 template< class DiscreteFunction, int method = -1 >
29 class KrylovInverseOperator;
30
31 template< class DiscreteFunction, int method >
32 struct KrylovInverseOperatorTraits
33 {
34 typedef DiscreteFunction DiscreteFunctionType;
36 typedef OperatorType PreconditionerType;
37
38 typedef OperatorType AssembledOperatorType;
39 typedef DiscreteFunction SolverDiscreteFunctionType ;
40
41 typedef KrylovInverseOperator< DiscreteFunction, method > InverseOperatorType;
42 typedef SolverParameter SolverParameterType;
43 };
44
45
46 template< class DiscreteFunction, int method >
47 class KrylovInverseOperator
48 : public InverseOperatorInterface< KrylovInverseOperatorTraits< DiscreteFunction, method > >
49 {
50 typedef KrylovInverseOperatorTraits< DiscreteFunction, method > Traits;
51 typedef InverseOperatorInterface< Traits > BaseType;
52
53 friend class InverseOperatorInterface< Traits >;
54 public:
55 typedef typename BaseType::DomainFunctionType DomainFunctionType;
56 typedef typename BaseType::RangeFunctionType RangeFunctionType;
57 typedef typename BaseType::OperatorType OperatorType;
58 typedef typename BaseType::PreconditionerType PreconditionerType;
59 typedef typename BaseType::AssembledOperatorType AssembledOperatorType;
60
61 using BaseType :: bind;
62 using BaseType :: unbind;
63 using BaseType :: setMaxLinearIterations;
64 using BaseType :: setMaxIterations;
65
66 public:
67 KrylovInverseOperator ( const OperatorType &op, const SolverParameter &parameter = SolverParameter(Parameter::container()) )
68 : KrylovInverseOperator( parameter )
69 {
70 bind( op );
71 }
72
73 KrylovInverseOperator ( const OperatorType &op, const PreconditionerType& preconditioner,
74 const SolverParameter &parameter = SolverParameter(Parameter::container()) )
75 : KrylovInverseOperator( op, &preconditioner, parameter )
76 {}
77
79 KrylovInverseOperator ( const SolverParameter &parameter = SolverParameter(Parameter::container()) )
80 : BaseType( parameter ),
81 precondObj_(),
82 verbose_( parameter.verbose() ),
83 method_( method < 0 ? parameter.solverMethod( supportedSolverMethods() ) : method ),
84 precondMethod_( parameter.preconditionMethod( supportedPreconditionMethods() ) )
85 {
86 // assert( parameter_->errorMeasure() == 0 );
87 }
88
89 template <class Operator>
90 void bind ( const Operator &op )
91 {
92 if( precondMethod_ && std::is_base_of< AssembledOperator< DomainFunctionType, DomainFunctionType >, Operator > :: value )
93 {
94 createPreconditioner( op );
95 }
96
97 if( precondObj_ )
98 BaseType::bind( op, *precondObj_ );
99 else
100 BaseType::bind( op );
101 }
102
103 static std::vector< int > supportedSolverMethods() {
104 return std::vector< int > ({ SolverParameter::gmres, // default solver
105 SolverParameter::cg,
106 SolverParameter::bicgstab });
107 }
108
109 static std::vector< int > supportedPreconditionMethods() {
110 return std::vector< int > ({ SolverParameter::none,
111 SolverParameter::sor,
112 SolverParameter::ssor,
113 SolverParameter::gauss_seidel,
114 SolverParameter::jacobi } );
115 }
116
117 protected:
118 int apply( const DomainFunctionType &u, RangeFunctionType &w ) const
119 {
120 std::ostream* os = nullptr;
121 // only set output when general verbose mode is enabled
122 // (basically to avoid output on every rank)
123 if( verbose_ && Parameter :: verbose( Parameter::solverStatistics ) )
124 {
125 os = &std::cout;
126 if( method_ == SolverParameter::gmres )
127 std::cout << "Fem::GMRES";
128 else if ( method_ == SolverParameter::bicgstab )
129 std::cout << "Fem::BiCGstab";
130 else if ( method_ == SolverParameter::cg )
131 std::cout << "Fem::CG";
132
133 std::cout << " preconditioning=";
134 if ( precondMethod_ == SolverParameter::none && BaseType::preconditioner_ )
135 std::cout << "custom" << std::endl;
136 else
137 std::cout << SolverParameter::preconditionMethodTable( precondMethod_ ) << std::endl;
138 }
139
140 int numIter = 0;
141
142 if( method_ == SolverParameter::gmres )
143 {
144 if( v_.empty() )
145 {
146 const int extra = ( preconditioner_ ) ? 2 : 1 ;
147 v_.reserve( parameter_->gmresRestart()+extra );
148 for( int i=0; i<=parameter_->gmresRestart(); ++i )
149 {
150 v_.emplace_back( DiscreteFunction( "GMRes::v", u.space() ) );
151 }
152 if( preconditioner_ )
153 v_.emplace_back( DiscreteFunction( "GMRes::z", u.space() ) );
154 }
155
156 // if solver convergence failed numIter will be negative
157 numIter = LinearSolver::gmres( *operator_, preconditioner_,
158 v_, w, u, parameter_->gmresRestart(),
159 parameter_->tolerance(), parameter_->maxIterations(),
160 parameter_->errorMeasure(), os );
161 }
162 else if( method_ == SolverParameter::bicgstab )
163 {
164 if( v_.empty() )
165 {
166 v_.emplace_back( DomainFunctionType( "BiCGStab::r", u.space() ) );
167 v_.emplace_back( DomainFunctionType( "BiCGStab::r*", u.space() ) );
168 v_.emplace_back( DomainFunctionType( "BiCGStab::p", u.space() ) );
169 v_.emplace_back( DomainFunctionType( "BiCGStab::s", u.space() ) );
170 v_.emplace_back( DomainFunctionType( "BiCGStab::tmp", u.space() ) );
171 if( preconditioner_ )
172 v_.emplace_back( DomainFunctionType( "BiCGStab::z", u.space() ) );
173 }
174
175 // if solver convergence failed numIter will be negative
176 numIter = LinearSolver::bicgstab( *operator_, preconditioner_,
177 v_, w, u,
178 parameter_->tolerance(), parameter_->maxIterations(),
179 parameter_->errorMeasure(), os );
180 }
181 else if( method_ == SolverParameter::cg )
182 {
183 if( v_.empty() )
184 {
185 v_.emplace_back( DomainFunctionType( "CG::h", u.space() ) );
186 v_.emplace_back( DomainFunctionType( "CG::r", u.space() ) );
187 v_.emplace_back( DomainFunctionType( "CG::p", u.space() ) );
188
189 if( preconditioner_ )
190 {
191 v_.emplace_back( DomainFunctionType( "CG::s", u.space() ) );
192 v_.emplace_back( DomainFunctionType( "CG::q", u.space() ) );
193 }
194 }
195
196 // if solver convergence failed numIter will be negative
197 numIter = LinearSolver::cg( *operator_, preconditioner_,
198 v_, w, u,
199 parameter_->tolerance(), parameter_->maxIterations(),
200 parameter_->errorMeasure(), os );
201 }
202 else
203 {
204 DUNE_THROW(InvalidStateException,"KrylovInverseOperator: invalid method " << method_ );
205 }
206
207 return numIter;
208 }
209
210 template <class LinearOperator>
211 KrylovInverseOperator ( const LinearOperator &op,
212 const PreconditionerType *preconditioner,
213 const SolverParameter &parameter = SolverParameter(Parameter::container()) )
214 : KrylovInverseOperator( parameter )
215 {
216 if( preconditioner )
217 BaseType::bind( op, *preconditioner );
218 else
219 BaseType::bind( op );
220 }
221
222 protected:
223 template <class LinearOperator>
224 void createPreconditioner( const LinearOperator &op )
225 {
226 if( precondMethod_ > 0 && std::is_base_of< AssembledOperator< DomainFunctionType, DomainFunctionType >, LinearOperator > :: value )
227 {
228 const int n = parameter_->preconditionerIteration();
229 const double w = parameter_->relaxation();
230
231 if( precondMethod_ == SolverParameter::jacobi )
232 {
233 // create diagonal preconditioner
234 precondObj_.reset( new FemJacobiPreconditioning< DomainFunctionType, LinearOperator >( op, n, w ) );
235 }
236 else if( precondMethod_ == SolverParameter::gauss_seidel )
237 {
238 // create diagonal preconditioner
239 precondObj_.reset( new FemGaussSeidelPreconditioning< DomainFunctionType, LinearOperator >( op, n, w ) );
240 }
241 else if( precondMethod_ == SolverParameter::sor )
242 {
243 // create diagonal preconditioner
244 precondObj_.reset( new FemSORPreconditioning< DomainFunctionType, LinearOperator >( op, n, w ) );
245 }
246 else if( precondMethod_ == SolverParameter::ssor )
247 {
248 // create diagonal preconditioner
249 precondObj_.reset( new FemSSORPreconditioning< DomainFunctionType, LinearOperator >( op, n, w ) );
250 }
251
252 // if preconditioner was created, set pointer
253 if( precondObj_ )
254 {
255 preconditioner_ = precondObj_.operator->();
256 }
257 }
258 }
259
260 using BaseType :: operator_;
261 using BaseType :: preconditioner_;
262
263 using BaseType :: parameter_;
264 using BaseType :: iterations_;
265
266 std::unique_ptr< PreconditionerType > precondObj_;
267
268 mutable std::vector< DomainFunctionType > v_;
269
270 const bool verbose_;
271
272 const int method_;
273 const int precondMethod_;
274 };
275
276
277 // CgInverseOperator
278 // -----------------
279
280 template< class DiscreteFunction >
281 using CgInverseOperator = KrylovInverseOperator< DiscreteFunction, SolverParameter :: cg >;
282
283
284 // BicgstabInverseOperator
285 // -----------------------
286
287 template< class DiscreteFunction >
288 using BicgstabInverseOperator = KrylovInverseOperator< DiscreteFunction, SolverParameter :: bicgstab >;
289
290
291 // GmresInverseOperator
292 // --------------------
293
294 template< class DiscreteFunction >
295 using GmresInverseOperator = KrylovInverseOperator< DiscreteFunction, SolverParameter :: gmres >;
296
297
298 // ParDGGeneralizedMinResInverseOperator
299 // -------------------------------------
300
301 template< class DiscreteFunction >
302 using ParDGGeneralizedMinResInverseOperator = GmresInverseOperator< DiscreteFunction >;
303
304 // ParDGBiCGStabInverseOperator
305 // ----------------------------
306
307 template< class DiscreteFunction >
308 using ParDGBiCGStabInverseOperator = BicgstabInverseOperator< DiscreteFunction >;
309
310 template <class DiscreteFunctionType, class OpType >
311 using OEMCGOp = CgInverseOperator< DiscreteFunctionType >;
312
313 template <class DiscreteFunctionType, class OpType >
314 using OEMBICGSTABOp = BicgstabInverseOperator< DiscreteFunctionType >;
315
316 template <class DiscreteFunctionType, class OpType >
317 using OEMBICGSQOp = BicgstabInverseOperator< DiscreteFunctionType >;
318
319 template <class DiscreteFunctionType, class OpType >
320 using OEMGMRESOp = GmresInverseOperator< DiscreteFunctionType >;
321
322 template <class DiscreteFunctionType, class OpType >
323 using GMRESOp = GmresInverseOperator< DiscreteFunctionType >;
324
325 } // namespace Fem
326
327} // namespace Dune
328
329#endif // DUNE_FEM_SOLVER_INVERSEOPERATORS_HH
Inverse operator base on CG method. Uses a runtime parameter fem.preconditioning which enables diagon...
Definition: cginverseoperator.hh:408
static bool verbose()
obtain the cached value for fem.verbose with default verbosity level 2
Definition: parameter.hh:466
forward declaration
Definition: discretefunction.hh:51
#define DUNE_THROW(E,...)
Definition: exceptions.hh:314
constexpr GeometryType none(unsigned int dim)
Returns a GeometryType representing a singular of dimension dim.
Definition: type.hh:471
Dune namespace
Definition: alignedallocator.hh:13
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden & Uni Heidelberg  |  generated with Hugo v0.111.3 (Feb 16, 23:40, 2026)