20#ifndef OPM_GASLIFT_SINGLE_WELL_IMPL_HEADER_INCLUDED
21#define OPM_GASLIFT_SINGLE_WELL_IMPL_HEADER_INCLUDED
24#ifndef OPM_GASLIFT_SINGLE_WELL_HEADER_INCLUDED
26#include <opm/simulators/wells/GasLiftSingleWell.hpp>
29#include <opm/common/TimingMacros.hpp>
31#include <opm/input/eclipse/Schedule/GasLiftOpt.hpp>
32#include <opm/input/eclipse/Schedule/Well/Well.hpp>
37#include <fmt/format.h>
41template<
typename TypeTag>
45 const SummaryState& summary_state,
50 GLiftSyncGroups &sync_groups,
51 const Parallel::Communication& comm,
61 simulator.vanguard().schedule(),
62 simulator.episodeIndex(),
66 , simulator_{simulator}
69 const auto& gl_well = *this->gl_well_;
70 if (this->useFixedAlq_(gl_well)) {
71 this->updateWellStateAlqFixedValue_(gl_well);
72 this->optimize_ =
false;
75 setAlqMaxRate_(gl_well);
76 this->optimize_ =
true;
79 setupPhaseVariables_();
84 this->orig_alq_ = this->well_state_.well(this->well_name_).alq_state.get();
85 if (this->optimize_) {
86 this->setAlqMinRate_(gl_well);
90 this->alpha_w_ = gl_well.weight_factor();
91 if (this->alpha_w_ <= 0 ) {
92 this->displayWarning_(
"Nonpositive value for alpha_w ignored");
102 this->alpha_g_ = gl_well.inc_weight_factor();
106 this->max_iterations_ = 1000;
114template<
typename TypeTag>
115typename GasLiftSingleWell<TypeTag>::RatesAndBhp
116GasLiftSingleWell<TypeTag>::
117computeWellRates_(Scalar bhp,
bool bhp_is_limited,
bool debug_output )
const
119 std::vector<Scalar> potentials(this->NUM_PHASES, 0.0);
120 this->well_.computeWellRatesWithBhp(this->simulator_,
123 this->deferred_logger_);
125 const std::string msg = fmt::format(
"computed well potentials given bhp {}, "
126 "oil: {}, gas: {}, water: {}", bhp,
127 -potentials[this->oil_pos_], -potentials[this->gas_pos_],
128 -potentials[this->water_pos_]);
129 this->displayDebugMessage_(msg);
132 std::ranges::transform(potentials, potentials.begin(),
133 [](
const auto& potential)
134 { return std::min(Scalar{0}, potential); });
135 return {-potentials[this->oil_pos_],
136 -potentials[this->gas_pos_],
137 -potentials[this->water_pos_],
143template<
typename TypeTag>
144std::optional<typename GasLiftSingleWell<TypeTag>::Scalar>
145GasLiftSingleWell<TypeTag>::
146computeBhpAtThpLimit_(Scalar alq, Scalar current_bhp,
bool debug_output)
const
152 auto bhp_at_thp_limit = this->well_.computeBhpAtThpLimitProdWithAlqUsingIPR(
156 this->summary_state_,
158 if (bhp_at_thp_limit) {
159 if (*bhp_at_thp_limit < this->controls_.bhp_limit) {
160 if (debug_output && this->debug) {
161 const std::string msg = fmt::format(
162 "Computed bhp ({}) from thp limit is below bhp limit ({}), (ALQ = {})."
163 " Using bhp limit instead",
164 *bhp_at_thp_limit, this->controls_.bhp_limit, alq
166 this->displayDebugMessage_(msg);
168 bhp_at_thp_limit = this->controls_.bhp_limit;
173 const std::string msg = fmt::format(
174 "Failed in getting converged bhp potential from thp limit (ALQ = {})", alq);
175 this->displayDebugMessage_(msg);
177 return bhp_at_thp_limit;
180template<
typename TypeTag>
182GasLiftSingleWell<TypeTag>::
183setupPhaseVariables_()
186 bool num_phases_ok = (FluidSystem::numActivePhases()== 3);
188 if (FluidSystem::numActivePhases()== 2) {
199 if ( FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)
200 && FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)
201 && !FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx) )
204 num_phases_ok =
true;
208 throw std::logic_error(
"Two-phase gas lift optimization only supported"
209 " for oil and water");
212 assert(num_phases_ok);
213 this->oil_pos_ = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::oilPhaseIdx);
214 this->gas_pos_ = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::gasPhaseIdx);
215 this->water_pos_ = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::waterPhaseIdx);
218template<
typename TypeTag>
220GasLiftSingleWell<TypeTag>::
221setAlqMaxRate_(
const GasLiftWell& well)
223 const auto& max_alq_optional = well.max_rate();
224 if (max_alq_optional) {
227 this->max_alq_ = *max_alq_optional;
233 const auto& table = well_.vfpProperties()->getProd()->getTable(
234 this->controls_.vfp_table_number);
235 const auto& alq_values = table.getALQAxis();
238 this->max_alq_ = alq_values.back();
242template<
typename TypeTag>
244GasLiftSingleWell<TypeTag>::
245checkThpControl_()
const
247 const int well_index = this->well_state_.index(this->well_name_).value();
248 const Well::ProducerCMode& control_mode =
249 this->well_state_.well(well_index).production_cmode;
250 bool thp_control = control_mode == Well::ProducerCMode::THP;
251 const auto& well = getWell();
252 thp_control = thp_control || well.thpLimitViolatedButNotSwitched();
255 this->displayDebugMessage_(
"Well is not under THP control, skipping iteration..");
Definition DeferredLogger.hpp:57
Definition GasLiftGroupInfo.hpp:47
Definition GasLiftSingleWellGeneric.hpp:48
Definition GasLiftSingleWell.hpp:39
Definition GroupState.hpp:41
Manages the initializing and running of time dependent problems.
Definition simulator.hh:84
Definition WellInterface.hpp:77
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