1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65
|
/*
Copyright 2019 Equinor ASA.
This file is part of the Open Porous Media project (OPM).
OPM is free software: you can redistribute it and/or modify
it under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OPM is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
GNU General Public License for more details.
You should have received a copy of the GNU General Public License
along with OPM. If not, see <http://www.gnu.org/licenses/>.
*/
#include <chrono>
#include <opm/input/eclipse/Schedule/SummaryState.hpp>
#include <opm/common/utility/TimeService.hpp>
#include <pybind11/stl.h>
#include <pybind11/chrono.h>
#include "export.hpp"
#include <python/cxx/OpmCommonPythonDoc.hpp>
namespace {
std::vector<std::string> groups(const SummaryState * st) {
return st->groups();
}
std::vector<std::string> wells(const SummaryState * st) {
return st->wells();
}
}
void python::common::export_SummaryState(py::module& module) {
using namespace Opm::Common::DocStrings;
py::class_<SummaryState, std::shared_ptr<SummaryState>>(module, "SummaryState", SummaryStateClass_docstring)
.def(py::init<std::time_t>())
.def("update", &SummaryState::update, py::arg("variable_name"), py::arg("value"), SummaryState_update_docstring)
.def("update_well_var", &SummaryState::update_well_var, py::arg("well_name"), py::arg("variable_name"), py::arg("new_value"), SummaryState_update_well_var_docstring)
.def("update_group_var", &SummaryState::update_group_var, py::arg("group_name"), py::arg("variable_name"), py::arg("new_value"), SummaryState_update_group_var_docstring)
.def("well_var", py::overload_cast<const std::string&, const std::string&>(&SummaryState::get_well_var, py::const_), py::arg("well_name"), py::arg("variable_name"), SummaryState_well_var_docstring)
.def("group_var", py::overload_cast<const std::string&, const std::string&>(&SummaryState::get_group_var, py::const_), py::arg("group_name"), py::arg("variable_name"), SummaryState_group_var_docstring)
.def("elapsed", &SummaryState::get_elapsed, SummaryState_elapsed_docstring)
.def_property_readonly("groups", groups, SummaryState_groups_docstring)
.def_property_readonly("wells", wells, SummaryState_wells_docstring)
.def("__contains__", &SummaryState::has, py::arg("variable_name"), SummaryState_contains_docstring)
.def("has_well_var", py::overload_cast<const std::string&, const std::string&>(&SummaryState::has_well_var, py::const_), py::arg("well_name"), py::arg("variable_name"), SummaryState_has_well_var_docstring)
.def("has_group_var", py::overload_cast<const std::string&, const std::string&>(&SummaryState::has_group_var, py::const_), py::arg("group_name"), py::arg("variable_name"), SummaryState_has_group_var_docstring)
.def("__setitem__", &SummaryState::set, py::arg("variable_name"), py::arg("new_value"), SummaryState_setitem_docstring)
.def("__getitem__", py::overload_cast<const std::string&>(&SummaryState::get, py::const_), py::arg("variable_name"), SummaryState_getitem_docstring)
;
}
|