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>
33
34#if HAVE_MPI
36#endif
37
38namespace Opm {
39
40template<typename TypeTag>
41bool
43maybeDoGasLiftOptimize(const Simulator& simulator,
44 const std::vector<WellInterfacePtr>& well_container,
45 const std::map<std::string, Scalar>& node_pressures,
46 const bool updatePotentials,
47 WellStateType& wellState,
48 GroupState<Scalar>& groupState,
49 DeferredLogger& deferred_logger)
50{
51 OPM_TIMEFUNCTION();
52 const auto& glo = simulator.vanguard().schedule().glo(simulator.episodeIndex());
53 if(!glo.active()) {
54 return false;
55 }
56 bool do_glift_optimization = false;
57 int num_wells_changed = 0;
58 const double simulation_time = simulator.time();
59 const Scalar min_wait = glo.min_wait();
60 // We only optimize if a min_wait time has past.
61 // If all_newton is true we still want to optimize several times pr timestep
62 // i.e. we also optimize if check simulation_time == last_glift_opt_time_
63 // that is when the last_glift_opt_time is already updated with the current time step
64 if (simulation_time == this->last_glift_opt_time_ ||
65 simulation_time >= (this->last_glift_opt_time_ + min_wait))
66 {
67 do_glift_optimization = true;
68 this->last_glift_opt_time_ = simulation_time;
69 }
70
71 // update potentials if changed by network (needed by the gaslift algorithm)
72 if (updatePotentials) {
73 updateWellPotentials(simulator, well_container, node_pressures, wellState, deferred_logger);
74 }
75
76 if (do_glift_optimization) {
77 GLiftOptWells glift_wells;
78 GLiftProdWells prod_wells;
79 GLiftWellStateMap state_map;
80 // NOTE: To make GasLiftGroupInfo (see below) independent of the TypeTag
81 // associated with *this (i.e. BlackoilWellModel<TypeTag>) we observe
82 // that GasLiftGroupInfo's only dependence on *this is that it needs to
83 // access the eclipse Wells in the well container (the eclipse Wells
84 // themselves are independent of the TypeTag).
85 // Hence, we extract them from the well container such that we can pass
86 // them to the GasLiftGroupInfo constructor.
87 GLiftEclWells ecl_well_map;
88 initGliftEclWellMap(well_container, ecl_well_map);
89 GasLiftGroupInfo group_info {
90 ecl_well_map,
91 simulator.vanguard().schedule(),
92 simulator.vanguard().summaryState(),
93 simulator.episodeIndex(),
94 simulator.model().newtonMethod().numIterations(),
95 deferred_logger,
96 wellState,
97 groupState,
98 simulator.vanguard().grid().comm(),
99 this->glift_debug
100 };
101 group_info.initialize();
102
103 gasLiftOptimizationStage1(simulator,
104 well_container,
105 wellState,
106 groupState,
107 prod_wells,
108 glift_wells,
109 group_info,
110 state_map,
111 deferred_logger);
112
113 this->gasLiftOptimizationStage2(simulator.vanguard().gridView().comm(),
114 simulator.vanguard().schedule(),
115 simulator.vanguard().summaryState(),
116 wellState,
117 groupState,
118 prod_wells,
119 glift_wells,
120 group_info,
121 state_map,
122 simulator.episodeIndex(),
123 deferred_logger);
124
125 if constexpr (glift_debug) {
126 std::vector<WellInterfaceGeneric<Scalar, IndexTraits>*> wc;
127 wc.reserve(well_container.size());
128 std::ranges::transform(well_container, std::back_inserter(wc),
129 [](const auto& w)
130 { return static_cast<WellInterfaceGeneric<Scalar, IndexTraits>*>(w.get()); });
131 this->gliftDebugShowALQ(wc,
132 wellState,
133 deferred_logger);
134 }
135 num_wells_changed = glift_wells.size();
136 }
137 num_wells_changed = simulator.vanguard().gridView().comm().sum(num_wells_changed);
138 return num_wells_changed > 0;
139}
140
141template<typename TypeTag>
142void
144gasLiftOptimizationStage1(const Simulator& simulator,
145 const std::vector<WellInterfacePtr>& well_container,
146 WellStateType& wellState,
147 GroupState<Scalar>& groupState,
148 GLiftProdWells& prod_wells,
149 GLiftOptWells &glift_wells,
151 GLiftWellStateMap& state_map,
152 DeferredLogger& deferred_logger)
153{
154 OPM_TIMEFUNCTION();
155 auto comm = simulator.vanguard().grid().comm();
156 int num_procs = comm.size();
157 // NOTE: Gas lift optimization stage 1 seems to be difficult
158 // to do in parallel since the wells are optimized on different
159 // processes and each process needs to know the current ALQ allocated
160 // to each group it is a memeber of in order to check group limits and avoid
161 // allocating more ALQ than necessary. (Surplus ALQ is removed in
162 // stage 2). In stage1, as each well is adding ALQ, the current group ALQ needs
163 // to be communicated to the other processes. But there is no common
164 // synchronization point that all process will reach in the
165 // runOptimizeLoop_() in GasLiftSingleWell.cpp.
166 //
167 // TODO: Maybe a better solution could be invented by distributing
168 // wells according to certain parent groups. Then updated group rates
169 // might not have to be communicated to the other processors.
170
171 // Currently, the best option seems to be to run this part sequentially
172 // (not in parallel).
173 //
174 // TODO: The simplest approach seems to be if a) one process could take
175 // ownership of all the wells (the union of all the wells in the
176 // well_container_ of each process) then this process could do the
177 // optimization, while the other processes could wait for it to
178 // finish (e.g. comm.barrier()), or alternatively, b) if all
179 // processes could take ownership of all the wells. Then there
180 // would be no need for synchronization here..
181 //
182 for (int i = 0; i< num_procs; i++) {
183 int num_rates_to_sync = 0; // communication variable
184 GLiftSyncGroups groups_to_sync;
185 if (comm.rank() == i) {
186 // Run stage1: Optimize single wells while also checking group limits
187 for (const auto& well : well_container) {
188 // NOTE: Only the wells in "group_info" needs to be optimized
189 if (group_info.hasWell(well->name())) {
190 gasLiftOptimizationStage1SingleWell(well.get(),
191 simulator,
192 wellState,
193 groupState,
194 prod_wells,
195 glift_wells,
196 group_info,
197 state_map,
198 groups_to_sync,
199 deferred_logger);
200 }
201 }
202 num_rates_to_sync = groups_to_sync.size();
203 }
204 {
205 OPM_TIMEBLOCK(WaitForGasLiftSyncGroups);
206 num_rates_to_sync = comm.sum(num_rates_to_sync);
207 }
208 if (num_rates_to_sync > 0) {
209 OPM_TIMEBLOCK(GasLiftSyncGroups);
210 std::vector<int> group_indexes;
211 group_indexes.reserve(num_rates_to_sync);
212 std::vector<Scalar> group_alq_rates;
213 group_alq_rates.reserve(num_rates_to_sync);
214 std::vector<Scalar> group_oil_rates;
215 group_oil_rates.reserve(num_rates_to_sync);
216 std::vector<Scalar> group_gas_rates;
217 group_gas_rates.reserve(num_rates_to_sync);
218 std::vector<Scalar> group_water_rates;
219 group_water_rates.reserve(num_rates_to_sync);
220 if (comm.rank() == i) {
221 for (auto idx : groups_to_sync) {
222 auto [oil_rate, gas_rate, water_rate, alq] = group_info.getRates(idx);
223 group_indexes.push_back(idx);
224 group_oil_rates.push_back(oil_rate);
225 group_gas_rates.push_back(gas_rate);
226 group_water_rates.push_back(water_rate);
227 group_alq_rates.push_back(alq);
228 }
229 } else {
230 group_indexes.resize(num_rates_to_sync);
231 group_oil_rates.resize(num_rates_to_sync);
232 group_gas_rates.resize(num_rates_to_sync);
233 group_water_rates.resize(num_rates_to_sync);
234 group_alq_rates.resize(num_rates_to_sync);
235 }
236#if HAVE_MPI
237 Parallel::MpiSerializer ser(comm);
238 ser.broadcast(Parallel::RootRank{i}, group_indexes, group_oil_rates,
239 group_gas_rates, group_water_rates, group_alq_rates);
240#endif
241 if (comm.rank() != i) {
242 for (int j = 0; j < num_rates_to_sync; ++j) {
243 group_info.updateRate(group_indexes[j],
244 group_oil_rates[j],
245 group_gas_rates[j],
246 group_water_rates[j],
247 group_alq_rates[j]);
248 }
249 }
250 }
251 }
252}
253
254// NOTE: this method cannot be const since it passes this->wellState()
255// (see below) to the GasLiftSingleWell constructor which accepts WellState
256// as a non-const reference..
257template<typename TypeTag>
258void
259BlackoilWellModelGasLift<TypeTag>::
260gasLiftOptimizationStage1SingleWell(WellInterface<TypeTag>* well,
261 const Simulator& simulator,
262 WellStateType& wellState,
263 GroupState<Scalar>& groupState,
264 GLiftProdWells& prod_wells,
265 GLiftOptWells& glift_wells,
266 GasLiftGroupInfo<Scalar, IndexTraits>& group_info,
267 GLiftWellStateMap& state_map,
268 GLiftSyncGroups& sync_groups,
269 DeferredLogger& deferred_logger)
270{
271 OPM_TIMEFUNCTION();
272 const auto& summary_state = simulator.vanguard().summaryState();
273 auto glift = std::make_unique<GasLiftSingleWell<TypeTag>>(*well,
274 simulator,
275 summary_state,
276 deferred_logger,
277 wellState,
278 groupState,
279 group_info,
280 sync_groups,
281 simulator.vanguard().gridView().comm(),
282 this->glift_debug);
283 auto state = glift->runOptimize(simulator.model().newtonMethod().numIterations());
284 if (state) {
285 state_map.emplace(well->name(), std::move(state));
286 glift_wells.emplace(well->name(), std::move(glift));
287 return;
288 }
289 prod_wells.insert({well->name(), well});
290}
291
292template<typename TypeTag>
293void
295initGliftEclWellMap(const std::vector<WellInterfacePtr>& well_container,
296 GLiftEclWells& ecl_well_map)
297{
298 for (const auto& well : well_container) {
299 ecl_well_map.try_emplace(well->name(), &well->wellEcl(), well->indexOfWell());
300 }
301}
302
303template<typename TypeTag>
304void
306updateWellPotentials(const Simulator& simulator,
307 const std::vector<WellInterfacePtr>& well_container,
308 const std::map<std::string, Scalar>& node_pressures,
309 WellStateType& wellState,
310 DeferredLogger& deferred_logger)
311{
312 auto well_state_copy = wellState;
313 const int np = wellState.numPhases();
314 std::string exc_msg;
315 auto exc_type = ExceptionType::NONE;
316 for (const auto& well : well_container) {
317 if (well->isInjector() || !well->wellEcl().predictionMode())
318 continue;
319
320 const auto it = node_pressures.find(well->wellEcl().groupName());
321 if (it != node_pressures.end()) {
322 std::string cur_exc_msg;
323 auto cur_exc_type = ExceptionType::NONE;
324 try {
325 std::vector<Scalar> potentials;
326 const auto& groupStateHelper = simulator.problem().wellModel().groupStateHelper();
327 well->computeWellPotentials(simulator, well_state_copy, groupStateHelper, potentials);
328 auto& ws = wellState.well(well->indexOfWell());
329 for (int p = 0; p < np; ++p) {
330 // make sure the potentials are positive
331 ws.well_potentials[p] = std::max(Scalar{0.0}, potentials[p]);
332 }
333 }
334 // catch all possible exception and store type and message.
335 OPM_PARALLEL_CATCH_CLAUSE(cur_exc_type, cur_exc_msg);
336 if (cur_exc_type != ExceptionType::NONE) {
337 exc_msg += fmt::format("\nFor well {}: {}", well->name(), cur_exc_msg);
338 }
339 exc_type = std::max(exc_type, cur_exc_type);
340 }
341 }
342 if (exc_type != ExceptionType::NONE) {
343 const std::string msg = "Updating well potentials after network balancing failed. Continue with current potentials";
344 deferred_logger.warning("WELL_POT_SOLVE_AFTER_NETWORK_FAILED", msg + exc_msg);
345 }
346}
347
348} // namespace Opm
#define OPM_PARALLEL_CATCH_CLAUSE(obptc_exc_type, obptc_exc_msg)
Inserts catch classes for the parallel try-catch.
Definition: DeferredLoggingErrorHelpers.hpp:166
Class for handling the gaslift in the blackoil well model.
Definition: BlackoilWellModelGasLift.hpp:96
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: BlackoilWellModelGasLift.hpp:108
bool maybeDoGasLiftOptimize(const Simulator &simulator, const std::vector< WellInterfacePtr > &well_container, const std::map< std::string, Scalar > &node_pressures, const bool updatePotentials, WellStateType &wellState, GroupState< Scalar > &groupState, DeferredLogger &deferred_logger)
Definition: BlackoilWellModelGasLift_impl.hpp:43
static void initGliftEclWellMap(const std::vector< WellInterfacePtr > &well_container, GLiftEclWells &ecl_well_map)
Definition: BlackoilWellModelGasLift_impl.hpp:295
typename Base::GLiftWellStateMap GLiftWellStateMap
Definition: BlackoilWellModelGasLift.hpp:107
typename Base::GLiftOptWells GLiftOptWells
Definition: BlackoilWellModelGasLift.hpp:104
typename GasLiftGroupInfo< Scalar, IndexTraits >::GLiftEclWells GLiftEclWells
Definition: BlackoilWellModelGasLift.hpp:103
typename Base::GLiftProdWells GLiftProdWells
Definition: BlackoilWellModelGasLift.hpp:105
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: BlackoilWellModelGasLift.hpp:98
Definition: DeferredLogger.hpp:57
void warning(const std::string &tag, const std::string &message)
Definition: GasLiftGroupInfo.hpp:46
bool hasWell(const std::string &well_name)
std::tuple< Scalar, Scalar, Scalar, Scalar > getRates(const int group_idx) const
void updateRate(int idx, Scalar oil_rate, Scalar gas_rate, Scalar water_rate, Scalar alq)
Definition: GroupState.hpp:41
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
Problem & problem()
Return the object which specifies the pysical setup of the simulation.
Definition: simulator.hh:265
Definition: WellInterfaceGeneric.hpp:53
Definition: WellState.hpp:66
@ NONE
Definition: DeferredLogger.hpp:46
Definition: blackoilbioeffectsmodules.hh:43