Consider Node Group Size When Calculating Maximum Group Size

This commit extends the 'maxGroupSize' function to also consider the
sizes of node groups (i.e., the number of child nodes/child groups)
when determining what the largest group size is on a particular step.

Add a unit test to demonstrate.
This commit is contained in:
Bård Skaflestad 2019-09-05 13:36:35 +02:00
parent 705be62be6
commit f813a04269
4 changed files with 202 additions and 17 deletions

View File

@ -318,6 +318,7 @@ list (APPEND TEST_DATA_FILES
) )
if(ENABLE_ECL_OUTPUT) if(ENABLE_ECL_OUTPUT)
list (APPEND TEST_DATA_FILES list (APPEND TEST_DATA_FILES
tests/expect-wdims.chldg.err.out
tests/expect-wdims.err.out tests/expect-wdims.err.out
tests/FIRST_SIM.DATA tests/FIRST_SIM.DATA
tests/FIRST_SIM_THPRES.DATA tests/FIRST_SIM_THPRES.DATA

View File

@ -27,7 +27,6 @@
#include <opm/parser/eclipse/EclipseState/EclipseState.hpp> #include <opm/parser/eclipse/EclipseState/EclipseState.hpp>
#include <opm/parser/eclipse/EclipseState/Runspec.hpp> #include <opm/parser/eclipse/EclipseState/Runspec.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/Group/Group2.hpp> #include <opm/parser/eclipse/EclipseState/Schedule/Group/Group2.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/Group/GTNode.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/Schedule.hpp> #include <opm/parser/eclipse/EclipseState/Schedule/Schedule.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/Well/WellConnections.hpp> #include <opm/parser/eclipse/EclipseState/Schedule/Well/WellConnections.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/ArrayDimChecker.hpp> #include <opm/parser/eclipse/EclipseState/Schedule/ArrayDimChecker.hpp>
@ -138,20 +137,6 @@ namespace {
WellDims::checkNumGroups (wdims, sched, ctxt, guard); WellDims::checkNumGroups (wdims, sched, ctxt, guard);
WellDims::checkGroupSize (wdims, sched, ctxt, guard); WellDims::checkGroupSize (wdims, sched, ctxt, guard);
} }
std::size_t groupSize(const Opm::GTNode& tree) {
auto nwgmax = std::size_t{0};
for (const auto& node : tree.groups()) {
const auto& group = node.group();
if (group.wellgroup())
nwgmax = std::max(nwgmax, group.wells().size());
else
nwgmax = std::max(nwgmax, groupSize(node));
}
return nwgmax;
}
} // Anonymous } // Anonymous
void void
@ -170,6 +155,15 @@ int
Opm::maxGroupSize(const Opm::Schedule& sched, Opm::maxGroupSize(const Opm::Schedule& sched,
const std::size_t step) const std::size_t step)
{ {
const auto& gt = sched.groupTree(step); int nwgmax = 0;
return static_cast<int>(groupSize(gt));
for (const auto& gnm : sched.groupNames(step)) {
const auto& grp = sched.getGroup2(gnm, step);
const auto gsz = grp.wellgroup()
? grp.numWells() : grp.groups().size();
nwgmax = std::max(nwgmax, static_cast<int>(gsz));
}
return nwgmax;
} }

View File

@ -0,0 +1,6 @@
Errors:
RUNSPEC_GROUPSIZE_TOO_LARGE: Run uses maximum group size of 6, but allocates at most 4 in RUNSPEC section. Increase item 4 of WELLDIMS accordingly.

View File

