additional python bindings for eclipse_grid

This commit is contained in:
Torbjørn Skille
2020-01-05 18:10:09 +01:00
parent c2b4d30be5
commit 4136868592
2 changed files with 127 additions and 2 deletions

View File

@@ -4,6 +4,8 @@
#include <opm/parser/eclipse/EclipseState/Grid/Fault.hpp>
#include <opm/parser/eclipse/EclipseState/Grid/FaceDir.hpp>
#include <pybind11/stl.h>
#include "converters.hpp"
#include "export.hpp"
@@ -39,6 +41,73 @@ namespace {
return grid.getCellVolume(i_idx, j_idx, k_idx);
}
py::array cellVolumeAll( const EclipseGrid& grid)
{
std::vector<double> cellVol;
std::array<int, 3> dims = grid.getNXYZ();
size_t nCells = dims[0]*dims[1]*dims[2];
cellVol.reserve(nCells);
for (size_t n = 0; n < nCells; n++)
cellVol.push_back(grid.getCellVolume(n));
return convert::numpy_array(cellVol);
}
py::array cellVolumeMask( const EclipseGrid& grid, std::vector<int>& mask)
{
std::array<int, 3> dims = grid.getNXYZ();
size_t nCells = dims[0]*dims[1]*dims[2];
if (nCells != mask.size())
throw std::logic_error("size of input mask doesn't match size of grid");
std::vector<double> cellVol(nCells, 0.0);
for (size_t n = 0; n < nCells; n++)
if (mask[n]==1)
cellVol[n] = grid.getCellVolume(n);
return convert::numpy_array(cellVol);
}
double cellDepth1G( const EclipseGrid& grid, size_t glob_idx) {
return grid.getCellDepth(glob_idx);
}
double cellDepth3( const EclipseGrid& grid, size_t i_idx, size_t j_idx, size_t k_idx) {
return grid.getCellDepth(i_idx, j_idx, k_idx);
}
py::array cellDepthAll( const EclipseGrid& grid)
{
std::vector<double> cellDepth;
std::array<int, 3> dims = grid.getNXYZ();
size_t nCells = dims[0]*dims[1]*dims[2];
cellDepth.reserve(nCells);
for (size_t n = 0; n < nCells; n++)
cellDepth.push_back(grid.getCellDepth(n));
return convert::numpy_array(cellDepth);
}
py::array cellDepthMask( const EclipseGrid& grid, std::vector<int>& mask)
{
std::array<int, 3> dims = grid.getNXYZ();
size_t nCells = dims[0]*dims[1]*dims[2];
if (nCells != mask.size())
throw std::logic_error("size of input mask doesn't match size of grid");
std::vector<double> cellDepth(nCells, 0.0);
for (size_t n = 0; n < nCells; n++)
if (mask[n]==1)
cellDepth[n] = grid.getCellDepth(n);
return convert::numpy_array(cellDepth);
}
}
void python::common::export_EclipseGrid(py::module& module) {
@@ -54,6 +123,12 @@ void python::common::export_EclipseGrid(py::module& module) {
.def( "getIJK", &getIJK )
.def( "getCellVolume", &cellVolume1G)
.def( "getCellVolume", &cellVolume3)
.def( "getCellVolume", &cellVolumeAll)
.def( "getCellVolume", &cellVolumeMask)
.def( "getCellDepth", &cellDepth1G)
.def( "getCellDepth", &cellDepth3)
.def( "getCellDepth", &cellDepthAll)
.def( "getCellDepth", &cellDepthMask)
;
}