Add Support for Reading MPI Partitioning From File

This commit introduces new, experimental support for loading a
partitioning of the cells from a text file.  The name of the file is
passed into the simulator using the new, hidden, command line option

    --external-partition=filename

and we perform some basic checking that the number of elements in the
partition matches the number of cells in the CpGrid object.
This commit is contained in:
Bård Skaflestad 2023-09-18 10:21:31 +02:00
parent 960663e7b8
commit 610c45aa77
5 changed files with 79 additions and 2 deletions

View File

@ -119,6 +119,12 @@ struct ZoltanParams {
using type = UndefinedProperty;
};
template <class TypeTag, class MyTypeTag>
struct ExternalPartition
{
using type = UndefinedProperty;
};
template<class TypeTag, class MyTypeTag>
struct AllowDistributedWells {
using type = UndefinedProperty;
@ -179,6 +185,12 @@ struct ZoltanParams<TypeTag,TTag::EclBaseVanguard> {
static constexpr auto value = "graph";
};
template <class TypeTag>
struct ExternalPartition<TypeTag, TTag::EclBaseVanguard>
{
static constexpr auto* value = "";
};
template<class TypeTag>
struct AllowDistributedWells<TypeTag, TTag::EclBaseVanguard> {
static constexpr bool value = false;
@ -268,6 +280,12 @@ public:
"See https://sandialabs.github.io/Zoltan/ug_html/ug.html "
"for available Zoltan options.");
EWOMS_HIDE_PARAM(TypeTag, ZoltanParams);
EWOMS_REGISTER_PARAM(TypeTag, std::string, ExternalPartition,
"Name of file from which to load an externally generated "
"partitioning of the model's active cells for MPI "
"distribution purposes. If empty, the built-in partitioning "
"method will be employed.");
EWOMS_HIDE_PARAM(TypeTag, ExternalPartition);
#endif
EWOMS_REGISTER_PARAM(TypeTag, bool, AllowDistributedWells,
"Allow the perforations of a well to be distributed to interior of multiple processes");
@ -297,6 +315,7 @@ public:
serialPartitioning_ = EWOMS_GET_PARAM(TypeTag, bool, SerialPartitioning);
zoltanImbalanceTol_ = EWOMS_GET_PARAM(TypeTag, double, ZoltanImbalanceTol);
zoltanParams_ = EWOMS_GET_PARAM(TypeTag, std::string, ZoltanParams);
externalPartitionFile_ = EWOMS_GET_PARAM(TypeTag, std::string, ExternalPartition);
#endif
enableDistributedWells_ = EWOMS_GET_PARAM(TypeTag, bool, AllowDistributedWells);
ignoredKeywords_ = EWOMS_GET_PARAM(TypeTag, std::string, IgnoreKeywords);

View File

@ -38,11 +38,58 @@
#include <opm/models/blackoil/blackoilproperties.hh>
#include <array>
#include <filesystem>
#include <fstream>
#include <functional>
#include <iterator>
#include <memory>
#include <stdexcept>
#include <tuple>
#include <vector>
#include <fmt/format.h>
#if HAVE_MPI
namespace Opm { namespace details {
class MPIPartitionFromFile
{
public:
explicit MPIPartitionFromFile(const std::filesystem::path& partitionFile)
: partitionFile_(partitionFile)
{}
std::vector<int> operator()(const Dune::CpGrid& grid) const;
private:
std::filesystem::path partitionFile_{};
};
inline std::vector<int>
MPIPartitionFromFile::operator()(const Dune::CpGrid& grid) const
{
std::ifstream pfile { this->partitionFile_ };
auto partition = std::vector<int> {
std::istream_iterator<int> { pfile },
std::istream_iterator<int> {}
};
if (partition.size() != static_cast<std::vector<int>::size_type>(grid.size(0))) {
throw std::invalid_argument {
fmt::format("Partition file '{}' with {} values does "
"not match CpGrid instance with {} cells",
this->partitionFile_.generic_string(),
partition.size(), grid.size(0))
};
}
return partition;
}
}} // namespace Opm::details
#endif // HAVE_MPI
namespace Opm {
template <class TypeTag>
class EclCpGridVanguard;
@ -224,6 +271,12 @@ public:
void loadBalance()
{
#if HAVE_MPI
if (const auto& extPFile = this->externalPartitionFile();
!extPFile.empty() && (extPFile != "none"))
{
this->setExternalLoadBalancer(details::MPIPartitionFromFile { extPFile });
}
this->doLoadBalance_(this->edgeWeightsMethod(), this->ownersFirst(),
this->serialPartitioning(), this->enableDistributedWells(),
this->zoltanImbalanceTol(), this->gridView(),

View File

@ -52,6 +52,7 @@
#include <ebos/femcpgridcompat.hh>
#endif //HAVE_DUNE_FEM
#include <algorithm>
#include <cassert>
#include <numeric>
#include <optional>

View File

@ -228,6 +228,11 @@ public:
*/
double zoltanImbalanceTol() const
{ return zoltanImbalanceTol_; }
const std::string& externalPartitionFile() const
{
return this->externalPartitionFile_;
}
#endif
/*!
@ -292,6 +297,7 @@ protected:
bool serialPartitioning_;
double zoltanImbalanceTol_;
std::string zoltanParams_;
std::string externalPartitionFile_{};
#endif
bool enableDistributedWells_;
std::string ignoredKeywords_;

View File

@ -97,7 +97,6 @@ struct LoadFile
using type = UndefinedProperty;
};
template<class TypeTag>
struct EnableTerminalOutput<TypeTag, TTag::EclFlowProblem> {
static constexpr bool value = true;
@ -135,7 +134,6 @@ struct LoadFile<TypeTag, TTag::EclFlowProblem>
static constexpr auto* value = "";
};
template <class TypeTag>
struct LoadStep<TypeTag, TTag::EclFlowProblem>
{