@ -185,6 +185,147 @@ WCONINJE
TSTEP TSTEP
100*30 / 100*30 /
END
)" };
return Opm::Parser{}.parseString(input);
}
Opm::Deck simCaseNodeGroupSizeFailure()
{
const auto input = std::string{ R"(RUNSPEC
TITLE
'Check Well Dimensions' /
DIMENS
20 20 15 /
OIL
WATER
METRIC
EQLDIMS
-- Defaulted
/
TABDIMS
-- Defaulted
/
WELLDIMS
-- NWMAX NCMAX NGMAX NWGMAX (too small)
2 20 16 4
/
-- ====================================================================
GRID
SPECGRID
20 20 15 1 F /
DXV
20*100.0 /
DYV
20*100.0 /
DZV
15*0.1 /
DEPTHZ
441*2000 /
PORO
6000*0.3 /
PERMX
6000*100.0 /
COPY
'PERMX' 'PERMY' /
'PERMX' 'PERMZ' /
/
MULTIPLY
'PERMZ' 0.1 /
/
-- ====================================================================
PROPS
SWOF
0 0 1 0
1 1 0 0 /
PVDO
1 1.0 0.5
800 0.99 0.51 /
PVTW
300 0.99 1.0e-6 0.25 0 /
DENSITY
850.0 1014.0 1.05 /
-- ====================================================================
SOLUTION
EQUIL
2000 300 2010 0.0 2000 10 /
-- ====================================================================
SUMMARY
ALL
-- ====================================================================
SCHEDULE
RPTRST
BASIC=5 FREQ=6 /
GRUPTREE
'G1' 'FIELD' /
'PLAT1' 'G1' /
'PLAT2' 'G1' /
'I-NORTH' 'PLAT1' /
'P-NORTH' 'PLAT1' /
'O-WEST' 'PLAT2' /
'I-SOUTH' 'PLAT2' /
'P-EAST' 'PLAT2' /
'G2' 'FIELD' /
'PLAT3' 'G2' /
'I-2' 'PLAT3' /
'I-21' 'PLAT3' /
'I-22' 'PLAT3' /
'I-23' 'PLAT3' /
'I-24' 'PLAT3' /
'I-25' 'PLAT3' /
/
WELSPECS
'I-N-1' 'I-NORTH' 1 1 2000.15 'WATER' /
'P-N-0' 'P-NORTH' 1 10 2000.15 'OIL' /
/
COMPDAT
'I-N-1' 0 0 2 10 'OPEN' 1* 1* 1.0 /
'P-N-0' 0 0 2 3 'OPEN' 1* 1* 1.0 /
/
WCONPROD
-- Well O/S Mode ORAT WRAT GRAT LRAT RESV BHP
'P-N-*' 'OPEN' 'LRAT' 1* 1* 1* 5E3 1* 100 /
/
WCONINJE
-- Well Type O/S Mode RATE RESV BHP
'I-N-*' 'WATER' 'OPEN' 'RATE' 25E3 1* 500 /
/
TSTEP
100*30 /
END END
)" }; )" };
@ -264,6 +405,17 @@ BOOST_AUTO_TEST_CASE(MaxGroupSize)
BOOST_CHECK_EQUAL(Opm::maxGroupSize(cse.sched, 1), 10); BOOST_CHECK_EQUAL(Opm::maxGroupSize(cse.sched, 1), 10);
} }
BOOST_AUTO_TEST_CASE(ManyChildGroups)
{
auto cse = CaseObjects {
simCaseNodeGroupSizeFailure(),
Opm::ParseContext{}
};
// 6 Child groups of 'PLAT3'
BOOST_CHECK_EQUAL(Opm::maxGroupSize(cse.sched, 1), 6);
}
BOOST_AUTO_TEST_CASE(WellDims) BOOST_AUTO_TEST_CASE(WellDims)
{ {
Opm::ParseContext parseContext; Opm::ParseContext parseContext;
@ -296,5 +448,37 @@ BOOST_AUTO_TEST_CASE(WellDims)
} }
} }
BOOST_AUTO_TEST_CASE(WellDims_ManyChildGroups)
{
Opm::ParseContext parseContext;
setWellDimsContext(Opm::InputError::THROW_EXCEPTION, parseContext);
auto cse = CaseObjects{ simCaseNodeGroupSizeFailure(), parseContext};
// There should be no failures in basic input layer
BOOST_CHECK(!cse.guard);
// There *should* be errors from dimension checking
BOOST_CHECK_THROW( Opm::checkConsistentArrayDimensions(cse.es , cse.sched,
parseContext, cse.guard),
std::invalid_argument);
setWellDimsContext(Opm::InputError::DELAYED_EXIT1, parseContext);
Opm::checkConsistentArrayDimensions(cse.es , cse.sched,
parseContext, cse.guard);
// There *should* be errors from dimension checking
BOOST_CHECK(cse.guard);
// Verify that we get expected output from ErrorGuard::dump()
boost::test_tools::output_test_stream output{"expect-wdims.chldg.err.out", true};
{
RedirectCERR stream(output.rdbuf());
cse.guard.dump();
BOOST_CHECK(output.match_pattern());
}
}
BOOST_AUTO_TEST_SUITE_END() BOOST_AUTO_TEST_SUITE_END()