version 3.8
istlsolvers.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_LINEAR_ISTL_SOLVERS_HH
13#define DUMUX_LINEAR_ISTL_SOLVERS_HH
14
15#include <memory>
16#include <variant>
17
18#include <dune/common/exceptions.hh>
19#include <dune/common/shared_ptr.hh>
20#include <dune/common/parallel/indexset.hh>
21#include <dune/common/parallel/mpicommunication.hh>
22#include <dune/grid/common/capabilities.hh>
23#include <dune/istl/solvers.hh>
24#include <dune/istl/solverfactory.hh>
25#include <dune/istl/owneroverlapcopy.hh>
26#include <dune/istl/scalarproducts.hh>
27#include <dune/istl/paamg/amg.hh>
28#include <dune/istl/paamg/pinfo.hh>
29
39
40#include <dune/istl/foreach.hh>
41
43
50template<class M>
51constexpr std::size_t preconditionerBlockLevel() noexcept
52{
54}
55
56template<template<class,class,class,int> class Preconditioner, int blockLevel = 1>
58{
59public:
60 template<class TL, class M>
61 auto operator() (TL typeList, const M& matrix, const Dune::ParameterTree& config)
62 {
63 using Matrix = typename Dune::TypeListElement<0, decltype(typeList)>::type;
64 using Domain = typename Dune::TypeListElement<1, decltype(typeList)>::type;
65 using Range = typename Dune::TypeListElement<2, decltype(typeList)>::type;
66 std::shared_ptr<Dune::Preconditioner<Domain, Range>> preconditioner
67 = std::make_shared<Preconditioner<Matrix, Domain, Range, blockLevel>>(matrix, config);
68 return preconditioner;
69 }
70};
71
72template<template<class,class,class> class Preconditioner>
74{
75 template<class TL, class M>
76 auto operator() (TL typeList, const M& matrix, const Dune::ParameterTree& config)
77 {
78 using Matrix = typename Dune::TypeListElement<0, decltype(typeList)>::type;
79 using Domain = typename Dune::TypeListElement<1, decltype(typeList)>::type;
80 using Range = typename Dune::TypeListElement<2, decltype(typeList)>::type;
81 std::shared_ptr<Dune::Preconditioner<Domain, Range>> preconditioner
82 = std::make_shared<Preconditioner<Matrix, Domain, Range>>(matrix, config);
83 return preconditioner;
84 }
85};
86
87using IstlAmgPreconditionerFactory = Dune::AMGCreator;
88
89template<class M, bool convert = false>
90struct MatrixForSolver { using type = M; };
91
92template<class M>
93struct MatrixForSolver<M, true>
94{ using type = std::decay_t<decltype(MatrixConverter<M>::multiTypeToBCRSMatrix(std::declval<M>()))>; };
95
96template<class V, bool convert = false>
97struct VectorForSolver { using type = V; };
98
99template<class V>
100struct VectorForSolver<V, true>
101{ using type = std::decay_t<decltype(VectorConverter<V>::multiTypeToBlockVector(std::declval<V>()))>; };
102
103template<class LSTraits, class LATraits, bool convert, bool parallel = LSTraits::canCommunicate>
105
106template<class LSTraits, class LATraits, bool convert>
107struct MatrixOperator<LSTraits, LATraits, convert, true>
108{
111#if HAVE_MPI
112 using type = std::variant<
113 std::shared_ptr<typename LSTraits::template Sequential<M, V>::LinearOperator>,
114 std::shared_ptr<typename LSTraits::template ParallelOverlapping<M, V>::LinearOperator>,
115 std::shared_ptr<typename LSTraits::template ParallelNonoverlapping<M, V>::LinearOperator>
116 >;
117#else
118 using type = std::variant<
119 std::shared_ptr<typename LSTraits::template Sequential<M, V>::LinearOperator>
120 >;
121#endif
122};
123
124template<class LSTraits, class LATraits, bool convert>
125struct MatrixOperator<LSTraits, LATraits, convert, false>
126{
129 using type = std::variant<
130 std::shared_ptr<typename LSTraits::template Sequential<M, V>::LinearOperator>
131 >;
132};
133
134} // end namespace Dumux::Detail::IstlSolvers
135
136namespace Dumux::Detail {
137
138struct IstlSolverResult : public Dune::InverseOperatorResult
139{
140 IstlSolverResult() = default;
143
144 IstlSolverResult(const Dune::InverseOperatorResult& o) : InverseOperatorResult(o) {}
145 IstlSolverResult(Dune::InverseOperatorResult&& o) : InverseOperatorResult(std::move(o)) {}
146
147 operator bool() const { return this->converged; }
148};
149
154template<class LinearSolverTraits, class LinearAlgebraTraits,
155 class InverseOperator, class PreconditionerFactory,
156 bool convertMultiTypeLATypes = false>
158{
159 using Matrix = typename LinearAlgebraTraits::Matrix;
160 using XVector = typename LinearAlgebraTraits::Vector;
161 using BVector = typename LinearAlgebraTraits::Vector;
162 using Scalar = typename InverseOperator::real_type;
163
164 using ScalarProduct = Dune::ScalarProduct<typename InverseOperator::domain_type>;
165
166 static constexpr bool convertMultiTypeVectorAndMatrix
167 = convertMultiTypeLATypes && isMultiTypeBlockVector<XVector>::value;
171 // a variant type that can hold sequential, overlapping, and non-overlapping operators
172 using MatrixOperatorHolder = typename Detail::IstlSolvers::MatrixOperator<
173 LinearSolverTraits, LinearAlgebraTraits, convertMultiTypeVectorAndMatrix
174 >::type;
175
176#if HAVE_MPI
177 using Comm = Dune::OwnerOverlapCopyCommunication<Dune::bigunsignedint<96>, int>;
179#endif
180
181 using ParameterInitializer = std::variant<std::string, Dune::ParameterTree>;
182public:
183
187 IstlIterativeLinearSolver(const ParameterInitializer& params = "")
188 {
189 if (Dune::MPIHelper::getCommunication().size() > 1)
190 DUNE_THROW(Dune::InvalidStateException, "Using sequential constructor for parallel run. Use signature with gridView and dofMapper!");
191
192 initializeParameters_(params);
193 solverCategory_ = Dune::SolverCategory::sequential;
194 scalarProduct_ = std::make_shared<ScalarProduct>();
195 }
196
200 template <class GridView, class DofMapper>
201 IstlIterativeLinearSolver(const GridView& gridView,
202 const DofMapper& dofMapper,
203 const ParameterInitializer& params = "")
204 {
205 initializeParameters_(params);
206#if HAVE_MPI
207 solverCategory_ = Detail::solverCategory<LinearSolverTraits>(gridView);
209 {
210
211 if (solverCategory_ != Dune::SolverCategory::sequential)
212 {
213 parallelHelper_ = std::make_shared<ParallelISTLHelper<LinearSolverTraits>>(gridView, dofMapper);
214 communication_ = std::make_shared<Comm>(gridView.comm(), solverCategory_);
215 scalarProduct_ = Dune::createScalarProduct<XVector>(*communication_, solverCategory_);
216 parallelHelper_->createParallelIndexSet(*communication_);
217 }
218 else
219 scalarProduct_ = std::make_shared<ScalarProduct>();
220 }
221 else
222 scalarProduct_ = std::make_shared<ScalarProduct>();
223#else
224 solverCategory_ = Dune::SolverCategory::sequential;
225 scalarProduct_ = std::make_shared<ScalarProduct>();
226#endif
227 }
228
229#if HAVE_MPI
233 template <class GridView, class DofMapper>
234 IstlIterativeLinearSolver(std::shared_ptr<Comm> communication,
235 std::shared_ptr<ScalarProduct> scalarProduct,
236 const GridView& gridView,
237 const DofMapper& dofMapper,
238 const ParameterInitializer& params = "")
239 {
240 initializeParameters_(params);
241 solverCategory_ = Detail::solverCategory(gridView);
242 scalarProduct_ = scalarProduct;
243 communication_ = communication;
245 {
246 if (solverCategory_ != Dune::SolverCategory::sequential)
247 {
248 parallelHelper_ = std::make_shared<ParallelISTLHelper<LinearSolverTraits>>(gridView, dofMapper);
249 parallelHelper_->createParallelIndexSet(communication);
250 }
251 }
252 }
253#endif
254
258 IstlSolverResult solve(Matrix& A, XVector& x, BVector& b)
259 { return solveSequentialOrParallel_(A, x, b); }
260
264 void setMatrix(std::shared_ptr<Matrix> A)
265 {
266 linearOperator_ = makeParallelOrSequentialLinearOperator_(std::move(A));
267 solver_ = constructPreconditionedSolver_(linearOperator_);
268 }
269
274 void setMatrix(Matrix& A)
275 { setMatrix(Dune::stackobject_to_shared_ptr(A)); }
276
280 IstlSolverResult solve(XVector& x, BVector& b) const
281 {
282 if (!solver_)
283 DUNE_THROW(Dune::InvalidStateException, "Called solve(x, b) but no linear operator has been set");
284
285 return solveSequentialOrParallel_(x, b, *solver_);
286 }
287
291 Scalar norm(const XVector& x) const
292 {
293#if HAVE_MPI
295 {
296 if (solverCategory_ == Dune::SolverCategory::nonoverlapping)
297 {
298 auto y(x); // make a copy because the vector needs to be made consistent
299 using GV = typename LinearSolverTraits::GridView;
300 using DM = typename LinearSolverTraits::DofMapper;
301 ParallelVectorHelper<GV, DM, LinearSolverTraits::dofCodim> vectorHelper(parallelHelper_->gridView(), parallelHelper_->dofMapper());
302 vectorHelper.makeNonOverlappingConsistent(y);
303 return scalarProduct_->norm(y);
304 }
305 }
306#endif
307 if constexpr (convertMultiTypeVectorAndMatrix)
308 {
310 return scalarProduct_->norm(y);
311 }
312 else
313 return scalarProduct_->norm(x);
314 }
315
319 const std::string& name() const
320 {
321 return name_;
322 }
323
327 void setResidualReduction(double residReduction)
328 {
329 params_["reduction"] = std::to_string(residReduction);
330
331 // reconstruct the solver with new parameters
332 if (solver_)
333 solver_ = constructPreconditionedSolver_(linearOperator_);
334 }
335
339 void setMaxIter(std::size_t maxIter)
340 {
341 params_["maxit"] = std::to_string(maxIter);
342
343 // reconstruct the solver with new parameters
344 if (solver_)
345 solver_ = constructPreconditionedSolver_(linearOperator_);
346 }
347
353 void setParams(const ParameterInitializer& params)
354 {
355 initializeParameters_(params);
356
357 // reconstruct the solver with new parameters
358 if (solver_)
359 solver_ = constructPreconditionedSolver_(linearOperator_);
360 }
361
362private:
363
364 void initializeParameters_(const ParameterInitializer& params)
365 {
366 if (std::holds_alternative<std::string>(params))
367 params_ = Dumux::LinearSolverParameters<LinearSolverTraits>::createParameterTree(std::get<std::string>(params));
368 else
369 params_ = std::get<Dune::ParameterTree>(params);
370 }
371
372 MatrixOperatorHolder makeSequentialLinearOperator_(std::shared_ptr<Matrix> A)
373 {
374 using SequentialTraits = typename LinearSolverTraits::template Sequential<MatrixForSolver, XVectorForSolver>;
375 if constexpr (convertMultiTypeVectorAndMatrix)
376 {
377 // create the BCRS matrix the IterativeSolver backend can handle
378 auto M = std::make_shared<MatrixForSolver>(MatrixConverter<Matrix>::multiTypeToBCRSMatrix(*A));
379 return std::make_shared<typename SequentialTraits::LinearOperator>(M);
380 }
381 else
382 {
383 return std::make_shared<typename SequentialTraits::LinearOperator>(A);
384 }
385 }
386
387 template<class ParallelTraits>
388 MatrixOperatorHolder makeParallelLinearOperator_(std::shared_ptr<Matrix> A, ParallelTraits = {})
389 {
390#if HAVE_MPI
391 // make matrix consistent
392 prepareMatrixParallel<LinearSolverTraits, ParallelTraits>(*A, *parallelHelper_);
393 return std::make_shared<typename ParallelTraits::LinearOperator>(std::move(A), *communication_);
394#else
395 DUNE_THROW(Dune::InvalidStateException, "Calling makeParallelLinearOperator for sequential run");
396#endif
397 }
398
399 MatrixOperatorHolder makeParallelOrSequentialLinearOperator_(std::shared_ptr<Matrix> A)
400 {
401 return executeSequentialOrParallel_(
402 [&]{ return makeSequentialLinearOperator_(std::move(A)); },
403 [&](auto traits){ return makeParallelLinearOperator_(std::move(A), traits); }
404 );
405 }
406
407 MatrixOperatorHolder makeSequentialLinearOperator_(Matrix& A)
408 { return makeSequentialLinearOperator_(Dune::stackobject_to_shared_ptr<Matrix>(A)); }
409
410 MatrixOperatorHolder makeParallelOrSequentialLinearOperator_(Matrix& A)
411 { return makeParallelOrSequentialLinearOperator_(Dune::stackobject_to_shared_ptr<Matrix>(A)); }
412
413 template<class ParallelTraits>
414 MatrixOperatorHolder makeParallelLinearOperator_(Matrix& A, ParallelTraits = {})
415 { return makeParallelLinearOperator_<ParallelTraits>(Dune::stackobject_to_shared_ptr<Matrix>(A)); }
416
417 IstlSolverResult solveSequential_(Matrix& A, XVector& x, BVector& b)
418 {
419 // construct solver from linear operator
420 auto linearOperatorHolder = makeSequentialLinearOperator_(A);
421 auto solver = constructPreconditionedSolver_(linearOperatorHolder);
422
423 return solveSequential_(x, b, *solver);
424 }
425
426 IstlSolverResult solveSequential_(XVector& x, BVector& b, InverseOperator& solver) const
427 {
428 Dune::InverseOperatorResult result;
429 if constexpr (convertMultiTypeVectorAndMatrix)
430 {
431 // create the vector the IterativeSolver backend can handle
432 BVectorForSolver bTmp = VectorConverter<BVector>::multiTypeToBlockVector(b);
433
434 // create a block vector to which the linear solver writes the solution
435 XVectorForSolver y(bTmp.size());
436
437 // solve linear system
438 solver.apply(y, bTmp, result);
439
440 // copy back the result y into x
441 if (result.converged)
443 }
444 else
445 {
446 // solve linear system
447 solver.apply(x, b, result);
448 }
449
450 return result;
451 }
452
453 IstlSolverResult solveSequentialOrParallel_(Matrix& A, XVector& x, BVector& b)
454 {
455 return executeSequentialOrParallel_(
456 [&]{ return solveSequential_(A, x, b); },
457 [&](auto traits){ return solveParallel_(A, x, b, traits); }
458 );
459 }
460
461 IstlSolverResult solveSequentialOrParallel_(XVector& x, BVector& b, InverseOperator& solver) const
462 {
463 return executeSequentialOrParallel_(
464 [&]{ return solveSequential_(x, b, solver); },
465 [&](auto traits){ return solveParallel_(x, b, solver, traits); }
466 );
467 }
468
469 template<class ParallelTraits>
470 IstlSolverResult solveParallel_(Matrix& A, XVector& x, BVector& b, ParallelTraits = {})
471 {
472 // construct solver from linear operator
473 auto linearOperatorHolder = makeParallelLinearOperator_<ParallelTraits>(A);
474 auto solver = constructPreconditionedSolver_(linearOperatorHolder);
475 return solveParallel_<ParallelTraits>(x, b, *solver);
476 }
477
478 template<class ParallelTraits>
479 IstlSolverResult solveParallel_(XVector& x, BVector& b, InverseOperator& solver, ParallelTraits = {}) const
480 {
481#if HAVE_MPI
482 // make right hand side consistent
483 prepareVectorParallel<LinearSolverTraits, ParallelTraits>(b, *parallelHelper_);
484
485 // solve linear system
486 Dune::InverseOperatorResult result;
487 solver.apply(x, b, result);
488 return result;
489#else
490 DUNE_THROW(Dune::InvalidStateException, "Calling makeParallelLinearOperator for sequential run");
491#endif
492 }
493
494
495 std::shared_ptr<InverseOperator> constructPreconditionedSolver_(MatrixOperatorHolder& ops)
496 {
497 return std::visit([&](auto&& op)
498 {
499 using LinearOperator = typename std::decay_t<decltype(op)>::element_type;
500 const auto& params = params_.sub("preconditioner");
501 using Prec = Dune::Preconditioner<typename LinearOperator::domain_type, typename LinearOperator::range_type>;
502 using TL = Dune::TypeList<typename LinearOperator::matrix_type, typename LinearOperator::domain_type, typename LinearOperator::range_type>;
503 std::shared_ptr<Prec> prec = PreconditionerFactory{}(TL{}, op, params);
504
505#if HAVE_MPI
506 if (prec->category() != op->category() && prec->category() == Dune::SolverCategory::sequential)
507 prec = Dune::wrapPreconditioner4Parallel(prec, op);
508#endif
509 return std::make_shared<InverseOperator>(op, scalarProduct_, prec, params_);
510 }, ops);
511 }
512
513 template<class Seq, class Par>
514 decltype(auto) executeSequentialOrParallel_(Seq&& sequentialAction, Par&& parallelAction) const
515 {
516#if HAVE_MPI
517 // For Dune::MultiTypeBlockMatrix there is currently no generic way
518 // of handling parallelism, we therefore can only solve these types of systems sequentially
519 if constexpr (isMultiTypeBlockMatrix<Matrix>::value || !LinearSolverTraits::canCommunicate)
520 return sequentialAction();
521 else
522 {
523 switch (solverCategory_)
524 {
525 case Dune::SolverCategory::sequential:
526 return sequentialAction();
527 case Dune::SolverCategory::nonoverlapping:
528 using NOTraits = typename LinearSolverTraits::template ParallelNonoverlapping<Matrix, XVector>;
529 return parallelAction(NOTraits{});
530 case Dune::SolverCategory::overlapping:
531 using OTraits = typename LinearSolverTraits::template ParallelOverlapping<Matrix, XVector>;
532 return parallelAction(OTraits{});
533 default: DUNE_THROW(Dune::InvalidStateException, "Unknown solver category");
534 }
535 }
536#else
537 return sequentialAction();
538#endif
539 }
540
541#if HAVE_MPI
542 std::shared_ptr<const ParallelHelper> parallelHelper_;
543 std::shared_ptr<Comm> communication_;
544#endif
545
546 Dune::SolverCategory::Category solverCategory_;
547 std::shared_ptr<ScalarProduct> scalarProduct_;
548
549 // for stored solvers (reuse matrix)
550 MatrixOperatorHolder linearOperator_;
551 // for stored solvers (reuse matrix)
552 std::shared_ptr<InverseOperator> solver_;
553
554 Dune::ParameterTree params_;
555 std::string name_;
556};
557
558} // end namespace Dumux::Detail
559
560namespace Dumux {
561
578template<class LSTraits, class LATraits>
580 Detail::IstlIterativeLinearSolver<LSTraits, LATraits,
581 Dune::BiCGSTABSolver<typename LATraits::SingleTypeVector>,
583 // the Dune::ILU preconditioners don't accept multi-type matrices
584 /*convert multi-type istl types?*/ true
585 >;
586
603template<class LSTraits, class LATraits>
605 Detail::IstlIterativeLinearSolver<LSTraits, LATraits,
606 Dune::RestartedGMResSolver<typename LATraits::SingleTypeVector>,
608 // the Dune::ILU preconditioners don't accept multi-type matrices
609 /*convert multi-type istl types?*/ true
610 >;
611
629template<class LSTraits, class LATraits>
631 Detail::IstlIterativeLinearSolver<LSTraits, LATraits,
632 Dune::BiCGSTABSolver<typename LATraits::Vector>,
634 >;
635
652template<class LSTraits, class LATraits>
654 Detail::IstlIterativeLinearSolver<LSTraits, LATraits,
655 Dune::CGSolver<typename LATraits::Vector>,
657 >;
658
672template<class LSTraits, class LATraits>
674 Detail::IstlIterativeLinearSolver<LSTraits, LATraits,
675 Dune::BiCGSTABSolver<typename LATraits::SingleTypeVector>,
677 // the AMG preconditioner doesn't accept multi-type matrices
678 /*convert multi-type istl types?*/ true
679 >;
680
693template<class LSTraits, class LATraits>
695 Detail::IstlIterativeLinearSolver<LSTraits, LATraits,
696 Dune::CGSolver<typename LATraits::SingleTypeVector>,
698 // the AMG preconditioner doesn't accept multi-type matrices
699 /*convert multi-type istl types?*/ true
700 >;
701
716template<class LSTraits, class LATraits>
718 Detail::IstlIterativeLinearSolver<LSTraits, LATraits,
719 Dune::BiCGSTABSolver<typename LATraits::Vector>,
721 >;
722
723} // end namespace Dumux
724
725namespace Dumux::Detail {
726
731template<class LSTraits, class LATraits, template<class M> class Solver,
732 bool convertMultiTypeVectorAndMatrix = isMultiTypeBlockVector<typename LATraits::Vector>::value>
734{
735 using Matrix = typename LATraits::Matrix;
736 using XVector = typename LATraits::Vector;
737 using BVector = typename LATraits::Vector;
738
742 using InverseOperator = Dune::InverseOperator<XVectorForSolver, BVectorForSolver>;
743public:
745
749 IstlSolverResult solve(const Matrix& A, XVector& x, const BVector& b)
750 {
751 return solve_(A, x, b);
752 }
753
757 IstlSolverResult solve(XVector& x, const BVector& b)
758 {
759 if (!solver_)
760 DUNE_THROW(Dune::InvalidStateException, "Called solve(x, b) but no linear operator has been set");
761
762 return solve_(x, b, *solver_);
763 }
764
768 void setMatrix(std::shared_ptr<Matrix> A)
769 {
770 if constexpr (convertMultiTypeVectorAndMatrix)
771 matrix_ = std::make_shared<MatrixForSolver>(MatrixConverter<Matrix>::multiTypeToBCRSMatrix(A));
772 else
773 matrix_ = A;
774
775 solver_ = std::make_shared<Solver<MatrixForSolver>>(*matrix_);
776 }
777
782 void setMatrix(Matrix& A)
783 { setMatrix(Dune::stackobject_to_shared_ptr(A)); }
784
788 std::string name() const
789 {
790 return "Direct solver";
791 }
792
793private:
794 IstlSolverResult solve_(const Matrix& A, XVector& x, const BVector& b)
795 {
796 // support dune-istl multi-type block vector/matrix by copying
797 if constexpr (convertMultiTypeVectorAndMatrix)
798 {
800 Solver<MatrixForSolver> solver(AA, this->verbosity() > 0);
801 return solve_(x, b, solver);
802 }
803 else
804 {
805 Solver<MatrixForSolver> solver(A, this->verbosity() > 0);
806 return solve_(x, b, solver);
807 }
808 }
809
810 IstlSolverResult solve_(XVector& x, const BVector& b, InverseOperator& solver) const
811 {
812 Dune::InverseOperatorResult result;
813
814 if constexpr (convertMultiTypeVectorAndMatrix)
815 {
817 XVectorForSolver xx(bb.size());
818 solver.apply(xx, bb, result);
819 checkResult_(xx, result);
820 if (result.converged)
822 return result;
823 }
824 else
825 {
826 BVectorForSolver bTmp(b);
827 solver.apply(x, bTmp, result);
828 checkResult_(x, result);
829 return result;
830 }
831 }
832
833
834 void checkResult_(XVectorForSolver& x, Dune::InverseOperatorResult& result) const
835 {
836 flatVectorForEach(x, [&](auto&& entry, std::size_t){
837 using std::isnan, std::isinf;
838 if (isnan(entry) || isinf(entry))
839 result.converged = false;
840 });
841 }
842
844 std::shared_ptr<MatrixForSolver> matrix_;
846 std::shared_ptr<InverseOperator> solver_;
847};
848
849} // end namespace Dumux::Detail
850
851#if HAVE_SUPERLU
852#include <dune/istl/superlu.hh>
853
854namespace Dumux {
855
864template<class LSTraits, class LATraits>
865using SuperLUIstlSolver = Detail::DirectIstlSolver<LSTraits, LATraits, Dune::SuperLU>;
866
867} // end namespace Dumux
868
869#endif // HAVE_SUPERLU
870
871#if HAVE_UMFPACK
872#include <dune/istl/umfpack.hh>
873#include <dune/common/version.hh>
874
875namespace Dumux {
876
885template<class LSTraits, class LATraits>
886using UMFPackIstlSolver = Detail::DirectIstlSolver<
887 LSTraits, LATraits, Dune::UMFPack
888#if DUNE_VERSION_GTE(DUNE_ISTL,2,10)
889 , false // no need to convert multi-type matrix anymore
890#endif
891>;
892
893} // end namespace Dumux
894
895#endif // HAVE_UMFPACK
896
897#endif
Direct dune-istl linear solvers.
Definition: istlsolvers.hh:734
void setMatrix(std::shared_ptr< Matrix > A)
Set the matrix A of the linear system Ax = b for reuse.
Definition: istlsolvers.hh:768
IstlSolverResult solve(const Matrix &A, XVector &x, const BVector &b)
Solve the linear system Ax = b.
Definition: istlsolvers.hh:749
std::string name() const
name of the linear solver
Definition: istlsolvers.hh:788
void setMatrix(Matrix &A)
Set the matrix A of the linear system Ax = b for reuse.
Definition: istlsolvers.hh:782
IstlSolverResult solve(XVector &x, const BVector &b)
Solve the linear system Ax = b using the matrix set with setMatrix.
Definition: istlsolvers.hh:757
Standard dune-istl iterative linear solvers.
Definition: istlsolvers.hh:158
IstlIterativeLinearSolver(const ParameterInitializer &params="")
Constructor for sequential solvers.
Definition: istlsolvers.hh:187
void setMatrix(Matrix &A)
Set the matrix A of the linear system Ax = b for reuse.
Definition: istlsolvers.hh:274
void setResidualReduction(double residReduction)
Set the residual reduction tolerance.
Definition: istlsolvers.hh:327
void setMaxIter(std::size_t maxIter)
Set the maximum number of linear solver iterations.
Definition: istlsolvers.hh:339
IstlIterativeLinearSolver(const GridView &gridView, const DofMapper &dofMapper, const ParameterInitializer &params="")
Constructor for parallel and sequential solvers.
Definition: istlsolvers.hh:201
const std::string & name() const
The name of the linear solver.
Definition: istlsolvers.hh:319
void setParams(const ParameterInitializer &params)
Set the linear solver parameters.
Definition: istlsolvers.hh:353
IstlSolverResult solve(Matrix &A, XVector &x, BVector &b)
Solve the linear system Ax = b.
Definition: istlsolvers.hh:258
void setMatrix(std::shared_ptr< Matrix > A)
Set the matrix A of the linear system Ax = b for reuse.
Definition: istlsolvers.hh:264
IstlSolverResult solve(XVector &x, BVector &b) const
Solve the linear system Ax = b where A has been set with setMatrix.
Definition: istlsolvers.hh:280
IstlIterativeLinearSolver(std::shared_ptr< Comm > communication, std::shared_ptr< ScalarProduct > scalarProduct, const GridView &gridView, const DofMapper &dofMapper, const ParameterInitializer &params="")
Constructor with custom scalar product and communication.
Definition: istlsolvers.hh:234
Scalar norm(const XVector &x) const
Compute the 2-norm of vector x.
Definition: istlsolvers.hh:291
auto operator()(TL typeList, const M &matrix, const Dune::ParameterTree &config)
Definition: istlsolvers.hh:61
Definition: parallelhelpers.hh:29
Base class for linear solvers.
Definition: solver.hh:27
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
A helper class that converts a Dune::MultiTypeBlockMatrix into a plain Dune::BCRSMatrix TODO: allow b...
Definition: matrixconverter.hh:35
static auto multiTypeToBCRSMatrix(const MultiTypeBlockMatrix &A)
Converts the matrix to a type the IterativeSolverBackend can handle.
Definition: matrixconverter.hh:46
Definition: parallelhelpers.hh:476
void makeNonOverlappingConsistent(Dune::BlockVector< Block, Alloc > &v) const
Make a vector consistent for non-overlapping domain decomposition methods.
Definition: parallelhelpers.hh:484
static void retrieveValues(MultiTypeBlockVector &x, const BlockVector &y)
Copies the entries of a Dune::BlockVector to a Dune::MultiTypeBlockVector.
Definition: matrixconverter.hh:229
static auto multiTypeToBlockVector(const MultiTypeBlockVector &b)
Converts a Dune::MultiTypeBlockVector to a plain 1x1 Dune::BlockVector.
Definition: matrixconverter.hh:203
constexpr std::size_t preconditionerBlockLevel() noexcept
Returns the block level for the preconditioner for a given matrix.
Definition: istlsolvers.hh:51
Define traits for linear algebra backends.
Generates a parameter tree required for the linear solvers and precondioners of the Dune ISTL.
Type traits to be used with matrix types.
A helper class that converts a Dune::MultiTypeBlockMatrix into a plain Dune::BCRSMatrix.
Definition: istlsolvers.hh:42
Dune::AMGCreator IstlAmgPreconditionerFactory
Definition: istlsolvers.hh:87
Distance implementation details.
Definition: cvfelocalresidual.hh:25
static constexpr bool canCommunicate
Definition: gridcapabilities.hh:51
Dune::SolverCategory::Category solverCategory(const GridView &gridView)
Definition: solvercategory.hh:20
Definition: adapt.hh:17
LinearSolverTraitsImpl< GridGeometry, typename GridGeometry::DiscretizationMethod > LinearSolverTraits
The type traits required for using the IstlFactoryBackend.
Definition: linearsolvertraits.hh:37
Provides a helper class for nonoverlapping decomposition.
Dumux preconditioners for iterative solvers.
Base class for linear solvers.
Solver category.
Definition: istlsolvers.hh:139
IstlSolverResult(IstlSolverResult &&)=default
IstlSolverResult(const Dune::InverseOperatorResult &o)
Definition: istlsolvers.hh:144
IstlSolverResult(Dune::InverseOperatorResult &&o)
Definition: istlsolvers.hh:145
IstlSolverResult(const IstlSolverResult &)=default
std::decay_t< decltype(MatrixConverter< M >::multiTypeToBCRSMatrix(std::declval< M >()))> type
Definition: istlsolvers.hh:94
Definition: istlsolvers.hh:90
M type
Definition: istlsolvers.hh:90
typename VectorForSolver< typename LATraits::Vector, convert >::type V
Definition: istlsolvers.hh:128
std::variant< std::shared_ptr< typename LSTraits::template Sequential< M, V >::LinearOperator > > type
Definition: istlsolvers.hh:131
typename MatrixForSolver< typename LATraits::Matrix, convert >::type M
Definition: istlsolvers.hh:127
std::variant< std::shared_ptr< typename LSTraits::template Sequential< M, V >::LinearOperator >, std::shared_ptr< typename LSTraits::template ParallelOverlapping< M, V >::LinearOperator >, std::shared_ptr< typename LSTraits::template ParallelNonoverlapping< M, V >::LinearOperator > > type
Definition: istlsolvers.hh:116
typename VectorForSolver< typename LATraits::Vector, convert >::type V
Definition: istlsolvers.hh:110
typename MatrixForSolver< typename LATraits::Matrix, convert >::type M
Definition: istlsolvers.hh:109
Definition: istlsolvers.hh:104
std::decay_t< decltype(VectorConverter< V >::multiTypeToBlockVector(std::declval< V >()))> type
Definition: istlsolvers.hh:101
Definition: istlsolvers.hh:97
V type
Definition: istlsolvers.hh:97
Definition: linearalgebratraits.hh:51
V Vector
Definition: linearalgebratraits.hh:53
M Matrix
Definition: linearalgebratraits.hh:52
Helper type to determine whether a given type is a Dune::MultiTypeBlockMatrix.
Definition: matrix.hh:37
Helper type to determine whether a given type is a Dune::MultiTypeBlockVector.
Definition: vector.hh:22
Type traits to be used with vector types.