ActionHandler.hpp
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
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 2 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 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
23
24#ifndef OPM_ACTION_HANDLER_HPP
25#define OPM_ACTION_HANDLER_HPP
26
28
29#include <functional>
30#include <string>
31#include <unordered_map>
32#include <vector>
33
34namespace Opm::Action {
35 class State;
36} // namespace Opm::Action
37
38namespace Opm {
39template<class Scalar> class BlackoilWellModelGeneric;
40class EclipseState;
41class Schedule;
42struct SimulatorUpdate;
43class SummaryState;
44class UDQState;
45} // namespace Opm
46
47namespace Opm {
48
50template<class Scalar>
52{
53public:
55 using TransFunc = std::function<void(bool)>;
56
73 ActionHandler(EclipseState& ecl_state,
74 Schedule& schedule,
75 Action::State& actionState,
76 SummaryState& summaryState,
79
89 void applyActions(int reportStep,
90 double sim_time,
91 const TransFunc& updateTrans);
92
99 void evalUDQAssignments(const unsigned episodeIdx,
100 UDQState& udq_state);
101
122 void applySimulatorUpdate(int report_step,
123 const SimulatorUpdate& sim_update,
124 const TransFunc& updateTrans,
125 bool& commit_wellstate);
126
127private:
129 EclipseState& ecl_state_;
130
132 Schedule& schedule_;
133
135 Action::State& actionState_;
136
138 SummaryState& summaryState_;
139
142
145};
146
147} // namespace Opm
148
149#endif // OPM_ACTION_HANDLER_HPP
Class handling Action support in simulator.
Definition: ActionHandler.hpp:52
ActionHandler(EclipseState &ecl_state, Schedule &schedule, Action::State &actionState, SummaryState &summaryState, BlackoilWellModelGeneric< Scalar > &wellModel, Parallel::Communication comm)
void evalUDQAssignments(const unsigned episodeIdx, UDQState &udq_state)
Evaluates UDQ assign statements.
std::function< void(bool)> TransFunc
Function handle to update transmissiblities.
Definition: ActionHandler.hpp:55
void applySimulatorUpdate(int report_step, const SimulatorUpdate &sim_update, const TransFunc &updateTrans, bool &commit_wellstate)
void applyActions(int reportStep, double sim_time, const TransFunc &updateTrans)
Class for handling the blackoil well model.
Definition: BlackoilWellModelGeneric.hpp:94
Definition: ActionHandler.hpp:34
Dune::Communication< MPIComm > Communication
Definition: ParallelCommunication.hpp:30
Definition: blackoilboundaryratevector.hh:39