Line | Branch | Exec | Source |
---|---|---|---|
1 | // -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*- | ||
2 | // vi: set et ts=4 sw=4 sts=4: | ||
3 | // | ||
4 | // SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder | ||
5 | // SPDX-License-Identifier: GPL-3.0-or-later | ||
6 | // | ||
7 | /*! | ||
8 | * \file | ||
9 | * \ingroup Linear | ||
10 | * \brief Dumux sequential linear solver backends | ||
11 | */ | ||
12 | #ifndef DUMUX_SEQ_SOLVER_BACKEND_HH | ||
13 | #define DUMUX_SEQ_SOLVER_BACKEND_HH | ||
14 | |||
15 | #include <type_traits> | ||
16 | #include <tuple> | ||
17 | #include <utility> | ||
18 | |||
19 | #include <dune/istl/preconditioners.hh> | ||
20 | #include <dune/istl/solvers.hh> | ||
21 | #include <dune/istl/io.hh> | ||
22 | #include <dune/common/indices.hh> | ||
23 | #include <dune/common/hybridutilities.hh> | ||
24 | |||
25 | #include <dumux/common/parameters.hh> | ||
26 | #include <dumux/common/typetraits/matrix.hh> | ||
27 | #include <dumux/common/typetraits/utility.hh> | ||
28 | #include <dumux/linear/solver.hh> | ||
29 | #include <dumux/linear/preconditioners.hh> | ||
30 | #include <dumux/linear/linearsolverparameters.hh> | ||
31 | #include <dumux/linear/parallelmatrixadapter.hh> | ||
32 | |||
33 | namespace Dumux { | ||
34 | |||
35 | /*! | ||
36 | * \ingroup Linear | ||
37 | * \brief A general solver backend allowing arbitrary preconditioners and solvers. | ||
38 | * | ||
39 | * This class is used as a base class for specific solver-preconditioner | ||
40 | * combinations. Several parameters from the group LinearSolver are read to | ||
41 | * customize the solver and preconditioner: | ||
42 | * | ||
43 | * - Verbosity: determines how verbose the linear solver should print output. | ||
44 | * - MaxIterations: the maximum number of iterations for the linear solver. | ||
45 | * - ResidualReduction: the threshold for declaration of convergence. | ||
46 | * - PreconditionerRelaxation: relaxation parameter for the preconditioner. | ||
47 | * - PreconditionerIterations: usually specifies the number of times the | ||
48 | * preconditioner is applied. In case of ILU(n), | ||
49 | * it specifies the order of the applied ILU. | ||
50 | */ | ||
51 | class IterativePreconditionedSolverImpl | ||
52 | { | ||
53 | public: | ||
54 | |||
55 | template<class Preconditioner, class Solver, class SolverInterface, class Matrix, class Vector> | ||
56 | [[deprecated("Removed after 3.8. Use solver from istlsolvers.hh")]] | ||
57 | static bool solve(const SolverInterface& s, const Matrix& A, Vector& x, const Vector& b, | ||
58 | const std::string& modelParamGroup = "") | ||
59 | { | ||
60 | Preconditioner precond(A, s.precondIter(), s.relaxation()); | ||
61 | |||
62 | // make a linear operator from a matrix | ||
63 | using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>; | ||
64 | MatrixAdapter linearOperator(A); | ||
65 | |||
66 | Solver solver(linearOperator, precond, s.residReduction(), s.maxIter(), s.verbosity()); | ||
67 | |||
68 | Vector bTmp(b); | ||
69 | |||
70 | Dune::InverseOperatorResult result; | ||
71 | solver.apply(x, bTmp, result); | ||
72 | |||
73 | return result.converged; | ||
74 | } | ||
75 | |||
76 | template<class Preconditioner, class Solver, class SolverInterface, class Matrix, class Vector> | ||
77 | [[deprecated("Removed after 3.8. Use solver from istlsolvers.hh")]] | ||
78 | static bool solveWithGMRes(const SolverInterface& s, const Matrix& A, Vector& x, const Vector& b, | ||
79 | const std::string& modelParamGroup = "") | ||
80 | { | ||
81 | // get the restart threshold | ||
82 | const int restartGMRes = getParamFromGroup<int>(modelParamGroup, "LinearSolver.GMResRestart", 10); | ||
83 | |||
84 | Preconditioner precond(A, s.precondIter(), s.relaxation()); | ||
85 | |||
86 | // make a linear operator from a matrix | ||
87 | using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>; | ||
88 | MatrixAdapter linearOperator(A); | ||
89 | |||
90 | Solver solver(linearOperator, precond, s.residReduction(), restartGMRes, s.maxIter(), s.verbosity()); | ||
91 | |||
92 | Vector bTmp(b); | ||
93 | |||
94 | Dune::InverseOperatorResult result; | ||
95 | solver.apply(x, bTmp, result); | ||
96 | |||
97 | return result.converged; | ||
98 | } | ||
99 | |||
100 | template<class Preconditioner, class Solver, class SolverInterface, class Matrix, class Vector> | ||
101 | [[deprecated("Removed after 3.8. Use solver from istlsolvers.hh")]] | ||
102 | static bool solveWithILU0Prec(const SolverInterface& s, const Matrix& A, Vector& x, const Vector& b, | ||
103 | const std::string& modelParamGroup = "") | ||
104 | { | ||
105 | Preconditioner precond(A, s.relaxation()); | ||
106 | |||
107 | using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>; | ||
108 | MatrixAdapter operatorA(A); | ||
109 | |||
110 | Solver solver(operatorA, precond, s.residReduction(), s.maxIter(), s.verbosity()); | ||
111 | |||
112 | Vector bTmp(b); | ||
113 | |||
114 | Dune::InverseOperatorResult result; | ||
115 | solver.apply(x, bTmp, result); | ||
116 | |||
117 | return result.converged; | ||
118 | } | ||
119 | |||
120 | // solve with RestartedGMRes (needs restartGMRes as additional argument) | ||
121 | template<class Preconditioner, class Solver, class SolverInterface, class Matrix, class Vector> | ||
122 | [[deprecated("Removed after 3.8. Use solver from istlsolvers.hh")]] | ||
123 | static bool solveWithILU0PrecGMRes(const SolverInterface& s, const Matrix& A, Vector& x, const Vector& b, | ||
124 | const std::string& modelParamGroup = "") | ||
125 | { | ||
126 | // get the restart threshold | ||
127 | const int restartGMRes = getParamFromGroup<int>(modelParamGroup, "LinearSolver.GMResRestart", 10); | ||
128 | |||
129 | Preconditioner precond(A, s.relaxation()); | ||
130 | |||
131 | using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>; | ||
132 | MatrixAdapter operatorA(A); | ||
133 | |||
134 | Solver solver(operatorA, precond, s.residReduction(), restartGMRes, s.maxIter(), s.verbosity()); | ||
135 | |||
136 | Vector bTmp(b); | ||
137 | |||
138 | Dune::InverseOperatorResult result; | ||
139 | solver.apply(x, bTmp, result); | ||
140 | |||
141 | return result.converged; | ||
142 | } | ||
143 | |||
144 | // solve with generic parameter tree | ||
145 | template<class Preconditioner, class Solver, class Matrix, class Vector> | ||
146 | static bool solveWithParamTree(const Matrix& A, Vector& x, const Vector& b, | ||
147 | const Dune::ParameterTree& params) | ||
148 | { | ||
149 | // make a linear operator from a matrix | ||
150 | using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>; | ||
151 | const auto linearOperator = std::make_shared<MatrixAdapter>(A); | ||
152 | |||
153 | auto precond = std::make_shared<Preconditioner>(linearOperator, params.sub("preconditioner")); | ||
154 | Solver solver(linearOperator, precond, params); | ||
155 | |||
156 | Vector bTmp(b); | ||
157 | |||
158 | Dune::InverseOperatorResult result; | ||
159 | solver.apply(x, bTmp, result); | ||
160 | |||
161 | return result.converged; | ||
162 | } | ||
163 | }; | ||
164 | |||
165 | /*! | ||
166 | * \ingroup Linear | ||
167 | * \brief Returns the block level for the preconditioner for a given matrix | ||
168 | * | ||
169 | * \tparam M The matrix. | ||
170 | */ | ||
171 | template<class M> | ||
172 | constexpr std::size_t preconditionerBlockLevel() noexcept | ||
173 | { | ||
174 | return isMultiTypeBlockMatrix<M>::value ? 2 : 1; | ||
175 | } | ||
176 | |||
177 | /*! | ||
178 | * \ingroup Linear | ||
179 | * \brief Solver for simple block-diagonal matrices (e.g. from explicit time stepping schemes) | ||
180 | * | ||
181 | * Solver: Single Jacobi iteration | ||
182 | * Preconditioner: Unity | ||
183 | */ | ||
184 |
1/6✗ Branch 0 not taken.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
2 | class ExplicitDiagonalSolver : public LinearSolver |
185 | { | ||
186 | public: | ||
187 | using LinearSolver::LinearSolver; | ||
188 | |||
189 | template<class Matrix, class Vector> | ||
190 | 2000 | bool solve(const Matrix& A, Vector& x, const Vector& b) | |
191 | { | ||
192 |
0/2✗ Branch 1 not taken.
✗ Branch 2 not taken.
|
2000 | Vector rhs(b); |
193 | 2000 | constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>(); | |
194 |
1/2✓ Branch 0 taken 2000 times.
✗ Branch 1 not taken.
|
4000 | Dune::SeqJac<Matrix, Vector, Vector, precondBlockLevel> jac(A, 1, 1.0); |
195 | jac.pre(x, rhs); | ||
196 | 2000 | jac.apply(x, rhs); | |
197 | 2000 | jac.post(x); | |
198 |
1/2✓ Branch 0 taken 2000 times.
✗ Branch 1 not taken.
|
4000 | return true; |
199 | } | ||
200 | |||
201 | std::string name() const | ||
202 | { | ||
203 | return "Explicit diagonal matrix solver"; | ||
204 | } | ||
205 | }; | ||
206 | |||
207 | /*! | ||
208 | * \name Solver for MultiTypeBlockMatrix's | ||
209 | */ | ||
210 | // \{ | ||
211 | |||
212 | /*! | ||
213 | * \ingroup Linear | ||
214 | * \brief A Uzawa preconditioned BiCGSTAB solver for saddle-point problems | ||
215 | */ | ||
216 | template <class LinearSolverTraits> | ||
217 | class UzawaBiCGSTABBackend : public LinearSolver | ||
218 | { | ||
219 | public: | ||
220 | using LinearSolver::LinearSolver; | ||
221 | |||
222 | template<class Matrix, class Vector> | ||
223 | bool solve(const Matrix& A, Vector& x, const Vector& b) | ||
224 | { | ||
225 | using Preconditioner = SeqUzawa<Matrix, Vector, Vector>; | ||
226 | using Solver = Dune::BiCGSTABSolver<Vector>; | ||
227 | static const auto solverParams = LinearSolverParameters<LinearSolverTraits>::createParameterTree(this->paramGroup()); | ||
228 | return IterativePreconditionedSolverImpl::template solveWithParamTree<Preconditioner, Solver>(A, x, b, solverParams); | ||
229 | } | ||
230 | |||
231 | std::string name() const | ||
232 | { | ||
233 | return "Uzawa preconditioned BiCGSTAB solver"; | ||
234 | } | ||
235 | }; | ||
236 | |||
237 | /*! | ||
238 | * \ingroup Linear | ||
239 | * \brief A simple ilu0 block diagonal preconditioner | ||
240 | */ | ||
241 | template<class M, class X, class Y, int blockLevel = 2> | ||
242 | 2608 | class BlockDiagILU0Preconditioner : public Dune::Preconditioner<X, Y> | |
243 | { | ||
244 | template<std::size_t i> | ||
245 | using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>; | ||
246 | |||
247 | template<std::size_t i> | ||
248 | using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>; | ||
249 | |||
250 | template<std::size_t i> | ||
251 | using BlockILU = Dune::SeqILU<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>, blockLevel-1>; | ||
252 | |||
253 | using ILUTuple = typename makeFromIndexedType<std::tuple, BlockILU, std::make_index_sequence<M::N()> >::type; | ||
254 | |||
255 | public: | ||
256 | //! \brief The matrix type the preconditioner is for. | ||
257 | using matrix_type = typename std::decay_t<M>; | ||
258 | //! \brief The domain type of the preconditioner. | ||
259 | using domain_type = X; | ||
260 | //! \brief The range type of the preconditioner. | ||
261 | using range_type = Y; | ||
262 | //! \brief The field type of the preconditioner. | ||
263 | using field_type = typename X::field_type; | ||
264 | |||
265 | /*! \brief Constructor. | ||
266 | |||
267 | Constructor gets all parameters to operate the prec. | ||
268 | \param m The (multi type block) matrix to operate on | ||
269 | \param w The relaxation factor | ||
270 | */ | ||
271 | 1304 | BlockDiagILU0Preconditioner(const M& m, double w = 1.0) | |
272 | 1304 | : BlockDiagILU0Preconditioner(m, w, std::make_index_sequence<M::N()>{}) | |
273 | { | ||
274 | static_assert(blockLevel >= 2, "Only makes sense for MultiTypeBlockMatrix!"); | ||
275 | } | ||
276 | |||
277 | 1304 | void pre (X& v, Y& d) final {} | |
278 | |||
279 | 153487 | void apply (X& v, const Y& d) final | |
280 | { | ||
281 | using namespace Dune::Hybrid; | ||
282 | 153487 | forEach(integralRange(Dune::Hybrid::size(ilu_)), [&](const auto i) | |
283 | { | ||
284 | ✗ | std::get<decltype(i)::value>(ilu_).apply(v[i], d[i]); | |
285 | }); | ||
286 | 153487 | } | |
287 | |||
288 | 1304 | void post (X&) final {} | |
289 | |||
290 | //! Category of the preconditioner (see SolverCategory::Category) | ||
291 | 1304 | Dune::SolverCategory::Category category() const final | |
292 | { | ||
293 | 1304 | return Dune::SolverCategory::sequential; | |
294 | } | ||
295 | |||
296 | private: | ||
297 | template<std::size_t... Is> | ||
298 | 1304 | BlockDiagILU0Preconditioner (const M& m, double w, std::index_sequence<Is...> is) | |
299 |
3/6✓ Branch 4 taken 1304 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 1304 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 1304 times.
✗ Branch 11 not taken.
|
5216 | : ilu_(std::make_tuple(BlockILU<Is>(m[Dune::index_constant<Is>{}][Dune::index_constant<Is>{}], w)...)) |
300 | 1304 | {} | |
301 | |||
302 | ILUTuple ilu_; | ||
303 | }; | ||
304 | |||
305 | |||
306 | /*! | ||
307 | * \ingroup Linear | ||
308 | * \brief A simple ilu0 block diagonal preconditioned BiCGSTABSolver | ||
309 | * \note expects a system as a multi-type block-matrix | ||
310 | * | A B | | ||
311 | * | C D | | ||
312 | */ | ||
313 |
1/6✗ Branch 0 not taken.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 27 times.
|
27 | class BlockDiagILU0BiCGSTABSolver : public LinearSolver |
314 | { | ||
315 | |||
316 | public: | ||
317 | using LinearSolver::LinearSolver; | ||
318 | |||
319 | template<class Matrix, class Vector> | ||
320 | 1304 | bool solve(const Matrix& M, Vector& x, const Vector& b) | |
321 | { | ||
322 | 1304 | BlockDiagILU0Preconditioner<Matrix, Vector, Vector> preconditioner(M); | |
323 |
1/4✓ Branch 1 taken 1304 times.
✗ Branch 2 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
2608 | Dumux::ParallelMultiTypeMatrixAdapter<Matrix, Vector, Vector> op(M); |
324 |
5/10✓ Branch 1 taken 1304 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1304 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 1304 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 1304 times.
✗ Branch 11 not taken.
✓ Branch 12 taken 1304 times.
✗ Branch 13 not taken.
|
6520 | Dune::BiCGSTABSolver<Vector> solver(op, preconditioner, this->residReduction(), |
325 | this->maxIter(), this->verbosity()); | ||
326 |
1/2✓ Branch 1 taken 1304 times.
✗ Branch 2 not taken.
|
2608 | auto bTmp(b); |
327 |
1/2✓ Branch 1 taken 1304 times.
✗ Branch 2 not taken.
|
1304 | solver.apply(x, bTmp, result_); |
328 | |||
329 | 2608 | return result_.converged; | |
330 | } | ||
331 | |||
332 | const Dune::InverseOperatorResult& result() const | ||
333 | { | ||
334 | return result_; | ||
335 | } | ||
336 | |||
337 | ✗ | std::string name() const | |
338 |
6/16✓ Branch 1 taken 19 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 1 times.
✓ Branch 4 taken 19 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 1 times.
✓ Branch 7 taken 19 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✓ Branch 10 taken 19 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 16 not taken.
✗ Branch 17 not taken.
|
40 | { return "block-diagonal ILU0-preconditioned BiCGSTAB solver"; } |
339 | |||
340 | private: | ||
341 | Dune::InverseOperatorResult result_; | ||
342 | }; | ||
343 | |||
344 | /*! | ||
345 | * \ingroup Linear | ||
346 | * \brief A simple ilu0 block diagonal preconditioned RestartedGMResSolver | ||
347 | * \note expects a system as a multi-type block-matrix | ||
348 | * | A B | | ||
349 | * | C D | | ||
350 | */ | ||
351 | class BlockDiagILU0RestartedGMResSolver : public LinearSolver | ||
352 | { | ||
353 | |||
354 | public: | ||
355 | using LinearSolver::LinearSolver; | ||
356 | |||
357 | template<int precondBlockLevel = 2, class Matrix, class Vector> | ||
358 | bool solve(const Matrix& M, Vector& x, const Vector& b) | ||
359 | { | ||
360 | BlockDiagILU0Preconditioner<Matrix, Vector, Vector> preconditioner(M); | ||
361 | Dumux::ParallelMultiTypeMatrixAdapter<Matrix, Vector, Vector> op(M); | ||
362 | static const int restartGMRes = getParamFromGroup<int>(this->paramGroup(), "LinearSolver.GMResRestart"); | ||
363 | Dune::RestartedGMResSolver<Vector> solver(op, preconditioner, this->residReduction(), restartGMRes, | ||
364 | this->maxIter(), this->verbosity()); | ||
365 | auto bTmp(b); | ||
366 | solver.apply(x, bTmp, result_); | ||
367 | |||
368 | return result_.converged; | ||
369 | } | ||
370 | |||
371 | const Dune::InverseOperatorResult& result() const | ||
372 | { | ||
373 | return result_; | ||
374 | } | ||
375 | |||
376 | std::string name() const | ||
377 | { return "block-diagonal ILU0-preconditioned restarted GMRes solver"; } | ||
378 | |||
379 | private: | ||
380 | Dune::InverseOperatorResult result_; | ||
381 | }; | ||
382 | |||
383 | /*! | ||
384 | * \ingroup Linear | ||
385 | * \brief A simple ilu0 block diagonal preconditioner | ||
386 | */ | ||
387 | template<class M, class X, class Y, int blockLevel = 2> | ||
388 | 3080 | class BlockDiagAMGPreconditioner : public Dune::Preconditioner<X, Y> | |
389 | { | ||
390 | template<std::size_t i> | ||
391 | using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>; | ||
392 | |||
393 | template<std::size_t i> | ||
394 | using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>; | ||
395 | |||
396 | template<std::size_t i> | ||
397 | using Smoother = Dune::SeqSSOR<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>>; | ||
398 | |||
399 | template<std::size_t i> | ||
400 | using LinearOperator = Dune::MatrixAdapter<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>>; | ||
401 | |||
402 | template<std::size_t i> | ||
403 | using ScalarProduct = Dune::SeqScalarProduct<VecBlockType<i>>; | ||
404 | |||
405 | template<std::size_t i> | ||
406 | using BlockAMG = Dune::Amg::AMG<LinearOperator<i>, VecBlockType<i>, Smoother<i>>; | ||
407 | |||
408 | using AMGTuple = typename makeFromIndexedType<std::tuple, BlockAMG, std::make_index_sequence<M::N()> >::type; | ||
409 | |||
410 | public: | ||
411 | //! \brief The matrix type the preconditioner is for. | ||
412 | using matrix_type = typename std::decay_t<M>; | ||
413 | //! \brief The domain type of the preconditioner. | ||
414 | using domain_type = X; | ||
415 | //! \brief The range type of the preconditioner. | ||
416 | using range_type = Y; | ||
417 | //! \brief The field type of the preconditioner. | ||
418 | using field_type = typename X::field_type; | ||
419 | |||
420 | /*! \brief Constructor. | ||
421 | |||
422 | Constructor gets all parameters to operate the prec. | ||
423 | \param lop The linear operator | ||
424 | \param c The criterion | ||
425 | \param sa The smoother arguments | ||
426 | */ | ||
427 | template<class LOP, class Criterion, class SmootherArgs> | ||
428 | 1540 | BlockDiagAMGPreconditioner(const LOP& lop, const Criterion& c, const SmootherArgs& sa) | |
429 | 1540 | : BlockDiagAMGPreconditioner(lop, c, sa, std::make_index_sequence<M::N()>{}) | |
430 | { | ||
431 | static_assert(blockLevel >= 2, "Only makes sense for MultiTypeBlockMatrix!"); | ||
432 | } | ||
433 | |||
434 | /*! | ||
435 | * \brief Prepare the preconditioner. | ||
436 | * | ||
437 | * A solver solves a linear operator equation A(v)=d by applying | ||
438 | * one or several steps of the preconditioner. The method pre() | ||
439 | * is called before the first apply operation. | ||
440 | * d and v are right hand side and solution vector of the linear | ||
441 | * system respectively. It may. e.g., scale the system, allocate memory or | ||
442 | * compute a (I)LU decomposition. | ||
443 | * Note: The ILU decomposition could also be computed in the constructor | ||
444 | * or with a separate method of the derived method if several | ||
445 | * linear systems with the same matrix are to be solved. | ||
446 | * | ||
447 | * \param v The left hand side of the equation. | ||
448 | * \param d The right hand side of the equation. | ||
449 | */ | ||
450 | 1540 | void pre (X& v, Y& d) final | |
451 | { | ||
452 | using namespace Dune::Hybrid; | ||
453 | 1540 | forEach(integralRange(Dune::Hybrid::size(amg_)), [&](const auto i) | |
454 | { | ||
455 | ✗ | std::get<decltype(i)::value>(amg_).pre(v[i], d[i]); | |
456 | }); | ||
457 | 1540 | } | |
458 | |||
459 | /*! | ||
460 | * \brief Apply one step of the preconditioner to the system A(v)=d. | ||
461 | * | ||
462 | * On entry v=0 and d=b-A(x) (although this might not be | ||
463 | * computed in that way. On exit v contains the update, i.e | ||
464 | * one step computes \f$ v = M^{-1} d \f$ where \f$ M \f$ is the | ||
465 | * approximate inverse of the operator \f$ A \f$ characterizing | ||
466 | * the preconditioner. | ||
467 | * \param v The update to be computed | ||
468 | * \param d The current defect. | ||
469 | */ | ||
470 | 2509 | void apply (X& v, const Y& d) final | |
471 | { | ||
472 | using namespace Dune::Hybrid; | ||
473 | 2509 | forEach(integralRange(Dune::Hybrid::size(amg_)), [&](const auto i) | |
474 | { | ||
475 | ✗ | std::get<decltype(i)::value>(amg_).apply(v[i], d[i]); | |
476 | }); | ||
477 | 2509 | } | |
478 | |||
479 | /*! | ||
480 | * \brief Clean up. | ||
481 | * | ||
482 | * This method is called after the last apply call for the | ||
483 | * linear system to be solved. Memory may be deallocated safely | ||
484 | * here. v is the solution of the linear equation. | ||
485 | * | ||
486 | * \param v The right hand side of the equation. | ||
487 | */ | ||
488 | 1540 | void post (X& v) final | |
489 | { | ||
490 | using namespace Dune::Hybrid; | ||
491 |
1/2✓ Branch 1 taken 1540 times.
✗ Branch 2 not taken.
|
1540 | forEach(integralRange(Dune::Hybrid::size(amg_)), [&](const auto i) |
492 | { | ||
493 |
6/36✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 16 not taken.
✗ Branch 17 not taken.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✗ Branch 28 not taken.
✗ Branch 29 not taken.
✗ Branch 31 not taken.
✗ Branch 32 not taken.
✗ Branch 34 not taken.
✗ Branch 35 not taken.
✓ Branch 49 taken 1540 times.
✗ Branch 50 not taken.
✓ Branch 52 taken 1540 times.
✗ Branch 53 not taken.
✓ Branch 55 taken 1540 times.
✗ Branch 56 not taken.
✓ Branch 58 taken 1540 times.
✗ Branch 59 not taken.
✓ Branch 61 taken 1540 times.
✗ Branch 62 not taken.
✓ Branch 64 taken 1540 times.
✗ Branch 65 not taken.
|
4620 | std::get<decltype(i)::value>(amg_).post(v[i]); |
494 | }); | ||
495 | 1540 | } | |
496 | |||
497 | //! Category of the preconditioner (see SolverCategory::Category) | ||
498 | 1540 | Dune::SolverCategory::Category category() const final | |
499 | { | ||
500 | 1540 | return Dune::SolverCategory::sequential; | |
501 | } | ||
502 | |||
503 | private: | ||
504 | template<class LOP, class Criterion, class SmootherArgs, std::size_t... Is> | ||
505 | 1540 | BlockDiagAMGPreconditioner (const LOP& lop, const Criterion& c, const SmootherArgs& sa, std::index_sequence<Is...> is) | |
506 |
8/16✓ Branch 9 taken 1540 times.
✗ Branch 10 not taken.
✓ Branch 12 taken 1540 times.
✗ Branch 13 not taken.
✓ Branch 15 taken 1540 times.
✗ Branch 16 not taken.
✓ Branch 18 taken 1540 times.
✗ Branch 19 not taken.
✓ Branch 21 taken 1540 times.
✗ Branch 22 not taken.
✓ Branch 24 taken 1540 times.
✗ Branch 25 not taken.
✓ Branch 27 taken 1540 times.
✗ Branch 28 not taken.
✓ Branch 30 taken 1540 times.
✗ Branch 31 not taken.
|
13860 | : amg_(std::make_tuple(BlockAMG<Is>(*std::get<Is>(lop), *std::get<Is>(c), *std::get<Is>(sa))...)) |
507 | 1540 | {} | |
508 | |||
509 | AMGTuple amg_; | ||
510 | }; | ||
511 | |||
512 | /*! | ||
513 | * \ingroup Linear | ||
514 | * \brief A simple ilu0 block diagonal preconditioned BiCGSTABSolver | ||
515 | * \note expects a system as a multi-type block-matrix | ||
516 | * | A B | | ||
517 | * | C D | | ||
518 | */ | ||
519 |
1/6✗ Branch 0 not taken.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 23 times.
|
23 | class BlockDiagAMGBiCGSTABSolver : public LinearSolver |
520 | { | ||
521 | template<class M, std::size_t i> | ||
522 | using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>; | ||
523 | |||
524 | template<class X, std::size_t i> | ||
525 | using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>; | ||
526 | |||
527 | template<class M, class X, std::size_t i> | ||
528 | using Smoother = Dune::SeqSSOR<DiagBlockType<M, i>, VecBlockType<X, i>, VecBlockType<X, i>>; | ||
529 | |||
530 | template<class M, class X, std::size_t i> | ||
531 | using SmootherArgs = typename Dune::Amg::SmootherTraits<Smoother<M, X, i>>::Arguments; | ||
532 | |||
533 | template<class M, std::size_t i> | ||
534 | using Criterion = Dune::Amg::CoarsenCriterion<Dune::Amg::SymmetricCriterion<DiagBlockType<M, i>, Dune::Amg::FirstDiagonal>>; | ||
535 | |||
536 | template<class M, class X, std::size_t i> | ||
537 | using LinearOperator = Dune::MatrixAdapter<DiagBlockType<M, i>, VecBlockType<X, i>, VecBlockType<X, i>>; | ||
538 | |||
539 | public: | ||
540 | using LinearSolver::LinearSolver; | ||
541 | |||
542 | // Solve saddle-point problem using a Schur complement based preconditioner | ||
543 | template<class Matrix, class Vector> | ||
544 | 1540 | bool solve(const Matrix& m, Vector& x, const Vector& b) | |
545 | { | ||
546 | //! \todo Check whether the default accumulation mode atOnceAccu is needed. | ||
547 | //! \todo make parameters changeable at runtime from input file / parameter tree | ||
548 | 1540 | Dune::Amg::Parameters params(15, 2000, 1.2, 1.6, Dune::Amg::atOnceAccu); | |
549 | 1540 | params.setDefaultValuesIsotropic(3); | |
550 | 3080 | params.setDebugLevel(this->verbosity()); | |
551 | |||
552 | 1540 | auto criterion = makeCriterion_<Criterion, Matrix>(params, std::make_index_sequence<Matrix::N()>{}); | |
553 |
1/2✓ Branch 1 taken 1540 times.
✗ Branch 2 not taken.
|
3080 | auto smootherArgs = makeSmootherArgs_<SmootherArgs, Matrix, Vector>(std::make_index_sequence<Matrix::N()>{}); |
554 | |||
555 | using namespace Dune::Hybrid; | ||
556 |
1/2✓ Branch 1 taken 1540 times.
✗ Branch 2 not taken.
|
1540 | forEach(std::make_index_sequence<Matrix::N()>{}, [&](const auto i) |
557 | { | ||
558 |
2/4✓ Branch 1 taken 1540 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1540 times.
✗ Branch 5 not taken.
|
3080 | auto& args = std::get<decltype(i)::value>(smootherArgs); |
559 |
2/4✓ Branch 1 taken 1540 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1540 times.
✗ Branch 5 not taken.
|
3080 | args->iterations = 1; |
560 |
4/8✓ Branch 1 taken 1540 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1540 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 1540 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 1540 times.
✗ Branch 11 not taken.
|
6160 | args->relaxationFactor = 1; |
561 | }); | ||
562 | |||
563 |
1/2✓ Branch 1 taken 1540 times.
✗ Branch 2 not taken.
|
3080 | auto linearOperator = makeLinearOperator_<LinearOperator, Matrix, Vector>(m, std::make_index_sequence<Matrix::N()>{}); |
564 | |||
565 |
1/2✓ Branch 1 taken 1540 times.
✗ Branch 2 not taken.
|
3080 | BlockDiagAMGPreconditioner<Matrix, Vector, Vector> preconditioner(linearOperator, criterion, smootherArgs); |
566 | |||
567 |
1/4✓ Branch 1 taken 1540 times.
✗ Branch 2 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
3080 | Dumux::ParallelMultiTypeMatrixAdapter<Matrix, Vector, Vector> op(m); |
568 |
5/10✓ Branch 1 taken 1540 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1540 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 1540 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 1540 times.
✗ Branch 11 not taken.
✓ Branch 12 taken 1540 times.
✗ Branch 13 not taken.
|
7700 | Dune::BiCGSTABSolver<Vector> solver(op, preconditioner, this->residReduction(), |
569 | this->maxIter(), this->verbosity()); | ||
570 |
1/2✓ Branch 1 taken 1540 times.
✗ Branch 2 not taken.
|
3080 | auto bTmp(b); |
571 |
1/2✓ Branch 1 taken 1540 times.
✗ Branch 2 not taken.
|
1540 | solver.apply(x, bTmp, result_); |
572 | |||
573 | 3080 | return result_.converged; | |
574 | } | ||
575 | |||
576 | const Dune::InverseOperatorResult& result() const | ||
577 | { | ||
578 | return result_; | ||
579 | } | ||
580 | |||
581 | std::string name() const | ||
582 | { return "block-diagonal AMG-preconditioned BiCGSTAB solver"; } | ||
583 | |||
584 | private: | ||
585 | template<template<class M, std::size_t i> class Criterion, class Matrix, class Params, std::size_t... Is> | ||
586 | 1540 | auto makeCriterion_(const Params& p, std::index_sequence<Is...>) | |
587 | { | ||
588 |
1/4✓ Branch 2 taken 1540 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
1540 | return std::make_tuple(std::make_shared<Criterion<Matrix, Is>>(p)...); |
589 | } | ||
590 | |||
591 | template<template<class M, class X, std::size_t i> class SmootherArgs, class Matrix, class Vector, std::size_t... Is> | ||
592 | 1540 | auto makeSmootherArgs_(std::index_sequence<Is...>) | |
593 | { | ||
594 |
1/4✓ Branch 2 taken 1540 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
1540 | return std::make_tuple(std::make_shared<SmootherArgs<Matrix, Vector, Is>>()...); |
595 | } | ||
596 | |||
597 | template<template<class M, class X, std::size_t i> class LinearOperator, class Matrix, class Vector, std::size_t... Is> | ||
598 | 1540 | auto makeLinearOperator_(const Matrix& m, std::index_sequence<Is...>) | |
599 | { | ||
600 |
3/8✓ Branch 4 taken 1540 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 1540 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 1540 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
|
4620 | return std::make_tuple(std::make_shared<LinearOperator<Matrix, Vector, Is>>(m[Dune::index_constant<Is>{}][Dune::index_constant<Is>{}])...); |
601 | } | ||
602 | |||
603 | Dune::InverseOperatorResult result_; | ||
604 | }; | ||
605 | |||
606 | // \} | ||
607 | |||
608 | } // end namespace Dumux | ||
609 | |||
610 | #endif | ||
611 |