My Project
Loading...
Searching...
No Matches
WellState.hpp
1/*
2 Copyright 2014 SINTEF ICT, Applied Mathematics.
3 Copyright 2017 IRIS AS
4
5 This file is part of the Open Porous Media project (OPM).
6
7 OPM is free software: you can redistribute it and/or modify
8 it under the terms of the GNU General Public License as published by
9 the Free Software Foundation, either version 3 of the License, or
10 (at your option) any later version.
11
12 OPM is distributed in the hope that it will be useful,
13 but WITHOUT ANY WARRANTY; without even the implied warranty of
14 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 GNU General Public License for more details.
16
17 You should have received a copy of the GNU General Public License
18 along with OPM. If not, see <http://www.gnu.org/licenses/>.
19*/
20
21#ifndef OPM_WELLSTATEFULLYIMPLICITBLACKOIL_HEADER_INCLUDED
22#define OPM_WELLSTATEFULLYIMPLICITBLACKOIL_HEADER_INCLUDED
23
24#include <opm/common/ErrorMacros.hpp>
25
26#include <opm/core/props/BlackoilPhases.hpp>
27
28#include <opm/simulators/wells/ALQState.hpp>
29#include <opm/simulators/wells/GlobalWellInfo.hpp>
30#include <opm/simulators/wells/PerfData.hpp>
31#include <opm/simulators/wells/PerforationData.hpp>
32#include <opm/simulators/wells/SegmentState.hpp>
33#include <opm/simulators/wells/SingleWellState.hpp>
34#include <opm/simulators/wells/WellContainer.hpp>
35
36#include <opm/output/data/Wells.hpp>
37
38#include <opm/input/eclipse/Schedule/Events.hpp>
39
40#include <dune/common/version.hh>
41#include <dune/common/parallel/mpihelper.hh>
42
43#include <functional>
44#include <map>
45#include <optional>
46#include <string>
47#include <utility>
48#include <vector>
49
50namespace Opm
51{
52
53class ParallelWellInfo;
54class Schedule;
55enum class WellStatus;
56
60{
61public:
62 static const uint64_t event_mask = ScheduleEvents::WELL_STATUS_CHANGE + ScheduleEvents::PRODUCTION_UPDATE + ScheduleEvents::INJECTION_UPDATE;
63 // TODO: same definition with WellInterface, eventually they should go to a common header file.
64 static const int Water = BlackoilPhases::Aqua;
65 static const int Oil = BlackoilPhases::Liquid;
66 static const int Gas = BlackoilPhases::Vapour;
67
68 // Only usable for testing purposes
69 explicit WellState(const ParallelWellInfo& pinfo);
70
71 explicit WellState(const PhaseUsage& pu)
72 : phase_usage_(pu)
73 {}
74
75 static WellState serializationTestObject(const ParallelWellInfo& pinfo);
76
77 std::size_t size() const {
78 return this->wells_.size();
79 }
80
81 std::vector<std::string> wells() const {
82 return this->wells_.wells();
83 }
84
85
86 int numWells() const
87 {
88 return this->size();
89 }
90
91 const ParallelWellInfo& parallelWellInfo(std::size_t well_index) const;
92
93
94
98 void init(const std::vector<double>& cellPressures,
99 const Schedule& schedule,
100 const std::vector<Well>& wells_ecl,
101 const std::vector<std::reference_wrapper<ParallelWellInfo>>& parallel_well_info,
102 const int report_step,
103 const WellState* prevState,
104 const std::vector<std::vector<PerforationData>>& well_perf_data,
106
107 void resize(const std::vector<Well>& wells_ecl,
108 const std::vector<std::reference_wrapper<ParallelWellInfo>>& parallel_well_info,
109 const Schedule& schedule,
110 const bool handle_ms_well,
111 const std::size_t numCells,
112 const std::vector<std::vector<PerforationData>>& well_perf_data,
114
115 void setCurrentWellRates(const std::string& wellName, const std::vector<double>& new_rates ) {
116 auto& [owner, rates] = this->well_rates.at(wellName);
117 if (owner)
118 rates = new_rates;
119 }
120
121 const std::vector<double>& currentWellRates(const std::string& wellName) const;
122
123 bool hasWellRates(const std::string& wellName) const {
124 return this->well_rates.find(wellName) != this->well_rates.end();
125 }
126
127 void clearWellRates()
128 {
129 this->well_rates.clear();
130 }
131
132 template<class Communication>
133 void gatherVectorsOnRoot(const std::vector< data::Connection >& from_connections,
134 std::vector< data::Connection >& to_connections,
135 const Communication& comm) const;
136
137 data::Wells
138 report(const int* globalCellIdxMap,
139 const std::function<bool(const int)>& wasDynamicallyClosed) const;
140
141 void reportConnections(std::vector<data::Connection>& connections, const PhaseUsage &pu,
142 std::size_t well_index,
143 const int* globalCellIdxMap) const;
144
146 void initWellStateMSWell(const std::vector<Well>& wells_ecl,
148
149 static void calculateSegmentRates(const std::vector<std::vector<int>>& segment_inlets, const std::vector<std::vector<int>>&segment_perforations,
150 const std::vector<double>& perforation_rates, const int np, const int segment, std::vector<double>& segment_rates);
151
152
153 template<class Comm>
154 void communicateGroupRates(const Comm& comm);
155
156 template<class Comm>
157 void updateGlobalIsGrup(const Comm& comm);
158
159 bool isInjectionGrup(const std::string& name) const {
160 return this->global_well_info.value().in_injecting_group(name);
161 }
162
163 bool isProductionGrup(const std::string& name) const {
164 return this->global_well_info.value().in_producing_group(name);
165 }
166
167 double getALQ( const std::string& name) const
168 {
169 return this->alq_state.get(name);
170 }
171
172 void setALQ( const std::string& name, double value)
173 {
174 this->alq_state.set(name, value);
175 }
176
177 int gliftGetDebugCounter() {
178 return this->alq_state.get_debug_counter();
179 }
180
181 void gliftSetDebugCounter(int value) {
182 return this->alq_state.set_debug_counter(value);
183 }
184
185 int gliftUpdateDebugCounter() {
186 return this->alq_state.update_debug_counter();
187 }
188
189 bool gliftCheckAlqOscillation(const std::string &name) const {
190 return this->alq_state.oscillation(name);
191 }
192
193 int gliftGetAlqDecreaseCount(const std::string &name) {
194 return this->alq_state.get_decrement_count(name);
195 }
196
197 int gliftGetAlqIncreaseCount(const std::string &name) {
198 return this->alq_state.get_increment_count(name);
199 }
200
201 void gliftUpdateAlqIncreaseCount(const std::string &name, bool increase) {
202 this->alq_state.update_count(name, increase);
203 }
204
205 void gliftTimeStepInit() {
206 this->alq_state.reset_count();
207 }
208
209 int wellNameToGlobalIdx(const std::string &name) {
210 return this->global_well_info.value().well_index(name);
211 }
212
213 std::string globalIdxToWellName(const int index) {
214 return this->global_well_info.value().well_name(index);
215 }
216
217 bool wellIsOwned(std::size_t well_index,
218 const std::string& wellName) const;
219
220 bool wellIsOwned(const std::string& wellName) const;
221
222 void updateStatus(int well_index, WellStatus status);
223
224 void openWell(int well_index);
225 void shutWell(int well_index);
226 void stopWell(int well_index);
227
229 int numPhases() const
230 {
231 return this->phase_usage_.num_phases;
232 }
233
234 const PhaseUsage& phaseUsage() const {
235 return this->phase_usage_;
236 }
237
239 std::vector<double>& wellRates(std::size_t well_index) { return this->wells_[well_index].surface_rates; }
240 const std::vector<double>& wellRates(std::size_t well_index) const { return this->wells_[well_index].surface_rates; }
241
242 const std::string& name(std::size_t well_index) const {
243 return this->wells_.well_name(well_index);
244 }
245
246 std::optional<std::size_t> index(const std::string& well_name) const {
247 return this->wells_.well_index(well_name);
248 }
249
250 const SingleWellState& operator[](std::size_t well_index) const {
251 return this->wells_[well_index];
252 }
253
254 const SingleWellState& operator[](const std::string& well_name) const {
255 return this->wells_[well_name];
256 }
257
258 SingleWellState& operator[](std::size_t well_index) {
259 return this->wells_[well_index];
260 }
261
262 SingleWellState& operator[](const std::string& well_name) {
263 return this->wells_[well_name];
264 }
265
266 const SingleWellState& well(std::size_t well_index) const {
267 return this->operator[](well_index);
268 }
269
270 const SingleWellState& well(const std::string& well_name) const {
271 return this->operator[](well_name);
272 }
273
274 SingleWellState& well(std::size_t well_index) {
275 return this->operator[](well_index);
276 }
277
278 SingleWellState& well(const std::string& well_name) {
279 return this->operator[](well_name);
280 }
281
282 bool has(const std::string& well_name) const {
283 return this->wells_.has(well_name);
284 }
285
286 bool operator==(const WellState&) const;
287
288 template<class Serializer>
289 void serializeOp(Serializer& serializer)
290 {
291 serializer(alq_state);
292 serializer(well_rates);
293 if (serializer.isSerializing()) {
294 serializer(wells_.size());
295 } else {
296 std::size_t size = 0;
297 serializer(size);
298 if (size != wells_.size()) {
299 OPM_THROW(std::runtime_error, "Error deserializing WellState: size mismatch");
300 }
301 }
302 for (auto& w : wells_) {
303 serializer(w);
304 }
305 }
306
307private:
308 PhaseUsage phase_usage_;
309
310 // The wells_ variable is essentially a map of all the wells on the current
311 // process. Observe that since a well can be split over several processes a
312 // well might appear in the WellContainer on different processes.
313 WellContainer<SingleWellState> wells_;
314
315 // The members alq_state, global_well_info and well_rates are map like
316 // structures which will have entries for *all* the wells in the system.
317
318 // Use of std::optional<> here is a technical crutch, the
319 // WellStateFullyImplicitBlackoil class should be default constructible,
320 // whereas the GlobalWellInfo is not.
321 std::optional<GlobalWellInfo> global_well_info;
322 ALQState alq_state;
323
324 // The well_rates variable is defined for all wells on all processors. The
325 // bool in the value pair is whether the current process owns the well or
326 // not.
327 std::map<std::string, std::pair<bool, std::vector<double>>> well_rates;
328
329 data::Segment
330 reportSegmentResults(const int well_id,
331 const int seg_ix,
332 const int seg_no) const;
333
334 // If the ALQ has changed since the previous report step,
335 // reset current_alq and update default_alq. ALQ is used for
336 // constant lift gas injection and for gas lift optimization
337 // (THP controlled wells).
338
339 void updateWellsDefaultALQ(const std::vector<Well>& wells_ecl);
340
346 void base_init(const std::vector<double>& cellPressures,
347 const std::vector<Well>& wells_ecl,
348 const std::vector<std::reference_wrapper<ParallelWellInfo>>& parallel_well_info,
349 const std::vector<std::vector<PerforationData>>& well_perf_data,
350 const SummaryState& summary_state);
351
352 void initSingleWell(const std::vector<double>& cellPressures,
353 const Well& well,
354 const std::vector<PerforationData>& well_perf_data,
355 const ParallelWellInfo& well_info,
356 const SummaryState& summary_state);
357
358 void initSingleProducer(const Well& well,
359 const ParallelWellInfo& well_info,
360 double pressure_first_connection,
361 const std::vector<PerforationData>& well_perf_data,
362 const SummaryState& summary_state);
363
364 void initSingleInjector(const Well& well,
365 const ParallelWellInfo& well_info,
366 double pressure_first_connection,
367 const std::vector<PerforationData>& well_perf_data,
368 const SummaryState& summary_state);
369
370};
371
372} // namespace Opm
373
374
375#endif // OPM_WELLSTATEFULLYIMPLICITBLACKOIL_HEADER_INCLUDED
Definition AquiferInterface.hpp:35
Class encapsulating some information about parallel wells.
Definition ParallelWellInfo.hpp:184
The state of a set of wells, tailored for use by the fully implicit blackoil simulator.
Definition WellState.hpp:60
void init(const std::vector< double > &cellPressures, const Schedule &schedule, const std::vector< Well > &wells_ecl, const std::vector< std::reference_wrapper< ParallelWellInfo > > &parallel_well_info, const int report_step, const WellState *prevState, const std::vector< std::vector< PerforationData > > &well_perf_data, const SummaryState &summary_state)
Allocate and initialize if wells is non-null.
Definition WellState.cpp:241
void initWellStateMSWell(const std::vector< Well > &wells_ecl, const WellState *prev_well_state)
init the MS well related.
Definition WellState.cpp:654
int numPhases() const
The number of phases present.
Definition WellState.hpp:229
std::vector< double > & wellRates(std::size_t well_index)
One rate per well and phase.
Definition WellState.hpp:239
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition BlackoilPhases.hpp:27
Definition BlackoilPhases.hpp:46