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 preconditioners for iterative solvers | ||
11 | */ | ||
12 | #ifndef DUMUX_LINEAR_PRECONDITIONERS_HH | ||
13 | #define DUMUX_LINEAR_PRECONDITIONERS_HH | ||
14 | |||
15 | #include <dune/common/exceptions.hh> | ||
16 | #include <dune/common/float_cmp.hh> | ||
17 | #include <dune/common/indices.hh> | ||
18 | #include <dune/istl/preconditioners.hh> | ||
19 | #include <dune/istl/paamg/amg.hh> | ||
20 | #include <dune/istl/gsetc.hh> | ||
21 | |||
22 | #if HAVE_UMFPACK | ||
23 | #include <dune/istl/umfpack.hh> | ||
24 | #endif | ||
25 | |||
26 | #include <dumux/assembly/coloring.hh> | ||
27 | #include <dumux/common/parameters.hh> | ||
28 | #include <dumux/common/typetraits/matrix.hh> | ||
29 | #include <dumux/linear/istlsolverregistry.hh> | ||
30 | #include <dumux/parallel/parallel_for.hh> | ||
31 | |||
32 | namespace Dumux { | ||
33 | |||
34 | /*! | ||
35 | * \ingroup Linear | ||
36 | * \brief A preconditioner based on the Uzawa algorithm for saddle-point problems of the form | ||
37 | * \f$ | ||
38 | \begin{pmatrix} | ||
39 | A & B \\ | ||
40 | C & D | ||
41 | \end{pmatrix} | ||
42 | |||
43 | \begin{pmatrix} | ||
44 | u\\ | ||
45 | p | ||
46 | \end{pmatrix} | ||
47 | |||
48 | = | ||
49 | |||
50 | \begin{pmatrix} | ||
51 | f\\ | ||
52 | g | ||
53 | \end{pmatrix} | ||
54 | * \f$ | ||
55 | * | ||
56 | * This preconditioner is especially suited for solving the incompressible (Navier-)Stokes equations. | ||
57 | * Here, \f$D = 0\f$ and \f$B = C^T\f$ if \f$\rho = 1\f$. | ||
58 | * We do not expect good convergence if energy or mass transport is considered. | ||
59 | * | ||
60 | * See: Benzi, M., Golub, G. H., & Liesen, J. (2005). Numerical solution of saddle point problems. Acta numerica, 14, 1-137 \cite benzi2005 and <BR> | ||
61 | * Ho, N., Olson, S. D., & Walker, H. F. (2017). Accelerating the Uzawa algorithm. SIAM Journal on Scientific Computing, 39(5), S461-S476 \cite ho2017 | ||
62 | * | ||
63 | * \tparam M Type of the matrix. | ||
64 | * \tparam X Type of the update. | ||
65 | * \tparam Y Type of the defect. | ||
66 | * \tparam l Preconditioner block level (for compatibility reasons, unused). | ||
67 | */ | ||
68 | template<class M, class X, class Y, int l = 1> | ||
69 | class SeqUzawa : public Dune::Preconditioner<X,Y> | ||
70 | { | ||
71 | static_assert(Dumux::isMultiTypeBlockMatrix<M>::value && M::M() == 2 && M::N() == 2, "SeqUzawa expects a 2x2 MultiTypeBlockMatrix."); | ||
72 | static_assert(l== 1, "SeqUzawa expects a block level of 1."); | ||
73 | |||
74 | using A = std::decay_t<decltype(std::declval<M>()[Dune::Indices::_0][Dune::Indices::_0])>; | ||
75 | using U = std::decay_t<decltype(std::declval<X>()[Dune::Indices::_0])>; | ||
76 | |||
77 | using Comm = Dune::Amg::SequentialInformation; | ||
78 | using LinearOperator = Dune::MatrixAdapter<A, U, U>; | ||
79 | using Smoother = Dune::SeqSSOR<A, U, U>; | ||
80 | using AMGSolverForA = Dune::Amg::AMG<LinearOperator, U, Smoother, Comm>; | ||
81 | |||
82 | public: | ||
83 | //! \brief The matrix type the preconditioner is for. | ||
84 | using matrix_type = M; | ||
85 | //! \brief The domain type of the preconditioner. | ||
86 | using domain_type = X; | ||
87 | //! \brief The range type of the preconditioner. | ||
88 | using range_type = Y; | ||
89 | //! \brief The field type of the preconditioner. | ||
90 | using field_type = typename X::field_type; | ||
91 | //! \brief Scalar type underlying the field_type. | ||
92 | using scalar_field_type = Dune::Simd::Scalar<field_type>; | ||
93 | |||
94 | /*! | ||
95 | * \brief Constructor | ||
96 | * | ||
97 | * \param op The linear operator. | ||
98 | * \param params Collection of parameters. | ||
99 | */ | ||
100 | 71 | SeqUzawa(const std::shared_ptr<const Dune::AssembledLinearOperator<M,X,Y>>& op, const Dune::ParameterTree& params) | |
101 | 142 | : matrix_(op->getmat()) | |
102 |
3/8✓ Branch 1 taken 71 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 71 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 71 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
|
142 | , numIterations_(params.get<std::size_t>("iterations")) |
103 | , relaxationFactor_(params.get<scalar_field_type>("relaxation")) | ||
104 |
3/8✓ Branch 1 taken 71 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 71 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 71 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
|
142 | , verbosity_(params.get<int>("verbosity")) |
105 | , paramGroup_(params.get<std::string>("ParameterGroup")) | ||
106 |
17/42✓ Branch 3 taken 71 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✓ Branch 6 taken 71 times.
✓ Branch 8 taken 71 times.
✗ Branch 9 not taken.
✓ Branch 11 taken 71 times.
✗ Branch 12 not taken.
✓ Branch 14 taken 71 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 71 times.
✓ Branch 19 taken 71 times.
✗ Branch 20 not taken.
✓ Branch 22 taken 71 times.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✓ Branch 25 taken 71 times.
✓ Branch 27 taken 71 times.
✗ Branch 28 not taken.
✓ Branch 30 taken 71 times.
✗ Branch 31 not taken.
✓ Branch 33 taken 71 times.
✗ Branch 34 not taken.
✓ Branch 36 taken 71 times.
✗ Branch 37 not taken.
✓ Branch 39 taken 71 times.
✗ Branch 40 not taken.
✗ Branch 41 not taken.
✓ Branch 42 taken 71 times.
✓ Branch 44 taken 71 times.
✗ Branch 45 not taken.
✓ Branch 47 taken 71 times.
✗ Branch 48 not taken.
✗ Branch 49 not taken.
✗ Branch 50 not taken.
✗ Branch 53 not taken.
✗ Branch 54 not taken.
✗ Branch 55 not taken.
✗ Branch 56 not taken.
✗ Branch 57 not taken.
✗ Branch 58 not taken.
|
284 | , useDirectVelocitySolverForA_(getParamFromGroup<bool>(paramGroup_, "LinearSolver.Preconditioner.DirectSolverForA", false)) |
107 | { | ||
108 |
1/4✓ Branch 1 taken 71 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
|
71 | const bool determineRelaxationFactor = getParamFromGroup<bool>(paramGroup_, "LinearSolver.Preconditioner.DetermineRelaxationFactor", true); |
109 | |||
110 | // AMG is needed for determination of omega | ||
111 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 71 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
71 | if (determineRelaxationFactor || !useDirectVelocitySolverForA_) |
112 |
1/2✓ Branch 1 taken 71 times.
✗ Branch 2 not taken.
|
71 | initAMG_(params); |
113 | |||
114 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 71 times.
|
71 | if (useDirectVelocitySolverForA_) |
115 | ✗ | initUMFPack_(); | |
116 | |||
117 |
1/2✓ Branch 0 taken 71 times.
✗ Branch 1 not taken.
|
71 | if (determineRelaxationFactor) |
118 |
1/2✓ Branch 1 taken 71 times.
✗ Branch 2 not taken.
|
71 | relaxationFactor_ = estimateOmega_(); |
119 | 71 | } | |
120 | |||
121 | /*! | ||
122 | * \brief Prepare the preconditioner. | ||
123 | */ | ||
124 | 71 | virtual void pre(X& x, Y& b) {} | |
125 | |||
126 | /*! | ||
127 | * \brief Apply the preconditioner | ||
128 | * | ||
129 | * \param update The update to be computed. | ||
130 | * \param currentDefect The current defect. | ||
131 | */ | ||
132 | 2440 | virtual void apply(X& update, const Y& currentDefect) | |
133 | { | ||
134 | using namespace Dune::Indices; | ||
135 | |||
136 | 4880 | auto& A = matrix_[_0][_0]; | |
137 | 4880 | auto& B = matrix_[_0][_1]; | |
138 | 4880 | auto& C = matrix_[_1][_0]; | |
139 | 4880 | auto& D = matrix_[_1][_1]; | |
140 | |||
141 | 2440 | const auto& f = currentDefect[_0]; | |
142 | 2440 | const auto& g = currentDefect[_1]; | |
143 | 2440 | auto& u = update[_0]; | |
144 | 2440 | auto& p = update[_1]; | |
145 | |||
146 | // incorporate Dirichlet cell values | ||
147 | // TODO: pass Dirichlet constraint handler from outside | ||
148 |
2/2✓ Branch 0 taken 5884500 times.
✓ Branch 1 taken 2440 times.
|
5886940 | for (std::size_t i = 0; i < D.N(); ++i) |
149 | { | ||
150 | 5884500 | const auto& block = D[i][i]; | |
151 |
4/4✓ Branch 1 taken 5884500 times.
✓ Branch 2 taken 5884500 times.
✓ Branch 3 taken 5884500 times.
✓ Branch 4 taken 5884500 times.
|
29422500 | for (auto rowIt = block.begin(); rowIt != block.end(); ++rowIt) |
152 |
8/8✓ Branch 0 taken 5811255 times.
✓ Branch 1 taken 73245 times.
✓ Branch 2 taken 5811255 times.
✓ Branch 3 taken 73245 times.
✓ Branch 4 taken 5811255 times.
✓ Branch 5 taken 73245 times.
✓ Branch 6 taken 5811255 times.
✓ Branch 7 taken 73245 times.
|
23464755 | if (Dune::FloatCmp::eq<scalar_field_type>(rowIt->one_norm(), 1.0)) |
153 | 219735 | p[i][rowIt.index()] = g[i][rowIt.index()]; | |
154 | } | ||
155 | |||
156 | // the actual Uzawa iteration | ||
157 |
2/2✓ Branch 0 taken 9310 times.
✓ Branch 1 taken 2440 times.
|
11750 | for (std::size_t k = 0; k < numIterations_; ++k) |
158 | { | ||
159 | // u_k+1 = u_k + Q_A^−1*(f − (A*u_k + B*p_k)), | ||
160 |
0/2✗ Branch 1 not taken.
✗ Branch 2 not taken.
|
18620 | auto uRhs = f; |
161 | 9310 | A.mmv(u, uRhs); | |
162 | 9310 | B.mmv(p, uRhs); | |
163 |
2/6✓ Branch 1 taken 9310 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 9310 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
|
18620 | auto uIncrement = u; |
164 |
1/2✓ Branch 1 taken 9310 times.
✗ Branch 2 not taken.
|
9310 | applySolverForA_(uIncrement, uRhs); |
165 | 9310 | u += uIncrement; | |
166 | |||
167 | // p_k+1 = p_k + omega*(g - C*u_k+1 - D*p_k) | ||
168 |
3/8✓ Branch 1 taken 9310 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 9310 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 9310 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
|
27930 | auto pIncrement = g; |
169 | 9310 | C.mmv(u, pIncrement); | |
170 | 9310 | D.mmv(p, pIncrement); | |
171 | 9310 | pIncrement *= relaxationFactor_; | |
172 | 9310 | p += pIncrement; | |
173 | |||
174 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 9310 times.
|
9310 | if (verbosity_ > 1) |
175 | { | ||
176 | ✗ | std::cout << "Uzawa iteration " << k | |
177 | ✗ | << ", residual: " << uRhs.two_norm() + pIncrement.two_norm()/relaxationFactor_ << std::endl; | |
178 | } | ||
179 | } | ||
180 | 2440 | } | |
181 | |||
182 | /*! | ||
183 | * \brief Clean up. | ||
184 | */ | ||
185 | 71 | virtual void post(X& x) {} | |
186 | |||
187 | //! Category of the preconditioner (see SolverCategory::Category) | ||
188 | 142 | virtual Dune::SolverCategory::Category category() const | |
189 | { | ||
190 | 142 | return Dune::SolverCategory::sequential; | |
191 | } | ||
192 | |||
193 | private: | ||
194 | |||
195 | 71 | void initAMG_(const Dune::ParameterTree& params) | |
196 | { | ||
197 | using namespace Dune::Indices; | ||
198 |
0/2✗ Branch 3 not taken.
✗ Branch 4 not taken.
|
213 | auto linearOperator = std::make_shared<LinearOperator>(matrix_[_0][_0]); |
199 |
4/8✓ Branch 1 taken 71 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 71 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 71 times.
✓ Branch 7 taken 71 times.
✗ Branch 8 not taken.
|
71 | amgSolverForA_ = std::make_unique<AMGSolverForA>(linearOperator, params); |
200 | 71 | } | |
201 | |||
202 | ✗ | void initUMFPack_() | |
203 | { | ||
204 | #if HAVE_UMFPACK | ||
205 | using namespace Dune::Indices; | ||
206 | ✗ | umfPackSolverForA_ = std::make_unique<Dune::UMFPack<A>>(matrix_[_0][_0]); | |
207 | #else | ||
208 | DUNE_THROW(Dune::InvalidStateException, "UMFPack not available. Use LinearSolver.Preconditioner.DirectVelocitySolver = false."); | ||
209 | #endif | ||
210 | ✗ | } | |
211 | |||
212 | /*! | ||
213 | * \brief Estimate the relaxation factor omega | ||
214 | * | ||
215 | * The optimal relaxation factor is omega = 2/(lambdaMin + lambdaMax), where lambdaMin and lambdaMax are the smallest and largest | ||
216 | * eigenvalues of the Schur complement -C*Ainv*B (assuming D = 0). | ||
217 | * lambdaMax can be easily determined using the power iteration algorithm (https://en.wikipedia.org/wiki/Power_iteration) and lambdaMin | ||
218 | * could be estimated in a similar manner. We do not consider lambdaMin because for certain cases, e.g., when C contains some rows of zeroes only, | ||
219 | * this estimate will fail. | ||
220 | * | ||
221 | * Instead we assume that lambdaMin is sufficiently close to lambdaMax such that omega = 1/lambdaMax. | ||
222 | * This seems to work rather well for various applications. | ||
223 | * We will underestimate omega by a factor of 2 in the worst case (i.e, lambdaMin = 0). | ||
224 | * | ||
225 | * When facing convergence issues, you may set LinearSolver.Preconditioner.Verbosity = 1 to see the estimate of lambdaMax. | ||
226 | * In a new simulation run, you can then set LinearSolver.Preconditioner.DetermineRelaxationFactor = false and set some other value | ||
227 | * for LinearSolver.Preconditioner.Relaxation based on the estimate of lambdaMax. | ||
228 | * | ||
229 | * See: Benzi, M., Golub, G. H., & Liesen, J. (2005). Numerical solution of saddle point problems. Acta numerica, 14, 1-137. | ||
230 | */ | ||
231 | 71 | scalar_field_type estimateOmega_() const | |
232 | { | ||
233 | using namespace Dune::Indices; | ||
234 | 142 | auto& A = matrix_[_0][_0]; | |
235 | 142 | auto& B = matrix_[_0][_1]; | |
236 | 142 | auto& C = matrix_[_1][_0]; | |
237 | |||
238 |
0/2✗ Branch 1 not taken.
✗ Branch 2 not taken.
|
71 | U x(A.M()); |
239 | 71 | x = 1.0; | |
240 | |||
241 | 71 | scalar_field_type omega = 0.0; | |
242 | 71 | scalar_field_type lambdaMax = 0.0; | |
243 | |||
244 |
4/6✓ Branch 0 taken 4 times.
✓ Branch 1 taken 67 times.
✓ Branch 3 taken 4 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 4 times.
✗ Branch 7 not taken.
|
71 | static const auto iterations = Dumux::getParamFromGroup<std::size_t>(paramGroup_, "LinearSolver.Preconditioner.PowerLawIterations", 5); |
245 | |||
246 | // apply power iteration x_k+1 = M*x_k/|M*x_k| for the matrix M = -C*Ainv*B | ||
247 |
2/2✓ Branch 0 taken 355 times.
✓ Branch 1 taken 71 times.
|
426 | for (std::size_t i = 0; i < iterations; ++i) |
248 | { | ||
249 | // bx = B*x | ||
250 |
1/4✓ Branch 1 taken 355 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
|
710 | U bx(x.size()); |
251 | 355 | B.mv(x, bx); | |
252 | |||
253 | // ainvbx = Ainv*(B*x) | ||
254 |
2/6✓ Branch 1 taken 355 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 355 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
|
710 | auto ainvbx = x; |
255 |
1/2✓ Branch 1 taken 355 times.
✗ Branch 2 not taken.
|
355 | applySolverForA_(ainvbx, bx); |
256 | |||
257 | // v = M*x = -C*(Ainv*B*x) | ||
258 |
2/6✓ Branch 1 taken 355 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 355 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
|
710 | U v(x.size()); |
259 | 355 | C.mv(ainvbx, v); | |
260 | 355 | v *= -1.0; | |
261 | |||
262 | // eigenvalue lambdaMax = xt*M*x/(xt*x) = xt*v/(xt*x); | ||
263 | 355 | lambdaMax = x.dot(v)/(x.dot(x)); | |
264 | |||
265 | // relaxation factor omega = 1/lambda; | ||
266 | 355 | omega = 1.0/lambdaMax; | |
267 | |||
268 | // new iterate x = M*x/|M*x| = v/|v| | ||
269 |
1/2✓ Branch 1 taken 355 times.
✗ Branch 2 not taken.
|
355 | x = v; |
270 |
1/2✓ Branch 1 taken 355 times.
✗ Branch 2 not taken.
|
710 | x /= v.two_norm(); |
271 | } | ||
272 | |||
273 |
2/2✓ Branch 0 taken 51 times.
✓ Branch 1 taken 20 times.
|
71 | if (verbosity_ > 0) |
274 | { | ||
275 |
1/2✓ Branch 2 taken 51 times.
✗ Branch 3 not taken.
|
102 | std::cout << "\n*** Uzawa Preconditioner ***" << std::endl; |
276 |
1/2✓ Branch 2 taken 51 times.
✗ Branch 3 not taken.
|
102 | std::cout << "Estimating relaxation factor based on Schur complement" << std::endl; |
277 |
2/4✓ Branch 2 taken 51 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 51 times.
✗ Branch 6 not taken.
|
102 | std::cout << "Largest estimated eigenvalue lambdaMax = " << lambdaMax << std::endl; |
278 |
2/4✓ Branch 2 taken 51 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 51 times.
✗ Branch 6 not taken.
|
102 | std::cout << "Relaxation factor omega = 1/lambdaMax = " << omega << std::endl; |
279 | } | ||
280 | |||
281 |
1/2✓ Branch 0 taken 71 times.
✗ Branch 1 not taken.
|
142 | return omega; |
282 | } | ||
283 | |||
284 | template<class Sol, class Rhs> | ||
285 | 9665 | void applySolverForA_(Sol& sol, Rhs& rhs) const | |
286 | { | ||
287 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 9665 times.
|
9665 | if (useDirectVelocitySolverForA_) |
288 | { | ||
289 | #if HAVE_UMFPACK | ||
290 | ✗ | Dune::InverseOperatorResult res; | |
291 | ✗ | umfPackSolverForA_->apply(sol, rhs, res); | |
292 | #endif | ||
293 | } | ||
294 | else | ||
295 | { | ||
296 | 19330 | amgSolverForA_->pre(sol, rhs); | |
297 | 19330 | amgSolverForA_->apply(sol, rhs); | |
298 | 19330 | amgSolverForA_->post(sol); | |
299 | } | ||
300 | 9665 | } | |
301 | |||
302 | //! \brief The matrix we operate on. | ||
303 | const M& matrix_; | ||
304 | //! \brief The number of steps to do in apply | ||
305 | const std::size_t numIterations_; | ||
306 | //! \brief The relaxation factor to use | ||
307 | scalar_field_type relaxationFactor_; | ||
308 | //! \brief The verbosity level | ||
309 | const int verbosity_; | ||
310 | |||
311 | std::unique_ptr<AMGSolverForA> amgSolverForA_; | ||
312 | #if HAVE_UMFPACK | ||
313 | std::unique_ptr<Dune::UMFPack<A>> umfPackSolverForA_; | ||
314 | #endif | ||
315 | const std::string paramGroup_; | ||
316 | const bool useDirectVelocitySolverForA_; | ||
317 | }; | ||
318 | |||
319 |
2/14✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 1 times.
✗ 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.
|
3 | DUMUX_REGISTER_PRECONDITIONER("uzawa", Dumux::MultiTypeBlockMatrixPreconditionerTag, Dune::defaultPreconditionerBlockLevelCreator<Dumux::SeqUzawa, 1>()); |
320 | |||
321 | |||
322 | /*! \brief Multi-threaded Jacobi preconditioner | ||
323 | * | ||
324 | * \tparam M The matrix type to operate on | ||
325 | * \tparam X Type of the update | ||
326 | * \tparam Y Type of the defect | ||
327 | * \tparam l The block level to invert. Default is 1 | ||
328 | */ | ||
329 | template<class M, class X, class Y, int l=1> | ||
330 | class ParMTJac : public Dune::Preconditioner<X,Y> | ||
331 | { | ||
332 | public: | ||
333 | //! \brief The matrix type the preconditioner is for. | ||
334 | typedef M matrix_type; | ||
335 | //! \brief The domain type of the preconditioner. | ||
336 | typedef X domain_type; | ||
337 | //! \brief The range type of the preconditioner. | ||
338 | typedef Y range_type; | ||
339 | //! \brief The field type of the preconditioner. | ||
340 | typedef typename X::field_type field_type; | ||
341 | //! \brief scalar type underlying the field_type | ||
342 | typedef Dune::Simd::Scalar<field_type> scalar_field_type; | ||
343 | //! \brief real scalar type underlying the field_type | ||
344 | typedef typename Dune::FieldTraits<scalar_field_type>::real_type real_field_type; | ||
345 | |||
346 | //! \brief Constructor. | ||
347 | 37 | ParMTJac (const M& A, int n, real_field_type w) | |
348 | 37 | : _A_(A), numSteps_(n), relaxationFactor_(w) | |
349 | 37 | { Dune::CheckIfDiagonalPresent<M,l>::check(_A_); } | |
350 | |||
351 | 37 | ParMTJac (const std::shared_ptr<const Dune::AssembledLinearOperator<M,X,Y>>& A, const Dune::ParameterTree& configuration) | |
352 | 74 | : ParMTJac(A->getmat(), configuration) | |
353 | 37 | {} | |
354 | |||
355 | 37 | ParMTJac (const M& A, const Dune::ParameterTree& configuration) | |
356 |
9/22✓ Branch 1 taken 37 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 37 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 37 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 37 times.
✗ Branch 11 not taken.
✓ Branch 13 taken 37 times.
✗ Branch 14 not taken.
✓ Branch 16 taken 37 times.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✓ Branch 19 taken 37 times.
✗ Branch 20 not taken.
✓ Branch 21 taken 37 times.
✗ Branch 22 not taken.
✓ Branch 23 taken 37 times.
✗ Branch 24 not taken.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✗ Branch 27 not taken.
|
74 | : ParMTJac(A, configuration.get<int>("iterations",1), configuration.get<real_field_type>("relaxation",1.0)) |
357 | 37 | {} | |
358 | |||
359 | 40 | void pre (X&, Y&) override {} | |
360 | |||
361 | 543 | void apply (X& update, const Y& defect) override | |
362 | { | ||
363 |
1/4✓ Branch 1 taken 543 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
|
1086 | X xOld(update); |
364 | 543 | const auto& A = _A_; | |
365 | |||
366 |
2/2✓ Branch 0 taken 795 times.
✓ Branch 1 taken 543 times.
|
1338 | for (int k=0; k<numSteps_; k++) |
367 | { | ||
368 |
2/2✓ Branch 0 taken 252 times.
✓ Branch 1 taken 543 times.
|
795 | if (k > 0) |
369 |
1/2✓ Branch 1 taken 252 times.
✗ Branch 2 not taken.
|
252 | xOld = update; |
370 | |||
371 | 795 | Dumux::parallelFor(A.N(), [&](const std::size_t i) | |
372 | { | ||
373 | 19956744 | auto& row = A[i]; | |
374 | 19956744 | auto v = update[i]; | |
375 | 19956744 | auto rhs = defect[i]; | |
376 | 19956744 | const auto endColIt = row.end(); | |
377 | 19956744 | auto colIt = row.begin(); | |
378 | |||
379 |
4/4✓ Branch 0 taken 39725448 times.
✓ Branch 1 taken 19956744 times.
✓ Branch 2 taken 39725448 times.
✓ Branch 3 taken 19956744 times.
|
119364384 | for (; colIt.index()<i; ++colIt) |
380 | 158901792 | colIt->mmv(xOld[colIt.index()], rhs); | |
381 | 19956744 | const auto diag = colIt; | |
382 |
4/4✓ Branch 0 taken 59682192 times.
✓ Branch 1 taken 19956744 times.
✓ Branch 2 taken 59682192 times.
✓ Branch 3 taken 19956744 times.
|
159277872 | for (; colIt != endColIt; ++colIt) |
383 | 298410960 | colIt->mmv(xOld[colIt.index()], rhs); | |
384 | |||
385 | 39913488 | Dune::algmeta_itsteps<l-1, typename M::block_type>::dbjac( | |
386 | 39913488 | *diag, v, rhs, relaxationFactor_ | |
387 | ); | ||
388 | |||
389 | 39913488 | update[i].axpy(relaxationFactor_, v); | |
390 | }); | ||
391 | } | ||
392 | 543 | } | |
393 | |||
394 | 40 | void post (X&) override {} | |
395 | |||
396 | //! Category of the preconditioner (see SolverCategory::Category) | ||
397 | 74 | Dune::SolverCategory::Category category() const override | |
398 | 74 | { return Dune::SolverCategory::sequential; } | |
399 | |||
400 | private: | ||
401 | const M& _A_; | ||
402 | int numSteps_; | ||
403 | real_field_type relaxationFactor_; | ||
404 | }; | ||
405 | |||
406 |
6/42✓ Branch 4 taken 14 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 14 times.
✗ Branch 8 not taken.
✓ Branch 13 taken 14 times.
✗ Branch 14 not taken.
✓ Branch 16 taken 14 times.
✗ Branch 17 not taken.
✓ Branch 22 taken 14 times.
✗ Branch 23 not taken.
✓ Branch 25 taken 14 times.
✗ 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 37 not taken.
✗ Branch 38 not taken.
✗ Branch 40 not taken.
✗ Branch 41 not taken.
✗ Branch 43 not taken.
✗ Branch 44 not taken.
✗ Branch 46 not taken.
✗ Branch 47 not taken.
✗ Branch 49 not taken.
✗ Branch 50 not taken.
✗ Branch 52 not taken.
✗ Branch 53 not taken.
✗ Branch 55 not taken.
✗ Branch 56 not taken.
✗ Branch 58 not taken.
✗ Branch 59 not taken.
✗ Branch 61 not taken.
✗ Branch 62 not taken.
✗ Branch 64 not taken.
✗ Branch 65 not taken.
✗ Branch 67 not taken.
✗ Branch 68 not taken.
✗ Branch 70 not taken.
✗ Branch 71 not taken.
|
126 | DUMUX_REGISTER_PRECONDITIONER("par_mt_jac", Dune::PreconditionerTag, Dune::defaultPreconditionerBlockLevelCreator<Dumux::ParMTJac>()); |
407 | |||
408 | |||
409 | namespace Detail::Preconditioners { | ||
410 | |||
411 | // compute coloring for parallel sweep | ||
412 | template<class M> | ||
413 | 82 | void computeColorsForMatrixSweep_(const M& matrix, std::deque<std::vector<std::size_t>>& coloredIndices) | |
414 | { | ||
415 | // allocate some temporary memory | ||
416 |
1/4✓ Branch 3 taken 82 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
|
164 | std::vector<int> colors(matrix.N(), -1); |
417 |
3/8✓ Branch 1 taken 82 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 82 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 82 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
|
246 | std::vector<int> neighborColors; neighborColors.reserve(30); |
418 |
3/6✓ Branch 1 taken 82 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 82 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 82 times.
✗ Branch 8 not taken.
|
246 | std::vector<bool> colorUsed; colorUsed.reserve(30); |
419 | |||
420 | // a row that has my row index in their column set cannot have the same color | ||
421 | // TODO this assumes a symmetric matrix pattern | ||
422 |
4/4✓ Branch 0 taken 267478 times.
✓ Branch 1 taken 82 times.
✓ Branch 2 taken 267478 times.
✓ Branch 3 taken 82 times.
|
267560 | for (std::size_t i = 0; i < colors.size(); ++i) |
423 | { | ||
424 |
2/2✓ Branch 0 taken 267396 times.
✓ Branch 1 taken 82 times.
|
267478 | neighborColors.clear(); |
425 | 267478 | auto& row = matrix[i]; | |
426 | 267478 | const auto endColIt = row.end(); | |
427 |
4/4✓ Branch 0 taken 1327118 times.
✓ Branch 1 taken 267478 times.
✓ Branch 2 taken 1327118 times.
✓ Branch 3 taken 267478 times.
|
1594596 | for (auto colIt = row.begin(); colIt != endColIt; ++colIt) |
428 |
3/6✓ Branch 1 taken 1327118 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1327118 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 1327118 times.
✗ Branch 8 not taken.
|
3981354 | neighborColors.push_back(colors[colIt.index()]); |
429 | |||
430 |
1/2✓ Branch 1 taken 267478 times.
✗ Branch 2 not taken.
|
267478 | const auto c = Detail::smallestAvailableColor(neighborColors, colorUsed); |
431 |
2/2✓ Branch 0 taken 267310 times.
✓ Branch 1 taken 168 times.
|
267478 | colors[i] = c; |
432 | |||
433 | // add element to the set of elements with the same color | ||
434 |
4/4✓ Branch 0 taken 267310 times.
✓ Branch 1 taken 168 times.
✓ Branch 2 taken 267310 times.
✓ Branch 3 taken 168 times.
|
534956 | if (c < coloredIndices.size()) |
435 |
1/2✓ Branch 2 taken 267310 times.
✗ Branch 3 not taken.
|
267310 | coloredIndices[c].push_back(i); |
436 | else | ||
437 |
1/2✓ Branch 1 taken 168 times.
✗ Branch 2 not taken.
|
168 | coloredIndices.emplace_back(1, i); |
438 | } | ||
439 | 82 | } | |
440 | |||
441 | // parallel SOR kernel (relaxed Gauss-Seidel) | ||
442 | template<bool forward, int l, class M, class X, class Y, class K> | ||
443 | 4040 | void parallelBlockSOR_(const M& A, X& update, const Y& b, const K& relaxationFactor, | |
444 | const std::deque<std::vector<std::size_t>>& coloredIndices) | ||
445 | { | ||
446 |
4/4✓ Branch 0 taken 4040 times.
✓ Branch 1 taken 2020 times.
✓ Branch 2 taken 4040 times.
✓ Branch 3 taken 2020 times.
|
20200 | for (int color = 0; color < coloredIndices.size(); ++color) |
447 | { | ||
448 | 8080 | const auto c = forward ? color : coloredIndices.size()-1-color; | |
449 | 8080 | const auto& rowIndices = coloredIndices[c]; | |
450 | 117604736 | Dumux::parallelFor(rowIndices.size(), [&](const std::size_t k) | |
451 | { | ||
452 | 58794288 | const auto i = rowIndices[k]; | |
453 | 58794288 | auto& row = A[i]; | |
454 | 58794288 | auto v = update[i]; | |
455 | 58794288 | auto rhs = b[i]; | |
456 | 58794288 | const auto endColIt = row.end(); | |
457 | 58794288 | auto colIt = row.begin(); | |
458 | |||
459 |
8/8✓ Branch 0 taken 58454088 times.
✓ Branch 1 taken 29360664 times.
✓ Branch 2 taken 58454088 times.
✓ Branch 3 taken 29360664 times.
✓ Branch 4 taken 58592408 times.
✓ Branch 5 taken 29433624 times.
✓ Branch 6 taken 58592408 times.
✓ Branch 7 taken 29433624 times.
|
351681568 | for (; colIt.index()<i; ++colIt) |
460 | 468185984 | colIt->mmv(update[colIt.index()], rhs); | |
461 | 58794288 | const auto diag = colIt; | |
462 |
8/8✓ Branch 0 taken 87814752 times.
✓ Branch 1 taken 29360664 times.
✓ Branch 2 taken 87814752 times.
✓ Branch 3 taken 29360664 times.
✓ Branch 4 taken 88026032 times.
✓ Branch 5 taken 29433624 times.
✓ Branch 6 taken 88026032 times.
✓ Branch 7 taken 29433624 times.
|
469270144 | for (; colIt != endColIt; ++colIt) |
463 | 879203920 | colIt->mmv(update[colIt.index()], rhs); | |
464 | |||
465 | if constexpr (forward) | ||
466 | 58867248 | Dune::algmeta_itsteps<l-1,typename M::block_type>::bsorf( | |
467 | 88300872 | *diag, v, rhs, relaxationFactor | |
468 | ); | ||
469 | else | ||
470 | 58721328 | Dune::algmeta_itsteps<l-1,typename M::block_type>::bsorb( | |
471 | 88081992 | *diag, v, rhs, relaxationFactor | |
472 | ); | ||
473 | |||
474 | 176382864 | update[i].axpy(relaxationFactor, v); | |
475 | }); | ||
476 | } | ||
477 | 4040 | } | |
478 | |||
479 | } // end namespace Detail | ||
480 | |||
481 | /*! \brief Multi-threaded SOR preconditioner using coloring | ||
482 | * | ||
483 | * \tparam M The matrix type to operate on | ||
484 | * \tparam X Type of the update | ||
485 | * \tparam Y Type of the defect | ||
486 | * \tparam l The block level to invert. Default is 1 | ||
487 | */ | ||
488 | template<class M, class X, class Y, int l=1> | ||
489 | class ParMTSOR : public Dune::Preconditioner<X,Y> | ||
490 | { | ||
491 | public: | ||
492 | //! \brief The matrix type the preconditioner is for. | ||
493 | typedef M matrix_type; | ||
494 | //! \brief The domain type of the preconditioner. | ||
495 | typedef X domain_type; | ||
496 | //! \brief The range type of the preconditioner. | ||
497 | typedef Y range_type; | ||
498 | //! \brief The field type of the preconditioner. | ||
499 | typedef typename X::field_type field_type; | ||
500 | //! \brief scalar type underlying the field_type | ||
501 | typedef Dune::Simd::Scalar<field_type> scalar_field_type; | ||
502 | //! \brief real scalar type underlying the field_type | ||
503 | typedef typename Dune::FieldTraits<scalar_field_type>::real_type real_field_type; | ||
504 | |||
505 | //! \brief Constructor. | ||
506 | 41 | ParMTSOR (const M& A, int n, real_field_type w) | |
507 | 41 | : _A_(A), numSteps_(n), relaxationFactor_(w) | |
508 | { | ||
509 | 41 | Dune::CheckIfDiagonalPresent<M,l>::check(_A_); | |
510 |
1/2✓ Branch 1 taken 41 times.
✗ Branch 2 not taken.
|
41 | Detail::Preconditioners::computeColorsForMatrixSweep_(_A_, colors_); |
511 | 41 | } | |
512 | |||
513 | 37 | ParMTSOR (const std::shared_ptr<const Dune::AssembledLinearOperator<M,X,Y>>& A, const Dune::ParameterTree& configuration) | |
514 | 74 | : ParMTSOR(A->getmat(), configuration) | |
515 | 37 | {} | |
516 | |||
517 | 37 | ParMTSOR (const M& A, const Dune::ParameterTree& configuration) | |
518 |
10/24✓ Branch 1 taken 37 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 37 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 37 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 37 times.
✗ Branch 11 not taken.
✓ Branch 13 taken 37 times.
✗ Branch 14 not taken.
✓ Branch 16 taken 37 times.
✗ Branch 17 not taken.
✓ Branch 19 taken 37 times.
✗ Branch 20 not taken.
✗ Branch 21 not taken.
✓ Branch 22 taken 37 times.
✗ Branch 23 not taken.
✓ Branch 24 taken 37 times.
✗ Branch 25 not taken.
✓ Branch 26 taken 37 times.
✗ Branch 27 not taken.
✗ Branch 28 not taken.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
|
74 | : ParMTSOR(A, configuration.get<int>("iterations",1), configuration.get<real_field_type>("relaxation",1.0)) |
519 | 37 | {} | |
520 | |||
521 | 40 | void pre (X&, Y&) override {} | |
522 | |||
523 | 190 | void apply (X& v, const Y& d) override | |
524 | { | ||
525 | 380 | this->template apply<true>(v,d); | |
526 | 190 | } | |
527 | |||
528 | template<bool forward> | ||
529 | void apply(X& v, const Y& d) | ||
530 | { | ||
531 |
4/6✓ Branch 0 taken 430 times.
✓ Branch 1 taken 310 times.
✓ Branch 2 taken 240 times.
✓ Branch 3 taken 120 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
1100 | for (int i=0; i<numSteps_; i++) |
532 | 670 | Detail::Preconditioners::parallelBlockSOR_<forward, l>(_A_, v, d, relaxationFactor_, colors_); | |
533 | } | ||
534 | |||
535 | 40 | void post (X&) override {} | |
536 | |||
537 | //! Category of the preconditioner (see SolverCategory::Category) | ||
538 | 74 | Dune::SolverCategory::Category category() const override | |
539 | 74 | { return Dune::SolverCategory::sequential; } | |
540 | |||
541 | private: | ||
542 | const M& _A_; | ||
543 | int numSteps_; | ||
544 | real_field_type relaxationFactor_; | ||
545 | |||
546 | //! for each color a vector of row indices that can be dealt with in parallel | ||
547 | std::deque<std::vector<std::size_t>> colors_; | ||
548 | }; | ||
549 | |||
550 |
6/42✓ Branch 4 taken 14 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 14 times.
✗ Branch 8 not taken.
✓ Branch 13 taken 14 times.
✗ Branch 14 not taken.
✓ Branch 16 taken 14 times.
✗ Branch 17 not taken.
✓ Branch 22 taken 14 times.
✗ Branch 23 not taken.
✓ Branch 25 taken 14 times.
✗ 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 37 not taken.
✗ Branch 38 not taken.
✗ Branch 40 not taken.
✗ Branch 41 not taken.
✗ Branch 43 not taken.
✗ Branch 44 not taken.
✗ Branch 46 not taken.
✗ Branch 47 not taken.
✗ Branch 49 not taken.
✗ Branch 50 not taken.
✗ Branch 52 not taken.
✗ Branch 53 not taken.
✗ Branch 55 not taken.
✗ Branch 56 not taken.
✗ Branch 58 not taken.
✗ Branch 59 not taken.
✗ Branch 61 not taken.
✗ Branch 62 not taken.
✗ Branch 64 not taken.
✗ Branch 65 not taken.
✗ Branch 67 not taken.
✗ Branch 68 not taken.
✗ Branch 70 not taken.
✗ Branch 71 not taken.
|
126 | DUMUX_REGISTER_PRECONDITIONER("par_mt_sor", Dune::PreconditionerTag, Dune::defaultPreconditionerBlockLevelCreator<Dumux::ParMTSOR>()); |
551 | |||
552 | |||
553 | /*! \brief Multithreaded SSOR preconditioner using coloring | ||
554 | * | ||
555 | * \tparam M The matrix type to operate on | ||
556 | * \tparam X Type of the update | ||
557 | * \tparam Y Type of the defect | ||
558 | * \tparam l The block level to invert. Default is 1 | ||
559 | */ | ||
560 | template<class M, class X, class Y, int l=1> | ||
561 | class ParMTSSOR : public Dune::Preconditioner<X,Y> | ||
562 | { | ||
563 | public: | ||
564 | //! \brief The matrix type the preconditioner is for. | ||
565 | typedef M matrix_type; | ||
566 | //! \brief The domain type of the preconditioner. | ||
567 | typedef X domain_type; | ||
568 | //! \brief The range type of the preconditioner. | ||
569 | typedef Y range_type; | ||
570 | //! \brief The field type of the preconditioner. | ||
571 | typedef typename X::field_type field_type; | ||
572 | //! \brief scalar type underlying the field_type | ||
573 | typedef Dune::Simd::Scalar<field_type> scalar_field_type; | ||
574 | //! \brief real scalar type underlying the field_type | ||
575 | typedef typename Dune::FieldTraits<scalar_field_type>::real_type real_field_type; | ||
576 | |||
577 | //! \brief Constructor. | ||
578 | 41 | ParMTSSOR (const M& A, int n, real_field_type w) | |
579 | 41 | : _A_(A), numSteps_(n), relaxationFactor_(w) | |
580 | { | ||
581 | 41 | Dune::CheckIfDiagonalPresent<M,l>::check(_A_); | |
582 |
1/2✓ Branch 1 taken 41 times.
✗ Branch 2 not taken.
|
41 | Detail::Preconditioners::computeColorsForMatrixSweep_(_A_, colors_); |
583 | 41 | } | |
584 | |||
585 | 37 | ParMTSSOR (const std::shared_ptr<const Dune::AssembledLinearOperator<M,X,Y>>& A, const Dune::ParameterTree& configuration) | |
586 | 74 | : ParMTSSOR(A->getmat(), configuration) | |
587 | 37 | {} | |
588 | |||
589 | 37 | ParMTSSOR (const M& A, const Dune::ParameterTree& configuration) | |
590 |
10/24✓ Branch 1 taken 37 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 37 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 37 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 37 times.
✗ Branch 11 not taken.
✓ Branch 13 taken 37 times.
✗ Branch 14 not taken.
✓ Branch 16 taken 37 times.
✗ Branch 17 not taken.
✓ Branch 19 taken 37 times.
✗ Branch 20 not taken.
✗ Branch 21 not taken.
✓ Branch 22 taken 37 times.
✗ Branch 23 not taken.
✓ Branch 24 taken 37 times.
✗ Branch 25 not taken.
✓ Branch 26 taken 37 times.
✗ Branch 27 not taken.
✗ Branch 28 not taken.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
|
74 | : ParMTSSOR(A, configuration.get<int>("iterations",1), configuration.get<real_field_type>("relaxation",1.0)) |
591 | 37 | {} | |
592 | |||
593 | 40 | void pre (X&, Y&) override {} | |
594 | |||
595 | 423 | void apply (X& v, const Y& d) override | |
596 | { | ||
597 |
2/2✓ Branch 0 taken 675 times.
✓ Branch 1 taken 423 times.
|
1098 | for (int i=0; i<numSteps_; i++) |
598 | { | ||
599 | 675 | Detail::Preconditioners::parallelBlockSOR_<true, l>(_A_, v, d, relaxationFactor_, colors_); | |
600 | 675 | Detail::Preconditioners::parallelBlockSOR_<false, l>(_A_, v, d, relaxationFactor_, colors_); | |
601 | } | ||
602 | 423 | } | |
603 | |||
604 | 40 | void post (X&) override {} | |
605 | |||
606 | //! Category of the preconditioner (see SolverCategory::Category) | ||
607 | 74 | Dune::SolverCategory::Category category() const override | |
608 | 74 | { return Dune::SolverCategory::sequential; } | |
609 | |||
610 | private: | ||
611 | const M& _A_; | ||
612 | int numSteps_; | ||
613 | real_field_type relaxationFactor_; | ||
614 | |||
615 | //! for each color a vector of row indices that can be dealt with in parallel | ||
616 | std::deque<std::vector<std::size_t>> colors_; | ||
617 | }; | ||
618 | |||
619 |
6/42✓ Branch 4 taken 14 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 14 times.
✗ Branch 8 not taken.
✓ Branch 13 taken 14 times.
✗ Branch 14 not taken.
✓ Branch 16 taken 14 times.
✗ Branch 17 not taken.
✓ Branch 22 taken 14 times.
✗ Branch 23 not taken.
✓ Branch 25 taken 14 times.
✗ 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 37 not taken.
✗ Branch 38 not taken.
✗ Branch 40 not taken.
✗ Branch 41 not taken.
✗ Branch 43 not taken.
✗ Branch 44 not taken.
✗ Branch 46 not taken.
✗ Branch 47 not taken.
✗ Branch 49 not taken.
✗ Branch 50 not taken.
✗ Branch 52 not taken.
✗ Branch 53 not taken.
✗ Branch 55 not taken.
✗ Branch 56 not taken.
✗ Branch 58 not taken.
✗ Branch 59 not taken.
✗ Branch 61 not taken.
✗ Branch 62 not taken.
✗ Branch 64 not taken.
✗ Branch 65 not taken.
✗ Branch 67 not taken.
✗ Branch 68 not taken.
✗ Branch 70 not taken.
✗ Branch 71 not taken.
|
126 | DUMUX_REGISTER_PRECONDITIONER("par_mt_ssor", Dune::PreconditionerTag, Dune::defaultPreconditionerBlockLevelCreator<Dumux::ParMTSSOR>()); |
620 | |||
621 | } // end namespace Dumux | ||
622 | |||
623 | namespace Dune::Amg { | ||
624 | |||
625 | // make it possible to use Dumux::ParMTJac as AMG smoother | ||
626 | template<class M, class X, class Y, int l> | ||
627 | struct ConstructionTraits<Dumux::ParMTJac<M,X,Y,l> > | ||
628 | { | ||
629 | using Arguments = DefaultConstructionArgs<SeqJac<M,X,Y,l> >; | ||
630 | 4 | static inline std::shared_ptr<Dumux::ParMTJac<M,X,Y,l>> construct(Arguments& args) | |
631 | { | ||
632 | 4 | return std::make_shared<Dumux::ParMTJac<M,X,Y,l>>( | |
633 | 8 | args.getMatrix(), args.getArgs().iterations, args.getArgs().relaxationFactor | |
634 | 4 | ); | |
635 | } | ||
636 | }; | ||
637 | |||
638 | // make it possible to use Dumux::ParMTSOR as AMG smoother | ||
639 | template<class M, class X, class Y, int l> | ||
640 | struct ConstructionTraits<Dumux::ParMTSOR<M,X,Y,l> > | ||
641 | { | ||
642 | using Arguments = DefaultConstructionArgs<SeqSOR<M,X,Y,l> >; | ||
643 | 4 | static inline std::shared_ptr<Dumux::ParMTSOR<M,X,Y,l>> construct(Arguments& args) | |
644 | { | ||
645 | 4 | return std::make_shared<Dumux::ParMTSOR<M,X,Y,l>>( | |
646 | 8 | args.getMatrix(), args.getArgs().iterations, args.getArgs().relaxationFactor | |
647 | 4 | ); | |
648 | } | ||
649 | }; | ||
650 | |||
651 | template<class M, class X, class Y, int l> | ||
652 | struct SmootherApplier<Dumux::ParMTSOR<M,X,Y,l> > | ||
653 | { | ||
654 | typedef Dumux::ParMTSOR<M,X,Y,l> Smoother; | ||
655 | typedef typename Smoother::range_type Range; | ||
656 | typedef typename Smoother::domain_type Domain; | ||
657 | |||
658 | static void preSmooth(Smoother& smoother, Domain& v, Range& d) | ||
659 | 120 | { smoother.template apply<true>(v,d); } | |
660 | |||
661 | static void postSmooth(Smoother& smoother, Domain& v, Range& d) | ||
662 | 120 | { smoother.template apply<false>(v,d); } | |
663 | }; | ||
664 | |||
665 | // make it possible to use Dumux::ParMTSSOR as AMG smoother | ||
666 | template<class M, class X, class Y, int l> | ||
667 | struct ConstructionTraits<Dumux::ParMTSSOR<M,X,Y,l> > | ||
668 | { | ||
669 | using Arguments = DefaultConstructionArgs<SeqSSOR<M,X,Y,l> >; | ||
670 | 4 | static inline std::shared_ptr<Dumux::ParMTSSOR<M,X,Y,l>> construct(Arguments& args) | |
671 | { | ||
672 | 4 | return std::make_shared<Dumux::ParMTSSOR<M,X,Y,l>>( | |
673 | 8 | args.getMatrix(), args.getArgs().iterations, args.getArgs().relaxationFactor | |
674 | 4 | ); | |
675 | } | ||
676 | }; | ||
677 | |||
678 | } // end namespace Dune::Amg | ||
679 | |||
680 | #endif | ||
681 |