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 Preconditioned iterative solver for the incompressible Stokes problem | ||
11 | */ | ||
12 | #ifndef DUMUX_LINEAR_STOKES_SOLVER_HH | ||
13 | #define DUMUX_LINEAR_STOKES_SOLVER_HH | ||
14 | |||
15 | #include <type_traits> | ||
16 | #include <memory> | ||
17 | #include <tuple> | ||
18 | |||
19 | #include <dune/common/parametertree.hh> | ||
20 | #include <dune/common/hybridutilities.hh> | ||
21 | #include <dune/common/exceptions.hh> | ||
22 | |||
23 | #include <dune/istl/matrixindexset.hh> | ||
24 | #include <dune/istl/preconditioner.hh> | ||
25 | #include <dune/istl/preconditioners.hh> | ||
26 | #include <dune/istl/solvers.hh> | ||
27 | #include <dune/istl/paamg/amg.hh> | ||
28 | |||
29 | #include <dumux/common/math.hh> | ||
30 | #include <dumux/linear/solver.hh> | ||
31 | #include <dumux/linear/preconditioners.hh> | ||
32 | #include <dumux/linear/linearsolverparameters.hh> | ||
33 | #include <dumux/linear/parallelhelpers.hh> | ||
34 | #include <dumux/linear/linearsolvertraits.hh> | ||
35 | #include <dumux/linear/parallelmatrixadapter.hh> | ||
36 | #include <dumux/discretization/extrusion.hh> | ||
37 | #include <dumux/linear/symmetrize_constraints.hh> | ||
38 | #include <dumux/assembly/jacobianpattern.hh> | ||
39 | |||
40 | namespace Dumux::Detail { | ||
41 | |||
42 | /*! | ||
43 | * \ingroup Linear | ||
44 | * \brief A Stokes preconditioner (saddle-point problem) for the problem | ||
45 | * \f$ | ||
46 | \begin{pmatrix} A & B \\ C & 0 \end{pmatrix} | ||
47 | \begin{pmatrix} u \\ p \end{pmatrix} = | ||
48 | \begin{pmatrix} f \\ g \end{pmatrix}, | ||
49 | * \f$ | ||
50 | * | ||
51 | * where A is the discrete velocity operator and B and C are discrete gradient and divergence operators. | ||
52 | * This preconditioner is especially suited for solving the incompressible Stokes equations. | ||
53 | * | ||
54 | * \tparam M Type of the matrix. | ||
55 | * \tparam X Type of the update. | ||
56 | * \tparam Y Type of the defect. | ||
57 | * \tparam l Preconditioner block level | ||
58 | */ | ||
59 | template<class M, class X, class Y, int l = 2> | ||
60 | class StokesPreconditioner : public Dune::Preconditioner<X,Y> | ||
61 | { | ||
62 | static_assert(Dumux::isMultiTypeBlockMatrix<M>::value && M::M() == 2 && M::N() == 2, "Expects a 2x2 MultiTypeBlockMatrix."); | ||
63 | static_assert(l == 2, "StokesPreconditioner expects a block level of 2 for Matrix type M."); | ||
64 | |||
65 | using A = std::decay_t<decltype(std::declval<M>()[Dune::Indices::_0][Dune::Indices::_0])>; | ||
66 | using U = std::decay_t<decltype(std::declval<X>()[Dune::Indices::_0])>; | ||
67 | |||
68 | using P = std::decay_t<decltype(std::declval<M>()[Dune::Indices::_1][Dune::Indices::_1])>; | ||
69 | using V = std::decay_t<decltype(std::declval<X>()[Dune::Indices::_1])>; | ||
70 | |||
71 | enum class Mode { symmetric, triangular, diagonal }; | ||
72 | |||
73 | public: | ||
74 | //! \brief The matrix type the preconditioner is for. | ||
75 | using matrix_type = M; | ||
76 | //! \brief The domain type of the preconditioner. | ||
77 | using domain_type = X; | ||
78 | //! \brief The range type of the preconditioner. | ||
79 | using range_type = Y; | ||
80 | //! \brief The field type of the preconditioner. | ||
81 | using field_type = typename X::field_type; | ||
82 | //! \brief Scalar type underlying the field_type. | ||
83 | using scalar_field_type = Dune::Simd::Scalar<field_type>; | ||
84 | //! \brief the type of the pressure operator | ||
85 | using PressureLinearOperator = Dune::MatrixAdapter<P,V,V>; | ||
86 | |||
87 | /*! | ||
88 | * \brief Constructor | ||
89 | * \param fullLinearOperator the Stokes linear operator | ||
90 | * \param pressureLinearOperator the linear operator to be used for the pressure space preconditioner | ||
91 | * \param params a parameter tree for the preconditioner configuration | ||
92 | */ | ||
93 | ✗ | StokesPreconditioner( | |
94 | const std::shared_ptr<const Dune::AssembledLinearOperator<M,X,Y>>& fullLinearOperator, | ||
95 | const std::shared_ptr<const PressureLinearOperator>& pressureLinearOperator, | ||
96 | const Dune::ParameterTree& params | ||
97 | ) | ||
98 | ✗ | : matrix_(fullLinearOperator->getmat()) | |
99 | ✗ | , pmatrix_(pressureLinearOperator->getmat()) | |
100 | ✗ | , verbosity_(params.get<int>("verbosity")) | |
101 | ✗ | , paramGroup_(params.get<std::string>("ParameterGroup")) | |
102 | { | ||
103 | ✗ | const auto mode = getParamFromGroup<std::string>( | |
104 | ✗ | paramGroup_, "LinearSolver.Preconditioner.Mode", "Diagonal" | |
105 | ); | ||
106 | |||
107 | ✗ | if (mode == "Symmetric") | |
108 | ✗ | mode_ = Mode::symmetric; | |
109 | ✗ | else if (mode == "Triangular") | |
110 | ✗ | mode_ = Mode::triangular; | |
111 | else | ||
112 | ✗ | mode_ = Mode::diagonal; | |
113 | |||
114 | ✗ | initPreconditioner_(params); | |
115 | ✗ | } | |
116 | |||
117 | /*! | ||
118 | * \brief Prepare the preconditioner. | ||
119 | */ | ||
120 | ✗ | void pre(X& update, Y& currentDefect) override {} | |
121 | |||
122 | /*! | ||
123 | * \brief Apply the preconditioner | ||
124 | * | ||
125 | * \param update The update to be computed. | ||
126 | * \param currentDefect The current defect. | ||
127 | * | ||
128 | * The currentDefect has be be in a consistent representation, | ||
129 | * Definition 2.3 Blatt and Bastian (2009) https://doi.org/10.1504/IJCSE.2008.021112 | ||
130 | * The update is initially zero. At exit the update has to be | ||
131 | * in a consistent representation. This usually requires communication. | ||
132 | */ | ||
133 | ✗ | void apply(X& update, const Y& currentDefect) override | |
134 | { | ||
135 | using namespace Dune::Indices; | ||
136 | |||
137 | ✗ | if (mode_ == Mode::symmetric) | |
138 | { | ||
139 | ✗ | update = applyRightBlock_(currentDefect); | |
140 | ✗ | update = applyDiagBlock_(update); | |
141 | ✗ | update = applyLeftBlock_(update); | |
142 | } | ||
143 | ✗ | else if (mode_ == Mode::triangular) | |
144 | { | ||
145 | ✗ | update = applyRightBlock_(currentDefect); | |
146 | ✗ | update = applyDiagBlock_(update); | |
147 | } | ||
148 | else | ||
149 | ✗ | update = applyDiagBlock_(currentDefect); | |
150 | ✗ | } | |
151 | |||
152 | /*! | ||
153 | * \brief Clean up. | ||
154 | */ | ||
155 | ✗ | void post(X& update) override {} | |
156 | |||
157 | //! Category of the preconditioner (see SolverCategory::Category) | ||
158 | ✗ | Dune::SolverCategory::Category category() const override | |
159 | { | ||
160 | ✗ | return Dune::SolverCategory::sequential; | |
161 | } | ||
162 | |||
163 | private: | ||
164 | ✗ | X applyRightBlock_(const Y& d) | |
165 | { | ||
166 | using namespace Dune::Indices; | ||
167 | |||
168 | // right bit of LDU decomposition | ||
169 | // applied to d | ||
170 | // | ||
171 | // | I 0 | | ||
172 | // | -C*inv(A) I | | ||
173 | // | ||
174 | |||
175 | ✗ | auto dTmp0 = d[_0]; | |
176 | ✗ | auto vTmp = d; vTmp = 0.0; | |
177 | |||
178 | // invert velocity block (apply to first block of d) | ||
179 | ✗ | applyPreconditionerForA_(vTmp[_0], dTmp0); | |
180 | |||
181 | // then multiply with C | ||
182 | ✗ | matrix_[_1][_0].mv(vTmp[_0], vTmp[_1]); | |
183 | |||
184 | // and subtract from d | ||
185 | ✗ | auto v = d; | |
186 | ✗ | v[_0] = vTmp[_0]; // already do A^-1 d of the diagonal block because we already computed it here | |
187 | ✗ | v[_1] -= vTmp[_1]; | |
188 | ✗ | return v; | |
189 | } | ||
190 | |||
191 | ✗ | X applyDiagBlock_(const Y& d) | |
192 | { | ||
193 | using namespace Dune::Indices; | ||
194 | |||
195 | // diagonal middle bit of LDU decomposition | ||
196 | ✗ | auto dTmp = d; | |
197 | ✗ | auto v = d; v = 0.0; | |
198 | |||
199 | // invert velocity block | ||
200 | ✗ | if (mode_ == Mode::diagonal) | |
201 | ✗ | applyPreconditionerForA_(v[_0], dTmp[_0]); | |
202 | |||
203 | // reuse the already computed A^-1 d (see applyRightBlock_) | ||
204 | else | ||
205 | ✗ | v[_0] = dTmp[_0]; | |
206 | |||
207 | // invert pressure block | ||
208 | ✗ | applyPreconditionerForP_(v[_1], dTmp[_1]); | |
209 | |||
210 | ✗ | return v; | |
211 | } | ||
212 | |||
213 | ✗ | X applyLeftBlock_(const Y& d) | |
214 | { | ||
215 | using namespace Dune::Indices; | ||
216 | |||
217 | // left bit of LDU decomposition | ||
218 | // applied to d | ||
219 | // | ||
220 | // | I -inv(A)*B | | ||
221 | // | 0 I | | ||
222 | // | ||
223 | |||
224 | ✗ | auto dTmp = d; | |
225 | ✗ | auto vTmp = d; vTmp = 0.0; | |
226 | |||
227 | // multiply with B | ||
228 | ✗ | matrix_[_0][_1].umv(dTmp[_1], vTmp[_0]); | |
229 | |||
230 | // invert velocity block (apply to first block of d) | ||
231 | ✗ | auto vTmp0 = vTmp[_0]; vTmp0 = 0.0; | |
232 | ✗ | applyPreconditionerForA_(vTmp0, vTmp[_0]); | |
233 | |||
234 | // and subtract from d | ||
235 | ✗ | auto v = d; | |
236 | ✗ | v[_0] -= vTmp0; | |
237 | |||
238 | ✗ | return v; | |
239 | } | ||
240 | |||
241 | ✗ | void initPreconditioner_(const Dune::ParameterTree& params) | |
242 | { | ||
243 | using namespace Dune::Indices; | ||
244 | |||
245 | ✗ | if (getParamFromGroup<bool>(paramGroup_, "LinearSolver.DirectSolverForVelocity", false)) | |
246 | { | ||
247 | #if HAVE_UMFPACK | ||
248 | ✗ | directSolverForA_ = std::make_shared<Dune::UMFPack<A>>(matrix_[_0][_0], verbosity_); | |
249 | using Wrap = Dune::InverseOperator2Preconditioner<Dune::InverseOperator<U, U>>; | ||
250 | ✗ | preconditionerForA_ = std::make_shared<Wrap>(*directSolverForA_); | |
251 | #else | ||
252 | DUNE_THROW(Dune::InvalidStateException, "Selected direct solver but UMFPack is not available."); | ||
253 | #endif | ||
254 | } | ||
255 | else | ||
256 | { | ||
257 | using VelLinearOperator = Dune::MatrixAdapter<A, U, U>; | ||
258 | ✗ | auto lopV = std::make_shared<VelLinearOperator>(matrix_[_0][_0]); | |
259 | ✗ | preconditionerForA_ = std::make_shared< | |
260 | Dune::Amg::AMG<VelLinearOperator, U, Dumux::ParMTSSOR<A,U,U>> | ||
261 | >(lopV, params); | ||
262 | } | ||
263 | |||
264 | ✗ | if (getParamFromGroup<bool>(paramGroup_, "LinearSolver.DirectSolverForPressure", false)) | |
265 | { | ||
266 | #if HAVE_UMFPACK | ||
267 | ✗ | directSolverForP_ = std::make_shared<Dune::UMFPack<P>>(pmatrix_, verbosity_); | |
268 | using Wrap = Dune::InverseOperator2Preconditioner<Dune::InverseOperator<V, V>>; | ||
269 | ✗ | preconditionerForP_ = std::make_shared<Wrap>(*directSolverForP_); | |
270 | #else | ||
271 | DUNE_THROW(Dune::InvalidStateException, "Selected direct solver but UMFPack is not available."); | ||
272 | #endif | ||
273 | } | ||
274 | else | ||
275 | { | ||
276 | using PressJacobi = Dumux::ParMTJac<P, V, V>; | ||
277 | ✗ | std::size_t numIterations = pmatrix_.nonzeroes() == pmatrix_.N() ? 1 : 10; | |
278 | ✗ | preconditionerForP_ = std::make_shared<PressJacobi>(pmatrix_, numIterations, 1.0); | |
279 | } | ||
280 | ✗ | } | |
281 | |||
282 | template<class Sol, class Rhs> | ||
283 | ✗ | void applyPreconditionerForA_(Sol& sol, Rhs& rhs) const | |
284 | { | ||
285 | ✗ | preconditionerForA_->pre(sol, rhs); | |
286 | ✗ | preconditionerForA_->apply(sol, rhs); | |
287 | ✗ | preconditionerForA_->post(sol); | |
288 | ✗ | } | |
289 | |||
290 | template<class Sol, class Rhs> | ||
291 | ✗ | void applyPreconditionerForP_(Sol& sol, Rhs& rhs) const | |
292 | { | ||
293 | ✗ | preconditionerForP_->pre(sol, rhs); | |
294 | ✗ | preconditionerForP_->apply(sol, rhs); | |
295 | ✗ | preconditionerForP_->post(sol); | |
296 | ✗ | } | |
297 | |||
298 | //! \brief The matrix we operate on. | ||
299 | const M& matrix_; | ||
300 | //! \brief The matrix we operate on. | ||
301 | const P& pmatrix_; | ||
302 | //! \brief The verbosity level | ||
303 | const int verbosity_; | ||
304 | |||
305 | std::shared_ptr<Dune::Preconditioner<U, U>> preconditionerForA_; | ||
306 | std::shared_ptr<Dune::Preconditioner<V, V>> preconditionerForP_; | ||
307 | std::shared_ptr<Dune::InverseOperator<U, U>> directSolverForA_; | ||
308 | std::shared_ptr<Dune::InverseOperator<V, V>> directSolverForP_; | ||
309 | |||
310 | const std::string paramGroup_; | ||
311 | Mode mode_; | ||
312 | }; | ||
313 | |||
314 | } // end namespace Detail | ||
315 | |||
316 | namespace Dumux { | ||
317 | |||
318 | /*! | ||
319 | * \ingroup Linear | ||
320 | * \brief Preconditioned iterative solver for the incompressible Stokes problem | ||
321 | * \note Uses StokesPreconditioner as preconditioner (tailored to the incompressible Stokes problem) | ||
322 | * \note No MPI parallelization implemented, some shared-memory parallelism is enabled | ||
323 | */ | ||
324 | template<class Matrix, class Vector, class VelocityGG, class PressureGG> | ||
325 | class StokesSolver | ||
326 | : public LinearSolver | ||
327 | { | ||
328 | using Preconditioner = Detail::StokesPreconditioner<Matrix, Vector, Vector>; | ||
329 | public: | ||
330 | /*! | ||
331 | * \brief Constructor | ||
332 | * \param vGridGeometry grid geometry of the velocity discretization | ||
333 | * \param pGridGeometry grid geometry of the pressure discretization | ||
334 | * \param dirichletDofs a vector (same size and shape as right hand side) where the dirichlet dofs are marked with 1.0 | ||
335 | * and all other entries are 0.0. | ||
336 | * \param paramGroup group prefix when looking up keys in the parameter tree | ||
337 | */ | ||
338 | 5 | StokesSolver(std::shared_ptr<const VelocityGG> vGridGeometry, | |
339 | std::shared_ptr<const PressureGG> pGridGeometry, | ||
340 | const Vector& dirichletDofs, | ||
341 | const std::string& paramGroup = "") | ||
342 | : LinearSolver(paramGroup) | ||
343 | 5 | , vGridGeometry_(std::move(vGridGeometry)) | |
344 | 5 | , pGridGeometry_(std::move(pGridGeometry)) | |
345 |
6/20✓ Branch 2 taken 5 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 5 times.
✗ Branch 6 not taken.
✓ Branch 8 taken 5 times.
✗ Branch 9 not taken.
✓ Branch 11 taken 5 times.
✗ Branch 12 not taken.
✓ Branch 14 taken 5 times.
✗ Branch 15 not taken.
✓ Branch 17 taken 5 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✗ Branch 26 not taken.
✗ Branch 27 not taken.
|
5 | , dirichletDofs_(dirichletDofs) |
346 | { | ||
347 |
2/4✓ Branch 1 taken 5 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 5 times.
✗ Branch 5 not taken.
|
10 | params_ = LinearSolverParameters<LinearSolverTraits<VelocityGG>>::createParameterTree(this->paramGroup()); |
348 |
2/4✓ Branch 1 taken 5 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 5 times.
✗ Branch 5 not taken.
|
10 | density_ = getParamFromGroup<double>(this->paramGroup(), "Component.LiquidDensity"); |
349 |
2/4✓ Branch 1 taken 5 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 5 times.
✗ Branch 5 not taken.
|
10 | viscosity_ = getParamFromGroup<double>(this->paramGroup(), "Component.LiquidDynamicViscosity"); |
350 |
2/6✓ Branch 1 taken 5 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 5 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
10 | weight_ = getParamFromGroup<double>(this->paramGroup(), "LinearSolver.Preconditioner.MassMatrixWeight", 1.0); |
351 |
3/6✓ Branch 1 taken 5 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 5 times.
✗ Branch 5 not taken.
✗ Branch 7 not taken.
✓ Branch 8 taken 5 times.
|
10 | solverType_ = getParamFromGroup<std::string>(this->paramGroup(), "LinearSolver.Type", "gmres"); |
352 |
2/4✓ Branch 1 taken 5 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 5 times.
|
5 | scalarProduct_ = std::make_shared<Dune::ScalarProduct<Vector>>(); |
353 | 5 | } | |
354 | |||
355 | ✗ | bool solve(const Matrix& A, Vector& x, const Vector& b) | |
356 | { | ||
357 | ✗ | auto bTmp = b; | |
358 | ✗ | auto ATmp = A; | |
359 | |||
360 | ✗ | return applyIterativeSolver_(ATmp, x, bTmp); | |
361 | } | ||
362 | |||
363 | Scalar norm(const Vector& b) const | ||
364 | { | ||
365 | ✗ | return scalarProduct_->norm(b); | |
366 | } | ||
367 | |||
368 | std::string name() const | ||
369 | { | ||
370 | return "Block-preconditioned Stokes solver"; | ||
371 | } | ||
372 | |||
373 | const Dune::InverseOperatorResult& result() const | ||
374 | { | ||
375 | return result_; | ||
376 | } | ||
377 | |||
378 | private: | ||
379 | ✗ | bool applyIterativeSolver_(Matrix& A, Vector& x, Vector& b) | |
380 | { | ||
381 | // make Dirichlet boundary conditions symmetric | ||
382 | ✗ | if (getParamFromGroup<bool>(this->paramGroup(), "LinearSolver.SymmetrizeDirichlet", true)) | |
383 | ✗ | symmetrizeConstraints(A, b, dirichletDofs_); | |
384 | |||
385 | // make Matrix symmetric on the block-scale | ||
386 | using namespace Dune::Indices; | ||
387 | ✗ | A[_1] *= -1.0/density_; | |
388 | ✗ | b[_1] *= -1.0/density_; | |
389 | |||
390 | ✗ | auto op = std::make_shared<Dumux::ParallelMultiTypeMatrixAdapter<Matrix, Vector, Vector>>(A); | |
391 | ✗ | auto pop = makePressureLinearOperator_<typename Preconditioner::PressureLinearOperator>(); | |
392 | ✗ | auto preconditioner = std::make_shared<Preconditioner>(op, pop, params_.sub("preconditioner")); | |
393 | ✗ | params_["verbose"] = pGridGeometry_->gridView().comm().rank() == 0 ? params_["verbose"] : "0"; | |
394 | |||
395 | // defaults to restarted GMRes | ||
396 | ✗ | std::unique_ptr<Dune::InverseOperator<Vector, Vector>> solver; | |
397 | ✗ | if (solverType_ == "minres") | |
398 | ✗ | solver = std::make_unique<Dune::MINRESSolver<Vector>>(op, scalarProduct_, preconditioner, params_); | |
399 | ✗ | else if (solverType_ == "bicgstab") | |
400 | ✗ | solver = std::make_unique<Dune::BiCGSTABSolver<Vector>>(op, scalarProduct_, preconditioner, params_); | |
401 | ✗ | else if (solverType_ == "gmres") | |
402 | ✗ | solver = std::make_unique<Dune::RestartedGMResSolver<Vector>>(op, scalarProduct_, preconditioner, params_); | |
403 | else | ||
404 | ✗ | DUNE_THROW(Dune::NotImplemented, "Solver choice " << solverType_ << " is not implemented"); | |
405 | |||
406 | ✗ | solver->apply(x, b, result_); | |
407 | |||
408 | ✗ | return result_.converged; | |
409 | } | ||
410 | |||
411 | template<class LinearOperator> | ||
412 | ✗ | std::shared_ptr<LinearOperator> makePressureLinearOperator_() | |
413 | { | ||
414 | using M = typename LinearOperator::matrix_type; | ||
415 | ✗ | auto massMatrix = createMassMatrix_<M>(); | |
416 | ✗ | return std::make_shared<LinearOperator>(massMatrix); | |
417 | } | ||
418 | |||
419 | template<class M, class DM = typename PressureGG::DiscretizationMethod> | ||
420 | ✗ | std::shared_ptr<M> createMassMatrix_() | |
421 | { | ||
422 | ✗ | auto massMatrix = std::make_shared<M>(); | |
423 | ✗ | massMatrix->setBuildMode(M::random); | |
424 | ✗ | const auto numDofs = pGridGeometry_->numDofs(); | |
425 | |||
426 | if constexpr (DM{} == DiscretizationMethods::cctpfa || DM{} == DiscretizationMethods::ccmpfa) | ||
427 | { | ||
428 | ✗ | Dune::MatrixIndexSet pattern; | |
429 | ✗ | pattern.resize(numDofs, numDofs); | |
430 | ✗ | for (unsigned int globalI = 0; globalI < numDofs; ++globalI) | |
431 | ✗ | pattern.add(globalI, globalI); | |
432 | ✗ | pattern.exportIdx(*massMatrix); | |
433 | |||
434 | ✗ | const auto& gv = pGridGeometry_->gridView(); | |
435 | ✗ | auto fvGeometry = localView(*pGridGeometry_); | |
436 | ✗ | for (const auto& element : elements(gv)) | |
437 | { | ||
438 | ✗ | fvGeometry.bindElement(element); | |
439 | ✗ | for (const auto& scv : scvs(fvGeometry)) | |
440 | { | ||
441 | using Extrusion = Extrusion_t<PressureGG>; | ||
442 | ✗ | const auto dofIndex = scv.dofIndex(); | |
443 | ✗ | if (element.partitionType() == Dune::GhostEntity) // do not modify ghosts | |
444 | ✗ | (*massMatrix)[dofIndex][dofIndex] = 1.0; | |
445 | else | ||
446 | ✗ | (*massMatrix)[dofIndex][dofIndex] += weight_*Extrusion::volume(fvGeometry, scv)/(2.0*viscosity_); | |
447 | } | ||
448 | } | ||
449 | |||
450 | ✗ | return massMatrix; | |
451 | } | ||
452 | else if constexpr (DM{} == DiscretizationMethods::box || DM{} == DiscretizationMethods::fcdiamond) | ||
453 | { | ||
454 | ✗ | getJacobianPattern<true>(*pGridGeometry_).exportIdx(*massMatrix); | |
455 | ✗ | (*massMatrix) = 0.0; | |
456 | |||
457 | ✗ | const auto& gv = pGridGeometry_->gridView(); | |
458 | ✗ | auto fvGeometry = localView(*pGridGeometry_); | |
459 | ✗ | std::vector<Dune::FieldVector<double, 1>> values; | |
460 | |||
461 | ✗ | for (const auto& element : elements(gv)) | |
462 | { | ||
463 | ✗ | const auto geometry = element.geometry(); | |
464 | ✗ | const auto& localBasis = pGridGeometry_->feCache().get(element.type()).localBasis(); | |
465 | ✗ | const auto numLocalDofs = localBasis.size(); | |
466 | ✗ | values.resize(numLocalDofs); | |
467 | |||
468 | ✗ | fvGeometry.bindElement(element); | |
469 | ✗ | for (const auto& scvJ : scvs(fvGeometry)) | |
470 | { | ||
471 | // Use mid-point rule (only works for linear ansatz functions) | ||
472 | ✗ | const auto globalJ = scvJ.dofIndex(); | |
473 | ✗ | const auto qWeightJ = PressureGG::Extrusion::volume(fvGeometry, scvJ); | |
474 | ✗ | const auto quadPos = geometry.local(scvJ.center()); | |
475 | ✗ | localBasis.evaluateFunction(quadPos, values); | |
476 | |||
477 | ✗ | for (const auto& scvI : scvs(fvGeometry)) | |
478 | { | ||
479 | ✗ | const auto valueIJ = values[scvI.localDofIndex()]*qWeightJ/(2.0*viscosity_); | |
480 | ✗ | (*massMatrix)[scvI.dofIndex()][globalJ][0][0] += valueIJ; | |
481 | } | ||
482 | } | ||
483 | } | ||
484 | |||
485 | ✗ | return massMatrix; | |
486 | } | ||
487 | |||
488 | DUNE_THROW(Dune::NotImplemented, "Mass matrix for discretization method not implemented"); | ||
489 | } | ||
490 | |||
491 | double density_, viscosity_, weight_; | ||
492 | Dune::InverseOperatorResult result_; | ||
493 | Dune::ParameterTree params_; | ||
494 | std::shared_ptr<const VelocityGG> vGridGeometry_; | ||
495 | std::shared_ptr<const PressureGG> pGridGeometry_; | ||
496 | const Vector& dirichletDofs_; | ||
497 | std::string solverType_; | ||
498 | |||
499 | std::shared_ptr<Dune::ScalarProduct<Vector>> scalarProduct_; | ||
500 | }; | ||
501 | |||
502 | } // end namespace Dumux | ||
503 | |||
504 | #endif | ||
505 |