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 Simple ilu0 block diagonal preconditioner that can store decomposition | ||
11 | */ | ||
12 | #ifndef DUMUX_EXAMPLE_1D3D_MULTIDOMAIN_LINEAR_SOLVER_HH | ||
13 | #define DUMUX_EXAMPLE_1D3D_MULTIDOMAIN_LINEAR_SOLVER_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/common/indices.hh> | ||
22 | #include <dune/common/hybridutilities.hh> | ||
23 | |||
24 | #include <dumux/common/parameters.hh> | ||
25 | #include <dumux/common/typetraits/matrix.hh> | ||
26 | #include <dumux/common/typetraits/utility.hh> | ||
27 | #include <dumux/linear/solver.hh> | ||
28 | #include <dumux/linear/preconditioners.hh> | ||
29 | #include <dumux/linear/linearsolverparameters.hh> | ||
30 | #include <dumux/linear/parallelmatrixadapter.hh> | ||
31 | |||
32 | namespace Dumux::Example { | ||
33 | |||
34 | /*! | ||
35 | * \ingroup Linear | ||
36 | * \brief A simple ilu0 block diagonal preconditioner | ||
37 | */ | ||
38 | template<class M, class X, class Y, int blockLevel = 2> | ||
39 | class BlockDiagILU0Preconditioner : public Dune::Preconditioner<X, Y> | ||
40 | { | ||
41 | template<std::size_t i> | ||
42 | using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>; | ||
43 | |||
44 | template<std::size_t i> | ||
45 | using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>; | ||
46 | |||
47 | template<std::size_t i> | ||
48 | using BlockILU = Dune::SeqILU<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>, blockLevel-1>; | ||
49 | |||
50 | using ILUTuple = typename makeFromIndexedType<std::tuple, BlockILU, std::make_index_sequence<M::N()> >::type; | ||
51 | |||
52 | public: | ||
53 | //! \brief The matrix type the preconditioner is for. | ||
54 | using matrix_type = typename std::decay_t<M>; | ||
55 | //! \brief The domain type of the preconditioner. | ||
56 | using domain_type = X; | ||
57 | //! \brief The range type of the preconditioner. | ||
58 | using range_type = Y; | ||
59 | //! \brief The field type of the preconditioner. | ||
60 | using field_type = typename X::field_type; | ||
61 | |||
62 | /*! \brief Constructor. | ||
63 | Constructor gets all parameters to operate the prec. | ||
64 | \param m The (multi type block) matrix to operate on | ||
65 | \param w The relaxation factor | ||
66 | */ | ||
67 | 1 | BlockDiagILU0Preconditioner(const M& m, double w = 1.0) | |
68 |
1/2✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
|
1 | : BlockDiagILU0Preconditioner(m, w, std::make_index_sequence<M::N()>{}) |
69 | { | ||
70 | static_assert(blockLevel >= 2, "Only makes sense for MultiTypeBlockMatrix!"); | ||
71 | } | ||
72 | |||
73 | 100 | void pre (X& v, Y& d) final {} | |
74 | |||
75 | 4817 | void apply (X& v, const Y& d) final | |
76 | { | ||
77 | using namespace Dune::Hybrid; | ||
78 | 4817 | forEach(integralRange(Dune::Hybrid::size(ilu_)), [&](const auto i) | |
79 | { | ||
80 | ✗ | std::get<decltype(i)::value>(ilu_).apply(v[i], d[i]); | |
81 | }); | ||
82 | 4817 | } | |
83 | |||
84 | 100 | void post (X&) final {} | |
85 | |||
86 | //! Category of the preconditioner (see SolverCategory::Category) | ||
87 | 1 | Dune::SolverCategory::Category category() const final | |
88 | { | ||
89 | 1 | return Dune::SolverCategory::sequential; | |
90 | } | ||
91 | |||
92 | private: | ||
93 | template<std::size_t... Is> | ||
94 | 1 | BlockDiagILU0Preconditioner (const M& m, double w, std::index_sequence<Is...> is) | |
95 |
3/6✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 1 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 1 times.
✗ Branch 11 not taken.
|
4 | : ilu_(std::make_tuple(BlockILU<Is>(m[Dune::index_constant<Is>{}][Dune::index_constant<Is>{}], w)...)) |
96 | 1 | {} | |
97 | |||
98 | ILUTuple ilu_; | ||
99 | }; | ||
100 | |||
101 | /*! | ||
102 | * \ingroup Linear | ||
103 | * \brief A simple ilu0 block diagonal preconditioned BiCGSTABSolver | ||
104 | * \note expects a system as a multi-type block-matrix | ||
105 | * | A B | | ||
106 | * | C D | | ||
107 | */ | ||
108 | template<class MatrixType, class VectorType> | ||
109 | class BlockDiagILU0BiCGSTABSolver : public LinearSolver | ||
110 | { | ||
111 | public: | ||
112 | using LinearSolver::LinearSolver; | ||
113 | |||
114 | // Solve saddle-point problem using a Schur complement based preconditioner | ||
115 | template<class Matrix, class Vector> | ||
116 | bool solve(const Matrix& m, Vector& x, const Vector& b) | ||
117 | { | ||
118 | setup(m); | ||
119 | return solve(x, b); | ||
120 | } | ||
121 | |||
122 | 1 | void setup(const MatrixType& m) | |
123 | { | ||
124 |
2/4✗ Branch 1 not taken.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 1 times.
|
1 | preconditioner_ = std::make_unique<BlockDiagILU0Preconditioner<MatrixType, VectorType, VectorType>>(m); |
125 |
2/4✗ Branch 1 not taken.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 1 times.
|
1 | linearOperator_ = std::make_shared<Dumux::ParallelMultiTypeMatrixAdapter<MatrixType, VectorType, VectorType>>(m); |
126 |
2/4✗ Branch 6 not taken.
✓ Branch 7 taken 1 times.
✗ Branch 8 not taken.
✓ Branch 9 taken 1 times.
|
6 | solver_ = std::make_unique<Dune::BiCGSTABSolver<VectorType>>(*linearOperator_, *preconditioner_, this->residReduction(), this->maxIter(), this->verbosity()); |
127 | |||
128 | 1 | isSetup_ = true; | |
129 | 1 | } | |
130 | |||
131 | 100 | bool solve(VectorType& x, const VectorType& b) | |
132 | { | ||
133 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 100 times.
|
100 | assert(isSetup_); |
134 | |||
135 | 100 | auto bTmp(b); | |
136 |
2/4✓ Branch 1 taken 100 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 100 times.
✗ Branch 5 not taken.
|
200 | solver_->apply(x, bTmp, result_); |
137 | |||
138 | 100 | return result_.converged; | |
139 | } | ||
140 | |||
141 | const Dune::InverseOperatorResult& result() const | ||
142 | { | ||
143 | return result_; | ||
144 | } | ||
145 | |||
146 | ✗ | std::string name() const | |
147 |
4/12✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✓ Branch 7 taken 100 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 100 times.
✗ Branch 11 not taken.
✓ Branch 13 taken 100 times.
✗ Branch 14 not taken.
✓ Branch 16 taken 100 times.
✗ Branch 17 not taken.
|
200 | { return "block-diagonal ILU0-preconditioned BiCGSTAB solver"; } |
148 | |||
149 | private: | ||
150 | bool isSetup_ = false; | ||
151 | |||
152 | std::shared_ptr<Dune::MatrixAdapter<MatrixType, VectorType, VectorType>> linearOperator_; | ||
153 | std::unique_ptr<Dune::Preconditioner<VectorType, VectorType>> preconditioner_; | ||
154 | std::unique_ptr<Dune::InverseOperator<VectorType, VectorType>> solver_; | ||
155 | |||
156 | Dune::InverseOperatorResult result_; | ||
157 | }; | ||
158 | |||
159 | } // end namespace Dumux::Example | ||
160 | |||
161 | #endif | ||
162 |