WellGroupHelper.hpp
Go to the documentation of this file.
1/*
2 Copyright 2025 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#ifndef OPM_WELLGROUPHELPER_HEADER_INCLUDED
20#define OPM_WELLGROUPHELPER_HEADER_INCLUDED
21
22#include <opm/common/TimingMacros.hpp>
23#include <opm/input/eclipse/EclipseState/Grid/FieldPropsManager.hpp>
24#include <opm/input/eclipse/Schedule/Group/GPMaint.hpp>
25#include <opm/input/eclipse/Schedule/Group/GSatProd.hpp>
26#include <opm/input/eclipse/Schedule/Group/GuideRate.hpp>
27#include <opm/input/eclipse/Schedule/Schedule.hpp>
28#include <opm/input/eclipse/Schedule/ScheduleState.hpp>
29#include <opm/input/eclipse/Schedule/SummaryState.hpp>
30#include <opm/material/fluidsystems/PhaseUsageInfo.hpp>
35
37
38#include <algorithm>
39#include <map>
40#include <memory>
41#include <optional>
42#include <stdexcept>
43#include <string>
44#include <vector>
45
46namespace Opm
47{
48
49template <typename Scalar, typename IndexTraits>
51{
52public:
53 // RAII guard for temporarily setting wellstate pointer
55 {
56 public:
58 : wgHelper_ {wgHelper}
59 , previous_state_ptr_ {wgHelper_.well_state_}
60 {
61 // Set the new state directly
62 wgHelper_.well_state_ = &well_state;
63 }
64
66 {
67 // Restore the previous state
68 wgHelper_.well_state_ = previous_state_ptr_;
69 }
70
71 // Delete copy and move operations
76
77 private:
78 WellGroupHelper& wgHelper_;
79 const WellState<Scalar, IndexTraits>* previous_state_ptr_;
80 };
81
82 // RAII guard for temporarily setting groupstate pointer
84 {
85 public:
87 : wgHelper_ {wgHelper}
88 , previous_state_ptr_ {wgHelper_.group_state_}
89 {
90 // Set the new state directly
91 wgHelper_.group_state_ = &group_state;
92 }
93
95 {
96 // Restore the previous state
97 wgHelper_.group_state_ = previous_state_ptr_;
98 }
99
100 // Delete copy and move operations
105
106 private:
107 WellGroupHelper& wgHelper_;
108 const GroupState<Scalar>* previous_state_ptr_;
109 };
110
112 GroupState<Scalar>& group_state,
113 const Schedule& schedule,
114 const SummaryState& summary_state,
115 const GuideRate& guide_rate,
116 const PhaseUsageInfo<IndexTraits>& phase_usage_info);
117
118 void accumulateGroupEfficiencyFactor(const Group& group, Scalar& factor) const;
119
120 std::pair<bool, Scalar> checkGroupConstraintsInj(const std::string& name,
121 const std::string& parent,
122 const Group& group,
123 const Scalar* rates,
124 const Phase injection_phase,
125 const Scalar efficiency_factor,
126 const std::vector<Scalar>& resv_coeff,
127 const bool check_guide_rate,
128 DeferredLogger& deferred_logger) const;
129
130 std::pair<bool, Scalar> checkGroupConstraintsProd(const std::string& name,
131 const std::string& parent,
132 const Group& group,
133 const Scalar* rates,
134 const Scalar efficiency_factor,
135 const std::vector<Scalar>& resv_coeff,
136 const bool check_guide_rate,
137 DeferredLogger& deferred_logger) const;
138
139 std::map<std::string, Scalar> computeNetworkPressures(const Network::ExtNetwork& network,
140 const VFPProdProperties<Scalar>& vfp_prod_props,
141 const Parallel::Communication& comm) const;
142
143 Scalar getGuideRate(const std::string& name, const GuideRateModel::Target target) const;
144
145 GuideRate::RateVector getProductionGroupRateVector(const std::string& group_name) const;
146
147 Scalar getWellGroupTargetInjector(const std::string& name,
148 const std::string& parent,
149 const Group& group,
150 const Scalar* rates,
151 const Phase injection_phase,
152 const Scalar efficiency_factor,
153 const std::vector<Scalar>& resv_coeff,
154 DeferredLogger& deferred_logger) const;
155
156 Scalar getWellGroupTargetProducer(const std::string& name,
157 const std::string& parent,
158 const Group& group,
159 const Scalar* rates,
160 const Scalar efficiency_factor,
161 const std::vector<Scalar>& resv_coeff,
162 DeferredLogger& deferred_logger) const;
163
164 GuideRate::RateVector getWellRateVector(const std::string& name) const;
165
166 std::vector<std::string> groupChainTopBot(const std::string& bottom, const std::string& top) const;
167
170 int groupControlledWells(const std::string& group_name,
171 const std::string& always_included_child,
172 const bool is_production_group,
173 const Phase injection_phase) const;
174
176 {
177 return *this->group_state_;
178 }
179
181 {
182 return this->phase_usage_info_;
183 }
184
186 {
187 return GroupStateGuard(*this, group_state);
188 }
189
191 {
192 return WellStateGuard(*this, well_state);
193 }
194
195 void setCmodeGroup(const Group& group, GroupState<Scalar>& group_state) const;
196
197 template <class AverageRegionalPressureType>
198 void
200 const FieldPropsManager& fp,
201 std::map<std::string, std::unique_ptr<AverageRegionalPressureType>>&
202 regional_average_pressure_calculator) const;
203
204 void setReportStep(int report_step)
205 {
206 report_step_ = report_step;
207 }
208
209 Scalar sumSolventRates(const Group& group, const bool is_injector) const;
210
211 Scalar sumWellResRates(const Group& group, const int phase_pos, const bool injector) const;
212
213 Scalar sumWellSurfaceRates(const Group& group, const int phase_pos, const bool injector) const;
214
215 Scalar sumWellPhaseRates(bool res_rates,
216 const Group& group,
217 const int phase_pos,
218 const bool injector,
219 const bool network = false) const;
220
221 template <class RegionalValues>
222 void updateGpMaintTargetForGroups(const Group& group,
223 const RegionalValues& regional_values,
224 const double dt,
225 GroupState<Scalar>& group_state) const;
226
229 int updateGroupControlledWells(const bool is_production_group,
230 const Phase injection_phase,
231 GroupState<Scalar>& group_state,
232 DeferredLogger& deferred_logger) const;
233
234 void updateGroupProductionRates(const Group& group, GroupState<Scalar>& group_state) const;
235
236 void updateGroupTargetReduction(const Group& group,
237 const bool is_injector,
238 GroupState<Scalar>& group_state) const;
239
241
242 void updateREINForGroups(const Group& group, bool sum_rank, GroupState<Scalar>& group_state) const;
243
244 void updateReservoirRatesInjectionGroups(const Group& group, GroupState<Scalar>& group_state) const;
245
246 void updateVREPForGroups(const Group& group, GroupState<Scalar>& group_state) const;
247
249
250 void updateSurfaceRatesInjectionGroups(const Group& group, GroupState<Scalar>& group_state) const;
251
252 void updateWellRates(const Group& group,
253 const WellState<Scalar, IndexTraits>& well_state_nupcol,
254 WellState<Scalar, IndexTraits>& well_state) const;
255
257 {
258 return *this->well_state_;
259 }
260
261 void updateWellRatesFromGroupTargetScale(const Scalar scale,
262 const Group& group,
263 bool is_injector,
264 WellState<Scalar, IndexTraits>& well_state) const;
265
267 std::pair<std::optional<std::string>, Scalar>
268 worstOffendingWell(const Group& group,
269 const Group::ProductionCMode& offended_control,
270 const Parallel::Communication& comm,
271 DeferredLogger& deferred_logger) const;
272
273private:
274 std::string controlGroup_(const Group& group) const;
275
276 GuideRate::RateVector getGuideRateVector_(const std::vector<Scalar>& rates) const;
277
279 bool isInGroupChainTopBot_(const std::string& bottom, const std::string& top) const;
280
281 int phaseToActivePhaseIdx_(const Phase phase) const;
282
283 Scalar satelliteInjectionRate_(const ScheduleState& sched,
284 const Group& group,
285 const int phase_pos,
286 bool res_rates) const;
287
288 Scalar satelliteProductionRate_(const ScheduleState& sched,
289 const Group& group,
290 const GSatProd::GSatProdGroupProp::Rate rate_comp,
291 bool res_rates) const;
292
293 std::optional<GSatProd::GSatProdGroupProp::Rate> selectRateComponent_(const int phase_pos) const;
294
295 int updateGroupControlledWellsRecursive_(const std::string& group_name,
296 const bool is_production_group,
297 const Phase injection_phase,
298 GroupState<Scalar>& group_state,
299 DeferredLogger& deferred_logger) const;
300
301 void updateGroupTargetReductionRecursive_(const Group& group,
302 const bool is_injector,
303 std::vector<Scalar>& group_target_reduction,
304 GroupState<Scalar>& group_state) const;
305
306 const WellState<Scalar, IndexTraits>* well_state_ {nullptr};
307 const GroupState<Scalar>* group_state_ {nullptr};
308 const Schedule& schedule_;
309 const SummaryState& summary_state_;
310 const GuideRate& guide_rate_;
311 int report_step_ {0};
312 // NOTE: The phase usage info seems to be read-only throughout the simulation, so it should be safe
313 // to store a reference to it here.
314 const PhaseUsageInfo<IndexTraits>& phase_usage_info_;
315};
316
317// -----------------------------------------------------------------------------
318// Template member function implementations
319// -----------------------------------------------------------------------------
320
321// NOTE: This template member function is defined here in the header because the
322// AverageRegionalPressureType type parameter depends on derived class types that are
323// not available when WellGroupHelper.cpp is compiled. Template functions
324// must be visible at their instantiation point.
325// See WellGroupHelper.cpp for detailed rationale.
326template <typename Scalar, typename IndexTraits>
327template <class AverageRegionalPressureType>
328void
330 const Group& group,
331 const FieldPropsManager& fp,
332 std::map<std::string, std::unique_ptr<AverageRegionalPressureType>>& regional_average_pressure_calculator)
333 const
334{
335 for (const std::string& groupName : group.groups()) {
336 this->setRegionAveragePressureCalculator(this->schedule_.getGroup(groupName, this->report_step_),
337 fp,
338 regional_average_pressure_calculator);
339 }
340 const auto& gpm = group.gpmaint();
341 if (!gpm)
342 return;
343
344 const auto& reg = gpm->region();
345 if (!reg)
346 return;
347
348 if (regional_average_pressure_calculator.count(reg->first) == 0) {
349 const std::string name = (reg->first.rfind("FIP", 0) == 0) ? reg->first : "FIP" + reg->first;
350 const auto& fipnum = fp.get_int(name);
351 regional_average_pressure_calculator[reg->first]
352 = std::make_unique<AverageRegionalPressureType>(fipnum);
353 }
354}
355
356// NOTE: This template member function is defined in the header because the
357// RegionalValues type parameter depends on derived class types that are
358// not available when WellGroupHelper.cpp is compiled. Template functions
359// must be visible at their instantiation point.
360// See WellGroupHelper.cpp for detailed rationale.
361template <typename Scalar, typename IndexTraits>
362template <class RegionalValues>
363void
365 const RegionalValues& regional_values,
366 const double dt,
367 GroupState<Scalar>& group_state) const
368{
369 OPM_TIMEFUNCTION();
370 for (const std::string& group_name : group.groups()) {
371 const Group& group_tmp = this->schedule_.getGroup(group_name, this->report_step_);
372 this->updateGpMaintTargetForGroups(group_tmp, regional_values, dt, group_state);
373 }
374 const auto& gpm = group.gpmaint();
375 if (!gpm)
376 return;
377
378 const auto& region = gpm->region();
379 if (!region)
380 return;
381
382 const auto [name, number] = *region;
383 const Scalar error = gpm->pressure_target() - regional_values.at(name)->pressure(number);
384 Scalar current_rate = 0.0;
385 const auto& pu = this->phase_usage_info_;
386 bool injection = true;
387 Scalar sign = 1.0;
388 switch (gpm->flow_target()) {
389 case GPMaint::FlowTarget::RESV_PROD: {
390 current_rate = -this->groupState().injection_vrep_rate(group.name());
391 injection = false;
392 sign = -1.0;
393 break;
394 }
395 case GPMaint::FlowTarget::RESV_OINJ: {
396 if (pu.phaseIsActive(IndexTraits::oilPhaseIdx)) {
397 const auto io = pu.canonicalToActivePhaseIdx(IndexTraits::oilPhaseIdx);
398 current_rate = group_state.injection_reservoir_rates(group.name())[io];
399 }
400 break;
401 }
402 case GPMaint::FlowTarget::RESV_WINJ: {
403 if (pu.phaseIsActive(IndexTraits::waterPhaseIdx)) {
404 const auto iw = pu.canonicalToActivePhaseIdx(IndexTraits::waterPhaseIdx);
405 current_rate = group_state.injection_reservoir_rates(group.name())[iw];
406 }
407 break;
408 }
409 case GPMaint::FlowTarget::RESV_GINJ: {
410 if (pu.phaseIsActive(IndexTraits::gasPhaseIdx)) {
411 const auto ig = pu.canonicalToActivePhaseIdx(IndexTraits::gasPhaseIdx);
412 current_rate = group_state.injection_reservoir_rates(group.name())[ig];
413 }
414 break;
415 }
416 case GPMaint::FlowTarget::SURF_OINJ: {
417 if (pu.phaseIsActive(IndexTraits::oilPhaseIdx)) {
418 const auto io = pu.canonicalToActivePhaseIdx(IndexTraits::oilPhaseIdx);
419 current_rate = group_state.injection_surface_rates(group.name())[io];
420 }
421 break;
422 }
423 case GPMaint::FlowTarget::SURF_WINJ: {
424 if (pu.phaseIsActive(IndexTraits::waterPhaseIdx)) {
425 const auto iw = pu.canonicalToActivePhaseIdx(IndexTraits::waterPhaseIdx);
426 current_rate = group_state.injection_surface_rates(group.name())[iw];
427 }
428 break;
429 }
430 case GPMaint::FlowTarget::SURF_GINJ: {
431 if (pu.phaseIsActive(IndexTraits::gasPhaseIdx)) {
432 const auto ig = pu.canonicalToActivePhaseIdx(IndexTraits::gasPhaseIdx);
433 current_rate = group_state.injection_surface_rates(group.name())[ig];
434 }
435 break;
436 }
437 default:
438 throw std::invalid_argument("Invalid Flow target type in GPMAINT");
439 }
440 auto& gpmaint_state = group_state.gpmaint(group.name());
441 // we only activate gpmaint if pressure is lower than the target regional pressure for injectors
442 // (i.e. error > 0) and higher for producers.
443 bool activate = (injection && error > 0) || (!injection && error < 0);
444 Scalar rate = 0.0;
445 if (activate) {
446 rate = gpm->rate(gpmaint_state, current_rate, error, dt);
447 } else {
448 gpm->resetState(gpmaint_state);
449 }
450 group_state.update_gpmaint_target(group.name(), std::max(Scalar {0.0}, sign * rate));
451}
452
453} // namespace Opm
454
455#endif
Definition: DeferredLogger.hpp:57
Definition: GroupState.hpp:41
const std::vector< Scalar > & injection_reservoir_rates(const std::string &gname) const
const std::vector< Scalar > & injection_surface_rates(const std::string &gname) const
GPMaint::State & gpmaint(const std::string &gname)
void update_gpmaint_target(const std::string &gname, Scalar target)
Definition: GasLiftGroupInfo.hpp:37
Definition: VFPProdProperties.hpp:38
Definition: WellGroupHelper.hpp:84
GroupStateGuard & operator=(GroupStateGuard &&)=delete
GroupStateGuard(WellGroupHelper &wgHelper, GroupState< Scalar > &group_state)
Definition: WellGroupHelper.hpp:86
GroupStateGuard & operator=(const GroupStateGuard &)=delete
GroupStateGuard(const GroupStateGuard &)=delete
GroupStateGuard(GroupStateGuard &&)=delete
~GroupStateGuard()
Definition: WellGroupHelper.hpp:94
Definition: WellGroupHelper.hpp:55
WellStateGuard(WellGroupHelper &wgHelper, WellState< Scalar, IndexTraits > &well_state)
Definition: WellGroupHelper.hpp:57
WellStateGuard(WellStateGuard &&)=delete
~WellStateGuard()
Definition: WellGroupHelper.hpp:65
WellStateGuard(const WellStateGuard &)=delete
WellStateGuard & operator=(const WellStateGuard &)=delete
WellStateGuard & operator=(WellStateGuard &&)=delete
Definition: WellGroupHelper.hpp:51
const GroupState< Scalar > & groupState() const
Definition: WellGroupHelper.hpp:175
WellGroupHelper(WellState< Scalar, IndexTraits > &well_state, GroupState< Scalar > &group_state, const Schedule &schedule, const SummaryState &summary_state, const GuideRate &guide_rate, const PhaseUsageInfo< IndexTraits > &phase_usage_info)
Scalar sumWellPhaseRates(bool res_rates, const Group &group, const int phase_pos, const bool injector, const bool network=false) const
void updateVREPForGroups(const Group &group, GroupState< Scalar > &group_state) const
Scalar getWellGroupTargetProducer(const std::string &name, const std::string &parent, const Group &group, const Scalar *rates, const Scalar efficiency_factor, const std::vector< Scalar > &resv_coeff, DeferredLogger &deferred_logger) const
std::vector< std::string > groupChainTopBot(const std::string &bottom, const std::string &top) const
Scalar sumSolventRates(const Group &group, const bool is_injector) const
std::pair< bool, Scalar > checkGroupConstraintsProd(const std::string &name, const std::string &parent, const Group &group, const Scalar *rates, const Scalar efficiency_factor, const std::vector< Scalar > &resv_coeff, const bool check_guide_rate, DeferredLogger &deferred_logger) const
Scalar sumWellSurfaceRates(const Group &group, const int phase_pos, const bool injector) const
void updateGroupTargetReduction(const Group &group, const bool is_injector, GroupState< Scalar > &group_state) const
void updateREINForGroups(const Group &group, bool sum_rank, GroupState< Scalar > &group_state) const
GuideRate::RateVector getProductionGroupRateVector(const std::string &group_name) const
int groupControlledWells(const std::string &group_name, const std::string &always_included_child, const bool is_production_group, const Phase injection_phase) const
void setReportStep(int report_step)
Definition: WellGroupHelper.hpp:204
void updateWellRatesFromGroupTargetScale(const Scalar scale, const Group &group, bool is_injector, WellState< Scalar, IndexTraits > &well_state) const
int updateGroupControlledWells(const bool is_production_group, const Phase injection_phase, GroupState< Scalar > &group_state, DeferredLogger &deferred_logger) const
const PhaseUsageInfo< IndexTraits > & phaseUsageInfo() const
Definition: WellGroupHelper.hpp:180
void setRegionAveragePressureCalculator(const Group &group, const FieldPropsManager &fp, std::map< std::string, std::unique_ptr< AverageRegionalPressureType > > &regional_average_pressure_calculator) const
Definition: WellGroupHelper.hpp:329
void updateGroupProductionRates(const Group &group, GroupState< Scalar > &group_state) const
Scalar sumWellResRates(const Group &group, const int phase_pos, const bool injector) const
void updateWellRates(const Group &group, const WellState< Scalar, IndexTraits > &well_state_nupcol, WellState< Scalar, IndexTraits > &well_state) const
void setCmodeGroup(const Group &group, GroupState< Scalar > &group_state) const
GuideRate::RateVector getWellRateVector(const std::string &name) const
void accumulateGroupEfficiencyFactor(const Group &group, Scalar &factor) const
GroupStateGuard pushGroupState(GroupState< Scalar > &group_state)
Definition: WellGroupHelper.hpp:185
void updateNetworkLeafNodeProductionRates(GroupState< Scalar > &group_state) const
std::map< std::string, Scalar > computeNetworkPressures(const Network::ExtNetwork &network, const VFPProdProperties< Scalar > &vfp_prod_props, const Parallel::Communication &comm) const
std::pair< bool, Scalar > checkGroupConstraintsInj(const std::string &name, const std::string &parent, const Group &group, const Scalar *rates, const Phase injection_phase, const Scalar efficiency_factor, const std::vector< Scalar > &resv_coeff, const bool check_guide_rate, DeferredLogger &deferred_logger) const
const WellState< Scalar, IndexTraits > & wellState() const
Definition: WellGroupHelper.hpp:256
Scalar getGuideRate(const std::string &name, const GuideRateModel::Target target) const
std::pair< std::optional< std::string >, Scalar > worstOffendingWell(const Group &group, const Group::ProductionCMode &offended_control, const Parallel::Communication &comm, DeferredLogger &deferred_logger) const
Returns the name of the worst offending well and its fraction (i.e. violated_phase / preferred_phase)
void updateState(WellState< Scalar, IndexTraits > &well_state, GroupState< Scalar > &group_state)
WellStateGuard pushWellState(WellState< Scalar, IndexTraits > &well_state)
Definition: WellGroupHelper.hpp:190
void updateGpMaintTargetForGroups(const Group &group, const RegionalValues &regional_values, const double dt, GroupState< Scalar > &group_state) const
Definition: WellGroupHelper.hpp:364
void updateReservoirRatesInjectionGroups(const Group &group, GroupState< Scalar > &group_state) const
void updateSurfaceRatesInjectionGroups(const Group &group, GroupState< Scalar > &group_state) const
Scalar getWellGroupTargetInjector(const std::string &name, const std::string &parent, const Group &group, const Scalar *rates, const Phase injection_phase, const Scalar efficiency_factor, const std::vector< Scalar > &resv_coeff, DeferredLogger &deferred_logger) const
Definition: WellState.hpp:66
Dune::Communication< MPIComm > Communication
Definition: ParallelCommunication.hpp:30
Definition: blackoilbioeffectsmodules.hh:43