version 3.7
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_SEQ_SOLVER_BACKEND_HH
13#define DUMUX_SEQ_SOLVER_BACKEND_HH
14
15#include <type_traits>
16#include <tuple>
17#include <utility>
18
19#include <dune/istl/preconditioners.hh>
20#include <dune/istl/solvers.hh>
21#include <dune/istl/superlu.hh>
22#include <dune/istl/umfpack.hh>
23#include <dune/istl/io.hh>
24#include <dune/common/indices.hh>
25#include <dune/common/hybridutilities.hh>
26
34
35namespace Dumux {
36
54{
55public:
56
57 template<class Preconditioner, class Solver, class SolverInterface, class Matrix, class Vector>
58 static bool solve(const SolverInterface& s, const Matrix& A, Vector& x, const Vector& b,
59 const std::string& modelParamGroup = "")
60 {
61 Preconditioner precond(A, s.precondIter(), s.relaxation());
62
63 // make a linear operator from a matrix
64 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
65 MatrixAdapter linearOperator(A);
66
67 Solver solver(linearOperator, precond, s.residReduction(), s.maxIter(), s.verbosity());
68
69 Vector bTmp(b);
70
71 Dune::InverseOperatorResult result;
72 solver.apply(x, bTmp, result);
73
74 return result.converged;
75 }
76
77 template<class Preconditioner, class Solver, class SolverInterface, class Matrix, class Vector>
78 static bool solveWithGMRes(const SolverInterface& s, const Matrix& A, Vector& x, const Vector& b,
79 const std::string& modelParamGroup = "")
80 {
81 // get the restart threshold
82 const int restartGMRes = getParamFromGroup<int>(modelParamGroup, "LinearSolver.GMResRestart", 10);
83
84 Preconditioner precond(A, s.precondIter(), s.relaxation());
85
86 // make a linear operator from a matrix
87 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
88 MatrixAdapter linearOperator(A);
89
90 Solver solver(linearOperator, precond, s.residReduction(), restartGMRes, s.maxIter(), s.verbosity());
91
92 Vector bTmp(b);
93
94 Dune::InverseOperatorResult result;
95 solver.apply(x, bTmp, result);
96
97 return result.converged;
98 }
99
100 template<class Preconditioner, class Solver, class SolverInterface, class Matrix, class Vector>
101 static bool solveWithILU0Prec(const SolverInterface& s, const Matrix& A, Vector& x, const Vector& b,
102 const std::string& modelParamGroup = "")
103 {
104 Preconditioner precond(A, s.relaxation());
105
106 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
107 MatrixAdapter operatorA(A);
108
109 Solver solver(operatorA, precond, s.residReduction(), s.maxIter(), s.verbosity());
110
111 Vector bTmp(b);
112
113 Dune::InverseOperatorResult result;
114 solver.apply(x, bTmp, result);
115
116 return result.converged;
117 }
118
119 // solve with RestartedGMRes (needs restartGMRes as additional argument)
120 template<class Preconditioner, class Solver, class SolverInterface, class Matrix, class Vector>
121 static bool solveWithILU0PrecGMRes(const SolverInterface& s, const Matrix& A, Vector& x, const Vector& b,
122 const std::string& modelParamGroup = "")
123 {
124 // get the restart threshold
125 const int restartGMRes = getParamFromGroup<int>(modelParamGroup, "LinearSolver.GMResRestart", 10);
126
127 Preconditioner precond(A, s.relaxation());
128
129 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
130 MatrixAdapter operatorA(A);
131
132 Solver solver(operatorA, precond, s.residReduction(), restartGMRes, s.maxIter(), s.verbosity());
133
134 Vector bTmp(b);
135
136 Dune::InverseOperatorResult result;
137 solver.apply(x, bTmp, result);
138
139 return result.converged;
140 }
141
142 // solve with generic parameter tree
143 template<class Preconditioner, class Solver, class Matrix, class Vector>
144 static bool solveWithParamTree(const Matrix& A, Vector& x, const Vector& b,
145 const Dune::ParameterTree& params)
146 {
147 // make a linear operator from a matrix
148 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
149 const auto linearOperator = std::make_shared<MatrixAdapter>(A);
150
151 auto precond = std::make_shared<Preconditioner>(linearOperator, params.sub("preconditioner"));
152 Solver solver(linearOperator, precond, params);
153
154 Vector bTmp(b);
155
156 Dune::InverseOperatorResult result;
157 solver.apply(x, bTmp, result);
158
159 return result.converged;
160 }
161};
162
169template<class M>
170constexpr std::size_t preconditionerBlockLevel() noexcept
171{
173}
174
192class [[deprecated("Removed after 3.7. Use solver from istlsolvers.hh")]] ILUnBiCGSTABBackend : public LinearSolver
193{
194public:
196
197 template<class Matrix, class Vector>
198 bool solve(const Matrix& A, Vector& x, const Vector& b)
199 {
200 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
201 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
202 using Solver = Dune::BiCGSTABSolver<Vector>;
203
204 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
205 }
206
207 std::string name() const
208 {
209 return "ILUn preconditioned BiCGSTAB solver";
210 }
211};
212
230class [[deprecated("Removed after 3.7. Use solver from istlsolvers.hh")]] SORBiCGSTABBackend : public LinearSolver
231{
232public:
234
235 template<class Matrix, class Vector>
236 bool solve(const Matrix& A, Vector& x, const Vector& b)
237 {
238 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
239 using Preconditioner = Dune::SeqSOR<Matrix, Vector, Vector, precondBlockLevel>;
240 using Solver = Dune::BiCGSTABSolver<Vector>;
241
242 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
243 }
244
245 std::string name() const
246 {
247 return "SOR preconditioned BiCGSTAB solver";
248 }
249};
250
268class [[deprecated("Removed after 3.7. Use solver from istlsolvers.hh")]] SSORBiCGSTABBackend : public LinearSolver
269{
270public:
272
273 template<class Matrix, class Vector>
274 bool solve(const Matrix& A, Vector& x, const Vector& b)
275 {
276 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
277 using Preconditioner = Dune::SeqSSOR<Matrix, Vector, Vector, precondBlockLevel>;
278 using Solver = Dune::BiCGSTABSolver<Vector>;
279
280 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
281 }
282
283 std::string name() const
284 {
285 return "SSOR preconditioned BiCGSTAB solver";
286 }
287};
288
306class [[deprecated("Removed after 3.7. Use solver from istlsolvers.hh")]] GSBiCGSTABBackend : public LinearSolver
307{
308public:
310
311 template<class Matrix, class Vector>
312 bool solve(const Matrix& A, Vector& x, const Vector& b)
313 {
314 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
315 using Preconditioner = Dune::SeqGS<Matrix, Vector, Vector, precondBlockLevel>;
316 using Solver = Dune::BiCGSTABSolver<Vector>;
317
318 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
319 }
320
321 std::string name() const
322 {
323 return "SSOR preconditioned BiCGSTAB solver";
324 }
325};
326
343class [[deprecated("Removed after 3.7. Use solver from istlsolvers.hh")]] JacBiCGSTABBackend : public LinearSolver
344{
345public:
347
348 template<class Matrix, class Vector>
349 bool solve(const Matrix& A, Vector& x, const Vector& b)
350 {
351 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
352 using Preconditioner = Dune::SeqJac<Matrix, Vector, Vector, precondBlockLevel>;
353 using Solver = Dune::BiCGSTABSolver<Vector>;
354
355 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
356 }
357
358 std::string name() const
359 {
360 return "Jac preconditioned BiCGSTAB solver";
361 }
362};
363
380class [[deprecated("Removed after 3.7. Use solver from istlsolvers.hh")]] ILUnCGBackend : public LinearSolver
381{
382public:
384
385 template<class Matrix, class Vector>
386 bool solve(const Matrix& A, Vector& x, const Vector& b)
387 {
388 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
389 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
390 using Solver = Dune::CGSolver<Vector>;
391
392 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
393 }
394
395 std::string name() const
396 {
397 return "ILUn preconditioned CG solver";
398 }
399};
400
417class [[deprecated("Removed after 3.7. Use solver from istlsolvers.hh")]] SORCGBackend : public LinearSolver
418{
419public:
421
422 template<class Matrix, class Vector>
423 bool solve(const Matrix& A, Vector& x, const Vector& b)
424 {
425 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
426 using Preconditioner = Dune::SeqSOR<Matrix, Vector, Vector, precondBlockLevel>;
427 using Solver = Dune::CGSolver<Vector>;
428
429 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
430 }
431
432 std::string name() const
433 {
434 return "SOR preconditioned CG solver";
435 }
436};
437
454class [[deprecated("Removed after 3.7. Use solver from istlsolvers.hh")]] SSORCGBackend : public LinearSolver
455{
456public:
458
459 template<class Matrix, class Vector>
460 bool solve(const Matrix& A, Vector& x, const Vector& b)
461 {
462 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
463 using Preconditioner = Dune::SeqSSOR<Matrix, Vector, Vector, precondBlockLevel>;
464 using Solver = Dune::CGSolver<Vector>;
465
466 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
467 }
468
469 std::string name() const
470 {
471 return "SSOR preconditioned CG solver";
472 }
473};
474
491class [[deprecated("Removed after 3.7. Use solver from istlsolvers.hh")]] GSCGBackend : public LinearSolver
492{
493public:
495
496 template<class Matrix, class Vector>
497 bool solve(const Matrix& A, Vector& x, const Vector& b)
498 {
499 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
500 using Preconditioner = Dune::SeqGS<Matrix, Vector, Vector, precondBlockLevel>;
501 using Solver = Dune::CGSolver<Vector>;
502
503 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
504 }
505
506 std::string name() const
507 {
508 return "GS preconditioned CG solver";
509 }
510};
511
527class [[deprecated("Removed after 3.7. Use solver from istlsolvers.hh")]] JacCGBackend : public LinearSolver
528{
529public:
531
532 template<class Matrix, class Vector>
533 bool solve(const Matrix& A, Vector& x, const Vector& b)
534 {
535 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
536 using Preconditioner = Dune::SeqJac<Matrix, Vector, Vector, precondBlockLevel>;
537 using Solver = Dune::CGSolver<Vector>;
538
539 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
540 }
541
542 std::string name() const
543 {
544 return "GS preconditioned CG solver";
545 }
546};
547
565class [[deprecated("Removed after 3.7. Use solver from istlsolvers.hh")]] SSORRestartedGMResBackend : public LinearSolver
566{
567public:
569
570 template<class Matrix, class Vector>
571 bool solve(const Matrix& A, Vector& x, const Vector& b)
572 {
573 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
574 using Preconditioner = Dune::SeqSSOR<Matrix, Vector, Vector, precondBlockLevel>;
575 using Solver = Dune::RestartedGMResSolver<Vector>;
576
577 return IterativePreconditionedSolverImpl::template solveWithGMRes<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
578 }
579
580 std::string name() const
581 {
582 return "SSOR preconditioned GMRes solver";
583 }
584};
585
602class [[deprecated("Removed after 3.7. Use solver from istlsolvers.hh")]] ILU0BiCGSTABBackend : public LinearSolver
603{
604public:
606
607 template<class Matrix, class Vector>
608 bool solve(const Matrix& A, Vector& x, const Vector& b)
609 {
610 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
611 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
612 using Solver = Dune::BiCGSTABSolver<Vector>;
613
614 return IterativePreconditionedSolverImpl::template solveWithILU0Prec<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
615 }
616
617 std::string name() const
618 {
619 return "ILU0 preconditioned BiCGSTAB solver";
620 }
621};
622
638class [[deprecated("Removed after 3.7. Use solver from istlsolvers.hh")]] ILU0CGBackend : public LinearSolver
639{
640public:
642
643 template<class Matrix, class Vector>
644 bool solve(const Matrix& A, Vector& x, const Vector& b)
645 {
646 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
647 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
648 using Solver = Dune::CGSolver<Vector>;
649
650 return IterativePreconditionedSolverImpl::template solveWithILU0Prec<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
651 }
652
653 std::string name() const
654 {
655 return "ILU0 preconditioned BiCGSTAB solver";
656 }
657};
658
675class [[deprecated("Removed after 3.7. Use solver from istlsolvers.hh")]] ILU0RestartedGMResBackend : public LinearSolver
676{
677public:
679
680 template<class Matrix, class Vector>
681 bool solve(const Matrix& A, Vector& x, const Vector& b)
682 {
683 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
684 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
685 using Solver = Dune::RestartedGMResSolver<Vector>;
686
687 return IterativePreconditionedSolverImpl::template solveWithILU0PrecGMRes<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
688 }
689
690 std::string name() const
691 {
692 return "ILU0 preconditioned BiCGSTAB solver";
693 }
694};
695
713class [[deprecated("Removed after 3.7. Use solver from istlsolvers.hh")]] ILUnRestartedGMResBackend : public LinearSolver
714{
715public:
717
718 template<class Matrix, class Vector>
719 bool solve(const Matrix& A, Vector& x, const Vector& b)
720 {
721 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
722 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
723 using Solver = Dune::RestartedGMResSolver<Vector>;
724
725 return IterativePreconditionedSolverImpl::template solveWithGMRes<Preconditioner, Solver>(*this, A, x, b, this->paramGroup());
726 }
727
728 std::string name() const
729 {
730 return "ILUn preconditioned GMRes solver";
731 }
732};
733
742{
743public:
745
746 template<class Matrix, class Vector>
747 bool solve(const Matrix& A, Vector& x, const Vector& b)
748 {
749 Vector rhs(b);
750 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
751 Dune::SeqJac<Matrix, Vector, Vector, precondBlockLevel> jac(A, 1, 1.0);
752 jac.pre(x, rhs);
753 jac.apply(x, rhs);
754 jac.post(x);
755 return true;
756 }
757
758 std::string name() const
759 {
760 return "Explicit diagonal matrix solver";
761 }
762};
763
764#if HAVE_SUPERLU
773class [[deprecated("Removed after 3.7. Use solver from istlsolvers.hh")]] SuperLUBackend : public LinearSolver
774{
775public:
777
778 template<class Matrix, class Vector>
779 bool solve(const Matrix& A, Vector& x, const Vector& b)
780 {
781 static_assert(isBCRSMatrix<Matrix>::value, "SuperLU only works with BCRS matrices!");
782 using BlockType = typename Matrix::block_type;
783 static_assert(BlockType::rows == BlockType::cols, "Matrix block must be quadratic!");
784 constexpr auto blockSize = BlockType::rows;
785
786 Dune::SuperLU<Matrix> solver(A, this->verbosity() > 0);
787
788 Vector bTmp(b);
789 solver.apply(x, bTmp, result_);
790
791 int size = x.size();
792 for (int i = 0; i < size; i++)
793 {
794 for (int j = 0; j < blockSize; j++)
795 {
796 using std::isnan;
797 using std::isinf;
798 if (isnan(x[i][j]) || isinf(x[i][j]))
799 {
800 result_.converged = false;
801 break;
802 }
803 }
804 }
805
806 return result_.converged;
807 }
808
809 std::string name() const
810 {
811 return "SuperLU solver";
812 }
813
814 const Dune::InverseOperatorResult& result() const
815 {
816 return result_;
817 }
818
819private:
820 Dune::InverseOperatorResult result_;
821};
822#endif // HAVE_SUPERLU
823
824#if HAVE_UMFPACK
847class [[deprecated("Removed after 3.7. Use solver from istlsolvers.hh")]] UMFPackBackend : public LinearSolver
848{
849public:
850
851 UMFPackBackend(const std::string& paramGroup = "")
852 : LinearSolver(paramGroup)
853 {
854 ordering_ = getParamFromGroup<int>(this->paramGroup(), "LinearSolver.UMFPackOrdering", 1);
855 }
856
858 void setOrdering(int i)
859 { ordering_ = i; }
860
862 int ordering() const
863 { return ordering_; }
864
865 template<class Matrix, class Vector>
866 bool solve(const Matrix& A, Vector& x, const Vector& b)
867 {
868 static_assert(isBCRSMatrix<Matrix>::value, "UMFPack only works with BCRS matrices!");
869 using BlockType = typename Matrix::block_type;
870 static_assert(BlockType::rows == BlockType::cols, "Matrix block must be quadratic!");
871 constexpr auto blockSize = BlockType::rows;
872
873 Dune::UMFPack<Matrix> solver;
874 solver.setVerbosity(this->verbosity() > 0);
875 solver.setOption(UMFPACK_ORDERING, ordering_);
876 solver.setMatrix(A);
877
878 Vector bTmp(b);
879 solver.apply(x, bTmp, result_);
880
881 int size = x.size();
882 for (int i = 0; i < size; i++)
883 {
884 for (int j = 0; j < blockSize; j++)
885 {
886 using std::isnan;
887 using std::isinf;
888 if (isnan(x[i][j]) || isinf(x[i][j]))
889 {
890 result_.converged = false;
891 break;
892 }
893 }
894 }
895
896 return result_.converged;
897 }
898
899 std::string name() const
900 {
901 return "UMFPack solver";
902 }
903
904 const Dune::InverseOperatorResult& result() const
905 {
906 return result_;
907 }
908
909private:
910 Dune::InverseOperatorResult result_;
911 int ordering_;
912};
913#endif // HAVE_UMFPACK
914
915
919// \{
920
925template <class LinearSolverTraits>
927{
928public:
930
931 template<class Matrix, class Vector>
932 bool solve(const Matrix& A, Vector& x, const Vector& b)
933 {
934 using Preconditioner = SeqUzawa<Matrix, Vector, Vector>;
935 using Solver = Dune::BiCGSTABSolver<Vector>;
936 static const auto solverParams = LinearSolverParameters<LinearSolverTraits>::createParameterTree(this->paramGroup());
937 return IterativePreconditionedSolverImpl::template solveWithParamTree<Preconditioner, Solver>(A, x, b, solverParams);
938 }
939
940 std::string name() const
941 {
942 return "Uzawa preconditioned BiCGSTAB solver";
943 }
944};
945
950template<class M, class X, class Y, int blockLevel = 2>
951class BlockDiagILU0Preconditioner : public Dune::Preconditioner<X, Y>
952{
953 template<std::size_t i>
954 using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>;
955
956 template<std::size_t i>
957 using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>;
958
959 template<std::size_t i>
960 using BlockILU = Dune::SeqILU<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>, blockLevel-1>;
961
962 using ILUTuple = typename makeFromIndexedType<std::tuple, BlockILU, std::make_index_sequence<M::N()> >::type;
963
964public:
966 using matrix_type = typename std::decay_t<M>;
968 using domain_type = X;
970 using range_type = Y;
972 using field_type = typename X::field_type;
973
980 BlockDiagILU0Preconditioner(const M& m, double w = 1.0)
981 : BlockDiagILU0Preconditioner(m, w, std::make_index_sequence<M::N()>{})
982 {
983 static_assert(blockLevel >= 2, "Only makes sense for MultiTypeBlockMatrix!");
984 }
985
986 void pre (X& v, Y& d) final {}
987
988 void apply (X& v, const Y& d) final
989 {
990 using namespace Dune::Hybrid;
991 forEach(integralRange(Dune::Hybrid::size(ilu_)), [&](const auto i)
992 {
993 std::get<decltype(i)::value>(ilu_).apply(v[i], d[i]);
994 });
995 }
996
997 void post (X&) final {}
998
1000 Dune::SolverCategory::Category category() const final
1001 {
1002 return Dune::SolverCategory::sequential;
1003 }
1004
1005private:
1006 template<std::size_t... Is>
1007 BlockDiagILU0Preconditioner (const M& m, double w, std::index_sequence<Is...> is)
1008 : ilu_(std::make_tuple(BlockILU<Is>(m[Dune::index_constant<Is>{}][Dune::index_constant<Is>{}], w)...))
1009 {}
1010
1011 ILUTuple ilu_;
1012};
1013
1014
1023{
1024
1025public:
1027
1028 template<class Matrix, class Vector>
1029 bool solve(const Matrix& M, Vector& x, const Vector& b)
1030 {
1033 Dune::BiCGSTABSolver<Vector> solver(op, preconditioner, this->residReduction(),
1034 this->maxIter(), this->verbosity());
1035 auto bTmp(b);
1036 solver.apply(x, bTmp, result_);
1037
1038 return result_.converged;
1039 }
1040
1041 const Dune::InverseOperatorResult& result() const
1042 {
1043 return result_;
1044 }
1045
1046 std::string name() const
1047 { return "block-diagonal ILU0-preconditioned BiCGSTAB solver"; }
1048
1049private:
1050 Dune::InverseOperatorResult result_;
1051};
1052
1061{
1062
1063public:
1065
1066 template<int precondBlockLevel = 2, class Matrix, class Vector>
1067 bool solve(const Matrix& M, Vector& x, const Vector& b)
1068 {
1071 static const int restartGMRes = getParamFromGroup<int>(this->paramGroup(), "LinearSolver.GMResRestart");
1072 Dune::RestartedGMResSolver<Vector> solver(op, preconditioner, this->residReduction(), restartGMRes,
1073 this->maxIter(), this->verbosity());
1074 auto bTmp(b);
1075 solver.apply(x, bTmp, result_);
1076
1077 return result_.converged;
1078 }
1079
1080 const Dune::InverseOperatorResult& result() const
1081 {
1082 return result_;
1083 }
1084
1085 std::string name() const
1086 { return "block-diagonal ILU0-preconditioned restarted GMRes solver"; }
1087
1088private:
1089 Dune::InverseOperatorResult result_;
1090};
1091
1096template<class M, class X, class Y, int blockLevel = 2>
1097class BlockDiagAMGPreconditioner : public Dune::Preconditioner<X, Y>
1098{
1099 template<std::size_t i>
1100 using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>;
1101
1102 template<std::size_t i>
1103 using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>;
1104
1105 template<std::size_t i>
1106 using Smoother = Dune::SeqSSOR<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>>;
1107
1108 template<std::size_t i>
1109 using LinearOperator = Dune::MatrixAdapter<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>>;
1110
1111 template<std::size_t i>
1112 using ScalarProduct = Dune::SeqScalarProduct<VecBlockType<i>>;
1113
1114 template<std::size_t i>
1115 using BlockAMG = Dune::Amg::AMG<LinearOperator<i>, VecBlockType<i>, Smoother<i>>;
1116
1117 using AMGTuple = typename makeFromIndexedType<std::tuple, BlockAMG, std::make_index_sequence<M::N()> >::type;
1118
1119public:
1121 using matrix_type = typename std::decay_t<M>;
1123 using domain_type = X;
1125 using range_type = Y;
1127 using field_type = typename X::field_type;
1128
1136 template<class LOP, class Criterion, class SmootherArgs>
1137 BlockDiagAMGPreconditioner(const LOP& lop, const Criterion& c, const SmootherArgs& sa)
1138 : BlockDiagAMGPreconditioner(lop, c, sa, std::make_index_sequence<M::N()>{})
1139 {
1140 static_assert(blockLevel >= 2, "Only makes sense for MultiTypeBlockMatrix!");
1141 }
1142
1159 void pre (X& v, Y& d) final
1160 {
1161 using namespace Dune::Hybrid;
1162 forEach(integralRange(Dune::Hybrid::size(amg_)), [&](const auto i)
1163 {
1164 std::get<decltype(i)::value>(amg_).pre(v[i], d[i]);
1165 });
1166 }
1167
1179 void apply (X& v, const Y& d) final
1180 {
1181 using namespace Dune::Hybrid;
1182 forEach(integralRange(Dune::Hybrid::size(amg_)), [&](const auto i)
1183 {
1184 std::get<decltype(i)::value>(amg_).apply(v[i], d[i]);
1185 });
1186 }
1187
1197 void post (X& v) final
1198 {
1199 using namespace Dune::Hybrid;
1200 forEach(integralRange(Dune::Hybrid::size(amg_)), [&](const auto i)
1201 {
1202 std::get<decltype(i)::value>(amg_).post(v[i]);
1203 });
1204 }
1205
1207 Dune::SolverCategory::Category category() const final
1208 {
1209 return Dune::SolverCategory::sequential;
1210 }
1211
1212private:
1213 template<class LOP, class Criterion, class SmootherArgs, std::size_t... Is>
1214 BlockDiagAMGPreconditioner (const LOP& lop, const Criterion& c, const SmootherArgs& sa, std::index_sequence<Is...> is)
1215 : amg_(std::make_tuple(BlockAMG<Is>(*std::get<Is>(lop), *std::get<Is>(c), *std::get<Is>(sa))...))
1216 {}
1217
1218 AMGTuple amg_;
1219};
1220
1229{
1230 template<class M, std::size_t i>
1231 using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>;
1232
1233 template<class X, std::size_t i>
1234 using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>;
1235
1236 template<class M, class X, std::size_t i>
1237 using Smoother = Dune::SeqSSOR<DiagBlockType<M, i>, VecBlockType<X, i>, VecBlockType<X, i>>;
1238
1239 template<class M, class X, std::size_t i>
1240 using SmootherArgs = typename Dune::Amg::SmootherTraits<Smoother<M, X, i>>::Arguments;
1241
1242 template<class M, std::size_t i>
1243 using Criterion = Dune::Amg::CoarsenCriterion<Dune::Amg::SymmetricCriterion<DiagBlockType<M, i>, Dune::Amg::FirstDiagonal>>;
1244
1245 template<class M, class X, std::size_t i>
1246 using LinearOperator = Dune::MatrixAdapter<DiagBlockType<M, i>, VecBlockType<X, i>, VecBlockType<X, i>>;
1247
1248public:
1250
1251 // Solve saddle-point problem using a Schur complement based preconditioner
1252 template<class Matrix, class Vector>
1253 bool solve(const Matrix& m, Vector& x, const Vector& b)
1254 {
1257 Dune::Amg::Parameters params(15, 2000, 1.2, 1.6, Dune::Amg::atOnceAccu);
1258 params.setDefaultValuesIsotropic(3);
1259 params.setDebugLevel(this->verbosity());
1260
1261 auto criterion = makeCriterion_<Criterion, Matrix>(params, std::make_index_sequence<Matrix::N()>{});
1262 auto smootherArgs = makeSmootherArgs_<SmootherArgs, Matrix, Vector>(std::make_index_sequence<Matrix::N()>{});
1263
1264 using namespace Dune::Hybrid;
1265 forEach(std::make_index_sequence<Matrix::N()>{}, [&](const auto i)
1266 {
1267 auto& args = std::get<decltype(i)::value>(smootherArgs);
1268 args->iterations = 1;
1269 args->relaxationFactor = 1;
1270 });
1271
1272 auto linearOperator = makeLinearOperator_<LinearOperator, Matrix, Vector>(m, std::make_index_sequence<Matrix::N()>{});
1273
1274 BlockDiagAMGPreconditioner<Matrix, Vector, Vector> preconditioner(linearOperator, criterion, smootherArgs);
1275
1277 Dune::BiCGSTABSolver<Vector> solver(op, preconditioner, this->residReduction(),
1278 this->maxIter(), this->verbosity());
1279 auto bTmp(b);
1280 solver.apply(x, bTmp, result_);
1281
1282 return result_.converged;
1283 }
1284
1285 const Dune::InverseOperatorResult& result() const
1286 {
1287 return result_;
1288 }
1289
1290 std::string name() const
1291 { return "block-diagonal AMG-preconditioned BiCGSTAB solver"; }
1292
1293private:
1294 template<template<class M, std::size_t i> class Criterion, class Matrix, class Params, std::size_t... Is>
1295 auto makeCriterion_(const Params& p, std::index_sequence<Is...>)
1296 {
1297 return std::make_tuple(std::make_shared<Criterion<Matrix, Is>>(p)...);
1298 }
1299
1300 template<template<class M, class X, std::size_t i> class SmootherArgs, class Matrix, class Vector, std::size_t... Is>
1301 auto makeSmootherArgs_(std::index_sequence<Is...>)
1302 {
1303 return std::make_tuple(std::make_shared<SmootherArgs<Matrix, Vector, Is>>()...);
1304 }
1305
1306 template<template<class M, class X, std::size_t i> class LinearOperator, class Matrix, class Vector, std::size_t... Is>
1307 auto makeLinearOperator_(const Matrix& m, std::index_sequence<Is...>)
1308 {
1309 return std::make_tuple(std::make_shared<LinearOperator<Matrix, Vector, Is>>(m[Dune::index_constant<Is>{}][Dune::index_constant<Is>{}])...);
1310 }
1311
1312 Dune::InverseOperatorResult result_;
1313};
1314
1315// \}
1316
1317} // end namespace Dumux
1318
1319#endif
A simple ilu0 block diagonal preconditioned BiCGSTABSolver.
Definition: seqsolverbackend.hh:1229
const Dune::InverseOperatorResult & result() const
Definition: seqsolverbackend.hh:1285
std::string name() const
Definition: seqsolverbackend.hh:1290
bool solve(const Matrix &m, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:1253
A simple ilu0 block diagonal preconditioner.
Definition: seqsolverbackend.hh:1098
void post(X &v) final
Clean up.
Definition: seqsolverbackend.hh:1197
typename std::decay_t< M > matrix_type
The matrix type the preconditioner is for.
Definition: seqsolverbackend.hh:1121
void apply(X &v, const Y &d) final
Apply one step of the preconditioner to the system A(v)=d.
Definition: seqsolverbackend.hh:1179
Dune::SolverCategory::Category category() const final
Category of the preconditioner (see SolverCategory::Category)
Definition: seqsolverbackend.hh:1207
Y range_type
The range type of the preconditioner.
Definition: seqsolverbackend.hh:1125
void pre(X &v, Y &d) final
Prepare the preconditioner.
Definition: seqsolverbackend.hh:1159
typename X::field_type field_type
The field type of the preconditioner.
Definition: seqsolverbackend.hh:1127
BlockDiagAMGPreconditioner(const LOP &lop, const Criterion &c, const SmootherArgs &sa)
Constructor.
Definition: seqsolverbackend.hh:1137
X domain_type
The domain type of the preconditioner.
Definition: seqsolverbackend.hh:1123
A simple ilu0 block diagonal preconditioned BiCGSTABSolver.
Definition: seqsolverbackend.hh:1023
std::string name() const
Definition: seqsolverbackend.hh:1046
const Dune::InverseOperatorResult & result() const
Definition: seqsolverbackend.hh:1041
bool solve(const Matrix &M, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:1029
A simple ilu0 block diagonal preconditioner.
Definition: seqsolverbackend.hh:952
typename X::field_type field_type
The field type of the preconditioner.
Definition: seqsolverbackend.hh:972
typename std::decay_t< M > matrix_type
The matrix type the preconditioner is for.
Definition: seqsolverbackend.hh:966
Y range_type
The range type of the preconditioner.
Definition: seqsolverbackend.hh:970
void pre(X &v, Y &d) final
Definition: seqsolverbackend.hh:986
X domain_type
The domain type of the preconditioner.
Definition: seqsolverbackend.hh:968
void post(X &) final
Definition: seqsolverbackend.hh:997
Dune::SolverCategory::Category category() const final
Category of the preconditioner (see SolverCategory::Category)
Definition: seqsolverbackend.hh:1000
void apply(X &v, const Y &d) final
Definition: seqsolverbackend.hh:988
BlockDiagILU0Preconditioner(const M &m, double w=1.0)
Constructor.
Definition: seqsolverbackend.hh:980
A simple ilu0 block diagonal preconditioned RestartedGMResSolver.
Definition: seqsolverbackend.hh:1061
std::string name() const
Definition: seqsolverbackend.hh:1085
const Dune::InverseOperatorResult & result() const
Definition: seqsolverbackend.hh:1080
bool solve(const Matrix &M, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:1067
Solver for simple block-diagonal matrices (e.g. from explicit time stepping schemes)
Definition: seqsolverbackend.hh:742
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:747
std::string name() const
Definition: seqsolverbackend.hh:758
Sequential GS-preconditioned BiCGSTAB solver.
Definition: seqsolverbackend.hh:307
std::string name() const
Definition: seqsolverbackend.hh:321
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:312
Sequential GS-preconditioned CG solver.
Definition: seqsolverbackend.hh:492
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:497
std::string name() const
Definition: seqsolverbackend.hh:506
Sequential ILU(0)-preconditioned BiCGSTAB solver.
Definition: seqsolverbackend.hh:603
std::string name() const
Definition: seqsolverbackend.hh:617
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:608
Sequential ILU(0)-preconditioned CG solver.
Definition: seqsolverbackend.hh:639
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:644
std::string name() const
Definition: seqsolverbackend.hh:653
Sequential ILU0-preconditioned GMRes solver.
Definition: seqsolverbackend.hh:676
std::string name() const
Definition: seqsolverbackend.hh:690
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:681
Sequential ILU(n)-preconditioned BiCSTAB solver.
Definition: seqsolverbackend.hh:193
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:198
std::string name() const
Definition: seqsolverbackend.hh:207
Sequential ILU(n)-preconditioned CG solver.
Definition: seqsolverbackend.hh:381
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:386
std::string name() const
Definition: seqsolverbackend.hh:395
Sequential ILU(n)-preconditioned GMRes solver.
Definition: seqsolverbackend.hh:714
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:719
std::string name() const
Definition: seqsolverbackend.hh:728
A general solver backend allowing arbitrary preconditioners and solvers.
Definition: seqsolverbackend.hh:54
static bool solveWithILU0PrecGMRes(const SolverInterface &s, const Matrix &A, Vector &x, const Vector &b, const std::string &modelParamGroup="")
Definition: seqsolverbackend.hh:121
static bool solveWithGMRes(const SolverInterface &s, const Matrix &A, Vector &x, const Vector &b, const std::string &modelParamGroup="")
Definition: seqsolverbackend.hh:78
static bool solveWithILU0Prec(const SolverInterface &s, const Matrix &A, Vector &x, const Vector &b, const std::string &modelParamGroup="")
Definition: seqsolverbackend.hh:101
static bool solve(const SolverInterface &s, const Matrix &A, Vector &x, const Vector &b, const std::string &modelParamGroup="")
Definition: seqsolverbackend.hh:58
static bool solveWithParamTree(const Matrix &A, Vector &x, const Vector &b, const Dune::ParameterTree &params)
Definition: seqsolverbackend.hh:144
Sequential Jacobi-preconditioned BiCSTAB solver.
Definition: seqsolverbackend.hh:344
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:349
std::string name() const
Definition: seqsolverbackend.hh:358
Sequential Jacobi-preconditioned CG solver.
Definition: seqsolverbackend.hh:528
std::string name() const
Definition: seqsolverbackend.hh:542
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:533
Base class for linear solvers.
Definition: solver.hh:27
Scalar residReduction() const
the linear solver residual reduction
Definition: solver.hh:98
int maxIter() const
the maximum number of linear solver iterations
Definition: solver.hh:90
const std::string & paramGroup() const
the parameter group for getting parameter from the parameter tree
Definition: solver.hh:78
LinearSolver(const std::string &paramGroup="")
Construct the solver.
Definition: solver.hh:43
int verbosity() const
the verbosity level
Definition: solver.hh:82
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:47
Adapter to turn a multi-type matrix into a thread-parallel linear operator. Adapts a matrix to the as...
Definition: parallelmatrixadapter.hh:28
Sequential SOR-preconditioned BiCSTAB solver.
Definition: seqsolverbackend.hh:231
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:236
std::string name() const
Definition: seqsolverbackend.hh:245
Sequential SOR-preconditioned CG solver.
Definition: seqsolverbackend.hh:418
std::string name() const
Definition: seqsolverbackend.hh:432
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:423
Sequential SSOR-preconditioned BiCGSTAB solver.
Definition: seqsolverbackend.hh:269
std::string name() const
Definition: seqsolverbackend.hh:283
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:274
Sequential SSOR-preconditioned CG solver.
Definition: seqsolverbackend.hh:455
std::string name() const
Definition: seqsolverbackend.hh:469
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:460
Sequential SSOR-preconditioned GMRes solver.
Definition: seqsolverbackend.hh:566
std::string name() const
Definition: seqsolverbackend.hh:580
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:571
A preconditioner based on the Uzawa algorithm for saddle-point problems of the form .
Definition: preconditioners.hh:70
A Uzawa preconditioned BiCGSTAB solver for saddle-point problems.
Definition: seqsolverbackend.hh:927
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:932
std::string name() const
Definition: seqsolverbackend.hh:940
constexpr std::size_t preconditionerBlockLevel() noexcept
Returns the block level for the preconditioner for a given matrix.
Definition: seqsolverbackend.hh:170
Generates a parameter tree required for the linear solvers and precondioners of the Dune ISTL.
Type traits to be used with matrix types.
constexpr std::size_t blockSize()
Definition: dunevectors.hh:27
Definition: adapt.hh:17
Definition: common/pdesolver.hh:24
A parallel version of a linear operator.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Dumux preconditioners for iterative solvers.
Base class for linear solvers.
Helper type to determine whether a given type is a Dune::MultiTypeBlockMatrix.
Definition: matrix.hh:37
Definition: utility.hh:28
Utilities for template meta programming.