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/Schedule/ResCoup/GrupSlav.hpp>
26#include <opm/input/eclipse/EclipseState/Grid/FieldPropsManager.hpp>
27#include <opm/input/eclipse/Schedule/Group/GPMaint.hpp>
28#include <opm/input/eclipse/Schedule/Group/GSatProd.hpp>
29#include <opm/input/eclipse/Schedule/Group/GuideRate.hpp>
30#include <opm/input/eclipse/Schedule/Schedule.hpp>
31#include <opm/input/eclipse/Schedule/ScheduleState.hpp>
32#include <opm/input/eclipse/Schedule/SummaryState.hpp>
33#include <opm/material/fluidsystems/PhaseUsageInfo.hpp>
40
42
43#include <algorithm>
44#include <map>
45#include <memory>
46#include <optional>
47#include <stdexcept>
48#include <string>
49#include <vector>
50
51namespace Opm
52{
53
54template <typename Scalar, typename IndexTraits>
56{
57public:
58 // RAII guard for temporarily setting wellstate pointer
60 {
61 public:
63 : groupStateHelper_ {groupStateHelper}
64 , previous_state_ptr_ {groupStateHelper_.well_state_}
65 {
66 // Set the new state directly
67 groupStateHelper_.well_state_ = &well_state;
68 }
69
71 {
72 // Restore the previous state
73 groupStateHelper_.well_state_ = previous_state_ptr_;
74 }
75
76 // Delete copy and move operations
81
82 private:
83 GroupStateHelper& groupStateHelper_;
84 const WellState<Scalar, IndexTraits>* previous_state_ptr_;
85 };
86
87 // RAII guard for temporarily setting groupstate pointer
89 {
90 public:
91 GroupStateGuard(GroupStateHelper& group_state_helper, GroupState<Scalar>& group_state)
92 : group_state_helper_ {group_state_helper}
93 , previous_state_ptr_ {group_state_helper.group_state_}
94 {
95 // Set the new state directly
96 group_state_helper_.group_state_ = &group_state;
97 }
98
100 {
101 // Restore the previous state
102 group_state_helper_.group_state_ = previous_state_ptr_;
103 }
104
105 // Delete copy and move operations
110
111 private:
112 GroupStateHelper& group_state_helper_;
113 GroupState<Scalar>* previous_state_ptr_ {nullptr};
114 };
115
129 {
130 public:
136 explicit ScopedLoggerGuard(const GroupStateHelper& helper, bool do_mpi_gather = true)
137 : helper_(&helper)
138 , previous_(helper.deferred_logger_)
139 , do_mpi_gather_(do_mpi_gather)
140 {
141 helper_->deferred_logger_ = &logger_;
142 }
143
145 {
146 if (helper_) {
147 if (do_mpi_gather_) {
148 // 1. Gather messages across MPI ranks
149 DeferredLogger global = gatherDeferredLogger(logger_, helper_->comm());
150
151 // 2. Log on rank 0 (if terminal_output enabled)
152 if (helper_->terminalOutput()) {
153 global.logMessages();
154 }
155 } else {
156 // Just log locally without MPI gather
157 if (helper_->terminalOutput()) {
158 logger_.logMessages();
159 }
160 }
161
162 // 3. Restore previous logger
163 helper_->deferred_logger_ = previous_;
164 }
165 }
166
167 // Delete copy operations and move assignment
171
172 // Move constructor required for pushLogger() return-by-value (must be
173 // available even if elided by RVO)
175 : helper_(other.helper_)
176 , logger_(std::move(other.logger_))
177 , previous_(other.previous_)
178 , do_mpi_gather_(other.do_mpi_gather_)
179 {
180 // Update the helper's pointer to our moved logger
181 if (helper_) {
182 helper_->deferred_logger_ = &logger_;
183 }
184 other.helper_ = nullptr;
185 other.previous_ = nullptr;
186 }
187
188 private:
189 // Pointer (not reference) to allow nulling in move constructor
190 const GroupStateHelper* helper_{nullptr};
191 DeferredLogger logger_; // Owned logger
192 DeferredLogger* previous_{nullptr}; // For restore
193 bool do_mpi_gather_{true}; // Whether to gather messages across MPI ranks
194 };
195
197
199 GroupState<Scalar>& group_state,
200 const Schedule& schedule,
201 const SummaryState& summary_state,
202 const GuideRate& guide_rate,
203 const PhaseUsageInfo<IndexTraits>& phase_usage_info,
205 bool terminal_output);
206
207 void accumulateGroupEfficiencyFactor(const Group& group, Scalar& factor) const;
208
209 std::pair<bool, Scalar> checkGroupConstraintsInj(const std::string& name,
210 const std::string& parent,
211 const Group& group,
212 const Scalar* rates,
213 const Phase injection_phase,
214 const Scalar efficiency_factor,
215 const std::vector<Scalar>& resv_coeff,
216 const bool check_guide_rate) const;
217
218 std::pair<bool, Scalar> checkGroupConstraintsProd(const std::string& name,
219 const std::string& parent,
220 const Group& group,
221 const Scalar* rates,
222 const Scalar efficiency_factor,
223 const std::vector<Scalar>& resv_coeff,
224 const bool check_guide_rate) const;
225
226 std::pair<Group::ProductionCMode, Scalar>
227 checkGroupProductionConstraints(const Group& group) const;
228
229 const Parallel::Communication& comm() const { return this->comm_; }
230
234 {
235 if (this->deferred_logger_ == nullptr) {
236 throw std::logic_error("DeferredLogger not set. Call pushLogger() first.");
237 }
238 return *this->deferred_logger_;
239 }
240
241 std::vector<Scalar> getGroupRatesAvailableForHigherLevelControl(const Group& group, const bool is_injector) const;
242
243 Scalar getInjectionGroupTarget(const Group& group,
244 const Phase& injection_phase,
245 const std::vector<Scalar>& resv_coeff) const;
246
250 GuideRateModel::Target getInjectionGuideTargetMode(Phase injection_phase) const;
251
252 Scalar getProductionGroupTarget(const Group& group) const;
253
255 Scalar getProductionGroupTargetForMode(const Group& group,
256 const Group::ProductionCMode cmode) const;
257
261 GuideRateModel::Target getProductionGuideTargetMode(const Group& group) const;
262
263 std::pair<Scalar, Group::ProductionCMode>
264 getAutoChokeGroupProductionTargetRate(const Group& bottom_group,
265 const Group& group,
266 const std::vector<Scalar>& resv_coeff,
267 Scalar efficiencyFactor) const;
268
269 GuideRate::RateVector getProductionGroupRateVector(const std::string& group_name) const;
270
271 std::optional<GroupTarget> getWellGroupTargetInjector(const std::string& name,
272 const std::string& parent,
273 const Group& group,
274 const Scalar* rates,
275 const Phase injection_phase,
276 const Scalar efficiency_factor,
277 const std::vector<Scalar>& resv_coeff) const;
278
279 std::optional<GroupTarget> getWellGroupTargetProducer(const std::string& name,
280 const std::string& parent,
281 const Group& group,
282 const Scalar* rates,
283 const Scalar efficiency_factor,
284 const std::vector<Scalar>& resv_coeff) const;
285
286 GuideRate::RateVector getWellRateVector(const std::string& name) const;
287
288 std::vector<std::string> groupChainTopBot(const std::string& bottom, const std::string& top) const;
289
292 int groupControlledWells(const std::string& group_name,
293 const std::string& always_included_child,
294 const bool is_production_group,
295 const Phase injection_phase) const;
296
298 {
299 return *this->group_state_;
300 }
301
302 const GuideRate& guideRate() const
303 {
304 return this->guide_rate_;
305 }
306
307 bool isRank0() const
308 {
309 return this->well_state_->isRank0();
310 }
311
312 bool isReservoirCouplingMaster() const { return rescoup_.isMaster(); }
313
314 bool isReservoirCouplingMasterGroup(const Group& group) const { return rescoup_.isMasterGroup(group.name()); }
315
316 bool isReservoirCouplingSlave() const { return rescoup_.isSlave(); }
317
318 bool isReservoirCouplingSlaveGroup(const Group& group) const { return rescoup_.isSlaveGroup(group.name()); }
319
320 constexpr int numPhases() const {
321 return this->wellState().numPhases();
322 }
323
324 int phaseToActivePhaseIdx(const Phase phase) const;
325
327 return this->phase_usage_info_;
328 }
329
331 {
332 return GroupStateGuard(*this, group_state);
333 }
334
348 ScopedLoggerGuard pushLogger(bool do_mpi_gather = true) const
349 {
350 return ScopedLoggerGuard(*this, do_mpi_gather);
351 }
352
354 {
355 return WellStateGuard(*this, well_state);
356 }
357
358 int reportStepIdx() const
359 {
360 return report_step_;
361 }
362
363 const Schedule& schedule() const
364 {
365 return this->schedule_;
366 }
367
369 return rescoup_;
370 }
371
373 return rescoup_;
374 }
375
377 return rescoup_.master();
378 }
379
381 return rescoup_.master();
382 }
383
385 return rescoup_.slave();
386 }
388 return rescoup_.slave();
389 }
390
391#ifdef RESERVOIR_COUPLING_ENABLED
392 void setReservoirCouplingMaster(ReservoirCouplingMaster<Scalar>* master) {
393 rescoup_.setMaster(master);
394 }
395 void setReservoirCouplingSlave(ReservoirCouplingSlave<Scalar>* slave) {
396 rescoup_.setSlave(slave);
397 }
398#endif
399
400 void setCmodeGroup(const Group& group);
401
402 template <class AverageRegionalPressureType>
403 void
405 const FieldPropsManager& fp,
406 std::map<std::string, std::unique_ptr<AverageRegionalPressureType>>&
407 regional_average_pressure_calculator) const;
408
409 void setReportStep(int report_step)
410 {
411 report_step_ = report_step;
412 }
413
414 const SummaryState& summaryState() const
415 {
416 return this->summary_state_;
417 }
418
419 Scalar sumSolventRates(const Group& group, const bool is_injector) const;
420
421 Scalar sumWellResRates(const Group& group, const int phase_pos, const bool injector) const;
422
423 Scalar sumWellSurfaceRates(const Group& group, const int phase_pos, const bool injector) const;
424
425 Scalar sumWellPhaseRates(bool res_rates,
426 const Group& group,
427 const int phase_pos,
428 const bool injector,
429 const bool network = false) const;
430
431 bool terminalOutput() const
432 {
433 return this->terminal_output_;
434 }
435
436 template <class RegionalValues>
437 void updateGpMaintTargetForGroups(const Group& group,
438 const RegionalValues& regional_values,
439 const double dt);
440
443 int updateGroupControlledWells(const bool is_production_group,
444 const Phase injection_phase);
445
446 void updateGroupProductionRates(const Group& group);
447
448 void updateGroupTargetReduction(const Group& group,
449 const bool is_injector);
450
452
453 void updateREINForGroups(const Group& group, bool sum_rank);
454
455 void updateReservoirRatesInjectionGroups(const Group& group);
456
457#ifdef RESERVOIR_COUPLING_ENABLED
467 void updateSlaveGroupCmodesFromMaster();
468#endif
469
471
472 void updateSurfaceRatesInjectionGroups(const Group& group);
473
474 void updateVREPForGroups(const Group& group);
475
476 void updateWellRates(const Group& group,
477 const WellState<Scalar, IndexTraits>& well_state_nupcol,
478 WellState<Scalar, IndexTraits>& well_state) const;
479
481 {
482 return *this->well_state_;
483 }
484
485 void updateWellRatesFromGroupTargetScale(const Scalar scale,
486 const Group& group,
487 bool is_injector,
488 WellState<Scalar, IndexTraits>& well_state) const;
489
491 std::pair<std::optional<std::string>, Scalar>
492 worstOffendingWell(const Group& group,
493 const Group::ProductionCMode& offended_control) const;
494
495private:
496#ifdef RESERVOIR_COUPLING_ENABLED
501 ReservoirCoupling::Phase activePhaseIdxToRescoupPhase_(int phase_pos) const;
502#endif
503
525 template<typename ReductionLambda, typename FractionLambda>
526 Scalar applyReductionsAndFractions_(const std::vector<std::string>& chain,
527 Scalar orig_target,
528 Scalar current_rate_available,
529 std::size_t local_reduction_level,
530 bool is_production_group,
531 Phase injection_phase,
532 ReductionLambda&& local_reduction_lambda,
533 FractionLambda&& local_fraction_lambda,
534 bool do_addback) const;
535
536 std::pair<Group::ProductionCMode, Scalar>
537 checkProductionRateConstraint_(const Group& group,
538 Group::ProductionCMode cmode,
539 Group::ProductionCMode currentControl,
540 Scalar target,
541 Scalar current_rate) const;
542
553 Scalar computeAddbackEfficiency_(const std::vector<std::string>& chain,
554 std::size_t local_reduction_level) const;
555
556 std::string controlGroup_(const Group& group) const;
557
558 #ifdef RESERVOIR_COUPLING_ENABLED
572 Scalar getEffectiveProductionLimit_(const std::string& gname,
573 Group::ProductionCMode rate_type,
574 Scalar slave_local_target) const;
575#endif // RESERVOIR_COUPLING_ENABLED
576
577 GuideRate::RateVector getGuideRateVector_(const std::vector<Scalar>& rates) const;
578
579 Scalar getInjectionGroupTargetForMode_(const Group& group,
580 const Phase& injection_phase,
581 const std::vector<Scalar>& resv_coeff,
582 const Group::InjectionCMode cmode) const;
583
594 std::size_t getLocalReductionLevel_(const std::vector<std::string>& chain,
595 bool is_production_group,
596 Phase injection_phase) const;
597
598#ifdef RESERVOIR_COUPLING_ENABLED
599 ReservoirCoupling::GrupSlav::FilterFlag getInjectionFilterFlag_(const std::string& group_name,
600 const Phase injection_phase) const;
601
602 ReservoirCoupling::GrupSlav::FilterFlag getProductionFilterFlag_(const std::string& group_name,
603 const Group::ProductionCMode cmode) const;
604
605 Scalar getReservoirCouplingMasterGroupRate_(const Group& group,
606 const int phase_pos,
607 const ReservoirCoupling::RateKind kind) const;
608#endif
609
610 Scalar getProductionConstraintTarget_(const Group& group,
611 Group::ProductionCMode cmode,
612 const Group::ProductionControls& controls) const;
613
614 Scalar getProductionGroupTargetForMode_(const Group& group, const Group::ProductionCMode cmode) const;
615
616 Scalar getSatelliteRate_(const Group& group,
617 const int phase_pos,
618 const bool res_rates,
619 const bool is_injector) const;
620
624 bool isAutoChokeGroupUnderperforming_(const Group& group) const;
625
627 bool isInGroupChainTopBot_(const std::string& bottom, const std::string& top) const;
628
629 bool isSatelliteGroup_(const Group& group) const;
630
631 Scalar satelliteInjectionRate_(const ScheduleState& sched,
632 const Group& group,
633 const int phase_pos,
634 bool res_rates) const;
635
636 Scalar satelliteProductionRate_(const ScheduleState& sched,
637 const Group& group,
638 const GSatProd::GSatProdGroupProp::Rate rate_comp,
639 bool res_rates) const;
640
641 std::optional<GSatProd::GSatProdGroupProp::Rate> selectRateComponent_(const int phase_pos) const;
642
653 Scalar subtractOtherPhaseResvInjection_(
654 Phase injection_phase,
655 Scalar base_reservoir_rate,
656 const std::vector<Scalar>& group_injection_reservoir_rates) const;
657
658 Scalar sumProductionRateForControlMode_(const Group& group, Group::ProductionCMode cmode) const;
659
660 int updateGroupControlledWellsRecursive_(const std::string& group_name,
661 const bool is_production_group,
662 const Phase injection_phase);
663
664 void updateGroupTargetReductionRecursive_(const Group& group,
665 const bool is_injector,
666 std::vector<Scalar>& group_target_reduction);
667
668 const WellState<Scalar, IndexTraits>* well_state_ {nullptr};
669 GroupState<Scalar>* group_state_ {nullptr};
670 const Schedule& schedule_;
671 const SummaryState& summary_state_;
672 const GuideRate& guide_rate_;
673 // NOTE: The deferred logger does not change the object "meaningful" state, so it should be ok to
674 // make it mutable and store a pointer to it here.
675 mutable DeferredLogger* deferred_logger_ {nullptr};
676 // NOTE: The phase usage info seems to be read-only throughout the simulation, so it should be safe
677 // to store a reference to it here.
678 const PhaseUsageInfo<IndexTraits>& phase_usage_info_;
679 const Parallel::Communication& comm_;
680 bool terminal_output_ {false};
681 int report_step_ {0};
682 ReservoirCoupling::Proxy<Scalar> rescoup_{};
683};
684
685// -----------------------------------------------------------------------------
686// Template member function implementations
687// -----------------------------------------------------------------------------
688
689// NOTE: This template member function is defined here in the header because the
690// AverageRegionalPressureType type parameter depends on derived class types that are
691// not available when GroupStateHelper.cpp is compiled. Template functions
692// must be visible at their instantiation point.
693// See GroupStateHelper.cpp for detailed rationale.
694template <typename Scalar, typename IndexTraits>
695template <class AverageRegionalPressureType>
696void
698 const Group& group,
699 const FieldPropsManager& fp,
700 std::map<std::string, std::unique_ptr<AverageRegionalPressureType>>& regional_average_pressure_calculator)
701 const
702{
703 for (const std::string& groupName : group.groups()) {
704 this->setRegionAveragePressureCalculator(this->schedule_.getGroup(groupName, this->report_step_),
705 fp,
706 regional_average_pressure_calculator);
707 }
708 const auto& gpm = group.gpmaint();
709 if (!gpm)
710 return;
711
712 const auto& reg = gpm->region();
713 if (!reg)
714 return;
715
716 if (regional_average_pressure_calculator.count(reg->first) == 0) {
717 const std::string name = (reg->first.rfind("FIP", 0) == 0) ? reg->first : "FIP" + reg->first;
718 const auto& fipnum = fp.get_int(name);
719 regional_average_pressure_calculator[reg->first]
720 = std::make_unique<AverageRegionalPressureType>(fipnum);
721 }
722}
723
724// NOTE: This template member function is defined in the header because the
725// RegionalValues type parameter depends on derived class types that are
726// not available when GroupStateHelper.cpp is compiled. Template functions
727// must be visible at their instantiation point.
728// See GroupStateHelper.cpp for detailed rationale.
729template <typename Scalar, typename IndexTraits>
730template <class RegionalValues>
731void
733 const RegionalValues& regional_values,
734 const double dt)
735{
736 OPM_TIMEFUNCTION();
737 for (const std::string& group_name : group.groups()) {
738 const Group& group_tmp = this->schedule_.getGroup(group_name, this->report_step_);
739 this->updateGpMaintTargetForGroups(group_tmp, regional_values, dt);
740 }
741 const auto& gpm = group.gpmaint();
742 if (!gpm)
743 return;
744
745 const auto& region = gpm->region();
746 if (!region)
747 return;
748 if (this->isReservoirCouplingMasterGroup(group)) {
749 // GPMAINT is not supported for reservoir coupling master groups since master groups do not have
750 // subordinate wells in the master reservoir, so the slaves cannot influence the master reservoir's
751 // average pressure.
752 // Even specifying GPMAINT on a group superior to the master group might not make sense, since if the
753 // superior target is distributed down to the master group with guide rate fractions, adjusting
754 // the master group's target (that is sent to the slave) could only indirectly influence the master
755 // reservoir's average pressure by affecting the guide rate fractions distributed to actual wells
756 // in the master reservoir.
758 std::runtime_error,
759 "GPMAINT is not supported for reservoir coupling master groups.",
760 this->deferredLogger()
761 );
762 return;
763 }
764 else if (this->isReservoirCouplingSlaveGroup(group)) {
765 // GPMAINT is not supported for reservoir coupling slave groups since their targets will be overridden
766 // by the corresponding master group's target anyway.
768 std::runtime_error,
769 "GPMAINT is not supported for reservoir coupling slave groups.",
770 this->deferredLogger()
771 );
772 return;
773 }
774 const auto [name, number] = *region;
775 const Scalar error = gpm->pressure_target() - regional_values.at(name)->pressure(number);
776 Scalar current_rate = 0.0;
777 const auto& pu = this->phase_usage_info_;
778 bool injection = true;
779 Scalar sign = 1.0;
780 switch (gpm->flow_target()) {
781 case GPMaint::FlowTarget::RESV_PROD: {
782 current_rate = -this->groupState().injection_vrep_rate(group.name());
783 injection = false;
784 sign = -1.0;
785 break;
786 }
787 case GPMaint::FlowTarget::RESV_OINJ: {
788 if (pu.phaseIsActive(IndexTraits::oilPhaseIdx)) {
789 const auto io = pu.canonicalToActivePhaseIdx(IndexTraits::oilPhaseIdx);
790 current_rate = this->groupState().injection_reservoir_rates(group.name())[io];
791 }
792 break;
793 }
794 case GPMaint::FlowTarget::RESV_WINJ: {
795 if (pu.phaseIsActive(IndexTraits::waterPhaseIdx)) {
796 const auto iw = pu.canonicalToActivePhaseIdx(IndexTraits::waterPhaseIdx);
797 current_rate = this->groupState().injection_reservoir_rates(group.name())[iw];
798 }
799 break;
800 }
801 case GPMaint::FlowTarget::RESV_GINJ: {
802 if (pu.phaseIsActive(IndexTraits::gasPhaseIdx)) {
803 const auto ig = pu.canonicalToActivePhaseIdx(IndexTraits::gasPhaseIdx);
804 current_rate = this->groupState().injection_reservoir_rates(group.name())[ig];
805 }
806 break;
807 }
808 case GPMaint::FlowTarget::SURF_OINJ: {
809 if (pu.phaseIsActive(IndexTraits::oilPhaseIdx)) {
810 const auto io = pu.canonicalToActivePhaseIdx(IndexTraits::oilPhaseIdx);
811 current_rate = this->groupState().injection_surface_rates(group.name())[io];
812 }
813 break;
814 }
815 case GPMaint::FlowTarget::SURF_WINJ: {
816 if (pu.phaseIsActive(IndexTraits::waterPhaseIdx)) {
817 const auto iw = pu.canonicalToActivePhaseIdx(IndexTraits::waterPhaseIdx);
818 current_rate = this->groupState().injection_surface_rates(group.name())[iw];
819 }
820 break;
821 }
822 case GPMaint::FlowTarget::SURF_GINJ: {
823 if (pu.phaseIsActive(IndexTraits::gasPhaseIdx)) {
824 const auto ig = pu.canonicalToActivePhaseIdx(IndexTraits::gasPhaseIdx);
825 current_rate = this->groupState().injection_surface_rates(group.name())[ig];
826 }
827 break;
828 }
829 default:
830 throw std::invalid_argument("Invalid Flow target type in GPMAINT");
831 }
832 auto& gpmaint_state = this->groupState().gpmaint(group.name());
833 // we only activate gpmaint if pressure is lower than the target regional pressure for injectors
834 // (i.e. error > 0) and higher for producers.
835 bool activate = (injection && error > 0) || (!injection && error < 0);
836 Scalar rate = 0.0;
837 if (activate) {
838 rate = gpm->rate(gpmaint_state, current_rate, error, dt);
839 } else {
840 gpm->resetState(gpmaint_state);
841 }
842 this->groupState().update_gpmaint_target(group.name(), std::max(Scalar {0.0}, sign * rate));
843}
844
845} // namespace Opm
846
847#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:89
GroupStateGuard(GroupStateGuard &&)=delete
GroupStateGuard & operator=(GroupStateGuard &&)=delete
GroupStateGuard(GroupStateHelper &group_state_helper, GroupState< Scalar > &group_state)
Definition: GroupStateHelper.hpp:91
GroupStateGuard(const GroupStateGuard &)=delete
~GroupStateGuard()
Definition: GroupStateHelper.hpp:99
GroupStateGuard & operator=(const GroupStateGuard &)=delete
RAII guard that owns a DeferredLogger and auto-gathers on destruction.
Definition: GroupStateHelper.hpp:129
ScopedLoggerGuard & operator=(const ScopedLoggerGuard &)=delete
ScopedLoggerGuard(ScopedLoggerGuard &&other) noexcept
Definition: GroupStateHelper.hpp:174
ScopedLoggerGuard(const ScopedLoggerGuard &)=delete
ScopedLoggerGuard(const GroupStateHelper &helper, bool do_mpi_gather=true)
Constructor for scoped logger guard.
Definition: GroupStateHelper.hpp:136
~ScopedLoggerGuard()
Definition: GroupStateHelper.hpp:144
ScopedLoggerGuard & operator=(ScopedLoggerGuard &&)=delete
Definition: GroupStateHelper.hpp:60
WellStateGuard(WellStateGuard &&)=delete
WellStateGuard(const WellStateGuard &)=delete
~WellStateGuard()
Definition: GroupStateHelper.hpp:70
WellStateGuard & operator=(const WellStateGuard &)=delete
WellStateGuard & operator=(WellStateGuard &&)=delete
WellStateGuard(GroupStateHelper &groupStateHelper, WellState< Scalar, IndexTraits > &well_state)
Definition: GroupStateHelper.hpp:62
Definition: GroupStateHelper.hpp:56
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:697
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:376
GroupState< Scalar > & groupState() const
Definition: GroupStateHelper.hpp:297
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:326
const ReservoirCouplingMaster< Scalar > & reservoirCouplingMaster() const
Definition: GroupStateHelper.hpp:380
const SummaryState & summaryState() const
Definition: GroupStateHelper.hpp:414
const GuideRate & guideRate() const
Definition: GroupStateHelper.hpp:302
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:368
void updateReservoirRatesInjectionGroups(const Group &group)
bool isReservoirCouplingSlaveGroup(const Group &group) const
Definition: GroupStateHelper.hpp:318
const WellState< Scalar, IndexTraits > & wellState() const
Definition: GroupStateHelper.hpp:480
constexpr int numPhases() const
Definition: GroupStateHelper.hpp:320
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:348
const Parallel::Communication & comm() const
Definition: GroupStateHelper.hpp:229
GuideRateModel::Target getInjectionGuideTargetMode(Phase injection_phase) const
Get the guide rate target mode for an injection phase.
bool isRank0() const
Definition: GroupStateHelper.hpp:307
std::pair< Group::ProductionCMode, Scalar > checkGroupProductionConstraints(const Group &group) const
Scalar sumSolventRates(const Group &group, const bool is_injector) const
const Schedule & schedule() const
Definition: GroupStateHelper.hpp:363
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:431
typename SingleWellState< Scalar, IndexTraits >::GroupTarget GroupTarget
Definition: GroupStateHelper.hpp:196
void updateGpMaintTargetForGroups(const Group &group, const RegionalValues &regional_values, const double dt)
Definition: GroupStateHelper.hpp:732
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:233
Scalar sumWellSurfaceRates(const Group &group, const int phase_pos, const bool injector) const
int reportStepIdx() const
Definition: GroupStateHelper.hpp:358
void setReportStep(int report_step)
Definition: GroupStateHelper.hpp:409
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:314
const ReservoirCoupling::Proxy< Scalar > & rescoup() const
Definition: GroupStateHelper.hpp:372
bool isReservoirCouplingMaster() const
Definition: GroupStateHelper.hpp:312
const ReservoirCouplingSlave< Scalar > & reservoirCouplingSlave() const
Definition: GroupStateHelper.hpp:387
void updateNetworkLeafNodeProductionRates()
Scalar getProductionGroupTargetForMode(const Group &group, const Group::ProductionCMode cmode) const
Get the production target for a specific control mode (not necessarily the active one).
WellStateGuard pushWellState(WellState< Scalar, IndexTraits > &well_state)
Definition: GroupStateHelper.hpp:353
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:330
ReservoirCouplingSlave< Scalar > & reservoirCouplingSlave()
Definition: GroupStateHelper.hpp:384
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:316
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:38
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:45
Opm::DeferredLogger gatherDeferredLogger(const Opm::DeferredLogger &local_deferredlogger, Parallel::Communication communicator)
Create a global log combining local logs.
Definition: SingleWellState.hpp:121