From eff87db94257bd3508c9a27ca771f921184bb3ac Mon Sep 17 00:00:00 2001 From: Kilian Weishaupt <kilian.weishaupt@iws.uni-stuttgart.de> Date: Tue, 5 Dec 2017 16:03:00 +0100 Subject: [PATCH] [properties] Remove LinearSolverBlockSize property * can be extracted from BCRS matrix directly * check for Matrix == BCRS matrix in a more concise way --- dumux/discretization/staggered/properties.hh | 2 -- dumux/linear/seqsolverbackend.hh | 29 ++++++++++++-------- 2 files changed, 17 insertions(+), 14 deletions(-) diff --git a/dumux/discretization/staggered/properties.hh b/dumux/discretization/staggered/properties.hh index b50b154b0d..7c42f72210 100644 --- a/dumux/discretization/staggered/properties.hh +++ b/dumux/discretization/staggered/properties.hh @@ -248,8 +248,6 @@ public: SET_TYPE_PROP(StaggeredModel, ElementSolutionVector, Dune::BlockVector<typename GET_PROP_TYPE(TypeTag, CellCenterPrimaryVariables)>); -SET_INT_PROP(StaggeredModel, LinearSolverBlockSize, 1); - //! Boundary types at a single degree of freedom SET_PROP(StaggeredModel, BoundaryTypes) { diff --git a/dumux/linear/seqsolverbackend.hh b/dumux/linear/seqsolverbackend.hh index 101cfab0a9..6549f61b16 100644 --- a/dumux/linear/seqsolverbackend.hh +++ b/dumux/linear/seqsolverbackend.hh @@ -37,6 +37,13 @@ namespace Dumux { +//! Helper type to determine whether a given type is a Dune::BCRSMatrix +template<class T> struct is_BCRSMatrix : public std::false_type {}; + +//! Helper type to determine whether a given type is a Dune::BCRSMatrix +template<class T> +struct is_BCRSMatrix<Dune::BCRSMatrix<T> > : public std::true_type {}; + /*! * \ingroup Linear * \brief A general solver backend allowing arbitrary preconditioners and solvers. @@ -733,13 +740,12 @@ public: template<class Matrix, class Vector> bool solve(const Matrix& A, Vector& x, const Vector& b) { - static const std::string paramGroup = GET_PROP_VALUE(TypeTag, ModelParameterGroup); - constexpr auto blockSize = GET_PROP_VALUE(TypeTag, LinearSolverBlockSize); - using MatrixBlock = typename Dune::FieldMatrix<double, blockSize, blockSize>; - using ISTLMatrix = typename Dune::BCRSMatrix<MatrixBlock>; - static_assert(std::is_same<Matrix, ISTLMatrix>::value, "SuperLU only works with BCRS matrices!"); + static_assert(is_BCRSMatrix<Matrix>::value, "SuperLU only works with BCRS matrices!"); + using BlockType = typename Matrix::block_type; + static_assert(BlockType::rows == BlockType::cols, "Matrix block must be quadratic!"); + constexpr auto blockSize = BlockType::rows; - Dune::SuperLU<ISTLMatrix> solver(A, this->verbosity() > 0); + Dune::SuperLU<Matrix> solver(A, this->verbosity() > 0); Vector bTmp(b); solver.apply(x, bTmp, result_); @@ -794,13 +800,12 @@ public: template<class Matrix, class Vector> bool solve(const Matrix& A, Vector& x, const Vector& b) { - static const std::string paramGroup = GET_PROP_VALUE(TypeTag, ModelParameterGroup); - constexpr auto blockSize = GET_PROP_VALUE(TypeTag, LinearSolverBlockSize); - using MatrixBlock = typename Dune::FieldMatrix<double, blockSize, blockSize>; - using ISTLMatrix = typename Dune::BCRSMatrix<MatrixBlock>; - static_assert(std::is_same<Matrix, ISTLMatrix>::value, "UMFPack only works with BCRS matrices!"); + static_assert(is_BCRSMatrix<Matrix>::value, "UMFPack only works with BCRS matrices!"); + using BlockType = typename Matrix::block_type; + static_assert(BlockType::rows == BlockType::cols, "Matrix block must be quadratic!"); + constexpr auto blockSize = BlockType::rows; - Dune::UMFPack<ISTLMatrix> solver(A, this->verbosity() > 0); + Dune::UMFPack<Matrix> solver(A, this->verbosity() > 0); Vector bTmp(b); solver.apply(x, bTmp, result_); -- GitLab