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
30namespace Opm {
31
32template<typename TypeTag>
35 const Simulator& simulator,
36 const SummaryState &summary_state,
37 DeferredLogger &deferred_logger,
38 WellState<Scalar>& well_state,
39 const GroupState<Scalar>& group_state,
40 GasLiftGroupInfo &group_info,
41 GLiftSyncGroups &sync_groups,
42 const Parallel::Communication& comm,
43 bool glift_debug
44 )
45 // The parent class GasLiftSingleWellGeneric contains all stuff
46 // that is not dependent on TypeTag
48 deferred_logger,
49 well_state,
50 group_state,
51 well.wellEcl(),
52 summary_state,
53 group_info,
54 well.phaseUsage(),
55 simulator.vanguard().schedule(),
56 simulator.episodeIndex(),
57 sync_groups,
58 comm,
59 glift_debug
60 )
61 , simulator_{simulator}
62 , well_{well}
63{
64 const auto& gl_well = *gl_well_;
65 if(useFixedAlq_(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 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 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>
112computeWellRates_( double bhp, bool bhp_is_limited, bool debug_output ) const
113{
114 std::vector<double> potentials(NUM_PHASES, 0.0);
115 this->well_.computeWellRatesWithBhp(
116 this->simulator_, bhp, potentials, this->deferred_logger_);
117 if (debug_output) {
118 const std::string msg = fmt::format("computed well potentials given bhp {}, "
119 "oil: {}, gas: {}, water: {}", bhp,
120 -potentials[this->oil_pos_], -potentials[this->gas_pos_],
121 -potentials[this->water_pos_]);
122 displayDebugMessage_(msg);
123 }
124
125 for (auto& potential : potentials) {
126 potential = std::min(0.0, potential);
127 }
128 return {-potentials[this->oil_pos_],
129 -potentials[this->gas_pos_],
130 -potentials[this->water_pos_],
131 bhp_is_limited
132 };
133}
134
135template<typename TypeTag>
136std::optional<double>
137GasLiftSingleWell<TypeTag>::
138computeBhpAtThpLimit_(double alq, bool debug_output) const
139{
140 auto bhp_at_thp_limit = this->well_.computeBhpAtThpLimitProdWithAlq(
141 this->simulator_,
142 this->summary_state_,
143 alq,
144 this->deferred_logger_);
145 if (bhp_at_thp_limit) {
146 if (*bhp_at_thp_limit < this->controls_.bhp_limit) {
147 if (debug_output && this->debug) {
148 const std::string msg = fmt::format(
149 "Computed bhp ({}) from thp limit is below bhp limit ({}), (ALQ = {})."
150 " Using bhp limit instead",
151 *bhp_at_thp_limit, this->controls_.bhp_limit, alq
152 );
153 displayDebugMessage_(msg);
154 }
155 bhp_at_thp_limit = this->controls_.bhp_limit;
156 }
157 //bhp_at_thp_limit = std::max(*bhp_at_thp_limit, this->controls_.bhp_limit);
158 }
159 else {
160 const std::string msg = fmt::format(
161 "Failed in getting converged bhp potential from thp limit (ALQ = {})", alq);
162 displayDebugMessage_(msg);
163 }
164 return bhp_at_thp_limit;
165}
166
167template<typename TypeTag>
168void
169GasLiftSingleWell<TypeTag>::
170setupPhaseVariables_()
171{
172 const auto& pu = this->phase_usage_;
173#ifndef NDEBUG
174 bool num_phases_ok = (pu.num_phases == 3);
175#endif
176 if (pu.num_phases == 2) {
177 // NOTE: We support two-phase oil-water flow, by setting the gas flow rate
178 // to zero. This is done by initializing the potential vector to zero:
179 //
180 // std::vector<double> potentials(NUM_PHASES, 0.0);
181 //
182 // see e.g. runOptimizeLoop_() in GasLiftSingleWellGeneric.cpp
183 // In addition the VFP calculations, e.g. to calculate BHP from THP
184 // has been adapted to the two-phase oil-water case, see the comment
185 // in WellInterfaceGeneric.cpp for the method adaptRatesForVFP() for
186 // more information.
187 if ( pu.phase_used[BlackoilPhases::Aqua] == 1
188 && pu.phase_used[BlackoilPhases::Liquid] == 1
189 && pu.phase_used[BlackoilPhases::Vapour] == 0)
190 {
191#ifndef NDEBUG
192 num_phases_ok = true; // two-phase oil-water is also supported
193#endif
194 }
195 else {
196 throw std::logic_error("Two-phase gas lift optimization only supported"
197 " for oil and water");
198 }
199 }
200 assert(num_phases_ok);
201 this->oil_pos_ = pu.phase_pos[Oil];
202 this->gas_pos_ = pu.phase_pos[Gas];
203 this->water_pos_ = pu.phase_pos[Water];
204}
205
206template<typename TypeTag>
207void
208GasLiftSingleWell<TypeTag>::
209setAlqMaxRate_(const GasLiftWell &well)
210{
211 auto& max_alq_optional = well.max_rate();
212 if (max_alq_optional) {
213 // NOTE: To prevent extrapolation of the VFP tables, any value
214 // entered here must not exceed the largest ALQ value in the well's VFP table.
215 this->max_alq_ = *max_alq_optional;
216 }
217 else { // i.e. WLIFTOPT, item 3 has been defaulted
218 // According to the manual for WLIFTOPT, item 3:
219 // The default value should be set to the largest ALQ
220 // value in the well's VFP table
221 const auto& table = well_.vfpProperties()->getProd()->getTable(
222 this->controls_.vfp_table_number);
223 const auto& alq_values = table.getALQAxis();
224 // Assume the alq_values are sorted in ascending order, so
225 // the last item should be the largest value:
226 this->max_alq_ = alq_values.back();
227 }
228}
229
230template<typename TypeTag>
231bool
232GasLiftSingleWell<TypeTag>::
233checkThpControl_() const
234{
235 const int well_index = this->well_state_.index(this->well_name_).value();
236 const Well::ProducerCMode& control_mode =
237 this->well_state_.well(well_index).production_cmode;
238 bool thp_control = control_mode == Well::ProducerCMode::THP;
239 const WellInterfaceGeneric &well = getWell();
240 thp_control = thp_control || well.thpLimitViolatedButNotSwitched();
241 if (this->debug) {
242 if (!thp_control) {
243 displayDebugMessage_("Well is not under THP control, skipping iteration..");
244 }
245 }
246 return thp_control;
247}
248
249}
@ Liquid
Definition: BlackoilPhases.hpp:42
@ Aqua
Definition: BlackoilPhases.hpp:42
@ Vapour
Definition: BlackoilPhases.hpp:42
Definition: DeferredLogger.hpp:57
WellState< double > & well_state_
Definition: GasLiftCommon.hpp:57
Definition: GasLiftGroupInfo.hpp:45
Definition: GasLiftSingleWellGeneric.hpp:49
double alpha_w_
Definition: GasLiftSingleWellGeneric.hpp:345
int max_iterations_
Definition: GasLiftSingleWellGeneric.hpp:353
const GasLiftWell * gl_well_
Definition: GasLiftSingleWellGeneric.hpp:357
double orig_alq_
Definition: GasLiftSingleWellGeneric.hpp:343
void setAlqMinRate_(const GasLiftWell &well)
bool optimize_
Definition: GasLiftSingleWellGeneric.hpp:359
void displayWarning_(const std::string &warning)
double alpha_g_
Definition: GasLiftSingleWellGeneric.hpp:346
bool useFixedAlq_(const GasLiftWell &well)
std::string well_name_
Definition: GasLiftSingleWellGeneric.hpp:355
void updateWellStateAlqFixedValue_(const GasLiftWell &well)
Definition: GasLiftSingleWell.hpp:38
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 &group_info, GLiftSyncGroups &sync_groups, const Parallel::Communication &comm, bool glift_debug)
Definition: GasLiftSingleWell_impl.hpp:34
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.
Definition: GasLiftSingleWellGeneric.hpp:122