opm-simulators
Loading...
Searching...
No Matches
GasLiftSingleWellGeneric.hpp
1/*
2 Copyright 2020 Equinor 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
20#ifndef OPM_GASLIFT_SINGLE_WELL_GENERIC_HEADER_INCLUDED
21#define OPM_GASLIFT_SINGLE_WELL_GENERIC_HEADER_INCLUDED
22
23#include <opm/input/eclipse/Schedule/Well/WellProductionControls.hpp>
24
25#include <opm/simulators/wells/GasLiftGroupInfo.hpp>
26#include <opm/simulators/wells/GasLiftCommon.hpp>
27
28#include <optional>
29#include <set>
30#include <stdexcept>
31#include <string>
32#include <tuple>
33#include <utility>
34
35namespace Opm {
36
37class DeferredLogger;
38class GasLiftWell;
39template<class Scalar> class GasLiftWellState;
40class Schedule;
41class SummaryState;
42template<typename Scalar, typename IndexTraits> class WellInterfaceGeneric;
43template<typename Scalar, typename IndexTraits> class WellState;
44template<class Scalar> class GroupState;
45
46template<typename Scalar, typename IndexTraits>
47class GasLiftSingleWellGeneric : public GasLiftCommon<Scalar, IndexTraits>
48{
49protected:
50 static constexpr int Water = IndexTraits::waterPhaseIdx;
51 static constexpr int Oil = IndexTraits::oilPhaseIdx;
52 static constexpr int Gas = IndexTraits::gasPhaseIdx;
53 static constexpr int NUM_PHASES = 3;
54 static constexpr Scalar ALQ_EPSILON = 1e-8;
55
56public:
57 using GLiftSyncGroups = std::set<int>;
58 using Rate = typename GasLiftGroupInfo<Scalar, IndexTraits>::Rate;
59 using MessageType = typename GasLiftCommon<Scalar, IndexTraits>::MessageType;
60
61 struct GradInfo
62 {
63 GradInfo() = default;
64 GradInfo(Scalar grad_,
65 Scalar new_oil_rate_,
66 bool oil_is_limited_,
67 Scalar new_gas_rate_,
68 bool gas_is_limited_,
69 Scalar new_water_rate_,
70 bool water_is_limited_,
71 Scalar alq_,
72 bool alq_is_limited_)
73 : grad{grad_}
74 , new_oil_rate{new_oil_rate_}
75 , oil_is_limited{oil_is_limited_}
76 , new_gas_rate{new_gas_rate_}
77 , gas_is_limited{gas_is_limited_}
78 , new_water_rate{new_water_rate_}
79 , water_is_limited{water_is_limited_}
80 , alq{alq_}
81 , alq_is_limited{alq_is_limited_}
82 {}
83
84 Scalar grad;
85 Scalar new_oil_rate;
86 bool oil_is_limited;
87 Scalar new_gas_rate;
88 bool gas_is_limited;
89 Scalar new_water_rate;
90 bool water_is_limited;
91 Scalar alq;
92 bool alq_is_limited;
93 };
94
95 const std::string& name() const { return well_name_; }
96
97 std::optional<GradInfo> calcIncOrDecGradient(Scalar oil_rate,
98 Scalar gas_rate,
99 Scalar water_rate,
100 Scalar alq,
101 const std::string& gr_name_dont_limit,
102 bool increase,
103 bool debug_output = true) const;
104
105 std::unique_ptr<GasLiftWellState<Scalar>> runOptimize(const int iteration_idx);
106
107 std::pair<Scalar, bool> wellTestALQ();
108
109 virtual const WellInterfaceGeneric<Scalar, IndexTraits>& getWell() const = 0;
110
111protected:
114 const GroupState<Scalar>& group_state,
115 const Well& ecl_well,
116 const SummaryState& summary_state,
118 const Schedule& schedule,
119 const int report_step_idx,
120 GLiftSyncGroups& sync_groups,
121 const Parallel::Communication& comm,
122 bool glift_debug);
123
124 struct LimitedRates;
125 struct BasicRates
126 {
127 BasicRates(const BasicRates& rates) :
128 oil{rates.oil},
129 gas{rates.gas},
130 water{rates.water},
131 bhp_is_limited{rates.bhp_is_limited}
132 {}
133
134 BasicRates(Scalar oil_,
135 Scalar gas_,
136 Scalar water_,
137 bool bhp_is_limited_)
138 : oil{oil_}
139 , gas{gas_}
140 , water{water_}
141 , bhp_is_limited{bhp_is_limited_}
142 {}
143
144 BasicRates& operator=(const BasicRates& rates)
145 {
146 oil = rates.oil;
147 gas = rates.gas;
148 water = rates.water;
149 bhp_is_limited = rates.bhp_is_limited;
150 return *this;
151 }
152
153 // This copy constructor cannot be defined inline here since LimitedRates
154 // has not been defined yet (it is defined below). Instead it is defined in
155 // in the .cpp file
156 explicit BasicRates(const LimitedRates& rates);
157
158 Scalar operator[](Rate rate_type) const
159 {
160 switch (rate_type) {
161 case Rate::oil:
162 return this->oil;
163 case Rate::gas:
164 return this->gas;
165 case Rate::water:
166 return this->water;
167 case Rate::liquid:
168 return this->oil + this->water;
169 default:
170 throw std::runtime_error("This should not happen");
171 }
172 }
173
174 Scalar oil, gas, water;
175 bool bhp_is_limited;
176 };
177
178 struct LimitedRates : public BasicRates
179 {
180 enum class LimitType {well, group, none};
181 LimitedRates(Scalar oil_,
182 Scalar gas_,
183 Scalar water_,
184 bool oil_is_limited_,
185 bool gas_is_limited_,
186 bool water_is_limited_,
187 bool bhp_is_limited_,
188 std::optional<Rate> oil_limiting_target_,
189 std ::optional<Rate> water_limiting_target_)
190 : BasicRates(oil_, gas_, water_, bhp_is_limited_)
191 , oil_is_limited{oil_is_limited_}
192 , gas_is_limited{gas_is_limited_}
193 , water_is_limited{water_is_limited_}
194 , oil_limiting_target{oil_limiting_target_}
195 , water_limiting_target{water_limiting_target_}
196 {
197 set_initial_limit_type_();
198 }
199
200 LimitedRates(const BasicRates& rates,
201 bool oil_is_limited_,
202 bool gas_is_limited_,
203 bool water_is_limited_)
204 : BasicRates(rates)
205 , oil_is_limited{oil_is_limited_}
206 , gas_is_limited{gas_is_limited_}
207 , water_is_limited{water_is_limited_}
208 {
209 set_initial_limit_type_();
210 }
211
212 bool limited() const
213 {
214 return oil_is_limited || gas_is_limited || water_is_limited;
215 }
216
217 // For a given ALQ value, were the rates limited due to group targets
218 // or due to well targets?
219 LimitType limit_type;
220 bool oil_is_limited;
221 bool gas_is_limited;
222 bool water_is_limited;
223 std::optional<Rate> oil_limiting_target;
224 std::optional<Rate> water_limiting_target;
225
226 private:
227 void set_initial_limit_type_()
228 {
229 limit_type = limited() ? LimitType::well : LimitType::none;
230 }
231 };
232
233 struct OptimizeState
234 {
235 OptimizeState( GasLiftSingleWellGeneric& parent_, bool increase_ )
236 : parent{parent_}
237 , increase{increase_}
238 , it{0}
239 , stop_iteration{false}
240 , bhp{-1}
241 {}
242
243 GasLiftSingleWellGeneric& parent;
244 bool increase;
245 int it;
246 bool stop_iteration;
247 Scalar bhp;
248
249 std::pair<std::optional<Scalar>,bool> addOrSubtractAlqIncrement(Scalar alq);
250 Scalar calcEcoGradient(Scalar oil_rate,
251 Scalar new_oil_rate,
252 Scalar gas_rate,
253 Scalar new_gas_rate);
254
255 bool checkAlqOutsideLimits(Scalar alq, Scalar oil_rate);
256 bool checkEcoGradient(Scalar gradient);
257 bool checkOilRateExceedsTarget(Scalar oil_rate);
258 bool checkRatesViolated(const LimitedRates& rates) const;
259
260 void debugShowIterationInfo(Scalar alq);
261
262 Scalar getBhpWithLimit();
263
264 void warn_(const std::string& msg) { parent.displayWarning_(msg); }
265 };
266
267 bool checkGroupALQrateExceeded(Scalar delta_alq,
268 const std::string& gr_name_dont_limit = "") const;
269 bool checkGroupTotalRateExceeded(Scalar delta_alq,
270 Scalar delta_gas_rate) const;
271
272 std::pair<std::optional<Scalar>, bool>
273 addOrSubtractAlqIncrement_(Scalar alq, bool increase) const;
274
275 Scalar calcEcoGradient_(Scalar oil_rate, Scalar new_oil_rate,
276 Scalar gas_rate, Scalar new_gas_rate, bool increase) const;
277
278 bool checkALQequal_(Scalar alq1, Scalar alq2) const;
279
280 bool checkGroupTargetsViolated(const BasicRates& rates,
281 const BasicRates& new_rates) const;
282 bool checkInitialALQmodified_(Scalar alq, Scalar initial_alq) const;
283
284 virtual bool checkThpControl_() const = 0;
285 virtual std::optional<Scalar > computeBhpAtThpLimit_(Scalar alq,
286 bool debug_output = true) const = 0;
287
288 std::pair<std::optional<Scalar>,Scalar>
289 computeConvergedBhpAtThpLimitByMaybeIncreasingALQ_() const;
290
291 std::pair<std::optional<BasicRates>,Scalar>
292 computeInitialWellRates_() const;
293
294 std::optional<LimitedRates>
295 computeLimitedWellRatesWithALQ_(Scalar alq) const;
296
297 virtual BasicRates computeWellRates_(Scalar bhp,
298 bool bhp_is_limited,
299 bool debug_output = true) const = 0;
300
301 std::optional<BasicRates> computeWellRatesWithALQ_(Scalar alq) const;
302
303 void debugCheckNegativeGradient_(Scalar grad, Scalar alq, Scalar new_alq,
304 Scalar oil_rate, Scalar new_oil_rate,
305 Scalar gas_rate, Scalar new_gas_rate,
306 bool increase) const;
307
308 void debugPrintWellStateRates() const;
309 void debugShowAlqIncreaseDecreaseCounts_();
310 void debugShowBhpAlqTable_();
311 void debugShowLimitingTargets_(const LimitedRates& rates) const;
312 void debugShowProducerControlMode() const;
313 void debugShowStartIteration_(Scalar alq, bool increase, Scalar oil_rate);
314 void debugShowTargets_();
315 void displayDebugMessage_(const std::string& msg) const override;
316 void displayWarning_(const std::string& warning);
317
318 std::pair<Scalar, bool> getBhpWithLimit_(Scalar bhp) const;
319 std::pair<Scalar, bool> getGasRateWithLimit_(const BasicRates& rates) const;
320 std::pair<Scalar, bool> getGasRateWithGroupLimit_(Scalar new_gas_rate,
321 Scalar gas_rate,
322 const std::string& gr_name_dont_limit) const;
323
324 std::pair<std::optional<LimitedRates>,Scalar >
325 getInitialRatesWithLimit_() const;
326
327 LimitedRates getLimitedRatesFromRates_(const BasicRates& rates) const;
328
329 std::tuple<Scalar,Scalar,bool,bool>
330 getLiquidRateWithGroupLimit_(const Scalar new_oil_rate,
331 const Scalar oil_rate,
332 const Scalar new_water_rate,
333 const Scalar water_rate,
334 const std::string& gr_name_dont_limit) const;
335
336 std::pair<Scalar, bool>
337 getOilRateWithGroupLimit_(Scalar new_oil_rate,
338 Scalar oil_rate,
339 const std::string& gr_name_dont_limit) const;
340
341 std::pair<Scalar, bool> getOilRateWithLimit_(const BasicRates& rates) const;
342
343 std::pair<Scalar, std::optional<Rate>>
344 getOilRateWithLimit2_(const BasicRates& rates) const;
345
346 Scalar getProductionTarget_(Rate rate) const;
347 Scalar getRate_(Rate rate_type, const BasicRates& rates) const;
348
349 std::pair<Scalar, std::optional<Rate>>
350 getRateWithLimit_(Rate rate_type, const BasicRates& rates) const;
351
352 std::tuple<Scalar, const std::string*, Scalar>
353 getRateWithGroupLimit_(Rate rate_type,
354 const Scalar new_rate,
355 const Scalar old_rate,
356 const std::string& gr_name_dont_limit) const;
357
358 std::pair<Scalar, bool>
359 getWaterRateWithGroupLimit_(Scalar new_water_rate,
360 Scalar water_rate,
361 const std::string& gr_name_dont_limit) const;
362
363 std::pair<Scalar, bool> getWaterRateWithLimit_(const BasicRates& rates) const;
364
365 std::pair<Scalar, std::optional<Rate>>
366 getWaterRateWithLimit2_(const BasicRates& rates) const;
367
368 BasicRates getWellStateRates_() const;
369 bool hasProductionControl_(Rate rate) const;
370
371 std::pair<LimitedRates, Scalar>
372 increaseALQtoPositiveOilRate_(Scalar alq,
373 const LimitedRates& orig_rates) const;
374
375 std::pair<LimitedRates, Scalar>
376 increaseALQtoMinALQ_(Scalar alq,
377 const LimitedRates& orig_rates) const;
378
379 void logSuccess_(Scalar alq,
380 const int iteration_idx);
381
382 std::pair<LimitedRates, Scalar>
383 maybeAdjustALQbeforeOptimizeLoop_(const LimitedRates& rates,
384 Scalar alq,
385 bool increase) const;
386
387 std::pair<LimitedRates, Scalar>
388 reduceALQtoGroupAlqLimits_(Scalar alq,
389 const LimitedRates& rates) const;
390
391 std::pair<LimitedRates, Scalar>
392 reduceALQtoGroupTarget(Scalar alq,
393 const LimitedRates& rates) const;
394
395 std::pair<LimitedRates, Scalar>
396 reduceALQtoWellTarget_(Scalar alq,
397 const LimitedRates& rates) const;
398
399 std::unique_ptr<GasLiftWellState<Scalar>> runOptimize1_();
400 std::unique_ptr<GasLiftWellState<Scalar>> runOptimize2_();
401 std::unique_ptr<GasLiftWellState<Scalar>> runOptimizeLoop_(bool increase);
402
403 void setAlqMinRate_(const GasLiftWell& well);
404 std::unique_ptr<GasLiftWellState<Scalar>> tryIncreaseLiftGas_();
405 std::unique_ptr<GasLiftWellState<Scalar>> tryDecreaseLiftGas_();
406
407 void updateGroupRates_(const LimitedRates& rates,
408 const LimitedRates& new_rates,
409 Scalar delta_alq) const;
410
412 updateRatesToGroupLimits_(const BasicRates& old_rates,
413 const LimitedRates& rates,
414 const std::string& gr_name = "") const;
415
416 void updateWellStateAlqFixedValue_(const GasLiftWell& well);
417 bool useFixedAlq_(const GasLiftWell& well);
418
419 void debugInfoGroupRatesExceedTarget(Rate rate_type,
420 const std::string& gr_name,
421 Scalar rate,
422 Scalar target) const;
423 void warnMaxIterationsExceeded_();
424
425 const Well& ecl_well_;
426 const SummaryState& summary_state_;
428 GLiftSyncGroups& sync_groups_;
429 const WellProductionControls controls_;
430
431 Scalar increment_;
432 Scalar max_alq_;
433 Scalar min_alq_;
434 Scalar orig_alq_;
435
436 Scalar alpha_w_;
437 Scalar alpha_g_;
438 Scalar eco_grad_;
439
440 int gas_pos_;
441 int oil_pos_;
442 int water_pos_;
443
444 int max_iterations_;
445
446 std::string well_name_;
447
448 const GasLiftWell* gl_well_;
449
450 bool optimize_;
451 bool debug_limit_increase_decrease_;
452 bool debug_abort_if_decrease_and_oil_is_limited_ = false;
453 bool debug_abort_if_increase_and_gas_is_limited_ = false;
454};
455
456} // namespace Opm
457
458#endif // OPM_GASLIFT_SINGLE_WELL_GENERIC_HEADER_INCLUDED
Definition DeferredLogger.hpp:57
Definition GasLiftGroupInfo.hpp:46
Definition GasLiftSingleWellGeneric.hpp:48
Definition GasLiftWellState.hpp:30
Definition GroupState.hpp:41
Definition WellInterfaceGeneric.hpp:53
The state of a set of wells, tailored for use by the fully implicit blackoil simulator.
Definition WellState.hpp:66
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:43
Definition GasLiftSingleWellGeneric.hpp:126
Definition GasLiftSingleWellGeneric.hpp:179