elasticity_preconditioners.hpp
Go to the documentation of this file.
1//==============================================================================
11//==============================================================================
12#ifndef ELASTICITY_PRECONDITIONERS_HPP_
13#define ELASTICITY_PRECONDITIONERS_HPP_
14
15#include <opm/common/utility/platform_dependent/disable_warnings.h>
16
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>
29
30#include <opm/common/utility/platform_dependent/reenable_warnings.h>
31
32#include <opm/grid/CpGrid.hpp>
35
36
37namespace Opm {
38namespace Elasticity {
39
40#if defined(HAVE_UMFPACK)
41typedef Dune::UMFPack<Matrix> LUSolver;
42#elif defined(HAVE_SUPERLU)
43typedef Dune::SuperLU<Matrix> LUSolver;
44#else
45static_assert(false, "Enable either SuperLU or UMFPACK");
46#endif
47
49typedef Dune::MatrixAdapter<Matrix,Vector,Vector> Operator;
50
52typedef Dune::SeqSSOR<Matrix, Vector, Vector> SSORSmoother;
53
55typedef Dune::SeqJac<Matrix, Vector, Vector> JACSmoother;
56
58typedef Dune::SeqILU<Matrix, Vector, Vector> ILUSmoother;
59
61typedef Dune::SeqOverlappingSchwarz<Matrix,Vector,
62 Dune::SymmetricMultiplicativeSchwarzMode, LUSolver> SchwarzSmoother;
63
65struct Schwarz {
66 typedef Dune::SeqOverlappingSchwarz<Matrix, Vector,
67 Dune::SymmetricMultiplicativeSchwarzMode,
68 LUSolver> type;
77 static std::shared_ptr<type>
78 setup(int /* pre */, int /* post */, int /* target */, int /* zcells */,
79 std::shared_ptr<Operator>& op, const Dune::CpGrid& gv,
80 ASMHandler<Dune::CpGrid>& A, bool& copy)
81 {
82 return std::shared_ptr<type>(setup2(op, gv, A, copy));
83 }
84
90 static type* setup2(std::shared_ptr<Operator>& op, const Dune::CpGrid& gv,
91 ASMHandler<Dune::CpGrid>& A, bool& copy);
92};
93
95template<class Smoother>
96struct AMG1 {
98 typedef Dune::Amg::FirstDiagonal CouplingMetric;
99
101 typedef Dune::Amg::SymmetricCriterion<Matrix, CouplingMetric> CritBase;
102
104 typedef Dune::Amg::CoarsenCriterion<CritBase> Criterion;
105
106 typedef Dune::Amg::AMG<Operator, Vector, Smoother> type;
107
116 static std::shared_ptr<type>
117 setup(int pre, int post, int target, int zcells,
118 std::shared_ptr<Operator>& op, const Dune::CpGrid& /* gv */,
119 ASMHandler<Dune::CpGrid>& /* A */, bool& copy)
120 {
121 Criterion crit;
123 args.relaxationFactor = 1.0;
124 crit.setCoarsenTarget(target);
125 crit.setGamma(1);
126 crit.setNoPreSmoothSteps(pre);
127 crit.setNoPostSmoothSteps(post);
128 crit.setDefaultValuesIsotropic(3, zcells);
129
130 std::cout << "\t collapsing 2x2x" << zcells << " cells per level" << std::endl;
131 copy = true;
132 return std::shared_ptr<type>(new type(*op, crit, args));
133 }
134};
135
137struct FastAMG {
138 typedef Dune::Amg::FastAMG<Operator, Vector> type;
139
148 static std::shared_ptr<type>
149 setup(int pre, int post, int target, int zcells,
150 std::shared_ptr<Operator>& op, const Dune::CpGrid& gv,
151 ASMHandler<Dune::CpGrid>& A, bool& copy);
152};
153
154
156 template<class Smoother>
157struct AMG2Level {
159 typedef Dune::Amg::AggregationLevelTransferPolicy<Operator,
161
162 typedef Dune::Amg::LevelTransferPolicy<Operator, Operator> LevelTransferPolicy;
163
164 typedef Dune::Amg::OneStepAMGCoarseSolverPolicy<Operator, Smoother,
166
167 typedef typename Dune::Amg::SmootherTraits<Smoother>::Arguments SmootherArgs;
168
169 typedef Dune::Amg::TwoLevelMethod<Operator, CoarsePolicy, Schwarz::type> type;
170
179 static std::shared_ptr<type>
180 setup(int pre, int post, int target, int zcells,
181 std::shared_ptr<Operator>& op, const Dune::CpGrid& gv,
182 ASMHandler<Dune::CpGrid>& A, bool& copy)
183 {
184 typename AMG1<Smoother>::Criterion crit;
185 SmootherArgs args;
186 args.relaxationFactor = 1.0;
187 crit.setCoarsenTarget(target);
188 crit.setGamma(1);
189 crit.setNoPreSmoothSteps(pre);
190 crit.setNoPostSmoothSteps(post);
191 crit.setDefaultValuesIsotropic(3, zcells);
192 CoarsePolicy coarsePolicy(args, crit);
193 TransferPolicy policy(crit);
194 std::shared_ptr<Schwarz::type> fsp(Schwarz::setup2(op, gv, A, copy));
195 copy = true;
196 return std::shared_ptr<type>(new type(*op, fsp, policy, coarsePolicy, pre, post));
197 }
198};
199
200}
201}
202
203#endif
Class handling finite element assembly.
Class handling finite element assembly.
Definition: asmhandler.hpp:35
Helper class with some matrix operations.
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::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
Dune::SeqOverlappingSchwarz< Matrix, Vector, Dune::SymmetricMultiplicativeSchwarzMode, LUSolver > SchwarzSmoother
Schwarz + ILU0 AMG smoother.
Definition: elasticity_preconditioners.hpp:62
Smoother
Smoother used in the AMG.
Definition: elasticity_upscale.hpp:75
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
Definition: ImplicitAssembly.hpp:43
An AMG.
Definition: elasticity_preconditioners.hpp:96
Dune::Amg::CoarsenCriterion< CritBase > Criterion
The coarsening criterion used in the AMG.
Definition: elasticity_preconditioners.hpp:104
Dune::Amg::FirstDiagonal CouplingMetric
The coupling metric used in the AMG.
Definition: elasticity_preconditioners.hpp:98
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 &copy)
Setup preconditioner.
Definition: elasticity_preconditioners.hpp:117
Dune::Amg::AMG< Operator, Vector, Smoother > type
Definition: elasticity_preconditioners.hpp:106
Dune::Amg::SymmetricCriterion< Matrix, CouplingMetric > CritBase
The coupling criterion used in the AMG.
Definition: elasticity_preconditioners.hpp:101
A two-level method with a coarse AMG solver.
Definition: elasticity_preconditioners.hpp:157
Dune::Amg::SmootherTraits< Smoother >::Arguments SmootherArgs
Definition: elasticity_preconditioners.hpp:167
Dune::Amg::OneStepAMGCoarseSolverPolicy< Operator, Smoother, typename AMG1< Smoother >::Criterion > CoarsePolicy
Definition: elasticity_preconditioners.hpp:165
Dune::Amg::LevelTransferPolicy< Operator, Operator > LevelTransferPolicy
Definition: elasticity_preconditioners.hpp:162
Dune::Amg::AggregationLevelTransferPolicy< Operator, typename AMG1< Smoother >::Criterion > TransferPolicy
AMG transfer policy.
Definition: elasticity_preconditioners.hpp:160
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 &copy)
Setup preconditioner.
Definition: elasticity_preconditioners.hpp:180
Dune::Amg::TwoLevelMethod< Operator, CoarsePolicy, Schwarz::type > type
Definition: elasticity_preconditioners.hpp:169
A FastAMG.
Definition: elasticity_preconditioners.hpp:137
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 &copy)
Setup preconditioner.
Dune::Amg::FastAMG< Operator, Vector > type
Definition: elasticity_preconditioners.hpp:138
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 &copy)
Setup preconditioner.
Definition: elasticity_preconditioners.hpp:78
static type * setup2(std::shared_ptr< Operator > &op, const Dune::CpGrid &gv, ASMHandler< Dune::CpGrid > &A, bool &copy)
Setup preconditioner.
Dune::SeqOverlappingSchwarz< Matrix, Vector, Dune::SymmetricMultiplicativeSchwarzMode, LUSolver > type
Definition: elasticity_preconditioners.hpp:68