From d66b01a998d460559d13cd7e728e32d4ce2af309 Mon Sep 17 00:00:00 2001 From: Kilian <kilian.weishaupt@iws.uni-stuttgart.de> Date: Mon, 6 Apr 2020 20:34:47 +0200 Subject: [PATCH] [pdesolver] Put checkSizesOfSubMatrices in base class * rename checkMatrix_ to checkSizesOfSubMatrices --- dumux/common/pdesolver.hh | 24 ++++++++++++++++++++++++ dumux/linear/pdesolver.hh | 27 ++------------------------- dumux/nonlinear/newtonsolver.hh | 27 ++------------------------- 3 files changed, 28 insertions(+), 50 deletions(-) diff --git a/dumux/common/pdesolver.hh b/dumux/common/pdesolver.hh index 12ae1f9dd1..e06c5cc4c9 100644 --- a/dumux/common/pdesolver.hh +++ b/dumux/common/pdesolver.hh @@ -26,6 +26,7 @@ #include <memory> +#include <dumux/common/typetraits/matrix.hh> #include <dumux/common/timeloop.hh> namespace Dumux { @@ -101,6 +102,29 @@ protected: LinearSolver& linearSolver() { return *linearSolver_; } + /*! + * \brief Helper function to assure the MultiTypeBlockMatrix's sub-blocks have the correct sizes. + */ + template<class M> + bool checkSizesOfSubMatrices(const M& A) const + { + static_assert(isMultiTypeBlockMatrix<M>::value, "This function can only be used with MultiTypeBlockMatrix"); + bool matrixHasCorrectSize = true; + using namespace Dune::Hybrid; + using namespace Dune::Indices; + forEach(A, [&matrixHasCorrectSize](const auto& rowOfMultiTypeBlockMatrix) + { + const auto numRowsLeftMostBlock = rowOfMultiTypeBlockMatrix[_0].N(); + + forEach(rowOfMultiTypeBlockMatrix, [&matrixHasCorrectSize, &numRowsLeftMostBlock](const auto& subBlock) + { + if (subBlock.N() != numRowsLeftMostBlock) + matrixHasCorrectSize = false; + }); + }); + return matrixHasCorrectSize; + } + private: std::shared_ptr<Assembler> assembler_; std::shared_ptr<LinearSolver> linearSolver_; diff --git a/dumux/linear/pdesolver.hh b/dumux/linear/pdesolver.hh index 8214d01ed6..dca3534564 100644 --- a/dumux/linear/pdesolver.hh +++ b/dumux/linear/pdesolver.hh @@ -269,8 +269,7 @@ private: SolutionVector& x, SolutionVector& b) { - // check matrix sizes - assert(checkMatrix_(A) && "Sub blocks of MultiType matrix have wrong sizes!"); + assert(this->checkSizesOfSubMatrices(A) && "Sub-blocks of MultiTypeBlockMatrix have wrong sizes!"); // TODO: automatically derive the precondBlockLevel return ls.template solve</*precondBlockLevel=*/2>(A, x, b); @@ -294,8 +293,7 @@ private: SolutionVector& x, SolutionVector& b) { - // check matrix sizes - assert(checkMatrix_(A) && "Sub blocks of MultiType matrix have wrong sizes!"); + assert(this->checkSizesOfSubMatrices(A) && "Sub-blocks of MultiTypeBlockMatrix have wrong sizes!"); // create the bcrs matrix the IterativeSolver backend can handle const auto M = MatrixConverter<JacobianMatrix>::multiTypeToBCRSMatrix(A); @@ -323,27 +321,6 @@ private: return converged; } - //! helper method to assure the MultiType matrix's sub blocks have the correct sizes - template<class M = JacobianMatrix> - typename std::enable_if_t<!isBCRSMatrix<M>(), bool> - checkMatrix_(const JacobianMatrix& A) - { - bool matrixHasCorrectSize = true; - using namespace Dune::Hybrid; - using namespace Dune::Indices; - forEach(A, [&matrixHasCorrectSize](const auto& rowOfMultiTypeMatrix) - { - const auto numRowsLeftMostBlock = rowOfMultiTypeMatrix[_0].N(); - - forEach(rowOfMultiTypeMatrix, [&matrixHasCorrectSize, &numRowsLeftMostBlock](const auto& subBlock) - { - if (subBlock.N() != numRowsLeftMostBlock) - matrixHasCorrectSize = false; - }); - }); - return matrixHasCorrectSize; - } - //! initialize the parameters by reading from the parameter tree void initParams_(const std::string& group = "") { diff --git a/dumux/nonlinear/newtonsolver.hh b/dumux/nonlinear/newtonsolver.hh index e938fdb372..4d3c0f2446 100644 --- a/dumux/nonlinear/newtonsolver.hh +++ b/dumux/nonlinear/newtonsolver.hh @@ -1095,8 +1095,7 @@ private: SolutionVector& x, SolutionVector& b) { - // check matrix sizes - assert(checkMatrix_(A) && "Sub blocks of MultiType matrix have wrong sizes!"); + assert(this->checkSizesOfSubMatrices(A) && "Sub-blocks of MultiTypeBlockMatrix have wrong sizes!"); // TODO: automatically derive the precondBlockLevel return ls.template solve</*precondBlockLevel=*/2>(A, x, b); @@ -1120,8 +1119,7 @@ private: SolutionVector& x, SolutionVector& b) { - // check matrix sizes - assert(checkMatrix_(A) && "Sub blocks of MultiType matrix have wrong sizes!"); + assert(this->checkSizesOfSubMatrices(A) && "Sub-blocks of MultiTypeBlockMatrix have wrong sizes!"); // create the bcrs matrix the IterativeSolver backend can handle const auto M = MatrixConverter<JacobianMatrix>::multiTypeToBCRSMatrix(A); @@ -1149,27 +1147,6 @@ private: return converged; } - //! helper method to assure the MultiType matrix's sub blocks have the correct sizes - template<class M = JacobianMatrix> - typename std::enable_if_t<!isBCRSMatrix<M>(), bool> - checkMatrix_(const JacobianMatrix& A) - { - bool matrixHasCorrectSize = true; - using namespace Dune::Hybrid; - using namespace Dune::Indices; - forEach(A, [&matrixHasCorrectSize](const auto& rowOfMultiTypeMatrix) - { - const auto numRowsLeftMostBlock = rowOfMultiTypeMatrix[_0].N(); - - forEach(rowOfMultiTypeMatrix, [&matrixHasCorrectSize, &numRowsLeftMostBlock](const auto& subBlock) - { - if (subBlock.N() != numRowsLeftMostBlock) - matrixHasCorrectSize = false; - }); - }); - return matrixHasCorrectSize; - } - //! initialize the parameters by reading from the parameter tree void initParams_(const std::string& group = "") { -- GitLab