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 {
35
36namespace Action {
37class ActionX;
38class State;
39}
40
41template<class Scalar> class BlackoilWellModelGeneric;
42class EclipseState;
43class Schedule;
44struct SimulatorUpdate;
45class SummaryState;
46class UDQState;
47
49template<class Scalar>
51{
52public:
54 using TransFunc = std::function<void(bool)>;
55
56 ActionHandler(EclipseState& ecl_state,
57 Schedule& schedule,
58 Action::State& actionState,
59 SummaryState& summaryState,
62
63 void applyActions(int reportStep,
64 double sim_time,
65 const TransFunc& updateTrans);
66
68 void evalUDQAssignments(const unsigned episodeIdx,
69 UDQState& udq_state);
70
71 private:
72 /*
73 This function is run after applyAction has been completed in the Schedule
74 implementation. The sim_update argument should have members & flags for
75 the simulator properties which need to be updated. This functionality is
76 probably not complete.
77 */
78 void applySimulatorUpdate(int report_step,
79 const SimulatorUpdate& sim_update,
80 bool& commit_wellstate,
81 const TransFunc& updateTrans);
82
83 std::unordered_map<std::string, Scalar>
84 fetchWellPI(int reportStep,
85 const Action::ActionX& action,
86 const std::vector<std::string>& matching_wells) const;
87
88 EclipseState& ecl_state_;
89 Schedule& schedule_;
90 Action::State& actionState_;
91 SummaryState& summaryState_;
94};
95
96} // namespace Opm
97
98#endif // OPM_ACTION_HANDLER_HPP
Class handling Action support in simulator.
Definition: ActionHandler.hpp:51
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:54
void applyActions(int reportStep, double sim_time, const TransFunc &updateTrans)
Class for handling the blackoil well model.
Definition: BlackoilWellModelGeneric.hpp:83
Dune::Communication< MPIComm > Communication
Definition: ParallelCommunication.hpp:30
Definition: blackoilboundaryratevector.hh:37