DUNE-FEM (unstable)

istlinverseoperators.hh
1#ifndef DUNE_FEM_SOLVER_ISTLINVERSEOPERATORS_HH
2#define DUNE_FEM_SOLVER_ISTLINVERSEOPERATORS_HH
3
4#include <dune/fem/function/blockvectorfunction.hh>
5#include <dune/fem/function/common/scalarproducts.hh>
6#include <dune/fem/operator/common/operator.hh>
7#include <dune/fem/io/parameter.hh>
8#include <dune/fem/misc/mpimanager.hh>
9
10#include <dune/fem/solver/parameter.hh>
11#include <dune/fem/solver/inverseoperatorinterface.hh>
12
13#if HAVE_DUNE_ISTL
15
16#include <dune/fem/operator/linear/istladapter.hh>
17#include <dune/fem/operator/linear/istloperator.hh>
18
20#include <dune/istl/solvers.hh>
21#include <dune/istl/preconditioner.hh>
22
23namespace Dune
24{
25
26 namespace Fem
27 {
28
29 // wrapper for Fem Preconditioners (Operators acting as preconditioners) into ISTL preconditioners
30 template< class Preconditioner >
31 class ISTLPreconditionAdapter
32 : public Dune::Preconditioner< typename Preconditioner::RangeFunctionType::DofStorageType, typename Preconditioner::DomainFunctionType::DofStorageType >
33 {
34 typedef ISTLPreconditionAdapter< Preconditioner > ThisType;
36
37 typedef typename Preconditioner::DomainFunctionType DomainFunctionType;
38 typedef typename Preconditioner::RangeFunctionType RangeFunctionType;
39
40 public:
41 typedef typename BaseType::domain_type domain_type;
42 typedef typename BaseType::range_type range_type;
43 typedef typename BaseType::field_type field_type;
44
45 typedef typename DomainFunctionType::DiscreteFunctionSpaceType DomainFunctionSpaceType;
46 typedef typename RangeFunctionType::DiscreteFunctionSpaceType RangeFunctionSpaceType;
47
48 ISTLPreconditionAdapter ( const Preconditioner *precon, const DomainFunctionSpaceType &domainSpace, const RangeFunctionSpaceType &rangeSpace )
49 : precon_( precon ),
50 domainSpace_( domainSpace ),
51 rangeSpace_( rangeSpace )
52 {}
53
54 // pre and post do nothing here
55 virtual void pre ( domain_type &x, range_type &y ) override {}
56 virtual void post ( domain_type &x ) override {}
57
58 virtual void apply ( domain_type &x, const range_type &y ) override
59 {
60 // no precon
61 if( !precon_ )
62 {
63 x = y;
64 }
65 else
66 {
67 // note: ISTL switches the arguments !!!
68 // it is assumed that we have a left preconditioner
69 RangeFunctionType px( "ISTLPreconditionAdapter::apply::x", rangeSpace_, x );
70 DomainFunctionType py( "ISTLPreconditionAdapter::apply::y", domainSpace_, y );
71
72 (*precon_)( px, py );
73 }
74 }
75
76 SolverCategory::Category category () const override { return SolverCategory::sequential; }
77
78 protected:
79 const Preconditioner *precon_;
80 const DomainFunctionSpaceType &domainSpace_;
81 const RangeFunctionSpaceType &rangeSpace_;
82 };
83
84
85 template< class BlockVector >
86 struct ISTLSolverReduction
87 {
88 ISTLSolverReduction ( std::shared_ptr< ISTLSolverParameter > parameter )
89 : parameter_( parameter )
90 {}
91
92 double operator() ( const Dune::LinearOperator< BlockVector, BlockVector > &op,
94 const BlockVector &rhs, const BlockVector &x ) const
95 {
96
97 if( parameter_->errorMeasure() == 0)
98 {
99 BlockVector residuum( rhs );
100 op.applyscaleadd( -1., x, residuum );
101 const double res = scp.norm( residuum );
102 return (res > 0 ? parameter_->tolerance() / res : 1e-3);
103 }
104 else
105 return parameter_->tolerance();
106 }
107
108 private:
109 std::shared_ptr<ISTLSolverParameter> parameter_;
110 };
111
112
113 struct ISTLInverseOperatorMethods
114 {
115 static std::vector< int > supportedSolverMethods() {
116 return std::vector< int > ({
117 SolverParameter::gmres, // default solver
118 SolverParameter::cg,
119 SolverParameter::bicgstab,
120 SolverParameter::minres,
121 SolverParameter::gradient,
122 SolverParameter::loop,
123 SolverParameter::superlu,
124 SolverParameter::fgmres
125 });
126 }
127 };
128
129 template< int method,
130 class X,
131 class Reduction = ISTLSolverReduction< X > >
132 struct ISTLSolverAdapter
133 {
134 typedef Reduction ReductionType;
135
136 typedef X domain_type;
137 typedef X range_type;
138
139 ISTLSolverAdapter ( const ReductionType &reduction, std::shared_ptr<ISTLSolverParameter> parameter )
140 : reduction_( reduction ),
141 method_( method < 0 ? parameter->solverMethod( ISTLInverseOperatorMethods::supportedSolverMethods() ) : method ),
142 parameter_( parameter )
143 {}
144
145 template<class Op, class ScP, class PC >
146 void operator () ( Op& op, ScP &scp, PC &pc,
147 range_type &rhs, domain_type &x,
148 Dune::InverseOperatorResult &result ) const
149 {
150 int verbosity = (Parameter::verbose( Parameter::solverStatistics ) && parameter_->verbose()) ? 2 : 0;
151 if ( verbosity && Parameter::verbose( Parameter::extendedStatistics ) )
152 verbosity += 2;
153
154 int maxIterations = std::min( std::numeric_limits< int >::max(), parameter_->maxIterations() );
155 if( method_ == SolverParameter::cg )
156 {
157 typedef Dune::CGSolver< X > SolverType;
158 SolverType solver( op, scp, pc, reduction_( op, scp, rhs, x ), maxIterations, verbosity );
159 solver.apply( x, rhs, result );
160 return ;
161 }
162 else if( method_ == SolverParameter::bicgstab )
163 {
164 typedef Dune::BiCGSTABSolver< X > SolverType;
165 SolverType solver( op, scp, pc, reduction_( op, scp, rhs, x ), maxIterations, verbosity );
166 solver.apply( x, rhs, result );
167 return ;
168 }
169 else if( method_ == SolverParameter::gmres )
170 {
171 typedef Dune::RestartedGMResSolver< X > SolverType;
172 SolverType solver( op, scp, pc, reduction_( op, scp, rhs, x ), parameter_->gmresRestart(), maxIterations, verbosity );
173 solver.apply( x, rhs, result );
174 return ;
175 }
176 else if( method_ == SolverParameter::fgmres )
177 {
179 SolverType solver( op, scp, pc, reduction_( op, scp, rhs, x ), parameter_->gmresRestart(), maxIterations, verbosity );
180 solver.apply( x, rhs, result );
181 return ;
182 }
183 else if( method_ == SolverParameter::minres )
184 {
185 typedef Dune::MINRESSolver< X > SolverType;
186 SolverType solver( op, scp, pc, reduction_( op, scp, rhs, x ), maxIterations, verbosity );
187 solver.apply( x, rhs, result );
188 }
189 else if( method_ == SolverParameter::gradient )
190 {
191 typedef Dune::GradientSolver< X > SolverType;
192 SolverType solver( op, scp, pc, reduction_( op, scp, rhs, x ), maxIterations, verbosity );
193 solver.apply( x, rhs, result );
194 return ;
195 }
196 else if( method_ == SolverParameter::loop )
197 {
198 typedef Dune::LoopSolver< X > SolverType;
199 SolverType solver( op, scp, pc, reduction_( op, scp, rhs, x ), maxIterations, verbosity );
200 solver.apply( x, rhs, result );
201 return ;
202 }
203 else if( method_ == SolverParameter::superlu )
204 {
205 callSuperLU( op, rhs, x, result );
206 }
207 else
208 {
209 DUNE_THROW(NotImplemented,"ISTLSolverAdapter::operator(): wrong method solver identifier" << method_ );
210 }
211 }
212
213 void setMaxLinearIterations( unsigned int maxIterations ) { parameter_->setMaxIterations(maxIterations); }
214 void setMaxIterations( unsigned int maxIterations ) { parameter_->setMaxIterations(maxIterations); }
215
216 std::shared_ptr<ISTLSolverParameter> parameter () const { return parameter_; }
217
218 protected:
219 template< class ImprovedMatrix >
220 void callSuperLU ( ISTLParallelMatrixAdapterInterface< ImprovedMatrix >& op,
221 range_type &rhs, domain_type &x,
222 Dune::InverseOperatorResult &result ) const
223 {
224#if HAVE_SUPERLU
225 const int verbosity = (Dune::Fem::Parameter::verbose( Parameter::solverStatistics ) && parameter_->verbose()) ? 2 : 0;
226 typedef typename ImprovedMatrix :: BaseType Matrix;
227 const ImprovedMatrix& matrix = op.getmat();
228 SuperLU< Matrix > solver( matrix, verbosity );
229 solver.apply( x, rhs, result );
230#else
231 DUNE_THROW(NotImplemented,"ISTLSolverAdapter::callSuperLU: SuperLU solver selected but SuperLU not available!");
232#endif
233 }
234
235 template< class Op >
236 void callSuperLU ( Op& op, range_type &rhs, domain_type &x,
237 Dune::InverseOperatorResult &result ) const
238 {
239 DUNE_THROW(NotImplemented,"ISTLSolverAdapter::callSuperLU: SuperLU only works for AssembledLinearOperators!");
240 }
241
242 ReductionType reduction_;
243 const int method_;
244
245 std::shared_ptr<ISTLSolverParameter> parameter_;
246 };
247
248
249
250 // ISTLInverseOperator
251 // -------------------
252
253 template< class DiscreteFunction, int method = -1,
254 class Preconditioner = Fem::Operator< DiscreteFunction, DiscreteFunction > >
255 class ISTLInverseOperator;
256
257 template< class DiscreteFunction, int method, class Preconditioner >
258 struct ISTLInverseOperatorTraits
259 {
260 typedef DiscreteFunction DiscreteFunctionType;
262 typedef Preconditioner PreconditionerType;
263
264 typedef ISTLBlockVectorDiscreteFunction< typename DiscreteFunction::DiscreteFunctionSpaceType > SolverDiscreteFunctionType ;
265 typedef Dune::Fem::ISTLLinearOperator< DiscreteFunction, DiscreteFunction > AssembledOperatorType;
266
267 typedef ISTLInverseOperator< DiscreteFunction, method, Preconditioner > InverseOperatorType;
268 typedef ISTLSolverParameter SolverParameterType;
269 };
270
271 template< class DiscreteFunction, int method, class Preconditioner >
272 class ISTLInverseOperator : public InverseOperatorInterface<
273 ISTLInverseOperatorTraits< DiscreteFunction,
274 method, Preconditioner > >
275 {
276 typedef ISTLInverseOperatorTraits< DiscreteFunction, method, Preconditioner > Traits;
277 typedef InverseOperatorInterface< Traits > BaseType;
278
279 friend class InverseOperatorInterface< Traits >;
280 public:
281 typedef typename BaseType::DomainFunctionType DomainFunctionType;
282 typedef typename BaseType::RangeFunctionType RangeFunctionType;
283 typedef typename BaseType::OperatorType OperatorType;
284 typedef typename BaseType::PreconditionerType PreconditionerType;
285 typedef typename BaseType::AssembledOperatorType AssembledOperatorType;
286
287 using BaseType :: bind;
288 using BaseType :: unbind;
289 using BaseType :: setMaxIterations;
290 using BaseType :: setMaxLinearIterations;
291
292 protected:
293 typedef typename Traits::SolverDiscreteFunctionType SolverDiscreteFunctionType;
294
295 typedef typename DomainFunctionType :: DiscreteFunctionSpaceType
296 DiscreteFunctionSpaceType;
297
298
299 typedef typename SolverDiscreteFunctionType::ScalarProductType ParallelScalarProductType;
300 typedef typename SolverDiscreteFunctionType::DofStorageType BlockVectorType;
301
302 typedef ISTLSolverAdapter< method, BlockVectorType > SolverAdapterType;
303 typedef typename SolverAdapterType::ReductionType ReductionType;
304 public:
305
306 //non-deprecated constructors
307 ISTLInverseOperator ( const ISTLSolverParameter & parameter = ISTLSolverParameter() )
308 : BaseType( parameter ), solverAdapter_( ReductionType( parameter_ ), parameter_ )
309 {}
310
311 ISTLInverseOperator ( const OperatorType &op,
312 const ISTLSolverParameter & parameter = ISTLSolverParameter() )
313 : ISTLInverseOperator ( parameter )
314 {
315 bind( op );
316 }
317
318 ISTLInverseOperator ( const OperatorType &op, PreconditionerType &preconditioner,
319 const ISTLSolverParameter & parameter = ISTLSolverParameter() )
320 : ISTLInverseOperator( parameter )
321 {
322 bind( op, preconditioner );
323 }
324
325 protected:
326 // apply for arbitrary domain function type and matching range function type
327 template <class DomainFunction>
328 int apply( const DomainFunction& u, SolverDiscreteFunctionType& w ) const
329 {
330 auto& scp = w.scalarProduct();
331 // u may not be a discrete function, therefore use w.space()
332 const DiscreteFunctionSpaceType& space = w.space();
333
335 std::shared_ptr< ISTLPreconditionerType > istlPre;
336 if constexpr (std::is_same< SolverDiscreteFunctionType, RangeFunctionType >::value )
337 {
338 typedef ISTLPreconditionAdapter< PreconditionerType > ISTLPreconditionerAdapterType;
339 istlPre.reset( new ISTLPreconditionerAdapterType( preconditioner_, space, space ) );
340 }
341
342 if( assembledOperator_ )
343 {
344 auto& matrix = assembledOperator_->matrixAdapter( *(solverAdapter_.parameter()) );
345 // if preconditioner_ was set use that one, otherwise the one from the matrix object
346 ISTLPreconditionerType& matrixPre = matrix.preconditionAdapter();
347 ISTLPreconditionerType& precon = ( preconditioner_ ) ? (*istlPre) : matrixPre;
348 return solve( matrix, scp, precon, u, w );
349 }
350
351 if constexpr (std::is_same< SolverDiscreteFunctionType, RangeFunctionType >::value )
352 {
353 assert( operator_ );
354 assert( istlPre );
355 typedef ISTLLinearOperatorAdapter< OperatorType > ISTLOperatorType;
356 ISTLOperatorType istlOperator( *operator_, space, space );
357 return solve( istlOperator, scp, *istlPre, u, w );
358 }
359
360 DUNE_THROW(InvalidStateException,"ISTLInverseOperator::apply: No valid operator found!");
361 return -1;
362 }
363
365 template< class OperatorAdapter, class ISTLPreconditioner, class DomainFunction >
366 int solve ( OperatorAdapter &istlOperator, ParallelScalarProductType &scp,
367 ISTLPreconditioner &preconditioner,
368 const DomainFunction& u,
369 SolverDiscreteFunctionType& w ) const
370 {
371 if( ! rhs_ )
372 {
373 // u may not be a discrete function, therefore use w.space()
374 rhs_.reset( new SolverDiscreteFunctionType( "ISTLInvOp::rhs", w.space() ) );
375 rightHandSideCopied_ = false;
376 }
377
378 if( ! rightHandSideCopied_ )
379 {
380 // copy right hand side since ISTL solvers seem to modify it
381 rhs_->assign( u );
382 rightHandSideCopied_ = true;
383 }
384
386 solverAdapter_.setMaxIterations( parameter_->maxIterations() );
387 solverAdapter_( istlOperator, scp, preconditioner, rhs_->blockVector(), w.blockVector(), result );
388 return (result.converged) ? result.iterations : -(result.iterations);
389 }
390
391 using BaseType :: operator_;
392 using BaseType :: assembledOperator_;
393 using BaseType :: preconditioner_;
394
395 using BaseType :: rhs_;
396 using BaseType :: x_;
397
398 using BaseType :: rightHandSideCopied_;
399 using BaseType :: parameter_;
400
401 mutable SolverAdapterType solverAdapter_;
402 };
403
404 } // namespace Fem
405
406} // namespace Dune
407
408#endif // #if HAVE_ISTL
409
410#endif // #ifndef DUNE_FEM_SOLVER_ISTLINVERSEOPERATORS_HH
Bi-conjugate Gradient Stabilized (BiCG-STAB)
Definition: solvers.hh:420
conjugate gradient method
Definition: solvers.hh:193
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
gradient method
Definition: solvers.hh:124
A linear operator.
Definition: operators.hh:68
virtual void applyscaleadd(field_type alpha, const X &x, Y &y) const =0
apply operator to x, scale and add:
Preconditioned loop solver.
Definition: solvers.hh:59
Minimal Residual Method (MINRES)
Definition: solvers.hh:610
Base class for matrix free definition of preconditioners.
Definition: preconditioner.hh:33
implements the Flexible Generalized Minimal Residual (FGMRes) method (right preconditioned)
Definition: solvers.hh:1130
implements the Generalized Minimal Residual (GMRes) method
Definition: solvers.hh:828
Base class for scalar product and norm computation.
Definition: scalarproducts.hh:52
virtual real_type norm(const X &x) const
Norm of a right-hand side vector. The vector must be consistent on the interior+border partition.
Definition: scalarproducts.hh:71
Various macros to work with Dune module version numbers.
Define base class for scalar product and norm.
#define DUNE_THROW(E,...)
Definition: exceptions.hh:314
constexpr auto max
Function object that returns the greater of the given values.
Definition: hybridutilities.hh:489
Dune namespace
Definition: alignedallocator.hh:13
Implementations of the inverse operator interface.
Statistics about the application of an inverse operator.
Definition: solver.hh:50
int iterations
Number of iterations.
Definition: solver.hh:69
bool converged
True if convergence criterion has been met.
Definition: solver.hh:75
Category
Definition: solvercategory.hh:23
@ sequential
Category for sequential solvers.
Definition: solvercategory.hh:25
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden & Uni Heidelberg  |  generated with Hugo v0.111.3 (May 31, 22:31, 2026)