3.1-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/common/version.hh>
36#include <dune/common/hybridutilities.hh>
37
43
44namespace Dumux {
45
63{
64public:
65
66 template<class Preconditioner, class Solver, class SolverInterface, class Matrix, class Vector>
67 static bool solve(const SolverInterface& s, const Matrix& A, Vector& x, const Vector& b,
68 const std::string& modelParamGroup = "")
69 {
70 Preconditioner precond(A, s.precondIter(), s.relaxation());
71
72 // make a linear operator from a matrix
73 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
74 MatrixAdapter linearOperator(A);
75
76 Solver solver(linearOperator, precond, s.residReduction(), s.maxIter(), s.verbosity());
77
78 Vector bTmp(b);
79
80 Dune::InverseOperatorResult result;
81 solver.apply(x, bTmp, result);
82
83 return result.converged;
84 }
85
86 template<class Preconditioner, class Solver, class SolverInterface, class Matrix, class Vector>
87 static bool solveWithGMRes(const SolverInterface& s, const Matrix& A, Vector& x, const Vector& b,
88 const std::string& modelParamGroup = "")
89 {
90 // get the restart threshold
91 const int restartGMRes = getParamFromGroup<double>(modelParamGroup, "LinearSolver.GMResRestart");
92
93 Preconditioner precond(A, s.precondIter(), s.relaxation());
94
95 // make a linear operator from a matrix
96 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
97 MatrixAdapter linearOperator(A);
98
99 Solver solver(linearOperator, precond, s.residReduction(), restartGMRes, s.maxIter(), s.verbosity());
100
101 Vector bTmp(b);
102
103 Dune::InverseOperatorResult result;
104 solver.apply(x, bTmp, result);
105
106 return result.converged;
107 }
108
109 template<class Preconditioner, class Solver, class SolverInterface, class Matrix, class Vector>
110 static bool solveWithILU0Prec(const SolverInterface& s, const Matrix& A, Vector& x, const Vector& b,
111 const std::string& modelParamGroup = "")
112 {
113 Preconditioner precond(A, s.relaxation());
114
115 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
116 MatrixAdapter operatorA(A);
117
118 Solver solver(operatorA, precond, s.residReduction(), s.maxIter(), s.verbosity());
119
120 Vector bTmp(b);
121
122 Dune::InverseOperatorResult result;
123 solver.apply(x, bTmp, result);
124
125 return result.converged;
126 }
127
128 // solve with RestartedGMRes (needs restartGMRes as additional argument)
129 template<class Preconditioner, class Solver, class SolverInterface, class Matrix, class Vector>
130 static bool solveWithILU0PrecGMRes(const SolverInterface& s, const Matrix& A, Vector& x, const Vector& b,
131 const std::string& modelParamGroup = "")
132 {
133 // get the restart threshold
134 const int restartGMRes = getParamFromGroup<int>(modelParamGroup, "LinearSolver.GMResRestart");
135
136 Preconditioner precond(A, s.relaxation());
137
138 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
139 MatrixAdapter operatorA(A);
140
141 Solver solver(operatorA, precond, s.residReduction(), restartGMRes, s.maxIter(), s.verbosity());
142
143 Vector bTmp(b);
144
145 Dune::InverseOperatorResult result;
146 solver.apply(x, bTmp, result);
147
148 return result.converged;
149 }
150};
151
170{
171public:
173
174 // precondBlockLevel: set this to more than one if the matrix to solve is nested multiple times
175 template<int precondBlockLevel = 1, class Matrix, class Vector>
176 bool solve(const Matrix& A, Vector& x, const Vector& b)
177 {
178 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
179 using Solver = Dune::BiCGSTABSolver<Vector>;
180
181 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
182 }
183
184 std::string name() const
185 {
186 return "ILUn preconditioned BiCGSTAB solver";
187 }
188};
189
208{
209public:
211
212 // precondBlockLevel: set this to more than one if the matrix to solve is nested multiple times
213 template<int precondBlockLevel = 1, class Matrix, class Vector>
214 bool solve(const Matrix& A, Vector& x, const Vector& b)
215 {
216 using Preconditioner = Dune::SeqSOR<Matrix, Vector, Vector, precondBlockLevel>;
217 using Solver = Dune::BiCGSTABSolver<Vector>;
218
219 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
220 }
221
222 std::string name() const
223 {
224 return "SOR preconditioned BiCGSTAB solver";
225 }
226};
227
246{
247public:
249
250 // precondBlockLevel: set this to more than one if the matrix to solve is nested multiple times
251 template<int precondBlockLevel = 1, class Matrix, class Vector>
252 bool solve(const Matrix& A, Vector& x, const Vector& b)
253 {
254 using Preconditioner = Dune::SeqSSOR<Matrix, Vector, Vector, precondBlockLevel>;
255 using Solver = Dune::BiCGSTABSolver<Vector>;
256
257 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
258 }
259
260 std::string name() const
261 {
262 return "SSOR preconditioned BiCGSTAB solver";
263 }
264};
265
284{
285public:
287
288 // precondBlockLevel: set this to more than one if the matrix to solve is nested multiple times
289 template<int precondBlockLevel = 1, class Matrix, class Vector>
290 bool solve(const Matrix& A, Vector& x, const Vector& b)
291 {
292 using Preconditioner = Dune::SeqGS<Matrix, Vector, Vector, precondBlockLevel>;
293 using Solver = Dune::BiCGSTABSolver<Vector>;
294
295 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
296 }
297
298 std::string name() const
299 {
300 return "SSOR preconditioned BiCGSTAB solver";
301 }
302};
303
321{
322public:
324
325 // precondBlockLevel: set this to more than one if the matrix to solve is nested multiple times
326 template<int precondBlockLevel = 1, class Matrix, class Vector>
327 bool solve(const Matrix& A, Vector& x, const Vector& b)
328 {
329 using Preconditioner = Dune::SeqJac<Matrix, Vector, Vector, precondBlockLevel>;
330 using Solver = Dune::BiCGSTABSolver<Vector>;
331
332 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
333 }
334
335 std::string name() const
336 {
337 return "Jac preconditioned BiCGSTAB solver";
338 }
339};
340
358{
359public:
361
362 // precondBlockLevel: set this to more than one if the matrix to solve is nested multiple times
363 template<int precondBlockLevel = 1, class Matrix, class Vector>
364 bool solve(const Matrix& A, Vector& x, const Vector& b)
365 {
366 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
367 using Solver = Dune::CGSolver<Vector>;
368
369 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
370 }
371
372 std::string name() const
373 {
374 return "ILUn preconditioned CG solver";
375 }
376};
377
395{
396public:
398
399 // precondBlockLevel: set this to more than one if the matrix to solve is nested multiple times
400 template<int precondBlockLevel = 1, class Matrix, class Vector>
401 bool solve(const Matrix& A, Vector& x, const Vector& b)
402 {
403 using Preconditioner = Dune::SeqSOR<Matrix, Vector, Vector, precondBlockLevel>;
404 using Solver = Dune::CGSolver<Vector>;
405
406 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
407 }
408
409 std::string name() const
410 {
411 return "SOR preconditioned CG solver";
412 }
413};
414
432{
433public:
435
436 // precondBlockLevel: set this to more than one if the matrix to solve is nested multiple times
437 template<int precondBlockLevel = 1, class Matrix, class Vector>
438 bool solve(const Matrix& A, Vector& x, const Vector& b)
439 {
440 using Preconditioner = Dune::SeqSSOR<Matrix, Vector, Vector, precondBlockLevel>;
441 using Solver = Dune::CGSolver<Vector>;
442
443 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
444 }
445
446 std::string name() const
447 {
448 return "SSOR preconditioned CG solver";
449 }
450};
451
469{
470public:
472
473 // precondBlockLevel: set this to more than one if the matrix to solve is nested multiple times
474 template<int precondBlockLevel = 1, class Matrix, class Vector>
475 bool solve(const Matrix& A, Vector& x, const Vector& b)
476 {
477 using Preconditioner = Dune::SeqGS<Matrix, Vector, Vector, precondBlockLevel>;
478 using Solver = Dune::CGSolver<Vector>;
479
480 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
481 }
482
483 std::string name() const
484 {
485 return "GS preconditioned CG solver";
486 }
487};
488
505{
506public:
508
509 // precondBlockLevel: set this to more than one if the matrix to solve is nested multiple times
510 template<int precondBlockLevel = 1, class Matrix, class Vector>
511 bool solve(const Matrix& A, Vector& x, const Vector& b)
512 {
513 using Preconditioner = Dune::SeqJac<Matrix, Vector, Vector, precondBlockLevel>;
514 using Solver = Dune::CGSolver<Vector>;
515
516 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
517 }
518
519 std::string name() const
520 {
521 return "GS preconditioned CG solver";
522 }
523};
524
543{
544public:
546
547 // precondBlockLevel: set this to more than one if the matrix to solve is nested multiple times
548 template<int precondBlockLevel = 1, class Matrix, class Vector>
549 bool solve(const Matrix& A, Vector& x, const Vector& b)
550 {
551 using Preconditioner = Dune::SeqSSOR<Matrix, Vector, Vector, precondBlockLevel>;
552 using Solver = Dune::RestartedGMResSolver<Vector>;
553
554 return IterativePreconditionedSolverImpl::template solveWithGMRes<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
555 }
556
557 std::string name() const
558 {
559 return "SSOR preconditioned GMRes solver";
560 }
561};
562
580{
581public:
583
584 // precondBlockLevel: set this to more than one if the matrix to solve is nested multiple times
585 template<int precondBlockLevel = 1, class Matrix, class Vector>
586 bool solve(const Matrix& A, Vector& x, const Vector& b)
587 {
588 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
589 using Solver = Dune::BiCGSTABSolver<Vector>;
590
591 return IterativePreconditionedSolverImpl::template solveWithILU0Prec<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
592 }
593
594 std::string name() const
595 {
596 return "ILU0 preconditioned BiCGSTAB solver";
597 }
598};
599
616{
617public:
619
620 // precondBlockLevel: set this to more than one if the matrix to solve is nested multiple times
621 template<int precondBlockLevel = 1, class Matrix, class Vector>
622 bool solve(const Matrix& A, Vector& x, const Vector& b)
623 {
624 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
625 using Solver = Dune::CGSolver<Vector>;
626
627 return IterativePreconditionedSolverImpl::template solveWithILU0Prec<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
628 }
629
630 std::string name() const
631 {
632 return "ILU0 preconditioned BiCGSTAB solver";
633 }
634};
635
653{
654public:
656
657 // precondBlockLevel: set this to more than one if the matrix to solve is nested multiple times
658 template<int precondBlockLevel = 1, class Matrix, class Vector>
659 bool solve(const Matrix& A, Vector& x, const Vector& b)
660 {
661 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
662 using Solver = Dune::RestartedGMResSolver<Vector>;
663
664 return IterativePreconditionedSolverImpl::template solveWithILU0PrecGMRes<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
665 }
666
667 std::string name() const
668 {
669 return "ILU0 preconditioned BiCGSTAB solver";
670 }
671};
672
691{
692public:
694
695 // precondBlockLevel: set this to more than one if the matrix to solve is nested multiple times
696 template<int precondBlockLevel = 1, class Matrix, class Vector>
697 bool solve(const Matrix& A, Vector& x, const Vector& b)
698 {
699 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
700 using Solver = Dune::RestartedGMResSolver<Vector>;
701
702 return IterativePreconditionedSolverImpl::template solveWithGMRes<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
703 }
704
705 std::string name() const
706 {
707 return "ILUn preconditioned GMRes solver";
708 }
709};
710
719{
720public:
722
723 // precondBlockLevel: set this to more than one if the matrix to solve is nested multiple times
724 template<int precondBlockLevel = 1, class Matrix, class Vector>
725 bool solve(const Matrix& A, Vector& x, const Vector& b)
726 {
727 Vector rhs(b);
728 Dune::SeqJac<Matrix, Vector, Vector, precondBlockLevel> jac(A, 1, 1.0);
729 jac.pre(x, rhs);
730 jac.apply(x, rhs);
731 jac.post(x);
732 return true;
733 }
734
735 std::string name() const
736 {
737 return "Explicit diagonal matrix solver";
738 }
739};
740
741#if HAVE_SUPERLU
750class SuperLUBackend : public LinearSolver
751{
752public:
754
755 // precondBlockLevel is unused and just here for compatibility with iterative solvers
756 template<int precondBlockLevel = 1, class Matrix, class Vector>
757 bool solve(const Matrix& A, Vector& x, const Vector& b)
758 {
759 static_assert(isBCRSMatrix<Matrix>::value, "SuperLU only works with BCRS matrices!");
760 using BlockType = typename Matrix::block_type;
761 static_assert(BlockType::rows == BlockType::cols, "Matrix block must be quadratic!");
762 constexpr auto blockSize = BlockType::rows;
763
764 Dune::SuperLU<Matrix> solver(A, this->verbosity() > 0);
765
766 Vector bTmp(b);
767 solver.apply(x, bTmp, result_);
768
769 int size = x.size();
770 for (int i = 0; i < size; i++)
771 {
772 for (int j = 0; j < blockSize; j++)
773 {
774 using std::isnan;
775 using std::isinf;
776 if (isnan(x[i][j]) || isinf(x[i][j]))
777 {
778 result_.converged = false;
779 break;
780 }
781 }
782 }
783
784 return result_.converged;
785 }
786
787 std::string name() const
788 {
789 return "SuperLU solver";
790 }
791
792 const Dune::InverseOperatorResult& result() const
793 {
794 return result_;
795 }
796
797private:
798 Dune::InverseOperatorResult result_;
799};
800#endif // HAVE_SUPERLU
801
802#if HAVE_UMFPACK
811class UMFPackBackend : public LinearSolver
812{
813public:
815
816 // precondBlockLevel is unused and just here for compatibility with iterative solvers
817 template<int precondBlockLevel = 1, class Matrix, class Vector>
818 bool solve(const Matrix& A, Vector& x, const Vector& b)
819 {
820 static_assert(isBCRSMatrix<Matrix>::value, "UMFPack only works with BCRS matrices!");
821 using BlockType = typename Matrix::block_type;
822 static_assert(BlockType::rows == BlockType::cols, "Matrix block must be quadratic!");
823 constexpr auto blockSize = BlockType::rows;
824
825 Dune::UMFPack<Matrix> solver(A, this->verbosity() > 0);
826
827 Vector bTmp(b);
828 solver.apply(x, bTmp, result_);
829
830 int size = x.size();
831 for (int i = 0; i < size; i++)
832 {
833 for (int j = 0; j < blockSize; j++)
834 {
835 using std::isnan;
836 using std::isinf;
837 if (isnan(x[i][j]) || isinf(x[i][j]))
838 {
839 result_.converged = false;
840 break;
841 }
842 }
843 }
844
845 return result_.converged;
846 }
847
848 std::string name() const
849 {
850 return "UMFPack solver";
851 }
852
853 const Dune::InverseOperatorResult& result() const
854 {
855 return result_;
856 }
857
858private:
859 Dune::InverseOperatorResult result_;
860};
861#endif // HAVE_UMFPACK
862
863
867// \{
868
873template<class M, class X, class Y, int blockLevel = 2>
874class BlockDiagILU0Preconditioner : public Dune::Preconditioner<X, Y>
875{
876 template<std::size_t i>
877 using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>;
878
879 template<std::size_t i>
880 using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>;
881
882#if DUNE_VERSION_NEWER(DUNE_ISTL,2,6)
883 template<std::size_t i>
884 using BlockILU = Dune::SeqILU<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>, blockLevel-1>;
885#else
886 template<std::size_t i>
887 using BlockILU = Dune::SeqILU0<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>, blockLevel-1>;
888#endif
889
890 using ILUTuple = typename makeFromIndexedType<std::tuple, BlockILU, std::make_index_sequence<M::size()> >::type;
891
892public:
894 using matrix_type = typename std::decay_t<M>;
896 using domain_type = X;
898 using range_type = Y;
900 using field_type = typename X::field_type;
901
908 BlockDiagILU0Preconditioner(const M& m, double w = 1.0)
909 : BlockDiagILU0Preconditioner(m, w, std::make_index_sequence<M::size()>{})
910 {
911 static_assert(blockLevel >= 2, "Only makes sense for MultiTypeBlockMatrix!");
912 }
913
914 void pre (X& v, Y& d) final {}
915
916 void apply (X& v, const Y& d) final
917 {
918 using namespace Dune::Hybrid;
919 forEach(integralRange(Dune::Hybrid::size(ilu_)), [&](const auto i)
920 {
921 std::get<decltype(i)::value>(ilu_).apply(v[i], d[i]);
922 });
923 }
924
925 void post (X&) final {}
926
928 Dune::SolverCategory::Category category() const final
929 {
930 return Dune::SolverCategory::sequential;
931 }
932
933private:
934 template<std::size_t... Is>
935 BlockDiagILU0Preconditioner (const M& m, double w, std::index_sequence<Is...> is)
936 : ilu_(std::make_tuple(BlockILU<Is>(m[Dune::index_constant<Is>{}][Dune::index_constant<Is>{}], w)...))
937 {}
938
939 ILUTuple ilu_;
940};
941
942
951{
952
953public:
955
956 // Solve saddle-point problem using a Schur complement based preconditioner
957 template<int precondBlockLevel = 2, class Matrix, class Vector>
958 bool solve(const Matrix& M, Vector& x, const Vector& b)
959 {
961 Dune::MatrixAdapter<Matrix, Vector, Vector> op(M);
962 Dune::BiCGSTABSolver<Vector> solver(op, preconditioner, this->residReduction(),
963 this->maxIter(), this->verbosity());
964 auto bTmp(b);
965 solver.apply(x, bTmp, result_);
966
967 return result_.converged;
968 }
969
970 const Dune::InverseOperatorResult& result() const
971 {
972 return result_;
973 }
974
975 std::string name() const
976 { return "block-diagonal ILU0 preconditioned BiCGSTAB solver"; }
977
978private:
979 Dune::InverseOperatorResult result_;
980};
981
986template<class M, class X, class Y, int blockLevel = 2>
987class BlockDiagAMGPreconditioner : public Dune::Preconditioner<X, Y>
988{
989 template<std::size_t i>
990 using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>;
991
992 template<std::size_t i>
993 using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>;
994
995 template<std::size_t i>
996 using Smoother = Dune::SeqSSOR<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>>;
997
998 template<std::size_t i>
999 using LinearOperator = Dune::MatrixAdapter<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>>;
1000
1001 template<std::size_t i>
1002 using ScalarProduct = Dune::SeqScalarProduct<VecBlockType<i>>;
1003
1004 template<std::size_t i>
1005 using BlockAMG = Dune::Amg::AMG<LinearOperator<i>, VecBlockType<i>, Smoother<i>>;
1006
1007 using AMGTuple = typename makeFromIndexedType<std::tuple, BlockAMG, std::make_index_sequence<M::size()> >::type;
1008
1009public:
1011 using matrix_type = typename std::decay_t<M>;
1013 using domain_type = X;
1015 using range_type = Y;
1017 using field_type = typename X::field_type;
1018
1026 template<class LOP, class Criterion, class SmootherArgs>
1027 BlockDiagAMGPreconditioner(const LOP& lop, const Criterion& c, const SmootherArgs& sa)
1028 : BlockDiagAMGPreconditioner(lop, c, sa, std::make_index_sequence<M::size()>{})
1029 {
1030 static_assert(blockLevel >= 2, "Only makes sense for MultiTypeBlockMatrix!");
1031 }
1032
1049 void pre (X& v, Y& d) final
1050 {
1051 using namespace Dune::Hybrid;
1052 forEach(integralRange(Dune::Hybrid::size(amg_)), [&](const auto i)
1053 {
1054 std::get<decltype(i)::value>(amg_).pre(v[i], d[i]);
1055 });
1056 }
1057
1069 void apply (X& v, const Y& d) final
1070 {
1071 using namespace Dune::Hybrid;
1072 forEach(integralRange(Dune::Hybrid::size(amg_)), [&](const auto i)
1073 {
1074 std::get<decltype(i)::value>(amg_).apply(v[i], d[i]);
1075 });
1076 }
1077
1087 void post (X& v) final
1088 {
1089 using namespace Dune::Hybrid;
1090 forEach(integralRange(Dune::Hybrid::size(amg_)), [&](const auto i)
1091 {
1092 std::get<decltype(i)::value>(amg_).post(v[i]);
1093 });
1094 }
1095
1097 Dune::SolverCategory::Category category() const final
1098 {
1099 return Dune::SolverCategory::sequential;
1100 }
1101
1102private:
1103 template<class LOP, class Criterion, class SmootherArgs, std::size_t... Is>
1104 BlockDiagAMGPreconditioner (const LOP& lop, const Criterion& c, const SmootherArgs& sa, std::index_sequence<Is...> is)
1105 : amg_(std::make_tuple(BlockAMG<Is>(*std::get<Is>(lop), *std::get<Is>(c), *std::get<Is>(sa))...))
1106 {}
1107
1108 AMGTuple amg_;
1109};
1110
1119{
1120 template<class M, std::size_t i>
1121 using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>;
1122
1123 template<class X, std::size_t i>
1124 using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>;
1125
1126 template<class M, class X, std::size_t i>
1127 using Smoother = Dune::SeqSSOR<DiagBlockType<M, i>, VecBlockType<X, i>, VecBlockType<X, i>>;
1128
1129 template<class M, class X, std::size_t i>
1130 using SmootherArgs = typename Dune::Amg::SmootherTraits<Smoother<M, X, i>>::Arguments;
1131
1132 template<class M, std::size_t i>
1133 using Criterion = Dune::Amg::CoarsenCriterion<Dune::Amg::SymmetricCriterion<DiagBlockType<M, i>, Dune::Amg::FirstDiagonal>>;
1134
1135 template<class M, class X, std::size_t i>
1136 using LinearOperator = Dune::MatrixAdapter<DiagBlockType<M, i>, VecBlockType<X, i>, VecBlockType<X, i>>;
1137
1138public:
1140
1141 // Solve saddle-point problem using a Schur complement based preconditioner
1142 template<int precondBlockLevel = 2, class Matrix, class Vector>
1143 bool solve(const Matrix& m, Vector& x, const Vector& b)
1144 {
1147 Dune::Amg::Parameters params(15, 2000, 1.2, 1.6, Dune::Amg::atOnceAccu);
1148 params.setDefaultValuesIsotropic(3);
1149 params.setDebugLevel(this->verbosity());
1150
1151 auto criterion = makeCriterion_<Criterion, Matrix>(params, std::make_index_sequence<Matrix::size()>{});
1152 auto smootherArgs = makeSmootherArgs_<SmootherArgs, Matrix, Vector>(std::make_index_sequence<Matrix::size()>{});
1153
1154 using namespace Dune::Hybrid;
1155 forEach(integralRange(Dune::Hybrid::size(m)), [&](const auto i)
1156 {
1157 auto& args = std::get<decltype(i)::value>(smootherArgs);
1158 args->iterations = 1;
1159 args->relaxationFactor = 1;
1160 });
1161
1162 auto linearOperator = makeLinearOperator_<LinearOperator, Matrix, Vector>(m, std::make_index_sequence<Matrix::size()>{});
1163
1164 BlockDiagAMGPreconditioner<Matrix, Vector, Vector> preconditioner(linearOperator, criterion, smootherArgs);
1165
1166 Dune::MatrixAdapter<Matrix, Vector, Vector> op(m);
1167 Dune::BiCGSTABSolver<Vector> solver(op, preconditioner, this->residReduction(),
1168 this->maxIter(), this->verbosity());
1169 auto bTmp(b);
1170 solver.apply(x, bTmp, result_);
1171
1172 return result_.converged;
1173 }
1174
1175 const Dune::InverseOperatorResult& result() const
1176 {
1177 return result_;
1178 }
1179
1180 std::string name() const
1181 { return "block-diagonal AMG preconditioned BiCGSTAB solver"; }
1182
1183private:
1184 template<template<class M, std::size_t i> class Criterion, class Matrix, class Params, std::size_t... Is>
1185 auto makeCriterion_(const Params& p, std::index_sequence<Is...>)
1186 {
1187 return std::make_tuple(std::make_shared<Criterion<Matrix, Is>>(p)...);
1188 }
1189
1190 template<template<class M, class X, std::size_t i> class SmootherArgs, class Matrix, class Vector, std::size_t... Is>
1191 auto makeSmootherArgs_(std::index_sequence<Is...>)
1192 {
1193 return std::make_tuple(std::make_shared<SmootherArgs<Matrix, Vector, Is>>()...);
1194 }
1195
1196 template<template<class M, class X, std::size_t i> class LinearOperator, class Matrix, class Vector, std::size_t... Is>
1197 auto makeLinearOperator_(const Matrix& m, std::index_sequence<Is...>)
1198 {
1199 return std::make_tuple(std::make_shared<LinearOperator<Matrix, Vector, Is>>(m[Dune::index_constant<Is>{}][Dune::index_constant<Is>{}])...);
1200 }
1201
1202 Dune::InverseOperatorResult result_;
1203};
1204
1205// \}
1206
1207} // end namespace Dumux
1208
1209#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.
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
Definition: common/properties/model.hh:34
Definition: utility.hh:40
A general solver backend allowing arbitrary preconditioners and solvers.
Definition: seqsolverbackend.hh:63
static bool solveWithILU0PrecGMRes(const SolverInterface &s, const Matrix &A, Vector &x, const Vector &b, const std::string &modelParamGroup="")
Definition: seqsolverbackend.hh:130
static bool solveWithGMRes(const SolverInterface &s, const Matrix &A, Vector &x, const Vector &b, const std::string &modelParamGroup="")
Definition: seqsolverbackend.hh:87
static bool solveWithILU0Prec(const SolverInterface &s, const Matrix &A, Vector &x, const Vector &b, const std::string &modelParamGroup="")
Definition: seqsolverbackend.hh:110
static bool solve(const SolverInterface &s, const Matrix &A, Vector &x, const Vector &b, const std::string &modelParamGroup="")
Definition: seqsolverbackend.hh:67
Sequential ILU(n)-preconditioned BiCSTAB solver.
Definition: seqsolverbackend.hh:170
std::string name() const
Definition: seqsolverbackend.hh:184
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:176
Sequential SOR-preconditioned BiCSTAB solver.
Definition: seqsolverbackend.hh:208
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:214
std::string name() const
Definition: seqsolverbackend.hh:222
Sequential SSOR-preconditioned BiCGSTAB solver.
Definition: seqsolverbackend.hh:246
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:252
std::string name() const
Definition: seqsolverbackend.hh:260
Sequential GS-preconditioned BiCGSTAB solver.
Definition: seqsolverbackend.hh:284
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:290
std::string name() const
Definition: seqsolverbackend.hh:298
Sequential Jacobi-preconditioned BiCSTAB solver.
Definition: seqsolverbackend.hh:321
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:327
std::string name() const
Definition: seqsolverbackend.hh:335
Sequential ILU(n)-preconditioned CG solver.
Definition: seqsolverbackend.hh:358
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:364
std::string name() const
Definition: seqsolverbackend.hh:372
Sequential SOR-preconditioned CG solver.
Definition: seqsolverbackend.hh:395
std::string name() const
Definition: seqsolverbackend.hh:409
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:401
Sequential SSOR-preconditioned CG solver.
Definition: seqsolverbackend.hh:432
std::string name() const
Definition: seqsolverbackend.hh:446
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:438
Sequential GS-preconditioned CG solver.
Definition: seqsolverbackend.hh:469
std::string name() const
Definition: seqsolverbackend.hh:483
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:475
Sequential Jacobi-preconditioned CG solver.
Definition: seqsolverbackend.hh:505
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:511
std::string name() const
Definition: seqsolverbackend.hh:519
Sequential SSOR-preconditioned GMRes solver.
Definition: seqsolverbackend.hh:543
std::string name() const
Definition: seqsolverbackend.hh:557
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:549
Sequential ILU(0)-preconditioned BiCGSTAB solver.
Definition: seqsolverbackend.hh:580
std::string name() const
Definition: seqsolverbackend.hh:594
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:586
Sequential ILU(0)-preconditioned CG solver.
Definition: seqsolverbackend.hh:616
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:622
std::string name() const
Definition: seqsolverbackend.hh:630
Sequential ILU0-preconditioned GMRes solver.
Definition: seqsolverbackend.hh:653
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:659
std::string name() const
Definition: seqsolverbackend.hh:667
Sequential ILU(n)-preconditioned GMRes solver.
Definition: seqsolverbackend.hh:691
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:697
std::string name() const
Definition: seqsolverbackend.hh:705
Solver for simple block-diagonal matrices (e.g. from explicit time stepping schemes)
Definition: seqsolverbackend.hh:719
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:725
std::string name() const
Definition: seqsolverbackend.hh:735
A simple ilu0 block diagonal preconditioner.
Definition: seqsolverbackend.hh:875
typename X::field_type field_type
The field type of the preconditioner.
Definition: seqsolverbackend.hh:900
typename std::decay_t< M > matrix_type
The matrix type the preconditioner is for.
Definition: seqsolverbackend.hh:894
Y range_type
The range type of the preconditioner.
Definition: seqsolverbackend.hh:898
void pre(X &v, Y &d) final
Definition: seqsolverbackend.hh:914
X domain_type
The domain type of the preconditioner.
Definition: seqsolverbackend.hh:896
void post(X &) final
Definition: seqsolverbackend.hh:925
Dune::SolverCategory::Category category() const final
Category of the preconditioner (see SolverCategory::Category)
Definition: seqsolverbackend.hh:928
void apply(X &v, const Y &d) final
Definition: seqsolverbackend.hh:916
BlockDiagILU0Preconditioner(const M &m, double w=1.0)
Constructor.
Definition: seqsolverbackend.hh:908
A simple ilu0 block diagonal preconditioned BiCGSTABSolver.
Definition: seqsolverbackend.hh:951
std::string name() const
Definition: seqsolverbackend.hh:975
const Dune::InverseOperatorResult & result() const
Definition: seqsolverbackend.hh:970
bool solve(const Matrix &M, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:958
A simple ilu0 block diagonal preconditioner.
Definition: seqsolverbackend.hh:988
void post(X &v) final
Clean up.
Definition: seqsolverbackend.hh:1087
typename std::decay_t< M > matrix_type
The matrix type the preconditioner is for.
Definition: seqsolverbackend.hh:1011
void apply(X &v, const Y &d) final
Apply one step of the preconditioner to the system A(v)=d.
Definition: seqsolverbackend.hh:1069
Dune::SolverCategory::Category category() const final
Category of the preconditioner (see SolverCategory::Category)
Definition: seqsolverbackend.hh:1097
Y range_type
The range type of the preconditioner.
Definition: seqsolverbackend.hh:1015
void pre(X &v, Y &d) final
Prepare the preconditioner.
Definition: seqsolverbackend.hh:1049
typename X::field_type field_type
The field type of the preconditioner.
Definition: seqsolverbackend.hh:1017
BlockDiagAMGPreconditioner(const LOP &lop, const Criterion &c, const SmootherArgs &sa)
Constructor.
Definition: seqsolverbackend.hh:1027
X domain_type
The domain type of the preconditioner.
Definition: seqsolverbackend.hh:1013
A simple ilu0 block diagonal preconditioned BiCGSTABSolver.
Definition: seqsolverbackend.hh:1119
bool solve(const Matrix &m, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:1143
const Dune::InverseOperatorResult & result() const
Definition: seqsolverbackend.hh:1175
std::string name() const
Definition: seqsolverbackend.hh:1180
Base class for linear solvers.
Definition: dumux/linear/solver.hh:37
double residReduction() const
the linear solver residual reduction
Definition: dumux/linear/solver.hh:97
int maxIter() const
the maximum number of linear solver iterations
Definition: dumux/linear/solver.hh:89
const std::string & paramGroup() const
the parameter group for getting parameter from the parameter tree
Definition: dumux/linear/solver.hh:77
LinearSolver(const std::string &paramGroup="")
Contruct the solver.
Definition: dumux/linear/solver.hh:52
int verbosity() const
the verbosity level
Definition: dumux/linear/solver.hh:81
Base class for linear solvers.