3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
preconditioners.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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_LINEAR_PRECONDITIONERS_HH
25#define DUMUX_LINEAR_PRECONDITIONERS_HH
26
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>
32#include <dune/common/version.hh>
33
34#if HAVE_UMFPACK
35#include <dune/istl/umfpack.hh>
36#endif
37
40
41#if DUNE_VERSION_GTE(DUNE_ISTL,2,7)
43#endif
44
45namespace Dumux {
46
47#if DUNE_VERSION_GTE(DUNE_ISTL,2,7)
82template<class M, class X, class Y, int l = 1>
83class SeqUzawa : public Dune::Preconditioner<X,Y>
84{
85 static_assert(Dumux::isMultiTypeBlockMatrix<M>::value && M::M() == 2 && M::N() == 2, "SeqUzawa expects a 2x2 MultiTypeBlockMatrix.");
86 static_assert(l== 1, "SeqUzawa expects a block level of 1.");
87
88 using A = std::decay_t<decltype(std::declval<M>()[Dune::Indices::_0][Dune::Indices::_0])>;
89 using U = std::decay_t<decltype(std::declval<X>()[Dune::Indices::_0])>;
90
91 using Comm = Dune::Amg::SequentialInformation;
92 using LinearOperator = Dune::MatrixAdapter<A, U, U>;
93 using Smoother = Dune::SeqSSOR<A, U, U>;
94 using AMGSolverForA = Dune::Amg::AMG<LinearOperator, U, Smoother, Comm>;
95
96public:
98 using matrix_type = M;
100 using domain_type = X;
102 using range_type = Y;
104 using field_type = typename X::field_type;
106 using scalar_field_type = Dune::Simd::Scalar<field_type>;
107
114#if DUNE_VERSION_GT(DUNE_ISTL,2,7)
115 SeqUzawa(const std::shared_ptr<const Dune::AssembledLinearOperator<M,X,Y>>& op, const Dune::ParameterTree& params)
116 : matrix_(op->getmat())
117#else
118 SeqUzawa(const M& mat, const Dune::ParameterTree& params)
119 : matrix_(mat)
120#endif
121 , numIterations_(params.get<std::size_t>("iterations"))
122 , relaxationFactor_(params.get<scalar_field_type>("relaxation"))
123 , verbosity_(params.get<int>("verbosity"))
124 , paramGroup_(params.get<std::string>("ParameterGroup"))
125 , useDirectVelocitySolverForA_(getParamFromGroup<bool>(paramGroup_, "LinearSolver.Preconditioner.DirectSolverForA", false))
126 {
127 const bool determineRelaxationFactor = getParamFromGroup<bool>(paramGroup_, "LinearSolver.Preconditioner.DetermineRelaxationFactor", true);
128
129 // AMG is needed for determination of omega
130 if (determineRelaxationFactor || !useDirectVelocitySolverForA_)
131 initAMG_(params);
132
133 if (useDirectVelocitySolverForA_)
134 initUMFPack_();
135
136 if (determineRelaxationFactor)
137 relaxationFactor_ = estimateOmega_();
138 }
139
143 virtual void pre(X& x, Y& b) {}
144
151 virtual void apply(X& update, const Y& currentDefect)
152 {
153 using namespace Dune::Indices;
154
155 auto& A = matrix_[_0][_0];
156 auto& B = matrix_[_0][_1];
157 auto& C = matrix_[_1][_0];
158 auto& D = matrix_[_1][_1];
159
160 const auto& f = currentDefect[_0];
161 const auto& g = currentDefect[_1];
162 auto& u = update[_0];
163 auto& p = update[_1];
164
165 // incorporate Dirichlet cell values
166 // TODO: pass Dirichlet constraint handler from outside
167 for (std::size_t i = 0; i < D.N(); ++i)
168 {
169 const auto& block = D[i][i];
170 for (auto rowIt = block.begin(); rowIt != block.end(); ++rowIt)
171 if (Dune::FloatCmp::eq<scalar_field_type>(rowIt->one_norm(), 1.0))
172 p[i][rowIt.index()] = g[i][rowIt.index()];
173 }
174
175 // the actual Uzawa iteration
176 for (std::size_t k = 0; k < numIterations_; ++k)
177 {
178 // u_k+1 = u_k + Q_A^−1*(f − (A*u_k + B*p_k)),
179 auto uRhs = f;
180 A.mmv(u, uRhs);
181 B.mmv(p, uRhs);
182 auto uIncrement = u;
183 applySolverForA_(uIncrement, uRhs);
184 u += uIncrement;
185
186 // p_k+1 = p_k + omega*(g - C*u_k+1 - D*p_k)
187 auto pIncrement = g;
188 C.mmv(u, pIncrement);
189 D.mmv(p, pIncrement);
190 pIncrement *= relaxationFactor_;
191 p += pIncrement;
192
193 if (verbosity_ > 1)
194 {
195 std::cout << "Uzawa iteration " << k
196 << ", residual: " << uRhs.two_norm() + pIncrement.two_norm()/relaxationFactor_ << std::endl;
197 }
198 }
199 }
200
204 virtual void post(X& x) {}
205
207 virtual Dune::SolverCategory::Category category() const
208 {
209 return Dune::SolverCategory::sequential;
210 }
211
212private:
213
214 void initAMG_(const Dune::ParameterTree& params)
215 {
216 using namespace Dune::Indices;
217 auto linearOperator = std::make_shared<LinearOperator>(matrix_[_0][_0]);
218 amgSolverForA_ = std::make_unique<AMGSolverForA>(linearOperator, params);
219 }
220
221 void initUMFPack_()
222 {
223#if HAVE_UMFPACK
224 using namespace Dune::Indices;
225 umfPackSolverForA_ = std::make_unique<Dune::UMFPack<A>>(matrix_[_0][_0]);
226#else
227 DUNE_THROW(Dune::InvalidStateException, "UMFPack not available. Use LinearSolver.Preconditioner.DirectVelocitySolver = false.");
228#endif
229 }
230
250 scalar_field_type estimateOmega_() const
251 {
252 using namespace Dune::Indices;
253 auto& A = matrix_[_0][_0];
254 auto& B = matrix_[_0][_1];
255 auto& C = matrix_[_1][_0];
256
257 U x(A.M());
258 x = 1.0;
259
260 scalar_field_type omega = 0.0;
261 scalar_field_type lambdaMax = 0.0;
262
263 static const auto iterations = Dumux::getParamFromGroup<std::size_t>(paramGroup_, "LinearSolver.Preconditioner.PowerLawIterations", 5);
264
265 // apply power iteration x_k+1 = M*x_k/|M*x_k| for the matrix M = -C*Ainv*B
266 for (std::size_t i = 0; i < iterations; ++i)
267 {
268 // bx = B*x
269 U bx(x.size());
270 B.mv(x, bx);
271
272 // ainvbx = Ainv*(B*x)
273 auto ainvbx = x;
274 applySolverForA_(ainvbx, bx);
275
276 // v = M*x = -C*(Ainv*B*x)
277 U v(x.size());
278 C.mv(ainvbx, v);
279 v *= -1.0;
280
281 // eigenvalue lambdaMax = xt*M*x/(xt*x) = xt*v/(xt*x);
282 lambdaMax = x.dot(v)/(x.dot(x));
283
284 // relaxation factor omega = 1/lambda;
285 omega = 1.0/lambdaMax;
286
287 // new iterate x = M*x/|M*x| = v/|v|
288 x = v;
289 x /= v.two_norm();
290 }
291
292 if (verbosity_ > 0)
293 {
294 std::cout << "\n*** Uzawa Preconditioner ***" << std::endl;
295 std::cout << "Estimating relaxation factor based on Schur complement" << std::endl;
296 std::cout << "Largest estimated eigenvalue lambdaMax = " << lambdaMax << std::endl;
297 std::cout << "Relaxation factor omega = 1/lambdaMax = " << omega << std::endl;
298 }
299
300 return omega;
301 }
302
303 template<class Sol, class Rhs>
304 void applySolverForA_(Sol& sol, Rhs& rhs) const
305 {
306 if (useDirectVelocitySolverForA_)
307 {
308#if HAVE_UMFPACK
309 Dune::InverseOperatorResult res;
310 umfPackSolverForA_->apply(sol, rhs, res);
311#endif
312 }
313 else
314 {
315 amgSolverForA_->pre(sol, rhs);
316 amgSolverForA_->apply(sol, rhs);
317 amgSolverForA_->post(sol);
318 }
319 }
320
322 const M& matrix_;
324 const std::size_t numIterations_;
326 scalar_field_type relaxationFactor_;
328 const int verbosity_;
329
330 std::unique_ptr<AMGSolverForA> amgSolverForA_;
331#if HAVE_UMFPACK
332 std::unique_ptr<Dune::UMFPack<A>> umfPackSolverForA_;
333#endif
334 const std::string paramGroup_;
335 const bool useDirectVelocitySolverForA_;
336};
337#endif // DUNE_VERSION_GTE(DUNE_ISTL,2,7)
338
339#if DUNE_VERSION_GTE(DUNE_ISTL,2,7)
340DUMUX_REGISTER_PRECONDITIONER("uzawa", Dumux::MultiTypeBlockMatrixPreconditionerTag, Dune::defaultPreconditionerBlockLevelCreator<Dumux::SeqUzawa, 1>());
341#endif
342
343} // end namespace Dumux
344
345#endif // DUMUX_LINEAR_PRECONDITIONERS_HH
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 ...
#define DUMUX_REGISTER_PRECONDITIONER(name, tag,...)
Register a Dumux preconditioner.
Definition: istlsolverregistry.hh:44
T getParamFromGroup(Args &&... args)
A free function to get a parameter from the parameter tree singleton with a model group.
Definition: parameters.hh:375
Definition: adapt.hh:29
Definition: common/pdesolver.hh:35
Helper type to determine whether a given type is a Dune::MultiTypeBlockMatrix.
Definition: matrix.hh:49