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