GroupStateHelper.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_GROUPSTATEHELPER_HEADER_INCLUDED
20#define OPM_GROUPSTATEHELPER_HEADER_INCLUDED
21
23
24#include <opm/common/TimingMacros.hpp>
25#include <opm/input/eclipse/EclipseState/Grid/FieldPropsManager.hpp>
26#include <opm/input/eclipse/Schedule/Group/GPMaint.hpp>
27#include <opm/input/eclipse/Schedule/Group/GSatProd.hpp>
28#include <opm/input/eclipse/Schedule/Group/GuideRate.hpp>
29#include <opm/input/eclipse/Schedule/Schedule.hpp>
30#include <opm/input/eclipse/Schedule/ScheduleState.hpp>
31#include <opm/input/eclipse/Schedule/SummaryState.hpp>
32#include <opm/material/fluidsystems/PhaseUsageInfo.hpp>
37
39
40#include <algorithm>
41#include <map>
42#include <memory>
43#include <optional>
44#include <stdexcept>
45#include <string>
46#include <vector>
47
48namespace Opm
49{
50
51template <typename Scalar, typename IndexTraits>
53{
54public:
55 // RAII guard for temporarily setting wellstate pointer
57 {
58 public:
60 : groupStateHelper_ {groupStateHelper}
61 , previous_state_ptr_ {groupStateHelper_.well_state_}
62 {
63 // Set the new state directly
64 groupStateHelper_.well_state_ = &well_state;
65 }
66
68 {
69 // Restore the previous state
70 groupStateHelper_.well_state_ = previous_state_ptr_;
71 }
72
73 // Delete copy and move operations
78
79 private:
80 GroupStateHelper& groupStateHelper_;
81 const WellState<Scalar, IndexTraits>* previous_state_ptr_;
82 };
83
84 // RAII guard for temporarily setting groupstate pointer
86 {
87 public:
88 GroupStateGuard(GroupStateHelper& group_state_helper, GroupState<Scalar>& group_state)
89 : group_state_helper_ {group_state_helper}
90 , previous_state_ptr_ {group_state_helper.group_state_}
91 {
92 // Set the new state directly
93 group_state_helper_.group_state_ = &group_state;
94 }
95
97 {
98 // Restore the previous state
99 group_state_helper_.group_state_ = previous_state_ptr_;
100 }
101
102 // Delete copy and move operations
107
108 private:
109 GroupStateHelper& group_state_helper_;
110 GroupState<Scalar>* previous_state_ptr_ {nullptr};
111 };
112
114 GroupState<Scalar>& group_state,
115 const Schedule& schedule,
116 const SummaryState& summary_state,
117 const GuideRate& guide_rate,
118 const PhaseUsageInfo<IndexTraits>& phase_usage_info);
119
120 void accumulateGroupEfficiencyFactor(const Group& group, Scalar& factor) const;
121
122 std::pair<bool, Scalar> checkGroupConstraintsInj(const std::string& name,
123 const std::string& parent,
124 const Group& group,
125 const Scalar* rates,
126 const Phase injection_phase,
127 const Scalar efficiency_factor,
128 const std::vector<Scalar>& resv_coeff,
129 const bool check_guide_rate,
130 DeferredLogger& deferred_logger) const;
131
132 std::pair<bool, Scalar> checkGroupConstraintsProd(const std::string& name,
133 const std::string& parent,
134 const Group& group,
135 const Scalar* rates,
136 const Scalar efficiency_factor,
137 const std::vector<Scalar>& resv_coeff,
138 const bool check_guide_rate,
139 DeferredLogger& deferred_logger) const;
140
141 Scalar getGuideRate(const std::string& name, const GuideRateModel::Target target) const;
142
143 Scalar getInjectionGroupTarget(const Group& group,
144 const Phase& injection_phase,
145 const std::vector<Scalar>& resv_coeff,
146 DeferredLogger& deferred_logger) const;
147
148 Scalar getProductionGroupTarget(const Group& group, DeferredLogger& deferred_logger) const;
149
150 GuideRate::RateVector getProductionGroupRateVector(const std::string& group_name) const;
151
152 std::optional<Scalar> getWellGroupTargetInjector(const std::string& name,
153 const std::string& parent,
154 const Group& group,
155 const Scalar* rates,
156 const Phase injection_phase,
157 const Scalar efficiency_factor,
158 const std::vector<Scalar>& resv_coeff,
159 DeferredLogger& deferred_logger) const;
160
161 std::optional<Scalar> getWellGroupTargetProducer(const std::string& name,
162 const std::string& parent,
163 const Group& group,
164 const Scalar* rates,
165 const Scalar efficiency_factor,
166 const std::vector<Scalar>& resv_coeff,
167 DeferredLogger& deferred_logger) const;
168
169 GuideRate::RateVector getWellRateVector(const std::string& name) const;
170
171 std::vector<std::string> groupChainTopBot(const std::string& bottom, const std::string& top) const;
172
175 int groupControlledWells(const std::string& group_name,
176 const std::string& always_included_child,
177 const bool is_production_group,
178 const Phase injection_phase) const;
179
181 {
182 return *this->group_state_;
183 }
184
186 {
187 return this->phase_usage_info_;
188 }
189
191 {
192 return GroupStateGuard(*this, group_state);
193 }
194
196 {
197 return WellStateGuard(*this, well_state);
198 }
199
200 int reportStepIdx() const
201 {
202 return report_step_;
203 }
204
205 const Schedule& schedule() const
206 {
207 return this->schedule_;
208 }
209
210 const GuideRate& guideRate() const
211 {
212 return this->guide_rate_;
213 }
214
215 int phaseToActivePhaseIdx(const Phase phase) const;
216
218 return this->phase_usage_info_;
219 }
220
221 // === Reservoir Coupling ===
222
225 const ReservoirCoupling::Proxy<Scalar>& rescoup() const { return rescoup_; }
226
227 bool isReservoirCouplingMaster() const { return rescoup_.isMaster(); }
228 bool isReservoirCouplingSlave() const { return rescoup_.isSlave(); }
229
232
233#ifdef RESERVOIR_COUPLING_ENABLED
234 void setReservoirCouplingMaster(ReservoirCouplingMaster<Scalar>* master) {
235 rescoup_.setMaster(master);
236 }
237 void setReservoirCouplingSlave(ReservoirCouplingSlave<Scalar>* slave) {
238 rescoup_.setSlave(slave);
239 }
240#endif
241
242 void setCmodeGroup(const Group& group);
243
244 template <class AverageRegionalPressureType>
245 void
247 const FieldPropsManager& fp,
248 std::map<std::string, std::unique_ptr<AverageRegionalPressureType>>&
249 regional_average_pressure_calculator) const;
250
251 void setReportStep(int report_step)
252 {
253 report_step_ = report_step;
254 }
255
256
257 const SummaryState& summaryState() const
258 {
259 return this->summary_state_;
260 }
261
262 Scalar sumSolventRates(const Group& group, const bool is_injector) const;
263
264 Scalar sumWellResRates(const Group& group, const int phase_pos, const bool injector) const;
265
266 Scalar sumWellSurfaceRates(const Group& group, const int phase_pos, const bool injector) const;
267
268 Scalar sumWellPhaseRates(bool res_rates,
269 const Group& group,
270 const int phase_pos,
271 const bool injector,
272 const bool network = false) const;
273
274 template <class RegionalValues>
275 void updateGpMaintTargetForGroups(const Group& group,
276 const RegionalValues& regional_values,
277 const double dt);
278
281 int updateGroupControlledWells(const bool is_production_group,
282 const Phase injection_phase,
283 DeferredLogger& deferred_logger);
284
285 void updateGroupProductionRates(const Group& group);
286
287 void updateGroupTargetReduction(const Group& group,
288 const bool is_injector);
289
291
292 void updateREINForGroups(const Group& group, bool sum_rank);
293
294 void updateReservoirRatesInjectionGroups(const Group& group);
295
297
298 void updateSurfaceRatesInjectionGroups(const Group& group);
299
300 void updateVREPForGroups(const Group& group);
301
302 void updateWellRates(const Group& group,
303 const WellState<Scalar, IndexTraits>& well_state_nupcol,
304 WellState<Scalar, IndexTraits>& well_state) const;
305
307 {
308 return *this->well_state_;
309 }
310
311 void updateWellRatesFromGroupTargetScale(const Scalar scale,
312 const Group& group,
313 bool is_injector,
314 WellState<Scalar, IndexTraits>& well_state) const;
315
317 std::pair<std::optional<std::string>, Scalar>
318 worstOffendingWell(const Group& group,
319 const Group::ProductionCMode& offended_control,
320 const Parallel::Communication& comm,
321 DeferredLogger& deferred_logger) const;
322
323private:
334 Scalar computeAddbackEfficiency_(const std::vector<std::string>& chain,
335 std::size_t local_reduction_level) const;
336
337 std::string controlGroup_(const Group& group) const;
338
339 GuideRate::RateVector getGuideRateVector_(const std::vector<Scalar>& rates) const;
340
342 bool isInGroupChainTopBot_(const std::string& bottom, const std::string& top) const;
343
344 Scalar satelliteInjectionRate_(const ScheduleState& sched,
345 const Group& group,
346 const int phase_pos,
347 bool res_rates) const;
348
349 Scalar satelliteProductionRate_(const ScheduleState& sched,
350 const Group& group,
351 const GSatProd::GSatProdGroupProp::Rate rate_comp,
352 bool res_rates) const;
353
354 std::optional<GSatProd::GSatProdGroupProp::Rate> selectRateComponent_(const int phase_pos) const;
355
356 int updateGroupControlledWellsRecursive_(const std::string& group_name,
357 const bool is_production_group,
358 const Phase injection_phase,
359 DeferredLogger& deferred_logger);
360
361 void updateGroupTargetReductionRecursive_(const Group& group,
362 const bool is_injector,
363 std::vector<Scalar>& group_target_reduction);
364
365 const WellState<Scalar, IndexTraits>* well_state_ {nullptr};
366 GroupState<Scalar>* group_state_ {nullptr};
367 const Schedule& schedule_;
368 const SummaryState& summary_state_;
369 const GuideRate& guide_rate_;
370 int report_step_ {0};
371 // NOTE: The phase usage info seems to be read-only throughout the simulation, so it should be safe
372 // to store a reference to it here.
373 const PhaseUsageInfo<IndexTraits>& phase_usage_info_;
374 ReservoirCoupling::Proxy<Scalar> rescoup_{};
375};
376
377// -----------------------------------------------------------------------------
378// Template member function implementations
379// -----------------------------------------------------------------------------
380
381// NOTE: This template member function is defined here in the header because the
382// AverageRegionalPressureType type parameter depends on derived class types that are
383// not available when GroupStateHelper.cpp is compiled. Template functions
384// must be visible at their instantiation point.
385// See GroupStateHelper.cpp for detailed rationale.
386template <typename Scalar, typename IndexTraits>
387template <class AverageRegionalPressureType>
388void
390 const Group& group,
391 const FieldPropsManager& fp,
392 std::map<std::string, std::unique_ptr<AverageRegionalPressureType>>& regional_average_pressure_calculator)
393 const
394{
395 for (const std::string& groupName : group.groups()) {
396 this->setRegionAveragePressureCalculator(this->schedule_.getGroup(groupName, this->report_step_),
397 fp,
398 regional_average_pressure_calculator);
399 }
400 const auto& gpm = group.gpmaint();
401 if (!gpm)
402 return;
403
404 const auto& reg = gpm->region();
405 if (!reg)
406 return;
407
408 if (regional_average_pressure_calculator.count(reg->first) == 0) {
409 const std::string name = (reg->first.rfind("FIP", 0) == 0) ? reg->first : "FIP" + reg->first;
410 const auto& fipnum = fp.get_int(name);
411 regional_average_pressure_calculator[reg->first]
412 = std::make_unique<AverageRegionalPressureType>(fipnum);
413 }
414}
415
416// NOTE: This template member function is defined in the header because the
417// RegionalValues type parameter depends on derived class types that are
418// not available when GroupStateHelper.cpp is compiled. Template functions
419// must be visible at their instantiation point.
420// See GroupStateHelper.cpp for detailed rationale.
421template <typename Scalar, typename IndexTraits>
422template <class RegionalValues>
423void
425 const RegionalValues& regional_values,
426 const double dt)
427{
428 OPM_TIMEFUNCTION();
429 for (const std::string& group_name : group.groups()) {
430 const Group& group_tmp = this->schedule_.getGroup(group_name, this->report_step_);
431 this->updateGpMaintTargetForGroups(group_tmp, regional_values, dt);
432 }
433 const auto& gpm = group.gpmaint();
434 if (!gpm)
435 return;
436
437 const auto& region = gpm->region();
438 if (!region)
439 return;
440
441 const auto [name, number] = *region;
442 const Scalar error = gpm->pressure_target() - regional_values.at(name)->pressure(number);
443 Scalar current_rate = 0.0;
444 const auto& pu = this->phase_usage_info_;
445 bool injection = true;
446 Scalar sign = 1.0;
447 switch (gpm->flow_target()) {
448 case GPMaint::FlowTarget::RESV_PROD: {
449 current_rate = -this->groupState().injection_vrep_rate(group.name());
450 injection = false;
451 sign = -1.0;
452 break;
453 }
454 case GPMaint::FlowTarget::RESV_OINJ: {
455 if (pu.phaseIsActive(IndexTraits::oilPhaseIdx)) {
456 const auto io = pu.canonicalToActivePhaseIdx(IndexTraits::oilPhaseIdx);
457 current_rate = this->groupState().injection_reservoir_rates(group.name())[io];
458 }
459 break;
460 }
461 case GPMaint::FlowTarget::RESV_WINJ: {
462 if (pu.phaseIsActive(IndexTraits::waterPhaseIdx)) {
463 const auto iw = pu.canonicalToActivePhaseIdx(IndexTraits::waterPhaseIdx);
464 current_rate = this->groupState().injection_reservoir_rates(group.name())[iw];
465 }
466 break;
467 }
468 case GPMaint::FlowTarget::RESV_GINJ: {
469 if (pu.phaseIsActive(IndexTraits::gasPhaseIdx)) {
470 const auto ig = pu.canonicalToActivePhaseIdx(IndexTraits::gasPhaseIdx);
471 current_rate = this->groupState().injection_reservoir_rates(group.name())[ig];
472 }
473 break;
474 }
475 case GPMaint::FlowTarget::SURF_OINJ: {
476 if (pu.phaseIsActive(IndexTraits::oilPhaseIdx)) {
477 const auto io = pu.canonicalToActivePhaseIdx(IndexTraits::oilPhaseIdx);
478 current_rate = this->groupState().injection_surface_rates(group.name())[io];
479 }
480 break;
481 }
482 case GPMaint::FlowTarget::SURF_WINJ: {
483 if (pu.phaseIsActive(IndexTraits::waterPhaseIdx)) {
484 const auto iw = pu.canonicalToActivePhaseIdx(IndexTraits::waterPhaseIdx);
485 current_rate = this->groupState().injection_surface_rates(group.name())[iw];
486 }
487 break;
488 }
489 case GPMaint::FlowTarget::SURF_GINJ: {
490 if (pu.phaseIsActive(IndexTraits::gasPhaseIdx)) {
491 const auto ig = pu.canonicalToActivePhaseIdx(IndexTraits::gasPhaseIdx);
492 current_rate = this->groupState().injection_surface_rates(group.name())[ig];
493 }
494 break;
495 }
496 default:
497 throw std::invalid_argument("Invalid Flow target type in GPMAINT");
498 }
499 auto& gpmaint_state = this->groupState().gpmaint(group.name());
500 // we only activate gpmaint if pressure is lower than the target regional pressure for injectors
501 // (i.e. error > 0) and higher for producers.
502 bool activate = (injection && error > 0) || (!injection && error < 0);
503 Scalar rate = 0.0;
504 if (activate) {
505 rate = gpm->rate(gpmaint_state, current_rate, error, dt);
506 } else {
507 gpm->resetState(gpmaint_state);
508 }
509 this->groupState().update_gpmaint_target(group.name(), std::max(Scalar {0.0}, sign * rate));
510}
511
512} // namespace Opm
513
514#endif // OPM_GROUPSTATE_HELPER_HEADER_INCLUDED
Definition: DeferredLogger.hpp:57
Definition: GroupStateHelper.hpp:86
GroupStateGuard(GroupStateGuard &&)=delete
GroupStateGuard & operator=(GroupStateGuard &&)=delete
GroupStateGuard(GroupStateHelper &group_state_helper, GroupState< Scalar > &group_state)
Definition: GroupStateHelper.hpp:88
GroupStateGuard(const GroupStateGuard &)=delete
~GroupStateGuard()
Definition: GroupStateHelper.hpp:96
GroupStateGuard & operator=(const GroupStateGuard &)=delete
Definition: GroupStateHelper.hpp:57
WellStateGuard(WellStateGuard &&)=delete
WellStateGuard(const WellStateGuard &)=delete
~WellStateGuard()
Definition: GroupStateHelper.hpp:67
WellStateGuard & operator=(const WellStateGuard &)=delete
WellStateGuard & operator=(WellStateGuard &&)=delete
WellStateGuard(GroupStateHelper &groupStateHelper, WellState< Scalar, IndexTraits > &well_state)
Definition: GroupStateHelper.hpp:59
Definition: GroupStateHelper.hpp:53
void updateWellRatesFromGroupTargetScale(const Scalar scale, const Group &group, bool is_injector, WellState< Scalar, IndexTraits > &well_state) const
void setRegionAveragePressureCalculator(const Group &group, const FieldPropsManager &fp, std::map< std::string, std::unique_ptr< AverageRegionalPressureType > > &regional_average_pressure_calculator) const
Definition: GroupStateHelper.hpp:389
void setCmodeGroup(const Group &group)
Scalar getInjectionGroupTarget(const Group &group, const Phase &injection_phase, const std::vector< Scalar > &resv_coeff, DeferredLogger &deferred_logger) const
std::vector< std::string > groupChainTopBot(const std::string &bottom, const std::string &top) const
void updateVREPForGroups(const Group &group)
void updateWellRates(const Group &group, const WellState< Scalar, IndexTraits > &well_state_nupcol, WellState< Scalar, IndexTraits > &well_state) const
ReservoirCouplingMaster< Scalar > & reservoirCouplingMaster()
Definition: GroupStateHelper.hpp:230
GroupState< Scalar > & groupState() const
Definition: GroupStateHelper.hpp:180
Scalar sumWellResRates(const Group &group, const int phase_pos, const bool injector) const
void accumulateGroupEfficiencyFactor(const Group &group, Scalar &factor) const
int groupControlledWells(const std::string &group_name, const std::string &always_included_child, const bool is_production_group, const Phase injection_phase) const
const PhaseUsageInfo< IndexTraits > & phaseUsage() const
Definition: GroupStateHelper.hpp:217
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
const SummaryState & summaryState() const
Definition: GroupStateHelper.hpp:257
const GuideRate & guideRate() const
Definition: GroupStateHelper.hpp:210
ReservoirCoupling::Proxy< Scalar > & rescoup()
Get the reservoir coupling proxy.
Definition: GroupStateHelper.hpp:224
void updateReservoirRatesInjectionGroups(const Group &group)
const WellState< Scalar, IndexTraits > & wellState() const
Definition: GroupStateHelper.hpp:306
const PhaseUsageInfo< IndexTraits > & phaseUsageInfo() const
Definition: GroupStateHelper.hpp:185
GuideRate::RateVector getWellRateVector(const std::string &name) const
int updateGroupControlledWells(const bool is_production_group, const Phase injection_phase, DeferredLogger &deferred_logger)
Scalar sumSolventRates(const Group &group, const bool is_injector) const
const Schedule & schedule() const
Definition: GroupStateHelper.hpp:205
void updateGroupTargetReduction(const Group &group, const bool is_injector)
std::optional< 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
GuideRate::RateVector getProductionGroupRateVector(const std::string &group_name) const
void updateGroupProductionRates(const Group &group)
GroupStateHelper(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)
void updateGpMaintTargetForGroups(const Group &group, const RegionalValues &regional_values, const double dt)
Definition: GroupStateHelper.hpp:424
int phaseToActivePhaseIdx(const Phase phase) const
Scalar sumWellSurfaceRates(const Group &group, const int phase_pos, const bool injector) const
int reportStepIdx() const
Definition: GroupStateHelper.hpp:200
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)
Scalar getGuideRate(const std::string &name, const GuideRateModel::Target target) const
void setReportStep(int report_step)
Definition: GroupStateHelper.hpp:251
void updateState(WellState< Scalar, IndexTraits > &well_state, GroupState< Scalar > &group_state)
const ReservoirCoupling::Proxy< Scalar > & rescoup() const
Definition: GroupStateHelper.hpp:225
bool isReservoirCouplingMaster() const
Definition: GroupStateHelper.hpp:227
void updateNetworkLeafNodeProductionRates()
WellStateGuard pushWellState(WellState< Scalar, IndexTraits > &well_state)
Definition: GroupStateHelper.hpp:195
GroupStateGuard pushGroupState(GroupState< Scalar > &group_state)
Definition: GroupStateHelper.hpp:190
ReservoirCouplingSlave< Scalar > & reservoirCouplingSlave()
Definition: GroupStateHelper.hpp:231
std::optional< 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
Scalar sumWellPhaseRates(bool res_rates, const Group &group, const int phase_pos, const bool injector, const bool network=false) const
void updateREINForGroups(const Group &group, bool sum_rank)
bool isReservoirCouplingSlave() const
Definition: GroupStateHelper.hpp:228
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
Scalar getProductionGroupTarget(const Group &group, DeferredLogger &deferred_logger) const
void updateSurfaceRatesInjectionGroups(const Group &group)
Definition: GroupState.hpp:41
Definition: GasLiftGroupInfo.hpp:37
Thin proxy for reservoir coupling master/slave pointers.
Definition: RescoupProxy.hpp:53
Definition: ReservoirCouplingMaster.hpp:38
Definition: ReservoirCouplingSlave.hpp:40
Definition: WellState.hpp:66
Dune::Communication< MPIComm > Communication
Definition: ParallelCommunication.hpp:30
Phase
Phase indices for reservoir coupling, we currently only support black-oil phases (oil,...
Definition: ReservoirCoupling.hpp:141
Definition: blackoilbioeffectsmodules.hh:43