3.6-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/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 // solve with generic parameter tree
155 template<class Preconditioner, class Solver, class Matrix, class Vector>
156 static bool solveWithParamTree(const Matrix& A, Vector& x, const Vector& b,
157 const Dune::ParameterTree& params)
158 {
159 // make a linear operator from a matrix
160 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
161 const auto linearOperator = std::make_shared<MatrixAdapter>(A);
162
163 auto precond = std::make_shared<Preconditioner>(linearOperator, params.sub("preconditioner"));
164 Solver solver(linearOperator, precond, params);
165
166 Vector bTmp(b);
167
168 Dune::InverseOperatorResult result;
169 solver.apply(x, bTmp, result);
170
171 return result.converged;
172 }
173};
174
181template<class M>
182constexpr std::size_t preconditionerBlockLevel() noexcept
183{
185}
186
205{
206public:
208
209 template<class Matrix, class Vector>
210 bool solve(const Matrix& A, Vector& x, const Vector& b)
211 {
212 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
213 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
214 using Solver = Dune::BiCGSTABSolver<Vector>;
215
216 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
217 }
218
219 std::string name() const
220 {
221 return "ILUn preconditioned BiCGSTAB solver";
222 }
223};
224
243{
244public:
246
247 template<class Matrix, class Vector>
248 bool solve(const Matrix& A, Vector& x, const Vector& b)
249 {
250 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
251 using Preconditioner = Dune::SeqSOR<Matrix, Vector, Vector, precondBlockLevel>;
252 using Solver = Dune::BiCGSTABSolver<Vector>;
253
254 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
255 }
256
257 std::string name() const
258 {
259 return "SOR preconditioned BiCGSTAB solver";
260 }
261};
262
281{
282public:
284
285 template<class Matrix, class Vector>
286 bool solve(const Matrix& A, Vector& x, const Vector& b)
287 {
288 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
289 using Preconditioner = Dune::SeqSSOR<Matrix, Vector, Vector, precondBlockLevel>;
290 using Solver = Dune::BiCGSTABSolver<Vector>;
291
292 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
293 }
294
295 std::string name() const
296 {
297 return "SSOR preconditioned BiCGSTAB solver";
298 }
299};
300
319{
320public:
322
323 template<class Matrix, class Vector>
324 bool solve(const Matrix& A, Vector& x, const Vector& b)
325 {
326 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
327 using Preconditioner = Dune::SeqGS<Matrix, Vector, Vector, precondBlockLevel>;
328 using Solver = Dune::BiCGSTABSolver<Vector>;
329
330 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
331 }
332
333 std::string name() const
334 {
335 return "SSOR preconditioned BiCGSTAB solver";
336 }
337};
338
356{
357public:
359
360 template<class Matrix, class Vector>
361 bool solve(const Matrix& A, Vector& x, const Vector& b)
362 {
363 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
364 using Preconditioner = Dune::SeqJac<Matrix, Vector, Vector, precondBlockLevel>;
365 using Solver = Dune::BiCGSTABSolver<Vector>;
366
367 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
368 }
369
370 std::string name() const
371 {
372 return "Jac preconditioned BiCGSTAB solver";
373 }
374};
375
393{
394public:
396
397 template<class Matrix, class Vector>
398 bool solve(const Matrix& A, Vector& x, const Vector& b)
399 {
400 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
401 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
402 using Solver = Dune::CGSolver<Vector>;
403
404 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
405 }
406
407 std::string name() const
408 {
409 return "ILUn preconditioned CG solver";
410 }
411};
412
430{
431public:
433
434 template<class Matrix, class Vector>
435 bool solve(const Matrix& A, Vector& x, const Vector& b)
436 {
437 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
438 using Preconditioner = Dune::SeqSOR<Matrix, Vector, Vector, precondBlockLevel>;
439 using Solver = Dune::CGSolver<Vector>;
440
441 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
442 }
443
444 std::string name() const
445 {
446 return "SOR preconditioned CG solver";
447 }
448};
449
467{
468public:
470
471 template<class Matrix, class Vector>
472 bool solve(const Matrix& A, Vector& x, const Vector& b)
473 {
474 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
475 using Preconditioner = Dune::SeqSSOR<Matrix, Vector, Vector, precondBlockLevel>;
476 using Solver = Dune::CGSolver<Vector>;
477
478 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
479 }
480
481 std::string name() const
482 {
483 return "SSOR preconditioned CG solver";
484 }
485};
486
504{
505public:
507
508 template<class Matrix, class Vector>
509 bool solve(const Matrix& A, Vector& x, const Vector& b)
510 {
511 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
512 using Preconditioner = Dune::SeqGS<Matrix, Vector, Vector, precondBlockLevel>;
513 using Solver = Dune::CGSolver<Vector>;
514
515 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
516 }
517
518 std::string name() const
519 {
520 return "GS preconditioned CG solver";
521 }
522};
523
540{
541public:
543
544 template<class Matrix, class Vector>
545 bool solve(const Matrix& A, Vector& x, const Vector& b)
546 {
547 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
548 using Preconditioner = Dune::SeqJac<Matrix, Vector, Vector, precondBlockLevel>;
549 using Solver = Dune::CGSolver<Vector>;
550
551 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
552 }
553
554 std::string name() const
555 {
556 return "GS preconditioned CG solver";
557 }
558};
559
578{
579public:
581
582 template<class Matrix, class Vector>
583 bool solve(const Matrix& A, Vector& x, const Vector& b)
584 {
585 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
586 using Preconditioner = Dune::SeqSSOR<Matrix, Vector, Vector, precondBlockLevel>;
587 using Solver = Dune::RestartedGMResSolver<Vector>;
588
589 return IterativePreconditionedSolverImpl::template solveWithGMRes<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
590 }
591
592 std::string name() const
593 {
594 return "SSOR preconditioned GMRes solver";
595 }
596};
597
615{
616public:
618
619 template<class Matrix, class Vector>
620 bool solve(const Matrix& A, Vector& x, const Vector& b)
621 {
622 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
623 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
624 using Solver = Dune::BiCGSTABSolver<Vector>;
625
626 return IterativePreconditionedSolverImpl::template solveWithILU0Prec<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
627 }
628
629 std::string name() const
630 {
631 return "ILU0 preconditioned BiCGSTAB solver";
632 }
633};
634
651{
652public:
654
655 template<class Matrix, class Vector>
656 bool solve(const Matrix& A, Vector& x, const Vector& b)
657 {
658 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
659 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
660 using Solver = Dune::CGSolver<Vector>;
661
662 return IterativePreconditionedSolverImpl::template solveWithILU0Prec<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
663 }
664
665 std::string name() const
666 {
667 return "ILU0 preconditioned BiCGSTAB solver";
668 }
669};
670
688{
689public:
691
692 template<class Matrix, class Vector>
693 bool solve(const Matrix& A, Vector& x, const Vector& b)
694 {
695 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
696 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
697 using Solver = Dune::RestartedGMResSolver<Vector>;
698
699 return IterativePreconditionedSolverImpl::template solveWithILU0PrecGMRes<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
700 }
701
702 std::string name() const
703 {
704 return "ILU0 preconditioned BiCGSTAB solver";
705 }
706};
707
726{
727public:
729
730 template<class Matrix, class Vector>
731 bool solve(const Matrix& A, Vector& x, const Vector& b)
732 {
733 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
734 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
735 using Solver = Dune::RestartedGMResSolver<Vector>;
736
737 return IterativePreconditionedSolverImpl::template solveWithGMRes<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
738 }
739
740 std::string name() const
741 {
742 return "ILUn preconditioned GMRes solver";
743 }
744};
745
754{
755public:
757
758 template<class Matrix, class Vector>
759 bool solve(const Matrix& A, Vector& x, const Vector& b)
760 {
761 Vector rhs(b);
762 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
763 Dune::SeqJac<Matrix, Vector, Vector, precondBlockLevel> jac(A, 1, 1.0);
764 jac.pre(x, rhs);
765 jac.apply(x, rhs);
766 jac.post(x);
767 return true;
768 }
769
770 std::string name() const
771 {
772 return "Explicit diagonal matrix solver";
773 }
774};
775
776#if HAVE_SUPERLU
785class SuperLUBackend : public LinearSolver
786{
787public:
789
790 template<class Matrix, class Vector>
791 bool solve(const Matrix& A, Vector& x, const Vector& b)
792 {
793 static_assert(isBCRSMatrix<Matrix>::value, "SuperLU only works with BCRS matrices!");
794 using BlockType = typename Matrix::block_type;
795 static_assert(BlockType::rows == BlockType::cols, "Matrix block must be quadratic!");
796 constexpr auto blockSize = BlockType::rows;
797
798 Dune::SuperLU<Matrix> solver(A, this->verbosity() > 0);
799
800 Vector bTmp(b);
801 solver.apply(x, bTmp, result_);
802
803 int size = x.size();
804 for (int i = 0; i < size; i++)
805 {
806 for (int j = 0; j < blockSize; j++)
807 {
808 using std::isnan;
809 using std::isinf;
810 if (isnan(x[i][j]) || isinf(x[i][j]))
811 {
812 result_.converged = false;
813 break;
814 }
815 }
816 }
817
818 return result_.converged;
819 }
820
821 std::string name() const
822 {
823 return "SuperLU solver";
824 }
825
826 const Dune::InverseOperatorResult& result() const
827 {
828 return result_;
829 }
830
831private:
832 Dune::InverseOperatorResult result_;
833};
834#endif // HAVE_SUPERLU
835
836#if HAVE_UMFPACK
859class UMFPackBackend : public LinearSolver
860{
861public:
862
863 UMFPackBackend(const std::string& paramGroup = "")
864 : LinearSolver(paramGroup)
865 {
866 ordering_ = getParamFromGroup<int>(this->paramGroup(), "LinearSolver.UMFPackOrdering", 1);
867 }
868
870 void setOrdering(int i)
871 { ordering_ = i; }
872
874 int ordering() const
875 { return ordering_; }
876
877 template<class Matrix, class Vector>
878 bool solve(const Matrix& A, Vector& x, const Vector& b)
879 {
880 static_assert(isBCRSMatrix<Matrix>::value, "UMFPack only works with BCRS matrices!");
881 using BlockType = typename Matrix::block_type;
882 static_assert(BlockType::rows == BlockType::cols, "Matrix block must be quadratic!");
883 constexpr auto blockSize = BlockType::rows;
884
885 Dune::UMFPack<Matrix> solver;
886 solver.setVerbosity(this->verbosity() > 0);
887 solver.setOption(UMFPACK_ORDERING, ordering_);
888 solver.setMatrix(A);
889
890 Vector bTmp(b);
891 solver.apply(x, bTmp, result_);
892
893 int size = x.size();
894 for (int i = 0; i < size; i++)
895 {
896 for (int j = 0; j < blockSize; j++)
897 {
898 using std::isnan;
899 using std::isinf;
900 if (isnan(x[i][j]) || isinf(x[i][j]))
901 {
902 result_.converged = false;
903 break;
904 }
905 }
906 }
907
908 return result_.converged;
909 }
910
911 std::string name() const
912 {
913 return "UMFPack solver";
914 }
915
916 const Dune::InverseOperatorResult& result() const
917 {
918 return result_;
919 }
920
921private:
922 Dune::InverseOperatorResult result_;
923 int ordering_;
924};
925#endif // HAVE_UMFPACK
926
927
931// \{
932
937template <class LinearSolverTraits>
939{
940public:
942
943 template<class Matrix, class Vector>
944 bool solve(const Matrix& A, Vector& x, const Vector& b)
945 {
946 using Preconditioner = SeqUzawa<Matrix, Vector, Vector>;
947 using Solver = Dune::BiCGSTABSolver<Vector>;
948 static const auto solverParams = LinearSolverParameters<LinearSolverTraits>::createParameterTree(this->paramGroup());
949 return IterativePreconditionedSolverImpl::template solveWithParamTree<Preconditioner, Solver>(A, x, b, solverParams);
950 }
951
952 std::string name() const
953 {
954 return "Uzawa preconditioned BiCGSTAB solver";
955 }
956};
957
962template<class M, class X, class Y, int blockLevel = 2>
963class BlockDiagILU0Preconditioner : public Dune::Preconditioner<X, Y>
964{
965 template<std::size_t i>
966 using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>;
967
968 template<std::size_t i>
969 using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>;
970
971 template<std::size_t i>
972 using BlockILU = Dune::SeqILU<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>, blockLevel-1>;
973
974 using ILUTuple = typename makeFromIndexedType<std::tuple, BlockILU, std::make_index_sequence<M::N()> >::type;
975
976public:
978 using matrix_type = typename std::decay_t<M>;
980 using domain_type = X;
982 using range_type = Y;
984 using field_type = typename X::field_type;
985
992 BlockDiagILU0Preconditioner(const M& m, double w = 1.0)
993 : BlockDiagILU0Preconditioner(m, w, std::make_index_sequence<M::N()>{})
994 {
995 static_assert(blockLevel >= 2, "Only makes sense for MultiTypeBlockMatrix!");
996 }
997
998 void pre (X& v, Y& d) final {}
999
1000 void apply (X& v, const Y& d) final
1001 {
1002 using namespace Dune::Hybrid;
1003 forEach(integralRange(Dune::Hybrid::size(ilu_)), [&](const auto i)
1004 {
1005 std::get<decltype(i)::value>(ilu_).apply(v[i], d[i]);
1006 });
1007 }
1008
1009 void post (X&) final {}
1010
1012 Dune::SolverCategory::Category category() const final
1013 {
1014 return Dune::SolverCategory::sequential;
1015 }
1016
1017private:
1018 template<std::size_t... Is>
1019 BlockDiagILU0Preconditioner (const M& m, double w, std::index_sequence<Is...> is)
1020 : ilu_(std::make_tuple(BlockILU<Is>(m[Dune::index_constant<Is>{}][Dune::index_constant<Is>{}], w)...))
1021 {}
1022
1023 ILUTuple ilu_;
1024};
1025
1026
1035{
1036
1037public:
1039
1040 template<class Matrix, class Vector>
1041 bool solve(const Matrix& M, Vector& x, const Vector& b)
1042 {
1044 Dune::MatrixAdapter<Matrix, Vector, Vector> op(M);
1045 Dune::BiCGSTABSolver<Vector> solver(op, preconditioner, this->residReduction(),
1046 this->maxIter(), this->verbosity());
1047 auto bTmp(b);
1048 solver.apply(x, bTmp, result_);
1049
1050 return result_.converged;
1051 }
1052
1053 const Dune::InverseOperatorResult& result() const
1054 {
1055 return result_;
1056 }
1057
1058 std::string name() const
1059 { return "block-diagonal ILU0-preconditioned BiCGSTAB solver"; }
1060
1061private:
1062 Dune::InverseOperatorResult result_;
1063};
1064
1073{
1074
1075public:
1077
1078 template<int precondBlockLevel = 2, class Matrix, class Vector>
1079 bool solve(const Matrix& M, Vector& x, const Vector& b)
1080 {
1082 Dune::MatrixAdapter<Matrix, Vector, Vector> op(M);
1083 static const int restartGMRes = getParamFromGroup<int>(this->paramGroup(), "LinearSolver.GMResRestart");
1084 Dune::RestartedGMResSolver<Vector> solver(op, preconditioner, this->residReduction(), restartGMRes,
1085 this->maxIter(), this->verbosity());
1086 auto bTmp(b);
1087 solver.apply(x, bTmp, result_);
1088
1089 return result_.converged;
1090 }
1091
1092 const Dune::InverseOperatorResult& result() const
1093 {
1094 return result_;
1095 }
1096
1097 std::string name() const
1098 { return "block-diagonal ILU0-preconditioned restarted GMRes solver"; }
1099
1100private:
1101 Dune::InverseOperatorResult result_;
1102};
1103
1108template<class M, class X, class Y, int blockLevel = 2>
1109class BlockDiagAMGPreconditioner : public Dune::Preconditioner<X, Y>
1110{
1111 template<std::size_t i>
1112 using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>;
1113
1114 template<std::size_t i>
1115 using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>;
1116
1117 template<std::size_t i>
1118 using Smoother = Dune::SeqSSOR<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>>;
1119
1120 template<std::size_t i>
1121 using LinearOperator = Dune::MatrixAdapter<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>>;
1122
1123 template<std::size_t i>
1124 using ScalarProduct = Dune::SeqScalarProduct<VecBlockType<i>>;
1125
1126 template<std::size_t i>
1127 using BlockAMG = Dune::Amg::AMG<LinearOperator<i>, VecBlockType<i>, Smoother<i>>;
1128
1129 using AMGTuple = typename makeFromIndexedType<std::tuple, BlockAMG, std::make_index_sequence<M::N()> >::type;
1130
1131public:
1133 using matrix_type = typename std::decay_t<M>;
1135 using domain_type = X;
1137 using range_type = Y;
1139 using field_type = typename X::field_type;
1140
1148 template<class LOP, class Criterion, class SmootherArgs>
1149 BlockDiagAMGPreconditioner(const LOP& lop, const Criterion& c, const SmootherArgs& sa)
1150 : BlockDiagAMGPreconditioner(lop, c, sa, std::make_index_sequence<M::N()>{})
1151 {
1152 static_assert(blockLevel >= 2, "Only makes sense for MultiTypeBlockMatrix!");
1153 }
1154
1171 void pre (X& v, Y& d) final
1172 {
1173 using namespace Dune::Hybrid;
1174 forEach(integralRange(Dune::Hybrid::size(amg_)), [&](const auto i)
1175 {
1176 std::get<decltype(i)::value>(amg_).pre(v[i], d[i]);
1177 });
1178 }
1179
1191 void apply (X& v, const Y& d) final
1192 {
1193 using namespace Dune::Hybrid;
1194 forEach(integralRange(Dune::Hybrid::size(amg_)), [&](const auto i)
1195 {
1196 std::get<decltype(i)::value>(amg_).apply(v[i], d[i]);
1197 });
1198 }
1199
1209 void post (X& v) final
1210 {
1211 using namespace Dune::Hybrid;
1212 forEach(integralRange(Dune::Hybrid::size(amg_)), [&](const auto i)
1213 {
1214 std::get<decltype(i)::value>(amg_).post(v[i]);
1215 });
1216 }
1217
1219 Dune::SolverCategory::Category category() const final
1220 {
1221 return Dune::SolverCategory::sequential;
1222 }
1223
1224private:
1225 template<class LOP, class Criterion, class SmootherArgs, std::size_t... Is>
1226 BlockDiagAMGPreconditioner (const LOP& lop, const Criterion& c, const SmootherArgs& sa, std::index_sequence<Is...> is)
1227 : amg_(std::make_tuple(BlockAMG<Is>(*std::get<Is>(lop), *std::get<Is>(c), *std::get<Is>(sa))...))
1228 {}
1229
1230 AMGTuple amg_;
1231};
1232
1241{
1242 template<class M, std::size_t i>
1243 using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>;
1244
1245 template<class X, std::size_t i>
1246 using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>;
1247
1248 template<class M, class X, std::size_t i>
1249 using Smoother = Dune::SeqSSOR<DiagBlockType<M, i>, VecBlockType<X, i>, VecBlockType<X, i>>;
1250
1251 template<class M, class X, std::size_t i>
1252 using SmootherArgs = typename Dune::Amg::SmootherTraits<Smoother<M, X, i>>::Arguments;
1253
1254 template<class M, std::size_t i>
1255 using Criterion = Dune::Amg::CoarsenCriterion<Dune::Amg::SymmetricCriterion<DiagBlockType<M, i>, Dune::Amg::FirstDiagonal>>;
1256
1257 template<class M, class X, std::size_t i>
1258 using LinearOperator = Dune::MatrixAdapter<DiagBlockType<M, i>, VecBlockType<X, i>, VecBlockType<X, i>>;
1259
1260public:
1262
1263 // Solve saddle-point problem using a Schur complement based preconditioner
1264 template<class Matrix, class Vector>
1265 bool solve(const Matrix& m, Vector& x, const Vector& b)
1266 {
1269 Dune::Amg::Parameters params(15, 2000, 1.2, 1.6, Dune::Amg::atOnceAccu);
1270 params.setDefaultValuesIsotropic(3);
1271 params.setDebugLevel(this->verbosity());
1272
1273 auto criterion = makeCriterion_<Criterion, Matrix>(params, std::make_index_sequence<Matrix::N()>{});
1274 auto smootherArgs = makeSmootherArgs_<SmootherArgs, Matrix, Vector>(std::make_index_sequence<Matrix::N()>{});
1275
1276 using namespace Dune::Hybrid;
1277 forEach(std::make_index_sequence<Matrix::N()>{}, [&](const auto i)
1278 {
1279 auto& args = std::get<decltype(i)::value>(smootherArgs);
1280 args->iterations = 1;
1281 args->relaxationFactor = 1;
1282 });
1283
1284 auto linearOperator = makeLinearOperator_<LinearOperator, Matrix, Vector>(m, std::make_index_sequence<Matrix::N()>{});
1285
1286 BlockDiagAMGPreconditioner<Matrix, Vector, Vector> preconditioner(linearOperator, criterion, smootherArgs);
1287
1288 Dune::MatrixAdapter<Matrix, Vector, Vector> op(m);
1289 Dune::BiCGSTABSolver<Vector> solver(op, preconditioner, this->residReduction(),
1290 this->maxIter(), this->verbosity());
1291 auto bTmp(b);
1292 solver.apply(x, bTmp, result_);
1293
1294 return result_.converged;
1295 }
1296
1297 const Dune::InverseOperatorResult& result() const
1298 {
1299 return result_;
1300 }
1301
1302 std::string name() const
1303 { return "block-diagonal AMG-preconditioned BiCGSTAB solver"; }
1304
1305private:
1306 template<template<class M, std::size_t i> class Criterion, class Matrix, class Params, std::size_t... Is>
1307 auto makeCriterion_(const Params& p, std::index_sequence<Is...>)
1308 {
1309 return std::make_tuple(std::make_shared<Criterion<Matrix, Is>>(p)...);
1310 }
1311
1312 template<template<class M, class X, std::size_t i> class SmootherArgs, class Matrix, class Vector, std::size_t... Is>
1313 auto makeSmootherArgs_(std::index_sequence<Is...>)
1314 {
1315 return std::make_tuple(std::make_shared<SmootherArgs<Matrix, Vector, Is>>()...);
1316 }
1317
1318 template<template<class M, class X, std::size_t i> class LinearOperator, class Matrix, class Vector, std::size_t... Is>
1319 auto makeLinearOperator_(const Matrix& m, std::index_sequence<Is...>)
1320 {
1321 return std::make_tuple(std::make_shared<LinearOperator<Matrix, Vector, Is>>(m[Dune::index_constant<Is>{}][Dune::index_constant<Is>{}])...);
1322 }
1323
1324 Dune::InverseOperatorResult result_;
1325};
1326
1327// \}
1328
1329} // end namespace Dumux
1330
1331#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:182
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
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: deprecated.hh:149
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:49
A preconditioner based on the Uzawa algorithm for saddle-point problems of the form .
Definition: preconditioners.hh:79
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
static bool solveWithParamTree(const Matrix &A, Vector &x, const Vector &b, const Dune::ParameterTree &params)
Definition: seqsolverbackend.hh:156
Sequential ILU(n)-preconditioned BiCSTAB solver.
Definition: seqsolverbackend.hh:205
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:210
std::string name() const
Definition: seqsolverbackend.hh:219
Sequential SOR-preconditioned BiCSTAB solver.
Definition: seqsolverbackend.hh:243
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:248
std::string name() const
Definition: seqsolverbackend.hh:257
Sequential SSOR-preconditioned BiCGSTAB solver.
Definition: seqsolverbackend.hh:281
std::string name() const
Definition: seqsolverbackend.hh:295
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:286
Sequential GS-preconditioned BiCGSTAB solver.
Definition: seqsolverbackend.hh:319
std::string name() const
Definition: seqsolverbackend.hh:333
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:324
Sequential Jacobi-preconditioned BiCSTAB solver.
Definition: seqsolverbackend.hh:356
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:361
std::string name() const
Definition: seqsolverbackend.hh:370
Sequential ILU(n)-preconditioned CG solver.
Definition: seqsolverbackend.hh:393
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:398
std::string name() const
Definition: seqsolverbackend.hh:407
Sequential SOR-preconditioned CG solver.
Definition: seqsolverbackend.hh:430
std::string name() const
Definition: seqsolverbackend.hh:444
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:435
Sequential SSOR-preconditioned CG solver.
Definition: seqsolverbackend.hh:467
std::string name() const
Definition: seqsolverbackend.hh:481
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:472
Sequential GS-preconditioned CG solver.
Definition: seqsolverbackend.hh:504
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:509
std::string name() const
Definition: seqsolverbackend.hh:518
Sequential Jacobi-preconditioned CG solver.
Definition: seqsolverbackend.hh:540
std::string name() const
Definition: seqsolverbackend.hh:554
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:545
Sequential SSOR-preconditioned GMRes solver.
Definition: seqsolverbackend.hh:578
std::string name() const
Definition: seqsolverbackend.hh:592
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:583
Sequential ILU(0)-preconditioned BiCGSTAB solver.
Definition: seqsolverbackend.hh:615
std::string name() const
Definition: seqsolverbackend.hh:629
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:620
Sequential ILU(0)-preconditioned CG solver.
Definition: seqsolverbackend.hh:651
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:656
std::string name() const
Definition: seqsolverbackend.hh:665
Sequential ILU0-preconditioned GMRes solver.
Definition: seqsolverbackend.hh:688
std::string name() const
Definition: seqsolverbackend.hh:702
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:693
Sequential ILU(n)-preconditioned GMRes solver.
Definition: seqsolverbackend.hh:726
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:731
std::string name() const
Definition: seqsolverbackend.hh:740
Solver for simple block-diagonal matrices (e.g. from explicit time stepping schemes)
Definition: seqsolverbackend.hh:754
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:759
std::string name() const
Definition: seqsolverbackend.hh:770
A Uzawa preconditioned BiCGSTAB solver for saddle-point problems.
Definition: seqsolverbackend.hh:939
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:944
std::string name() const
Definition: seqsolverbackend.hh:952
A simple ilu0 block diagonal preconditioner.
Definition: seqsolverbackend.hh:964
typename X::field_type field_type
The field type of the preconditioner.
Definition: seqsolverbackend.hh:984
typename std::decay_t< M > matrix_type
The matrix type the preconditioner is for.
Definition: seqsolverbackend.hh:978
Y range_type
The range type of the preconditioner.
Definition: seqsolverbackend.hh:982
void pre(X &v, Y &d) final
Definition: seqsolverbackend.hh:998
X domain_type
The domain type of the preconditioner.
Definition: seqsolverbackend.hh:980
void post(X &) final
Definition: seqsolverbackend.hh:1009
Dune::SolverCategory::Category category() const final
Category of the preconditioner (see SolverCategory::Category)
Definition: seqsolverbackend.hh:1012
void apply(X &v, const Y &d) final
Definition: seqsolverbackend.hh:1000
BlockDiagILU0Preconditioner(const M &m, double w=1.0)
Constructor.
Definition: seqsolverbackend.hh:992
A simple ilu0 block diagonal preconditioned BiCGSTABSolver.
Definition: seqsolverbackend.hh:1035
std::string name() const
Definition: seqsolverbackend.hh:1058
const Dune::InverseOperatorResult & result() const
Definition: seqsolverbackend.hh:1053
bool solve(const Matrix &M, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:1041
A simple ilu0 block diagonal preconditioned RestartedGMResSolver.
Definition: seqsolverbackend.hh:1073
std::string name() const
Definition: seqsolverbackend.hh:1097
const Dune::InverseOperatorResult & result() const
Definition: seqsolverbackend.hh:1092
bool solve(const Matrix &M, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:1079
A simple ilu0 block diagonal preconditioner.
Definition: seqsolverbackend.hh:1110
void post(X &v) final
Clean up.
Definition: seqsolverbackend.hh:1209
typename std::decay_t< M > matrix_type
The matrix type the preconditioner is for.
Definition: seqsolverbackend.hh:1133
void apply(X &v, const Y &d) final
Apply one step of the preconditioner to the system A(v)=d.
Definition: seqsolverbackend.hh:1191
Dune::SolverCategory::Category category() const final
Category of the preconditioner (see SolverCategory::Category)
Definition: seqsolverbackend.hh:1219
Y range_type
The range type of the preconditioner.
Definition: seqsolverbackend.hh:1137
void pre(X &v, Y &d) final
Prepare the preconditioner.
Definition: seqsolverbackend.hh:1171
typename X::field_type field_type
The field type of the preconditioner.
Definition: seqsolverbackend.hh:1139
BlockDiagAMGPreconditioner(const LOP &lop, const Criterion &c, const SmootherArgs &sa)
Constructor.
Definition: seqsolverbackend.hh:1149
X domain_type
The domain type of the preconditioner.
Definition: seqsolverbackend.hh:1135
A simple ilu0 block diagonal preconditioned BiCGSTABSolver.
Definition: seqsolverbackend.hh:1241
const Dune::InverseOperatorResult & result() const
Definition: seqsolverbackend.hh:1297
std::string name() const
Definition: seqsolverbackend.hh:1302
bool solve(const Matrix &m, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:1265
Base class for linear solvers.
Definition: solver.hh:37
Scalar 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="")
Construct the solver.
Definition: solver.hh:53
int verbosity() const
the verbosity level
Definition: solver.hh:83