12#ifndef ELASTICITY_PRECONDITIONERS_HPP_
13#define ELASTICITY_PRECONDITIONERS_HPP_
15#include <opm/common/utility/platform_dependent/disable_warnings.h>
17#include <dune/common/fmatrix.hh>
18#include <dune/istl/bcrsmatrix.hh>
19#include <dune/istl/matrixmatrix.hh>
20#include <dune/istl/ilu.hh>
21#include <dune/istl/solvers.hh>
22#include <dune/istl/preconditioners.hh>
23#include <dune/istl/superlu.hh>
24#include <dune/istl/umfpack.hh>
25#include <dune/istl/paamg/amg.hh>
26#include <dune/istl/paamg/fastamg.hh>
27#include <dune/istl/paamg/twolevelmethod.hh>
28#include <dune/istl/overlappingschwarz.hh>
30#include <opm/common/utility/platform_dependent/reenable_warnings.h>
32#include <opm/grid/CpGrid.hpp>
40#if defined(HAVE_SUITESPARSE_UMFPACK)
41typedef Dune::UMFPack<Matrix> LUSolver;
42#elif defined(HAVE_SUPERLU)
43typedef Dune::SuperLU<Matrix> LUSolver;
45static_assert(
false,
"Enable either SuperLU or UMFPACK");
49typedef Dune::MatrixAdapter<Matrix,Vector,Vector>
Operator;
67 Dune::SymmetricMultiplicativeSchwarzMode,
78 static std::shared_ptr<type>
80 std::shared_ptr<Operator>& op,
const Dune::CpGrid& gv,
83 return std::shared_ptr<type>(
setup2(op, gv, A, copy));
91 static type*
setup2(std::shared_ptr<Operator>& op,
const Dune::CpGrid& gv,
96template<
class Smoother>
102 typedef Dune::Amg::SymmetricCriterion<Matrix, CouplingMetric>
CritBase;
105 typedef Dune::Amg::CoarsenCriterion<CritBase>
Criterion;
107 typedef Dune::Amg::AMG<Operator, Vector, Smoother> type;
117 static std::shared_ptr<type>
118 setup(
int pre,
int post,
int target,
int zcells,
119 std::shared_ptr<Operator>& op,
const Dune::CpGrid& ,
124 args.relaxationFactor = 1.0;
125 crit.setCoarsenTarget(target);
127 crit.setNoPreSmoothSteps(pre);
128 crit.setNoPostSmoothSteps(post);
129 crit.setDefaultValuesIsotropic(3, zcells);
131 std::cout <<
"\t collapsing 2x2x" << zcells <<
" cells per level" << std::endl;
133 return std::shared_ptr<type>(
new type(*op, crit, args));
139 typedef Dune::Amg::FastAMG<Operator, Vector> type;
150 static std::shared_ptr<type>
151 setup(
int pre,
int post,
int target,
int zcells,
152 std::shared_ptr<Operator>& op,
const Dune::CpGrid& gv,
158 template<
class Smoother>
161 typedef Dune::Amg::AggregationLevelTransferPolicy<
Operator,
164 typedef Dune::Amg::LevelTransferPolicy<Operator, Operator> LevelTransferPolicy;
169 typedef typename Dune::Amg::SmootherTraits<Smoother>::Arguments SmootherArgs;
171 typedef Dune::Amg::TwoLevelMethod<Operator, CoarsePolicy, Schwarz::type> type;
182 static std::shared_ptr<type>
183 setup(
int pre,
int post,
int target,
int zcells,
184 std::shared_ptr<Operator>& op,
const Dune::CpGrid& gv,
189 args.relaxationFactor = 1.0;
190 crit.setCoarsenTarget(target);
192 crit.setNoPreSmoothSteps(pre);
193 crit.setNoPostSmoothSteps(post);
194 crit.setDefaultValuesIsotropic(3, zcells);
195 CoarsePolicy coarsePolicy(args, crit);
199 return std::shared_ptr<type>(
new type(*op, fsp, policy, coarsePolicy, pre, post));
Class handling finite element assembly.
Class handling finite element assembly.
Definition asmhandler.hpp:35
Dune::MatrixAdapter< Matrix, Vector, Vector > Operator
A linear operator.
Definition elasticity_preconditioners.hpp:49
Dune::SeqJac< Matrix, Vector, Vector > JACSmoother
GJ AMG smoother.
Definition elasticity_preconditioners.hpp:55
Dune::SeqOverlappingSchwarz< Matrix, Vector, Dune::SymmetricMultiplicativeSchwarzMode, LUSolver > SchwarzSmoother
Schwarz + ILU0 AMG smoother.
Definition elasticity_preconditioners.hpp:62
Dune::SeqSSOR< Matrix, Vector, Vector > SSORSmoother
SSOR AMG smoother.
Definition elasticity_preconditioners.hpp:52
Dune::SeqILU< Matrix, Vector, Vector > ILUSmoother
ILU0 AMG smoother.
Definition elasticity_preconditioners.hpp:58
Smoother
Smoother used in the AMG.
Definition elasticity_upscale.hpp:75
Helper class with some matrix operations.
Dune::BCRSMatrix< Dune::FieldMatrix< double, 1, 1 > > Matrix
A sparse matrix holding our operator.
Definition matrixops.hpp:27
Dune::BlockVector< Dune::FieldVector< double, 1 > > Vector
A vector holding our RHS.
Definition matrixops.hpp:33
Inverting small matrices.
Definition ImplicitAssembly.hpp:43
An AMG.
Definition elasticity_preconditioners.hpp:97
Dune::Amg::CoarsenCriterion< CritBase > Criterion
The coarsening criterion used in the AMG.
Definition elasticity_preconditioners.hpp:105
Dune::Amg::FirstDiagonal CouplingMetric
The coupling metric used in the AMG.
Definition elasticity_preconditioners.hpp:99
static std::shared_ptr< type > setup(int pre, int post, int target, int zcells, std::shared_ptr< Operator > &op, const Dune::CpGrid &, ASMHandler< Dune::CpGrid > &, bool ©)
Setup preconditioner.
Definition elasticity_preconditioners.hpp:118
Dune::Amg::SymmetricCriterion< Matrix, CouplingMetric > CritBase
The coupling criterion used in the AMG.
Definition elasticity_preconditioners.hpp:102
A two-level method with a coarse AMG solver.
Definition elasticity_preconditioners.hpp:159
Dune::Amg::AggregationLevelTransferPolicy< Operator, typename AMG1< Smoother >::Criterion > TransferPolicy
AMG transfer policy.
Definition elasticity_preconditioners.hpp:162
static std::shared_ptr< type > setup(int pre, int post, int target, int zcells, std::shared_ptr< Operator > &op, const Dune::CpGrid &gv, ASMHandler< Dune::CpGrid > &A, bool ©)
Setup preconditioner.
Definition elasticity_preconditioners.hpp:183
A FastAMG.
Definition elasticity_preconditioners.hpp:138
static std::shared_ptr< type > setup(int pre, int post, int target, int zcells, std::shared_ptr< Operator > &op, const Dune::CpGrid &gv, ASMHandler< Dune::CpGrid > &A, bool ©)
Setup preconditioner.
Definition elasticity_preconditioners.cpp:20
Overlapping Schwarz preconditioner.
Definition elasticity_preconditioners.hpp:65
static std::shared_ptr< type > setup(int, int, int, int, std::shared_ptr< Operator > &op, const Dune::CpGrid &gv, ASMHandler< Dune::CpGrid > &A, bool ©)
Setup preconditioner.
Definition elasticity_preconditioners.hpp:79
static type * setup2(std::shared_ptr< Operator > &op, const Dune::CpGrid &gv, ASMHandler< Dune::CpGrid > &A, bool ©)
Setup preconditioner.
Definition elasticity_preconditioners.cpp:37