BlackoilWellModelGasLift_impl.hpp
Go to the documentation of this file.
1/*
2 Copyright 2016 - 2019 SINTEF Digital, Mathematics & Cybernetics.
3 Copyright 2016 - 2018 Equinor ASA.
4 Copyright 2017 Dr. Blatt - HPC-Simulation-Software & Services
5 Copyright 2016 - 2018 Norce AS
6
7 This file is part of the Open Porous Media project (OPM).
8
9 OPM is free software: you can redistribute it and/or modify
10 it under the terms of the GNU General Public License as published by
11 the Free Software Foundation, either version 3 of the License, or
12 (at your option) any later version.
13
14 OPM is distributed in the hope that it will be useful,
15 but WITHOUT ANY WARRANTY; without even the implied warranty of
16 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17 GNU General Public License for more details.
18
19 You should have received a copy of the GNU General Public License
20 along with OPM. If not, see <http://www.gnu.org/licenses/>.
21*/
22#define OPM_BLACKOILWELLMODEL_GASLIFT_IMPL_HEADER_INCLUDED
23#define OPM_BLACKOILWELLMODEL_GASLIFT_IMPL_HEADER_INCLUDED
24
25// Improve IDE experience
26#ifndef OPM_BLACKOILWELLMODEL_GASLIFT_HEADER_INCLUDED
27#include <config.h>
29#endif
30#include <opm/common/TimingMacros.hpp>
32
33#if HAVE_MPI
35#endif
36
37namespace Opm {
38
39template<typename TypeTag>
40bool
42maybeDoGasLiftOptimize(const Simulator& simulator,
43 const std::vector<WellInterfacePtr>& well_container,
44 WellState<Scalar>& wellState,
45 GroupState<Scalar>& groupState,
46 DeferredLogger& deferred_logger)
47{
48 OPM_TIMEFUNCTION();
49 const auto& glo = simulator.vanguard().schedule().glo(simulator.episodeIndex());
50 if(!glo.active()) {
51 return false;
52 }
53 bool do_glift_optimization = false;
54 int num_wells_changed = 0;
55 const double simulation_time = simulator.time();
56 const Scalar min_wait = glo.min_wait();
57 // We only optimize if a min_wait time has past.
58 // If all_newton is true we still want to optimize several times pr timestep
59 // i.e. we also optimize if check simulation_time == last_glift_opt_time_
60 // that is when the last_glift_opt_time is already updated with the current time step
61 if (simulation_time == this->last_glift_opt_time_ ||
62 simulation_time >= (this->last_glift_opt_time_ + min_wait))
63 {
64 do_glift_optimization = true;
65 this->last_glift_opt_time_ = simulation_time;
66 }
67
68 if (do_glift_optimization) {
69 GLiftOptWells glift_wells;
70 GLiftProdWells prod_wells;
71 GLiftWellStateMap state_map;
72 // NOTE: To make GasLiftGroupInfo (see below) independent of the TypeTag
73 // associated with *this (i.e. BlackoilWellModel<TypeTag>) we observe
74 // that GasLiftGroupInfo's only dependence on *this is that it needs to
75 // access the eclipse Wells in the well container (the eclipse Wells
76 // themselves are independent of the TypeTag).
77 // Hence, we extract them from the well container such that we can pass
78 // them to the GasLiftGroupInfo constructor.
79 GLiftEclWells ecl_well_map;
80 initGliftEclWellMap(well_container, ecl_well_map);
81 GasLiftGroupInfo group_info {
82 ecl_well_map,
83 simulator.vanguard().schedule(),
84 simulator.vanguard().summaryState(),
85 simulator.episodeIndex(),
86 simulator.model().newtonMethod().numIterations(),
87 phase_usage_,
88 deferred_logger,
89 wellState,
90 groupState,
91 simulator.vanguard().grid().comm(),
92 this->glift_debug
93 };
94 group_info.initialize();
95
96 gasLiftOptimizationStage1(simulator,
97 well_container,
98 wellState,
99 groupState,
100 prod_wells,
101 glift_wells,
102 group_info,
103 state_map,
104 deferred_logger);
105
106 this->gasLiftOptimizationStage2(simulator.vanguard().gridView().comm(),
107 simulator.vanguard().schedule(),
108 simulator.vanguard().summaryState(),
109 wellState,
110 groupState,
111 prod_wells,
112 glift_wells,
113 group_info,
114 state_map,
115 simulator.episodeIndex(),
116 deferred_logger);
117
118 if constexpr (glift_debug) {
119 std::vector<WellInterfaceGeneric<Scalar>*> wc;
120 wc.reserve(well_container.size());
121 std::transform(well_container.begin(), well_container.end(),
122 std::back_inserter(wc),
123 [](const auto& w)
124 { return static_cast<WellInterfaceGeneric<Scalar>*>(w.get()); });
125 this->gliftDebugShowALQ(wc,
126 wellState,
127 deferred_logger);
128 }
129 num_wells_changed = glift_wells.size();
130 }
131 num_wells_changed = simulator.vanguard().gridView().comm().sum(num_wells_changed);
132 return num_wells_changed > 0;
133}
134
135template<typename TypeTag>
136void
138gasLiftOptimizationStage1(const Simulator& simulator,
139 const std::vector<WellInterfacePtr>& well_container,
140 WellState<Scalar>& wellState,
141 GroupState<Scalar>& groupState,
142 GLiftProdWells& prod_wells,
143 GLiftOptWells &glift_wells,
144 GasLiftGroupInfo<Scalar>& group_info,
145 GLiftWellStateMap& state_map,
146 DeferredLogger& deferred_logger)
147{
148 OPM_TIMEFUNCTION();
149 auto comm = simulator.vanguard().grid().comm();
150 int num_procs = comm.size();
151 // NOTE: Gas lift optimization stage 1 seems to be difficult
152 // to do in parallel since the wells are optimized on different
153 // processes and each process needs to know the current ALQ allocated
154 // to each group it is a memeber of in order to check group limits and avoid
155 // allocating more ALQ than necessary. (Surplus ALQ is removed in
156 // stage 2). In stage1, as each well is adding ALQ, the current group ALQ needs
157 // to be communicated to the other processes. But there is no common
158 // synchronization point that all process will reach in the
159 // runOptimizeLoop_() in GasLiftSingleWell.cpp.
160 //
161 // TODO: Maybe a better solution could be invented by distributing
162 // wells according to certain parent groups. Then updated group rates
163 // might not have to be communicated to the other processors.
164
165 // Currently, the best option seems to be to run this part sequentially
166 // (not in parallel).
167 //
168 // TODO: The simplest approach seems to be if a) one process could take
169 // ownership of all the wells (the union of all the wells in the
170 // well_container_ of each process) then this process could do the
171 // optimization, while the other processes could wait for it to
172 // finish (e.g. comm.barrier()), or alternatively, b) if all
173 // processes could take ownership of all the wells. Then there
174 // would be no need for synchronization here..
175 //
176 for (int i = 0; i< num_procs; i++) {
177 int num_rates_to_sync = 0; // communication variable
178 GLiftSyncGroups groups_to_sync;
179 if (comm.rank() == i) {
180 // Run stage1: Optimize single wells while also checking group limits
181 for (const auto& well : well_container) {
182 // NOTE: Only the wells in "group_info" needs to be optimized
183 if (group_info.hasWell(well->name())) {
184 gasLiftOptimizationStage1SingleWell(well.get(),
185 simulator,
186 wellState,
187 groupState,
188 prod_wells,
189 glift_wells,
190 group_info,
191 state_map,
192 groups_to_sync,
193 deferred_logger);
194 }
195 }
196 num_rates_to_sync = groups_to_sync.size();
197 }
198 {
199 OPM_TIMEBLOCK(WaitForGasLiftSyncGroups);
200 num_rates_to_sync = comm.sum(num_rates_to_sync);
201 }
202 if (num_rates_to_sync > 0) {
203 OPM_TIMEBLOCK(GasLiftSyncGroups);
204 std::vector<int> group_indexes;
205 group_indexes.reserve(num_rates_to_sync);
206 std::vector<Scalar> group_alq_rates;
207 group_alq_rates.reserve(num_rates_to_sync);
208 std::vector<Scalar> group_oil_rates;
209 group_oil_rates.reserve(num_rates_to_sync);
210 std::vector<Scalar> group_gas_rates;
211 group_gas_rates.reserve(num_rates_to_sync);
212 std::vector<Scalar> group_water_rates;
213 group_water_rates.reserve(num_rates_to_sync);
214 if (comm.rank() == i) {
215 for (auto idx : groups_to_sync) {
216 auto [oil_rate, gas_rate, water_rate, alq] = group_info.getRates(idx);
217 group_indexes.push_back(idx);
218 group_oil_rates.push_back(oil_rate);
219 group_gas_rates.push_back(gas_rate);
220 group_water_rates.push_back(water_rate);
221 group_alq_rates.push_back(alq);
222 }
223 } else {
224 group_indexes.resize(num_rates_to_sync);
225 group_oil_rates.resize(num_rates_to_sync);
226 group_gas_rates.resize(num_rates_to_sync);
227 group_water_rates.resize(num_rates_to_sync);
228 group_alq_rates.resize(num_rates_to_sync);
229 }
230#if HAVE_MPI
231 Parallel::MpiSerializer ser(comm);
232 ser.broadcast(Parallel::RootRank{i}, group_indexes, group_oil_rates,
233 group_gas_rates, group_water_rates, group_alq_rates);
234#endif
235 if (comm.rank() != i) {
236 for (int j = 0; j < num_rates_to_sync; ++j) {
237 group_info.updateRate(group_indexes[j],
238 group_oil_rates[j],
239 group_gas_rates[j],
240 group_water_rates[j],
241 group_alq_rates[j]);
242 }
243 }
244 if constexpr (glift_debug) {
245 int counter = 0;
246 if (comm.rank() == i) {
247 counter = wellState.gliftGetDebugCounter();
248 }
249 counter = comm.sum(counter);
250 if (comm.rank() != i) {
251 wellState.gliftSetDebugCounter(counter);
252 }
253 }
254 }
255 }
256}
257
258// NOTE: this method cannot be const since it passes this->wellState()
259// (see below) to the GasLiftSingleWell constructor which accepts WellState
260// as a non-const reference..
261template<typename TypeTag>
262void
263BlackoilWellModelGasLift<TypeTag>::
264gasLiftOptimizationStage1SingleWell(WellInterface<TypeTag>* well,
265 const Simulator& simulator,
266 WellState<Scalar>& wellState,
267 GroupState<Scalar>& groupState,
268 GLiftProdWells& prod_wells,
269 GLiftOptWells& glift_wells,
270 GasLiftGroupInfo<Scalar>& group_info,
271 GLiftWellStateMap& state_map,
272 GLiftSyncGroups& sync_groups,
273 DeferredLogger& deferred_logger)
274{
275 OPM_TIMEFUNCTION();
276 const auto& summary_state = simulator.vanguard().summaryState();
277 auto glift = std::make_unique<GasLiftSingleWell<TypeTag>>(*well,
278 simulator,
279 summary_state,
280 deferred_logger,
281 wellState,
282 groupState,
283 group_info,
284 sync_groups,
285 simulator.vanguard().gridView().comm(),
286 this->glift_debug);
287 auto state = glift->runOptimize(simulator.model().newtonMethod().numIterations());
288 if (state) {
289 state_map.emplace(well->name(), std::move(state));
290 glift_wells.emplace(well->name(), std::move(glift));
291 return;
292 }
293 prod_wells.insert({well->name(), well});
294}
295
296template<typename TypeTag>
297void
299initGliftEclWellMap(const std::vector<WellInterfacePtr>& well_container,
300 GLiftEclWells& ecl_well_map)
301{
302 for (const auto& well : well_container) {
303 ecl_well_map.try_emplace(well->name(), &well->wellEcl(), well->indexOfWell());
304 }
305}
306
307} // namespace Opm
Class for handling the gaslift in the blackoil well model.
Definition: BlackoilWellModelGasLift.hpp:94
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: BlackoilWellModelGasLift.hpp:105
static void initGliftEclWellMap(const std::vector< WellInterfacePtr > &well_container, GLiftEclWells &ecl_well_map)
Definition: BlackoilWellModelGasLift_impl.hpp:299
bool maybeDoGasLiftOptimize(const Simulator &simulator, const std::vector< WellInterfacePtr > &well_container, WellState< Scalar > &wellState, GroupState< Scalar > &groupState, DeferredLogger &deferred_logger)
Definition: BlackoilWellModelGasLift_impl.hpp:42
typename Base::GLiftWellStateMap GLiftWellStateMap
Definition: BlackoilWellModelGasLift.hpp:104
typename GasLiftGroupInfo< Scalar >::GLiftEclWells GLiftEclWells
Definition: BlackoilWellModelGasLift.hpp:100
typename Base::GLiftOptWells GLiftOptWells
Definition: BlackoilWellModelGasLift.hpp:101
typename Base::GLiftProdWells GLiftProdWells
Definition: BlackoilWellModelGasLift.hpp:102
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: BlackoilWellModelGasLift.hpp:99
Definition: DeferredLogger.hpp:57
Definition: GasLiftGroupInfo.hpp:46
std::tuple< Scalar, Scalar, Scalar, Scalar > getRates(const int group_idx) const
bool hasWell(const std::string &well_name)
void updateRate(int idx, Scalar oil_rate, Scalar gas_rate, Scalar water_rate, Scalar alq)
Definition: GroupState.hpp:43
Manages the initializing and running of time dependent problems.
Definition: simulator.hh:84
Vanguard & vanguard()
Return a reference to the grid manager of simulation.
Definition: simulator.hh:234
Definition: WellState.hpp:65
Definition: blackoilboundaryratevector.hh:39