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:
72 const ScheduleBlock& block_,
73 const DeckKeyword& keyword_,
74 const ScheduleGrid& grid_,
75 const std::size_t currentStep_,
77 bool action_mode_,
78 const ParseContext& parseContext_,
79 ErrorGuard& errors_,
80 SimulatorUpdate* sim_update_,
81 const std::unordered_map<std::string, double>* target_wellpi_,
82 std::unordered_map<std::string, double>& wpimult_global_factor_,
83 WelSegsSet* welsegs_wells_,
84 std::set<std::string>* compsegs_wells_)
85 : block(block_)
86 , keyword(keyword_)
87 , currentStep(currentStep_)
88 , matches(matches_)
89 , action_mode(action_mode_)
90 , parseContext(parseContext_)
91 , errors(errors_)
92 , wpimult_global_factor(wpimult_global_factor_)
93 , grid(grid_)
94 , target_wellpi(target_wellpi_)
95 , welsegs_wells(welsegs_wells_)
96 , compsegs_wells(compsegs_wells_)
97 , sim_update(sim_update_)
98 , schedule_(schedule)
99 {}
100
102 void affected_well(const std::string& well_name);
103
105 void welpi_well(const std::string& well_name);
106
108 void record_tran_change();
109
112
115
117 const ScheduleStatic& static_schedule() const;
118
120 void welsegs_handled(const std::string& well_name);
121
123 void compsegs_handled(const std::string& well_name);
124
126 void setExitCode(int code);
127
129 bool updateWellStatus(const std::string& well,
130 WellStatus status,
131 const std::optional<KeywordLocation>& location = {});
132
134 WellStatus getWellStatus(const std::string& well) const;
135
137 void addGroup(const std::string& groupName);
139 void addGroupToGroup(const std::string& parent_group,
140 const std::string& child_group);
141
143 void welspecsCreateNewWell(const DeckRecord& record,
144 const std::string& wellName,
145 const std::string& groupName);
146
148 void welspecsUpdateExistingWells(const DeckRecord& record,
149 const std::vector<std::string>& wellNames,
150 const std::string& groupName);
151
153 double getWellPI(const std::string& well_name) const;
154
156 double elapsed_seconds() const;
157
159 void invalidNamePattern(const std::string& namePattern) const;
160
162 const Action::WGNames& action_wgnames() const;
163
165 std::vector<std::string> groupNames(const std::string& pattern) const;
166
177 bool hasWell(const std::string& pattern) const;
178
184 bool hasGroup(const std::string& pattern) const;
185
188 std::vector<std::string>
189 wellNames(const std::string& pattern) const;
190
194 std::vector<std::string>
195 wellNames(const std::string& pattern, bool allowEmpty) const;
196
197 const ScheduleBlock& block;
198 const DeckKeyword& keyword;
199 const std::size_t currentStep;
200 const Action::Result::MatchingEntities& matches;
201 const bool action_mode;
202 const ParseContext& parseContext;
203 ErrorGuard& errors;
204 std::unordered_map<std::string, double>& wpimult_global_factor;
205 const ScheduleGrid& grid;
206
207private:
208 const std::unordered_map<std::string, double>* target_wellpi{nullptr};
209 WelSegsSet* welsegs_wells{nullptr};
210 std::set<std::string>* compsegs_wells{nullptr};
211 SimulatorUpdate* sim_update{nullptr};
212 Schedule& schedule_;
213};
214
215} // end namespace Opm
216
217#endif // HANDLER_CONTEXT_HPP
Container of matching entities.
Definition ActionResult.hpp:174
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:194
void setExitCode(int code)
Set exit code.
Definition HandlerContext.cpp:91
void welspecsCreateNewWell(const DeckRecord &record, const std::string &wellName, const std::string &groupName)
Create a new Well from a WELSPECS record.
Definition HandlerContext.cpp:200
ScheduleState & state()
Returns a reference to current state.
Definition HandlerContext.cpp:86
void addGroup(const std::string &groupName)
Adds a group to the schedule.
Definition HandlerContext.cpp:189
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:171
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:229
bool hasWell(const std::string &pattern) const
Whether or not any existing well matches a name pattern.
Definition HandlerContext.cpp:155
WellStatus getWellStatus(const std::string &well) const
Get status of a well.
Definition HandlerContext.cpp:104
bool updateWellStatus(const std::string &well, WellStatus status, const std::optional< KeywordLocation > &location={})
Update status of a well.
Definition HandlerContext.cpp:96
const Action::WGNames & action_wgnames() const
Obtain action well group names.
Definition HandlerContext.cpp:150
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:114
bool hasGroup(const std::string &pattern) const
Whether or not any existing group matches a name pattern.
Definition HandlerContext.cpp:165
void invalidNamePattern(const std::string &namePattern) const
Adds parse error for an invalid name pattern.
Definition HandlerContext.cpp:133
std::vector< std::string > wellNames(const std::string &pattern) const
Obtain well names from a pattern.
Definition HandlerContext.cpp:183
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_)
Definition HandlerContext.hpp:71
double elapsed_seconds() const
Returns elapsed time since simulation start in seconds.
Definition HandlerContext.cpp:128
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:109
void affected_well(const std::string &well_name)
Mark that a well has changed.
Definition HandlerContext.cpp:44
Definition ParseContext.hpp:84
Definition ScheduleBlock.hpp:52
Collection of intersected cells and associate properties for all simulation grids,...
Definition ScheduleGrid.hpp:50
Definition ScheduleState.hpp:101
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