22#define OPM_BLACKOILWELLMODEL_GASLIFT_IMPL_HEADER_INCLUDED
23#define OPM_BLACKOILWELLMODEL_GASLIFT_IMPL_HEADER_INCLUDED
26#ifndef OPM_BLACKOILWELLMODEL_GASLIFT_HEADER_INCLUDED
28#include <opm/simulators/wells/BlackoilWellModelGasLift.hpp>
30#include <opm/common/TimingMacros.hpp>
31#include <opm/simulators/wells/GasLiftSingleWell.hpp>
32#include <opm/simulators/utils/DeferredLoggingErrorHelpers.hpp>
35#include <opm/simulators/utils/MPISerializer.hpp>
40template<
typename TypeTag>
44 const std::vector<WellInterfacePtr>& well_container,
45 const std::map<std::string, Scalar>& node_pressures,
46 const bool updatePotentials,
47 WellStateType& wellState,
52 const auto& glo = simulator.vanguard().schedule().glo(simulator.episodeIndex());
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();
64 if (simulation_time == this->last_glift_opt_time_ ||
65 simulation_time >= (this->last_glift_opt_time_ + min_wait))
67 do_glift_optimization =
true;
68 this->last_glift_opt_time_ = simulation_time;
72 if (updatePotentials) {
73 updateWellPotentials(simulator, well_container, node_pressures, wellState, deferred_logger);
76 if (do_glift_optimization) {
77 GLiftOptWells glift_wells;
78 GLiftProdWells prod_wells;
79 GLiftWellStateMap state_map;
87 GLiftEclWells ecl_well_map;
88 initGliftEclWellMap(well_container, ecl_well_map);
89 const auto& iterCtx = simulator.problem().iterationContext();
92 simulator.vanguard().schedule(),
93 simulator.vanguard().summaryState(),
94 simulator.episodeIndex(),
99 simulator.vanguard().grid().comm(),
102 group_info.initialize();
104 gasLiftOptimizationStage1(simulator,
114 this->gasLiftOptimizationStage2(simulator.vanguard().gridView().comm(),
115 simulator.vanguard().schedule(),
116 simulator.vanguard().summaryState(),
123 simulator.episodeIndex(),
126 if constexpr (glift_debug) {
127 std::vector<WellInterfaceGeneric<Scalar, IndexTraits>*> wc;
128 wc.reserve(well_container.size());
129 std::ranges::transform(well_container, std::back_inserter(wc),
132 this->gliftDebugShowALQ(wc,
136 num_wells_changed = glift_wells.size();
138 num_wells_changed = simulator.vanguard().gridView().comm().sum(num_wells_changed);
139 return num_wells_changed > 0;
142template<
typename TypeTag>
146 const std::vector<WellInterfacePtr>& well_container,
147 WellStateType& wellState,
149 GLiftProdWells& prod_wells,
150 GLiftOptWells &glift_wells,
152 GLiftWellStateMap& state_map,
156 auto comm = simulator.vanguard().grid().comm();
157 int num_procs = comm.size();
183 for (
int i = 0; i< num_procs; i++) {
184 int num_rates_to_sync = 0;
185 GLiftSyncGroups groups_to_sync;
186 if (comm.rank() == i) {
188 for (
const auto& well : well_container) {
190 if (group_info.hasWell(well->name())) {
191 gasLiftOptimizationStage1SingleWell(well.get(),
203 num_rates_to_sync = groups_to_sync.size();
206 OPM_TIMEBLOCK(WaitForGasLiftSyncGroups);
207 num_rates_to_sync = comm.sum(num_rates_to_sync);
209 if (num_rates_to_sync > 0) {
210 OPM_TIMEBLOCK(GasLiftSyncGroups);
211 std::vector<int> group_indexes;
212 group_indexes.reserve(num_rates_to_sync);
213 std::vector<Scalar> group_alq_rates;
214 group_alq_rates.reserve(num_rates_to_sync);
215 std::vector<Scalar> group_oil_rates;
216 group_oil_rates.reserve(num_rates_to_sync);
217 std::vector<Scalar> group_gas_rates;
218 group_gas_rates.reserve(num_rates_to_sync);
219 std::vector<Scalar> group_water_rates;
220 group_water_rates.reserve(num_rates_to_sync);
221 if (comm.rank() == i) {
222 for (
auto idx : groups_to_sync) {
223 auto [oil_rate, gas_rate, water_rate, alq] = group_info.getRates(idx);
224 group_indexes.push_back(idx);
225 group_oil_rates.push_back(oil_rate);
226 group_gas_rates.push_back(gas_rate);
227 group_water_rates.push_back(water_rate);
228 group_alq_rates.push_back(alq);
231 group_indexes.resize(num_rates_to_sync);
232 group_oil_rates.resize(num_rates_to_sync);
233 group_gas_rates.resize(num_rates_to_sync);
234 group_water_rates.resize(num_rates_to_sync);
235 group_alq_rates.resize(num_rates_to_sync);
240 group_gas_rates, group_water_rates, group_alq_rates);
242 if (comm.rank() != i) {
243 for (
int j = 0; j < num_rates_to_sync; ++j) {
244 group_info.updateRate(group_indexes[j],
247 group_water_rates[j],
258template<
typename TypeTag>
263 WellStateType& wellState,
265 GLiftProdWells& prod_wells,
266 GLiftOptWells& glift_wells,
268 GLiftWellStateMap& state_map,
269 GLiftSyncGroups& sync_groups,
273 const auto& summary_state = simulator.vanguard().summaryState();
274 auto glift = std::make_unique<GasLiftSingleWell<TypeTag>>(*well,
282 simulator.vanguard().gridView().comm(),
284 const auto& iterCtx = simulator.problem().iterationContext();
285 auto state = glift->runOptimize(iterCtx.iteration());
287 state_map.emplace(well->name(), std::move(state));
288 glift_wells.emplace(well->name(), std::move(glift));
291 prod_wells.insert({well->name(), well});
294template<
typename TypeTag>
298 GLiftEclWells& ecl_well_map)
300 for (
const auto& well : well_container) {
301 ecl_well_map.try_emplace(well->name(), &well->wellEcl(), well->indexOfWell());
305template<
typename TypeTag>
309 const std::vector<WellInterfacePtr>& well_container,
310 const std::map<std::string, Scalar>& node_pressures,
311 WellStateType& wellState,
314 auto well_state_copy = wellState;
315 const int np = wellState.numPhases();
317 auto exc_type = ExceptionType::NONE;
318 for (
const auto& well : well_container) {
319 if (well->isInjector() || !well->wellEcl().predictionMode())
322 const auto it = node_pressures.find(well->wellEcl().groupName());
323 if (it != node_pressures.end()) {
324 std::string cur_exc_msg;
325 auto cur_exc_type = ExceptionType::NONE;
327 std::vector<Scalar> potentials;
328 const auto& groupStateHelper = simulator.problem().wellModel().groupStateHelper();
329 well->computeWellPotentials(simulator, well_state_copy, groupStateHelper, potentials);
330 auto& ws = wellState.well(well->indexOfWell());
331 for (
int p = 0; p < np; ++p) {
333 ws.well_potentials[p] = std::max(Scalar{0.0}, potentials[p]);
337 OPM_PARALLEL_CATCH_CLAUSE(cur_exc_type, cur_exc_msg);
338 if (cur_exc_type != ExceptionType::NONE) {
339 exc_msg += fmt::format(
"\nFor well {}: {}", well->name(), cur_exc_msg);
341 exc_type = std::max(exc_type, cur_exc_type);
344 if (exc_type != ExceptionType::NONE) {
345 const std::string msg =
"Updating well potentials after network balancing failed. Continue with current potentials";
346 deferred_logger.warning(
"WELL_POT_SOLVE_AFTER_NETWORK_FAILED", msg + exc_msg);
Class for handling the gaslift in the blackoil well model.
Definition BlackoilWellModelGasLift.hpp:96
Definition DeferredLogger.hpp:57
Definition GasLiftGroupInfo.hpp:47
Definition GroupState.hpp:41
Class for serializing and broadcasting data using MPI.
Definition MPISerializer.hpp:38
Manages the initializing and running of time dependent problems.
Definition simulator.hh:84
Definition WellInterfaceGeneric.hpp:53
Definition WellInterface.hpp:77
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:45
Avoid mistakes in calls to broadcast() by wrapping the root argument in an explicit type.
Definition MPISerializer.hpp:33