From 148e9cae9e214d0cf80730ec9c8ca7c4d34a5018 Mon Sep 17 00:00:00 2001 From: Timo Koch <timo.koch@iws.uni-stuttgart.de> Date: Tue, 28 Feb 2023 15:11:30 +0100 Subject: [PATCH] [solver][fix] Make sure to only call communication code if we are able to communicate --- dumux/linear/istlsolvers.hh | 29 ++++++++++++++++++++--------- 1 file changed, 20 insertions(+), 9 deletions(-) diff --git a/dumux/linear/istlsolvers.hh b/dumux/linear/istlsolvers.hh index 42dfbf1914..97f5a8316a 100644 --- a/dumux/linear/istlsolvers.hh +++ b/dumux/linear/istlsolvers.hh @@ -143,16 +143,24 @@ public: { initializeParameters_(paramGroup); #if HAVE_MPI - solverCategory_ = Detail::solverCategory<LinearSolverTraits>(gridView); - if (solverCategory_ != Dune::SolverCategory::sequential) + if constexpr (LinearSolverTraits::canCommunicate) { - parallelHelper_ = std::make_unique<ParallelISTLHelper<LinearSolverTraits>>(gridView, dofMapper); - communication_ = std::make_shared<Comm>(gridView.comm(), solverCategory_); - scalarProduct_ = Dune::createScalarProduct<XVector>(*communication_, solverCategory_); - parallelHelper_->createParallelIndexSet(*communication_); + solverCategory_ = Detail::solverCategory<LinearSolverTraits>(gridView); + if (solverCategory_ != Dune::SolverCategory::sequential) + { + parallelHelper_ = std::make_unique<ParallelISTLHelper<LinearSolverTraits>>(gridView, dofMapper); + communication_ = std::make_shared<Comm>(gridView.comm(), solverCategory_); + scalarProduct_ = Dune::createScalarProduct<XVector>(*communication_, solverCategory_); + parallelHelper_->createParallelIndexSet(*communication_); + } + else + scalarProduct_ = std::make_shared<ScalarProduct>(); } else + { + solverCategory_ = Dune::SolverCategory::sequential; scalarProduct_ = std::make_shared<ScalarProduct>(); + } #else solverCategory_ = Dune::SolverCategory::sequential; scalarProduct_ = std::make_shared<ScalarProduct>(); @@ -174,10 +182,13 @@ public: solverCategory_ = Detail::solverCategory(gridView); scalarProduct_ = scalarProduct; communication_ = communication; - if (solverCategory_ != Dune::SolverCategory::sequential) + if constexpr (LinearSolverTraits::canCommunicate) { - parallelHelper_ = std::make_unique<ParallelISTLHelper<LinearSolverTraits>>(gridView, dofMapper); - parallelHelper_->createParallelIndexSet(communication); + if (solverCategory_ != Dune::SolverCategory::sequential) + { + parallelHelper_ = std::make_unique<ParallelISTLHelper<LinearSolverTraits>>(gridView, dofMapper); + parallelHelper_->createParallelIndexSet(communication); + } } } #endif -- GitLab