version 3.8
stokes_solver.hh
Go to the documentation of this file.
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//
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>
38
39namespace Dumux::Detail {
40
58template<class M, class X, class Y, int l = 2>
59class StokesPreconditioner : public Dune::Preconditioner<X,Y>
60{
61 static_assert(Dumux::isMultiTypeBlockMatrix<M>::value && M::M() == 2 && M::N() == 2, "Expects a 2x2 MultiTypeBlockMatrix.");
62 static_assert(l == 2, "StokesPreconditioner expects a block level of 2 for Matrix type M.");
63
64 using A = std::decay_t<decltype(std::declval<M>()[Dune::Indices::_0][Dune::Indices::_0])>;
65 using U = std::decay_t<decltype(std::declval<X>()[Dune::Indices::_0])>;
66
67 using P = std::decay_t<decltype(std::declval<M>()[Dune::Indices::_1][Dune::Indices::_1])>;
68 using V = std::decay_t<decltype(std::declval<X>()[Dune::Indices::_1])>;
69
70 enum class Mode { symmetric, triangular, diagonal };
71
72public:
74 using matrix_type = M;
76 using domain_type = X;
78 using range_type = Y;
80 using field_type = typename X::field_type;
82 using scalar_field_type = Dune::Simd::Scalar<field_type>;
84 using PressureLinearOperator = Dune::MatrixAdapter<P,V,V>;
85
93 const std::shared_ptr<const Dune::AssembledLinearOperator<M,X,Y>>& fullLinearOperator,
94 const std::shared_ptr<const PressureLinearOperator>& pressureLinearOperator,
95 const Dune::ParameterTree& params
96 )
97 : matrix_(fullLinearOperator->getmat())
98 , pmatrix_(pressureLinearOperator->getmat())
99 , verbosity_(params.get<int>("verbosity"))
100 , paramGroup_(params.get<std::string>("ParameterGroup"))
101 {
102 const auto mode = getParamFromGroup<std::string>(
103 paramGroup_, "LinearSolver.Preconditioner.Mode", "Diagonal"
104 );
105
106 if (mode == "Symmetric")
107 mode_ = Mode::symmetric;
108 else if (mode == "Triangular")
109 mode_ = Mode::triangular;
110 else
111 mode_ = Mode::diagonal;
112
113 initPreconditioner_(params);
114 }
115
119 void pre(X& update, Y& currentDefect) override {}
120
132 void apply(X& update, const Y& currentDefect) override
133 {
134 using namespace Dune::Indices;
135
136 if (mode_ == Mode::symmetric)
137 {
138 update = applyRightBlock_(currentDefect);
139 update = applyDiagBlock_(update);
140 update = applyLeftBlock_(update);
141 }
142 else if (mode_ == Mode::triangular)
143 {
144 update = applyRightBlock_(currentDefect);
145 update = applyDiagBlock_(update);
146 }
147 else
148 update = applyDiagBlock_(currentDefect);
149 }
150
154 void post(X& update) override {}
155
157 Dune::SolverCategory::Category category() const override
158 {
159 return Dune::SolverCategory::sequential;
160 }
161
162private:
163 X applyRightBlock_(const Y& d)
164 {
165 using namespace Dune::Indices;
166
167 // right bit of LDU decomposition
168 // applied to d
169 //
170 // | I 0 |
171 // | -C*inv(A) I |
172 //
173
174 auto dTmp0 = d[_0];
175 auto vTmp = d; vTmp = 0.0;
176
177 // invert velocity block (apply to first block of d)
178 applyPreconditionerForA_(vTmp[_0], dTmp0);
179
180 // then multiply with C
181 matrix_[_1][_0].mv(vTmp[_0], vTmp[_1]);
182
183 // and subtract from d
184 auto v = d;
185 v[_0] = vTmp[_0]; // already do A^-1 d of the diagonal block because we already computed it here
186 v[_1] -= vTmp[_1];
187 return v;
188 }
189
190 X applyDiagBlock_(const Y& d)
191 {
192 using namespace Dune::Indices;
193
194 // diagonal middle bit of LDU decomposition
195 auto dTmp = d;
196 auto v = d; v = 0.0;
197
198 // invert velocity block
199 if (mode_ == Mode::diagonal)
200 applyPreconditionerForA_(v[_0], dTmp[_0]);
201
202 // reuse the already computed A^-1 d (see applyRightBlock_)
203 else
204 v[_0] = dTmp[_0];
205
206 // invert pressure block
207 applyPreconditionerForP_(v[_1], dTmp[_1]);
208
209 return v;
210 }
211
212 X applyLeftBlock_(const Y& d)
213 {
214 using namespace Dune::Indices;
215
216 // left bit of LDU decomposition
217 // applied to d
218 //
219 // | I -inv(A)*B |
220 // | 0 I |
221 //
222
223 auto dTmp = d;
224 auto vTmp = d; vTmp = 0.0;
225
226 // multiply with B
227 matrix_[_0][_1].umv(dTmp[_1], vTmp[_0]);
228
229 // invert velocity block (apply to first block of d)
230 auto vTmp0 = vTmp[_0]; vTmp0 = 0.0;
231 applyPreconditionerForA_(vTmp0, vTmp[_0]);
232
233 // and subtract from d
234 auto v = d;
235 v[_0] -= vTmp0;
236
237 return v;
238 }
239
240 void initPreconditioner_(const Dune::ParameterTree& params)
241 {
242 using namespace Dune::Indices;
243
244 if (getParamFromGroup<bool>(paramGroup_, "LinearSolver.DirectSolverForVelocity", false))
245 {
246#if HAVE_UMFPACK
247 directSolver_ = std::make_shared<Dune::UMFPack<A>>(matrix_[_0][_0], verbosity_);
248 using Wrap = Dune::InverseOperator2Preconditioner<Dune::InverseOperator<U, U>>;
249 preconditionerForA_ = std::make_shared<Wrap>(*directSolver_);
250#else
251 DUNE_THROW(Dune::InvalidStateException, "Selected direct solver but UMFPack is not available.");
252#endif
253 }
254 else
255 {
256 using VelLinearOperator = Dune::MatrixAdapter<A, U, U>;
257 auto lopV = std::make_shared<VelLinearOperator>(matrix_[_0][_0]);
258 preconditionerForA_ = std::make_shared<
259 Dune::Amg::AMG<VelLinearOperator, U, Dumux::ParMTSSOR<A,U,U>>
260 >(lopV, params);
261 }
262
263 using PressJacobi = Dumux::ParMTJac<P, V, V>;
264 preconditionerForP_ = std::make_shared<PressJacobi>(pmatrix_, 1, 1.0);
265 }
266
267 template<class Sol, class Rhs>
268 void applyPreconditionerForA_(Sol& sol, Rhs& rhs) const
269 {
270 preconditionerForA_->pre(sol, rhs);
271 preconditionerForA_->apply(sol, rhs);
272 preconditionerForA_->post(sol);
273 }
274
275 template<class Sol, class Rhs>
276 void applyPreconditionerForP_(Sol& sol, Rhs& rhs) const
277 {
278 preconditionerForP_->pre(sol, rhs);
279 preconditionerForP_->apply(sol, rhs);
280 preconditionerForP_->post(sol);
281 }
282
284 const M& matrix_;
286 const P& pmatrix_;
288 const int verbosity_;
289
290 std::shared_ptr<Dune::Preconditioner<U, U>> preconditionerForA_;
291 std::shared_ptr<Dune::Preconditioner<V, V>> preconditionerForP_;
292 std::shared_ptr<Dune::InverseOperator<U, U>> directSolver_;
293
294 const std::string paramGroup_;
295 Mode mode_;
296};
297
298} // end namespace Detail
299
300namespace Dumux {
301
308template<class Matrix, class Vector, class VelocityGG, class PressureGG>
310: public LinearSolver
311{
313public:
322 StokesSolver(std::shared_ptr<const VelocityGG> vGridGeometry,
323 std::shared_ptr<const PressureGG> pGridGeometry,
324 const Vector& dirichletDofs,
325 const std::string& paramGroup = "")
327 , vGridGeometry_(std::move(vGridGeometry))
328 , pGridGeometry_(std::move(pGridGeometry))
329 , dirichletDofs_(dirichletDofs)
330 {
331 params_ = LinearSolverParameters<LinearSolverTraits<VelocityGG>>::createParameterTree(this->paramGroup());
332 density_ = getParamFromGroup<double>(this->paramGroup(), "Component.LiquidDensity");
333 viscosity_ = getParamFromGroup<double>(this->paramGroup(), "Component.LiquidDynamicViscosity");
334 weight_ = getParamFromGroup<double>(this->paramGroup(), "LinearSolver.Preconditioner.MassMatrixWeight", 1.0);
335 solverType_ = getParamFromGroup<std::string>(this->paramGroup(), "LinearSolver.Type", "gmres");
336 scalarProduct_ = std::make_shared<Dune::ScalarProduct<Vector>>();
337 }
338
339 bool solve(const Matrix& A, Vector& x, const Vector& b)
340 {
341 auto bTmp = b;
342 auto ATmp = A;
343
344 return applyIterativeSolver_(ATmp, x, bTmp);
345 }
346
347 Scalar norm(const Vector& b) const
348 {
349 return scalarProduct_->norm(b);
350 }
351
352 std::string name() const
353 {
354 return "Block-preconditioned Stokes solver";
355 }
356
357 const Dune::InverseOperatorResult& result() const
358 {
359 return result_;
360 }
361
362private:
363 bool applyIterativeSolver_(Matrix& A, Vector& x, Vector& b)
364 {
365 // make Dirichlet boundary conditions symmetric
366 if (getParamFromGroup<bool>(this->paramGroup(), "LinearSolver.SymmetrizeDirichlet", true))
367 symmetrizeConstraints(A, b, dirichletDofs_);
368
369 // make Matrix symmetric on the block-scale
370 using namespace Dune::Indices;
371 A[_1] *= -1.0/density_;
372 b[_1] *= -1.0/density_;
373
374 auto op = std::make_shared<Dumux::ParallelMultiTypeMatrixAdapter<Matrix, Vector, Vector>>(A);
375 auto pop = makePressureLinearOperator_<typename Preconditioner::PressureLinearOperator>();
376 auto preconditioner = std::make_shared<Preconditioner>(op, pop, params_.sub("preconditioner"));
377 params_["verbose"] = pGridGeometry_->gridView().comm().rank() == 0 ? params_["verbose"] : "0";
378
379 // defaults to restarted GMRes
380 std::unique_ptr<Dune::InverseOperator<Vector, Vector>> solver;
381 if (solverType_ == "minres")
382 solver = std::make_unique<Dune::MINRESSolver<Vector>>(op, scalarProduct_, preconditioner, params_);
383 else if (solverType_ == "bicgstab")
384 solver = std::make_unique<Dune::BiCGSTABSolver<Vector>>(op, scalarProduct_, preconditioner, params_);
385 else if (solverType_ == "gmres")
386 solver = std::make_unique<Dune::RestartedGMResSolver<Vector>>(op, scalarProduct_, preconditioner, params_);
387 else
388 DUNE_THROW(Dune::NotImplemented, "Solver choice " << solverType_ << " is not implemented");
389
390 solver->apply(x, b, result_);
391
392 return result_.converged;
393 }
394
395 template<class LinearOperator>
396 std::shared_ptr<LinearOperator> makePressureLinearOperator_()
397 {
398 using M = typename LinearOperator::matrix_type;
399 auto massMatrix = createMassMatrix_<M>();
400 return std::make_shared<LinearOperator>(massMatrix);
401 }
402
403 template<class M>
404 std::shared_ptr<M> createMassMatrix_()
405 {
406 auto massMatrix = std::make_shared<M>();
407 massMatrix->setBuildMode(M::random);
408 const auto numDofs = pGridGeometry_->numDofs();
409
410 Dune::MatrixIndexSet pattern;
411 pattern.resize(numDofs, numDofs);
412 for (unsigned int globalI = 0; globalI < numDofs; ++globalI)
413 pattern.add(globalI, globalI);
414 pattern.exportIdx(*massMatrix);
415
416 const auto& gv = pGridGeometry_->gridView();
417 auto fvGeometry = localView(*pGridGeometry_);
418 for (const auto& element : elements(gv))
419 {
420 fvGeometry.bindElement(element);
421 for (const auto& scv : scvs(fvGeometry))
422 {
423 using Extrusion = Extrusion_t<PressureGG>;
424 const auto dofIndex = scv.dofIndex();
425 if (element.partitionType() == Dune::GhostEntity) // do not modify ghosts
426 (*massMatrix)[dofIndex][dofIndex] = 1.0;
427 else
428 (*massMatrix)[dofIndex][dofIndex] += weight_*Extrusion::volume(fvGeometry, scv)/(2.0*viscosity_);
429 }
430 }
431
432 return massMatrix;
433 }
434
435 double density_, viscosity_, weight_;
436 Dune::InverseOperatorResult result_;
437 Dune::ParameterTree params_;
438 std::shared_ptr<const VelocityGG> vGridGeometry_;
439 std::shared_ptr<const PressureGG> pGridGeometry_;
440 const Vector& dirichletDofs_;
441 std::string solverType_;
442
443 std::shared_ptr<Dune::ScalarProduct<Vector>> scalarProduct_;
444};
445
446} // end namespace Dumux
447
448#endif
A Stokes preconditioner (saddle-point problem) for the problem .
Definition: stokes_solver.hh:60
Dune::SolverCategory::Category category() const override
Category of the preconditioner (see SolverCategory::Category)
Definition: stokes_solver.hh:157
X domain_type
The domain type of the preconditioner.
Definition: stokes_solver.hh:76
void post(X &update) override
Clean up.
Definition: stokes_solver.hh:154
typename X::field_type field_type
The field type of the preconditioner.
Definition: stokes_solver.hh:80
Dune::Simd::Scalar< field_type > scalar_field_type
Scalar type underlying the field_type.
Definition: stokes_solver.hh:82
M matrix_type
The matrix type the preconditioner is for.
Definition: stokes_solver.hh:74
Dune::MatrixAdapter< P, V, V > PressureLinearOperator
the type of the pressure operator
Definition: stokes_solver.hh:84
Y range_type
The range type of the preconditioner.
Definition: stokes_solver.hh:78
void pre(X &update, Y &currentDefect) override
Prepare the preconditioner.
Definition: stokes_solver.hh:119
StokesPreconditioner(const std::shared_ptr< const Dune::AssembledLinearOperator< M, X, Y > > &fullLinearOperator, const std::shared_ptr< const PressureLinearOperator > &pressureLinearOperator, const Dune::ParameterTree &params)
Constructor.
Definition: stokes_solver.hh:92
void apply(X &update, const Y &currentDefect) override
Apply the preconditioner.
Definition: stokes_solver.hh:132
Base class for linear solvers.
Definition: solver.hh:27
const std::string & paramGroup() const
the parameter group for getting parameter from the parameter tree
Definition: solver.hh:78
double Scalar
Definition: solver.hh:31
Generates a parameter tree required for the linear solvers and precondioners of the Dune ISTL.
Definition: linearsolverparameters.hh:40
Multi-threaded Jacobi preconditioner.
Definition: preconditioners.hh:331
Preconditioned iterative solver for the incompressible Stokes problem.
Definition: stokes_solver.hh:311
const Dune::InverseOperatorResult & result() const
Definition: stokes_solver.hh:357
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: stokes_solver.hh:339
StokesSolver(std::shared_ptr< const VelocityGG > vGridGeometry, std::shared_ptr< const PressureGG > pGridGeometry, const Vector &dirichletDofs, const std::string &paramGroup="")
Constructor.
Definition: stokes_solver.hh:322
std::string name() const
Definition: stokes_solver.hh:352
Scalar norm(const Vector &b) const
Definition: stokes_solver.hh:347
Helper classes to compute the integration elements.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
auto volume(const Geometry &geo, unsigned int integrationOrder=4)
The volume of a given geometry.
Definition: volume.hh:159
void symmetrizeConstraints(Dune::BCRSMatrix< MBlock > &A, Vector &b, const Vector &constrainedRows)
Symmetrize the constrained system Ax=b.
Definition: symmetrize_constraints.hh:117
Generates a parameter tree required for the linear solvers and precondioners of the Dune ISTL.
Define traits for linear solvers.
Define some often used mathematical functions.
Distance implementation details.
Definition: cvfelocalresidual.hh:25
Definition: adapt.hh:17
Provides a helper class for nonoverlapping decomposition.
A parallel version of a linear operator.
Dumux preconditioners for iterative solvers.
Base class for linear solvers.
Helper type to determine whether a given type is a Dune::MultiTypeBlockMatrix.
Definition: matrix.hh:37
Helper function to symmetrize row-constraints in constrained linear systems.