From b9b6f7defd458203813243707d2833913335baa8 Mon Sep 17 00:00:00 2001 From: Markus Blatt Date: Mon, 10 Oct 2022 21:58:46 +0200 Subject: [PATCH] Use Communication instead of CollectiveCommunication for DUNE>=2.7. The later is removed in 2.9. --- opm/models/nonlinear/newtonmethod.hh | 8 ++++++++ opm/models/utils/simulator.hh | 19 +++++++++++++++++-- .../linalg/overlappingscalarproduct.hh | 13 +++++++++++-- 3 files changed, 36 insertions(+), 4 deletions(-) diff --git a/opm/models/nonlinear/newtonmethod.hh b/opm/models/nonlinear/newtonmethod.hh index 68a051b9e..11f4414ef 100644 --- a/opm/models/nonlinear/newtonmethod.hh +++ b/opm/models/nonlinear/newtonmethod.hh @@ -181,14 +181,22 @@ class NewtonMethod using ConvergenceWriter = GetPropType; using Communicator = typename Dune::MPIHelper::MPICommunicator; +#if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 7) + using CollectiveCommunication = typename Dune::Communication; +#else using CollectiveCommunication = Dune::CollectiveCommunication; +#endif public: NewtonMethod(Simulator& simulator) : simulator_(simulator) , endIterMsgStream_(std::ostringstream::out) , linearSolver_(simulator) +#if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 7) , comm_(Dune::MPIHelper::getCommunicator()) +#else + , comm_(Dune::MPIHelper::getCollectiveCommunicator()) +#endif , convergenceWriter_(asImp_()) { lastError_ = 1e100; diff --git a/opm/models/utils/simulator.hh b/opm/models/utils/simulator.hh index cd5fc6607..48da9b7d0 100644 --- a/opm/models/utils/simulator.hh +++ b/opm/models/utils/simulator.hh @@ -47,9 +47,24 @@ #include #include -#define EWOMS_CATCH_PARALLEL_EXCEPTIONS_FATAL(code) \ +namespace Opm +{ +namespace detail +{ +inline auto getMPIHelperCommunication() +{ +#if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 7) + return Dune::MPIHelper::getCommunication(); +#else + return Dune::MPIHelper::getCollectiveCommunication(); +#endif +} +} // end namespace detail +} // end namespace Opm + +#define EWOMS_CATCH_PARALLEL_EXCEPTIONS_FATAL(code) \ { \ - const auto& comm = Dune::MPIHelper::getCollectiveCommunication(); \ + const auto& comm = ::Opm::detail::getMPIHelperCommunication(); \ bool exceptionThrown = false; \ try { code; } \ catch (const Dune::Exception& e) { \ diff --git a/opm/simulators/linalg/overlappingscalarproduct.hh b/opm/simulators/linalg/overlappingscalarproduct.hh index eff0a8af3..a04325386 100644 --- a/opm/simulators/linalg/overlappingscalarproduct.hh +++ b/opm/simulators/linalg/overlappingscalarproduct.hh @@ -43,8 +43,12 @@ class OverlappingScalarProduct { public: using field_type = typename OverlappingBlockVector::field_type; - using CollectiveCommunication = typename Dune::CollectiveCommunication; +#if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 7) + using CollectiveCommunication = typename Dune::Communication; +#else + using CollectiveCommunication = typename Dune::CollectiveCommunication; +#endif using real_type = typename Dune::ScalarProduct::real_type; //! the kind of computations supported by the operator. Either overlapping or non-overlapping @@ -52,7 +56,12 @@ public: { return Dune::SolverCategory::overlapping; } OverlappingScalarProduct(const Overlap& overlap) - : overlap_(overlap), comm_( Dune::MPIHelper::getCollectiveCommunication() ) + : overlap_(overlap), +#if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 7) + comm_( Dune::MPIHelper::getCommunication() ) +#else + comm_( Dune::MPIHelper::getCollectiveCommunication() ) +#endif {} #if DUNE_VERSION_NEWER(DUNE_ISTL, 2,7)