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