PyBlackOilSimulator.hpp
Go to the documentation of this file.
1/*
2 Copyright 2020 Equinor ASA.
3
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 3 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18*/
19
20#ifndef OPM_PY_BLACKOIL_SIMULATOR_HEADER_INCLUDED
21#define OPM_PY_BLACKOIL_SIMULATOR_HEADER_INCLUDED
22
23#include <python/simulators/PyMain.hpp>
24
27
33
34#include <opm/input/eclipse/EclipseState/EclipseState.hpp>
35#include <opm/input/eclipse/EclipseState/SummaryConfig/SummaryConfig.hpp>
36#include <opm/input/eclipse/Schedule/Schedule.hpp>
37
38#include <map>
39#include <memory>
40#include <string>
41#include <vector>
42
43namespace Opm::Pybind {
44
46{
47private:
50
51public:
52 PyBlackOilSimulator(const std::string& deckFilename,
53 const std::vector<std::string>& args);
55 std::shared_ptr<Opm::Deck> deck,
56 std::shared_ptr<Opm::EclipseState> state,
57 std::shared_ptr<Opm::Schedule> schedule,
58 std::shared_ptr<Opm::SummaryConfig> summary_config);
59 void advance(int report_step);
62 py::array_t<double> getFluidStateVariable(const std::string &name) const;
63 py::array_t<double> getCellVolumes();
64 double getDT();
65 py::array_t<double> getPorosity();
66 py::array_t<double> getPrimaryVariable(const std::string &variable) const;
67 py::array_t<int> getPrimaryVarMeaning(const std::string &variable) const;
68 std::map<std::string, int> getPrimaryVarMeaningMap(const std::string &variable) const;
69 int run();
71 py::array_t<double, py::array::c_style | py::array::forcecast> array);
73 const std::string &idx_name,
74 py::array_t<double,
75 py::array::c_style | py::array::forcecast> array);
76 void setupMpi(bool init_mpi, bool finalize_mpi);
77 int step();
79 int stepInit();
80
81private:
82 Opm::FlowMain<TypeTag>& getFlowMain() const;
83 PyFluidState<TypeTag>& getFluidState() const;
84 PyMaterialState<TypeTag>& getMaterialState() const;
85
86 const std::string deck_filename_;
87 bool has_run_init_ = false;
88 bool has_run_cleanup_ = false;
89 bool mpi_init_ = true;
90 bool mpi_finalize_ = true;
91 //bool debug_ = false;
92 // This *must* be declared before other pointers
93 // to simulator objects. This in order to deinitialize
94 // MPI at the correct time (ie after the other objects).
95 std::unique_ptr<Opm::PyMain> main_;
96
97 std::unique_ptr<Opm::FlowMain<TypeTag>> flow_main_;
98 Simulator* simulator_;
99 std::unique_ptr<PyFluidState<TypeTag>> fluid_state_;
100 std::unique_ptr<PyMaterialState<TypeTag>> material_state_;
101 std::shared_ptr<Opm::Deck> deck_;
102 std::shared_ptr<Opm::EclipseState> eclipse_state_;
103 std::shared_ptr<Opm::Schedule> schedule_;
104 std::shared_ptr<Opm::SummaryConfig> summary_config_;
105 std::vector<std::string> args_;
106};
107
108} // namespace Opm::Pybind
109#endif // OPM_PY_BLACKOIL_SIMULATOR_HEADER_INCLUDED
Definition: PyBlackOilSimulator.hpp:46
void setPorosity(py::array_t< double, py::array::c_style|py::array::forcecast > array)
PyBlackOilSimulator(const std::string &deckFilename, const std::vector< std::string > &args)
PyBlackOilSimulator(std::shared_ptr< Opm::Deck > deck, std::shared_ptr< Opm::EclipseState > state, std::shared_ptr< Opm::Schedule > schedule, std::shared_ptr< Opm::SummaryConfig > summary_config)
py::array_t< double > getPorosity()
void setPrimaryVariable(const std::string &idx_name, py::array_t< double, py::array::c_style|py::array::forcecast > array)
void setupMpi(bool init_mpi, bool finalize_mpi)
py::array_t< double > getPrimaryVariable(const std::string &variable) const
py::array_t< double > getCellVolumes()
py::array_t< double > getFluidStateVariable(const std::string &name) const
py::array_t< int > getPrimaryVarMeaning(const std::string &variable) const
void advance(int report_step)
std::map< std::string, int > getPrimaryVarMeaningMap(const std::string &variable) const
Definition: PyFluidState.hpp:35
Definition: PyMaterialState.hpp:34
Definition: Pybind11Exporter.hpp:11
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:233
This file provides the infrastructure to retrieve run-time parameters.
The Opm property system, traits with inheritance.
Definition: TTagFlowProblemTPFA.hpp:34