1 #ifndef DUNE_FEM_SOLVER_ISTLINVERSEOPERATORS_HH 2 #define DUNE_FEM_SOLVER_ISTLINVERSEOPERATORS_HH 4 #include <dune/common/nullptr.hh> 14 #include <dune/istl/scalarproducts.hh> 15 #include <dune/istl/solvers.hh> 16 #include <dune/istl/preconditioner.hh> 25 template<
class Preconditioner >
26 class ISTLPreconditionAdapter
27 :
public Dune::Preconditioner< typename Preconditioner::RangeFunctionType::DofStorageType, typename Preconditioner::DomainFunctionType::DofStorageType >
29 typedef ISTLPreconditionAdapter< Preconditioner > ThisType;
30 typedef Dune::Preconditioner< typename Preconditioner::RangeFunctionType::DofStorageType, typename Preconditioner::DomainFunctionType::DofStorageType > BaseType;
32 typedef typename Preconditioner::DomainFunctionType DomainFunctionType;
33 typedef typename Preconditioner::RangeFunctionType RangeFunctionType;
36 enum {category=SolverCategory::sequential};
38 typedef typename BaseType::domain_type domain_type;
39 typedef typename BaseType::range_type range_type;
40 typedef typename BaseType::field_type field_type;
42 typedef typename DomainFunctionType::DiscreteFunctionSpaceType DomainFunctionSpaceType;
43 typedef typename RangeFunctionType::DiscreteFunctionSpaceType RangeFunctionSpaceType;
45 ISTLPreconditionAdapter (
const Preconditioner *precon,
const DomainFunctionSpaceType &domainSpace,
const RangeFunctionSpaceType &rangeSpace )
47 domainSpace_( domainSpace ),
48 rangeSpace_( rangeSpace )
52 virtual void pre ( domain_type &x, range_type &y ) {}
53 virtual void post ( domain_type &x ) {}
55 virtual void apply ( domain_type &x,
const range_type &y )
66 RangeFunctionType px(
"ISTLPreconditionAdapter::apply::x", rangeSpace_, x );
67 DomainFunctionType py(
"ISTLPreconditionAdapter::apply::y", domainSpace_, y );
74 const Preconditioner *precon_;
75 const DomainFunctionSpaceType &domainSpace_;
76 const RangeFunctionSpaceType &rangeSpace_;
80 template<
class BlockVector >
81 struct ISTLSolverReduction
83 ISTLSolverReduction (
double redEps,
double absLimit )
88 double operator() (
const Dune::LinearOperator< BlockVector, BlockVector > &op,
89 Dune::ScalarProduct< BlockVector > &scp,
90 const BlockVector &rhs,
const BlockVector &x )
const 94 BlockVector residuum( rhs );
95 op.applyscaleadd( -1., x, residuum );
96 const double res = scp.norm( residuum );
97 return (res > 0 ? absLimit_ / res : 1e-3);
110 template<
class Solver,
class Reduction = ISTLSolverReduction<
typename Solver::range_type > >
111 struct ISTLSolverAdapter
113 typedef Solver SolverType;
114 typedef Reduction ReductionType;
116 typedef typename SolverType::domain_type domain_type;
117 typedef typename SolverType::range_type range_type;
119 ISTLSolverAdapter (
const ReductionType &reduction,
unsigned int maxIterations,
int verbose,
121 : reduction_( reduction ),
122 maxIterations_( maxIterations ),
126 template<
class Op,
class ScP,
class PC >
127 void operator () ( Op& op, ScP &scp, PC &pc,
128 range_type &rhs, domain_type &x,
129 Dune::InverseOperatorResult &result )
const 132 SolverType solver( op, scp, pc, reduction_( op, scp, rhs, x ), maxIterations, verbose_ );
133 solver.apply( x, rhs, result );
137 ReductionType reduction_;
138 unsigned int maxIterations_;
143 template<
class X,
class Y,
class F,
class Reduction >
144 struct ISTLSolverAdapter<
Dune::RestartedGMResSolver< X, Y, F>, Reduction >
146 typedef Dune::RestartedGMResSolver< X, Y, F> SolverType;
147 typedef Reduction ReductionType;
149 typedef typename SolverType::domain_type domain_type;
150 typedef typename SolverType::range_type range_type;
152 ISTLSolverAdapter (
const ReductionType &reduction,
unsigned int restart,
unsigned int maxIterations,
int verbose,
154 : reduction_( reduction ),
156 maxIterations_( maxIterations ),
160 ISTLSolverAdapter (
const Reduction &reduction,
unsigned int maxIterations,
int verbose,
162 : reduction_( reduction ),
163 restart_( parameter.getValue< int >(
"fem.solver.gmres.restart", 20 ) ),
164 maxIterations_( maxIterations ),
168 template<
class Op,
class ScP,
class PC >
169 void operator () ( Op& op, ScP &scp, PC &pc,
170 range_type &rhs, domain_type &x,
171 Dune::InverseOperatorResult &result )
const 174 SolverType solver( op, scp, pc, reduction_( op, scp, rhs, x ), restart_, maxIterations, verbose_ );
175 solver.apply( x, rhs, result );
179 ReductionType reduction_;
180 unsigned int restart_;
181 unsigned int maxIterations_;
187 struct ISTLLoopSolver {
typedef LoopSolver< X > Type; };
190 struct ISTLGradientSolver {
typedef GradientSolver< X > Type; };
193 struct ISTLCGSolver {
typedef CGSolver< X > Type; };
196 struct ISTLBiCGSTABSolver {
typedef BiCGSTABSolver< X > Type; };
199 struct ISTLMINRESSolver {
typedef MINRESSolver< X > Type; };
202 struct ISTLRestartedGMRes {
typedef RestartedGMResSolver< X > Type; };
208 template<
class DiscreteFunction,
template<
class >
class Solver,
210 class ISTLInverseOperator
211 :
public Operator< DiscreteFunction, DiscreteFunction >
216 typedef typename BaseType::DomainFunctionType DomainFunctionType;
217 typedef typename BaseType::RangeFunctionType RangeFunctionType;
220 typedef Preconditioner PreconditionerType;
223 typedef ISTLLinearOperatorAdapter< OperatorType > ISTLOperatorType;
224 typedef ISTLPreconditionAdapter< OperatorType > ISTLPreconditionerAdapterType;
226 typedef Fem::ParallelScalarProduct< RangeFunctionType > ParallelScalarProductType;
227 typedef typename DomainFunctionType::DofStorageType BlockVectorType;
229 typedef ISTLSolverAdapter< typename Solver< BlockVectorType >::Type > SolverAdapterType;
230 typedef typename SolverAdapterType::ReductionType ReductionType;
233 typedef typename SolverAdapterType::SolverType SolverType;
235 ISTLInverseOperator (
const OperatorType &op,
236 double redEps,
double absLimit,
unsigned int maxIterations,
bool verbose,
238 : ISTLInverseOperator( op, nullptr, redEps, absLimit, maxIterations, verbose, parameter ) {}
240 ISTLInverseOperator (
const OperatorType &op,
241 double redEps,
double absLimit,
243 : ISTLInverseOperator( op, nullptr, redEps, absLimit,
std::numeric_limits< unsigned int >::
max(), false, parameter ) {}
245 ISTLInverseOperator (
const OperatorType &op,
246 double redEps,
double absLimit,
unsigned int maxIterations,
248 : ISTLInverseOperator( op, nullptr, redEps, absLimit, maxIterations, false, parameter ) {}
251 ISTLInverseOperator (
const OperatorType &op, PreconditionerType &preconditioner,
252 double redEps,
double absLimit,
unsigned int maxIterations,
bool verbose,
254 : ISTLInverseOperator( op, &preconditioner, redEps, absLimit, maxIterations, verbose, parameter ) {}
256 ISTLInverseOperator (
const OperatorType &op, PreconditionerType &preconditioner,
257 double redEps,
double absLimit,
259 : ISTLInverseOperator( op, &preconditioner, redEps, absLimit,
std::numeric_limits< unsigned int >::
max(), false, parameter ) {}
261 ISTLInverseOperator (
const OperatorType &op, PreconditionerType &preconditioner,
262 double redEps,
double absLimit,
unsigned int maxIterations,
264 : ISTLInverseOperator( op, &preconditioner, redEps, absLimit, maxIterations, false, parameter ) {}
267 virtual void operator() (
const DomainFunctionType &u, RangeFunctionType &w )
const 269 ISTLOperatorType istlOperator( operator_, w.space(), u.space() );
270 ParallelScalarProductType scp( u.space() );
272 if( !preconditioner_ )
274 ISTLPreconditionerAdapterType istlPreconditioner(
nullptr, w.space(), u.space() );
275 solve( istlOperator, scp, istlPreconditioner, u, w );
278 solve( istlOperator, scp, *preconditioner_, u, w );
281 unsigned int iterations ()
const {
return result_.iterations; }
284 ISTLInverseOperator (
const OperatorType &op, PreconditionerType *preconditioner,
285 double redEps,
double absLimit,
unsigned int maxIterations,
bool verbose,
288 preconditioner_( preconditioner ),
289 solverAdapter_( ReductionType( redEps, absLimit ), maxIterations, (Parameter::verbose() && verbose) ? 2 : 0, parameter )
292 void solve ( ISTLOperatorType &istlOperator, ParallelScalarProductType &scp,
293 const OperatorType &preconditioner,
294 const DomainFunctionType &u, RangeFunctionType &w )
const 296 ISTLPreconditionerAdapterType istlPreconditioner( &preconditioner, w.space(), u.space() );
297 solve( istlOperator, scp, istlPreconditioner, u, w );
300 template<
class ISTLPreconditioner >
301 void solve ( ISTLOperatorType &istlOperator, ParallelScalarProductType &scp,
302 ISTLPreconditioner &preconditioner,
303 const DomainFunctionType &u, RangeFunctionType &w )
const 305 BlockVectorType rhs( u.blockVector() );
306 solverAdapter_( istlOperator, scp, preconditioner, rhs, w.blockVector(), result_ );
309 const OperatorType &operator_;
310 PreconditionerType *preconditioner_;
311 SolverAdapterType solverAdapter_;
312 mutable Dune::InverseOperatorResult result_;
319 #endif // #if HAVE_ISTL 321 #endif // #ifndef DUNE_FEM_SOLVER_ISTLINVERSEOPERATORS_HH
static constexpr T min(T a)
Definition: utility.hh:81
BasicParameterReader< std::function< const std::string *(const std::string &, const std::string *) > > ParameterReader
Definition: reader.hh:264
static constexpr T max(T a)
Definition: utility.hh:65
static double max(const Double &v, const double p)
Definition: double.hh:387
Definition: coordinate.hh:4
static ParameterContainer & container()
Definition: io/parameter.hh:190