12#ifndef DUMUX_LINEAR_ISTL_SOLVERS_HH
13#define DUMUX_LINEAR_ISTL_SOLVERS_HH
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>
40#include <dune/istl/foreach.hh>
56template<
template<
class,
class,
class,
int>
class Preconditioner,
int blockLevel = 1>
60 template<
class TL,
class M>
61 auto operator() (TL typeList,
const M& matrix,
const Dune::ParameterTree& config)
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;
72template<
template<
class,
class,
class>
class Preconditioner>
75 template<
class TL,
class M>
76 auto operator() (TL typeList,
const M& matrix,
const Dune::ParameterTree& config)
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;
89template<
class M,
bool convert = false>
96template<
class V,
bool convert = false>
103template<
class LSTraits,
class LATraits,
bool convert,
bool parallel = LSTraits::canCommunicate>
106template<
class LSTraits,
class LATraits,
bool convert>
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>
118 using type = std::variant<
119 std::shared_ptr<typename LSTraits::template Sequential<M, V>::LinearOperator>
124template<
class LSTraits,
class LATraits,
bool convert>
130 std::shared_ptr<typename LSTraits::template Sequential<M, V>::LinearOperator>
147 operator bool()
const {
return this->converged; }
155 class InverseOperator,
class PreconditionerFactory,
156 bool convertMultiTypeLATypes =
false>
162 using Scalar =
typename InverseOperator::real_type;
164 using ScalarProduct = Dune::ScalarProduct<typename InverseOperator::domain_type>;
166 static constexpr bool convertMultiTypeVectorAndMatrix
177 using Comm = Dune::OwnerOverlapCopyCommunication<Dune::bigunsignedint<96>,
int>;
181 using ParameterInitializer = std::variant<std::string, Dune::ParameterTree>;
189 if (Dune::MPIHelper::getCommunication().size() > 1)
190 DUNE_THROW(Dune::InvalidStateException,
"Using sequential constructor for parallel run. Use signature with gridView and dofMapper!");
192 initializeParameters_(params);
193 solverCategory_ = Dune::SolverCategory::sequential;
194 scalarProduct_ = std::make_shared<ScalarProduct>();
200 template <
class Gr
idView,
class DofMapper>
202 const DofMapper& dofMapper,
203 const ParameterInitializer& params =
"")
205 initializeParameters_(params);
207 solverCategory_ = Detail::solverCategory<LinearSolverTraits>(gridView);
211 if (solverCategory_ != Dune::SolverCategory::sequential)
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_);
219 scalarProduct_ = std::make_shared<ScalarProduct>();
222 scalarProduct_ = std::make_shared<ScalarProduct>();
224 solverCategory_ = Dune::SolverCategory::sequential;
225 scalarProduct_ = std::make_shared<ScalarProduct>();
233 template <
class Gr
idView,
class DofMapper>
235 std::shared_ptr<ScalarProduct> scalarProduct,
236 const GridView& gridView,
237 const DofMapper& dofMapper,
238 const ParameterInitializer& params =
"")
240 initializeParameters_(params);
242 scalarProduct_ = scalarProduct;
243 communication_ = communication;
246 if (solverCategory_ != Dune::SolverCategory::sequential)
248 parallelHelper_ = std::make_shared<ParallelISTLHelper<LinearSolverTraits>>(gridView, dofMapper);
249 parallelHelper_->createParallelIndexSet(communication);
259 {
return solveSequentialOrParallel_(A, x, b); }
266 linearOperator_ = makeParallelOrSequentialLinearOperator_(std::move(A));
267 solver_ = constructPreconditionedSolver_(linearOperator_);
275 {
setMatrix(Dune::stackobject_to_shared_ptr(A)); }
283 DUNE_THROW(Dune::InvalidStateException,
"Called solve(x, b) but no linear operator has been set");
285 return solveSequentialOrParallel_(x, b, *solver_);
291 Scalar
norm(
const XVector& x)
const
296 if (solverCategory_ == Dune::SolverCategory::nonoverlapping)
299 using GV =
typename LinearSolverTraits::GridView;
300 using DM =
typename LinearSolverTraits::DofMapper;
303 return scalarProduct_->norm(y);
307 if constexpr (convertMultiTypeVectorAndMatrix)
310 return scalarProduct_->norm(y);
313 return scalarProduct_->norm(x);
319 const std::string&
name()
const
329 params_[
"reduction"] = std::to_string(residReduction);
333 solver_ = constructPreconditionedSolver_(linearOperator_);
341 params_[
"maxit"] = std::to_string(maxIter);
345 solver_ = constructPreconditionedSolver_(linearOperator_);
355 initializeParameters_(params);
359 solver_ = constructPreconditionedSolver_(linearOperator_);
364 void initializeParameters_(
const ParameterInitializer& params)
366 if (std::holds_alternative<std::string>(params))
369 params_ = std::get<Dune::ParameterTree>(params);
372 MatrixOperatorHolder makeSequentialLinearOperator_(std::shared_ptr<Matrix> A)
374 using SequentialTraits =
typename LinearSolverTraits::template Sequential<MatrixForSolver, XVectorForSolver>;
375 if constexpr (convertMultiTypeVectorAndMatrix)
379 return std::make_shared<typename SequentialTraits::LinearOperator>(M);
383 return std::make_shared<typename SequentialTraits::LinearOperator>(A);
387 template<
class ParallelTraits>
388 MatrixOperatorHolder makeParallelLinearOperator_(std::shared_ptr<Matrix> A, ParallelTraits = {})
392 prepareMatrixParallel<LinearSolverTraits, ParallelTraits>(*A, *parallelHelper_);
393 return std::make_shared<typename ParallelTraits::LinearOperator>(std::move(A), *communication_);
395 DUNE_THROW(Dune::InvalidStateException,
"Calling makeParallelLinearOperator for sequential run");
399 MatrixOperatorHolder makeParallelOrSequentialLinearOperator_(std::shared_ptr<Matrix> A)
401 return executeSequentialOrParallel_(
402 [&]{
return makeSequentialLinearOperator_(std::move(A)); },
403 [&](
auto traits){
return makeParallelLinearOperator_(std::move(A), traits); }
407 MatrixOperatorHolder makeSequentialLinearOperator_(Matrix& A)
408 {
return makeSequentialLinearOperator_(Dune::stackobject_to_shared_ptr<Matrix>(A)); }
410 MatrixOperatorHolder makeParallelOrSequentialLinearOperator_(Matrix& A)
411 {
return makeParallelOrSequentialLinearOperator_(Dune::stackobject_to_shared_ptr<Matrix>(A)); }
413 template<
class ParallelTraits>
414 MatrixOperatorHolder makeParallelLinearOperator_(Matrix& A, ParallelTraits = {})
415 {
return makeParallelLinearOperator_<ParallelTraits>(Dune::stackobject_to_shared_ptr<Matrix>(A)); }
417 IstlSolverResult solveSequential_(Matrix& A, XVector& x, BVector& b)
420 auto linearOperatorHolder = makeSequentialLinearOperator_(A);
421 auto solver = constructPreconditionedSolver_(linearOperatorHolder);
423 return solveSequential_(x, b, *solver);
426 IstlSolverResult solveSequential_(XVector& x, BVector& b, InverseOperator& solver)
const
428 Dune::InverseOperatorResult result;
429 if constexpr (convertMultiTypeVectorAndMatrix)
435 XVectorForSolver y(bTmp.size());
438 solver.apply(y, bTmp, result);
441 if (result.converged)
447 solver.apply(x, b, result);
453 IstlSolverResult solveSequentialOrParallel_(Matrix& A, XVector& x, BVector& b)
455 return executeSequentialOrParallel_(
456 [&]{
return solveSequential_(A, x, b); },
457 [&](
auto traits){
return solveParallel_(A, x, b, traits); }
461 IstlSolverResult solveSequentialOrParallel_(XVector& x, BVector& b, InverseOperator& solver)
const
463 return executeSequentialOrParallel_(
464 [&]{
return solveSequential_(x, b, solver); },
465 [&](
auto traits){
return solveParallel_(x, b, solver, traits); }
469 template<
class ParallelTraits>
470 IstlSolverResult solveParallel_(Matrix& A, XVector& x, BVector& b, ParallelTraits = {})
473 auto linearOperatorHolder = makeParallelLinearOperator_<ParallelTraits>(A);
474 auto solver = constructPreconditionedSolver_(linearOperatorHolder);
475 return solveParallel_<ParallelTraits>(x, b, *solver);
478 template<
class ParallelTraits>
479 IstlSolverResult solveParallel_(XVector& x, BVector& b, InverseOperator& solver, ParallelTraits = {})
const
483 prepareVectorParallel<LinearSolverTraits, ParallelTraits>(b, *parallelHelper_);
486 Dune::InverseOperatorResult result;
487 solver.apply(x, b, result);
490 DUNE_THROW(Dune::InvalidStateException,
"Calling makeParallelLinearOperator for sequential run");
495 std::shared_ptr<InverseOperator> constructPreconditionedSolver_(MatrixOperatorHolder& ops)
497 return std::visit([&](
auto&& op)
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);
506 if (prec->category() != op->category() && prec->category() == Dune::SolverCategory::sequential)
507 prec = Dune::wrapPreconditioner4Parallel(prec, op);
509 return std::make_shared<InverseOperator>(op, scalarProduct_, prec, params_);
513 template<
class Seq,
class Par>
514 decltype(
auto) executeSequentialOrParallel_(Seq&& sequentialAction, Par&& parallelAction)
const
520 return sequentialAction();
523 switch (solverCategory_)
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");
537 return sequentialAction();
542 std::shared_ptr<const ParallelHelper> parallelHelper_;
543 std::shared_ptr<Comm> communication_;
546 Dune::SolverCategory::Category solverCategory_;
547 std::shared_ptr<ScalarProduct> scalarProduct_;
550 MatrixOperatorHolder linearOperator_;
552 std::shared_ptr<InverseOperator> solver_;
554 Dune::ParameterTree params_;
578template<
class LSTraits,
class LATraits>
581 Dune::BiCGSTABSolver<typename LATraits::SingleTypeVector>,
603template<
class LSTraits,
class LATraits>
606 Dune::RestartedGMResSolver<typename LATraits::SingleTypeVector>,
629template<
class LSTraits,
class LATraits>
632 Dune::BiCGSTABSolver<typename LATraits::Vector>,
652template<
class LSTraits,
class LATraits>
655 Dune::CGSolver<typename LATraits::Vector>,
672template<
class LSTraits,
class LATraits>
675 Dune::BiCGSTABSolver<typename LATraits::SingleTypeVector>,
693template<
class LSTraits,
class LATraits>
696 Dune::CGSolver<typename LATraits::SingleTypeVector>,
716template<
class LSTraits,
class LATraits>
719 Dune::BiCGSTABSolver<typename LATraits::Vector>,
731template<
class LSTraits,
class LATraits,
template<
class M>
class Solver,
732 bool convertMultiTypeVectorAndMatrix = isMultiTypeBlockVector<typename LATraits::Vector>::value>
735 using Matrix =
typename LATraits::Matrix;
736 using XVector =
typename LATraits::Vector;
737 using BVector =
typename LATraits::Vector;
742 using InverseOperator = Dune::InverseOperator<XVectorForSolver, BVectorForSolver>;
751 return solve_(A, x, b);
760 DUNE_THROW(Dune::InvalidStateException,
"Called solve(x, b) but no linear operator has been set");
762 return solve_(x, b, *solver_);
770 if constexpr (convertMultiTypeVectorAndMatrix)
775 solver_ = std::make_shared<Solver<MatrixForSolver>>(*matrix_);
783 {
setMatrix(Dune::stackobject_to_shared_ptr(A)); }
790 return "Direct solver";
797 if constexpr (convertMultiTypeVectorAndMatrix)
800 Solver<MatrixForSolver> solver(AA, this->
verbosity() > 0);
801 return solve_(x, b, solver);
805 Solver<MatrixForSolver> solver(A, this->
verbosity() > 0);
806 return solve_(x, b, solver);
810 IstlSolverResult solve_(XVector& x,
const BVector& b, InverseOperator& solver)
const
812 Dune::InverseOperatorResult result;
814 if constexpr (convertMultiTypeVectorAndMatrix)
817 XVectorForSolver xx(bb.size());
818 solver.apply(xx, bb, result);
819 checkResult_(xx, result);
820 if (result.converged)
826 BVectorForSolver bTmp(b);
827 solver.apply(x, bTmp, result);
828 checkResult_(x, result);
834 void checkResult_(XVectorForSolver& x, Dune::InverseOperatorResult& result)
const
836 flatVectorForEach(x, [&](
auto&& entry, std::size_t){
837 using std::isnan, std::isinf;
838 if (isnan(entry) || isinf(entry))
839 result.converged =
false;
844 std::shared_ptr<MatrixForSolver> matrix_;
846 std::shared_ptr<InverseOperator> solver_;
852#include <dune/istl/superlu.hh>
864template<
class LSTraits,
class LATraits>
865using SuperLUIstlSolver = Detail::DirectIstlSolver<LSTraits, LATraits, Dune::SuperLU>;
872#include <dune/istl/umfpack.hh>
873#include <dune/common/version.hh>
885template<
class LSTraits,
class LATraits>
886using UMFPackIstlSolver = Detail::DirectIstlSolver<
887 LSTraits, LATraits, Dune::UMFPack
888#if DUNE_VERSION_GTE(DUNE_ISTL,2,10)
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 ¶ms="")
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 ¶ms="")
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 ¶ms)
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 ¶ms="")
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
Definition: istlsolvers.hh:58
auto operator()(TL typeList, const M &matrix, const Dune::ParameterTree &config)
Definition: istlsolvers.hh:61
Definition: istlsolvers.hh:74
Definition: parallelhelpers.hh:29
Base class for linear solvers.
Definition: solver.hh:27
LinearSolver(const std::string ¶mGroup="")
Construct the solver.
Definition: solver.hh:43
int verbosity() const
the verbosity level
Definition: solver.hh:82
static Dune::ParameterTree createParameterTree(const std::string ¶mGroup="")
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
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.
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()=default
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.