Get cell volumes from Python

Adds a new method get_cell_volumes() to the opm.simulators Python module
that returns a python list of the cell volumes in the black oil
simulator.
This commit is contained in:
Håkon Hægland
2023-09-08 09:39:24 +02:00
parent 8bf45ea86b
commit e2f62644ae
7 changed files with 72 additions and 245 deletions

View File

@@ -44,6 +44,7 @@ public:
std::shared_ptr<Opm::Schedule> schedule,
std::shared_ptr<Opm::SummaryConfig> summary_config);
bool checkSimulationFinished();
py::array_t<double> getCellVolumes();
py::array_t<double> getPorosity();
int run();
void setPorosity(
@@ -56,18 +57,18 @@ public:
const Opm::FlowMainEbos<TypeTag>& getFlowMainEbos() const;
private:
const std::string deckFilename_;
bool hasRunInit_ = false;
bool hasRunCleanup_ = false;
const std::string deck_filename_;
bool has_run_init_ = false;
bool has_run_cleanup_ = false;
//bool debug_ = false;
// This *must* be declared before other pointers
// to simulator objects. This in order to deinitialize
// MPI at the correct time (ie after the other objects).
std::unique_ptr<Opm::Main> main_;
std::unique_ptr<Opm::FlowMainEbos<TypeTag>> mainEbos_;
Simulator *ebosSimulator_;
std::unique_ptr<PyMaterialState<TypeTag>> materialState_;
std::unique_ptr<Opm::FlowMainEbos<TypeTag>> main_ebos_;
Simulator *ebos_simulator_;
std::unique_ptr<PyMaterialState<TypeTag>> material_state_;
std::shared_ptr<Opm::Deck> deck_;
std::shared_ptr<Opm::EclipseState> eclipse_state_;
std::shared_ptr<Opm::Schedule> schedule_;

View File

@@ -42,14 +42,14 @@ namespace Opm::Pybind
using GridView = GetPropType<TypeTag, Opm::Properties::GridView>;
public:
PyMaterialState(Simulator *ebosSimulator)
: ebosSimulator_(ebosSimulator) { }
PyMaterialState(Simulator *ebos_simulator)
: ebos_simulator_(ebos_simulator) { }
std::unique_ptr<double []> getCellVolumes( std::size_t *size);
std::unique_ptr<double []> getPorosity( std::size_t *size);
void setPorosity(const double *poro, std::size_t size);
private:
Simulator *ebosSimulator_;
Simulator *ebos_simulator_;
};
}

View File

@@ -17,6 +17,8 @@
along with OPM. If not, see <http://www.gnu.org/licenses/>.
*/
#include <fmt/format.h>
namespace Opm::Pybind {
template <class TypeTag>
@@ -24,11 +26,11 @@ std::unique_ptr<double []>
PyMaterialState<TypeTag>::
getCellVolumes( std::size_t *size)
{
Model &model = ebosSimulator_->model();
Model &model = this->ebos_simulator_->model();
*size = model.numGridDof();
auto array = std::make_unique<double []>(*size);
for (unsigned dofIdx = 0; dofIdx < *size; ++dofIdx) {
array[dofIdx] = model.dofTotalVolume(dofIdx);
for (unsigned dof_idx = 0; dof_idx < *size; ++dof_idx) {
array[dof_idx] = model.dofTotalVolume(dof_idx);
}
return array;
}
@@ -38,12 +40,12 @@ std::unique_ptr<double []>
PyMaterialState<TypeTag>::
getPorosity( std::size_t *size)
{
Problem &problem = ebosSimulator_->problem();
Model &model = ebosSimulator_->model();
Problem &problem = this->ebos_simulator_->problem();
Model &model = this->ebos_simulator_->model();
*size = model.numGridDof();
auto array = std::make_unique<double []>(*size);
for (unsigned dofIdx = 0; dofIdx < *size; ++dofIdx) {
array[dofIdx] = problem.referencePorosity(dofIdx, /*timeIdx*/0);
for (unsigned dof_idx = 0; dof_idx < *size; ++dof_idx) {
array[dof_idx] = problem.referencePorosity(dof_idx, /*timeIdx*/0);
}
return array;
}
@@ -53,17 +55,17 @@ void
PyMaterialState<TypeTag>::
setPorosity(const double *poro, std::size_t size)
{
Problem &problem = ebosSimulator_->problem();
Model &model = ebosSimulator_->model();
Problem &problem = this->ebos_simulator_->problem();
Model &model = this->ebos_simulator_->model();
auto model_size = model.numGridDof();
if (model_size != size) {
std::ostringstream message;
message << "Cannot set porosity. Expected array of size: "
<< model_size << ", got array of size: " << size;
throw std::runtime_error(message.str());
const std::string msg = fmt::format(
"Cannot set porosity. Expected array of size: {}, got array of size: ",
model_size, size);
throw std::runtime_error(msg);
}
for (unsigned dofIdx = 0; dofIdx < size; ++dofIdx) {
problem.setPorosity(poro[dofIdx], dofIdx);
for (unsigned dof_idx = 0; dof_idx < size; ++dof_idx) {
problem.setPorosity(poro[dof_idx], dof_idx);
}
}
} //namespace Opm::Pybind

View File

@@ -1,64 +0,0 @@
/*
Copyright 2020 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/>.
*/
#ifndef OPM_SIMULATORS_HEADER_INCLUDED
#define OPM_SIMULATORS_HEADER_INCLUDED
#include <opm/simulators/flow/Main.hpp>
#include <opm/simulators/flow/FlowMainEbos.hpp>
#include <opm/models/utils/propertysystem.hh>
#include <pybind11/pybind11.h>
#include <pybind11/numpy.h>
namespace py = pybind11;
namespace Opm::Pybind {
class BlackOilSimulator
{
private:
using TypeTag = Opm::Properties::TTag::EclFlowProblemTpfa;
using Simulator = Opm::GetPropType<TypeTag, Opm::Properties::Simulator>;
public:
BlackOilSimulator( const std::string &deckFilename);
py::array_t<double> getPorosity();
int run();
void setPorosity(
py::array_t<double, py::array::c_style | py::array::forcecast> array);
int step();
int stepInit();
int stepCleanup();
private:
const std::string deckFilename_;
bool hasRunInit_ = false;
bool hasRunCleanup_ = false;
// This *must* be declared before other pointers
// to simulator objects. This in order to deinitialize
// MPI at the correct time (ie after the other objects).
std::unique_ptr<Opm::Main> main_;
std::unique_ptr<Opm::FlowMainEbos<TypeTag>> mainEbos_;
Simulator *ebosSimulator_;
std::unique_ptr<PyMaterialState<TypeTag>> materialState_;
};
} // namespace Opm::Pybind
#endif // OPM_SIMULATORS_HEADER_INCLUDED