Support at Least NTFIP Distinict Regions for Inter-Region Flow

This commit ensures that we have backing support for region IDs up
to and including NTFIP (TABDIMS(5), REGDIMS(1)).  The existing setup
would fail (segmentation violation) in the case of a summary vector
of the form

    ROFT
      36 31 /
    /

when the maximum FIPNUM value was 30.  We nevertheless support
maximum FIPNUM values exceeding NTFIP.

We add a new optional parameter to the EclInterRegionFlowMap
constructor.  The parameter allows client code to specifiy a
"minimum maximum" region ID that all ranks must support.  This value
will be enforced during parallel aggregation.
This commit is contained in:
Bård Skaflestad
2024-01-23 18:54:46 +01:00
parent 24ebb77257
commit eb9ead5577
4 changed files with 108 additions and 4 deletions

View File

@@ -31,12 +31,17 @@
#include <opm/output/data/Solution.hpp>
#include <opm/input/eclipse/EclipseState/EclipseState.hpp>
#include <opm/input/eclipse/EclipseState/Runspec.hpp>
#include <opm/input/eclipse/EclipseState/SummaryConfig/SummaryConfig.hpp>
#include <opm/input/eclipse/EclipseState/Tables/Regdims.hpp>
#include <opm/input/eclipse/EclipseState/Tables/Tabdims.hpp>
#include <opm/input/eclipse/Schedule/RFTConfig.hpp>
#include <opm/input/eclipse/Schedule/Schedule.hpp>
#include <opm/input/eclipse/Schedule/SummaryState.hpp>
#include <opm/input/eclipse/Schedule/Well/Well.hpp>
#include <opm/input/eclipse/Schedule/Well/WellConnections.hpp>
#include <opm/input/eclipse/Units/Units.hpp>
#include <opm/simulators/utils/PressureAverage.hpp>
@@ -136,6 +141,12 @@ std::string EclString(const Opm::Inplace::Phase phase)
return regions;
}
std::size_t declaredMaxRegionID(const Opm::Runspec& rspec)
{
return std::max(rspec.tabdims().getNumFIPRegions(),
rspec.regdims().getNTFIP());
}
}
namespace Opm {
@@ -160,7 +171,9 @@ EclGenericOutputBlackoilModule(const EclipseState& eclState,
, schedule_(schedule)
, summaryConfig_(summaryConfig)
, summaryState_(summaryState)
, interRegionFlows_(numCells(eclState), defineInterRegionFlowArrays(eclState, summaryConfig))
, interRegionFlows_(numCells(eclState),
defineInterRegionFlowArrays(eclState, summaryConfig),
declaredMaxRegionID(eclState.runspec()))
, logOutput_(eclState, schedule, summaryState)
, enableEnergy_(enableEnergy)
, enableTemperature_(enableTemperature)

View File

@@ -105,7 +105,10 @@ assignGlobalMaxRegionID(const std::size_t regID)
return false;
}
this->maxGlobalRegionID_ = regID;
if (regID > this->maxGlobalRegionID_) {
this->maxGlobalRegionID_ = regID;
}
return true;
}
@@ -128,7 +131,8 @@ Opm::EclInterRegFlowMap::createMapFromNames(std::vector<std::string> names)
Opm::EclInterRegFlowMap::
EclInterRegFlowMap(const std::size_t numCells,
const std::vector<SingleRegion>& regions)
const std::vector<SingleRegion>& regions,
const std::size_t declaredMaxRegID)
{
this->regionMaps_.reserve(regions.size());
this->names_.reserve(regions.size());
@@ -139,6 +143,12 @@ EclInterRegFlowMap(const std::size_t numCells,
this->regionMaps_.emplace_back(region.definition);
this->names_.push_back(region.name);
}
if (declaredMaxRegID > std::size_t{0}) {
for (auto& regionMap : this->regionMaps_) {
regionMap.assignGlobalMaxRegionID(declaredMaxRegID);
}
}
}
void

View File

@@ -210,8 +210,13 @@ namespace Opm {
/// overlap cells if applicable.
///
/// \param[in] regions All applicable region definition arrays.
///
/// \param[in] declaredMaxRegID Declared maximum region ID in the
/// run-typically from the TABDIMS and/or REGDIMS keywords. Used
/// for sizing internal data structures if greater than zero.
explicit EclInterRegFlowMap(const std::size_t numCells,
const std::vector<SingleRegion>& regions);
const std::vector<SingleRegion>& regions,
const std::size_t declaredMaxRegID = 0);
EclInterRegFlowMap(const EclInterRegFlowMap& rhs) = default;
EclInterRegFlowMap(EclInterRegFlowMap&& rhs) noexcept = default;

View File

@@ -3382,4 +3382,80 @@ BOOST_AUTO_TEST_CASE(Four_Processes)
}
}
BOOST_AUTO_TEST_CASE(Four_Processes_Declared_MaxID)
{
using Map = Opm::EclInterRegFlowMap;
const auto fipnum = all_same_region();
const auto fipspl = left_right_split_region();
const auto fipchk = checker_board_region();
const auto fipsep = all_separate_region();
auto rank = std::vector<Map>(4, Map {
fipnum.size(), {
{ "FIPNUM", std::cref(fipnum) },
{ "FIPSPL", std::cref(fipspl) },
{ "FIPCHK", std::cref(fipchk) },
{ "FIPSEP", std::cref(fipsep) },
}, 42});
addConnections(FourProc::P1::isInterior, rank[0]);
rank[0].assignGlobalMaxRegionID({3, 2, 50, 4});
addConnections(FourProc::P2::isInterior, rank[1]);
addConnections(FourProc::P3::isInterior, rank[2]);
addConnections(FourProc::P4::isInterior, rank[3]);
rank[3].assignGlobalMaxRegionID({5, 2, 2, 4});
for (auto& r : rank) {
r.compress();
}
auto buffer = makeMessageBuffer();
rank[1].write(buffer);
rank[0].write(buffer);
rank[3].write(buffer);
rank[2].read(buffer);
rank[2].read(buffer);
rank[2].read(buffer);
BOOST_CHECK_MESSAGE(rank[2].readIsConsistent(),
"Consistent write() must yield consistent read()");
BOOST_CHECK_THROW(rank[2].addConnection({ 0, TwoProc::P1::isInterior(0) },
{ 1, TwoProc::P1::isInterior(1) }, conn_01()),
std::logic_error);
rank[2].compress();
const auto iregFlows = rank[2].getInterRegFlows();
BOOST_REQUIRE_EQUAL(iregFlows.size(), std::size_t{4});
// FIPNUM
{
const auto& map = iregFlows[0];
BOOST_CHECK_EQUAL(map.numRegions(), 42);
for (auto r1 = 0*map.numRegions(); r1 < map.numRegions(); ++r1) {
for (auto r2 = r1 + 1; r2 < map.numRegions(); ++r2) {
auto flow = map.getInterRegFlows(r1, r2);
BOOST_CHECK_MESSAGE(! flow.has_value(),
"There must not be inter-regional flow "
"between regions " << r1 << " and " << r2);
}
}
}
// FIPSPL
BOOST_CHECK_EQUAL(iregFlows[1].numRegions(), 42);
// FIPCHK
BOOST_CHECK_EQUAL(iregFlows[2].numRegions(), 50);
// FIPSEP
BOOST_CHECK_EQUAL(iregFlows[3].numRegions(), 42);
}
BOOST_AUTO_TEST_SUITE_END() // MultiArray_Wrapper