add collective communication to generic vanguard

This commit is contained in:
Arne Morten Kvarving 2021-05-12 10:24:53 +02:00
parent ac4ee51a5f
commit d99f642dff
3 changed files with 36 additions and 12 deletions

View File

@ -61,6 +61,7 @@ std::unique_ptr<EclipseState> EclGenericVanguard::externalEclState_;
std::unique_ptr<Schedule> EclGenericVanguard::externalEclSchedule_;
std::unique_ptr<SummaryConfig> EclGenericVanguard::externalEclSummaryConfig_;
std::unique_ptr<UDQState> EclGenericVanguard::externalUDQState_;
std::unique_ptr<EclGenericVanguard::CommunicationType> EclGenericVanguard::comm_;
EclGenericVanguard::EclGenericVanguard()
: python(std::make_shared<Python>())

View File

@ -29,6 +29,10 @@
#include <opm/grid/common/GridEnums.hpp>
#include <dune/common/version.hh>
#include <dune/common/parallel/collectivecommunication.hh>
#include <dune/common/parallel/mpihelper.hh>
#include <array>
#include <memory>
#include <optional>
@ -57,6 +61,12 @@ class EclGenericVanguard {
public:
using ParallelWellStruct = std::vector<std::pair<std::string,bool>>;
#if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 7)
using CommunicationType = Dune::Communication<Dune::MPIHelper::MPICommunicator>;
#else
using CommunicationType = Dune::CollectiveCommunication<Dune::MPIHelper::MPICommunicator>;
#endif
/*!
* \brief Constructor.
* \details Needs to be in compile unit.
@ -260,6 +270,17 @@ public:
const ParallelWellStruct& parallelWells() const
{ return parallelWells_; }
//! \brief Set global communication.
static void setCommunication(std::unique_ptr<CommunicationType> comm)
{ comm_ = std::move(comm); }
//! \brief Obtain global communicator.
static CommunicationType& comm()
{
assert(comm_);
return *comm_;
}
protected:
void updateOutputDir_(std::string outputDir,
bool enableEclCompatFile);
@ -279,7 +300,7 @@ protected:
static std::unique_ptr<Schedule> externalEclSchedule_;
static std::unique_ptr<SummaryConfig> externalEclSummaryConfig_;
static std::unique_ptr<UDQState> externalUDQState_;
static std::unique_ptr<CommunicationType> comm_;
std::string caseName_;
std::string fileName_;

View File

@ -98,11 +98,6 @@ namespace Opm {
// with incorrect locale settings.
resetLocale();
# if HAVE_DUNE_FEM
Dune::Fem::MPIManager::initialize(argc, argv);
# else
Dune::MPIHelper::instance(argc, argv);
# endif
FlowMainEbos<TypeTag> mainfunc(argc, argv, outputCout, outputFiles);
return mainfunc.execute();
}
@ -152,6 +147,14 @@ namespace Opm
{
}
~Main()
{
#if HAVE_MPI && !HAVE_DUNE_FEM
EclGenericVanguard::setCommunication(nullptr);
MPI_Finalize();
#endif
}
int runDynamic()
{
int exitCode = EXIT_SUCCESS;
@ -361,12 +364,11 @@ namespace Opm
Dune::Fem::MPIManager::initialize(argc_, argv_);
int mpiRank = Dune::Fem::MPIManager::rank();
#else
// the design of the plain dune MPIHelper class is quite flawed: there is no way to
// get the instance without having the argc and argv parameters available and it is
// not possible to determine the MPI rank and size without an instance. (IOW: the
// rank() and size() methods are supposed to be static.)
const auto& mpiHelper = Dune::MPIHelper::instance(argc_, argv_);
int mpiRank = mpiHelper.rank();
#if HAVE_MPI
MPI_Init(&argc_, &argv_);
#endif
EclGenericVanguard::setCommunication(std::make_unique<EclGenericVanguard::CommunicationType>());
int mpiRank = EclGenericVanguard::comm().rank();
#endif
// we always want to use the default locale, and thus spare us the trouble