Expose "REGDIMS" Keyword from TableManager

Some of these values are needed for restart purposes (go into the
INTEHEAD) vector.
This commit is contained in:
Bård Skaflestad
2018-03-22 00:00:48 +01:00
parent b4d69b929d
commit 15eb84d3e9
3 changed files with 50 additions and 0 deletions

View File

@@ -64,6 +64,7 @@ namespace Opm {
const Tabdims& getTabdims() const;
const Eqldims& getEqldims() const;
const Aqudims& getAqudims() const;
const Regdims& getRegdims() const;
/*
WIll return max{ Tabdims::NTFIP , Regdims::NTFIP }.

View File

@@ -566,6 +566,10 @@ namespace Opm {
return m_aqudims;
}
const Regdims& TableManager::getRegdims() const {
return *this->m_regdims;
}
/*
const std::vector<SwofTable>& TableManager::getSwofTables() const {
return m_swofTables;

View File

@@ -31,6 +31,7 @@
// keyword specific table classes
#include <opm/parser/eclipse/EclipseState/Tables/PlyrockTable.hpp>
#include <opm/parser/eclipse/EclipseState/Tables/Regdims.hpp>
#include <opm/parser/eclipse/EclipseState/Tables/SwofTable.hpp>
#include <opm/parser/eclipse/EclipseState/Tables/SgwfnTable.hpp>
#include <opm/parser/eclipse/EclipseState/Tables/SwfnTable.hpp>
@@ -1357,3 +1358,47 @@ BOOST_AUTO_TEST_CASE( TestParseTABDIMS ) {
Opm::Parser parser;
BOOST_CHECK_NO_THROW( parser.parseString(data, Opm::ParseContext()));
}
BOOST_AUTO_TEST_CASE (Regdims_Entries) {
// All defaulted
{
const auto input = std::string {
R"~(
REGDIMS
/
)~" };
const auto tabMgr = ::Opm::TableManager {
::Opm::Parser{}.parseString(input)
};
const auto& rd = tabMgr.getRegdims();
BOOST_CHECK_EQUAL(rd.getNTFIP() , std::size_t{1});
BOOST_CHECK_EQUAL(rd.getNMFIPR(), std::size_t{1});
BOOST_CHECK_EQUAL(rd.getNRFREG(), std::size_t{0});
BOOST_CHECK_EQUAL(rd.getNTFREG(), std::size_t{0});
}
// All user-specified
{
const auto input = std::string {
R"~(
REGDIMS
11 22 33 44 55 66 77 88 99 110
/
)~" };
const auto tabMgr = ::Opm::TableManager {
::Opm::Parser{}.parseString(input)
};
const auto& rd = tabMgr.getRegdims();
BOOST_CHECK_EQUAL(rd.getNTFIP() , std::size_t{11});
BOOST_CHECK_EQUAL(rd.getNMFIPR(), std::size_t{22});
BOOST_CHECK_EQUAL(rd.getNRFREG(), std::size_t{33});
BOOST_CHECK_EQUAL(rd.getNTFREG(), std::size_t{44});
}
}