GasLiftSingleWell_impl.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// Improve IDE experience
21#ifndef OPM_GASLIFT_SINGLE_WELL_HEADER_INCLUDED
22#define OPM_GASLIFT_SINGLE_WELL_IMPL_HEADER_INCLUDED
23#include <config.h>
25#endif
26
27#include <opm/input/eclipse/Schedule/GasLiftOpt.hpp>
28#include <fmt/format.h>
29
30#include <string>
31#include <vector>
32
33namespace Opm {
34
35template<typename TypeTag>
38 const Simulator& simulator,
39 const SummaryState& summary_state,
40 DeferredLogger& deferred_logger,
41 WellState<Scalar>& well_state,
42 const GroupState<Scalar>& group_state,
43 GasLiftGroupInfo<Scalar>& group_info,
44 GLiftSyncGroups &sync_groups,
45 const Parallel::Communication& comm,
46 bool glift_debug)
47 // The parent class GasLiftSingleWellGeneric contains all stuff
48 // that is not dependent on TypeTag
49 : GasLiftSingleWellGeneric<Scalar>(deferred_logger,
50 well_state,
51 group_state,
52 well.wellEcl(),
53 summary_state,
54 group_info,
55 well.phaseUsage(),
56 simulator.vanguard().schedule(),
57 simulator.episodeIndex(),
58 sync_groups,
59 comm,
60 glift_debug)
61 , simulator_{simulator}
62 , well_{well}
63{
64 const auto& gl_well = *this->gl_well_;
65 if (this->useFixedAlq_(gl_well)) {
66 this->updateWellStateAlqFixedValue_(gl_well);
67 this->optimize_ = false; // lift gas supply is fixed
68 }
69 else {
70 setAlqMaxRate_(gl_well);
71 this->optimize_ = true;
72 }
73
74 setupPhaseVariables_();
75 // get the alq value used for this well for the previous iteration (a
76 // nonlinear iteration in assemble() in BlackoilWellModel).
77 // If gas lift optimization has not been applied to this well yet, the
78 // default value is used.
79 this->orig_alq_ = this->well_state_.getALQ(this->well_name_);
80 if (this->optimize_) {
81 this->setAlqMinRate_(gl_well);
82 // NOTE: According to item 4 in WLIFTOPT, this value does not
83 // have to be positive.
84 // TODO: Does it make sense to have a negative value?
85 this->alpha_w_ = gl_well.weight_factor();
86 if (this->alpha_w_ <= 0 ) {
87 this->displayWarning_("Nonpositive value for alpha_w ignored");
88 this->alpha_w_ = 1.0;
89 }
90
91 // NOTE: According to item 6 in WLIFTOPT:
92 // "If this value is greater than zero, the incremental gas rate will influence
93 // the calculation of the incremental gradient and may be used
94 // to discourage the allocation of lift gas to wells which produce more gas."
95 // TODO: Does this mean that we should ignore this value if it
96 // is negative?
97 this->alpha_g_ = gl_well.inc_weight_factor();
98
99 // TODO: adhoc value.. Should we keep max_iterations_ as a safety measure
100 // or does it not make sense to have it?
101 this->max_iterations_ = 1000;
102 }
103}
104
105/****************************************
106 * Private methods in alphabetical order
107 ****************************************/
108
109template<typename TypeTag>
110typename GasLiftSingleWell<TypeTag>::BasicRates
112computeWellRates_(Scalar bhp, bool bhp_is_limited, bool debug_output ) const
113{
114 std::vector<Scalar> potentials(this->NUM_PHASES, 0.0);
115 this->well_.computeWellRatesWithBhp(this->simulator_,
116 bhp,
117 potentials,
118 this->deferred_logger_);
119 if (debug_output) {
120 const std::string msg = fmt::format("computed well potentials given bhp {}, "
121 "oil: {}, gas: {}, water: {}", bhp,
122 -potentials[this->oil_pos_], -potentials[this->gas_pos_],
123 -potentials[this->water_pos_]);
124 this->displayDebugMessage_(msg);
125 }
126
127 for (auto& potential : potentials) {
128 potential = std::min(Scalar{0.0}, potential);
129 }
130 return {-potentials[this->oil_pos_],
131 -potentials[this->gas_pos_],
132 -potentials[this->water_pos_],
133 bhp_is_limited
134 };
135}
136
137template<typename TypeTag>
138std::optional<typename GasLiftSingleWell<TypeTag>::Scalar>
139GasLiftSingleWell<TypeTag>::
140computeBhpAtThpLimit_(Scalar alq, bool debug_output) const
141{
142 auto bhp_at_thp_limit = this->well_.computeBhpAtThpLimitProdWithAlq(
143 this->simulator_,
144 this->summary_state_,
145 alq,
146 this->deferred_logger_);
147 if (bhp_at_thp_limit) {
148 if (*bhp_at_thp_limit < this->controls_.bhp_limit) {
149 if (debug_output && this->debug) {
150 const std::string msg = fmt::format(
151 "Computed bhp ({}) from thp limit is below bhp limit ({}), (ALQ = {})."
152 " Using bhp limit instead",
153 *bhp_at_thp_limit, this->controls_.bhp_limit, alq
154 );
155 this->displayDebugMessage_(msg);
156 }
157 bhp_at_thp_limit = this->controls_.bhp_limit;
158 }
159 //bhp_at_thp_limit = std::max(*bhp_at_thp_limit, this->controls_.bhp_limit);
160 }
161 else {
162 const std::string msg = fmt::format(
163 "Failed in getting converged bhp potential from thp limit (ALQ = {})", alq);
164 this->displayDebugMessage_(msg);
165 }
166 return bhp_at_thp_limit;
167}
168
169template<typename TypeTag>
170void
171GasLiftSingleWell<TypeTag>::
172setupPhaseVariables_()
173{
174 const auto& pu = this->phase_usage_;
175#ifndef NDEBUG
176 bool num_phases_ok = (pu.num_phases == 3);
177#endif
178 if (pu.num_phases == 2) {
179 // NOTE: We support two-phase oil-water flow, by setting the gas flow rate
180 // to zero. This is done by initializing the potential vector to zero:
181 //
182 // std::vector<double> potentials(NUM_PHASES, 0.0);
183 //
184 // see e.g. runOptimizeLoop_() in GasLiftSingleWellGeneric.cpp
185 // In addition the VFP calculations, e.g. to calculate BHP from THP
186 // has been adapted to the two-phase oil-water case, see the comment
187 // in WellInterfaceGeneric.cpp for the method adaptRatesForVFP() for
188 // more information.
189 if ( pu.phase_used[BlackoilPhases::Aqua] == 1
190 && pu.phase_used[BlackoilPhases::Liquid] == 1
191 && pu.phase_used[BlackoilPhases::Vapour] == 0)
192 {
193#ifndef NDEBUG
194 num_phases_ok = true; // two-phase oil-water is also supported
195#endif
196 }
197 else {
198 throw std::logic_error("Two-phase gas lift optimization only supported"
199 " for oil and water");
200 }
201 }
202 assert(num_phases_ok);
203 this->oil_pos_ = pu.phase_pos[this->Oil];
204 this->gas_pos_ = pu.phase_pos[this->Gas];
205 this->water_pos_ = pu.phase_pos[this->Water];
206}
207
208template<typename TypeTag>
209void
210GasLiftSingleWell<TypeTag>::
211setAlqMaxRate_(const GasLiftWell& well)
212{
213 auto& max_alq_optional = well.max_rate();
214 if (max_alq_optional) {
215 // NOTE: To prevent extrapolation of the VFP tables, any value
216 // entered here must not exceed the largest ALQ value in the well's VFP table.
217 this->max_alq_ = *max_alq_optional;
218 }
219 else { // i.e. WLIFTOPT, item 3 has been defaulted
220 // According to the manual for WLIFTOPT, item 3:
221 // The default value should be set to the largest ALQ
222 // value in the well's VFP table
223 const auto& table = well_.vfpProperties()->getProd()->getTable(
224 this->controls_.vfp_table_number);
225 const auto& alq_values = table.getALQAxis();
226 // Assume the alq_values are sorted in ascending order, so
227 // the last item should be the largest value:
228 this->max_alq_ = alq_values.back();
229 }
230}
231
232template<typename TypeTag>
233bool
234GasLiftSingleWell<TypeTag>::
235checkThpControl_() const
236{
237 const int well_index = this->well_state_.index(this->well_name_).value();
238 const Well::ProducerCMode& control_mode =
239 this->well_state_.well(well_index).production_cmode;
240 bool thp_control = control_mode == Well::ProducerCMode::THP;
241 const auto& well = getWell();
242 thp_control = thp_control || well.thpLimitViolatedButNotSwitched();
243 if (this->debug) {
244 if (!thp_control) {
245 this->displayDebugMessage_("Well is not under THP control, skipping iteration..");
246 }
247 }
248 return thp_control;
249}
250
251}
@ Liquid
Definition: BlackoilPhases.hpp:42
@ Aqua
Definition: BlackoilPhases.hpp:42
@ Vapour
Definition: BlackoilPhases.hpp:42
Definition: DeferredLogger.hpp:57
WellState< GetPropType< TypeTag, Properties::Scalar > > & well_state_
Definition: GasLiftCommon.hpp:55
Definition: GasLiftGroupInfo.hpp:45
Definition: GasLiftSingleWellGeneric.hpp:49
bool optimize_
Definition: GasLiftSingleWellGeneric.hpp:453
std::string well_name_
Definition: GasLiftSingleWellGeneric.hpp:449
GetPropType< TypeTag, Properties::Scalar > alpha_w_
Definition: GasLiftSingleWellGeneric.hpp:439
GetPropType< TypeTag, Properties::Scalar > alpha_g_
Definition: GasLiftSingleWellGeneric.hpp:440
GetPropType< TypeTag, Properties::Scalar > orig_alq_
Definition: GasLiftSingleWellGeneric.hpp:437
int max_iterations_
Definition: GasLiftSingleWellGeneric.hpp:447
const GasLiftWell * gl_well_
Definition: GasLiftSingleWellGeneric.hpp:451
Definition: GasLiftSingleWell.hpp:36
GasLiftSingleWell(const WellInterface< TypeTag > &well, const Simulator &simulator, const SummaryState &summary_state, DeferredLogger &deferred_logger, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, GasLiftGroupInfo< Scalar > &group_info, GLiftSyncGroups &sync_groups, const Parallel::Communication &comm, bool glift_debug)
Definition: GasLiftSingleWell_impl.hpp:37
Definition: GroupState.hpp:35
Definition: WellInterface.hpp:75
Definition: WellState.hpp:62
Scalar getALQ(const std::string &name) const
Definition: WellState.hpp:175
Dune::Communication< MPIComm > Communication
Definition: ParallelCommunication.hpp:30
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
PhaseUsage phaseUsage(const Phases &phases)
Determine the active phases.