3.2-git
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/version.hh>
38#include <dune/common/hybridutilities.hh>
39
47
48namespace Dumux {
49
67{
68public:
69
70 template<class Preconditioner, class Solver, class SolverInterface, class Matrix, class Vector>
71 static bool solve(const SolverInterface& s, const Matrix& A, Vector& x, const Vector& b,
72 const std::string& modelParamGroup = "")
73 {
74 Preconditioner precond(A, s.precondIter(), s.relaxation());
75
76 // make a linear operator from a matrix
77 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
78 MatrixAdapter linearOperator(A);
79
80 Solver solver(linearOperator, precond, s.residReduction(), s.maxIter(), s.verbosity());
81
82 Vector bTmp(b);
83
84 Dune::InverseOperatorResult result;
85 solver.apply(x, bTmp, result);
86
87 return result.converged;
88 }
89
90 template<class Preconditioner, class Solver, class SolverInterface, class Matrix, class Vector>
91 static bool solveWithGMRes(const SolverInterface& s, const Matrix& A, Vector& x, const Vector& b,
92 const std::string& modelParamGroup = "")
93 {
94 // get the restart threshold
95 const int restartGMRes = getParamFromGroup<int>(modelParamGroup, "LinearSolver.GMResRestart", 10);
96
97 Preconditioner precond(A, s.precondIter(), s.relaxation());
98
99 // make a linear operator from a matrix
100 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
101 MatrixAdapter linearOperator(A);
102
103 Solver solver(linearOperator, precond, s.residReduction(), restartGMRes, s.maxIter(), s.verbosity());
104
105 Vector bTmp(b);
106
107 Dune::InverseOperatorResult result;
108 solver.apply(x, bTmp, result);
109
110 return result.converged;
111 }
112
113 template<class Preconditioner, class Solver, class SolverInterface, class Matrix, class Vector>
114 static bool solveWithILU0Prec(const SolverInterface& s, const Matrix& A, Vector& x, const Vector& b,
115 const std::string& modelParamGroup = "")
116 {
117 Preconditioner precond(A, s.relaxation());
118
119 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
120 MatrixAdapter operatorA(A);
121
122 Solver solver(operatorA, precond, s.residReduction(), s.maxIter(), s.verbosity());
123
124 Vector bTmp(b);
125
126 Dune::InverseOperatorResult result;
127 solver.apply(x, bTmp, result);
128
129 return result.converged;
130 }
131
132 // solve with RestartedGMRes (needs restartGMRes as additional argument)
133 template<class Preconditioner, class Solver, class SolverInterface, class Matrix, class Vector>
134 static bool solveWithILU0PrecGMRes(const SolverInterface& s, const Matrix& A, Vector& x, const Vector& b,
135 const std::string& modelParamGroup = "")
136 {
137 // get the restart threshold
138 const int restartGMRes = getParamFromGroup<int>(modelParamGroup, "LinearSolver.GMResRestart", 10);
139
140 Preconditioner precond(A, s.relaxation());
141
142 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
143 MatrixAdapter operatorA(A);
144
145 Solver solver(operatorA, precond, s.residReduction(), restartGMRes, s.maxIter(), s.verbosity());
146
147 Vector bTmp(b);
148
149 Dune::InverseOperatorResult result;
150 solver.apply(x, bTmp, result);
151
152 return result.converged;
153 }
154
155#if DUNE_VERSION_GTE(DUNE_ISTL,2,7)
156 // solve with generic parameter tree
157 template<class Preconditioner, class Solver, class Matrix, class Vector>
158 static bool solveWithParamTree(const Matrix& A, Vector& x, const Vector& b,
159 const Dune::ParameterTree& params)
160 {
161 // make a linear operator from a matrix
162 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
163 const auto linearOperator = std::make_shared<MatrixAdapter>(A);
164
165#if DUNE_VERSION_GT(DUNE_ISTL,2,7)
166 auto precond = std::make_shared<Preconditioner>(linearOperator, params.sub("preconditioner"));
167#else
168 auto precond = std::make_shared<Preconditioner>(A, params.sub("preconditioner"));
169#endif
170 Solver solver(linearOperator, precond, params);
171
172 Vector bTmp(b);
173
174 Dune::InverseOperatorResult result;
175 solver.apply(x, bTmp, result);
176
177 return result.converged;
178 }
179#endif
180};
181
188template<class M>
189constexpr std::size_t preconditionerBlockLevel() noexcept
190{
192}
193
212{
213public:
215
216 template<class Matrix, class Vector>
217 bool solve(const Matrix& A, Vector& x, const Vector& b)
218 {
219 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
220 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
221 using Solver = Dune::BiCGSTABSolver<Vector>;
222
223 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
224 }
225
226 std::string name() const
227 {
228 return "ILUn preconditioned BiCGSTAB solver";
229 }
230};
231
250{
251public:
253
254 template<class Matrix, class Vector>
255 bool solve(const Matrix& A, Vector& x, const Vector& b)
256 {
257 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
258 using Preconditioner = Dune::SeqSOR<Matrix, Vector, Vector, precondBlockLevel>;
259 using Solver = Dune::BiCGSTABSolver<Vector>;
260
261 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
262 }
263
264 std::string name() const
265 {
266 return "SOR preconditioned BiCGSTAB solver";
267 }
268};
269
288{
289public:
291
292 template<class Matrix, class Vector>
293 bool solve(const Matrix& A, Vector& x, const Vector& b)
294 {
295 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
296 using Preconditioner = Dune::SeqSSOR<Matrix, Vector, Vector, precondBlockLevel>;
297 using Solver = Dune::BiCGSTABSolver<Vector>;
298
299 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
300 }
301
302 std::string name() const
303 {
304 return "SSOR preconditioned BiCGSTAB solver";
305 }
306};
307
326{
327public:
329
330 template<class Matrix, class Vector>
331 bool solve(const Matrix& A, Vector& x, const Vector& b)
332 {
333 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
334 using Preconditioner = Dune::SeqGS<Matrix, Vector, Vector, precondBlockLevel>;
335 using Solver = Dune::BiCGSTABSolver<Vector>;
336
337 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
338 }
339
340 std::string name() const
341 {
342 return "SSOR preconditioned BiCGSTAB solver";
343 }
344};
345
363{
364public:
366
367 template<class Matrix, class Vector>
368 bool solve(const Matrix& A, Vector& x, const Vector& b)
369 {
370 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
371 using Preconditioner = Dune::SeqJac<Matrix, Vector, Vector, precondBlockLevel>;
372 using Solver = Dune::BiCGSTABSolver<Vector>;
373
374 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
375 }
376
377 std::string name() const
378 {
379 return "Jac preconditioned BiCGSTAB solver";
380 }
381};
382
400{
401public:
403
404 template<class Matrix, class Vector>
405 bool solve(const Matrix& A, Vector& x, const Vector& b)
406 {
407 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
408 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
409 using Solver = Dune::CGSolver<Vector>;
410
411 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
412 }
413
414 std::string name() const
415 {
416 return "ILUn preconditioned CG solver";
417 }
418};
419
437{
438public:
440
441 template<class Matrix, class Vector>
442 bool solve(const Matrix& A, Vector& x, const Vector& b)
443 {
444 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
445 using Preconditioner = Dune::SeqSOR<Matrix, Vector, Vector, precondBlockLevel>;
446 using Solver = Dune::CGSolver<Vector>;
447
448 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
449 }
450
451 std::string name() const
452 {
453 return "SOR preconditioned CG solver";
454 }
455};
456
474{
475public:
477
478 template<class Matrix, class Vector>
479 bool solve(const Matrix& A, Vector& x, const Vector& b)
480 {
481 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
482 using Preconditioner = Dune::SeqSSOR<Matrix, Vector, Vector, precondBlockLevel>;
483 using Solver = Dune::CGSolver<Vector>;
484
485 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
486 }
487
488 std::string name() const
489 {
490 return "SSOR preconditioned CG solver";
491 }
492};
493
511{
512public:
514
515 template<class Matrix, class Vector>
516 bool solve(const Matrix& A, Vector& x, const Vector& b)
517 {
518 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
519 using Preconditioner = Dune::SeqGS<Matrix, Vector, Vector, precondBlockLevel>;
520 using Solver = Dune::CGSolver<Vector>;
521
522 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
523 }
524
525 std::string name() const
526 {
527 return "GS preconditioned CG solver";
528 }
529};
530
547{
548public:
550
551 template<class Matrix, class Vector>
552 bool solve(const Matrix& A, Vector& x, const Vector& b)
553 {
554 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
555 using Preconditioner = Dune::SeqJac<Matrix, Vector, Vector, precondBlockLevel>;
556 using Solver = Dune::CGSolver<Vector>;
557
558 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
559 }
560
561 std::string name() const
562 {
563 return "GS preconditioned CG solver";
564 }
565};
566
585{
586public:
588
589 template<class Matrix, class Vector>
590 bool solve(const Matrix& A, Vector& x, const Vector& b)
591 {
592 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
593 using Preconditioner = Dune::SeqSSOR<Matrix, Vector, Vector, precondBlockLevel>;
594 using Solver = Dune::RestartedGMResSolver<Vector>;
595
596 return IterativePreconditionedSolverImpl::template solveWithGMRes<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
597 }
598
599 std::string name() const
600 {
601 return "SSOR preconditioned GMRes solver";
602 }
603};
604
622{
623public:
625
626 template<class Matrix, class Vector>
627 bool solve(const Matrix& A, Vector& x, const Vector& b)
628 {
629 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
630 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
631 using Solver = Dune::BiCGSTABSolver<Vector>;
632
633 return IterativePreconditionedSolverImpl::template solveWithILU0Prec<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
634 }
635
636 std::string name() const
637 {
638 return "ILU0 preconditioned BiCGSTAB solver";
639 }
640};
641
658{
659public:
661
662 template<class Matrix, class Vector>
663 bool solve(const Matrix& A, Vector& x, const Vector& b)
664 {
665 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
666 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
667 using Solver = Dune::CGSolver<Vector>;
668
669 return IterativePreconditionedSolverImpl::template solveWithILU0Prec<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
670 }
671
672 std::string name() const
673 {
674 return "ILU0 preconditioned BiCGSTAB solver";
675 }
676};
677
695{
696public:
698
699 template<class Matrix, class Vector>
700 bool solve(const Matrix& A, Vector& x, const Vector& b)
701 {
702 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
703 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
704 using Solver = Dune::RestartedGMResSolver<Vector>;
705
706 return IterativePreconditionedSolverImpl::template solveWithILU0PrecGMRes<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
707 }
708
709 std::string name() const
710 {
711 return "ILU0 preconditioned BiCGSTAB solver";
712 }
713};
714
733{
734public:
736
737 template<class Matrix, class Vector>
738 bool solve(const Matrix& A, Vector& x, const Vector& b)
739 {
740 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
741 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
742 using Solver = Dune::RestartedGMResSolver<Vector>;
743
744 return IterativePreconditionedSolverImpl::template solveWithGMRes<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
745 }
746
747 std::string name() const
748 {
749 return "ILUn preconditioned GMRes solver";
750 }
751};
752
761{
762public:
764
765 template<class Matrix, class Vector>
766 bool solve(const Matrix& A, Vector& x, const Vector& b)
767 {
768 Vector rhs(b);
769 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
770 Dune::SeqJac<Matrix, Vector, Vector, precondBlockLevel> jac(A, 1, 1.0);
771 jac.pre(x, rhs);
772 jac.apply(x, rhs);
773 jac.post(x);
774 return true;
775 }
776
777 std::string name() const
778 {
779 return "Explicit diagonal matrix solver";
780 }
781};
782
783#if HAVE_SUPERLU
792class SuperLUBackend : public LinearSolver
793{
794public:
796
797 template<class Matrix, class Vector>
798 bool solve(const Matrix& A, Vector& x, const Vector& b)
799 {
800 static_assert(isBCRSMatrix<Matrix>::value, "SuperLU only works with BCRS matrices!");
801 using BlockType = typename Matrix::block_type;
802 static_assert(BlockType::rows == BlockType::cols, "Matrix block must be quadratic!");
803 constexpr auto blockSize = BlockType::rows;
804
805 Dune::SuperLU<Matrix> solver(A, this->verbosity() > 0);
806
807 Vector bTmp(b);
808 solver.apply(x, bTmp, result_);
809
810 int size = x.size();
811 for (int i = 0; i < size; i++)
812 {
813 for (int j = 0; j < blockSize; j++)
814 {
815 using std::isnan;
816 using std::isinf;
817 if (isnan(x[i][j]) || isinf(x[i][j]))
818 {
819 result_.converged = false;
820 break;
821 }
822 }
823 }
824
825 return result_.converged;
826 }
827
828 std::string name() const
829 {
830 return "SuperLU solver";
831 }
832
833 const Dune::InverseOperatorResult& result() const
834 {
835 return result_;
836 }
837
838private:
839 Dune::InverseOperatorResult result_;
840};
841#endif // HAVE_SUPERLU
842
843#if HAVE_UMFPACK
852class UMFPackBackend : public LinearSolver
853{
854public:
856
857 template<class Matrix, class Vector>
858 bool solve(const Matrix& A, Vector& x, const Vector& b)
859 {
860 static_assert(isBCRSMatrix<Matrix>::value, "UMFPack only works with BCRS matrices!");
861 using BlockType = typename Matrix::block_type;
862 static_assert(BlockType::rows == BlockType::cols, "Matrix block must be quadratic!");
863 constexpr auto blockSize = BlockType::rows;
864
865 Dune::UMFPack<Matrix> solver(A, this->verbosity() > 0);
866
867 Vector bTmp(b);
868 solver.apply(x, bTmp, result_);
869
870 int size = x.size();
871 for (int i = 0; i < size; i++)
872 {
873 for (int j = 0; j < blockSize; j++)
874 {
875 using std::isnan;
876 using std::isinf;
877 if (isnan(x[i][j]) || isinf(x[i][j]))
878 {
879 result_.converged = false;
880 break;
881 }
882 }
883 }
884
885 return result_.converged;
886 }
887
888 std::string name() const
889 {
890 return "UMFPack solver";
891 }
892
893 const Dune::InverseOperatorResult& result() const
894 {
895 return result_;
896 }
897
898private:
899 Dune::InverseOperatorResult result_;
900};
901#endif // HAVE_UMFPACK
902
903
907// \{
908
909#if DUNE_VERSION_GTE(DUNE_ISTL,2,7)
914template <class LinearSolverTraits>
915class UzawaBiCGSTABBackend : public LinearSolver
916{
917public:
919
920 template<class Matrix, class Vector>
921 bool solve(const Matrix& A, Vector& x, const Vector& b)
922 {
923 using Preconditioner = SeqUzawa<Matrix, Vector, Vector>;
924 using Solver = Dune::BiCGSTABSolver<Vector>;
925 static const auto solverParams = LinearSolverParameters<LinearSolverTraits>::createParameterTree(this->paramGroup());
926 return IterativePreconditionedSolverImpl::template solveWithParamTree<Preconditioner, Solver>(A, x, b, solverParams);
927 }
928
929 std::string name() const
930 {
931 return "Uzawa preconditioned BiCGSTAB solver";
932 }
933};
934#endif
935
940template<class M, class X, class Y, int blockLevel = 2>
941class BlockDiagILU0Preconditioner : public Dune::Preconditioner<X, Y>
942{
943 template<std::size_t i>
944 using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>;
945
946 template<std::size_t i>
947 using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>;
948
949#if DUNE_VERSION_NEWER(DUNE_ISTL,2,6)
950 template<std::size_t i>
951 using BlockILU = Dune::SeqILU<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>, blockLevel-1>;
952#else
953 template<std::size_t i>
954 using BlockILU = Dune::SeqILU0<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>, blockLevel-1>;
955#endif
956
957 using ILUTuple = typename makeFromIndexedType<std::tuple, BlockILU, std::make_index_sequence<M::N()> >::type;
958
959public:
961 using matrix_type = typename std::decay_t<M>;
963 using domain_type = X;
965 using range_type = Y;
967 using field_type = typename X::field_type;
968
975 BlockDiagILU0Preconditioner(const M& m, double w = 1.0)
976 : BlockDiagILU0Preconditioner(m, w, std::make_index_sequence<M::N()>{})
977 {
978 static_assert(blockLevel >= 2, "Only makes sense for MultiTypeBlockMatrix!");
979 }
980
981 void pre (X& v, Y& d) final {}
982
983 void apply (X& v, const Y& d) final
984 {
985 using namespace Dune::Hybrid;
986 forEach(integralRange(Dune::Hybrid::size(ilu_)), [&](const auto i)
987 {
988 std::get<decltype(i)::value>(ilu_).apply(v[i], d[i]);
989 });
990 }
991
992 void post (X&) final {}
993
995 Dune::SolverCategory::Category category() const final
996 {
997 return Dune::SolverCategory::sequential;
998 }
999
1000private:
1001 template<std::size_t... Is>
1002 BlockDiagILU0Preconditioner (const M& m, double w, std::index_sequence<Is...> is)
1003 : ilu_(std::make_tuple(BlockILU<Is>(m[Dune::index_constant<Is>{}][Dune::index_constant<Is>{}], w)...))
1004 {}
1005
1006 ILUTuple ilu_;
1007};
1008
1009
1018{
1019
1020public:
1022
1023 // Solve saddle-point problem using a Schur complement based preconditioner
1024 template<class Matrix, class Vector>
1025 bool solve(const Matrix& M, Vector& x, const Vector& b)
1026 {
1028 Dune::MatrixAdapter<Matrix, Vector, Vector> op(M);
1029 Dune::BiCGSTABSolver<Vector> solver(op, preconditioner, this->residReduction(),
1030 this->maxIter(), this->verbosity());
1031 auto bTmp(b);
1032 solver.apply(x, bTmp, result_);
1033
1034 return result_.converged;
1035 }
1036
1037 const Dune::InverseOperatorResult& result() const
1038 {
1039 return result_;
1040 }
1041
1042 std::string name() const
1043 { return "block-diagonal ILU0 preconditioned BiCGSTAB solver"; }
1044
1045private:
1046 Dune::InverseOperatorResult result_;
1047};
1048
1053template<class M, class X, class Y, int blockLevel = 2>
1054class BlockDiagAMGPreconditioner : public Dune::Preconditioner<X, Y>
1055{
1056 template<std::size_t i>
1057 using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>;
1058
1059 template<std::size_t i>
1060 using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>;
1061
1062 template<std::size_t i>
1063 using Smoother = Dune::SeqSSOR<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>>;
1064
1065 template<std::size_t i>
1066 using LinearOperator = Dune::MatrixAdapter<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>>;
1067
1068 template<std::size_t i>
1069 using ScalarProduct = Dune::SeqScalarProduct<VecBlockType<i>>;
1070
1071 template<std::size_t i>
1072 using BlockAMG = Dune::Amg::AMG<LinearOperator<i>, VecBlockType<i>, Smoother<i>>;
1073
1074 using AMGTuple = typename makeFromIndexedType<std::tuple, BlockAMG, std::make_index_sequence<M::N()> >::type;
1075
1076public:
1078 using matrix_type = typename std::decay_t<M>;
1080 using domain_type = X;
1082 using range_type = Y;
1084 using field_type = typename X::field_type;
1085
1093 template<class LOP, class Criterion, class SmootherArgs>
1094 BlockDiagAMGPreconditioner(const LOP& lop, const Criterion& c, const SmootherArgs& sa)
1095 : BlockDiagAMGPreconditioner(lop, c, sa, std::make_index_sequence<M::N()>{})
1096 {
1097 static_assert(blockLevel >= 2, "Only makes sense for MultiTypeBlockMatrix!");
1098 }
1099
1116 void pre (X& v, Y& d) final
1117 {
1118 using namespace Dune::Hybrid;
1119 forEach(integralRange(Dune::Hybrid::size(amg_)), [&](const auto i)
1120 {
1121 std::get<decltype(i)::value>(amg_).pre(v[i], d[i]);
1122 });
1123 }
1124
1136 void apply (X& v, const Y& d) final
1137 {
1138 using namespace Dune::Hybrid;
1139 forEach(integralRange(Dune::Hybrid::size(amg_)), [&](const auto i)
1140 {
1141 std::get<decltype(i)::value>(amg_).apply(v[i], d[i]);
1142 });
1143 }
1144
1154 void post (X& v) final
1155 {
1156 using namespace Dune::Hybrid;
1157 forEach(integralRange(Dune::Hybrid::size(amg_)), [&](const auto i)
1158 {
1159 std::get<decltype(i)::value>(amg_).post(v[i]);
1160 });
1161 }
1162
1164 Dune::SolverCategory::Category category() const final
1165 {
1166 return Dune::SolverCategory::sequential;
1167 }
1168
1169private:
1170 template<class LOP, class Criterion, class SmootherArgs, std::size_t... Is>
1171 BlockDiagAMGPreconditioner (const LOP& lop, const Criterion& c, const SmootherArgs& sa, std::index_sequence<Is...> is)
1172 : amg_(std::make_tuple(BlockAMG<Is>(*std::get<Is>(lop), *std::get<Is>(c), *std::get<Is>(sa))...))
1173 {}
1174
1175 AMGTuple amg_;
1176};
1177
1186{
1187 template<class M, std::size_t i>
1188 using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>;
1189
1190 template<class X, std::size_t i>
1191 using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>;
1192
1193 template<class M, class X, std::size_t i>
1194 using Smoother = Dune::SeqSSOR<DiagBlockType<M, i>, VecBlockType<X, i>, VecBlockType<X, i>>;
1195
1196 template<class M, class X, std::size_t i>
1197 using SmootherArgs = typename Dune::Amg::SmootherTraits<Smoother<M, X, i>>::Arguments;
1198
1199 template<class M, std::size_t i>
1200 using Criterion = Dune::Amg::CoarsenCriterion<Dune::Amg::SymmetricCriterion<DiagBlockType<M, i>, Dune::Amg::FirstDiagonal>>;
1201
1202 template<class M, class X, std::size_t i>
1203 using LinearOperator = Dune::MatrixAdapter<DiagBlockType<M, i>, VecBlockType<X, i>, VecBlockType<X, i>>;
1204
1205public:
1207
1208 // Solve saddle-point problem using a Schur complement based preconditioner
1209 template<class Matrix, class Vector>
1210 bool solve(const Matrix& m, Vector& x, const Vector& b)
1211 {
1214 Dune::Amg::Parameters params(15, 2000, 1.2, 1.6, Dune::Amg::atOnceAccu);
1215 params.setDefaultValuesIsotropic(3);
1216 params.setDebugLevel(this->verbosity());
1217
1218 auto criterion = makeCriterion_<Criterion, Matrix>(params, std::make_index_sequence<Matrix::N()>{});
1219 auto smootherArgs = makeSmootherArgs_<SmootherArgs, Matrix, Vector>(std::make_index_sequence<Matrix::N()>{});
1220
1221 using namespace Dune::Hybrid;
1222 forEach(std::make_index_sequence<Matrix::N()>{}, [&](const auto i)
1223 {
1224 auto& args = std::get<decltype(i)::value>(smootherArgs);
1225 args->iterations = 1;
1226 args->relaxationFactor = 1;
1227 });
1228
1229 auto linearOperator = makeLinearOperator_<LinearOperator, Matrix, Vector>(m, std::make_index_sequence<Matrix::N()>{});
1230
1231 BlockDiagAMGPreconditioner<Matrix, Vector, Vector> preconditioner(linearOperator, criterion, smootherArgs);
1232
1233 Dune::MatrixAdapter<Matrix, Vector, Vector> op(m);
1234 Dune::BiCGSTABSolver<Vector> solver(op, preconditioner, this->residReduction(),
1235 this->maxIter(), this->verbosity());
1236 auto bTmp(b);
1237 solver.apply(x, bTmp, result_);
1238
1239 return result_.converged;
1240 }
1241
1242 const Dune::InverseOperatorResult& result() const
1243 {
1244 return result_;
1245 }
1246
1247 std::string name() const
1248 { return "block-diagonal AMG preconditioned BiCGSTAB solver"; }
1249
1250private:
1251 template<template<class M, std::size_t i> class Criterion, class Matrix, class Params, std::size_t... Is>
1252 auto makeCriterion_(const Params& p, std::index_sequence<Is...>)
1253 {
1254 return std::make_tuple(std::make_shared<Criterion<Matrix, Is>>(p)...);
1255 }
1256
1257 template<template<class M, class X, std::size_t i> class SmootherArgs, class Matrix, class Vector, std::size_t... Is>
1258 auto makeSmootherArgs_(std::index_sequence<Is...>)
1259 {
1260 return std::make_tuple(std::make_shared<SmootherArgs<Matrix, Vector, Is>>()...);
1261 }
1262
1263 template<template<class M, class X, std::size_t i> class LinearOperator, class Matrix, class Vector, std::size_t... Is>
1264 auto makeLinearOperator_(const Matrix& m, std::index_sequence<Is...>)
1265 {
1266 return std::make_tuple(std::make_shared<LinearOperator<Matrix, Vector, Is>>(m[Dune::index_constant<Is>{}][Dune::index_constant<Is>{}])...);
1267 }
1268
1269 Dune::InverseOperatorResult result_;
1270};
1271
1272// \}
1273
1274} // end namespace Dumux
1275
1276#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:189
Definition: adapt.hh:29
Definition: common/pdesolver.hh:35
Helper type to determine whether a given type is a Dune::MultiTypeBlockMatrix.
Definition: matrix.hh:49
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 general solver backend allowing arbitrary preconditioners and solvers.
Definition: seqsolverbackend.hh:67
static bool solveWithILU0PrecGMRes(const SolverInterface &s, const Matrix &A, Vector &x, const Vector &b, const std::string &modelParamGroup="")
Definition: seqsolverbackend.hh:134
static bool solveWithGMRes(const SolverInterface &s, const Matrix &A, Vector &x, const Vector &b, const std::string &modelParamGroup="")
Definition: seqsolverbackend.hh:91
static bool solveWithILU0Prec(const SolverInterface &s, const Matrix &A, Vector &x, const Vector &b, const std::string &modelParamGroup="")
Definition: seqsolverbackend.hh:114
static bool solve(const SolverInterface &s, const Matrix &A, Vector &x, const Vector &b, const std::string &modelParamGroup="")
Definition: seqsolverbackend.hh:71
Sequential ILU(n)-preconditioned BiCSTAB solver.
Definition: seqsolverbackend.hh:212
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:217
std::string name() const
Definition: seqsolverbackend.hh:226
Sequential SOR-preconditioned BiCSTAB solver.
Definition: seqsolverbackend.hh:250
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:255
std::string name() const
Definition: seqsolverbackend.hh:264
Sequential SSOR-preconditioned BiCGSTAB solver.
Definition: seqsolverbackend.hh:288
std::string name() const
Definition: seqsolverbackend.hh:302
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:293
Sequential GS-preconditioned BiCGSTAB solver.
Definition: seqsolverbackend.hh:326
std::string name() const
Definition: seqsolverbackend.hh:340
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:331
Sequential Jacobi-preconditioned BiCSTAB solver.
Definition: seqsolverbackend.hh:363
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:368
std::string name() const
Definition: seqsolverbackend.hh:377
Sequential ILU(n)-preconditioned CG solver.
Definition: seqsolverbackend.hh:400
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:405
std::string name() const
Definition: seqsolverbackend.hh:414
Sequential SOR-preconditioned CG solver.
Definition: seqsolverbackend.hh:437
std::string name() const
Definition: seqsolverbackend.hh:451
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:442
Sequential SSOR-preconditioned CG solver.
Definition: seqsolverbackend.hh:474
std::string name() const
Definition: seqsolverbackend.hh:488
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:479
Sequential GS-preconditioned CG solver.
Definition: seqsolverbackend.hh:511
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:516
std::string name() const
Definition: seqsolverbackend.hh:525
Sequential Jacobi-preconditioned CG solver.
Definition: seqsolverbackend.hh:547
std::string name() const
Definition: seqsolverbackend.hh:561
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:552
Sequential SSOR-preconditioned GMRes solver.
Definition: seqsolverbackend.hh:585
std::string name() const
Definition: seqsolverbackend.hh:599
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:590
Sequential ILU(0)-preconditioned BiCGSTAB solver.
Definition: seqsolverbackend.hh:622
std::string name() const
Definition: seqsolverbackend.hh:636
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:627
Sequential ILU(0)-preconditioned CG solver.
Definition: seqsolverbackend.hh:658
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:663
std::string name() const
Definition: seqsolverbackend.hh:672
Sequential ILU0-preconditioned GMRes solver.
Definition: seqsolverbackend.hh:695
std::string name() const
Definition: seqsolverbackend.hh:709
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:700
Sequential ILU(n)-preconditioned GMRes solver.
Definition: seqsolverbackend.hh:733
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:738
std::string name() const
Definition: seqsolverbackend.hh:747
Solver for simple block-diagonal matrices (e.g. from explicit time stepping schemes)
Definition: seqsolverbackend.hh:761
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:766
std::string name() const
Definition: seqsolverbackend.hh:777
A simple ilu0 block diagonal preconditioner.
Definition: seqsolverbackend.hh:942
typename X::field_type field_type
The field type of the preconditioner.
Definition: seqsolverbackend.hh:967
typename std::decay_t< M > matrix_type
The matrix type the preconditioner is for.
Definition: seqsolverbackend.hh:961
Y range_type
The range type of the preconditioner.
Definition: seqsolverbackend.hh:965
void pre(X &v, Y &d) final
Definition: seqsolverbackend.hh:981
X domain_type
The domain type of the preconditioner.
Definition: seqsolverbackend.hh:963
void post(X &) final
Definition: seqsolverbackend.hh:992
Dune::SolverCategory::Category category() const final
Category of the preconditioner (see SolverCategory::Category)
Definition: seqsolverbackend.hh:995
void apply(X &v, const Y &d) final
Definition: seqsolverbackend.hh:983
BlockDiagILU0Preconditioner(const M &m, double w=1.0)
Constructor.
Definition: seqsolverbackend.hh:975
A simple ilu0 block diagonal preconditioned BiCGSTABSolver.
Definition: seqsolverbackend.hh:1018
std::string name() const
Definition: seqsolverbackend.hh:1042
const Dune::InverseOperatorResult & result() const
Definition: seqsolverbackend.hh:1037
bool solve(const Matrix &M, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:1025
A simple ilu0 block diagonal preconditioner.
Definition: seqsolverbackend.hh:1055
void post(X &v) final
Clean up.
Definition: seqsolverbackend.hh:1154
typename std::decay_t< M > matrix_type
The matrix type the preconditioner is for.
Definition: seqsolverbackend.hh:1078
void apply(X &v, const Y &d) final
Apply one step of the preconditioner to the system A(v)=d.
Definition: seqsolverbackend.hh:1136
Dune::SolverCategory::Category category() const final
Category of the preconditioner (see SolverCategory::Category)
Definition: seqsolverbackend.hh:1164
Y range_type
The range type of the preconditioner.
Definition: seqsolverbackend.hh:1082
void pre(X &v, Y &d) final
Prepare the preconditioner.
Definition: seqsolverbackend.hh:1116
typename X::field_type field_type
The field type of the preconditioner.
Definition: seqsolverbackend.hh:1084
BlockDiagAMGPreconditioner(const LOP &lop, const Criterion &c, const SmootherArgs &sa)
Constructor.
Definition: seqsolverbackend.hh:1094
X domain_type
The domain type of the preconditioner.
Definition: seqsolverbackend.hh:1080
A simple ilu0 block diagonal preconditioned BiCGSTABSolver.
Definition: seqsolverbackend.hh:1186
const Dune::InverseOperatorResult & result() const
Definition: seqsolverbackend.hh:1242
std::string name() const
Definition: seqsolverbackend.hh:1247
bool solve(const Matrix &m, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:1210
Base class for linear solvers.
Definition: solver.hh:37
double residReduction() const
the linear solver residual reduction
Definition: solver.hh:131
int maxIter() const
the maximum number of linear solver iterations
Definition: solver.hh:123
const std::string & paramGroup() const
the parameter group for getting parameter from the parameter tree
Definition: solver.hh:111
LinearSolver(const std::string &paramGroup="")
Contruct the solver.
Definition: solver.hh:53
int verbosity() const
the verbosity level
Definition: solver.hh:115