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 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124
|
#include <opm/input/eclipse/Deck/Deck.hpp>
#include <opm/input/eclipse/EclipseState/EclipseState.hpp>
#include <opm/input/eclipse/EclipseState/Grid/FaultCollection.hpp>
#include "export.hpp"
#include <python/cxx/OpmCommonPythonDoc.hpp>
namespace {
py::list getNNC( const EclipseState& state ) {
py::list l;
for( const auto& x : state.getInputNNC().input() )
l.append( py::make_tuple( x.cell1, x.cell2, x.trans ) );
return l;
}
py::list faultNames( const EclipseState& state ) {
py::list l;
const auto& fc = state.getFaults();
for (size_t i = 0; i < fc.size(); i++) {
const auto& f = fc.getFault(i);
l.append(f.getName());
}
return l;
}
py::dict jfunc( const EclipseState& s) {
const auto& tm = s.getTableManager();
if (!tm.useJFunc())
return py::dict();
const auto& j = tm.getJFunc();
std::string flag = "BOTH";
std::string dir = "XY";
if (j.flag() == JFunc::Flag::WATER)
flag = "WATER";
else if (j.flag() == JFunc::Flag::GAS)
flag = "GAS";
if (j.direction() == JFunc::Direction::X)
dir = "X";
else if (j.direction() == JFunc::Direction::Y)
dir = "Y";
else if (j.direction() == JFunc::Direction::Z)
dir = "Z";
py::dict ret;
ret["FLAG"] = flag;
ret["DIRECTION"] = dir;
ret["ALPHA_FACTOR"] = j.alphaFactor();
ret["BETA_FACTOR"] = j.betaFactor();
if (j.flag() == JFunc::Flag::WATER || j.flag() == JFunc::Flag::BOTH)
ret["OIL_WATER"] = j.owSurfaceTension();
if (j.flag() == JFunc::Flag::GAS || j.flag() == JFunc::Flag::BOTH)
ret["GAS_OIL"] = j.goSurfaceTension();
return ret;
}
const std::string faceDir( FaceDir::DirEnum dir ) {
switch (dir) {
case FaceDir::DirEnum::Unknown:
break;
case FaceDir::DirEnum::XPlus: return "X+";
case FaceDir::DirEnum::XMinus: return "X-";
case FaceDir::DirEnum::YPlus: return "Y+";
case FaceDir::DirEnum::YMinus: return "Y-";
case FaceDir::DirEnum::ZPlus: return "Z+";
case FaceDir::DirEnum::ZMinus: return "Z-";
}
return "Unknown direction";
}
py::list faultFaces( const EclipseState& state, const std::string& name ) {
py::list l;
const auto& gr = state.getInputGrid(); // used for global -> IJK
const auto& fc = state.getFaults();
const Fault& f = fc.getFault(name);
for (const auto& ff : f) {
// for each fault face
for (size_t g : ff) {
// for global index g in ff
const auto ijk = gr.getIJK(g);
l.append(py::make_tuple(ijk[0], ijk[1], ijk[2], faceDir(ff.getDir())));
}
}
return l;
}
const FieldPropsManager& get_field_props(const EclipseState& state) {
return state.fieldProps();
}
}
void python::common::export_EclipseState(py::module& module) {
using namespace Opm::Common::DocStrings;
// Note: In the below class we std::shared_ptr as the holder type, see:
//
// https://pybind11.readthedocs.io/en/stable/advanced/smart_ptrs.html
//
// this makes it possible to share the returned object with e.g. and
// opm.simulators.BlackOilSimulator Python object
//
py::class_< EclipseState, std::shared_ptr<EclipseState> >( module, "EclipseState", EclipseStateClass_docstring)
.def(py::init<const Deck&>(), py::arg("deck"), EclipseState_init_docstring)
.def_property_readonly( "title", &EclipseState::getTitle, EclipseState_title_docstring)
.def( "field_props", &get_field_props, ref_internal, EclipseState_field_props_docstring)
.def( "grid", &EclipseState::getInputGrid, ref_internal, EclipseState_grid_docstring)
.def( "config", &EclipseState::cfg, ref_internal, EclipseState_config_docstring)
.def( "tables", &EclipseState::getTableManager, ref_internal, EclipseState_tables_docstring)
.def( "has_input_nnc", &EclipseState::hasInputNNC, EclipseState_has_input_nnc_docstring)
.def( "simulation", &EclipseState::getSimulationConfig, ref_internal, EclipseState_simulation_docstring)
.def( "input_nnc", &getNNC, EclipseState_input_nnc_docstring)
.def( "faultNames", &faultNames, EclipseState_faultNames_docstring)
.def( "faultFaces", &faultFaces, py::arg("fault_name"), EclipseState_faultFaces_docstring)
.def( "jfunc", &jfunc, EclipseState_jfunc_docstring)
;
}
|