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>
39
41
42#include <algorithm>
43#include <map>
44#include <memory>
45#include <optional>
46#include <stdexcept>
47#include <string>
48#include <vector>
49
50namespace Opm
51{
52
53template <typename Scalar, typename IndexTraits>
55{
56public:
57 // RAII guard for temporarily setting wellstate pointer
59 {
60 public:
62 : groupStateHelper_ {groupStateHelper}
63 , previous_state_ptr_ {groupStateHelper_.well_state_}
64 {
65 // Set the new state directly
66 groupStateHelper_.well_state_ = &well_state;
67 }
68
70 {
71 // Restore the previous state
72 groupStateHelper_.well_state_ = previous_state_ptr_;
73 }
74
75 // Delete copy and move operations
80
81 private:
82 GroupStateHelper& groupStateHelper_;
83 const WellState<Scalar, IndexTraits>* previous_state_ptr_;
84 };
85
86 // RAII guard for temporarily setting groupstate pointer
88 {
89 public:
90 GroupStateGuard(GroupStateHelper& group_state_helper, GroupState<Scalar>& group_state)
91 : group_state_helper_ {group_state_helper}
92 , previous_state_ptr_ {group_state_helper.group_state_}
93 {
94 // Set the new state directly
95 group_state_helper_.group_state_ = &group_state;
96 }
97
99 {
100 // Restore the previous state
101 group_state_helper_.group_state_ = previous_state_ptr_;
102 }
103
104 // Delete copy and move operations
109
110 private:
111 GroupStateHelper& group_state_helper_;
112 GroupState<Scalar>* previous_state_ptr_ {nullptr};
113 };
114
128 {
129 public:
135 explicit ScopedLoggerGuard(const GroupStateHelper& helper, bool do_mpi_gather = true)
136 : helper_(&helper)
137 , previous_(helper.deferred_logger_)
138 , do_mpi_gather_(do_mpi_gather)
139 {
140 helper_->deferred_logger_ = &logger_;
141 }
142
144 {
145 if (helper_) {
146 if (do_mpi_gather_) {
147 // 1. Gather messages across MPI ranks
148 DeferredLogger global = gatherDeferredLogger(logger_, helper_->comm());
149
150 // 2. Log on rank 0 (if terminal_output enabled)
151 if (helper_->terminalOutput()) {
152 global.logMessages();
153 }
154 } else {
155 // Just log locally without MPI gather
156 if (helper_->terminalOutput()) {
157 logger_.logMessages();
158 }
159 }
160
161 // 3. Restore previous logger
162 helper_->deferred_logger_ = previous_;
163 }
164 }
165
166 // Delete copy operations and move assignment
170
171 // Move constructor required for pushLogger() return-by-value (must be
172 // available even if elided by RVO)
174 : helper_(other.helper_)
175 , logger_(std::move(other.logger_))
176 , previous_(other.previous_)
177 , do_mpi_gather_(other.do_mpi_gather_)
178 {
179 // Update the helper's pointer to our moved logger
180 if (helper_) {
181 helper_->deferred_logger_ = &logger_;
182 }
183 other.helper_ = nullptr;
184 other.previous_ = nullptr;
185 }
186
187 private:
188 // Pointer (not reference) to allow nulling in move constructor
189 const GroupStateHelper* helper_{nullptr};
190 DeferredLogger logger_; // Owned logger
191 DeferredLogger* previous_{nullptr}; // For restore
192 bool do_mpi_gather_{true}; // Whether to gather messages across MPI ranks
193 };
194
196 GroupState<Scalar>& group_state,
197 const Schedule& schedule,
198 const SummaryState& summary_state,
199 const GuideRate& guide_rate,
200 const PhaseUsageInfo<IndexTraits>& phase_usage_info,
202 bool terminal_output);
203
204 void accumulateGroupEfficiencyFactor(const Group& group, Scalar& factor) const;
205
206 std::pair<bool, Scalar> checkGroupConstraintsInj(const std::string& name,
207 const std::string& parent,
208 const Group& group,
209 const Scalar* rates,
210 const Phase injection_phase,
211 const Scalar efficiency_factor,
212 const std::vector<Scalar>& resv_coeff,
213 const bool check_guide_rate) const;
214
215 std::pair<bool, Scalar> checkGroupConstraintsProd(const std::string& name,
216 const std::string& parent,
217 const Group& group,
218 const Scalar* rates,
219 const Scalar efficiency_factor,
220 const std::vector<Scalar>& resv_coeff,
221 const bool check_guide_rate) const;
222
223 const Parallel::Communication& comm() const { return this->comm_; }
224
228 {
229 if (this->deferred_logger_ == nullptr) {
230 throw std::logic_error("DeferredLogger not set. Call pushLogger() first.");
231 }
232 return *this->deferred_logger_;
233 }
234
235 std::vector<Scalar> getGroupRatesAvailableForHigherLevelControl(const Group& group, const bool is_injector) const;
236
237 Scalar getInjectionGroupTarget(const Group& group,
238 const Phase& injection_phase,
239 const std::vector<Scalar>& resv_coeff) const;
240
244 GuideRateModel::Target getInjectionGuideTargetMode(Phase injection_phase) const;
245
246 Scalar getProductionGroupTarget(const Group& group) const;
247
251 GuideRateModel::Target getProductionGuideTargetMode(const Group& group) const;
252
253 std::pair<Scalar, Group::ProductionCMode>
254 getAutoChokeGroupProductionTargetRate(const Group& bottom_group,
255 const Group& group,
256 const std::vector<Scalar>& resv_coeff,
257 Scalar efficiencyFactor) const;
258
259 GuideRate::RateVector getProductionGroupRateVector(const std::string& group_name) const;
260
262
263 std::optional<GroupTarget> getWellGroupTargetInjector(const std::string& name,
264 const std::string& parent,
265 const Group& group,
266 const Scalar* rates,
267 const Phase injection_phase,
268 const Scalar efficiency_factor,
269 const std::vector<Scalar>& resv_coeff) const;
270
271 std::optional<GroupTarget> getWellGroupTargetProducer(const std::string& name,
272 const std::string& parent,
273 const Group& group,
274 const Scalar* rates,
275 const Scalar efficiency_factor,
276 const std::vector<Scalar>& resv_coeff) const;
277
278 GuideRate::RateVector getWellRateVector(const std::string& name) const;
279
280 std::vector<std::string> groupChainTopBot(const std::string& bottom, const std::string& top) const;
281
284 int groupControlledWells(const std::string& group_name,
285 const std::string& always_included_child,
286 const bool is_production_group,
287 const Phase injection_phase) const;
288
290 {
291 return *this->group_state_;
292 }
293
294 const GuideRate& guideRate() const
295 {
296 return this->guide_rate_;
297 }
298
299 bool isRank0() const
300 {
301 return this->well_state_->isRank0();
302 }
303
304 bool isReservoirCouplingMaster() const { return rescoup_.isMaster(); }
305
306 bool isReservoirCouplingMasterGroup(const Group& group) const { return rescoup_.isMasterGroup(group.name()); }
307
308 bool isReservoirCouplingSlave() const { return rescoup_.isSlave(); }
309
310 bool isReservoirCouplingSlaveGroup(const Group& group) const { return rescoup_.isSlaveGroup(group.name()); }
311
312 constexpr int numPhases() const {
313 return this->wellState().numPhases();
314 }
315
316 int phaseToActivePhaseIdx(const Phase phase) const;
317
319 return this->phase_usage_info_;
320 }
321
323 {
324 return GroupStateGuard(*this, group_state);
325 }
326
340 ScopedLoggerGuard pushLogger(bool do_mpi_gather = true) const
341 {
342 return ScopedLoggerGuard(*this, do_mpi_gather);
343 }
344
346 {
347 return WellStateGuard(*this, well_state);
348 }
349
350 int reportStepIdx() const
351 {
352 return report_step_;
353 }
354
355 const Schedule& schedule() const
356 {
357 return this->schedule_;
358 }
359
361 return rescoup_;
362 }
363
365 return rescoup_;
366 }
367
369 return rescoup_.master();
370 }
371
373 return rescoup_.master();
374 }
375
377 return rescoup_.slave();
378 }
380 return rescoup_.slave();
381 }
382
383#ifdef RESERVOIR_COUPLING_ENABLED
384 void setReservoirCouplingMaster(ReservoirCouplingMaster<Scalar>* master) {
385 rescoup_.setMaster(master);
386 }
387 void setReservoirCouplingSlave(ReservoirCouplingSlave<Scalar>* slave) {
388 rescoup_.setSlave(slave);
389 }
390#endif
391
392 void setCmodeGroup(const Group& group);
393
394 template <class AverageRegionalPressureType>
395 void
397 const FieldPropsManager& fp,
398 std::map<std::string, std::unique_ptr<AverageRegionalPressureType>>&
399 regional_average_pressure_calculator) const;
400
401 void setReportStep(int report_step)
402 {
403 report_step_ = report_step;
404 }
405
406 const SummaryState& summaryState() const
407 {
408 return this->summary_state_;
409 }
410
411 Scalar sumSolventRates(const Group& group, const bool is_injector) const;
412
413 Scalar sumWellResRates(const Group& group, const int phase_pos, const bool injector) const;
414
415 Scalar sumWellSurfaceRates(const Group& group, const int phase_pos, const bool injector) const;
416
417 Scalar sumWellPhaseRates(bool res_rates,
418 const Group& group,
419 const int phase_pos,
420 const bool injector,
421 const bool network = false) const;
422
423 bool terminalOutput() const
424 {
425 return this->terminal_output_;
426 }
427
428 template <class RegionalValues>
429 void updateGpMaintTargetForGroups(const Group& group,
430 const RegionalValues& regional_values,
431 const double dt);
432
435 int updateGroupControlledWells(const bool is_production_group,
436 const Phase injection_phase);
437
438 void updateGroupProductionRates(const Group& group);
439
440 void updateGroupTargetReduction(const Group& group,
441 const bool is_injector);
442
444
445 void updateREINForGroups(const Group& group, bool sum_rank);
446
447 void updateReservoirRatesInjectionGroups(const Group& group);
448
450
451 void updateSurfaceRatesInjectionGroups(const Group& group);
452
453 void updateVREPForGroups(const Group& group);
454
455 void updateWellRates(const Group& group,
456 const WellState<Scalar, IndexTraits>& well_state_nupcol,
457 WellState<Scalar, IndexTraits>& well_state) const;
458
460 {
461 return *this->well_state_;
462 }
463
464 void updateWellRatesFromGroupTargetScale(const Scalar scale,
465 const Group& group,
466 bool is_injector,
467 WellState<Scalar, IndexTraits>& well_state) const;
468
470 std::pair<std::optional<std::string>, Scalar>
471 worstOffendingWell(const Group& group,
472 const Group::ProductionCMode& offended_control) const;
473
474private:
475#ifdef RESERVOIR_COUPLING_ENABLED
480 ReservoirCoupling::Phase activePhaseIdxToRescoupPhase_(int phase_pos) const;
481#endif
482
504 template<typename ReductionLambda, typename FractionLambda>
505 Scalar applyReductionsAndFractions_(const std::vector<std::string>& chain,
506 Scalar orig_target,
507 Scalar current_rate_available,
508 std::size_t local_reduction_level,
509 bool is_production_group,
510 Phase injection_phase,
511 ReductionLambda&& local_reduction_lambda,
512 FractionLambda&& local_fraction_lambda,
513 bool do_addback) const;
514
525 Scalar computeAddbackEfficiency_(const std::vector<std::string>& chain,
526 std::size_t local_reduction_level) const;
527
528 std::string controlGroup_(const Group& group) const;
529
530 GuideRate::RateVector getGuideRateVector_(const std::vector<Scalar>& rates) const;
531
542 std::size_t getLocalReductionLevel_(const std::vector<std::string>& chain,
543 bool is_production_group,
544 Phase injection_phase) const;
545
546#ifdef RESERVOIR_COUPLING_ENABLED
547 Scalar getReservoirCouplingMasterGroupRate_(const Group& group,
548 const int phase_pos,
549 ReservoirCoupling::RateKind kind) const;
550#endif
551
552 Scalar getSatelliteRate_(const Group& group,
553 const int phase_pos,
554 const bool res_rates,
555 const bool is_injector) const;
556
560 bool isAutoChokeGroupUnderperforming_(const Group& group) const;
561
563 bool isInGroupChainTopBot_(const std::string& bottom, const std::string& top) const;
564
565 bool isSatelliteGroup_(const Group& group) const;
566
567 Scalar satelliteInjectionRate_(const ScheduleState& sched,
568 const Group& group,
569 const int phase_pos,
570 bool res_rates) const;
571
572 Scalar satelliteProductionRate_(const ScheduleState& sched,
573 const Group& group,
574 const GSatProd::GSatProdGroupProp::Rate rate_comp,
575 bool res_rates) const;
576
577 std::optional<GSatProd::GSatProdGroupProp::Rate> selectRateComponent_(const int phase_pos) const;
578
579 int updateGroupControlledWellsRecursive_(const std::string& group_name,
580 const bool is_production_group,
581 const Phase injection_phase);
582
583 void updateGroupTargetReductionRecursive_(const Group& group,
584 const bool is_injector,
585 std::vector<Scalar>& group_target_reduction);
586
587 const WellState<Scalar, IndexTraits>* well_state_ {nullptr};
588 GroupState<Scalar>* group_state_ {nullptr};
589 const Schedule& schedule_;
590 const SummaryState& summary_state_;
591 const GuideRate& guide_rate_;
592 // NOTE: The deferred logger does not change the object "meaningful" state, so it should be ok to
593 // make it mutable and store a pointer to it here.
594 mutable DeferredLogger* deferred_logger_ {nullptr};
595 // NOTE: The phase usage info seems to be read-only throughout the simulation, so it should be safe
596 // to store a reference to it here.
597 const PhaseUsageInfo<IndexTraits>& phase_usage_info_;
598 const Parallel::Communication& comm_;
599 bool terminal_output_ {false};
600 int report_step_ {0};
601 ReservoirCoupling::Proxy<Scalar> rescoup_{};
602};
603
604// -----------------------------------------------------------------------------
605// Template member function implementations
606// -----------------------------------------------------------------------------
607
608// NOTE: This template member function is defined here in the header because the
609// AverageRegionalPressureType type parameter depends on derived class types that are
610// not available when GroupStateHelper.cpp is compiled. Template functions
611// must be visible at their instantiation point.
612// See GroupStateHelper.cpp for detailed rationale.
613template <typename Scalar, typename IndexTraits>
614template <class AverageRegionalPressureType>
615void
617 const Group& group,
618 const FieldPropsManager& fp,
619 std::map<std::string, std::unique_ptr<AverageRegionalPressureType>>& regional_average_pressure_calculator)
620 const
621{
622 for (const std::string& groupName : group.groups()) {
623 this->setRegionAveragePressureCalculator(this->schedule_.getGroup(groupName, this->report_step_),
624 fp,
625 regional_average_pressure_calculator);
626 }
627 const auto& gpm = group.gpmaint();
628 if (!gpm)
629 return;
630
631 const auto& reg = gpm->region();
632 if (!reg)
633 return;
634
635 if (regional_average_pressure_calculator.count(reg->first) == 0) {
636 const std::string name = (reg->first.rfind("FIP", 0) == 0) ? reg->first : "FIP" + reg->first;
637 const auto& fipnum = fp.get_int(name);
638 regional_average_pressure_calculator[reg->first]
639 = std::make_unique<AverageRegionalPressureType>(fipnum);
640 }
641}
642
643// NOTE: This template member function is defined in the header because the
644// RegionalValues type parameter depends on derived class types that are
645// not available when GroupStateHelper.cpp is compiled. Template functions
646// must be visible at their instantiation point.
647// See GroupStateHelper.cpp for detailed rationale.
648template <typename Scalar, typename IndexTraits>
649template <class RegionalValues>
650void
652 const RegionalValues& regional_values,
653 const double dt)
654{
655 OPM_TIMEFUNCTION();
656 for (const std::string& group_name : group.groups()) {
657 const Group& group_tmp = this->schedule_.getGroup(group_name, this->report_step_);
658 this->updateGpMaintTargetForGroups(group_tmp, regional_values, dt);
659 }
660 const auto& gpm = group.gpmaint();
661 if (!gpm)
662 return;
663
664 const auto& region = gpm->region();
665 if (!region)
666 return;
667 if (this->isReservoirCouplingMasterGroup(group)) {
668 // GPMAINT is not supported for reservoir coupling master groups since master groups do not have
669 // subordinate wells in the master reservoir, so the slaves cannot influence the master reservoir's
670 // average pressure.
671 // Even specifying GPMAINT on a group superior to the master group might not make sense, since if the
672 // superior target is distributed down to the master group with guide rate fractions, adjusting
673 // the master group's target (that is sent to the slave) could only indirectly influence the master
674 // reservoir's average pressure by affecting the guide rate fractions distributed to actual wells
675 // in the master reservoir.
677 std::runtime_error,
678 "GPMAINT is not supported for reservoir coupling master groups.",
679 this->deferredLogger()
680 );
681 return;
682 }
683 else if (this->isReservoirCouplingSlaveGroup(group)) {
684 // GPMAINT is not supported for reservoir coupling slave groups since their targets will be overridden
685 // by the corresponding master group's target anyway.
687 std::runtime_error,
688 "GPMAINT is not supported for reservoir coupling slave groups.",
689 this->deferredLogger()
690 );
691 return;
692 }
693 const auto [name, number] = *region;
694 const Scalar error = gpm->pressure_target() - regional_values.at(name)->pressure(number);
695 Scalar current_rate = 0.0;
696 const auto& pu = this->phase_usage_info_;
697 bool injection = true;
698 Scalar sign = 1.0;
699 switch (gpm->flow_target()) {
700 case GPMaint::FlowTarget::RESV_PROD: {
701 current_rate = -this->groupState().injection_vrep_rate(group.name());
702 injection = false;
703 sign = -1.0;
704 break;
705 }
706 case GPMaint::FlowTarget::RESV_OINJ: {
707 if (pu.phaseIsActive(IndexTraits::oilPhaseIdx)) {
708 const auto io = pu.canonicalToActivePhaseIdx(IndexTraits::oilPhaseIdx);
709 current_rate = this->groupState().injection_reservoir_rates(group.name())[io];
710 }
711 break;
712 }
713 case GPMaint::FlowTarget::RESV_WINJ: {
714 if (pu.phaseIsActive(IndexTraits::waterPhaseIdx)) {
715 const auto iw = pu.canonicalToActivePhaseIdx(IndexTraits::waterPhaseIdx);
716 current_rate = this->groupState().injection_reservoir_rates(group.name())[iw];
717 }
718 break;
719 }
720 case GPMaint::FlowTarget::RESV_GINJ: {
721 if (pu.phaseIsActive(IndexTraits::gasPhaseIdx)) {
722 const auto ig = pu.canonicalToActivePhaseIdx(IndexTraits::gasPhaseIdx);
723 current_rate = this->groupState().injection_reservoir_rates(group.name())[ig];
724 }
725 break;
726 }
727 case GPMaint::FlowTarget::SURF_OINJ: {
728 if (pu.phaseIsActive(IndexTraits::oilPhaseIdx)) {
729 const auto io = pu.canonicalToActivePhaseIdx(IndexTraits::oilPhaseIdx);
730 current_rate = this->groupState().injection_surface_rates(group.name())[io];
731 }
732 break;
733 }
734 case GPMaint::FlowTarget::SURF_WINJ: {
735 if (pu.phaseIsActive(IndexTraits::waterPhaseIdx)) {
736 const auto iw = pu.canonicalToActivePhaseIdx(IndexTraits::waterPhaseIdx);
737 current_rate = this->groupState().injection_surface_rates(group.name())[iw];
738 }
739 break;
740 }
741 case GPMaint::FlowTarget::SURF_GINJ: {
742 if (pu.phaseIsActive(IndexTraits::gasPhaseIdx)) {
743 const auto ig = pu.canonicalToActivePhaseIdx(IndexTraits::gasPhaseIdx);
744 current_rate = this->groupState().injection_surface_rates(group.name())[ig];
745 }
746 break;
747 }
748 default:
749 throw std::invalid_argument("Invalid Flow target type in GPMAINT");
750 }
751 auto& gpmaint_state = this->groupState().gpmaint(group.name());
752 // we only activate gpmaint if pressure is lower than the target regional pressure for injectors
753 // (i.e. error > 0) and higher for producers.
754 bool activate = (injection && error > 0) || (!injection && error < 0);
755 Scalar rate = 0.0;
756 if (activate) {
757 rate = gpm->rate(gpmaint_state, current_rate, error, dt);
758 } else {
759 gpm->resetState(gpmaint_state);
760 }
761 this->groupState().update_gpmaint_target(group.name(), std::max(Scalar {0.0}, sign * rate));
762}
763
764} // namespace Opm
765
766#endif // OPM_GROUPSTATE_HELPER_HEADER_INCLUDED
#define OPM_DEFLOG_THROW(Exception, message, deferred_logger)
Definition: DeferredLoggingErrorHelpers.hpp:45
Definition: DeferredLogger.hpp:57
Definition: GroupStateHelper.hpp:88
GroupStateGuard(GroupStateGuard &&)=delete
GroupStateGuard & operator=(GroupStateGuard &&)=delete
GroupStateGuard(GroupStateHelper &group_state_helper, GroupState< Scalar > &group_state)
Definition: GroupStateHelper.hpp:90
GroupStateGuard(const GroupStateGuard &)=delete
~GroupStateGuard()
Definition: GroupStateHelper.hpp:98
GroupStateGuard & operator=(const GroupStateGuard &)=delete
RAII guard that owns a DeferredLogger and auto-gathers on destruction.
Definition: GroupStateHelper.hpp:128
ScopedLoggerGuard & operator=(const ScopedLoggerGuard &)=delete
ScopedLoggerGuard(ScopedLoggerGuard &&other) noexcept
Definition: GroupStateHelper.hpp:173
ScopedLoggerGuard(const ScopedLoggerGuard &)=delete
ScopedLoggerGuard(const GroupStateHelper &helper, bool do_mpi_gather=true)
Constructor for scoped logger guard.
Definition: GroupStateHelper.hpp:135
~ScopedLoggerGuard()
Definition: GroupStateHelper.hpp:143
ScopedLoggerGuard & operator=(ScopedLoggerGuard &&)=delete
Definition: GroupStateHelper.hpp:59
WellStateGuard(WellStateGuard &&)=delete
WellStateGuard(const WellStateGuard &)=delete
~WellStateGuard()
Definition: GroupStateHelper.hpp:69
WellStateGuard & operator=(const WellStateGuard &)=delete
WellStateGuard & operator=(WellStateGuard &&)=delete
WellStateGuard(GroupStateHelper &groupStateHelper, WellState< Scalar, IndexTraits > &well_state)
Definition: GroupStateHelper.hpp:61
Definition: GroupStateHelper.hpp:55
void updateWellRatesFromGroupTargetScale(const Scalar scale, const Group &group, bool is_injector, WellState< Scalar, IndexTraits > &well_state) const
std::pair< Scalar, Group::ProductionCMode > getAutoChokeGroupProductionTargetRate(const Group &bottom_group, const Group &group, const std::vector< Scalar > &resv_coeff, Scalar efficiencyFactor) 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:616
void setCmodeGroup(const Group &group)
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
std::optional< GroupTarget > 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) const
ReservoirCouplingMaster< Scalar > & reservoirCouplingMaster()
Definition: GroupStateHelper.hpp:368
GroupState< Scalar > & groupState() const
Definition: GroupStateHelper.hpp:289
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:318
const ReservoirCouplingMaster< Scalar > & reservoirCouplingMaster() const
Definition: GroupStateHelper.hpp:372
const SummaryState & summaryState() const
Definition: GroupStateHelper.hpp:406
const GuideRate & guideRate() const
Definition: GroupStateHelper.hpp:294
std::pair< std::optional< std::string >, Scalar > worstOffendingWell(const Group &group, const Group::ProductionCMode &offended_control) const
Returns the name of the worst offending well and its fraction (i.e. violated_phase / preferred_phase)
ReservoirCoupling::Proxy< Scalar > & rescoup()
Definition: GroupStateHelper.hpp:360
void updateReservoirRatesInjectionGroups(const Group &group)
bool isReservoirCouplingSlaveGroup(const Group &group) const
Definition: GroupStateHelper.hpp:310
const WellState< Scalar, IndexTraits > & wellState() const
Definition: GroupStateHelper.hpp:459
constexpr int numPhases() const
Definition: GroupStateHelper.hpp:312
GuideRate::RateVector getWellRateVector(const std::string &name) const
ScopedLoggerGuard pushLogger(bool do_mpi_gather=true) const
Push a new logger onto the stack with auto-cleanup on destruction.
Definition: GroupStateHelper.hpp:340
const Parallel::Communication & comm() const
Definition: GroupStateHelper.hpp:223
GuideRateModel::Target getInjectionGuideTargetMode(Phase injection_phase) const
Get the guide rate target mode for an injection phase.
bool isRank0() const
Definition: GroupStateHelper.hpp:299
Scalar sumSolventRates(const Group &group, const bool is_injector) const
const Schedule & schedule() const
Definition: GroupStateHelper.hpp:355
void updateGroupTargetReduction(const Group &group, const bool is_injector)
GuideRate::RateVector getProductionGroupRateVector(const std::string &group_name) const
void updateGroupProductionRates(const Group &group)
bool terminalOutput() const
Definition: GroupStateHelper.hpp:423
typename SingleWellState< Scalar, IndexTraits >::GroupTarget GroupTarget
Definition: GroupStateHelper.hpp:261
void updateGpMaintTargetForGroups(const Group &group, const RegionalValues &regional_values, const double dt)
Definition: GroupStateHelper.hpp:651
int phaseToActivePhaseIdx(const Phase phase) 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) const
DeferredLogger & deferredLogger() const
Get the deferred logger.
Definition: GroupStateHelper.hpp:227
Scalar sumWellSurfaceRates(const Group &group, const int phase_pos, const bool injector) const
int reportStepIdx() const
Definition: GroupStateHelper.hpp:350
void setReportStep(int report_step)
Definition: GroupStateHelper.hpp:401
int updateGroupControlledWells(const bool is_production_group, const Phase injection_phase)
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:306
const ReservoirCoupling::Proxy< Scalar > & rescoup() const
Definition: GroupStateHelper.hpp:364
bool isReservoirCouplingMaster() const
Definition: GroupStateHelper.hpp:304
const ReservoirCouplingSlave< Scalar > & reservoirCouplingSlave() const
Definition: GroupStateHelper.hpp:379
void updateNetworkLeafNodeProductionRates()
WellStateGuard pushWellState(WellState< Scalar, IndexTraits > &well_state)
Definition: GroupStateHelper.hpp:345
Scalar getInjectionGroupTarget(const Group &group, const Phase &injection_phase, const std::vector< Scalar > &resv_coeff) const
GroupStateGuard pushGroupState(GroupState< Scalar > &group_state)
Definition: GroupStateHelper.hpp:322
ReservoirCouplingSlave< Scalar > & reservoirCouplingSlave()
Definition: GroupStateHelper.hpp:376
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:308
GuideRateModel::Target getProductionGuideTargetMode(const Group &group) const
Get the guide rate target mode for a production group.
Scalar getProductionGroupTarget(const Group &group) 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) const
std::optional< GroupTarget > 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) 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, bool terminal_output)
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
RateKind
Selects which kind of rate to retrieve from slave group data.
Definition: ReservoirCoupling.hpp:168
Phase
Phase indices for reservoir coupling, we currently only support black-oil phases (oil,...
Definition: ReservoirCoupling.hpp:156
Definition: blackoilbioeffectsmodules.hh:43
Opm::DeferredLogger gatherDeferredLogger(const Opm::DeferredLogger &local_deferredlogger, Parallel::Communication communicator)
Create a global log combining local logs.
Definition: SingleWellState.hpp:121