My Project
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 <cstddef>
25#include <optional>
26#include <set>
27#include <string>
28#include <unordered_map>
29#include <vector>
30
31namespace Opm {
32
33namespace Action { class WGNames; }
34class DeckKeyword;
35class DeckRecord;
36class ErrorGuard;
37class ParseContext;
38class Schedule;
39class ScheduleBlock;
40class ScheduleGrid;
41class ScheduleState;
42struct ScheduleStatic;
43struct SimulatorUpdate;
44enum class WellStatus;
45class WelSegsSet;
46
48{
49public:
53 const ScheduleBlock& block_,
54 const DeckKeyword& keyword_,
55 const ScheduleGrid& grid_,
56 const std::size_t currentStep_,
57 const std::vector<std::string>& matching_wells_,
58 bool actionx_mode_,
59 const ParseContext& parseContext_,
60 ErrorGuard& errors_,
61 SimulatorUpdate* sim_update_,
62 const std::unordered_map<std::string, double>* target_wellpi_,
63 std::unordered_map<std::string, double>& wpimult_global_factor_,
64 WelSegsSet* welsegs_wells_,
65 std::set<std::string>* compsegs_wells_)
66 : block(block_)
67 , keyword(keyword_)
68 , currentStep(currentStep_)
69 , matching_wells(matching_wells_)
70 , actionx_mode(actionx_mode_)
71 , parseContext(parseContext_)
72 , errors(errors_)
73 , wpimult_global_factor(wpimult_global_factor_)
74 , grid(grid_)
75 , target_wellpi(target_wellpi_)
76 , welsegs_wells(welsegs_wells_)
77 , compsegs_wells(compsegs_wells_)
78 , sim_update(sim_update_)
79 , schedule_(schedule)
80 {}
81
83 void affected_well(const std::string& well_name);
84
86 void record_tran_change();
87
90
93
95 const ScheduleStatic& static_schedule() const;
96
98 void welsegs_handled(const std::string& well_name);
99
101 void compsegs_handled(const std::string& well_name);
102
104 void setExitCode(int code);
105
107 bool updateWellStatus(const std::string& well,
108 WellStatus status,
109 std::optional<KeywordLocation> location = {});
110
112 WellStatus getWellStatus(const std::string& well) const;
113
115 void addGroup(const std::string& groupName);
117 void addGroupToGroup(const std::string& parent_group,
118 const std::string& child_group);
119
121 void welspecsCreateNewWell(const DeckRecord& record,
122 const std::string& wellName,
123 const std::string& groupName);
124
126 void welspecsUpdateExistingWells(const DeckRecord& record,
127 const std::vector<std::string>& wellNames,
128 const std::string& groupName);
129
131 double getWellPI(const std::string& well_name) const;
132
134 double elapsed_seconds() const;
135
137 void invalidNamePattern(const std::string& namePattern) const;
138
140 const Action::WGNames& action_wgnames() const;
141
143 std::vector<std::string> groupNames(const std::string& pattern) const;
144
147 std::vector<std::string>
148 wellNames(const std::string& pattern) const;
149
153 std::vector<std::string>
154 wellNames(const std::string& pattern, bool allowEmpty) const;
155
156 const ScheduleBlock& block;
157 const DeckKeyword& keyword;
158 const std::size_t currentStep;
159 const std::vector<std::string>& matching_wells;
160 const bool actionx_mode;
161 const ParseContext& parseContext;
162 ErrorGuard& errors;
163 std::unordered_map<std::string, double>& wpimult_global_factor;
164 const ScheduleGrid& grid;
165
166private:
167 const std::unordered_map<std::string, double>* target_wellpi{nullptr};
168 WelSegsSet* welsegs_wells{nullptr};
169 std::set<std::string>* compsegs_wells{nullptr};
170 SimulatorUpdate* sim_update{nullptr};
171 Schedule& schedule_;
172};
173
174} // end namespace Opm
175
176#endif // HANDLER_CONTEXT_HPP
Definition DeckKeyword.hpp:36
Definition ErrorGuard.hpp:29
Definition HandlerContext.hpp:48
void compsegs_handled(const std::string &well_name)
Mark that the well occured in a COMPSEGS keyword.
Definition HandlerContext.cpp:71
void addGroupToGroup(const std::string &parent_group, const std::string &child_group)
Adds a group to a group.
Definition HandlerContext.cpp:173
HandlerContext(Schedule &schedule, const ScheduleBlock &block_, const DeckKeyword &keyword_, const ScheduleGrid &grid_, const std::size_t currentStep_, const std::vector< std::string > &matching_wells_, bool actionx_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:52
void setExitCode(int code)
Set exit code.
Definition HandlerContext.cpp:83
void welspecsCreateNewWell(const DeckRecord &record, const std::string &wellName, const std::string &groupName)
Create a new Well from a WELSPECS record.
Definition HandlerContext.cpp:179
ScheduleState & state()
Returns a reference to current state.
Definition HandlerContext.cpp:78
void addGroup(const std::string &groupName)
Adds a group to the schedule.
Definition HandlerContext.cpp:168
void record_well_structure_change()
Mark that well structure has changed.
Definition HandlerContext.cpp:57
std::vector< std::string > groupNames(const std::string &pattern) const
Obtain well group names from a pattern.
Definition HandlerContext.cpp:150
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:208
WellStatus getWellStatus(const std::string &well) const
Get status of a well.
Definition HandlerContext.cpp:96
const Action::WGNames & action_wgnames() const
Obtain action well group names.
Definition HandlerContext.cpp:144
void record_tran_change()
Mark that transmissibilities must be recalculated.
Definition HandlerContext.cpp:50
double getWellPI(const std::string &well_name) const
Obtain PI for a well.
Definition HandlerContext.cpp:106
void invalidNamePattern(const std::string &namePattern) const
Adds parse error for an invalid name pattern.
Definition HandlerContext.cpp:125
std::vector< std::string > wellNames(const std::string &pattern) const
Obtain well names from a pattern.
Definition HandlerContext.cpp:162
double elapsed_seconds() const
Returns elapsed time since simulation start in seconds.
Definition HandlerContext.cpp:120
void welsegs_handled(const std::string &well_name)
Mark that the well occured in a WELSEGS keyword.
Definition HandlerContext.cpp:64
bool updateWellStatus(const std::string &well, WellStatus status, std::optional< KeywordLocation > location={})
Update status of a well.
Definition HandlerContext.cpp:88
const ScheduleStatic & static_schedule() const
Returns a const-ref to the static schedule.
Definition HandlerContext.cpp:101
void affected_well(const std::string &well_name)
Mark that a well has changed.
Definition HandlerContext.cpp:43
Definition ParseContext.hpp:84
Definition ScheduleBlock.hpp:52
Definition ScheduleGrid.hpp:37
Definition ScheduleState.hpp:91
Definition Schedule.hpp:88
Definition WelSegsSet.hpp:34
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30
Definition ScheduleStatic.hpp:40
This struct is used to communicate back from the Schedule::applyAction() what needs to be updated in ...
Definition SimulatorUpdate.hpp:33