version 3.8
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// 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_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
31
32namespace Dumux {
33
68template<class M, class X, class Y, int l = 1>
69class 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
82public:
84 using matrix_type = M;
86 using domain_type = X;
88 using range_type = Y;
90 using field_type = typename X::field_type;
92 using scalar_field_type = Dune::Simd::Scalar<field_type>;
93
100 SeqUzawa(const std::shared_ptr<const Dune::AssembledLinearOperator<M,X,Y>>& op, const Dune::ParameterTree& params)
101 : matrix_(op->getmat())
102 , numIterations_(params.get<std::size_t>("iterations"))
103 , relaxationFactor_(params.get<scalar_field_type>("relaxation"))
104 , verbosity_(params.get<int>("verbosity"))
105 , paramGroup_(params.get<std::string>("ParameterGroup"))
106 , useDirectVelocitySolverForA_(getParamFromGroup<bool>(paramGroup_, "LinearSolver.Preconditioner.DirectSolverForA", false))
107 {
108 const bool determineRelaxationFactor = getParamFromGroup<bool>(paramGroup_, "LinearSolver.Preconditioner.DetermineRelaxationFactor", true);
109
110 // AMG is needed for determination of omega
111 if (determineRelaxationFactor || !useDirectVelocitySolverForA_)
112 initAMG_(params);
113
114 if (useDirectVelocitySolverForA_)
115 initUMFPack_();
116
117 if (determineRelaxationFactor)
118 relaxationFactor_ = estimateOmega_();
119 }
120
124 virtual void pre(X& x, Y& b) {}
125
132 virtual void apply(X& update, const Y& currentDefect)
133 {
134 using namespace Dune::Indices;
135
136 auto& A = matrix_[_0][_0];
137 auto& B = matrix_[_0][_1];
138 auto& C = matrix_[_1][_0];
139 auto& D = matrix_[_1][_1];
140
141 const auto& f = currentDefect[_0];
142 const auto& g = currentDefect[_1];
143 auto& u = update[_0];
144 auto& p = update[_1];
145
146 // incorporate Dirichlet cell values
147 // TODO: pass Dirichlet constraint handler from outside
148 for (std::size_t i = 0; i < D.N(); ++i)
149 {
150 const auto& block = D[i][i];
151 for (auto rowIt = block.begin(); rowIt != block.end(); ++rowIt)
152 if (Dune::FloatCmp::eq<scalar_field_type>(rowIt->one_norm(), 1.0))
153 p[i][rowIt.index()] = g[i][rowIt.index()];
154 }
155
156 // the actual Uzawa iteration
157 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 auto uRhs = f;
161 A.mmv(u, uRhs);
162 B.mmv(p, uRhs);
163 auto uIncrement = u;
164 applySolverForA_(uIncrement, uRhs);
165 u += uIncrement;
166
167 // p_k+1 = p_k + omega*(g - C*u_k+1 - D*p_k)
168 auto pIncrement = g;
169 C.mmv(u, pIncrement);
170 D.mmv(p, pIncrement);
171 pIncrement *= relaxationFactor_;
172 p += pIncrement;
173
174 if (verbosity_ > 1)
175 {
176 std::cout << "Uzawa iteration " << k
177 << ", residual: " << uRhs.two_norm() + pIncrement.two_norm()/relaxationFactor_ << std::endl;
178 }
179 }
180 }
181
185 virtual void post(X& x) {}
186
188 virtual Dune::SolverCategory::Category category() const
189 {
190 return Dune::SolverCategory::sequential;
191 }
192
193private:
194
195 void initAMG_(const Dune::ParameterTree& params)
196 {
197 using namespace Dune::Indices;
198 auto linearOperator = std::make_shared<LinearOperator>(matrix_[_0][_0]);
199 amgSolverForA_ = std::make_unique<AMGSolverForA>(linearOperator, params);
200 }
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
231 scalar_field_type estimateOmega_() const
232 {
233 using namespace Dune::Indices;
234 auto& A = matrix_[_0][_0];
235 auto& B = matrix_[_0][_1];
236 auto& C = matrix_[_1][_0];
237
238 U x(A.M());
239 x = 1.0;
240
241 scalar_field_type omega = 0.0;
242 scalar_field_type lambdaMax = 0.0;
243
244 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 for (std::size_t i = 0; i < iterations; ++i)
248 {
249 // bx = B*x
250 U bx(x.size());
251 B.mv(x, bx);
252
253 // ainvbx = Ainv*(B*x)
254 auto ainvbx = x;
255 applySolverForA_(ainvbx, bx);
256
257 // v = M*x = -C*(Ainv*B*x)
258 U v(x.size());
259 C.mv(ainvbx, v);
260 v *= -1.0;
261
262 // eigenvalue lambdaMax = xt*M*x/(xt*x) = xt*v/(xt*x);
263 lambdaMax = x.dot(v)/(x.dot(x));
264
265 // relaxation factor omega = 1/lambda;
266 omega = 1.0/lambdaMax;
267
268 // new iterate x = M*x/|M*x| = v/|v|
269 x = v;
270 x /= v.two_norm();
271 }
272
273 if (verbosity_ > 0)
274 {
275 std::cout << "\n*** Uzawa Preconditioner ***" << std::endl;
276 std::cout << "Estimating relaxation factor based on Schur complement" << std::endl;
277 std::cout << "Largest estimated eigenvalue lambdaMax = " << lambdaMax << std::endl;
278 std::cout << "Relaxation factor omega = 1/lambdaMax = " << omega << std::endl;
279 }
280
281 return omega;
282 }
283
284 template<class Sol, class Rhs>
285 void applySolverForA_(Sol& sol, Rhs& rhs) const
286 {
287 if (useDirectVelocitySolverForA_)
288 {
289#if HAVE_UMFPACK
290 Dune::InverseOperatorResult res;
291 umfPackSolverForA_->apply(sol, rhs, res);
292#endif
293 }
294 else
295 {
296 amgSolverForA_->pre(sol, rhs);
297 amgSolverForA_->apply(sol, rhs);
298 amgSolverForA_->post(sol);
299 }
300 }
301
303 const M& matrix_;
305 const std::size_t numIterations_;
307 scalar_field_type relaxationFactor_;
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
319DUMUX_REGISTER_PRECONDITIONER("uzawa", Dumux::MultiTypeBlockMatrixPreconditionerTag, Dune::defaultPreconditionerBlockLevelCreator<Dumux::SeqUzawa, 1>());
320
321
329template<class M, class X, class Y, int l=1>
330class ParMTJac : public Dune::Preconditioner<X,Y>
331{
332public:
334 typedef M matrix_type;
336 typedef X domain_type;
338 typedef Y range_type;
340 typedef typename X::field_type field_type;
342 typedef Dune::Simd::Scalar<field_type> scalar_field_type;
344 typedef typename Dune::FieldTraits<scalar_field_type>::real_type real_field_type;
345
347 ParMTJac (const M& A, int n, real_field_type w)
348 : _A_(A), numSteps_(n), relaxationFactor_(w)
349 { Dune::CheckIfDiagonalPresent<M,l>::check(_A_); }
350
351 ParMTJac (const std::shared_ptr<const Dune::AssembledLinearOperator<M,X,Y>>& A, const Dune::ParameterTree& configuration)
352 : ParMTJac(A->getmat(), configuration)
353 {}
354
355 ParMTJac (const M& A, const Dune::ParameterTree& configuration)
356 : ParMTJac(A, configuration.get<int>("iterations",1), configuration.get<real_field_type>("relaxation",1.0))
357 {}
358
359 void pre (X&, Y&) override {}
360
361 void apply (X& update, const Y& defect) override
362 {
363 X xOld(update);
364 const auto& A = _A_;
365
366 for (int k=0; k<numSteps_; k++)
367 {
368 if (k > 0)
369 xOld = update;
370
371 Dumux::parallelFor(A.N(), [&](const std::size_t i)
372 {
373 auto& row = A[i];
374 auto v = update[i];
375 auto rhs = defect[i];
376 const auto endColIt = row.end();
377 auto colIt = row.begin();
378
379 for (; colIt.index()<i; ++colIt)
380 colIt->mmv(xOld[colIt.index()], rhs);
381 const auto diag = colIt;
382 for (; colIt != endColIt; ++colIt)
383 colIt->mmv(xOld[colIt.index()], rhs);
384
385 Dune::algmeta_itsteps<l-1, typename M::block_type>::dbjac(
386 *diag, v, rhs, relaxationFactor_
387 );
388
389 update[i].axpy(relaxationFactor_, v);
390 });
391 }
392 }
393
394 void post (X&) override {}
395
397 Dune::SolverCategory::Category category() const override
398 { return Dune::SolverCategory::sequential; }
399
400private:
401 const M& _A_;
402 int numSteps_;
403 real_field_type relaxationFactor_;
404};
405
406DUMUX_REGISTER_PRECONDITIONER("par_mt_jac", Dune::PreconditionerTag, Dune::defaultPreconditionerBlockLevelCreator<Dumux::ParMTJac>());
407
408
409namespace Detail::Preconditioners {
410
411// compute coloring for parallel sweep
412template<class M>
413void computeColorsForMatrixSweep_(const M& matrix, std::deque<std::vector<std::size_t>>& coloredIndices)
414{
415 // allocate some temporary memory
416 std::vector<int> colors(matrix.N(), -1);
417 std::vector<int> neighborColors; neighborColors.reserve(30);
418 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 for (std::size_t i = 0; i < colors.size(); ++i)
423 {
424 neighborColors.clear();
425 auto& row = matrix[i];
426 const auto endColIt = row.end();
427 for (auto colIt = row.begin(); colIt != endColIt; ++colIt)
428 neighborColors.push_back(colors[colIt.index()]);
429
430 const auto c = Detail::smallestAvailableColor(neighborColors, colorUsed);
431 colors[i] = c;
432
433 // add element to the set of elements with the same color
434 if (c < coloredIndices.size())
435 coloredIndices[c].push_back(i);
436 else
437 coloredIndices.emplace_back(1, i);
438 }
439}
440
441// parallel SOR kernel (relaxed Gauss-Seidel)
442template<bool forward, int l, class M, class X, class Y, class K>
443void parallelBlockSOR_(const M& A, X& update, const Y& b, const K& relaxationFactor,
444 const std::deque<std::vector<std::size_t>>& coloredIndices)
445{
446 for (int color = 0; color < coloredIndices.size(); ++color)
447 {
448 const auto c = forward ? color : coloredIndices.size()-1-color;
449 const auto& rowIndices = coloredIndices[c];
450 Dumux::parallelFor(rowIndices.size(), [&](const std::size_t k)
451 {
452 const auto i = rowIndices[k];
453 auto& row = A[i];
454 auto v = update[i];
455 auto rhs = b[i];
456 const auto endColIt = row.end();
457 auto colIt = row.begin();
458
459 for (; colIt.index()<i; ++colIt)
460 colIt->mmv(update[colIt.index()], rhs);
461 const auto diag = colIt;
462 for (; colIt != endColIt; ++colIt)
463 colIt->mmv(update[colIt.index()], rhs);
464
465 if constexpr (forward)
466 Dune::algmeta_itsteps<l-1,typename M::block_type>::bsorf(
467 *diag, v, rhs, relaxationFactor
468 );
469 else
470 Dune::algmeta_itsteps<l-1,typename M::block_type>::bsorb(
471 *diag, v, rhs, relaxationFactor
472 );
473
474 update[i].axpy(relaxationFactor, v);
475 });
476 }
477}
478
479} // end namespace Detail
480
488template<class M, class X, class Y, int l=1>
489class ParMTSOR : public Dune::Preconditioner<X,Y>
490{
491public:
493 typedef M matrix_type;
495 typedef X domain_type;
497 typedef Y range_type;
499 typedef typename X::field_type field_type;
501 typedef Dune::Simd::Scalar<field_type> scalar_field_type;
503 typedef typename Dune::FieldTraits<scalar_field_type>::real_type real_field_type;
504
506 ParMTSOR (const M& A, int n, real_field_type w)
507 : _A_(A), numSteps_(n), relaxationFactor_(w)
508 {
509 Dune::CheckIfDiagonalPresent<M,l>::check(_A_);
511 }
512
513 ParMTSOR (const std::shared_ptr<const Dune::AssembledLinearOperator<M,X,Y>>& A, const Dune::ParameterTree& configuration)
514 : ParMTSOR(A->getmat(), configuration)
515 {}
516
517 ParMTSOR (const M& A, const Dune::ParameterTree& configuration)
518 : ParMTSOR(A, configuration.get<int>("iterations",1), configuration.get<real_field_type>("relaxation",1.0))
519 {}
520
521 void pre (X&, Y&) override {}
522
523 void apply (X& v, const Y& d) override
524 {
525 this->template apply<true>(v,d);
526 }
527
528 template<bool forward>
529 void apply(X& v, const Y& d)
530 {
531 for (int i=0; i<numSteps_; i++)
532 Detail::Preconditioners::parallelBlockSOR_<forward, l>(_A_, v, d, relaxationFactor_, colors_);
533 }
534
535 void post (X&) override {}
536
538 Dune::SolverCategory::Category category() const override
539 { return Dune::SolverCategory::sequential; }
540
541private:
542 const M& _A_;
543 int numSteps_;
544 real_field_type relaxationFactor_;
545
547 std::deque<std::vector<std::size_t>> colors_;
548};
549
550DUMUX_REGISTER_PRECONDITIONER("par_mt_sor", Dune::PreconditionerTag, Dune::defaultPreconditionerBlockLevelCreator<Dumux::ParMTSOR>());
551
552
560template<class M, class X, class Y, int l=1>
561class ParMTSSOR : public Dune::Preconditioner<X,Y>
562{
563public:
565 typedef M matrix_type;
567 typedef X domain_type;
569 typedef Y range_type;
571 typedef typename X::field_type field_type;
573 typedef Dune::Simd::Scalar<field_type> scalar_field_type;
575 typedef typename Dune::FieldTraits<scalar_field_type>::real_type real_field_type;
576
578 ParMTSSOR (const M& A, int n, real_field_type w)
579 : _A_(A), numSteps_(n), relaxationFactor_(w)
580 {
581 Dune::CheckIfDiagonalPresent<M,l>::check(_A_);
583 }
584
585 ParMTSSOR (const std::shared_ptr<const Dune::AssembledLinearOperator<M,X,Y>>& A, const Dune::ParameterTree& configuration)
586 : ParMTSSOR(A->getmat(), configuration)
587 {}
588
589 ParMTSSOR (const M& A, const Dune::ParameterTree& configuration)
590 : ParMTSSOR(A, configuration.get<int>("iterations",1), configuration.get<real_field_type>("relaxation",1.0))
591 {}
592
593 void pre (X&, Y&) override {}
594
595 void apply (X& v, const Y& d) override
596 {
597 for (int i=0; i<numSteps_; i++)
598 {
599 Detail::Preconditioners::parallelBlockSOR_<true, l>(_A_, v, d, relaxationFactor_, colors_);
600 Detail::Preconditioners::parallelBlockSOR_<false, l>(_A_, v, d, relaxationFactor_, colors_);
601 }
602 }
603
604 void post (X&) override {}
605
607 Dune::SolverCategory::Category category() const override
608 { return Dune::SolverCategory::sequential; }
609
610private:
611 const M& _A_;
612 int numSteps_;
613 real_field_type relaxationFactor_;
614
616 std::deque<std::vector<std::size_t>> colors_;
617};
618
619DUMUX_REGISTER_PRECONDITIONER("par_mt_ssor", Dune::PreconditionerTag, Dune::defaultPreconditionerBlockLevelCreator<Dumux::ParMTSSOR>());
620
621} // end namespace Dumux
622
623namespace Dune::Amg {
624
625// make it possible to use Dumux::ParMTJac as AMG smoother
626template<class M, class X, class Y, int l>
627struct ConstructionTraits<Dumux::ParMTJac<M,X,Y,l> >
628{
629 using Arguments = DefaultConstructionArgs<SeqJac<M,X,Y,l> >;
630 static inline std::shared_ptr<Dumux::ParMTJac<M,X,Y,l>> construct(Arguments& args)
631 {
632 return std::make_shared<Dumux::ParMTJac<M,X,Y,l>>(
633 args.getMatrix(), args.getArgs().iterations, args.getArgs().relaxationFactor
634 );
635 }
636};
637
638// make it possible to use Dumux::ParMTSOR as AMG smoother
639template<class M, class X, class Y, int l>
640struct ConstructionTraits<Dumux::ParMTSOR<M,X,Y,l> >
641{
642 using Arguments = DefaultConstructionArgs<SeqSOR<M,X,Y,l> >;
643 static inline std::shared_ptr<Dumux::ParMTSOR<M,X,Y,l>> construct(Arguments& args)
644 {
645 return std::make_shared<Dumux::ParMTSOR<M,X,Y,l>>(
646 args.getMatrix(), args.getArgs().iterations, args.getArgs().relaxationFactor
647 );
648 }
649};
650
651template<class M, class X, class Y, int l>
652struct SmootherApplier<Dumux::ParMTSOR<M,X,Y,l> >
653{
655 typedef typename Smoother::range_type Range;
657
658 static void preSmooth(Smoother& smoother, Domain& v, Range& d)
659 { smoother.template apply<true>(v,d); }
660
661 static void postSmooth(Smoother& smoother, Domain& v, Range& d)
662 { smoother.template apply<false>(v,d); }
663};
664
665// make it possible to use Dumux::ParMTSSOR as AMG smoother
666template<class M, class X, class Y, int l>
667struct ConstructionTraits<Dumux::ParMTSSOR<M,X,Y,l> >
668{
669 using Arguments = DefaultConstructionArgs<SeqSSOR<M,X,Y,l> >;
670 static inline std::shared_ptr<Dumux::ParMTSSOR<M,X,Y,l>> construct(Arguments& args)
671 {
672 return std::make_shared<Dumux::ParMTSSOR<M,X,Y,l>>(
673 args.getMatrix(), args.getArgs().iterations, args.getArgs().relaxationFactor
674 );
675 }
676};
677
678} // end namespace Dune::Amg
679
680#endif
Multi-threaded Jacobi preconditioner.
Definition: preconditioners.hh:331
ParMTJac(const M &A, const Dune::ParameterTree &configuration)
Definition: preconditioners.hh:355
ParMTJac(const std::shared_ptr< const Dune::AssembledLinearOperator< M, X, Y > > &A, const Dune::ParameterTree &configuration)
Definition: preconditioners.hh:351
ParMTJac(const M &A, int n, real_field_type w)
Constructor.
Definition: preconditioners.hh:347
X::field_type field_type
The field type of the preconditioner.
Definition: preconditioners.hh:340
void pre(X &, Y &) override
Definition: preconditioners.hh:359
Dune::FieldTraits< scalar_field_type >::real_type real_field_type
real scalar type underlying the field_type
Definition: preconditioners.hh:344
void post(X &) override
Definition: preconditioners.hh:394
M matrix_type
The matrix type the preconditioner is for.
Definition: preconditioners.hh:334
Dune::Simd::Scalar< field_type > scalar_field_type
scalar type underlying the field_type
Definition: preconditioners.hh:342
Y range_type
The range type of the preconditioner.
Definition: preconditioners.hh:338
X domain_type
The domain type of the preconditioner.
Definition: preconditioners.hh:336
Dune::SolverCategory::Category category() const override
Category of the preconditioner (see SolverCategory::Category)
Definition: preconditioners.hh:397
void apply(X &update, const Y &defect) override
Definition: preconditioners.hh:361
Multi-threaded SOR preconditioner using coloring.
Definition: preconditioners.hh:490
ParMTSOR(const M &A, const Dune::ParameterTree &configuration)
Definition: preconditioners.hh:517
X domain_type
The domain type of the preconditioner.
Definition: preconditioners.hh:495
M matrix_type
The matrix type the preconditioner is for.
Definition: preconditioners.hh:493
void apply(X &v, const Y &d)
Definition: preconditioners.hh:529
ParMTSOR(const M &A, int n, real_field_type w)
Constructor.
Definition: preconditioners.hh:506
void apply(X &v, const Y &d) override
Definition: preconditioners.hh:523
Dune::Simd::Scalar< field_type > scalar_field_type
scalar type underlying the field_type
Definition: preconditioners.hh:501
void post(X &) override
Definition: preconditioners.hh:535
Dune::SolverCategory::Category category() const override
Category of the preconditioner (see SolverCategory::Category)
Definition: preconditioners.hh:538
Y range_type
The range type of the preconditioner.
Definition: preconditioners.hh:497
void pre(X &, Y &) override
Definition: preconditioners.hh:521
X::field_type field_type
The field type of the preconditioner.
Definition: preconditioners.hh:499
ParMTSOR(const std::shared_ptr< const Dune::AssembledLinearOperator< M, X, Y > > &A, const Dune::ParameterTree &configuration)
Definition: preconditioners.hh:513
Dune::FieldTraits< scalar_field_type >::real_type real_field_type
real scalar type underlying the field_type
Definition: preconditioners.hh:503
Multithreaded SSOR preconditioner using coloring.
Definition: preconditioners.hh:562
ParMTSSOR(const M &A, int n, real_field_type w)
Constructor.
Definition: preconditioners.hh:578
Dune::Simd::Scalar< field_type > scalar_field_type
scalar type underlying the field_type
Definition: preconditioners.hh:573
ParMTSSOR(const std::shared_ptr< const Dune::AssembledLinearOperator< M, X, Y > > &A, const Dune::ParameterTree &configuration)
Definition: preconditioners.hh:585
void pre(X &, Y &) override
Definition: preconditioners.hh:593
Dune::SolverCategory::Category category() const override
Category of the preconditioner (see SolverCategory::Category)
Definition: preconditioners.hh:607
void apply(X &v, const Y &d) override
Definition: preconditioners.hh:595
X domain_type
The domain type of the preconditioner.
Definition: preconditioners.hh:567
Dune::FieldTraits< scalar_field_type >::real_type real_field_type
real scalar type underlying the field_type
Definition: preconditioners.hh:575
ParMTSSOR(const M &A, const Dune::ParameterTree &configuration)
Definition: preconditioners.hh:589
M matrix_type
The matrix type the preconditioner is for.
Definition: preconditioners.hh:565
Y range_type
The range type of the preconditioner.
Definition: preconditioners.hh:569
void post(X &) override
Definition: preconditioners.hh:604
X::field_type field_type
The field type of the preconditioner.
Definition: preconditioners.hh:571
A preconditioner based on the Uzawa algorithm for saddle-point problems of the form .
Definition: preconditioners.hh:70
virtual Dune::SolverCategory::Category category() const
Category of the preconditioner (see SolverCategory::Category)
Definition: preconditioners.hh:188
M matrix_type
The matrix type the preconditioner is for.
Definition: preconditioners.hh:84
virtual void pre(X &x, Y &b)
Prepare the preconditioner.
Definition: preconditioners.hh:124
SeqUzawa(const std::shared_ptr< const Dune::AssembledLinearOperator< M, X, Y > > &op, const Dune::ParameterTree &params)
Constructor.
Definition: preconditioners.hh:100
virtual void post(X &x)
Clean up.
Definition: preconditioners.hh:185
Dune::Simd::Scalar< field_type > scalar_field_type
Scalar type underlying the field_type.
Definition: preconditioners.hh:92
typename X::field_type field_type
The field type of the preconditioner.
Definition: preconditioners.hh:90
X domain_type
The domain type of the preconditioner.
Definition: preconditioners.hh:86
Y range_type
The range type of the preconditioner.
Definition: preconditioners.hh:88
virtual void apply(X &update, const Y &currentDefect)
Apply the preconditioner.
Definition: preconditioners.hh:132
Coloring schemes for shared-memory-parallel assembly.
void parallelFor(const std::size_t count, const FunctorType &functor)
A parallel for loop (multithreading)
Definition: parallel_for.hh:160
T getParamFromGroup(Args &&... args)
A free function to get a parameter from the parameter tree singleton with a model group.
Definition: parameters.hh:149
The specialized Dumux macro and tag for the ISTL registry to choose the solver and preconditioner at ...
Type traits to be used with matrix types.
void parallelBlockSOR_(const M &A, X &update, const Y &b, const K &relaxationFactor, const std::deque< std::vector< std::size_t > > &coloredIndices)
Definition: preconditioners.hh:443
void computeColorsForMatrixSweep_(const M &matrix, std::deque< std::vector< std::size_t > > &coloredIndices)
Definition: preconditioners.hh:413
Definition: adapt.hh:17
DUMUX_REGISTER_PRECONDITIONER("uzawa", Dumux::MultiTypeBlockMatrixPreconditionerTag, Dune::defaultPreconditionerBlockLevelCreator< Dumux::SeqUzawa, 1 >())
Definition: preconditioners.hh:623
Parallel for loop (multithreading)
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Helper type to determine whether a given type is a Dune::MultiTypeBlockMatrix.
Definition: matrix.hh:37
DefaultConstructionArgs< SeqJac< M, X, Y, l > > Arguments
Definition: preconditioners.hh:629
static std::shared_ptr< Dumux::ParMTJac< M, X, Y, l > > construct(Arguments &args)
Definition: preconditioners.hh:630
static std::shared_ptr< Dumux::ParMTSOR< M, X, Y, l > > construct(Arguments &args)
Definition: preconditioners.hh:643
DefaultConstructionArgs< SeqSOR< M, X, Y, l > > Arguments
Definition: preconditioners.hh:642
DefaultConstructionArgs< SeqSSOR< M, X, Y, l > > Arguments
Definition: preconditioners.hh:669
static std::shared_ptr< Dumux::ParMTSSOR< M, X, Y, l > > construct(Arguments &args)
Definition: preconditioners.hh:670
static void preSmooth(Smoother &smoother, Domain &v, Range &d)
Definition: preconditioners.hh:658
Smoother::range_type Range
Definition: preconditioners.hh:655
Dumux::ParMTSOR< M, X, Y, l > Smoother
Definition: preconditioners.hh:654
Smoother::domain_type Domain
Definition: preconditioners.hh:656
static void postSmooth(Smoother &smoother, Domain &v, Range &d)
Definition: preconditioners.hh:661