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>::BasicRates
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::transform(potentials.begin(), potentials.end(),
134 [](
const auto& potential)
135 { return std::min(Scalar{0}, potential); });
136 return {-potentials[this->oil_pos_],
137 -potentials[this->gas_pos_],
138 -potentials[this->water_pos_],
143template<
typename TypeTag>
144std::optional<typename GasLiftSingleWell<TypeTag>::Scalar>
145GasLiftSingleWell<TypeTag>::
146computeBhpAtThpLimit_(Scalar alq,
bool debug_output)
const
149 auto bhp_at_thp_limit = this->well_.computeBhpAtThpLimitProdWithAlq(
151 this->summary_state_,
153 this->deferred_logger_,
155 if (bhp_at_thp_limit) {
156 if (*bhp_at_thp_limit < this->controls_.bhp_limit) {
157 if (debug_output && this->debug) {
158 const std::string msg = fmt::format(
159 "Computed bhp ({}) from thp limit is below bhp limit ({}), (ALQ = {})."
160 " Using bhp limit instead",
161 *bhp_at_thp_limit, this->controls_.bhp_limit, alq
163 this->displayDebugMessage_(msg);
165 bhp_at_thp_limit = this->controls_.bhp_limit;
170 const std::string msg = fmt::format(
171 "Failed in getting converged bhp potential from thp limit (ALQ = {})", alq);
172 this->displayDebugMessage_(msg);
174 return bhp_at_thp_limit;
177template<
typename TypeTag>
179GasLiftSingleWell<TypeTag>::
180setupPhaseVariables_()
183 bool num_phases_ok = (FluidSystem::numActivePhases()== 3);
185 if (FluidSystem::numActivePhases()== 2) {
196 if ( FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)
197 && FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)
198 && !FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx) )
201 num_phases_ok =
true;
205 throw std::logic_error(
"Two-phase gas lift optimization only supported"
206 " for oil and water");
209 assert(num_phases_ok);
210 this->oil_pos_ = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::oilPhaseIdx);
211 this->gas_pos_ = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::gasPhaseIdx);
212 this->water_pos_ = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::waterPhaseIdx);
215template<
typename TypeTag>
217GasLiftSingleWell<TypeTag>::
218setAlqMaxRate_(
const GasLiftWell& well)
220 const auto& max_alq_optional = well.max_rate();
221 if (max_alq_optional) {
224 this->max_alq_ = *max_alq_optional;
230 const auto& table = well_.vfpProperties()->getProd()->getTable(
231 this->controls_.vfp_table_number);
232 const auto& alq_values = table.getALQAxis();
235 this->max_alq_ = alq_values.back();
239template<
typename TypeTag>
241GasLiftSingleWell<TypeTag>::
242checkThpControl_()
const
244 const int well_index = this->well_state_.index(this->well_name_).value();
245 const Well::ProducerCMode& control_mode =
246 this->well_state_.well(well_index).production_cmode;
247 bool thp_control = control_mode == Well::ProducerCMode::THP;
248 const auto& well = getWell();
249 thp_control = thp_control || well.thpLimitViolatedButNotSwitched();
252 this->displayDebugMessage_(
"Well is not under THP control, skipping iteration..");
Definition DeferredLogger.hpp:57
Definition GasLiftGroupInfo.hpp:46
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:76
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