Merge pull request #2833 from hakonhagland/py_sched2

Implement support for adding keywords to the schedule from Python.
This commit is contained in:
Joakim Hove 2021-11-15 08:47:51 +01:00 committed by GitHub
commit 893363e225
No known key found for this signature in database
GPG Key ID: 4AEE18F83AFDEB23
4 changed files with 82 additions and 5 deletions

View File

@ -301,6 +301,7 @@ namespace Opm
for the schedule instances created by loading a restart file.
*/
static bool cmp(const Schedule& sched1, const Schedule& sched2, std::size_t report_step);
void applyKeywords(std::vector<DeckKeyword*>& keywords, std::size_t timeStep);
template<class Serializer>
void serializeOp(Serializer& serializer)

View File

@ -45,7 +45,7 @@ namespace {
void python::common::export_Deck(py::module &module) {
// Note: In the below class we std::shared_ptr as the holder type, see:
// Note: In the below class we use std::shared_ptr as the holder type, see:
//
// https://pybind11.readthedocs.io/en/stable/advanced/smart_ptrs.html
//

View File

@ -1,7 +1,9 @@
#include <ctime>
#include <chrono>
#include <opm/parser/eclipse/Parser/Parser.hpp>
#include <opm/parser/eclipse/Deck/Deck.hpp>
#include <opm/parser/eclipse/Deck/DeckKeyword.hpp>
#include <opm/parser/eclipse/EclipseState/EclipseState.hpp>
@ -92,13 +94,42 @@ namespace {
return sch[index];
}
void insert_keywords(
Schedule& sch,
const std::string& deck_string,
std::size_t index,
const UnitSystem& unit_system
)
{
Parser parser;
std::string str {unit_system.deck_name() + "\n\n" + deck_string};
auto deck = parser.parseString(str);
std::vector<DeckKeyword*> keywords;
for (auto &keyword : deck) {
keywords.push_back(&keyword);
}
sch.applyKeywords(keywords, index);
}
// NOTE: this overload does currently not work, see PR #2833. The plan
// is to fix this in a later commit. For now, the overload insert_keywords()
// above taking a deck_string (std::string) instead of a list of DeckKeywords
// has to be used instead.
void insert_keywords(
Schedule& sch, py::list& deck_keywords, std::size_t index)
{
Parser parser;
std::vector<DeckKeyword*> keywords;
for (py::handle item : deck_keywords) {
DeckKeyword &keyword = item.cast<DeckKeyword&>();
keywords.push_back(&keyword);
}
sch.applyKeywords(keywords, index);
}
}
void python::common::export_Schedule(py::module& module) {
@ -107,7 +138,7 @@ void python::common::export_Schedule(py::module& module) {
.def("group", &get_group, ref_internal);
// Note: In the below class we std::shared_ptr as the holder type, see:
// Note: In the below class we use std::shared_ptr as the holder type, see:
//
// https://pybind11.readthedocs.io/en/stable/advanced/smart_ptrs.html
//
@ -128,6 +159,14 @@ void python::common::export_Schedule(py::module& module) {
.def( "get_wells", &Schedule::getWells)
.def("well_names", py::overload_cast<const std::string&>(&Schedule::wellNames, py::const_))
.def( "get_well", &get_well)
.def( "insert_keywords",
py::overload_cast<Schedule&, py::list&, std::size_t>(&insert_keywords),
py::arg("keywords"), py::arg("step"))
.def( "insert_keywords",
py::overload_cast<
Schedule&, const std::string&, std::size_t, const UnitSystem&
>(&insert_keywords),
py::arg("data"), py::arg("step"), py::arg("unit_system"))
.def( "__contains__", &has_well );
}

View File

@ -1273,6 +1273,43 @@ void Schedule::iterateScheduleSection(std::size_t load_start, std::size_t load_e
return std::chrono::duration_cast<std::chrono::seconds>(elapsed).count();
}
void Schedule::applyKeywords(
std::vector<DeckKeyword*>& keywords, std::size_t reportStep) {
ParseContext parseContext;
ErrorGuard errors;
ScheduleGrid grid(this->completed_cells);
SimulatorUpdate sim_update;
std::unordered_map<std::string, double> target_wellpi;
std::vector<std::string> matching_wells;
const std::string prefix = "| "; /* logger prefix string */
this->snapshots.resize(reportStep + 1);
auto& input_block = this->m_sched_deck[reportStep];
for (auto keyword : keywords) {
input_block.push_back(*keyword);
this->handleKeyword(reportStep,
input_block,
*keyword,
parseContext,
errors,
grid,
matching_wells,
/*actionx_mode=*/false,
&sim_update,
&target_wellpi);
}
this->end_report(reportStep);
if (reportStep < this->m_sched_deck.size() - 1) {
iterateScheduleSection(
reportStep + 1,
this->m_sched_deck.size(),
parseContext,
errors,
grid,
&target_wellpi,
prefix);
}
}
SimulatorUpdate Schedule::applyAction(std::size_t reportStep, const time_point&, const Action::ActionX& action, const Action::Result& result, const std::unordered_map<std::string, double>& target_wellpi) {
const std::string prefix = "| ";