opm-common
Loading...
Searching...
No Matches
HandlerContext.hpp
1/*
2 Copyright 2013 Statoil 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#ifndef HANDLER_CONTEXT_HPP
20#define HANDLER_CONTEXT_HPP
21
22#include <opm/common/OpmLog/KeywordLocation.hpp>
23
24#include <opm/input/eclipse/Schedule/Action/ActionResult.hpp>
25
26#include <cstddef>
27#include <cstdint>
28#include <optional>
29#include <set>
30#include <string>
31#include <unordered_map>
32#include <vector>
33
34namespace Opm::Action {
35 class WGNames;
36}
37
38namespace Opm {
39class DeckKeyword;
40class DeckRecord;
41class ErrorGuard;
42class ParseContext;
43class Schedule;
44class ScheduleBlock;
45class ScheduleGrid;
46class ScheduleState;
47struct ScheduleStatic;
48struct SimulatorUpdate;
49enum class WellStatus : std::uint8_t;
50class WelSegsSet;
51}
52
53namespace Opm {
55{
56public:
73 const ScheduleBlock& block_,
74 const DeckKeyword& keyword_,
75 const ScheduleGrid& grid_,
76 const std::size_t currentStep_,
78 bool action_mode_,
79 const ParseContext& parseContext_,
80 ErrorGuard& errors_,
81 SimulatorUpdate* sim_update_,
82 const std::unordered_map<std::string, double>* target_wellpi_,
83 std::unordered_map<std::string, double>& wpimult_global_factor_,
84 WelSegsSet* welsegs_wells_,
85 std::set<std::string>* compsegs_wells_,
86 std::set<std::string>* comptraj_wells_)
87 : block(block_)
88 , keyword(keyword_)
89 , currentStep(currentStep_)
90 , matches(matches_)
91 , action_mode(action_mode_)
92 , parseContext(parseContext_)
93 , errors(errors_)
94 , wpimult_global_factor(wpimult_global_factor_)
95 , grid(grid_)
96 , target_wellpi(target_wellpi_)
97 , welsegs_wells(welsegs_wells_)
98 , compsegs_wells(compsegs_wells_)
99 , comptraj_wells(comptraj_wells_)
100 , sim_update(sim_update_)
101 , schedule_(schedule)
102 {}
103
105 void affected_well(const std::string& well_name);
106
108 void welpi_well(const std::string& well_name);
109
111 void record_tran_change();
112
115
118
120 const ScheduleStatic& static_schedule() const;
121
123 void welsegs_handled(const std::string& well_name);
124
126 void compsegs_handled(const std::string& well_name);
127
129 void comptraj_handled(const std::string& well_name);
130
132 void setExitCode(int code);
133
135 bool updateWellStatus(const std::string& well,
136 WellStatus status,
137 const std::optional<KeywordLocation>& location = {});
138
140 WellStatus getWellStatus(const std::string& well) const;
141
143 void addGroup(const std::string& groupName);
145 void addGroupToGroup(const std::string& parent_group,
146 const std::string& child_group);
147
149 void welspecsCreateNewWell(const DeckRecord& record,
150 const std::string& wellName,
151 const std::string& groupName);
152
154 void welspecsUpdateExistingWells(const DeckRecord& record,
155 const std::vector<std::string>& wellNames,
156 const std::string& groupName);
157
159 double getWellPI(const std::string& well_name) const;
160
162 double elapsed_seconds() const;
163
165 void invalidNamePattern(const std::string& namePattern) const;
166
168 const Action::WGNames& action_wgnames() const;
169
171 std::vector<std::string> groupNames(const std::string& pattern) const;
172
183 bool hasWell(const std::string& pattern) const;
184
190 bool hasGroup(const std::string& pattern) const;
191
194 std::vector<std::string>
195 wellNames(const std::string& pattern) const;
196
200 std::vector<std::string>
201 wellNames(const std::string& pattern, bool allowEmpty) const;
202
203 const ScheduleBlock& block;
204 const DeckKeyword& keyword;
205 const std::size_t currentStep;
206 const Action::Result::MatchingEntities& matches;
207 const bool action_mode;
208 const ParseContext& parseContext;
209 ErrorGuard& errors;
210 std::unordered_map<std::string, double>& wpimult_global_factor;
211 const ScheduleGrid& grid;
212
213private:
214 const std::unordered_map<std::string, double>* target_wellpi{nullptr};
215 WelSegsSet* welsegs_wells{nullptr};
216 std::set<std::string>* compsegs_wells{nullptr};
217 std::set<std::string>* comptraj_wells{nullptr};
218 SimulatorUpdate* sim_update{nullptr};
219 Schedule& schedule_;
220};
221
222} // end namespace Opm
223
224#endif // HANDLER_CONTEXT_HPP
Container of matching entities.
Definition ActionResult.hpp:173
Definition WGNames.hpp:29
Definition DeckKeyword.hpp:36
Definition DeckRecord.hpp:32
Definition ErrorGuard.hpp:30
void compsegs_handled(const std::string &well_name)
Mark that the well occured in a COMPSEGS keyword.
Definition HandlerContext.cpp:79
void addGroupToGroup(const std::string &parent_group, const std::string &child_group)
Adds a group to a group.
Definition HandlerContext.cpp:201
void setExitCode(int code)
Set exit code.
Definition HandlerContext.cpp:98
void welspecsCreateNewWell(const DeckRecord &record, const std::string &wellName, const std::string &groupName)
Create a new Well from a WELSPECS record.
Definition HandlerContext.cpp:207
ScheduleState & state()
Returns a reference to current state.
Definition HandlerContext.cpp:93
void addGroup(const std::string &groupName)
Adds a group to the schedule.
Definition HandlerContext.cpp:196
void comptraj_handled(const std::string &well_name)
Mark that the well occured in a COMPTRAJ keyword.
Definition HandlerContext.cpp:86
void record_well_structure_change()
Mark that well structure has changed.
Definition HandlerContext.cpp:65
void welpi_well(const std::string &well_name)
Mark that a well is affected by WELPI.
Definition HandlerContext.cpp:51
std::vector< std::string > groupNames(const std::string &pattern) const
Obtain well group names from a pattern.
Definition HandlerContext.cpp:178
void welspecsUpdateExistingWells(const DeckRecord &record, const std::vector< std::string > &wellNames, const std::string &groupName)
Update one or more existing wells from a WELSPECS record.
Definition HandlerContext.cpp:236
bool hasWell(const std::string &pattern) const
Whether or not any existing well matches a name pattern.
Definition HandlerContext.cpp:162
WellStatus getWellStatus(const std::string &well) const
Get status of a well.
Definition HandlerContext.cpp:111
bool updateWellStatus(const std::string &well, WellStatus status, const std::optional< KeywordLocation > &location={})
Update status of a well.
Definition HandlerContext.cpp:103
const Action::WGNames & action_wgnames() const
Obtain action well group names.
Definition HandlerContext.cpp:157
void record_tran_change()
Mark that transmissibilities must be recalculated.
Definition HandlerContext.cpp:58
double getWellPI(const std::string &well_name) const
Obtain PI for a well.
Definition HandlerContext.cpp:121
HandlerContext(Schedule &schedule, const ScheduleBlock &block_, const DeckKeyword &keyword_, const ScheduleGrid &grid_, const std::size_t currentStep_, const Action::Result::MatchingEntities &matches_, bool action_mode_, const ParseContext &parseContext_, ErrorGuard &errors_, SimulatorUpdate *sim_update_, const std::unordered_map< std::string, double > *target_wellpi_, std::unordered_map< std::string, double > &wpimult_global_factor_, WelSegsSet *welsegs_wells_, std::set< std::string > *compsegs_wells_, std::set< std::string > *comptraj_wells_)
Definition HandlerContext.hpp:72
bool hasGroup(const std::string &pattern) const
Whether or not any existing group matches a name pattern.
Definition HandlerContext.cpp:172
void invalidNamePattern(const std::string &namePattern) const
Adds parse error for an invalid name pattern.
Definition HandlerContext.cpp:140
std::vector< std::string > wellNames(const std::string &pattern) const
Obtain well names from a pattern.
Definition HandlerContext.cpp:190
double elapsed_seconds() const
Returns elapsed time since simulation start in seconds.
Definition HandlerContext.cpp:135
void welsegs_handled(const std::string &well_name)
Mark that the well occured in a WELSEGS keyword.
Definition HandlerContext.cpp:72
const ScheduleStatic & static_schedule() const
Returns a const-ref to the static schedule.
Definition HandlerContext.cpp:116
void affected_well(const std::string &well_name)
Mark that a well has changed.
Definition HandlerContext.cpp:44
Control parser behaviour in failure conditions.
Definition ParseContext.hpp:115
Definition ScheduleBlock.hpp:52
Collection of intersected cells and associate properties for all simulation grids,...
Definition ScheduleGrid.hpp:50
Definition ScheduleState.hpp:106
Definition Schedule.hpp:101
Definition WelSegsSet.hpp:34
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30
Initial state of Schedule object created from information in SOLUTION section.
Definition ScheduleStatic.hpp:48
This struct is used to communicate back from the Schedule::applyAction() what needs to be updated in ...
Definition SimulatorUpdate.hpp:33