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/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>
41#include <dune/istl/foreach.hh>
57template<
template<
class,
class,
class,
int>
class Preconditioner,
int blockLevel = 1>
61 template<
class TL,
class M>
62 auto operator() (TL typeList,
const M& matrix,
const Dune::ParameterTree& config)
64 using Matrix =
typename Dune::TypeListElement<0,
decltype(typeList)>::type;
65 using Domain =
typename Dune::TypeListElement<1,
decltype(typeList)>::type;
66 using Range =
typename Dune::TypeListElement<2,
decltype(typeList)>::type;
67 std::shared_ptr<Dune::Preconditioner<Domain, Range>> preconditioner
68 = std::make_shared<Preconditioner<Matrix, Domain, Range, blockLevel>>(matrix, config);
69 return preconditioner;
73template<
template<
class,
class,
class>
class Preconditioner>
76 template<
class TL,
class M>
77 auto operator() (TL typeList,
const M& matrix,
const Dune::ParameterTree& config)
79 using Matrix =
typename Dune::TypeListElement<0,
decltype(typeList)>::type;
80 using Domain =
typename Dune::TypeListElement<1,
decltype(typeList)>::type;
81 using Range =
typename Dune::TypeListElement<2,
decltype(typeList)>::type;
82 std::shared_ptr<Dune::Preconditioner<Domain, Range>> preconditioner
83 = std::make_shared<Preconditioner<Matrix, Domain, Range>>(matrix, config);
84 return preconditioner;
90template<
class M,
bool convert = false>
97template<
class V,
bool convert = false>
104template<
class LSTraits,
class LATraits,
bool convert,
bool parallel = LSTraits::canCommunicate>
107template<
class LSTraits,
class LATraits,
bool convert>
114 std::shared_ptr<typename LSTraits::template Sequential<M, V>::LinearOperator>,
115 std::shared_ptr<typename LSTraits::template ParallelOverlapping<M, V>::LinearOperator>,
116 std::shared_ptr<typename LSTraits::template ParallelNonoverlapping<M, V>::LinearOperator>
119 using type = std::variant<
120 std::shared_ptr<typename LSTraits::template Sequential<M, V>::LinearOperator>
125template<
class LSTraits,
class LATraits,
bool convert>
131 std::shared_ptr<typename LSTraits::template Sequential<M, V>::LinearOperator>
148 operator bool()
const {
return this->converged; }
156 class InverseOperator,
class PreconditionerFactory,
157 bool convertMultiTypeLATypes =
false>
163 using Scalar =
typename InverseOperator::real_type;
165 using ScalarProduct = Dune::ScalarProduct<typename InverseOperator::domain_type>;
167 static constexpr bool convertMultiTypeVectorAndMatrix
178 using Comm = Dune::OwnerOverlapCopyCommunication<Dune::bigunsignedint<96>,
int>;
182 using ParameterInitializer = std::variant<std::string, Dune::ParameterTree>;
190 if (Dune::MPIHelper::getCommunication().size() > 1)
191 DUNE_THROW(Dune::InvalidStateException,
"Using sequential constructor for parallel run. Use signature with gridView and dofMapper!");
193 initializeParameters_(params);
194 solverCategory_ = Dune::SolverCategory::sequential;
195 scalarProduct_ = std::make_shared<ScalarProduct>();
201 template <
class Gr
idView,
class DofMapper>
203 const DofMapper& dofMapper,
204 const ParameterInitializer& params =
"")
206 initializeParameters_(params);
208 solverCategory_ = Detail::solverCategory<LinearSolverTraits>(gridView);
212 if (solverCategory_ != Dune::SolverCategory::sequential)
214 parallelHelper_ = std::make_shared<ParallelISTLHelper<LinearSolverTraits>>(gridView, dofMapper);
215 communication_ = std::make_shared<Comm>(gridView.comm(), solverCategory_);
216 scalarProduct_ = Dune::createScalarProduct<XVector>(*communication_, solverCategory_);
217 parallelHelper_->createParallelIndexSet(*communication_);
220 scalarProduct_ = std::make_shared<ScalarProduct>();
223 scalarProduct_ = std::make_shared<ScalarProduct>();
225 solverCategory_ = Dune::SolverCategory::sequential;
226 scalarProduct_ = std::make_shared<ScalarProduct>();
234 template <
class Gr
idView,
class DofMapper>
236 std::shared_ptr<ScalarProduct> scalarProduct,
237 const GridView& gridView,
238 const DofMapper& dofMapper,
239 const ParameterInitializer& params =
"")
241 initializeParameters_(params);
243 scalarProduct_ = scalarProduct;
244 communication_ = communication;
247 if (solverCategory_ != Dune::SolverCategory::sequential)
249 parallelHelper_ = std::make_shared<ParallelISTLHelper<LinearSolverTraits>>(gridView, dofMapper);
250 parallelHelper_->createParallelIndexSet(communication);
260 {
return solveSequentialOrParallel_(A, x, b); }
267 linearOperator_ = makeParallelOrSequentialLinearOperator_(std::move(A));
268 solver_ = constructPreconditionedSolver_(linearOperator_);
276 {
setMatrix(Dune::stackobject_to_shared_ptr(A)); }
284 DUNE_THROW(Dune::InvalidStateException,
"Called solve(x, b) but no linear operator has been set");
286 return solveSequentialOrParallel_(x, b, *solver_);
292 Scalar
norm(
const XVector& x)
const
297 if (solverCategory_ == Dune::SolverCategory::nonoverlapping)
300 using GV =
typename LinearSolverTraits::GridView;
301 using DM =
typename LinearSolverTraits::DofMapper;
304 return scalarProduct_->norm(y);
308 if constexpr (convertMultiTypeVectorAndMatrix)
311 return scalarProduct_->norm(y);
314 return scalarProduct_->norm(x);
320 const std::string&
name()
const
330 params_[
"reduction"] = std::to_string(residReduction);
334 solver_ = constructPreconditionedSolver_(linearOperator_);
342 params_[
"maxit"] = std::to_string(maxIter);
346 solver_ = constructPreconditionedSolver_(linearOperator_);
356 initializeParameters_(params);
360 solver_ = constructPreconditionedSolver_(linearOperator_);
365 void initializeParameters_(
const ParameterInitializer& params)
367 if (std::holds_alternative<std::string>(params))
370 params_ = std::get<Dune::ParameterTree>(params);
373 MatrixOperatorHolder makeSequentialLinearOperator_(std::shared_ptr<Matrix> A)
375 using SequentialTraits =
typename LinearSolverTraits::template Sequential<MatrixForSolver, XVectorForSolver>;
376 if constexpr (convertMultiTypeVectorAndMatrix)
380 return std::make_shared<typename SequentialTraits::LinearOperator>(M);
384 return std::make_shared<typename SequentialTraits::LinearOperator>(A);
388 template<
class ParallelTraits>
389 MatrixOperatorHolder makeParallelLinearOperator_(std::shared_ptr<Matrix> A, ParallelTraits = {})
393 prepareMatrixParallel<LinearSolverTraits, ParallelTraits>(*A, *parallelHelper_);
394 return std::make_shared<typename ParallelTraits::LinearOperator>(std::move(A), *communication_);
396 DUNE_THROW(Dune::InvalidStateException,
"Calling makeParallelLinearOperator for sequential run");
400 MatrixOperatorHolder makeParallelOrSequentialLinearOperator_(std::shared_ptr<Matrix> A)
402 return executeSequentialOrParallel_(
403 [&]{
return makeSequentialLinearOperator_(std::move(A)); },
404 [&](
auto traits){
return makeParallelLinearOperator_(std::move(A), traits); }
408 MatrixOperatorHolder makeSequentialLinearOperator_(Matrix& A)
409 {
return makeSequentialLinearOperator_(Dune::stackobject_to_shared_ptr<Matrix>(A)); }
411 MatrixOperatorHolder makeParallelOrSequentialLinearOperator_(Matrix& A)
412 {
return makeParallelOrSequentialLinearOperator_(Dune::stackobject_to_shared_ptr<Matrix>(A)); }
414 template<
class ParallelTraits>
415 MatrixOperatorHolder makeParallelLinearOperator_(Matrix& A, ParallelTraits = {})
416 {
return makeParallelLinearOperator_<ParallelTraits>(Dune::stackobject_to_shared_ptr<Matrix>(A)); }
418 IstlSolverResult solveSequential_(Matrix& A, XVector& x, BVector& b)
421 auto linearOperatorHolder = makeSequentialLinearOperator_(A);
422 auto solver = constructPreconditionedSolver_(linearOperatorHolder);
424 return solveSequential_(x, b, *solver);
427 IstlSolverResult solveSequential_(XVector& x, BVector& b, InverseOperator& solver)
const
429 Dune::InverseOperatorResult result;
430 if constexpr (convertMultiTypeVectorAndMatrix)
436 XVectorForSolver y(bTmp.size());
439 solver.apply(y, bTmp, result);
442 if (result.converged)
448 solver.apply(x, b, result);
454 IstlSolverResult solveSequentialOrParallel_(Matrix& A, XVector& x, BVector& b)
456 return executeSequentialOrParallel_(
457 [&]{
return solveSequential_(A, x, b); },
458 [&](
auto traits){
return solveParallel_(A, x, b, traits); }
462 IstlSolverResult solveSequentialOrParallel_(XVector& x, BVector& b, InverseOperator& solver)
const
464 return executeSequentialOrParallel_(
465 [&]{
return solveSequential_(x, b, solver); },
466 [&](
auto traits){
return solveParallel_(x, b, solver, traits); }
470 template<
class ParallelTraits>
471 IstlSolverResult solveParallel_(Matrix& A, XVector& x, BVector& b, ParallelTraits = {})
474 auto linearOperatorHolder = makeParallelLinearOperator_<ParallelTraits>(A);
475 auto solver = constructPreconditionedSolver_(linearOperatorHolder);
476 return solveParallel_<ParallelTraits>(x, b, *solver);
479 template<
class ParallelTraits>
480 IstlSolverResult solveParallel_(XVector& x, BVector& b, InverseOperator& solver, ParallelTraits = {})
const
484 prepareVectorParallel<LinearSolverTraits, ParallelTraits>(b, *parallelHelper_);
487 Dune::InverseOperatorResult result;
488 solver.apply(x, b, result);
491 DUNE_THROW(Dune::InvalidStateException,
"Calling makeParallelLinearOperator for sequential run");
496 std::shared_ptr<InverseOperator> constructPreconditionedSolver_(MatrixOperatorHolder& ops)
498 return std::visit([&](
auto&& op)
500 using LinearOperator =
typename std::decay_t<
decltype(op)>::element_type;
501 const auto& params = params_.sub(
"preconditioner");
502 using Prec = Dune::Preconditioner<typename LinearOperator::domain_type, typename LinearOperator::range_type>;
503#if DUNE_VERSION_GTE(DUNE_ISTL,2,11)
504 using OpTraits = Dune::OperatorTraits<LinearOperator>;
505 std::shared_ptr<Prec> prec = PreconditionerFactory{}(OpTraits{}, op, params);
507 using TL = Dune::TypeList<typename LinearOperator::matrix_type, typename LinearOperator::domain_type, typename LinearOperator::range_type>;
508 std::shared_ptr<Prec> prec = PreconditionerFactory{}(TL{}, op, params);
512#if DUNE_VERSION_LT(DUNE_ISTL,2,11)
513 if (prec->category() != op->category() && prec->category() == Dune::SolverCategory::sequential)
514 prec = Dune::wrapPreconditioner4Parallel(prec, op);
516 if constexpr (OpTraits::isParallel)
518 using Comm =
typename OpTraits::comm_type;
519 const Comm& comm = OpTraits::getCommOrThrow(op);
520 if (op->category() == Dune::SolverCategory::overlapping && prec->category() == Dune::SolverCategory::sequential)
521 prec = std::make_shared<Dune::BlockPreconditioner<typename OpTraits::domain_type, typename OpTraits::range_type,Comm> >(prec, comm);
522 else if (op->category() == Dune::SolverCategory::nonoverlapping && prec->category() == Dune::SolverCategory::sequential)
523 prec = std::make_shared<Dune::NonoverlappingBlockPreconditioner<Comm, Prec> >(prec, comm);
527 return std::make_shared<InverseOperator>(op, scalarProduct_, prec, params_);
531 template<
class Seq,
class Par>
532 decltype(
auto) executeSequentialOrParallel_(Seq&& sequentialAction, Par&& parallelAction)
const
538 return sequentialAction();
541 switch (solverCategory_)
543 case Dune::SolverCategory::sequential:
544 return sequentialAction();
545 case Dune::SolverCategory::nonoverlapping:
546 using NOTraits =
typename LinearSolverTraits::template ParallelNonoverlapping<Matrix, XVector>;
547 return parallelAction(NOTraits{});
548 case Dune::SolverCategory::overlapping:
549 using OTraits =
typename LinearSolverTraits::template ParallelOverlapping<Matrix, XVector>;
550 return parallelAction(OTraits{});
551 default: DUNE_THROW(Dune::InvalidStateException,
"Unknown solver category");
555 return sequentialAction();
560 std::shared_ptr<const ParallelHelper> parallelHelper_;
561 std::shared_ptr<Comm> communication_;
564 Dune::SolverCategory::Category solverCategory_;
565 std::shared_ptr<ScalarProduct> scalarProduct_;
568 MatrixOperatorHolder linearOperator_;
570 std::shared_ptr<InverseOperator> solver_;
572 Dune::ParameterTree params_;
596template<
class LSTraits,
class LATraits>
599 Dune::BiCGSTABSolver<typename LATraits::SingleTypeVector>,
621template<
class LSTraits,
class LATraits>
624 Dune::RestartedGMResSolver<typename LATraits::SingleTypeVector>,
647template<
class LSTraits,
class LATraits>
650 Dune::BiCGSTABSolver<typename LATraits::Vector>,
670template<
class LSTraits,
class LATraits>
673 Dune::CGSolver<typename LATraits::Vector>,
690template<
class LSTraits,
class LATraits>
693 Dune::BiCGSTABSolver<typename LATraits::SingleTypeVector>,
711template<
class LSTraits,
class LATraits>
714 Dune::CGSolver<typename LATraits::SingleTypeVector>,
734template<
class LSTraits,
class LATraits>
737 Dune::BiCGSTABSolver<typename LATraits::Vector>,
749template<
class LSTraits,
class LATraits,
template<
class M>
class Solver,
750 bool convertMultiTypeVectorAndMatrix = isMultiTypeBlockVector<typename LATraits::Vector>::value>
753 using Matrix =
typename LATraits::Matrix;
754 using XVector =
typename LATraits::Vector;
755 using BVector =
typename LATraits::Vector;
760 using InverseOperator = Dune::InverseOperator<XVectorForSolver, BVectorForSolver>;
769 return solve_(A, x, b);
778 DUNE_THROW(Dune::InvalidStateException,
"Called solve(x, b) but no linear operator has been set");
780 return solve_(x, b, *solver_);
788 if constexpr (convertMultiTypeVectorAndMatrix)
793 solver_ = std::make_shared<Solver<MatrixForSolver>>(*matrix_);
801 {
setMatrix(Dune::stackobject_to_shared_ptr(A)); }
808 return "Direct solver";
815 if constexpr (convertMultiTypeVectorAndMatrix)
818 Solver<MatrixForSolver> solver(AA, this->
verbosity() > 0);
819 return solve_(x, b, solver);
823 Solver<MatrixForSolver> solver(A, this->
verbosity() > 0);
824 return solve_(x, b, solver);
828 IstlSolverResult solve_(XVector& x,
const BVector& b, InverseOperator& solver)
const
830 Dune::InverseOperatorResult result;
832 if constexpr (convertMultiTypeVectorAndMatrix)
835 XVectorForSolver xx(bb.size());
836 solver.apply(xx, bb, result);
837 checkResult_(xx, result);
838 if (result.converged)
844 BVectorForSolver bTmp(b);
845 solver.apply(x, bTmp, result);
846 checkResult_(x, result);
852 void checkResult_(XVectorForSolver& x, Dune::InverseOperatorResult& result)
const
854 flatVectorForEach(x, [&](
auto&& entry, std::size_t){
855 using std::isnan, std::isinf;
856 if (isnan(entry) || isinf(entry))
857 result.converged =
false;
862 std::shared_ptr<MatrixForSolver> matrix_;
864 std::shared_ptr<InverseOperator> solver_;
870#include <dune/istl/superlu.hh>
882template<
class LSTraits,
class LATraits>
883using SuperLUIstlSolver = Detail::DirectIstlSolver<LSTraits, LATraits, Dune::SuperLU>;
890#include <dune/istl/umfpack.hh>
902template<
class LSTraits,
class LATraits>
903using UMFPackIstlSolver = Detail::DirectIstlSolver<
904 LSTraits, LATraits, Dune::UMFPack
905#if DUNE_VERSION_GTE(DUNE_ISTL,2,10)
Direct dune-istl linear solvers.
Definition: istlsolvers.hh:752
void setMatrix(std::shared_ptr< Matrix > A)
Set the matrix A of the linear system Ax = b for reuse.
Definition: istlsolvers.hh:786
IstlSolverResult solve(const Matrix &A, XVector &x, const BVector &b)
Solve the linear system Ax = b.
Definition: istlsolvers.hh:767
std::string name() const
name of the linear solver
Definition: istlsolvers.hh:806
void setMatrix(Matrix &A)
Set the matrix A of the linear system Ax = b for reuse.
Definition: istlsolvers.hh:800
IstlSolverResult solve(XVector &x, const BVector &b)
Solve the linear system Ax = b using the matrix set with setMatrix.
Definition: istlsolvers.hh:775
Standard dune-istl iterative linear solvers.
Definition: istlsolvers.hh:159
IstlIterativeLinearSolver(const ParameterInitializer ¶ms="")
Constructor for sequential solvers.
Definition: istlsolvers.hh:188
void setMatrix(Matrix &A)
Set the matrix A of the linear system Ax = b for reuse.
Definition: istlsolvers.hh:275
void setResidualReduction(double residReduction)
Set the residual reduction tolerance.
Definition: istlsolvers.hh:328
void setMaxIter(std::size_t maxIter)
Set the maximum number of linear solver iterations.
Definition: istlsolvers.hh:340
IstlIterativeLinearSolver(const GridView &gridView, const DofMapper &dofMapper, const ParameterInitializer ¶ms="")
Constructor for parallel and sequential solvers.
Definition: istlsolvers.hh:202
const std::string & name() const
The name of the linear solver.
Definition: istlsolvers.hh:320
void setParams(const ParameterInitializer ¶ms)
Set the linear solver parameters.
Definition: istlsolvers.hh:354
IstlSolverResult solve(Matrix &A, XVector &x, BVector &b)
Solve the linear system Ax = b.
Definition: istlsolvers.hh:259
void setMatrix(std::shared_ptr< Matrix > A)
Set the matrix A of the linear system Ax = b for reuse.
Definition: istlsolvers.hh:265
IstlSolverResult solve(XVector &x, BVector &b) const
Solve the linear system Ax = b where A has been set with setMatrix.
Definition: istlsolvers.hh:281
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:235
Scalar norm(const XVector &x) const
Compute the 2-norm of vector x.
Definition: istlsolvers.hh:292
Definition: istlsolvers.hh:59
auto operator()(TL typeList, const M &matrix, const Dune::ParameterTree &config)
Definition: istlsolvers.hh:62
Definition: istlsolvers.hh:75
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:52
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:88
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:140
IstlSolverResult(IstlSolverResult &&)=default
IstlSolverResult(const Dune::InverseOperatorResult &o)
Definition: istlsolvers.hh:145
IstlSolverResult(Dune::InverseOperatorResult &&o)
Definition: istlsolvers.hh:146
IstlSolverResult()=default
IstlSolverResult(const IstlSolverResult &)=default
std::decay_t< decltype(MatrixConverter< M >::multiTypeToBCRSMatrix(std::declval< M >()))> type
Definition: istlsolvers.hh:95
Definition: istlsolvers.hh:91
M type
Definition: istlsolvers.hh:91
typename VectorForSolver< typename LATraits::Vector, convert >::type V
Definition: istlsolvers.hh:129
std::variant< std::shared_ptr< typename LSTraits::template Sequential< M, V >::LinearOperator > > type
Definition: istlsolvers.hh:132
typename MatrixForSolver< typename LATraits::Matrix, convert >::type M
Definition: istlsolvers.hh:128
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:117
typename VectorForSolver< typename LATraits::Vector, convert >::type V
Definition: istlsolvers.hh:111
typename MatrixForSolver< typename LATraits::Matrix, convert >::type M
Definition: istlsolvers.hh:110
Definition: istlsolvers.hh:105
std::decay_t< decltype(VectorConverter< V >::multiTypeToBlockVector(std::declval< V >()))> type
Definition: istlsolvers.hh:102
Definition: istlsolvers.hh:98
V type
Definition: istlsolvers.hh:98
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.