Merge pull request #372 from chflo/OPM_139_write_ecl_restart_wellinfo
OPM-139: Added methods to Schedule
This commit is contained in:
commit
bb9cf64776
@ -537,6 +537,12 @@ namespace Opm {
|
|||||||
return m_wells.size();
|
return m_wells.size();
|
||||||
}
|
}
|
||||||
|
|
||||||
|
size_t Schedule::numWells(size_t timestep) const {
|
||||||
|
std::vector<WellConstPtr> wells = getWells(timestep);
|
||||||
|
return wells.size();
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
bool Schedule::hasWell(const std::string& wellName) const {
|
bool Schedule::hasWell(const std::string& wellName) const {
|
||||||
return m_wells.hasKey( wellName );
|
return m_wells.hasKey( wellName );
|
||||||
}
|
}
|
||||||
@ -663,4 +669,19 @@ namespace Opm {
|
|||||||
}
|
}
|
||||||
else throw std::invalid_argument("String " + eclipseString + " not recognized as a boolean-convertible string.");
|
else throw std::invalid_argument("String " + eclipseString + " not recognized as a boolean-convertible string.");
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
|
size_t Schedule::getMaxNumCompletionsForWells(size_t timestep) const {
|
||||||
|
size_t ncwmax = 0;
|
||||||
|
const std::vector<WellConstPtr>& wells = getWells();
|
||||||
|
for (auto wellIter=wells.begin(); wellIter != wells.end(); ++wellIter) {
|
||||||
|
WellConstPtr wellPtr = *wellIter;
|
||||||
|
CompletionSetConstPtr completionsSetPtr = wellPtr->getCompletions(timestep);
|
||||||
|
|
||||||
|
if (completionsSetPtr->size() > ncwmax )
|
||||||
|
ncwmax = completionsSetPtr->size();
|
||||||
|
|
||||||
|
}
|
||||||
|
return ncwmax;
|
||||||
|
}
|
||||||
}
|
}
|
||||||
|
@ -46,6 +46,8 @@ namespace Opm
|
|||||||
TimeMapConstPtr getTimeMap() const;
|
TimeMapConstPtr getTimeMap() const;
|
||||||
|
|
||||||
size_t numWells() const;
|
size_t numWells() const;
|
||||||
|
size_t numWells(size_t timestep) const;
|
||||||
|
size_t getMaxNumCompletionsForWells(size_t timestep) const;
|
||||||
bool hasWell(const std::string& wellName) const;
|
bool hasWell(const std::string& wellName) const;
|
||||||
WellPtr getWell(const std::string& wellName) const;
|
WellPtr getWell(const std::string& wellName) const;
|
||||||
std::vector<WellConstPtr> getWells() const;
|
std::vector<WellConstPtr> getWells() const;
|
||||||
@ -56,6 +58,7 @@ namespace Opm
|
|||||||
size_t numGroups() const;
|
size_t numGroups() const;
|
||||||
bool hasGroup(const std::string& groupName) const;
|
bool hasGroup(const std::string& groupName) const;
|
||||||
GroupPtr getGroup(const std::string& groupName) const;
|
GroupPtr getGroup(const std::string& groupName) const;
|
||||||
|
|
||||||
|
|
||||||
private:
|
private:
|
||||||
TimeMapPtr m_timeMap;
|
TimeMapPtr m_timeMap;
|
||||||
|
@ -85,6 +85,43 @@ static DeckPtr createDeckWithWellsOrdered() {
|
|||||||
return parser.parseString(input);
|
return parser.parseString(input);
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
|
static DeckPtr createDeckWithWellsAndCompletionData() {
|
||||||
|
Opm::Parser parser;
|
||||||
|
std::string input =
|
||||||
|
"START -- 0 \n"
|
||||||
|
"1 NOV 1979 / \n"
|
||||||
|
"SCHEDULE\n"
|
||||||
|
"DATES -- 1\n"
|
||||||
|
" 1 DES 1979/ \n"
|
||||||
|
"/\n"
|
||||||
|
"WELSPECS\n"
|
||||||
|
" \'OP_1\' \'OP'\ 9 9 1* \'OIL\' 1* 1* 1* 1* 1* 1* 1* / \n"
|
||||||
|
" \'OP_2\' \'OP'\ 8 8 1* \'OIL\' 1* 1* 1* 1* 1* 1* 1* / \n"
|
||||||
|
" \'OP_3\' \'OP'\ 7 7 1* \'OIL\' 1* 1* 1* 1* 1* 1* 1* / \n"
|
||||||
|
"/\n"
|
||||||
|
"COMPDAT\n"
|
||||||
|
" \'OP_1\' 9 9 1 1 \'OPEN\' 1* 32.948 0.311 3047.839 1* 1* \'X\' 22.100 / \n"
|
||||||
|
" \'OP_1\' 9 9 2 2 \'OPEN\' 1* 46.825 0.311 4332.346 1* 1* \'X\' 22.123 / \n"
|
||||||
|
" \'OP_2\' 8 8 1 3 \'OPEN\' 1* 1.168 0.311 107.872 1* 1* \'Y\' 21.925 / \n"
|
||||||
|
" \'OP_2\' 8 7 3 3 \'OPEN\' 1* 15.071 0.311 1391.859 1* 1* \'Y\' 21.920 / \n"
|
||||||
|
" \'OP_2\' 8 7 3 6 \'OPEN\' 1* 6.242 0.311 576.458 1* 1* \'Y\' 21.915 / \n"
|
||||||
|
" \'OP_3\' 7 7 1 1 \'OPEN\' 1* 27.412 0.311 2445.337 1* 1* \'Y\' 18.521 / \n"
|
||||||
|
" \'OP_3\' 7 7 2 2 \'OPEN\' 1* 55.195 0.311 4923.842 1* 1* \'Y\' 18.524 / \n"
|
||||||
|
"/\n"
|
||||||
|
"DATES -- 2,3\n"
|
||||||
|
" 10 JUL 2007 / \n"
|
||||||
|
" 10 AUG 2007 / \n"
|
||||||
|
"/\n"
|
||||||
|
"COMPDAT\n"
|
||||||
|
" \'OP_1\' 9 9 3 9 \'OPEN\' 1* 32.948 0.311 3047.839 1* 1* \'X\' 22.100 / \n"
|
||||||
|
"/\n";
|
||||||
|
|
||||||
|
|
||||||
|
return parser.parseString(input);
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
BOOST_AUTO_TEST_CASE(CreateScheduleDeckMissingReturnsDefaults) {
|
BOOST_AUTO_TEST_CASE(CreateScheduleDeckMissingReturnsDefaults) {
|
||||||
DeckPtr deck(new Deck());
|
DeckPtr deck(new Deck());
|
||||||
DeckKeywordPtr keyword(new DeckKeyword("SCHEDULE"));
|
DeckKeywordPtr keyword(new DeckKeyword("SCHEDULE"));
|
||||||
@ -215,3 +252,22 @@ BOOST_AUTO_TEST_CASE(WellsIteratorWithRegex_HasWells_WellsReturned) {
|
|||||||
BOOST_CHECK_EQUAL(1U, wells.size());
|
BOOST_CHECK_EQUAL(1U, wells.size());
|
||||||
}
|
}
|
||||||
|
|
||||||
|
BOOST_AUTO_TEST_CASE(ReturnNumWellsTimestep) {
|
||||||
|
DeckPtr deck = createDeckWithWells();
|
||||||
|
Schedule schedule(deck);
|
||||||
|
|
||||||
|
BOOST_CHECK_EQUAL(schedule.numWells(0), 1);
|
||||||
|
BOOST_CHECK_EQUAL(schedule.numWells(1), 1);
|
||||||
|
BOOST_CHECK_EQUAL(schedule.numWells(2), 1);
|
||||||
|
BOOST_CHECK_EQUAL(schedule.numWells(3), 3);
|
||||||
|
}
|
||||||
|
|
||||||
|
BOOST_AUTO_TEST_CASE(ReturnMaxNumCompletionsForWellsInTimestep) {
|
||||||
|
DeckPtr deck = createDeckWithWellsAndCompletionData();
|
||||||
|
Schedule schedule(deck);
|
||||||
|
|
||||||
|
BOOST_CHECK_EQUAL(schedule.getMaxNumCompletionsForWells(1), 7);
|
||||||
|
BOOST_CHECK_EQUAL(schedule.getMaxNumCompletionsForWells(3), 9);
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
Loading…
Reference in New Issue
Block a user