3.3.0
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
851class UMFPackBackend : public LinearSolver
852{
853public:
855
856 template<class Matrix, class Vector>
857 bool solve(const Matrix& A, Vector& x, const Vector& b)
858 {
859 static_assert(isBCRSMatrix<Matrix>::value, "UMFPack only works with BCRS matrices!");
860 using BlockType = typename Matrix::block_type;
861 static_assert(BlockType::rows == BlockType::cols, "Matrix block must be quadratic!");
862 constexpr auto blockSize = BlockType::rows;
863
864 Dune::UMFPack<Matrix> solver(A, this->verbosity() > 0);
865
866 Vector bTmp(b);
867 solver.apply(x, bTmp, result_);
868
869 int size = x.size();
870 for (int i = 0; i < size; i++)
871 {
872 for (int j = 0; j < blockSize; j++)
873 {
874 using std::isnan;
875 using std::isinf;
876 if (isnan(x[i][j]) || isinf(x[i][j]))
877 {
878 result_.converged = false;
879 break;
880 }
881 }
882 }
883
884 return result_.converged;
885 }
886
887 std::string name() const
888 {
889 return "UMFPack solver";
890 }
891
892 const Dune::InverseOperatorResult& result() const
893 {
894 return result_;
895 }
896
897private:
898 Dune::InverseOperatorResult result_;
899};
900#endif // HAVE_UMFPACK
901
902
906// \{
907
912template <class LinearSolverTraits>
914{
915public:
917
918 template<class Matrix, class Vector>
919 bool solve(const Matrix& A, Vector& x, const Vector& b)
920 {
921 using Preconditioner = SeqUzawa<Matrix, Vector, Vector>;
922 using Solver = Dune::BiCGSTABSolver<Vector>;
923 static const auto solverParams = LinearSolverParameters<LinearSolverTraits>::createParameterTree(this->paramGroup());
924 return IterativePreconditionedSolverImpl::template solveWithParamTree<Preconditioner, Solver>(A, x, b, solverParams);
925 }
926
927 std::string name() const
928 {
929 return "Uzawa preconditioned BiCGSTAB solver";
930 }
931};
932
937template<class M, class X, class Y, int blockLevel = 2>
938class BlockDiagILU0Preconditioner : public Dune::Preconditioner<X, Y>
939{
940 template<std::size_t i>
941 using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>;
942
943 template<std::size_t i>
944 using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>;
945
946 template<std::size_t i>
947 using BlockILU = Dune::SeqILU<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>, blockLevel-1>;
948
949 using ILUTuple = typename makeFromIndexedType<std::tuple, BlockILU, std::make_index_sequence<M::N()> >::type;
950
951public:
953 using matrix_type = typename std::decay_t<M>;
955 using domain_type = X;
957 using range_type = Y;
959 using field_type = typename X::field_type;
960
967 BlockDiagILU0Preconditioner(const M& m, double w = 1.0)
968 : BlockDiagILU0Preconditioner(m, w, std::make_index_sequence<M::N()>{})
969 {
970 static_assert(blockLevel >= 2, "Only makes sense for MultiTypeBlockMatrix!");
971 }
972
973 void pre (X& v, Y& d) final {}
974
975 void apply (X& v, const Y& d) final
976 {
977 using namespace Dune::Hybrid;
978 forEach(integralRange(Dune::Hybrid::size(ilu_)), [&](const auto i)
979 {
980 std::get<decltype(i)::value>(ilu_).apply(v[i], d[i]);
981 });
982 }
983
984 void post (X&) final {}
985
987 Dune::SolverCategory::Category category() const final
988 {
989 return Dune::SolverCategory::sequential;
990 }
991
992private:
993 template<std::size_t... Is>
994 BlockDiagILU0Preconditioner (const M& m, double w, std::index_sequence<Is...> is)
995 : ilu_(std::make_tuple(BlockILU<Is>(m[Dune::index_constant<Is>{}][Dune::index_constant<Is>{}], w)...))
996 {}
997
998 ILUTuple ilu_;
999};
1000
1001
1010{
1011
1012public:
1014
1015 template<class Matrix, class Vector>
1016 bool solve(const Matrix& M, Vector& x, const Vector& b)
1017 {
1019 Dune::MatrixAdapter<Matrix, Vector, Vector> op(M);
1020 Dune::BiCGSTABSolver<Vector> solver(op, preconditioner, this->residReduction(),
1021 this->maxIter(), this->verbosity());
1022 auto bTmp(b);
1023 solver.apply(x, bTmp, result_);
1024
1025 return result_.converged;
1026 }
1027
1028 const Dune::InverseOperatorResult& result() const
1029 {
1030 return result_;
1031 }
1032
1033 std::string name() const
1034 { return "block-diagonal ILU0-preconditioned BiCGSTAB solver"; }
1035
1036private:
1037 Dune::InverseOperatorResult result_;
1038};
1039
1048{
1049
1050public:
1052
1053 template<int precondBlockLevel = 2, class Matrix, class Vector>
1054 bool solve(const Matrix& M, Vector& x, const Vector& b)
1055 {
1057 Dune::MatrixAdapter<Matrix, Vector, Vector> op(M);
1058 static const int restartGMRes = getParamFromGroup<double>(this->paramGroup(), "LinearSolver.GMResRestart");
1059 Dune::RestartedGMResSolver<Vector> solver(op, preconditioner, this->residReduction(), restartGMRes,
1060 this->maxIter(), this->verbosity());
1061 auto bTmp(b);
1062 solver.apply(x, bTmp, result_);
1063
1064 return result_.converged;
1065 }
1066
1067 const Dune::InverseOperatorResult& result() const
1068 {
1069 return result_;
1070 }
1071
1072 std::string name() const
1073 { return "block-diagonal ILU0-preconditioned restarted GMRes solver"; }
1074
1075private:
1076 Dune::InverseOperatorResult result_;
1077};
1078
1083template<class M, class X, class Y, int blockLevel = 2>
1084class BlockDiagAMGPreconditioner : public Dune::Preconditioner<X, Y>
1085{
1086 template<std::size_t i>
1087 using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>;
1088
1089 template<std::size_t i>
1090 using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>;
1091
1092 template<std::size_t i>
1093 using Smoother = Dune::SeqSSOR<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>>;
1094
1095 template<std::size_t i>
1096 using LinearOperator = Dune::MatrixAdapter<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>>;
1097
1098 template<std::size_t i>
1099 using ScalarProduct = Dune::SeqScalarProduct<VecBlockType<i>>;
1100
1101 template<std::size_t i>
1102 using BlockAMG = Dune::Amg::AMG<LinearOperator<i>, VecBlockType<i>, Smoother<i>>;
1103
1104 using AMGTuple = typename makeFromIndexedType<std::tuple, BlockAMG, std::make_index_sequence<M::N()> >::type;
1105
1106public:
1108 using matrix_type = typename std::decay_t<M>;
1110 using domain_type = X;
1112 using range_type = Y;
1114 using field_type = typename X::field_type;
1115
1123 template<class LOP, class Criterion, class SmootherArgs>
1124 BlockDiagAMGPreconditioner(const LOP& lop, const Criterion& c, const SmootherArgs& sa)
1125 : BlockDiagAMGPreconditioner(lop, c, sa, std::make_index_sequence<M::N()>{})
1126 {
1127 static_assert(blockLevel >= 2, "Only makes sense for MultiTypeBlockMatrix!");
1128 }
1129
1146 void pre (X& v, Y& d) final
1147 {
1148 using namespace Dune::Hybrid;
1149 forEach(integralRange(Dune::Hybrid::size(amg_)), [&](const auto i)
1150 {
1151 std::get<decltype(i)::value>(amg_).pre(v[i], d[i]);
1152 });
1153 }
1154
1166 void apply (X& v, const Y& d) final
1167 {
1168 using namespace Dune::Hybrid;
1169 forEach(integralRange(Dune::Hybrid::size(amg_)), [&](const auto i)
1170 {
1171 std::get<decltype(i)::value>(amg_).apply(v[i], d[i]);
1172 });
1173 }
1174
1184 void post (X& v) final
1185 {
1186 using namespace Dune::Hybrid;
1187 forEach(integralRange(Dune::Hybrid::size(amg_)), [&](const auto i)
1188 {
1189 std::get<decltype(i)::value>(amg_).post(v[i]);
1190 });
1191 }
1192
1194 Dune::SolverCategory::Category category() const final
1195 {
1196 return Dune::SolverCategory::sequential;
1197 }
1198
1199private:
1200 template<class LOP, class Criterion, class SmootherArgs, std::size_t... Is>
1201 BlockDiagAMGPreconditioner (const LOP& lop, const Criterion& c, const SmootherArgs& sa, std::index_sequence<Is...> is)
1202 : amg_(std::make_tuple(BlockAMG<Is>(*std::get<Is>(lop), *std::get<Is>(c), *std::get<Is>(sa))...))
1203 {}
1204
1205 AMGTuple amg_;
1206};
1207
1216{
1217 template<class M, std::size_t i>
1218 using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>;
1219
1220 template<class X, std::size_t i>
1221 using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>;
1222
1223 template<class M, class X, std::size_t i>
1224 using Smoother = Dune::SeqSSOR<DiagBlockType<M, i>, VecBlockType<X, i>, VecBlockType<X, i>>;
1225
1226 template<class M, class X, std::size_t i>
1227 using SmootherArgs = typename Dune::Amg::SmootherTraits<Smoother<M, X, i>>::Arguments;
1228
1229 template<class M, std::size_t i>
1230 using Criterion = Dune::Amg::CoarsenCriterion<Dune::Amg::SymmetricCriterion<DiagBlockType<M, i>, Dune::Amg::FirstDiagonal>>;
1231
1232 template<class M, class X, std::size_t i>
1233 using LinearOperator = Dune::MatrixAdapter<DiagBlockType<M, i>, VecBlockType<X, i>, VecBlockType<X, i>>;
1234
1235public:
1237
1238 // Solve saddle-point problem using a Schur complement based preconditioner
1239 template<class Matrix, class Vector>
1240 bool solve(const Matrix& m, Vector& x, const Vector& b)
1241 {
1244 Dune::Amg::Parameters params(15, 2000, 1.2, 1.6, Dune::Amg::atOnceAccu);
1245 params.setDefaultValuesIsotropic(3);
1246 params.setDebugLevel(this->verbosity());
1247
1248 auto criterion = makeCriterion_<Criterion, Matrix>(params, std::make_index_sequence<Matrix::N()>{});
1249 auto smootherArgs = makeSmootherArgs_<SmootherArgs, Matrix, Vector>(std::make_index_sequence<Matrix::N()>{});
1250
1251 using namespace Dune::Hybrid;
1252 forEach(std::make_index_sequence<Matrix::N()>{}, [&](const auto i)
1253 {
1254 auto& args = std::get<decltype(i)::value>(smootherArgs);
1255 args->iterations = 1;
1256 args->relaxationFactor = 1;
1257 });
1258
1259 auto linearOperator = makeLinearOperator_<LinearOperator, Matrix, Vector>(m, std::make_index_sequence<Matrix::N()>{});
1260
1261 BlockDiagAMGPreconditioner<Matrix, Vector, Vector> preconditioner(linearOperator, criterion, smootherArgs);
1262
1263 Dune::MatrixAdapter<Matrix, Vector, Vector> op(m);
1264 Dune::BiCGSTABSolver<Vector> solver(op, preconditioner, this->residReduction(),
1265 this->maxIter(), this->verbosity());
1266 auto bTmp(b);
1267 solver.apply(x, bTmp, result_);
1268
1269 return result_.converged;
1270 }
1271
1272 const Dune::InverseOperatorResult& result() const
1273 {
1274 return result_;
1275 }
1276
1277 std::string name() const
1278 { return "block-diagonal AMG-preconditioned BiCGSTAB solver"; }
1279
1280private:
1281 template<template<class M, std::size_t i> class Criterion, class Matrix, class Params, std::size_t... Is>
1282 auto makeCriterion_(const Params& p, std::index_sequence<Is...>)
1283 {
1284 return std::make_tuple(std::make_shared<Criterion<Matrix, Is>>(p)...);
1285 }
1286
1287 template<template<class M, class X, std::size_t i> class SmootherArgs, class Matrix, class Vector, std::size_t... Is>
1288 auto makeSmootherArgs_(std::index_sequence<Is...>)
1289 {
1290 return std::make_tuple(std::make_shared<SmootherArgs<Matrix, Vector, Is>>()...);
1291 }
1292
1293 template<template<class M, class X, std::size_t i> class LinearOperator, class Matrix, class Vector, std::size_t... Is>
1294 auto makeLinearOperator_(const Matrix& m, std::index_sequence<Is...>)
1295 {
1296 return std::make_tuple(std::make_shared<LinearOperator<Matrix, Vector, Is>>(m[Dune::index_constant<Is>{}][Dune::index_constant<Is>{}])...);
1297 }
1298
1299 Dune::InverseOperatorResult result_;
1300};
1301
1302// \}
1303
1304} // end namespace Dumux
1305
1306#endif
Type traits to be used with matrix types.
Utilities for template meta programming.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Provides a parallel linear solver based on the ISTL AMG preconditioner and the ISTL BiCGSTAB solver.
Generates a parameter tree required for the linear solvers and precondioners of the Dune ISTL.
Dumux preconditioners for iterative solvers.
Base class for linear solvers.
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
Definition: common/pdesolver.hh:35
constexpr std::size_t blockSize()
Definition: nonlinear/newtonsolver.hh:156
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:914
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:919
std::string name() const
Definition: seqsolverbackend.hh:927
A simple ilu0 block diagonal preconditioner.
Definition: seqsolverbackend.hh:939
typename X::field_type field_type
The field type of the preconditioner.
Definition: seqsolverbackend.hh:959
typename std::decay_t< M > matrix_type
The matrix type the preconditioner is for.
Definition: seqsolverbackend.hh:953
Y range_type
The range type of the preconditioner.
Definition: seqsolverbackend.hh:957
void pre(X &v, Y &d) final
Definition: seqsolverbackend.hh:973
X domain_type
The domain type of the preconditioner.
Definition: seqsolverbackend.hh:955
void post(X &) final
Definition: seqsolverbackend.hh:984
Dune::SolverCategory::Category category() const final
Category of the preconditioner (see SolverCategory::Category)
Definition: seqsolverbackend.hh:987
void apply(X &v, const Y &d) final
Definition: seqsolverbackend.hh:975
BlockDiagILU0Preconditioner(const M &m, double w=1.0)
Constructor.
Definition: seqsolverbackend.hh:967
A simple ilu0 block diagonal preconditioned BiCGSTABSolver.
Definition: seqsolverbackend.hh:1010
std::string name() const
Definition: seqsolverbackend.hh:1033
const Dune::InverseOperatorResult & result() const
Definition: seqsolverbackend.hh:1028
bool solve(const Matrix &M, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:1016
A simple ilu0 block diagonal preconditioned RestartedGMResSolver.
Definition: seqsolverbackend.hh:1048
std::string name() const
Definition: seqsolverbackend.hh:1072
const Dune::InverseOperatorResult & result() const
Definition: seqsolverbackend.hh:1067
bool solve(const Matrix &M, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:1054
A simple ilu0 block diagonal preconditioner.
Definition: seqsolverbackend.hh:1085
void post(X &v) final
Clean up.
Definition: seqsolverbackend.hh:1184
typename std::decay_t< M > matrix_type
The matrix type the preconditioner is for.
Definition: seqsolverbackend.hh:1108
void apply(X &v, const Y &d) final
Apply one step of the preconditioner to the system A(v)=d.
Definition: seqsolverbackend.hh:1166
Dune::SolverCategory::Category category() const final
Category of the preconditioner (see SolverCategory::Category)
Definition: seqsolverbackend.hh:1194
Y range_type
The range type of the preconditioner.
Definition: seqsolverbackend.hh:1112
void pre(X &v, Y &d) final
Prepare the preconditioner.
Definition: seqsolverbackend.hh:1146
typename X::field_type field_type
The field type of the preconditioner.
Definition: seqsolverbackend.hh:1114
BlockDiagAMGPreconditioner(const LOP &lop, const Criterion &c, const SmootherArgs &sa)
Constructor.
Definition: seqsolverbackend.hh:1124
X domain_type
The domain type of the preconditioner.
Definition: seqsolverbackend.hh:1110
A simple ilu0 block diagonal preconditioned BiCGSTABSolver.
Definition: seqsolverbackend.hh:1216
const Dune::InverseOperatorResult & result() const
Definition: seqsolverbackend.hh:1272
std::string name() const
Definition: seqsolverbackend.hh:1277
bool solve(const Matrix &m, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:1240
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