24#ifndef DUMUX_LINEAR_PRECONDITIONERS_HH
25#define DUMUX_LINEAR_PRECONDITIONERS_HH
27#include <dune/common/exceptions.hh>
28#include <dune/common/float_cmp.hh>
29#include <dune/common/indices.hh>
30#include <dune/istl/preconditioners.hh>
31#include <dune/istl/paamg/amg.hh>
34#include <dune/istl/umfpack.hh>
77template<
class M,
class X,
class Y,
int l = 1>
78class SeqUzawa :
public Dune::Preconditioner<X,Y>
81 static_assert(l== 1,
"SeqUzawa expects a block level of 1.");
83 using A = std::decay_t<decltype(std::declval<M>()[Dune::Indices::_0][Dune::Indices::_0])>;
84 using U = std::decay_t<decltype(std::declval<X>()[Dune::Indices::_0])>;
86 using Comm = Dune::Amg::SequentialInformation;
87 using LinearOperator = Dune::MatrixAdapter<A, U, U>;
88 using Smoother = Dune::SeqSSOR<A, U, U>;
89 using AMGSolverForA = Dune::Amg::AMG<LinearOperator, U, Smoother, Comm>;
109 SeqUzawa(
const std::shared_ptr<
const Dune::AssembledLinearOperator<M,X,Y>>& op,
const Dune::ParameterTree& params)
110 : matrix_(op->getmat())
111 , numIterations_(params.get<std::size_t>(
"iterations"))
113 , verbosity_(params.get<int>(
"verbosity"))
114 , paramGroup_(params.get<std::string>(
"ParameterGroup"))
115 , useDirectVelocitySolverForA_(
getParamFromGroup<bool>(paramGroup_,
"LinearSolver.Preconditioner.DirectSolverForA", false))
117 const bool determineRelaxationFactor = getParamFromGroup<bool>(paramGroup_,
"LinearSolver.Preconditioner.DetermineRelaxationFactor",
true);
120 if (determineRelaxationFactor || !useDirectVelocitySolverForA_)
123 if (useDirectVelocitySolverForA_)
126 if (determineRelaxationFactor)
127 relaxationFactor_ = estimateOmega_();
133 virtual void pre(X& x, Y& b) {}
141 virtual void apply(X& update,
const Y& currentDefect)
143 using namespace Dune::Indices;
145 auto& A = matrix_[_0][_0];
146 auto& B = matrix_[_0][_1];
147 auto& C = matrix_[_1][_0];
148 auto& D = matrix_[_1][_1];
150 const auto& f = currentDefect[_0];
151 const auto& g = currentDefect[_1];
152 auto& u = update[_0];
153 auto& p = update[_1];
157 for (std::size_t i = 0; i < D.N(); ++i)
159 const auto& block = D[i][i];
160 for (
auto rowIt = block.begin(); rowIt != block.end(); ++rowIt)
161 if (Dune::FloatCmp::eq<scalar_field_type>(rowIt->one_norm(), 1.0))
162 p[i][rowIt.index()] = g[i][rowIt.index()];
166 for (std::size_t k = 0; k < numIterations_; ++k)
173 applySolverForA_(uIncrement, uRhs);
178 C.mmv(u, pIncrement);
179 D.mmv(p, pIncrement);
180 pIncrement *= relaxationFactor_;
185 std::cout <<
"Uzawa iteration " << k
186 <<
", residual: " << uRhs.two_norm() + pIncrement.two_norm()/relaxationFactor_ << std::endl;
197 virtual Dune::SolverCategory::Category
category()
const
199 return Dune::SolverCategory::sequential;
204 void initAMG_(
const Dune::ParameterTree& params)
206 using namespace Dune::Indices;
207 auto linearOperator = std::make_shared<LinearOperator>(matrix_[_0][_0]);
208 amgSolverForA_ = std::make_unique<AMGSolverForA>(linearOperator, params);
214 using namespace Dune::Indices;
215 umfPackSolverForA_ = std::make_unique<Dune::UMFPack<A>>(matrix_[_0][_0]);
217 DUNE_THROW(Dune::InvalidStateException,
"UMFPack not available. Use LinearSolver.Preconditioner.DirectVelocitySolver = false.");
242 using namespace Dune::Indices;
243 auto& A = matrix_[_0][_0];
244 auto& B = matrix_[_0][_1];
245 auto& C = matrix_[_1][_0];
253 static const auto iterations = Dumux::getParamFromGroup<std::size_t>(paramGroup_,
"LinearSolver.Preconditioner.PowerLawIterations", 5);
256 for (std::size_t i = 0; i < iterations; ++i)
264 applySolverForA_(ainvbx, bx);
272 lambdaMax = x.dot(v)/(x.dot(x));
275 omega = 1.0/lambdaMax;
284 std::cout <<
"\n*** Uzawa Preconditioner ***" << std::endl;
285 std::cout <<
"Estimating relaxation factor based on Schur complement" << std::endl;
286 std::cout <<
"Largest estimated eigenvalue lambdaMax = " << lambdaMax << std::endl;
287 std::cout <<
"Relaxation factor omega = 1/lambdaMax = " << omega << std::endl;
293 template<
class Sol,
class Rhs>
294 void applySolverForA_(Sol& sol, Rhs& rhs)
const
296 if (useDirectVelocitySolverForA_)
299 Dune::InverseOperatorResult res;
300 umfPackSolverForA_->apply(sol, rhs, res);
305 amgSolverForA_->pre(sol, rhs);
306 amgSolverForA_->apply(sol, rhs);
307 amgSolverForA_->post(sol);
314 const std::size_t numIterations_;
318 const int verbosity_;
320 std::unique_ptr<AMGSolverForA> amgSolverForA_;
322 std::unique_ptr<Dune::UMFPack<A>> umfPackSolverForA_;
324 const std::string paramGroup_;
325 const bool useDirectVelocitySolverForA_;
Type traits to be used with matrix types.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The specialized Dumux macro and tag for the ISTL registry to choose the solver and preconditioner at ...
T getParamFromGroup(Args &&... args)
A free function to get a parameter from the parameter tree singleton with a model group.
Definition: parameters.hh:161
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
DUMUX_REGISTER_PRECONDITIONER("uzawa", Dumux::MultiTypeBlockMatrixPreconditionerTag, Dune::defaultPreconditionerBlockLevelCreator< Dumux::SeqUzawa, 1 >())
Helper type to determine whether a given type is a Dune::MultiTypeBlockMatrix.
Definition: matrix.hh:49
A preconditioner based on the Uzawa algorithm for saddle-point problems of the form .
Definition: preconditioners.hh:79
virtual Dune::SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: preconditioners.hh:197
M matrix_type
The matrix type the preconditioner is for.
Definition: preconditioners.hh:93
virtual void pre(X &x, Y &b)
Prepare the preconditioner.
Definition: preconditioners.hh:133
SeqUzawa(const std::shared_ptr< const Dune::AssembledLinearOperator< M, X, Y > > &op, const Dune::ParameterTree ¶ms)
Constructor.
Definition: preconditioners.hh:109
virtual void post(X &x)
Clean up.
Definition: preconditioners.hh:194
Dune::Simd::Scalar< field_type > scalar_field_type
Scalar type underlying the field_type.
Definition: preconditioners.hh:101
typename X::field_type field_type
The field type of the preconditioner.
Definition: preconditioners.hh:99
X domain_type
The domain type of the preconditioner.
Definition: preconditioners.hh:95
Y range_type
The range type of the preconditioner.
Definition: preconditioners.hh:97
virtual void apply(X &update, const Y ¤tDefect)
Apply the preconditioner.
Definition: preconditioners.hh:141