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>
54template<
template<
class,
class,
class,
int>
class Preconditioner,
int blockLevel = 1>
58 template<
class TL,
class M>
59 auto operator() (TL typeList,
const M& matrix,
const Dune::ParameterTree& config)
61 using Matrix =
typename Dune::TypeListElement<0,
decltype(typeList)>::type;
62 using Domain =
typename Dune::TypeListElement<1,
decltype(typeList)>::type;
63 using Range =
typename Dune::TypeListElement<2,
decltype(typeList)>::type;
64 std::shared_ptr<Dune::Preconditioner<Domain, Range>> preconditioner
65 = std::make_shared<Preconditioner<Matrix, Domain, Range, blockLevel>>(matrix, config);
66 return preconditioner;
70template<
template<
class,
class,
class>
class Preconditioner>
73 template<
class TL,
class M>
74 auto operator() (TL typeList,
const M& matrix,
const Dune::ParameterTree& config)
76 using Matrix =
typename Dune::TypeListElement<0,
decltype(typeList)>::type;
77 using Domain =
typename Dune::TypeListElement<1,
decltype(typeList)>::type;
78 using Range =
typename Dune::TypeListElement<2,
decltype(typeList)>::type;
79 std::shared_ptr<Dune::Preconditioner<Domain, Range>> preconditioner
80 = std::make_shared<Preconditioner<Matrix, Domain, Range>>(matrix, config);
81 return preconditioner;
87template<
class M,
bool convert = false>
94template<
class V,
bool convert = false>
101template<
class LSTraits,
class LATraits,
bool convert,
bool parallel = LSTraits::canCommunicate>
104template<
class LSTraits,
class LATraits,
bool convert>
111 std::shared_ptr<typename LSTraits::template Sequential<M, V>::LinearOperator>,
112 std::shared_ptr<typename LSTraits::template ParallelOverlapping<M, V>::LinearOperator>,
113 std::shared_ptr<typename LSTraits::template ParallelNonoverlapping<M, V>::LinearOperator>
116 using type = std::variant<
117 std::shared_ptr<typename LSTraits::template Sequential<M, V>::LinearOperator>
122template<
class LSTraits,
class LATraits,
bool convert>
128 std::shared_ptr<typename LSTraits::template Sequential<M, V>::LinearOperator>
145 operator bool()
const {
return this->converged; }
153 class InverseOperator,
class PreconditionerFactory,
154 bool convertMultiTypeLATypes =
false>
160 using Scalar =
typename InverseOperator::real_type;
162 using ScalarProduct = Dune::ScalarProduct<typename InverseOperator::domain_type>;
164 static constexpr bool convertMultiTypeVectorAndMatrix
175 using Comm = Dune::OwnerOverlapCopyCommunication<Dune::bigunsignedint<96>,
int>;
179 using ParameterInitializer = std::variant<std::string, Dune::ParameterTree>;
187 if (Dune::MPIHelper::getCommunication().size() > 1)
188 DUNE_THROW(Dune::InvalidStateException,
"Using sequential constructor for parallel run. Use signature with gridView and dofMapper!");
190 initializeParameters_(params);
191 solverCategory_ = Dune::SolverCategory::sequential;
192 scalarProduct_ = std::make_shared<ScalarProduct>();
198 template <
class Gr
idView,
class DofMapper>
200 const DofMapper& dofMapper,
201 const ParameterInitializer& params =
"")
203 initializeParameters_(params);
205 solverCategory_ = Detail::solverCategory<LinearSolverTraits>(gridView);
209 if (solverCategory_ != Dune::SolverCategory::sequential)
211 parallelHelper_ = std::make_shared<ParallelISTLHelper<LinearSolverTraits>>(gridView, dofMapper);
212 communication_ = std::make_shared<Comm>(gridView.comm(), solverCategory_);
213 scalarProduct_ = Dune::createScalarProduct<XVector>(*communication_, solverCategory_);
214 parallelHelper_->createParallelIndexSet(*communication_);
217 scalarProduct_ = std::make_shared<ScalarProduct>();
220 scalarProduct_ = std::make_shared<ScalarProduct>();
222 solverCategory_ = Dune::SolverCategory::sequential;
223 scalarProduct_ = std::make_shared<ScalarProduct>();
231 template <
class Gr
idView,
class DofMapper>
233 std::shared_ptr<ScalarProduct> scalarProduct,
234 const GridView& gridView,
235 const DofMapper& dofMapper,
236 const ParameterInitializer& params =
"")
238 initializeParameters_(params);
240 scalarProduct_ = scalarProduct;
241 communication_ = communication;
244 if (solverCategory_ != Dune::SolverCategory::sequential)
246 parallelHelper_ = std::make_shared<ParallelISTLHelper<LinearSolverTraits>>(gridView, dofMapper);
247 parallelHelper_->createParallelIndexSet(communication);
257 {
return solveSequentialOrParallel_(A, x, b); }
264 linearOperator_ = makeParallelOrSequentialLinearOperator_(std::move(A));
265 solver_ = constructPreconditionedSolver_(linearOperator_);
273 {
setMatrix(Dune::stackobject_to_shared_ptr(A)); }
281 DUNE_THROW(Dune::InvalidStateException,
"Called solve(x, b) but no linear operator has been set");
283 return solveSequentialOrParallel_(x, b, *solver_);
289 Scalar
norm(
const XVector& x)
const
294 if (solverCategory_ == Dune::SolverCategory::nonoverlapping)
297 using GV =
typename LinearSolverTraits::GridView;
298 using DM =
typename LinearSolverTraits::DofMapper;
301 return scalarProduct_->norm(y);
305 if constexpr (convertMultiTypeVectorAndMatrix)
308 return scalarProduct_->norm(y);
311 return scalarProduct_->norm(x);
317 const std::string&
name()
const
327 params_[
"reduction"] = std::to_string(residReduction);
331 solver_ = constructPreconditionedSolver_(linearOperator_);
339 params_[
"maxit"] = std::to_string(maxIter);
343 solver_ = constructPreconditionedSolver_(linearOperator_);
353 initializeParameters_(params);
357 solver_ = constructPreconditionedSolver_(linearOperator_);
362 void initializeParameters_(
const ParameterInitializer& params)
364 if (std::holds_alternative<std::string>(params))
367 params_ = std::get<Dune::ParameterTree>(params);
370 MatrixOperatorHolder makeSequentialLinearOperator_(std::shared_ptr<Matrix> A)
372 using SequentialTraits =
typename LinearSolverTraits::template Sequential<MatrixForSolver, XVectorForSolver>;
373 if constexpr (convertMultiTypeVectorAndMatrix)
377 return std::make_shared<typename SequentialTraits::LinearOperator>(M);
381 return std::make_shared<typename SequentialTraits::LinearOperator>(A);
385 template<
class ParallelTraits>
386 MatrixOperatorHolder makeParallelLinearOperator_(std::shared_ptr<Matrix> A, ParallelTraits = {})
390 prepareMatrixParallel<LinearSolverTraits, ParallelTraits>(*A, *parallelHelper_);
391 return std::make_shared<typename ParallelTraits::LinearOperator>(std::move(A), *communication_);
393 DUNE_THROW(Dune::InvalidStateException,
"Calling makeParallelLinearOperator for sequential run");
397 MatrixOperatorHolder makeParallelOrSequentialLinearOperator_(std::shared_ptr<Matrix> A)
399 return executeSequentialOrParallel_(
400 [&]{
return makeSequentialLinearOperator_(std::move(A)); },
401 [&](
auto traits){
return makeParallelLinearOperator_(std::move(A), traits); }
405 MatrixOperatorHolder makeSequentialLinearOperator_(Matrix& A)
406 {
return makeSequentialLinearOperator_(Dune::stackobject_to_shared_ptr<Matrix>(A)); }
408 MatrixOperatorHolder makeParallelOrSequentialLinearOperator_(Matrix& A)
409 {
return makeParallelOrSequentialLinearOperator_(Dune::stackobject_to_shared_ptr<Matrix>(A)); }
411 template<
class ParallelTraits>
412 MatrixOperatorHolder makeParallelLinearOperator_(Matrix& A, ParallelTraits = {})
413 {
return makeParallelLinearOperator_<ParallelTraits>(Dune::stackobject_to_shared_ptr<Matrix>(A)); }
415 IstlSolverResult solveSequential_(Matrix& A, XVector& x, BVector& b)
418 auto linearOperatorHolder = makeSequentialLinearOperator_(A);
419 auto solver = constructPreconditionedSolver_(linearOperatorHolder);
421 return solveSequential_(x, b, *solver);
424 IstlSolverResult solveSequential_(XVector& x, BVector& b, InverseOperator& solver)
const
426 Dune::InverseOperatorResult result;
427 if constexpr (convertMultiTypeVectorAndMatrix)
433 XVectorForSolver y(bTmp.size());
436 solver.apply(y, bTmp, result);
439 if (result.converged)
445 solver.apply(x, b, result);
451 IstlSolverResult solveSequentialOrParallel_(Matrix& A, XVector& x, BVector& b)
453 return executeSequentialOrParallel_(
454 [&]{
return solveSequential_(A, x, b); },
455 [&](
auto traits){
return solveParallel_(A, x, b, traits); }
459 IstlSolverResult solveSequentialOrParallel_(XVector& x, BVector& b, InverseOperator& solver)
const
461 return executeSequentialOrParallel_(
462 [&]{
return solveSequential_(x, b, solver); },
463 [&](
auto traits){
return solveParallel_(x, b, solver, traits); }
467 template<
class ParallelTraits>
468 IstlSolverResult solveParallel_(Matrix& A, XVector& x, BVector& b, ParallelTraits = {})
471 auto linearOperatorHolder = makeParallelLinearOperator_<ParallelTraits>(A);
472 auto solver = constructPreconditionedSolver_(linearOperatorHolder);
473 return solveParallel_<ParallelTraits>(x, b, *solver);
476 template<
class ParallelTraits>
477 IstlSolverResult solveParallel_(XVector& x, BVector& b, InverseOperator& solver, ParallelTraits = {})
const
481 prepareVectorParallel<LinearSolverTraits, ParallelTraits>(b, *parallelHelper_);
484 Dune::InverseOperatorResult result;
485 solver.apply(x, b, result);
488 DUNE_THROW(Dune::InvalidStateException,
"Calling makeParallelLinearOperator for sequential run");
493 std::shared_ptr<InverseOperator> constructPreconditionedSolver_(MatrixOperatorHolder& ops)
495 return std::visit([&](
auto&& op)
497 using LinearOperator =
typename std::decay_t<
decltype(op)>::element_type;
498 const auto& params = params_.sub(
"preconditioner");
499 using Prec = Dune::Preconditioner<typename LinearOperator::domain_type, typename LinearOperator::range_type>;
500 using TL = Dune::TypeList<typename LinearOperator::matrix_type, typename LinearOperator::domain_type, typename LinearOperator::range_type>;
501 std::shared_ptr<Prec> prec = PreconditionerFactory{}(TL{}, op, params);
504 if (prec->category() != op->category() && prec->category() == Dune::SolverCategory::sequential)
505 prec = Dune::wrapPreconditioner4Parallel(prec, op);
507 return std::make_shared<InverseOperator>(op, scalarProduct_, prec, params_);
511 template<
class Seq,
class Par>
512 decltype(
auto) executeSequentialOrParallel_(Seq&& sequentialAction, Par&& parallelAction)
const
518 return sequentialAction();
521 switch (solverCategory_)
523 case Dune::SolverCategory::sequential:
524 return sequentialAction();
525 case Dune::SolverCategory::nonoverlapping:
526 using NOTraits =
typename LinearSolverTraits::template ParallelNonoverlapping<Matrix, XVector>;
527 return parallelAction(NOTraits{});
528 case Dune::SolverCategory::overlapping:
529 using OTraits =
typename LinearSolverTraits::template ParallelOverlapping<Matrix, XVector>;
530 return parallelAction(OTraits{});
531 default: DUNE_THROW(Dune::InvalidStateException,
"Unknown solver category");
535 return sequentialAction();
540 std::shared_ptr<const ParallelHelper> parallelHelper_;
541 std::shared_ptr<Comm> communication_;
544 Dune::SolverCategory::Category solverCategory_;
545 std::shared_ptr<ScalarProduct> scalarProduct_;
548 MatrixOperatorHolder linearOperator_;
550 std::shared_ptr<InverseOperator> solver_;
552 Dune::ParameterTree params_;
576template<
class LSTraits,
class LATraits>
579 Dune::BiCGSTABSolver<typename LATraits::SingleTypeVector>,
601template<
class LSTraits,
class LATraits>
604 Dune::RestartedGMResSolver<typename LATraits::SingleTypeVector>,
627template<
class LSTraits,
class LATraits>
630 Dune::BiCGSTABSolver<typename LATraits::Vector>,
650template<
class LSTraits,
class LATraits>
653 Dune::CGSolver<typename LATraits::Vector>,
670template<
class LSTraits,
class LATraits>
673 Dune::BiCGSTABSolver<typename LATraits::SingleTypeVector>,
691template<
class LSTraits,
class LATraits>
694 Dune::CGSolver<typename LATraits::SingleTypeVector>,
714template<
class LSTraits,
class LATraits>
717 Dune::BiCGSTABSolver<typename LATraits::Vector>,
729template<
class LSTraits,
class LATraits,
template<
class M>
class Solver>
732 using Matrix =
typename LATraits::Matrix;
733 using XVector =
typename LATraits::Vector;
734 using BVector =
typename LATraits::Vector;
740 using InverseOperator = Dune::InverseOperator<XVectorForSolver, BVectorForSolver>;
749 return solve_(A, x, b);
758 DUNE_THROW(Dune::InvalidStateException,
"Called solve(x, b) but no linear operator has been set");
760 return solve_(x, b, *solver_);
768 if constexpr (convertMultiTypeVectorAndMatrix)
773 solver_ = std::make_shared<Solver<MatrixForSolver>>(*matrix_);
781 {
setMatrix(Dune::stackobject_to_shared_ptr(A)); }
788 return "Direct solver";
798 Solver<MatrixForSolver> solver(AA, this->
verbosity() > 0);
799 return solve_(x, b, solver);
803 Solver<MatrixForSolver> solver(A, this->
verbosity() > 0);
804 return solve_(x, b, solver);
808 IstlSolverResult solve_(XVector& x,
const BVector& b, InverseOperator& solver)
const
810 static_assert(isBCRSMatrix<MatrixForSolver>::value,
"Direct solver only works with BCRS matrices!");
811 static_assert(MatrixForSolver::block_type::rows == MatrixForSolver::block_type::cols,
"Matrix block must be quadratic!");
813 Dune::InverseOperatorResult result;
814 if constexpr (isMultiTypeBlockVector<BVector>())
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);
833 void checkResult_(
const XVectorForSolver& x, Dune::InverseOperatorResult& result)
const
836 for (
int i = 0; i < size; i++)
838 for (
int j = 0; j < x[i].size(); j++)
842 if (isnan(x[i][j]) || isinf(x[i][j]))
844 result.converged =
false;
852 std::shared_ptr<MatrixForSolver> matrix_;
854 std::shared_ptr<InverseOperator> solver_;
860#include <dune/istl/superlu.hh>
872template<
class LSTraits,
class LATraits>
873using SuperLUIstlSolver = Detail::DirectIstlSolver<LSTraits, LATraits, Dune::SuperLU>;
880#include <dune/istl/umfpack.hh>
892template<
class LSTraits,
class LATraits>
893using UMFPackIstlSolver = Detail::DirectIstlSolver<LSTraits, LATraits, Dune::UMFPack>;
Direct dune-istl linear solvers.
Definition: istlsolvers.hh:731
IstlSolverResult solve(const Matrix &A, XVector &x, const BVector &b)
Solve the linear system Ax = b.
Definition: istlsolvers.hh:747
void setMatrix(std::shared_ptr< Matrix > A)
Set the matrix A of the linear system Ax = b for reuse.
Definition: istlsolvers.hh:766
IstlSolverResult solve(XVector &x, const BVector &b)
Solve the linear system Ax = b using the matrix set with setMatrix.
Definition: istlsolvers.hh:755
void setMatrix(Matrix &A)
Set the matrix A of the linear system Ax = b for reuse.
Definition: istlsolvers.hh:780
std::string name() const
name of the linear solver
Definition: istlsolvers.hh:786
Standard dune-istl iterative linear solvers.
Definition: istlsolvers.hh:156
IstlIterativeLinearSolver(const ParameterInitializer ¶ms="")
Constructor for sequential solvers.
Definition: istlsolvers.hh:185
void setMatrix(Matrix &A)
Set the matrix A of the linear system Ax = b for reuse.
Definition: istlsolvers.hh:272
void setResidualReduction(double residReduction)
Set the residual reduction tolerance.
Definition: istlsolvers.hh:325
void setMaxIter(std::size_t maxIter)
Set the maximum number of linear solver iterations.
Definition: istlsolvers.hh:337
IstlIterativeLinearSolver(const GridView &gridView, const DofMapper &dofMapper, const ParameterInitializer ¶ms="")
Constructor for parallel and sequential solvers.
Definition: istlsolvers.hh:199
const std::string & name() const
The name of the linear solver.
Definition: istlsolvers.hh:317
void setParams(const ParameterInitializer ¶ms)
Set the linear solver parameters.
Definition: istlsolvers.hh:351
IstlSolverResult solve(Matrix &A, XVector &x, BVector &b)
Solve the linear system Ax = b.
Definition: istlsolvers.hh:256
void setMatrix(std::shared_ptr< Matrix > A)
Set the matrix A of the linear system Ax = b for reuse.
Definition: istlsolvers.hh:262
IstlSolverResult solve(XVector &x, BVector &b) const
Solve the linear system Ax = b where A has been set with setMatrix.
Definition: istlsolvers.hh:278
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:232
Scalar norm(const XVector &x) const
Compute the 2-norm of vector x.
Definition: istlsolvers.hh:289
Definition: istlsolvers.hh:56
auto operator()(TL typeList, const M &matrix, const Dune::ParameterTree &config)
Definition: istlsolvers.hh:59
Definition: istlsolvers.hh:72
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:49
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:40
Dune::AMGCreator IstlAmgPreconditionerFactory
Definition: istlsolvers.hh:85
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:15
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:137
IstlSolverResult(IstlSolverResult &&)=default
IstlSolverResult(const Dune::InverseOperatorResult &o)
Definition: istlsolvers.hh:142
IstlSolverResult(Dune::InverseOperatorResult &&o)
Definition: istlsolvers.hh:143
IstlSolverResult()=default
IstlSolverResult(const IstlSolverResult &)=default
std::decay_t< decltype(MatrixConverter< M >::multiTypeToBCRSMatrix(std::declval< M >()))> type
Definition: istlsolvers.hh:92
Definition: istlsolvers.hh:88
M type
Definition: istlsolvers.hh:88
typename VectorForSolver< typename LATraits::Vector, convert >::type V
Definition: istlsolvers.hh:126
std::variant< std::shared_ptr< typename LSTraits::template Sequential< M, V >::LinearOperator > > type
Definition: istlsolvers.hh:129
typename MatrixForSolver< typename LATraits::Matrix, convert >::type M
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:114
typename VectorForSolver< typename LATraits::Vector, convert >::type V
Definition: istlsolvers.hh:108
typename MatrixForSolver< typename LATraits::Matrix, convert >::type M
Definition: istlsolvers.hh:107
Definition: istlsolvers.hh:102
std::decay_t< decltype(VectorConverter< V >::multiTypeToBlockVector(std::declval< V >()))> type
Definition: istlsolvers.hh:99
Definition: istlsolvers.hh:95
V type
Definition: istlsolvers.hh:95
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.