3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
seqsolverbackend.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_SEQ_SOLVER_BACKEND_HH
25#define DUMUX_SEQ_SOLVER_BACKEND_HH
26
27#include <type_traits>
28#include <tuple>
29#include <utility>
30
31#include <dune/istl/preconditioners.hh>
32#include <dune/istl/solvers.hh>
33#include <dune/istl/superlu.hh>
34#include <dune/istl/umfpack.hh>
35#include <dune/istl/io.hh>
36#include <dune/common/indices.hh>
37#include <dune/common/hybridutilities.hh>
38
46
47namespace Dumux {
48
66{
67public:
68
69 template<class Preconditioner, class Solver, class SolverInterface, class Matrix, class Vector>
70 static bool solve(const SolverInterface& s, const Matrix& A, Vector& x, const Vector& b,
71 const std::string& modelParamGroup = "")
72 {
73 Preconditioner precond(A, s.precondIter(), s.relaxation());
74
75 // make a linear operator from a matrix
76 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
77 MatrixAdapter linearOperator(A);
78
79 Solver solver(linearOperator, precond, s.residReduction(), s.maxIter(), s.verbosity());
80
81 Vector bTmp(b);
82
83 Dune::InverseOperatorResult result;
84 solver.apply(x, bTmp, result);
85
86 return result.converged;
87 }
88
89 template<class Preconditioner, class Solver, class SolverInterface, class Matrix, class Vector>
90 static bool solveWithGMRes(const SolverInterface& s, const Matrix& A, Vector& x, const Vector& b,
91 const std::string& modelParamGroup = "")
92 {
93 // get the restart threshold
94 const int restartGMRes = getParamFromGroup<int>(modelParamGroup, "LinearSolver.GMResRestart", 10);
95
96 Preconditioner precond(A, s.precondIter(), s.relaxation());
97
98 // make a linear operator from a matrix
99 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
100 MatrixAdapter linearOperator(A);
101
102 Solver solver(linearOperator, precond, s.residReduction(), restartGMRes, s.maxIter(), s.verbosity());
103
104 Vector bTmp(b);
105
106 Dune::InverseOperatorResult result;
107 solver.apply(x, bTmp, result);
108
109 return result.converged;
110 }
111
112 template<class Preconditioner, class Solver, class SolverInterface, class Matrix, class Vector>
113 static bool solveWithILU0Prec(const SolverInterface& s, const Matrix& A, Vector& x, const Vector& b,
114 const std::string& modelParamGroup = "")
115 {
116 Preconditioner precond(A, s.relaxation());
117
118 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
119 MatrixAdapter operatorA(A);
120
121 Solver solver(operatorA, precond, s.residReduction(), s.maxIter(), s.verbosity());
122
123 Vector bTmp(b);
124
125 Dune::InverseOperatorResult result;
126 solver.apply(x, bTmp, result);
127
128 return result.converged;
129 }
130
131 // solve with RestartedGMRes (needs restartGMRes as additional argument)
132 template<class Preconditioner, class Solver, class SolverInterface, class Matrix, class Vector>
133 static bool solveWithILU0PrecGMRes(const SolverInterface& s, const Matrix& A, Vector& x, const Vector& b,
134 const std::string& modelParamGroup = "")
135 {
136 // get the restart threshold
137 const int restartGMRes = getParamFromGroup<int>(modelParamGroup, "LinearSolver.GMResRestart", 10);
138
139 Preconditioner precond(A, s.relaxation());
140
141 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
142 MatrixAdapter operatorA(A);
143
144 Solver solver(operatorA, precond, s.residReduction(), restartGMRes, s.maxIter(), s.verbosity());
145
146 Vector bTmp(b);
147
148 Dune::InverseOperatorResult result;
149 solver.apply(x, bTmp, result);
150
151 return result.converged;
152 }
153
154#if DUNE_VERSION_GTE(DUNE_ISTL,2,7)
155 // solve with generic parameter tree
156 template<class Preconditioner, class Solver, class Matrix, class Vector>
157 static bool solveWithParamTree(const Matrix& A, Vector& x, const Vector& b,
158 const Dune::ParameterTree& params)
159 {
160 // make a linear operator from a matrix
161 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
162 const auto linearOperator = std::make_shared<MatrixAdapter>(A);
163
164#if DUNE_VERSION_GTE(DUNE_ISTL,2,8)
165 auto precond = std::make_shared<Preconditioner>(linearOperator, params.sub("preconditioner"));
166#else
167 auto precond = std::make_shared<Preconditioner>(A, params.sub("preconditioner"));
168#endif
169 Solver solver(linearOperator, precond, params);
170
171 Vector bTmp(b);
172
173 Dune::InverseOperatorResult result;
174 solver.apply(x, bTmp, result);
175
176 return result.converged;
177 }
178#endif
179};
180
187template<class M>
188constexpr std::size_t preconditionerBlockLevel() noexcept
189{
191}
192
211{
212public:
214
215 template<class Matrix, class Vector>
216 bool solve(const Matrix& A, Vector& x, const Vector& b)
217 {
218 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
219 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
220 using Solver = Dune::BiCGSTABSolver<Vector>;
221
222 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
223 }
224
225 std::string name() const
226 {
227 return "ILUn preconditioned BiCGSTAB solver";
228 }
229};
230
249{
250public:
252
253 template<class Matrix, class Vector>
254 bool solve(const Matrix& A, Vector& x, const Vector& b)
255 {
256 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
257 using Preconditioner = Dune::SeqSOR<Matrix, Vector, Vector, precondBlockLevel>;
258 using Solver = Dune::BiCGSTABSolver<Vector>;
259
260 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
261 }
262
263 std::string name() const
264 {
265 return "SOR preconditioned BiCGSTAB solver";
266 }
267};
268
287{
288public:
290
291 template<class Matrix, class Vector>
292 bool solve(const Matrix& A, Vector& x, const Vector& b)
293 {
294 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
295 using Preconditioner = Dune::SeqSSOR<Matrix, Vector, Vector, precondBlockLevel>;
296 using Solver = Dune::BiCGSTABSolver<Vector>;
297
298 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
299 }
300
301 std::string name() const
302 {
303 return "SSOR preconditioned BiCGSTAB solver";
304 }
305};
306
325{
326public:
328
329 template<class Matrix, class Vector>
330 bool solve(const Matrix& A, Vector& x, const Vector& b)
331 {
332 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
333 using Preconditioner = Dune::SeqGS<Matrix, Vector, Vector, precondBlockLevel>;
334 using Solver = Dune::BiCGSTABSolver<Vector>;
335
336 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
337 }
338
339 std::string name() const
340 {
341 return "SSOR preconditioned BiCGSTAB solver";
342 }
343};
344
362{
363public:
365
366 template<class Matrix, class Vector>
367 bool solve(const Matrix& A, Vector& x, const Vector& b)
368 {
369 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
370 using Preconditioner = Dune::SeqJac<Matrix, Vector, Vector, precondBlockLevel>;
371 using Solver = Dune::BiCGSTABSolver<Vector>;
372
373 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
374 }
375
376 std::string name() const
377 {
378 return "Jac preconditioned BiCGSTAB solver";
379 }
380};
381
399{
400public:
402
403 template<class Matrix, class Vector>
404 bool solve(const Matrix& A, Vector& x, const Vector& b)
405 {
406 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
407 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
408 using Solver = Dune::CGSolver<Vector>;
409
410 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
411 }
412
413 std::string name() const
414 {
415 return "ILUn preconditioned CG solver";
416 }
417};
418
436{
437public:
439
440 template<class Matrix, class Vector>
441 bool solve(const Matrix& A, Vector& x, const Vector& b)
442 {
443 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
444 using Preconditioner = Dune::SeqSOR<Matrix, Vector, Vector, precondBlockLevel>;
445 using Solver = Dune::CGSolver<Vector>;
446
447 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
448 }
449
450 std::string name() const
451 {
452 return "SOR preconditioned CG solver";
453 }
454};
455
473{
474public:
476
477 template<class Matrix, class Vector>
478 bool solve(const Matrix& A, Vector& x, const Vector& b)
479 {
480 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
481 using Preconditioner = Dune::SeqSSOR<Matrix, Vector, Vector, precondBlockLevel>;
482 using Solver = Dune::CGSolver<Vector>;
483
484 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
485 }
486
487 std::string name() const
488 {
489 return "SSOR preconditioned CG solver";
490 }
491};
492
510{
511public:
513
514 template<class Matrix, class Vector>
515 bool solve(const Matrix& A, Vector& x, const Vector& b)
516 {
517 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
518 using Preconditioner = Dune::SeqGS<Matrix, Vector, Vector, precondBlockLevel>;
519 using Solver = Dune::CGSolver<Vector>;
520
521 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
522 }
523
524 std::string name() const
525 {
526 return "GS preconditioned CG solver";
527 }
528};
529
546{
547public:
549
550 template<class Matrix, class Vector>
551 bool solve(const Matrix& A, Vector& x, const Vector& b)
552 {
553 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
554 using Preconditioner = Dune::SeqJac<Matrix, Vector, Vector, precondBlockLevel>;
555 using Solver = Dune::CGSolver<Vector>;
556
557 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
558 }
559
560 std::string name() const
561 {
562 return "GS preconditioned CG solver";
563 }
564};
565
584{
585public:
587
588 template<class Matrix, class Vector>
589 bool solve(const Matrix& A, Vector& x, const Vector& b)
590 {
591 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
592 using Preconditioner = Dune::SeqSSOR<Matrix, Vector, Vector, precondBlockLevel>;
593 using Solver = Dune::RestartedGMResSolver<Vector>;
594
595 return IterativePreconditionedSolverImpl::template solveWithGMRes<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
596 }
597
598 std::string name() const
599 {
600 return "SSOR preconditioned GMRes solver";
601 }
602};
603
621{
622public:
624
625 template<class Matrix, class Vector>
626 bool solve(const Matrix& A, Vector& x, const Vector& b)
627 {
628 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
629 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
630 using Solver = Dune::BiCGSTABSolver<Vector>;
631
632 return IterativePreconditionedSolverImpl::template solveWithILU0Prec<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
633 }
634
635 std::string name() const
636 {
637 return "ILU0 preconditioned BiCGSTAB solver";
638 }
639};
640
657{
658public:
660
661 template<class Matrix, class Vector>
662 bool solve(const Matrix& A, Vector& x, const Vector& b)
663 {
664 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
665 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
666 using Solver = Dune::CGSolver<Vector>;
667
668 return IterativePreconditionedSolverImpl::template solveWithILU0Prec<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
669 }
670
671 std::string name() const
672 {
673 return "ILU0 preconditioned BiCGSTAB solver";
674 }
675};
676
694{
695public:
697
698 template<class Matrix, class Vector>
699 bool solve(const Matrix& A, Vector& x, const Vector& b)
700 {
701 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
702 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
703 using Solver = Dune::RestartedGMResSolver<Vector>;
704
705 return IterativePreconditionedSolverImpl::template solveWithILU0PrecGMRes<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
706 }
707
708 std::string name() const
709 {
710 return "ILU0 preconditioned BiCGSTAB solver";
711 }
712};
713
732{
733public:
735
736 template<class Matrix, class Vector>
737 bool solve(const Matrix& A, Vector& x, const Vector& b)
738 {
739 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
740 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
741 using Solver = Dune::RestartedGMResSolver<Vector>;
742
743 return IterativePreconditionedSolverImpl::template solveWithGMRes<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
744 }
745
746 std::string name() const
747 {
748 return "ILUn preconditioned GMRes solver";
749 }
750};
751
760{
761public:
763
764 template<class Matrix, class Vector>
765 bool solve(const Matrix& A, Vector& x, const Vector& b)
766 {
767 Vector rhs(b);
768 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
769 Dune::SeqJac<Matrix, Vector, Vector, precondBlockLevel> jac(A, 1, 1.0);
770 jac.pre(x, rhs);
771 jac.apply(x, rhs);
772 jac.post(x);
773 return true;
774 }
775
776 std::string name() const
777 {
778 return "Explicit diagonal matrix solver";
779 }
780};
781
782#if HAVE_SUPERLU
791class SuperLUBackend : public LinearSolver
792{
793public:
795
796 template<class Matrix, class Vector>
797 bool solve(const Matrix& A, Vector& x, const Vector& b)
798 {
799 static_assert(isBCRSMatrix<Matrix>::value, "SuperLU only works with BCRS matrices!");
800 using BlockType = typename Matrix::block_type;
801 static_assert(BlockType::rows == BlockType::cols, "Matrix block must be quadratic!");
802 constexpr auto blockSize = BlockType::rows;
803
804 Dune::SuperLU<Matrix> solver(A, this->verbosity() > 0);
805
806 Vector bTmp(b);
807 solver.apply(x, bTmp, result_);
808
809 int size = x.size();
810 for (int i = 0; i < size; i++)
811 {
812 for (int j = 0; j < blockSize; j++)
813 {
814 using std::isnan;
815 using std::isinf;
816 if (isnan(x[i][j]) || isinf(x[i][j]))
817 {
818 result_.converged = false;
819 break;
820 }
821 }
822 }
823
824 return result_.converged;
825 }
826
827 std::string name() const
828 {
829 return "SuperLU solver";
830 }
831
832 const Dune::InverseOperatorResult& result() const
833 {
834 return result_;
835 }
836
837private:
838 Dune::InverseOperatorResult result_;
839};
840#endif // HAVE_SUPERLU
841
842#if HAVE_UMFPACK
865class UMFPackBackend : public LinearSolver
866{
867public:
868
869 UMFPackBackend(const std::string& paramGroup = "")
870 : LinearSolver(paramGroup)
871 {
872 ordering_ = getParamFromGroup<int>(this->paramGroup(), "LinearSolver.UMFPackOrdering", 1);
873 }
874
876 void setOrdering(int i)
877 { ordering_ = i; }
878
880 int ordering() const
881 { return ordering_; }
882
883 template<class Matrix, class Vector>
884 bool solve(const Matrix& A, Vector& x, const Vector& b)
885 {
886 static_assert(isBCRSMatrix<Matrix>::value, "UMFPack only works with BCRS matrices!");
887 using BlockType = typename Matrix::block_type;
888 static_assert(BlockType::rows == BlockType::cols, "Matrix block must be quadratic!");
889 constexpr auto blockSize = BlockType::rows;
890
891 Dune::UMFPack<Matrix> solver;
892 solver.setVerbosity(this->verbosity() > 0);
893 solver.setOption(UMFPACK_ORDERING, ordering_);
894 solver.setMatrix(A);
895
896 Vector bTmp(b);
897 solver.apply(x, bTmp, result_);
898
899 int size = x.size();
900 for (int i = 0; i < size; i++)
901 {
902 for (int j = 0; j < blockSize; j++)
903 {
904 using std::isnan;
905 using std::isinf;
906 if (isnan(x[i][j]) || isinf(x[i][j]))
907 {
908 result_.converged = false;
909 break;
910 }
911 }
912 }
913
914 return result_.converged;
915 }
916
917 std::string name() const
918 {
919 return "UMFPack solver";
920 }
921
922 const Dune::InverseOperatorResult& result() const
923 {
924 return result_;
925 }
926
927private:
928 Dune::InverseOperatorResult result_;
929 int ordering_;
930};
931#endif // HAVE_UMFPACK
932
933
937// \{
938
943template <class LinearSolverTraits>
945{
946public:
948
949 template<class Matrix, class Vector>
950 bool solve(const Matrix& A, Vector& x, const Vector& b)
951 {
952 using Preconditioner = SeqUzawa<Matrix, Vector, Vector>;
953 using Solver = Dune::BiCGSTABSolver<Vector>;
954 static const auto solverParams = LinearSolverParameters<LinearSolverTraits>::createParameterTree(this->paramGroup());
955 return IterativePreconditionedSolverImpl::template solveWithParamTree<Preconditioner, Solver>(A, x, b, solverParams);
956 }
957
958 std::string name() const
959 {
960 return "Uzawa preconditioned BiCGSTAB solver";
961 }
962};
963
968template<class M, class X, class Y, int blockLevel = 2>
969class BlockDiagILU0Preconditioner : public Dune::Preconditioner<X, Y>
970{
971 template<std::size_t i>
972 using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>;
973
974 template<std::size_t i>
975 using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>;
976
977 template<std::size_t i>
978 using BlockILU = Dune::SeqILU<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>, blockLevel-1>;
979
980 using ILUTuple = typename makeFromIndexedType<std::tuple, BlockILU, std::make_index_sequence<M::N()> >::type;
981
982public:
984 using matrix_type = typename std::decay_t<M>;
986 using domain_type = X;
988 using range_type = Y;
990 using field_type = typename X::field_type;
991
998 BlockDiagILU0Preconditioner(const M& m, double w = 1.0)
999 : BlockDiagILU0Preconditioner(m, w, std::make_index_sequence<M::N()>{})
1000 {
1001 static_assert(blockLevel >= 2, "Only makes sense for MultiTypeBlockMatrix!");
1002 }
1003
1004 void pre (X& v, Y& d) final {}
1005
1006 void apply (X& v, const Y& d) final
1007 {
1008 using namespace Dune::Hybrid;
1009 forEach(integralRange(Dune::Hybrid::size(ilu_)), [&](const auto i)
1010 {
1011 std::get<decltype(i)::value>(ilu_).apply(v[i], d[i]);
1012 });
1013 }
1014
1015 void post (X&) final {}
1016
1018 Dune::SolverCategory::Category category() const final
1019 {
1020 return Dune::SolverCategory::sequential;
1021 }
1022
1023private:
1024 template<std::size_t... Is>
1025 BlockDiagILU0Preconditioner (const M& m, double w, std::index_sequence<Is...> is)
1026 : ilu_(std::make_tuple(BlockILU<Is>(m[Dune::index_constant<Is>{}][Dune::index_constant<Is>{}], w)...))
1027 {}
1028
1029 ILUTuple ilu_;
1030};
1031
1032
1041{
1042
1043public:
1045
1046 template<class Matrix, class Vector>
1047 bool solve(const Matrix& M, Vector& x, const Vector& b)
1048 {
1050 Dune::MatrixAdapter<Matrix, Vector, Vector> op(M);
1051 Dune::BiCGSTABSolver<Vector> solver(op, preconditioner, this->residReduction(),
1052 this->maxIter(), this->verbosity());
1053 auto bTmp(b);
1054 solver.apply(x, bTmp, result_);
1055
1056 return result_.converged;
1057 }
1058
1059 const Dune::InverseOperatorResult& result() const
1060 {
1061 return result_;
1062 }
1063
1064 std::string name() const
1065 { return "block-diagonal ILU0-preconditioned BiCGSTAB solver"; }
1066
1067private:
1068 Dune::InverseOperatorResult result_;
1069};
1070
1079{
1080
1081public:
1083
1084 template<int precondBlockLevel = 2, class Matrix, class Vector>
1085 bool solve(const Matrix& M, Vector& x, const Vector& b)
1086 {
1088 Dune::MatrixAdapter<Matrix, Vector, Vector> op(M);
1089 static const int restartGMRes = getParamFromGroup<int>(this->paramGroup(), "LinearSolver.GMResRestart");
1090 Dune::RestartedGMResSolver<Vector> solver(op, preconditioner, this->residReduction(), restartGMRes,
1091 this->maxIter(), this->verbosity());
1092 auto bTmp(b);
1093 solver.apply(x, bTmp, result_);
1094
1095 return result_.converged;
1096 }
1097
1098 const Dune::InverseOperatorResult& result() const
1099 {
1100 return result_;
1101 }
1102
1103 std::string name() const
1104 { return "block-diagonal ILU0-preconditioned restarted GMRes solver"; }
1105
1106private:
1107 Dune::InverseOperatorResult result_;
1108};
1109
1114template<class M, class X, class Y, int blockLevel = 2>
1115class BlockDiagAMGPreconditioner : public Dune::Preconditioner<X, Y>
1116{
1117 template<std::size_t i>
1118 using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>;
1119
1120 template<std::size_t i>
1121 using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>;
1122
1123 template<std::size_t i>
1124 using Smoother = Dune::SeqSSOR<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>>;
1125
1126 template<std::size_t i>
1127 using LinearOperator = Dune::MatrixAdapter<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>>;
1128
1129 template<std::size_t i>
1130 using ScalarProduct = Dune::SeqScalarProduct<VecBlockType<i>>;
1131
1132 template<std::size_t i>
1133 using BlockAMG = Dune::Amg::AMG<LinearOperator<i>, VecBlockType<i>, Smoother<i>>;
1134
1135 using AMGTuple = typename makeFromIndexedType<std::tuple, BlockAMG, std::make_index_sequence<M::N()> >::type;
1136
1137public:
1139 using matrix_type = typename std::decay_t<M>;
1141 using domain_type = X;
1143 using range_type = Y;
1145 using field_type = typename X::field_type;
1146
1154 template<class LOP, class Criterion, class SmootherArgs>
1155 BlockDiagAMGPreconditioner(const LOP& lop, const Criterion& c, const SmootherArgs& sa)
1156 : BlockDiagAMGPreconditioner(lop, c, sa, std::make_index_sequence<M::N()>{})
1157 {
1158 static_assert(blockLevel >= 2, "Only makes sense for MultiTypeBlockMatrix!");
1159 }
1160
1177 void pre (X& v, Y& d) final
1178 {
1179 using namespace Dune::Hybrid;
1180 forEach(integralRange(Dune::Hybrid::size(amg_)), [&](const auto i)
1181 {
1182 std::get<decltype(i)::value>(amg_).pre(v[i], d[i]);
1183 });
1184 }
1185
1197 void apply (X& v, const Y& d) final
1198 {
1199 using namespace Dune::Hybrid;
1200 forEach(integralRange(Dune::Hybrid::size(amg_)), [&](const auto i)
1201 {
1202 std::get<decltype(i)::value>(amg_).apply(v[i], d[i]);
1203 });
1204 }
1205
1215 void post (X& v) final
1216 {
1217 using namespace Dune::Hybrid;
1218 forEach(integralRange(Dune::Hybrid::size(amg_)), [&](const auto i)
1219 {
1220 std::get<decltype(i)::value>(amg_).post(v[i]);
1221 });
1222 }
1223
1225 Dune::SolverCategory::Category category() const final
1226 {
1227 return Dune::SolverCategory::sequential;
1228 }
1229
1230private:
1231 template<class LOP, class Criterion, class SmootherArgs, std::size_t... Is>
1232 BlockDiagAMGPreconditioner (const LOP& lop, const Criterion& c, const SmootherArgs& sa, std::index_sequence<Is...> is)
1233 : amg_(std::make_tuple(BlockAMG<Is>(*std::get<Is>(lop), *std::get<Is>(c), *std::get<Is>(sa))...))
1234 {}
1235
1236 AMGTuple amg_;
1237};
1238
1247{
1248 template<class M, std::size_t i>
1249 using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>;
1250
1251 template<class X, std::size_t i>
1252 using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>;
1253
1254 template<class M, class X, std::size_t i>
1255 using Smoother = Dune::SeqSSOR<DiagBlockType<M, i>, VecBlockType<X, i>, VecBlockType<X, i>>;
1256
1257 template<class M, class X, std::size_t i>
1258 using SmootherArgs = typename Dune::Amg::SmootherTraits<Smoother<M, X, i>>::Arguments;
1259
1260 template<class M, std::size_t i>
1261 using Criterion = Dune::Amg::CoarsenCriterion<Dune::Amg::SymmetricCriterion<DiagBlockType<M, i>, Dune::Amg::FirstDiagonal>>;
1262
1263 template<class M, class X, std::size_t i>
1264 using LinearOperator = Dune::MatrixAdapter<DiagBlockType<M, i>, VecBlockType<X, i>, VecBlockType<X, i>>;
1265
1266public:
1268
1269 // Solve saddle-point problem using a Schur complement based preconditioner
1270 template<class Matrix, class Vector>
1271 bool solve(const Matrix& m, Vector& x, const Vector& b)
1272 {
1275 Dune::Amg::Parameters params(15, 2000, 1.2, 1.6, Dune::Amg::atOnceAccu);
1276 params.setDefaultValuesIsotropic(3);
1277 params.setDebugLevel(this->verbosity());
1278
1279 auto criterion = makeCriterion_<Criterion, Matrix>(params, std::make_index_sequence<Matrix::N()>{});
1280 auto smootherArgs = makeSmootherArgs_<SmootherArgs, Matrix, Vector>(std::make_index_sequence<Matrix::N()>{});
1281
1282 using namespace Dune::Hybrid;
1283 forEach(std::make_index_sequence<Matrix::N()>{}, [&](const auto i)
1284 {
1285 auto& args = std::get<decltype(i)::value>(smootherArgs);
1286 args->iterations = 1;
1287 args->relaxationFactor = 1;
1288 });
1289
1290 auto linearOperator = makeLinearOperator_<LinearOperator, Matrix, Vector>(m, std::make_index_sequence<Matrix::N()>{});
1291
1292 BlockDiagAMGPreconditioner<Matrix, Vector, Vector> preconditioner(linearOperator, criterion, smootherArgs);
1293
1294 Dune::MatrixAdapter<Matrix, Vector, Vector> op(m);
1295 Dune::BiCGSTABSolver<Vector> solver(op, preconditioner, this->residReduction(),
1296 this->maxIter(), this->verbosity());
1297 auto bTmp(b);
1298 solver.apply(x, bTmp, result_);
1299
1300 return result_.converged;
1301 }
1302
1303 const Dune::InverseOperatorResult& result() const
1304 {
1305 return result_;
1306 }
1307
1308 std::string name() const
1309 { return "block-diagonal AMG-preconditioned BiCGSTAB solver"; }
1310
1311private:
1312 template<template<class M, std::size_t i> class Criterion, class Matrix, class Params, std::size_t... Is>
1313 auto makeCriterion_(const Params& p, std::index_sequence<Is...>)
1314 {
1315 return std::make_tuple(std::make_shared<Criterion<Matrix, Is>>(p)...);
1316 }
1317
1318 template<template<class M, class X, std::size_t i> class SmootherArgs, class Matrix, class Vector, std::size_t... Is>
1319 auto makeSmootherArgs_(std::index_sequence<Is...>)
1320 {
1321 return std::make_tuple(std::make_shared<SmootherArgs<Matrix, Vector, Is>>()...);
1322 }
1323
1324 template<template<class M, class X, std::size_t i> class LinearOperator, class Matrix, class Vector, std::size_t... Is>
1325 auto makeLinearOperator_(const Matrix& m, std::index_sequence<Is...>)
1326 {
1327 return std::make_tuple(std::make_shared<LinearOperator<Matrix, Vector, Is>>(m[Dune::index_constant<Is>{}][Dune::index_constant<Is>{}])...);
1328 }
1329
1330 Dune::InverseOperatorResult result_;
1331};
1332
1333// \}
1334
1335} // end namespace Dumux
1336
1337#endif
Dumux preconditioners for iterative solvers.
Generates a parameter tree required for the linear solvers and precondioners of the Dune ISTL.
Base class for linear solvers.
Provides a parallel linear solver based on the ISTL AMG preconditioner and the ISTL BiCGSTAB solver.
Utilities for template meta programming.
Type traits to be used with matrix types.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
constexpr std::size_t preconditionerBlockLevel() noexcept
Returns the block level for the preconditioner for a given matrix.
Definition: seqsolverbackend.hh:188
Definition: adapt.hh:29
constexpr std::size_t blockSize()
Definition: nonlinear/newtonsolver.hh:185
typename BlockTypeHelper< SolutionVector, Dune::IsNumber< SolutionVector >::value >::type BlockType
Definition: nonlinear/newtonsolver.hh:198
Definition: common/pdesolver.hh:36
Helper type to determine whether a given type is a Dune::MultiTypeBlockMatrix.
Definition: matrix.hh:49
Definition: utility.hh:40
static Dune::ParameterTree createParameterTree(const std::string &paramGroup="")
Create a tree containing parameters required for the linear solvers and precondioners of the Dune IST...
Definition: linearsolverparameters.hh:46
A preconditioner based on the Uzawa algorithm for saddle-point problems of the form .
Definition: preconditioners.hh:80
A general solver backend allowing arbitrary preconditioners and solvers.
Definition: seqsolverbackend.hh:66
static bool solveWithILU0PrecGMRes(const SolverInterface &s, const Matrix &A, Vector &x, const Vector &b, const std::string &modelParamGroup="")
Definition: seqsolverbackend.hh:133
static bool solveWithGMRes(const SolverInterface &s, const Matrix &A, Vector &x, const Vector &b, const std::string &modelParamGroup="")
Definition: seqsolverbackend.hh:90
static bool solveWithILU0Prec(const SolverInterface &s, const Matrix &A, Vector &x, const Vector &b, const std::string &modelParamGroup="")
Definition: seqsolverbackend.hh:113
static bool solve(const SolverInterface &s, const Matrix &A, Vector &x, const Vector &b, const std::string &modelParamGroup="")
Definition: seqsolverbackend.hh:70
Sequential ILU(n)-preconditioned BiCSTAB solver.
Definition: seqsolverbackend.hh:211
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:216
std::string name() const
Definition: seqsolverbackend.hh:225
Sequential SOR-preconditioned BiCSTAB solver.
Definition: seqsolverbackend.hh:249
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:254
std::string name() const
Definition: seqsolverbackend.hh:263
Sequential SSOR-preconditioned BiCGSTAB solver.
Definition: seqsolverbackend.hh:287
std::string name() const
Definition: seqsolverbackend.hh:301
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:292
Sequential GS-preconditioned BiCGSTAB solver.
Definition: seqsolverbackend.hh:325
std::string name() const
Definition: seqsolverbackend.hh:339
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:330
Sequential Jacobi-preconditioned BiCSTAB solver.
Definition: seqsolverbackend.hh:362
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:367
std::string name() const
Definition: seqsolverbackend.hh:376
Sequential ILU(n)-preconditioned CG solver.
Definition: seqsolverbackend.hh:399
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:404
std::string name() const
Definition: seqsolverbackend.hh:413
Sequential SOR-preconditioned CG solver.
Definition: seqsolverbackend.hh:436
std::string name() const
Definition: seqsolverbackend.hh:450
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:441
Sequential SSOR-preconditioned CG solver.
Definition: seqsolverbackend.hh:473
std::string name() const
Definition: seqsolverbackend.hh:487
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:478
Sequential GS-preconditioned CG solver.
Definition: seqsolverbackend.hh:510
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:515
std::string name() const
Definition: seqsolverbackend.hh:524
Sequential Jacobi-preconditioned CG solver.
Definition: seqsolverbackend.hh:546
std::string name() const
Definition: seqsolverbackend.hh:560
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:551
Sequential SSOR-preconditioned GMRes solver.
Definition: seqsolverbackend.hh:584
std::string name() const
Definition: seqsolverbackend.hh:598
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:589
Sequential ILU(0)-preconditioned BiCGSTAB solver.
Definition: seqsolverbackend.hh:621
std::string name() const
Definition: seqsolverbackend.hh:635
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:626
Sequential ILU(0)-preconditioned CG solver.
Definition: seqsolverbackend.hh:657
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:662
std::string name() const
Definition: seqsolverbackend.hh:671
Sequential ILU0-preconditioned GMRes solver.
Definition: seqsolverbackend.hh:694
std::string name() const
Definition: seqsolverbackend.hh:708
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:699
Sequential ILU(n)-preconditioned GMRes solver.
Definition: seqsolverbackend.hh:732
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:737
std::string name() const
Definition: seqsolverbackend.hh:746
Solver for simple block-diagonal matrices (e.g. from explicit time stepping schemes)
Definition: seqsolverbackend.hh:760
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:765
std::string name() const
Definition: seqsolverbackend.hh:776
A Uzawa preconditioned BiCGSTAB solver for saddle-point problems.
Definition: seqsolverbackend.hh:945
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:950
std::string name() const
Definition: seqsolverbackend.hh:958
A simple ilu0 block diagonal preconditioner.
Definition: seqsolverbackend.hh:970
typename X::field_type field_type
The field type of the preconditioner.
Definition: seqsolverbackend.hh:990
typename std::decay_t< M > matrix_type
The matrix type the preconditioner is for.
Definition: seqsolverbackend.hh:984
Y range_type
The range type of the preconditioner.
Definition: seqsolverbackend.hh:988
void pre(X &v, Y &d) final
Definition: seqsolverbackend.hh:1004
X domain_type
The domain type of the preconditioner.
Definition: seqsolverbackend.hh:986
void post(X &) final
Definition: seqsolverbackend.hh:1015
Dune::SolverCategory::Category category() const final
Category of the preconditioner (see SolverCategory::Category)
Definition: seqsolverbackend.hh:1018
void apply(X &v, const Y &d) final
Definition: seqsolverbackend.hh:1006
BlockDiagILU0Preconditioner(const M &m, double w=1.0)
Constructor.
Definition: seqsolverbackend.hh:998
A simple ilu0 block diagonal preconditioned BiCGSTABSolver.
Definition: seqsolverbackend.hh:1041
std::string name() const
Definition: seqsolverbackend.hh:1064
const Dune::InverseOperatorResult & result() const
Definition: seqsolverbackend.hh:1059
bool solve(const Matrix &M, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:1047
A simple ilu0 block diagonal preconditioned RestartedGMResSolver.
Definition: seqsolverbackend.hh:1079
std::string name() const
Definition: seqsolverbackend.hh:1103
const Dune::InverseOperatorResult & result() const
Definition: seqsolverbackend.hh:1098
bool solve(const Matrix &M, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:1085
A simple ilu0 block diagonal preconditioner.
Definition: seqsolverbackend.hh:1116
void post(X &v) final
Clean up.
Definition: seqsolverbackend.hh:1215
typename std::decay_t< M > matrix_type
The matrix type the preconditioner is for.
Definition: seqsolverbackend.hh:1139
void apply(X &v, const Y &d) final
Apply one step of the preconditioner to the system A(v)=d.
Definition: seqsolverbackend.hh:1197
Dune::SolverCategory::Category category() const final
Category of the preconditioner (see SolverCategory::Category)
Definition: seqsolverbackend.hh:1225
Y range_type
The range type of the preconditioner.
Definition: seqsolverbackend.hh:1143
void pre(X &v, Y &d) final
Prepare the preconditioner.
Definition: seqsolverbackend.hh:1177
typename X::field_type field_type
The field type of the preconditioner.
Definition: seqsolverbackend.hh:1145
BlockDiagAMGPreconditioner(const LOP &lop, const Criterion &c, const SmootherArgs &sa)
Constructor.
Definition: seqsolverbackend.hh:1155
X domain_type
The domain type of the preconditioner.
Definition: seqsolverbackend.hh:1141
A simple ilu0 block diagonal preconditioned BiCGSTABSolver.
Definition: seqsolverbackend.hh:1247
const Dune::InverseOperatorResult & result() const
Definition: seqsolverbackend.hh:1303
std::string name() const
Definition: seqsolverbackend.hh:1308
bool solve(const Matrix &m, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:1271
Base class for linear solvers.
Definition: solver.hh:37
double residReduction() const
the linear solver residual reduction
Definition: solver.hh:99
int maxIter() const
the maximum number of linear solver iterations
Definition: solver.hh:91
const std::string & paramGroup() const
the parameter group for getting parameter from the parameter tree
Definition: solver.hh:79
LinearSolver(const std::string &paramGroup="")
Contruct the solver.
Definition: solver.hh:53
int verbosity() const
the verbosity level
Definition: solver.hh:83