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