24#ifndef OPM_ACTION_HANDLER_HPP
25#define OPM_ACTION_HANDLER_HPP
27#include <opm/simulators/utils/ParallelCommunication.hpp>
31#include <unordered_map>
41class BlackoilWellModelGeneric;
44struct SimulatorUpdate;
57 Action::State& actionState,
58 SummaryState& summaryState,
60 Parallel::Communication comm);
62 void applyActions(
int reportStep,
77 void applySimulatorUpdate(
int report_step,
78 const SimulatorUpdate& sim_update,
79 bool& commit_wellstate,
82 std::unordered_map<std::string, double>
83 fetchWellPI(
int reportStep,
84 const Action::ActionX& action,
85 const std::vector<std::string>& matching_wells)
const;
87 EclipseState& ecl_state_;
89 Action::State& actionState_;
90 SummaryState& summaryState_;
92 Parallel::Communication comm_;
Class handling Action support in simulator.
Definition ActionHandler.hpp:50
void evalUDQAssignments(const unsigned episodeIdx, UDQState &udq_state)
Evaluates UDQ assign statements.
Definition ActionHandler.cpp:256
std::function< void(bool)> TransFunc
Function handle to update transmissiblities.
Definition ActionHandler.hpp:53
Class for handling the blackoil well model.
Definition BlackoilWellModelGeneric.hpp:83
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition BlackoilPhases.hpp:27