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,
120
121 void accumulateGroupEfficiencyFactor(const Group& group, Scalar& factor) const;
122
123 std::pair<bool, Scalar> checkGroupConstraintsInj(const std::string& name,
124 const std::string& parent,
125 const Group& group,
126 const Scalar* rates,
127 const Phase injection_phase,
128 const Scalar efficiency_factor,
129 const std::vector<Scalar>& resv_coeff,
130 const bool check_guide_rate,
131 DeferredLogger& deferred_logger) const;
132
133 std::pair<bool, Scalar> checkGroupConstraintsProd(const std::string& name,
134 const std::string& parent,
135 const Group& group,
136 const Scalar* rates,
137 const Scalar efficiency_factor,
138 const std::vector<Scalar>& resv_coeff,
139 const bool check_guide_rate,
140 DeferredLogger& deferred_logger) const;
141
142 const Parallel::Communication& comm() const { return this->comm_; }
143
144 std::vector<Scalar> getGroupRatesAvailableForHigherLevelControl(const Group& group, const bool is_injector) const;
145
146 Scalar getGuideRate(const std::string& name, const GuideRateModel::Target target) const;
147
148 Scalar getInjectionGroupTarget(const Group& group,
149 const Phase& injection_phase,
150 const std::vector<Scalar>& resv_coeff,
151 DeferredLogger& deferred_logger) const;
152
153 Scalar getProductionGroupTarget(const Group& group, DeferredLogger& deferred_logger) const;
154
155 GuideRate::RateVector getProductionGroupRateVector(const std::string& group_name) const;
156
157 std::optional<Scalar> getWellGroupTargetInjector(const std::string& name,
158 const std::string& parent,
159 const Group& group,
160 const Scalar* rates,
161 const Phase injection_phase,
162 const Scalar efficiency_factor,
163 const std::vector<Scalar>& resv_coeff,
164 DeferredLogger& deferred_logger) const;
165
166 std::optional<Scalar> getWellGroupTargetProducer(const std::string& name,
167 const std::string& parent,
168 const Group& group,
169 const Scalar* rates,
170 const Scalar efficiency_factor,
171 const std::vector<Scalar>& resv_coeff,
172 DeferredLogger& deferred_logger) const;
173
174 GuideRate::RateVector getWellRateVector(const std::string& name) const;
175
176 std::vector<std::string> groupChainTopBot(const std::string& bottom, const std::string& top) const;
177
180 int groupControlledWells(const std::string& group_name,
181 const std::string& always_included_child,
182 const bool is_production_group,
183 const Phase injection_phase) const;
184
186 {
187 return *this->group_state_;
188 }
189
190 const GuideRate& guideRate() const
191 {
192 return this->guide_rate_;
193 }
194
195 bool isRank0() const
196 {
197 return this->well_state_->isRank0();
198 }
199
200 bool isReservoirCouplingMaster() const { return rescoup_.isMaster(); }
201
202 bool isReservoirCouplingMasterGroup(const Group& group) const { return rescoup_.isMasterGroup(group.name()); }
203
204 bool isReservoirCouplingSlave() const { return rescoup_.isSlave(); }
205
206 constexpr int numPhases() const {
207 return this->wellState().numPhases();
208 }
209
210 int phaseToActivePhaseIdx(const Phase phase) const;
211
213 return this->phase_usage_info_;
214 }
215
217 {
218 return GroupStateGuard(*this, group_state);
219 }
220
222 {
223 return WellStateGuard(*this, well_state);
224 }
225
226 int reportStepIdx() const
227 {
228 return report_step_;
229 }
230
231 const Schedule& schedule() const
232 {
233 return this->schedule_;
234 }
235
237 return rescoup_;
238 }
239
241 return rescoup_;
242 }
243
245 return rescoup_.master();
246 }
247
249 return rescoup_.master();
250 }
251
253 return rescoup_.slave();
254 }
256 return rescoup_.slave();
257 }
258
259#ifdef RESERVOIR_COUPLING_ENABLED
260 void setReservoirCouplingMaster(ReservoirCouplingMaster<Scalar>* master) {
261 rescoup_.setMaster(master);
262 }
263 void setReservoirCouplingSlave(ReservoirCouplingSlave<Scalar>* slave) {
264 rescoup_.setSlave(slave);
265 }
266#endif
267
268 void setCmodeGroup(const Group& group);
269
270 template <class AverageRegionalPressureType>
271 void
273 const FieldPropsManager& fp,
274 std::map<std::string, std::unique_ptr<AverageRegionalPressureType>>&
275 regional_average_pressure_calculator) const;
276
277 void setReportStep(int report_step)
278 {
279 report_step_ = report_step;
280 }
281
282
283 const SummaryState& summaryState() const
284 {
285 return this->summary_state_;
286 }
287
288 Scalar sumSolventRates(const Group& group, const bool is_injector) const;
289
290 Scalar sumWellResRates(const Group& group, const int phase_pos, const bool injector) const;
291
292 Scalar sumWellSurfaceRates(const Group& group, const int phase_pos, const bool injector) const;
293
294 Scalar sumWellPhaseRates(bool res_rates,
295 const Group& group,
296 const int phase_pos,
297 const bool injector,
298 const bool network = false) const;
299
300 template <class RegionalValues>
301 void updateGpMaintTargetForGroups(const Group& group,
302 const RegionalValues& regional_values,
303 const double dt);
304
307 int updateGroupControlledWells(const bool is_production_group,
308 const Phase injection_phase,
309 DeferredLogger& deferred_logger);
310
311 void updateGroupProductionRates(const Group& group);
312
313 void updateGroupTargetReduction(const Group& group,
314 const bool is_injector);
315
317
318 void updateREINForGroups(const Group& group, bool sum_rank);
319
320 void updateReservoirRatesInjectionGroups(const Group& group);
321
323
324 void updateSurfaceRatesInjectionGroups(const Group& group);
325
326 void updateVREPForGroups(const Group& group);
327
328 void updateWellRates(const Group& group,
329 const WellState<Scalar, IndexTraits>& well_state_nupcol,
330 WellState<Scalar, IndexTraits>& well_state) const;
331
333 {
334 return *this->well_state_;
335 }
336
337 void updateWellRatesFromGroupTargetScale(const Scalar scale,
338 const Group& group,
339 bool is_injector,
340 WellState<Scalar, IndexTraits>& well_state) const;
341
343 std::pair<std::optional<std::string>, Scalar>
344 worstOffendingWell(const Group& group,
345 const Group::ProductionCMode& offended_control,
346 DeferredLogger& deferred_logger) const;
347
348private:
349#ifdef RESERVOIR_COUPLING_ENABLED
354 ReservoirCoupling::Phase activePhaseIdxToRescoupPhase_(int phase_pos) const;
355#endif
356
367 Scalar computeAddbackEfficiency_(const std::vector<std::string>& chain,
368 std::size_t local_reduction_level) const;
369
370 std::string controlGroup_(const Group& group) const;
371
372 GuideRate::RateVector getGuideRateVector_(const std::vector<Scalar>& rates) const;
373
374 Scalar getReservoirCouplingMasterGroupRate_(const Group& group,
375 const int phase_pos,
376 const bool res_rates,
377 const bool is_injector) const;
378
379 Scalar getSatelliteRate_(const Group& group,
380 const int phase_pos,
381 const bool res_rates,
382 const bool is_injector) const;
383
385 bool isInGroupChainTopBot_(const std::string& bottom, const std::string& top) const;
386
387 bool isSatelliteGroup_(const Group& group) const;
388
389 Scalar satelliteInjectionRate_(const ScheduleState& sched,
390 const Group& group,
391 const int phase_pos,
392 bool res_rates) const;
393
394 Scalar satelliteProductionRate_(const ScheduleState& sched,
395 const Group& group,
396 const GSatProd::GSatProdGroupProp::Rate rate_comp,
397 bool res_rates) const;
398
399 std::optional<GSatProd::GSatProdGroupProp::Rate> selectRateComponent_(const int phase_pos) const;
400
401 int updateGroupControlledWellsRecursive_(const std::string& group_name,
402 const bool is_production_group,
403 const Phase injection_phase,
404 DeferredLogger& deferred_logger);
405
406 void updateGroupTargetReductionRecursive_(const Group& group,
407 const bool is_injector,
408 std::vector<Scalar>& group_target_reduction);
409
410 const WellState<Scalar, IndexTraits>* well_state_ {nullptr};
411 GroupState<Scalar>* group_state_ {nullptr};
412 const Schedule& schedule_;
413 const SummaryState& summary_state_;
414 const GuideRate& guide_rate_;
415 // NOTE: The phase usage info seems to be read-only throughout the simulation, so it should be safe
416 // to store a reference to it here.
417 const PhaseUsageInfo<IndexTraits>& phase_usage_info_;
418 const Parallel::Communication& comm_;
419 int report_step_ {0};
420 ReservoirCoupling::Proxy<Scalar> rescoup_{};
421};
422
423// -----------------------------------------------------------------------------
424// Template member function implementations
425// -----------------------------------------------------------------------------
426
427// NOTE: This template member function is defined here in the header because the
428// AverageRegionalPressureType type parameter depends on derived class types that are
429// not available when GroupStateHelper.cpp is compiled. Template functions
430// must be visible at their instantiation point.
431// See GroupStateHelper.cpp for detailed rationale.
432template <typename Scalar, typename IndexTraits>
433template <class AverageRegionalPressureType>
434void
436 const Group& group,
437 const FieldPropsManager& fp,
438 std::map<std::string, std::unique_ptr<AverageRegionalPressureType>>& regional_average_pressure_calculator)
439 const
440{
441 for (const std::string& groupName : group.groups()) {
442 this->setRegionAveragePressureCalculator(this->schedule_.getGroup(groupName, this->report_step_),
443 fp,
444 regional_average_pressure_calculator);
445 }
446 const auto& gpm = group.gpmaint();
447 if (!gpm)
448 return;
449
450 const auto& reg = gpm->region();
451 if (!reg)
452 return;
453
454 if (regional_average_pressure_calculator.count(reg->first) == 0) {
455 const std::string name = (reg->first.rfind("FIP", 0) == 0) ? reg->first : "FIP" + reg->first;
456 const auto& fipnum = fp.get_int(name);
457 regional_average_pressure_calculator[reg->first]
458 = std::make_unique<AverageRegionalPressureType>(fipnum);
459 }
460}
461
462// NOTE: This template member function is defined in the header because the
463// RegionalValues type parameter depends on derived class types that are
464// not available when GroupStateHelper.cpp is compiled. Template functions
465// must be visible at their instantiation point.
466// See GroupStateHelper.cpp for detailed rationale.
467template <typename Scalar, typename IndexTraits>
468template <class RegionalValues>
469void
471 const RegionalValues& regional_values,
472 const double dt)
473{
474 OPM_TIMEFUNCTION();
475 for (const std::string& group_name : group.groups()) {
476 const Group& group_tmp = this->schedule_.getGroup(group_name, this->report_step_);
477 this->updateGpMaintTargetForGroups(group_tmp, regional_values, dt);
478 }
479 const auto& gpm = group.gpmaint();
480 if (!gpm)
481 return;
482
483 const auto& region = gpm->region();
484 if (!region)
485 return;
486
487 const auto [name, number] = *region;
488 const Scalar error = gpm->pressure_target() - regional_values.at(name)->pressure(number);
489 Scalar current_rate = 0.0;
490 const auto& pu = this->phase_usage_info_;
491 bool injection = true;
492 Scalar sign = 1.0;
493 switch (gpm->flow_target()) {
494 case GPMaint::FlowTarget::RESV_PROD: {
495 current_rate = -this->groupState().injection_vrep_rate(group.name());
496 injection = false;
497 sign = -1.0;
498 break;
499 }
500 case GPMaint::FlowTarget::RESV_OINJ: {
501 if (pu.phaseIsActive(IndexTraits::oilPhaseIdx)) {
502 const auto io = pu.canonicalToActivePhaseIdx(IndexTraits::oilPhaseIdx);
503 current_rate = this->groupState().injection_reservoir_rates(group.name())[io];
504 }
505 break;
506 }
507 case GPMaint::FlowTarget::RESV_WINJ: {
508 if (pu.phaseIsActive(IndexTraits::waterPhaseIdx)) {
509 const auto iw = pu.canonicalToActivePhaseIdx(IndexTraits::waterPhaseIdx);
510 current_rate = this->groupState().injection_reservoir_rates(group.name())[iw];
511 }
512 break;
513 }
514 case GPMaint::FlowTarget::RESV_GINJ: {
515 if (pu.phaseIsActive(IndexTraits::gasPhaseIdx)) {
516 const auto ig = pu.canonicalToActivePhaseIdx(IndexTraits::gasPhaseIdx);
517 current_rate = this->groupState().injection_reservoir_rates(group.name())[ig];
518 }
519 break;
520 }
521 case GPMaint::FlowTarget::SURF_OINJ: {
522 if (pu.phaseIsActive(IndexTraits::oilPhaseIdx)) {
523 const auto io = pu.canonicalToActivePhaseIdx(IndexTraits::oilPhaseIdx);
524 current_rate = this->groupState().injection_surface_rates(group.name())[io];
525 }
526 break;
527 }
528 case GPMaint::FlowTarget::SURF_WINJ: {
529 if (pu.phaseIsActive(IndexTraits::waterPhaseIdx)) {
530 const auto iw = pu.canonicalToActivePhaseIdx(IndexTraits::waterPhaseIdx);
531 current_rate = this->groupState().injection_surface_rates(group.name())[iw];
532 }
533 break;
534 }
535 case GPMaint::FlowTarget::SURF_GINJ: {
536 if (pu.phaseIsActive(IndexTraits::gasPhaseIdx)) {
537 const auto ig = pu.canonicalToActivePhaseIdx(IndexTraits::gasPhaseIdx);
538 current_rate = this->groupState().injection_surface_rates(group.name())[ig];
539 }
540 break;
541 }
542 default:
543 throw std::invalid_argument("Invalid Flow target type in GPMAINT");
544 }
545 auto& gpmaint_state = this->groupState().gpmaint(group.name());
546 // we only activate gpmaint if pressure is lower than the target regional pressure for injectors
547 // (i.e. error > 0) and higher for producers.
548 bool activate = (injection && error > 0) || (!injection && error < 0);
549 Scalar rate = 0.0;
550 if (activate) {
551 rate = gpm->rate(gpmaint_state, current_rate, error, dt);
552 } else {
553 gpm->resetState(gpmaint_state);
554 }
555 this->groupState().update_gpmaint_target(group.name(), std::max(Scalar {0.0}, sign * rate));
556}
557
558} // namespace Opm
559
560#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:435
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:244
GroupState< Scalar > & groupState() const
Definition: GroupStateHelper.hpp:185
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:212
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
std::pair< std::optional< std::string >, Scalar > worstOffendingWell(const Group &group, const Group::ProductionCMode &offended_control, DeferredLogger &deferred_logger) const
Returns the name of the worst offending well and its fraction (i.e. violated_phase / preferred_phase)
const ReservoirCouplingMaster< Scalar > & reservoirCouplingMaster() const
Definition: GroupStateHelper.hpp:248
const SummaryState & summaryState() const
Definition: GroupStateHelper.hpp:283
const GuideRate & guideRate() const
Definition: GroupStateHelper.hpp:190
ReservoirCoupling::Proxy< Scalar > & rescoup()
Definition: GroupStateHelper.hpp:236
void updateReservoirRatesInjectionGroups(const Group &group)
const WellState< Scalar, IndexTraits > & wellState() const
Definition: GroupStateHelper.hpp:332
constexpr int numPhases() const
Definition: GroupStateHelper.hpp:206
GuideRate::RateVector getWellRateVector(const std::string &name) const
const Parallel::Communication & comm() const
Definition: GroupStateHelper.hpp:142
bool isRank0() const
Definition: GroupStateHelper.hpp:195
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:231
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)
void updateGpMaintTargetForGroups(const Group &group, const RegionalValues &regional_values, const double dt)
Definition: GroupStateHelper.hpp:470
int phaseToActivePhaseIdx(const Phase phase) const
Scalar sumWellSurfaceRates(const Group &group, const int phase_pos, const bool injector) const
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, const Parallel::Communication &comm)
int reportStepIdx() const
Definition: GroupStateHelper.hpp:226
Scalar getGuideRate(const std::string &name, const GuideRateModel::Target target) const
void setReportStep(int report_step)
Definition: GroupStateHelper.hpp:277
void updateState(WellState< Scalar, IndexTraits > &well_state, GroupState< Scalar > &group_state)
std::vector< Scalar > getGroupRatesAvailableForHigherLevelControl(const Group &group, const bool is_injector) const
bool isReservoirCouplingMasterGroup(const Group &group) const
Definition: GroupStateHelper.hpp:202
const ReservoirCoupling::Proxy< Scalar > & rescoup() const
Definition: GroupStateHelper.hpp:240
bool isReservoirCouplingMaster() const
Definition: GroupStateHelper.hpp:200
const ReservoirCouplingSlave< Scalar > & reservoirCouplingSlave() const
Definition: GroupStateHelper.hpp:255
void updateNetworkLeafNodeProductionRates()
WellStateGuard pushWellState(WellState< Scalar, IndexTraits > &well_state)
Definition: GroupStateHelper.hpp:221
GroupStateGuard pushGroupState(GroupState< Scalar > &group_state)
Definition: GroupStateHelper.hpp:216
ReservoirCouplingSlave< Scalar > & reservoirCouplingSlave()
Definition: GroupStateHelper.hpp:252
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:204
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:54
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:147
Definition: blackoilbioeffectsmodules.hh:43