version 3.11-dev
Loading...
Searching...
No Matches
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-FileCopyrightText: 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/version.hh>
21#include <dune/common/parallel/indexset.hh>
22#include <dune/common/parallel/mpicommunication.hh>
23#include <dune/grid/common/capabilities.hh>
24#include <dune/istl/solvers.hh>
25#include <dune/istl/solverfactory.hh>
26#include <dune/istl/owneroverlapcopy.hh>
27#include <dune/istl/scalarproducts.hh>
28#include <dune/istl/paamg/amg.hh>
29#include <dune/istl/paamg/pinfo.hh>
30
40
41#include <dune/istl/foreach.hh>
42
44
51template<class M>
52constexpr std::size_t preconditionerBlockLevel() noexcept
53{
55}
56
57template<template<class,class,class,int> class Preconditioner, int blockLevel = 1>
59{
60public:
61 template<class OI, class M>
62 auto operator() (OI opInfo, const M& matrix, const Dune::ParameterTree& config)
63 {
64#if DUNE_VERSION_LT(DUNE_ISTL,2,11)
65 using Matrix = typename Dune::TypeListElement<0, decltype(opInfo)>::type;
66 using Domain = typename Dune::TypeListElement<1, decltype(opInfo)>::type;
67 using Range = typename Dune::TypeListElement<2, decltype(opInfo)>::type;
68#else
69 using OpInfo = std::decay_t<decltype(opInfo)>;
70 using Matrix = typename OpInfo::matrix_type;
71 using Domain = typename OpInfo::domain_type;
72 using Range = typename OpInfo::range_type;
73#endif
74 std::shared_ptr<Dune::Preconditioner<Domain, Range>> preconditioner
75 = std::make_shared<Preconditioner<Matrix, Domain, Range, blockLevel>>(matrix, config);
76 return preconditioner;
77 }
78};
79
80template<template<class,class,class> class Preconditioner>
82{
83 template<class OI, class M>
84 auto operator() (OI opInfo, const M& matrix, const Dune::ParameterTree& config)
85 {
86#if DUNE_VERSION_LT(DUNE_ISTL,2,11)
87 using Matrix = typename Dune::TypeListElement<0, decltype(opInfo)>::type;
88 using Domain = typename Dune::TypeListElement<1, decltype(opInfo)>::type;
89 using Range = typename Dune::TypeListElement<2, decltype(opInfo)>::type;
90#else
91 using OpInfo = std::decay_t<decltype(opInfo)>;
92 using Matrix = typename OpInfo::matrix_type;
93 using Domain = typename OpInfo::domain_type;
94 using Range = typename OpInfo::range_type;
95#endif
96 std::shared_ptr<Dune::Preconditioner<Domain, Range>> preconditioner
97 = std::make_shared<Preconditioner<Matrix, Domain, Range>>(matrix, config);
98 return preconditioner;
99 }
100};
101
102using IstlAmgPreconditionerFactory = Dune::AMGCreator;
103
104template<class M, bool convert = false>
105struct MatrixForSolver { using type = M; };
106
107template<class M>
108struct MatrixForSolver<M, true>
109{ using type = std::decay_t<decltype(MatrixConverter<M>::multiTypeToBCRSMatrix(std::declval<M>()))>; };
110
111template<class V, bool convert = false>
112struct VectorForSolver { using type = V; };
113
114template<class V>
115struct VectorForSolver<V, true>
116{ using type = std::decay_t<decltype(VectorConverter<V>::multiTypeToBlockVector(std::declval<V>()))>; };
117
118template<class LSTraits, class LATraits, bool convert, bool parallel = LSTraits::canCommunicate>
120
121template<class LSTraits, class LATraits, bool convert>
122struct MatrixOperator<LSTraits, LATraits, convert, true>
123{
126#if HAVE_MPI
127 using type = std::variant<
128 std::shared_ptr<typename LSTraits::template Sequential<M, V>::LinearOperator>,
129 std::shared_ptr<typename LSTraits::template ParallelOverlapping<M, V>::LinearOperator>,
130 std::shared_ptr<typename LSTraits::template ParallelNonoverlapping<M, V>::LinearOperator>
131 >;
132#else
133 using type = std::variant<
134 std::shared_ptr<typename LSTraits::template Sequential<M, V>::LinearOperator>
135 >;
136#endif
137};
138
139template<class LSTraits, class LATraits, bool convert>
140struct MatrixOperator<LSTraits, LATraits, convert, false>
141{
144 using type = std::variant<
145 std::shared_ptr<typename LSTraits::template Sequential<M, V>::LinearOperator>
146 >;
147};
148
149} // end namespace Dumux::Detail::IstlSolvers
150
151namespace Dumux::Detail {
152
153struct IstlSolverResult : public Dune::InverseOperatorResult
154{
155 IstlSolverResult() = default;
158
159 IstlSolverResult(const Dune::InverseOperatorResult& o) : InverseOperatorResult(o) {}
160 IstlSolverResult(Dune::InverseOperatorResult&& o) : InverseOperatorResult(std::move(o)) {}
161
162 operator bool() const { return this->converged; }
163};
164
169template<class LinearSolverTraits, class LinearAlgebraTraits,
170 class InverseOperator, class PreconditionerFactory,
171 bool convertMultiTypeLATypes = false>
173{
174 using Matrix = typename LinearAlgebraTraits::Matrix;
175 using XVector = typename LinearAlgebraTraits::Vector;
176 using BVector = typename LinearAlgebraTraits::Vector;
177 using Scalar = typename InverseOperator::real_type;
178
179 using ScalarProduct = Dune::ScalarProduct<typename InverseOperator::domain_type>;
180
181 static constexpr bool convertMultiTypeVectorAndMatrix
182 = convertMultiTypeLATypes && isMultiTypeBlockVector<XVector>::value;
186 // a variant type that can hold sequential, overlapping, and non-overlapping operators
187 using MatrixOperatorHolder = typename Detail::IstlSolvers::MatrixOperator<
188 LinearSolverTraits, LinearAlgebraTraits, convertMultiTypeVectorAndMatrix
189 >::type;
190
191#if HAVE_MPI
192 using Comm = Dune::OwnerOverlapCopyCommunication<Dune::bigunsignedint<96>, int>;
193 using ParallelHelper = ParallelISTLHelper<LinearSolverTraits>;
194#endif
195
196 using ParameterInitializer = std::variant<std::string, Dune::ParameterTree>;
197public:
198
202 IstlIterativeLinearSolver(const ParameterInitializer& params = "")
203 {
204 if (Dune::MPIHelper::getCommunication().size() > 1)
205 DUNE_THROW(Dune::InvalidStateException, "Using sequential constructor for parallel run. Use signature with gridView and dofMapper!");
206
207 initializeParameters_(params);
208 solverCategory_ = Dune::SolverCategory::sequential;
209 scalarProduct_ = std::make_shared<ScalarProduct>();
210 }
211
215 template <class GridView, class DofMapper>
216 IstlIterativeLinearSolver(const GridView& gridView,
217 const DofMapper& dofMapper,
218 const ParameterInitializer& params = "")
219 {
220 initializeParameters_(params, gridView.comm());
221#if HAVE_MPI
222 solverCategory_ = Detail::solverCategory<LinearSolverTraits>(gridView);
223 if constexpr (LinearSolverTraits::canCommunicate)
224 {
225
226 if (solverCategory_ != Dune::SolverCategory::sequential)
227 {
228 parallelHelper_ = std::make_shared<ParallelISTLHelper<LinearSolverTraits>>(gridView, dofMapper);
229 communication_ = std::make_shared<Comm>(gridView.comm(), solverCategory_);
230 scalarProduct_ = Dune::createScalarProduct<XVector>(*communication_, solverCategory_);
231 parallelHelper_->createParallelIndexSet(*communication_);
232 }
233 else
234 scalarProduct_ = std::make_shared<ScalarProduct>();
235 }
236 else
237 scalarProduct_ = std::make_shared<ScalarProduct>();
238#else
239 solverCategory_ = Dune::SolverCategory::sequential;
240 scalarProduct_ = std::make_shared<ScalarProduct>();
241#endif
242 }
243
244#if HAVE_MPI
248 template <class GridView, class DofMapper>
249 IstlIterativeLinearSolver(std::shared_ptr<Comm> communication,
250 std::shared_ptr<ScalarProduct> scalarProduct,
251 const GridView& gridView,
252 const DofMapper& dofMapper,
253 const ParameterInitializer& params = "")
254 {
255 initializeParameters_(params, gridView.comm());
256 solverCategory_ = Detail::solverCategory(gridView);
257 scalarProduct_ = scalarProduct;
258 communication_ = communication;
259 if constexpr (LinearSolverTraits::canCommunicate)
260 {
261 if (solverCategory_ != Dune::SolverCategory::sequential)
262 {
263 parallelHelper_ = std::make_shared<ParallelISTLHelper<LinearSolverTraits>>(gridView, dofMapper);
264 parallelHelper_->createParallelIndexSet(communication);
265 }
266 }
267 }
268#endif
269
273 IstlSolverResult solve(Matrix& A, XVector& x, BVector& b)
274 { return solveSequentialOrParallel_(A, x, b); }
275
279 void setMatrix(std::shared_ptr<Matrix> A)
280 {
281 linearOperator_ = makeParallelOrSequentialLinearOperator_(std::move(A));
282 solver_ = constructPreconditionedSolver_(linearOperator_);
283 }
284
289 void setMatrix(Matrix& A)
290 { setMatrix(Dune::stackobject_to_shared_ptr(A)); }
291
295 IstlSolverResult solve(XVector& x, BVector& b) const
296 {
297 if (!solver_)
298 DUNE_THROW(Dune::InvalidStateException, "Called solve(x, b) but no linear operator has been set");
299
300 return solveSequentialOrParallel_(x, b, *solver_);
301 }
302
306 Scalar norm(const XVector& x) const
307 {
308#if HAVE_MPI
309 if constexpr (LinearSolverTraits::canCommunicate)
310 {
311 if (solverCategory_ == Dune::SolverCategory::nonoverlapping)
312 {
313 auto y(x); // make a copy because the vector needs to be made consistent
314 using GV = typename LinearSolverTraits::GridView;
315 using DM = typename LinearSolverTraits::DofMapper;
316 if constexpr (requires { LinearSolverTraits::dofCodims; })
317 {
318 MultiCodimParallelVectorHelper<GV, DM> vectorHelper(parallelHelper_->gridView(), parallelHelper_->dofMapper());
319 vectorHelper.makeNonOverlappingConsistent(y, LinearSolverTraits::dofCodims);
320 }
321 else
322 {
323 ParallelVectorHelper<GV, DM, LinearSolverTraits::dofCodim> vectorHelper(parallelHelper_->gridView(), parallelHelper_->dofMapper());
324 vectorHelper.makeNonOverlappingConsistent(y);
325 }
326 return scalarProduct_->norm(y);
327 }
328 }
329#endif
330 if constexpr (convertMultiTypeVectorAndMatrix)
331 {
333 return scalarProduct_->norm(y);
334 }
335 else
336 return scalarProduct_->norm(x);
337 }
338
342 const std::string& name() const
343 {
344 return name_;
345 }
346
350 void setResidualReduction(double residReduction)
351 {
352 params_["reduction"] = std::to_string(residReduction);
353
354 // reconstruct the solver with new parameters
355 if (solver_)
356 solver_ = constructPreconditionedSolver_(linearOperator_);
357 }
358
362 void setMaxIter(std::size_t maxIter)
363 {
364 params_["maxit"] = std::to_string(maxIter);
365
366 // reconstruct the solver with new parameters
367 if (solver_)
368 solver_ = constructPreconditionedSolver_(linearOperator_);
369 }
370
376 void setParams(const ParameterInitializer& params)
377 {
378#if HAVE_MPI
379 if (communication_)
380 initializeParameters_(params, communication_->communicator());
381 else
382 initializeParameters_(params);
383#else
384 initializeParameters_(params);
385#endif
386
387 // reconstruct the solver with new parameters
388 if (solver_)
389 solver_ = constructPreconditionedSolver_(linearOperator_);
390 }
391
392private:
393
394 void initializeParameters_(const ParameterInitializer& params)
395 {
396 if (std::holds_alternative<std::string>(params))
397 params_ = Dumux::LinearSolverParameters<LinearSolverTraits>::createParameterTree(std::get<std::string>(params));
398 else
399 params_ = std::get<Dune::ParameterTree>(params);
400 }
401
402 template <class Comm>
403 void initializeParameters_(const ParameterInitializer& params, const Comm& comm)
404 {
405 initializeParameters_(params);
406
407 // disable verbose output on all ranks except rank 0
408 if (comm.rank() != 0)
410 }
411
412 MatrixOperatorHolder makeSequentialLinearOperator_(std::shared_ptr<Matrix> A)
413 {
414 using SequentialTraits = typename LinearSolverTraits::template Sequential<MatrixForSolver, XVectorForSolver>;
415 if constexpr (convertMultiTypeVectorAndMatrix)
416 {
417 // create the BCRS matrix the IterativeSolver backend can handle
418 auto M = std::make_shared<MatrixForSolver>(MatrixConverter<Matrix>::multiTypeToBCRSMatrix(*A));
419 return std::make_shared<typename SequentialTraits::LinearOperator>(M);
420 }
421 else
422 {
423 return std::make_shared<typename SequentialTraits::LinearOperator>(A);
424 }
425 }
426
427 template<class ParallelTraits>
428 MatrixOperatorHolder makeParallelLinearOperator_(std::shared_ptr<Matrix> A, ParallelTraits = {})
429 {
430#if HAVE_MPI
431 // make matrix consistent
433 return std::make_shared<typename ParallelTraits::LinearOperator>(std::move(A), *communication_);
434#else
435 DUNE_THROW(Dune::InvalidStateException, "Calling makeParallelLinearOperator for sequential run");
436#endif
437 }
438
439 MatrixOperatorHolder makeParallelOrSequentialLinearOperator_(std::shared_ptr<Matrix> A)
440 {
441 return executeSequentialOrParallel_(
442 [&]{ return makeSequentialLinearOperator_(std::move(A)); },
443 [&](auto traits){ return makeParallelLinearOperator_(std::move(A), traits); }
444 );
445 }
446
447 MatrixOperatorHolder makeSequentialLinearOperator_(Matrix& A)
448 { return makeSequentialLinearOperator_(Dune::stackobject_to_shared_ptr<Matrix>(A)); }
449
450 MatrixOperatorHolder makeParallelOrSequentialLinearOperator_(Matrix& A)
451 { return makeParallelOrSequentialLinearOperator_(Dune::stackobject_to_shared_ptr<Matrix>(A)); }
452
453 template<class ParallelTraits>
454 MatrixOperatorHolder makeParallelLinearOperator_(Matrix& A, ParallelTraits = {})
455 { return makeParallelLinearOperator_<ParallelTraits>(Dune::stackobject_to_shared_ptr<Matrix>(A)); }
456
457 IstlSolverResult solveSequential_(Matrix& A, XVector& x, BVector& b)
458 {
459 // construct solver from linear operator
460 auto linearOperatorHolder = makeSequentialLinearOperator_(A);
461 auto solver = constructPreconditionedSolver_(linearOperatorHolder);
462
463 return solveSequential_(x, b, *solver);
464 }
465
466 IstlSolverResult solveSequential_(XVector& x, BVector& b, InverseOperator& solver) const
467 {
468 Dune::InverseOperatorResult result;
469 if constexpr (convertMultiTypeVectorAndMatrix)
470 {
471 // create the vector the IterativeSolver backend can handle
472 BVectorForSolver bTmp = VectorConverter<BVector>::multiTypeToBlockVector(b);
473
474 // create a block vector to which the linear solver writes the solution
475 XVectorForSolver y(bTmp.size());
476
477 // solve linear system
478 solver.apply(y, bTmp, result);
479
480 // copy back the result y into x
481 if (result.converged)
483 }
484 else
485 {
486 // solve linear system
487 solver.apply(x, b, result);
488 }
489
490 return result;
491 }
492
493 IstlSolverResult solveSequentialOrParallel_(Matrix& A, XVector& x, BVector& b)
494 {
495 return executeSequentialOrParallel_(
496 [&]{ return solveSequential_(A, x, b); },
497 [&](auto traits){ return solveParallel_(A, x, b, traits); }
498 );
499 }
500
501 IstlSolverResult solveSequentialOrParallel_(XVector& x, BVector& b, InverseOperator& solver) const
502 {
503 return executeSequentialOrParallel_(
504 [&]{ return solveSequential_(x, b, solver); },
505 [&](auto traits){ return solveParallel_(x, b, solver, traits); }
506 );
507 }
508
509 template<class ParallelTraits>
510 IstlSolverResult solveParallel_(Matrix& A, XVector& x, BVector& b, ParallelTraits = {})
511 {
512 // construct solver from linear operator
513 auto linearOperatorHolder = makeParallelLinearOperator_<ParallelTraits>(A);
514 auto solver = constructPreconditionedSolver_(linearOperatorHolder);
515 return solveParallel_<ParallelTraits>(x, b, *solver);
516 }
517
518 template<class ParallelTraits>
519 IstlSolverResult solveParallel_(XVector& x, BVector& b, InverseOperator& solver, ParallelTraits = {}) const
520 {
521#if HAVE_MPI
522 // make right hand side consistent
524
525 // solve linear system
526 Dune::InverseOperatorResult result;
527 solver.apply(x, b, result);
528 return result;
529#else
530 DUNE_THROW(Dune::InvalidStateException, "Calling makeParallelLinearOperator for sequential run");
531#endif
532 }
533
534
535 std::shared_ptr<InverseOperator> constructPreconditionedSolver_(MatrixOperatorHolder& ops)
536 {
537 return std::visit([&](auto&& op)
538 {
539 using LinearOperator = typename std::decay_t<decltype(op)>::element_type;
540 const auto& params = params_.sub("preconditioner");
541 using Prec = Dune::Preconditioner<typename LinearOperator::domain_type, typename LinearOperator::range_type>;
542#if DUNE_VERSION_GTE(DUNE_ISTL,2,11)
543 using OpTraits = Dune::OperatorTraits<LinearOperator>;
544 std::shared_ptr<Prec> prec = PreconditionerFactory{}(OpTraits{}, op, params);
545#else
546 using TL = Dune::TypeList<typename LinearOperator::matrix_type, typename LinearOperator::domain_type, typename LinearOperator::range_type>;
547 std::shared_ptr<Prec> prec = PreconditionerFactory{}(TL{}, op, params);
548#endif
549
550#if HAVE_MPI
551#if DUNE_VERSION_LT(DUNE_ISTL,2,11)
552 if (prec->category() != op->category() && prec->category() == Dune::SolverCategory::sequential)
553 prec = Dune::wrapPreconditioner4Parallel(prec, op);
554#else
555 if constexpr (OpTraits::isParallel)
556 {
557 using Comm = typename OpTraits::comm_type;
558 const Comm& comm = OpTraits::getCommOrThrow(op);
559 if (op->category() == Dune::SolverCategory::overlapping && prec->category() == Dune::SolverCategory::sequential)
560 prec = std::make_shared<Dune::BlockPreconditioner<typename OpTraits::domain_type, typename OpTraits::range_type,Comm> >(prec, comm);
561 else if (op->category() == Dune::SolverCategory::nonoverlapping && prec->category() == Dune::SolverCategory::sequential)
562 prec = std::make_shared<Dune::NonoverlappingBlockPreconditioner<Comm, Prec> >(prec, comm);
563 }
564#endif
565#endif
566 return std::make_shared<InverseOperator>(op, scalarProduct_, prec, params_);
567 }, ops);
568 }
569
570 template<class Seq, class Par>
571 decltype(auto) executeSequentialOrParallel_(Seq&& sequentialAction, Par&& parallelAction) const
572 {
573#if HAVE_MPI
574 // For Dune::MultiTypeBlockMatrix there is currently no generic way
575 // of handling parallelism, we therefore can only solve these types of systems sequentially
576 if constexpr (isMultiTypeBlockMatrix<Matrix>::value || !LinearSolverTraits::canCommunicate)
577 return sequentialAction();
578 else
579 {
580 switch (solverCategory_)
581 {
582 case Dune::SolverCategory::sequential:
583 return sequentialAction();
584 case Dune::SolverCategory::nonoverlapping:
585 using NOTraits = typename LinearSolverTraits::template ParallelNonoverlapping<Matrix, XVector>;
586 return parallelAction(NOTraits{});
587 case Dune::SolverCategory::overlapping:
588 using OTraits = typename LinearSolverTraits::template ParallelOverlapping<Matrix, XVector>;
589 return parallelAction(OTraits{});
590 default: DUNE_THROW(Dune::InvalidStateException, "Unknown solver category");
591 }
592 }
593#else
594 return sequentialAction();
595#endif
596 }
597
598#if HAVE_MPI
599 std::shared_ptr<const ParallelHelper> parallelHelper_;
600 std::shared_ptr<Comm> communication_;
601#endif
602
603 Dune::SolverCategory::Category solverCategory_;
604 std::shared_ptr<ScalarProduct> scalarProduct_;
605
606 // for stored solvers (reuse matrix)
607 MatrixOperatorHolder linearOperator_;
608 // for stored solvers (reuse matrix)
609 std::shared_ptr<InverseOperator> solver_;
610
611 Dune::ParameterTree params_;
612 std::string name_;
613};
614
615} // end namespace Dumux::Detail
616
617namespace Dumux {
618
635template<class LSTraits, class LATraits>
637 Detail::IstlIterativeLinearSolver<LSTraits, LATraits,
638 Dune::BiCGSTABSolver<typename LATraits::SingleTypeVector>,
640 // the Dune::ILU preconditioners don't accept multi-type matrices
641 /*convert multi-type istl types?*/ true
642 >;
643
660template<class LSTraits, class LATraits>
662 Detail::IstlIterativeLinearSolver<LSTraits, LATraits,
663 Dune::RestartedGMResSolver<typename LATraits::SingleTypeVector>,
665 // the Dune::ILU preconditioners don't accept multi-type matrices
666 /*convert multi-type istl types?*/ true
667 >;
668
686template<class LSTraits, class LATraits>
688 Detail::IstlIterativeLinearSolver<LSTraits, LATraits,
689 Dune::BiCGSTABSolver<typename LATraits::Vector>,
691 >;
692
709template<class LSTraits, class LATraits>
711 Detail::IstlIterativeLinearSolver<LSTraits, LATraits,
712 Dune::CGSolver<typename LATraits::Vector>,
714 >;
715
729template<class LSTraits, class LATraits>
731 Detail::IstlIterativeLinearSolver<LSTraits, LATraits,
732 Dune::BiCGSTABSolver<typename LATraits::SingleTypeVector>,
734 // the AMG preconditioner doesn't accept multi-type matrices
735 /*convert multi-type istl types?*/ true
736 >;
737
750template<class LSTraits, class LATraits>
752 Detail::IstlIterativeLinearSolver<LSTraits, LATraits,
753 Dune::CGSolver<typename LATraits::SingleTypeVector>,
755 // the AMG preconditioner doesn't accept multi-type matrices
756 /*convert multi-type istl types?*/ true
757 >;
758
773template<class LSTraits, class LATraits>
775 Detail::IstlIterativeLinearSolver<LSTraits, LATraits,
776 Dune::BiCGSTABSolver<typename LATraits::Vector>,
778 >;
779
780} // end namespace Dumux
781
782namespace Dumux::Detail {
783
788template<class LSTraits, class LATraits, template<class M> class Solver,
789 bool convertMultiTypeVectorAndMatrix = isMultiTypeBlockVector<typename LATraits::Vector>::value>
791{
792 using Matrix = typename LATraits::Matrix;
793 using XVector = typename LATraits::Vector;
794 using BVector = typename LATraits::Vector;
795
799 using InverseOperator = Dune::InverseOperator<XVectorForSolver, BVectorForSolver>;
800public:
802
806 IstlSolverResult solve(const Matrix& A, XVector& x, const BVector& b)
807 {
808 return solve_(A, x, b);
809 }
810
814 IstlSolverResult solve(XVector& x, const BVector& b)
815 {
816 if (!solver_)
817 DUNE_THROW(Dune::InvalidStateException, "Called solve(x, b) but no linear operator has been set");
818
819 return solve_(x, b, *solver_);
820 }
821
825 void setMatrix(std::shared_ptr<Matrix> A)
826 {
827 if constexpr (convertMultiTypeVectorAndMatrix)
828 matrix_ = std::make_shared<MatrixForSolver>(MatrixConverter<Matrix>::multiTypeToBCRSMatrix(A));
829 else
830 matrix_ = A;
831
832 solver_ = std::make_shared<Solver<MatrixForSolver>>(*matrix_);
833 }
834
839 void setMatrix(Matrix& A)
840 { setMatrix(Dune::stackobject_to_shared_ptr(A)); }
841
845 std::string name() const
846 {
847 return "Direct solver";
848 }
849
850private:
851 IstlSolverResult solve_(const Matrix& A, XVector& x, const BVector& b)
852 {
853 // support dune-istl multi-type block vector/matrix by copying
854 if constexpr (convertMultiTypeVectorAndMatrix)
855 {
857 Solver<MatrixForSolver> solver(AA, this->verbosity() > 0);
858 return solve_(x, b, solver);
859 }
860 else
861 {
862 Solver<MatrixForSolver> solver(A, this->verbosity() > 0);
863 return solve_(x, b, solver);
864 }
865 }
866
867 IstlSolverResult solve_(XVector& x, const BVector& b, InverseOperator& solver) const
868 {
869 Dune::InverseOperatorResult result;
870
871 if constexpr (convertMultiTypeVectorAndMatrix)
872 {
874 XVectorForSolver xx(bb.size());
875 solver.apply(xx, bb, result);
876 checkResult_(xx, result);
877 if (result.converged)
879 return result;
880 }
881 else
882 {
883 BVectorForSolver bTmp(b);
884 solver.apply(x, bTmp, result);
885 checkResult_(x, result);
886 return result;
887 }
888 }
889
890
891 void checkResult_(XVectorForSolver& x, Dune::InverseOperatorResult& result) const
892 {
893 flatVectorForEach(x, [&](auto&& entry, std::size_t){
894 using std::isnan, std::isinf;
895 if (isnan(entry) || isinf(entry))
896 result.converged = false;
897 });
898 }
899
901 std::shared_ptr<MatrixForSolver> matrix_;
903 std::shared_ptr<InverseOperator> solver_;
904};
905
906} // end namespace Dumux::Detail
907
908#if HAVE_SUPERLU
909#include <dune/istl/superlu.hh>
910
911namespace Dumux {
912
921template<class LSTraits, class LATraits>
923
924} // end namespace Dumux
925
926#endif // HAVE_SUPERLU
927
928#if HAVE_UMFPACK
929#include <dune/istl/umfpack.hh>
930
931namespace Dumux {
932
941template<class LSTraits, class LATraits>
943
944} // end namespace Dumux
945
946#endif // HAVE_UMFPACK
947
948#endif
Direct dune-istl linear solvers.
Definition istlsolvers.hh:791
void setMatrix(std::shared_ptr< Matrix > A)
Set the matrix A of the linear system Ax = b for reuse.
Definition istlsolvers.hh:825
IstlSolverResult solve(const Matrix &A, XVector &x, const BVector &b)
Solve the linear system Ax = b.
Definition istlsolvers.hh:806
std::string name() const
name of the linear solver
Definition istlsolvers.hh:845
void setMatrix(Matrix &A)
Set the matrix A of the linear system Ax = b for reuse.
Definition istlsolvers.hh:839
LinearSolver(const std::string &paramGroup="")
Construct the solver.
Definition solver.hh:43
IstlSolverResult solve(XVector &x, const BVector &b)
Solve the linear system Ax = b using the matrix set with setMatrix.
Definition istlsolvers.hh:814
Standard dune-istl iterative linear solvers.
Definition istlsolvers.hh:173
IstlIterativeLinearSolver(const ParameterInitializer &params="")
Constructor for sequential solvers.
Definition istlsolvers.hh:202
void setMatrix(Matrix &A)
Set the matrix A of the linear system Ax = b for reuse.
Definition istlsolvers.hh:289
void setResidualReduction(double residReduction)
Set the residual reduction tolerance.
Definition istlsolvers.hh:350
void setMaxIter(std::size_t maxIter)
Set the maximum number of linear solver iterations.
Definition istlsolvers.hh:362
IstlIterativeLinearSolver(const GridView &gridView, const DofMapper &dofMapper, const ParameterInitializer &params="")
Constructor for parallel and sequential solvers.
Definition istlsolvers.hh:216
const std::string & name() const
The name of the linear solver.
Definition istlsolvers.hh:342
void setParams(const ParameterInitializer &params)
Set the linear solver parameters.
Definition istlsolvers.hh:376
IstlSolverResult solve(Matrix &A, XVector &x, BVector &b)
Solve the linear system Ax = b.
Definition istlsolvers.hh:273
void setMatrix(std::shared_ptr< Matrix > A)
Set the matrix A of the linear system Ax = b for reuse.
Definition istlsolvers.hh:279
IstlSolverResult solve(XVector &x, BVector &b) const
Solve the linear system Ax = b where A has been set with setMatrix.
Definition istlsolvers.hh:295
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:249
Scalar norm(const XVector &x) const
Compute the 2-norm of vector x.
Definition istlsolvers.hh:306
auto operator()(OI opInfo, const M &matrix, const Dune::ParameterTree &config)
Definition istlsolvers.hh:62
LinearSolver(const std::string &paramGroup="")
Construct the solver.
Definition solver.hh:43
int verbosity() const
the verbosity level
Definition solver.hh:82
static void disableVerbosity(Dune::ParameterTree &params)
Definition linearsolverparameters.hh:97
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
static auto multiTypeToBCRSMatrix(const MultiTypeBlockMatrix &A)
Converts the matrix to a type the IterativeSolverBackend can handle.
Definition matrixconverter.hh:46
Definition parallelhelpers.hh:580
void makeNonOverlappingConsistent(Dune::BlockVector< Block, Alloc > &v, const std::bitset< numCodims > &activeCodims) const
Make a vector consistent for non-overlapping domain decomposition methods.
Definition parallelhelpers.hh:588
Definition parallelhelpers.hh:522
void makeNonOverlappingConsistent(Dune::BlockVector< Block, Alloc > &v) const
Make a vector consistent for non-overlapping domain decomposition methods.
Definition parallelhelpers.hh:530
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
Detail::IstlIterativeLinearSolver< LSTraits, LATraits, Dune::BiCGSTABSolver< typename LATraits::Vector >, Detail::IstlSolvers::IstlDefaultBlockLevelPreconditionerFactory< Dumux::SeqUzawa > > UzawaBiCGSTABIstlSolver
An Uzawa preconditioned BiCGSTAB solver using dune-istl.
Definition istlsolvers.hh:774
constexpr std::size_t preconditionerBlockLevel() noexcept
Returns the block level for the preconditioner for a given matrix.
Definition istlsolvers.hh:52
Detail::IstlIterativeLinearSolver< LSTraits, LATraits, Dune::RestartedGMResSolver< typename LATraits::SingleTypeVector >, Detail::IstlSolvers::IstlDefaultBlockLevelPreconditionerFactory< Dune::SeqILU >, true > ILURestartedGMResIstlSolver
An ILU preconditioned GMres solver using dune-istl.
Definition istlsolvers.hh:661
Detail::IstlIterativeLinearSolver< LSTraits, LATraits, Dune::BiCGSTABSolver< typename LATraits::Vector >, Detail::IstlSolvers::IstlDefaultBlockLevelPreconditionerFactory< Dune::SeqSSOR > > SSORBiCGSTABIstlSolver
An SSOR-preconditioned BiCGSTAB solver using dune-istl.
Definition istlsolvers.hh:687
Detail::IstlIterativeLinearSolver< LSTraits, LATraits, Dune::BiCGSTABSolver< typename LATraits::SingleTypeVector >, Detail::IstlSolvers::IstlDefaultBlockLevelPreconditionerFactory< Dune::SeqILU >, true > ILUBiCGSTABIstlSolver
An ILU preconditioned BiCGSTAB solver using dune-istl.
Definition istlsolvers.hh:636
Detail::IstlIterativeLinearSolver< LSTraits, LATraits, Dune::BiCGSTABSolver< typename LATraits::SingleTypeVector >, Detail::IstlSolvers::IstlAmgPreconditionerFactory, true > AMGBiCGSTABIstlSolver
An AMG preconditioned BiCGSTAB solver using dune-istl.
Definition istlsolvers.hh:730
Detail::ParallelISTLHelperImpl< LinearSolverTraits, LinearSolverTraits::canCommunicate > ParallelISTLHelper
A parallel helper class providing a parallel decomposition of all degrees of freedom.
Definition parallelhelpers.hh:514
Detail::IstlIterativeLinearSolver< LSTraits, LATraits, Dune::CGSolver< typename LATraits::SingleTypeVector >, Detail::IstlSolvers::IstlAmgPreconditionerFactory, true > AMGCGIstlSolver
An AMG preconditioned CG solver using dune-istl.
Definition istlsolvers.hh:751
Detail::IstlIterativeLinearSolver< LSTraits, LATraits, Dune::CGSolver< typename LATraits::Vector >, Detail::IstlSolvers::IstlDefaultBlockLevelPreconditionerFactory< Dune::SeqSSOR > > SSORCGIstlSolver
An SSOR-preconditioned CG solver using dune-istl.
Definition istlsolvers.hh:710
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:43
Dune::AMGCreator IstlAmgPreconditionerFactory
Definition istlsolvers.hh:102
Definition linearalgebratraits.hh:21
Definition cvfelocalresidual.hh:25
Dune::SolverCategory::Category solverCategory(const GridView &gridView)
Definition solvercategory.hh:20
Definition adapt.hh:17
void prepareMatrixParallel(Matrix &A, ParallelHelper &pHelper)
Prepare a matrix for parallel solvers.
Definition parallelhelpers.hh:1420
LinearSolverTraitsImpl< GridGeometry, typename GridGeometry::DiscretizationMethod > LinearSolverTraits
The type traits required for using the IstlFactoryBackend.
Definition linearsolvertraits.hh:39
void prepareVectorParallel(Vector &b, ParallelHelper &pHelper)
Prepare a vector for parallel solvers.
Definition parallelhelpers.hh:1451
Provides a helper class for nonoverlapping decomposition.
Dumux preconditioners for iterative solvers.
Base class for linear solvers.
Solver category.
Definition istlsolvers.hh:154
IstlSolverResult(IstlSolverResult &&)=default
IstlSolverResult(const Dune::InverseOperatorResult &o)
Definition istlsolvers.hh:159
IstlSolverResult(Dune::InverseOperatorResult &&o)
Definition istlsolvers.hh:160
IstlSolverResult(const IstlSolverResult &)=default
std::decay_t< decltype(MatrixConverter< M >::multiTypeToBCRSMatrix(std::declval< M >()))> type
Definition istlsolvers.hh:109
Definition istlsolvers.hh:105
M type
Definition istlsolvers.hh:105
typename VectorForSolver< typename LATraits::Vector, convert >::type V
Definition istlsolvers.hh:143
std::variant< std::shared_ptr< typename LSTraits::template Sequential< M, V >::LinearOperator > > type
Definition istlsolvers.hh:144
typename MatrixForSolver< typename LATraits::Matrix, convert >::type M
Definition istlsolvers.hh:142
typename VectorForSolver< typename LATraits::Vector, convert >::type V
Definition istlsolvers.hh:125
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:127
typename MatrixForSolver< typename LATraits::Matrix, convert >::type M
Definition istlsolvers.hh:124
Definition istlsolvers.hh:119
std::decay_t< decltype(VectorConverter< V >::multiTypeToBlockVector(std::declval< V >()))> type
Definition istlsolvers.hh:116
Definition istlsolvers.hh:112
V type
Definition istlsolvers.hh:112
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.