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 <unordered_set>
50#include <vector>
51
52namespace Opm
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 GuideRateModel::Target
264 getProductionGuideTargetModeFromControlMode(const Group::ProductionCMode cmode) const;
265
266 std::pair<Scalar, Group::ProductionCMode>
267 getAutoChokeGroupProductionTargetRate(const Group& bottom_group,
268 const Group& group,
269 const std::vector<Scalar>& resv_coeff,
270 Scalar efficiencyFactor) const;
271
272 GuideRate::RateVector getProductionGroupRateVector(const std::string& group_name) const;
273
274 std::optional<GroupTarget> getWellGroupTargetInjector(const std::string& name,
275 const std::string& parent,
276 const Group& group,
277 const Scalar* rates,
278 const Phase injection_phase,
279 const Scalar efficiency_factor,
280 const std::vector<Scalar>& resv_coeff) const;
281
282 std::optional<GroupTarget> getWellGroupTargetProducer(const std::string& name,
283 const std::string& parent,
284 const Group& group,
285 const Scalar* rates,
286 const Scalar efficiency_factor,
287 const std::vector<Scalar>& resv_coeff,
288 std::optional<Group::ProductionCMode> target_cmode = std::nullopt) const;
289
290 GuideRate::RateVector getWellRateVector(const std::string& name) const;
291
292 std::vector<std::string> groupChainTopBot(const std::string& bottom, const std::string& top) const;
293
296 int groupControlledWells(const std::string& group_name,
297 const std::string& always_included_child,
298 const bool is_production_group,
299 const Phase injection_phase) const;
300
302 {
303 return *this->group_state_;
304 }
305
306 const GuideRate& guideRate() const
307 {
308 return this->guide_rate_;
309 }
310
320 bool isMasterGroupEligibleForGuideRate(const std::string& group_name) const;
321
322 bool isRank0() const
323 {
324 return this->well_state_->isRank0();
325 }
326
327 bool isReservoirCouplingMaster() const { return rescoup_.isMaster(); }
328
329 bool isReservoirCouplingMasterGroup(const Group& group) const { return rescoup_.isMasterGroup(group.name()); }
330
331 bool isReservoirCouplingSlave() const { return rescoup_.isSlave(); }
332
333 bool isReservoirCouplingSlaveGroup(const Group& group) const { return rescoup_.isSlaveGroup(group.name()); }
334
335 constexpr int numPhases() const {
336 return this->wellState().numPhases();
337 }
338
339 int phaseToActivePhaseIdx(const Phase phase) const;
340
342 return this->phase_usage_info_;
343 }
344
346 {
347 return GroupStateGuard(*this, group_state);
348 }
349
363 ScopedLoggerGuard pushLogger(bool do_mpi_gather = true) const
364 {
365 return ScopedLoggerGuard(*this, do_mpi_gather);
366 }
367
369 {
370 return WellStateGuard(*this, well_state);
371 }
372
373 int reportStepIdx() const
374 {
375 return report_step_;
376 }
377
378 const Schedule& schedule() const
379 {
380 return this->schedule_;
381 }
382
384 return rescoup_;
385 }
386
388 return rescoup_;
389 }
390
392 return rescoup_.master();
393 }
394
396 return rescoup_.master();
397 }
398
400 return rescoup_.slave();
401 }
403 return rescoup_.slave();
404 }
405
406#ifdef RESERVOIR_COUPLING_ENABLED
407 void setReservoirCouplingMaster(ReservoirCouplingMaster<Scalar>* master) {
408 rescoup_.setMaster(master);
409 }
410 void setReservoirCouplingSlave(ReservoirCouplingSlave<Scalar>* slave) {
411 rescoup_.setSlave(slave);
412 }
413#endif
414
415 void setCmodeGroup(const Group& group);
416
417 template <class AverageRegionalPressureType>
418 void
420 const FieldPropsManager& fp,
421 std::map<std::string, std::unique_ptr<AverageRegionalPressureType>>&
422 regional_average_pressure_calculator) const;
423
424 void setReportStep(int report_step)
425 {
426 report_step_ = report_step;
427 }
428
429 const SummaryState& summaryState() const
430 {
431 return this->summary_state_;
432 }
433
434 Scalar sumSolventRates(const Group& group, const bool is_injector) const;
435
436 Scalar sumWellResRates(const Group& group, const int phase_pos, const bool injector) const;
437
438 Scalar sumWellSurfaceRates(const Group& group, const int phase_pos, const bool injector) const;
439
440 Scalar sumWellPhaseRates(bool res_rates,
441 const Group& group,
442 const int phase_pos,
443 const bool injector,
444 const bool network = false) const;
445
446 bool terminalOutput() const
447 {
448 return this->terminal_output_;
449 }
450
451 template <class RegionalValues>
452 void updateGpMaintTargetForGroups(const Group& group,
453 const RegionalValues& regional_values,
454 const double dt);
455
458 int updateGroupControlledWells(const bool is_production_group,
459 const Phase injection_phase);
460
461 void updateGroupProductionRates(const Group& group);
462
463 void updatePreviousGroupProductionRates(const Group& group);
464
465 void updateGroupTargetReduction(const Group& group,
466 const bool is_injector);
467
469
482
483 void updateREINForGroups(const Group& group, bool sum_rank);
484
485 void updateReservoirRatesInjectionGroups(const Group& group);
486
487#ifdef RESERVOIR_COUPLING_ENABLED
497 void updateSlaveGroupCmodesFromMaster();
498#endif
499
501
502 void updateSurfaceRatesInjectionGroups(const Group& group);
503
504 void updateVREPForGroups(const Group& group);
505
506 void updateWellRates(const Group& group,
507 const WellState<Scalar, IndexTraits>& well_state_nupcol,
508 WellState<Scalar, IndexTraits>& well_state) const;
509
511 {
512 return *this->well_state_;
513 }
514
515 void updateWellRatesFromGroupTargetScale(const Scalar scale,
516 const Group& group,
517 bool is_injector,
518 WellState<Scalar, IndexTraits>& well_state) const;
519
521 std::pair<std::optional<std::string>, Scalar>
522 worstOffendingWell(const Group& group,
523 const Group::ProductionCMode& offended_control) const;
524
525private:
526
548 template<typename ReductionLambda, typename FractionLambda>
549 Scalar applyReductionsAndFractions_(const std::vector<std::string>& chain,
550 Scalar orig_target,
551 Scalar current_rate_available,
552 std::size_t local_reduction_level,
553 bool is_production_group,
554 Phase injection_phase,
555 ReductionLambda&& local_reduction_lambda,
556 FractionLambda&& local_fraction_lambda,
557 bool do_addback) const;
558
559 std::pair<Group::ProductionCMode, Scalar>
560 checkProductionRateConstraint_(const Group& group,
561 Group::ProductionCMode cmode,
562 Group::ProductionCMode currentControl,
563 Scalar target,
564 Scalar current_rate) const;
565
573 std::unordered_set<std::string> collectTargetedProductionGroups_() const;
574
585 Scalar computeAddbackEfficiency_(const std::vector<std::string>& chain,
586 std::size_t local_reduction_level) const;
587
588 std::string controlGroup_(const Group& group) const;
589
590 GuideRate::RateVector getGuideRateVector_(const std::vector<Scalar>& rates) const;
591
592 Scalar getInjectionGroupTargetForMode_(const Group& group,
593 const Phase& injection_phase,
594 const std::vector<Scalar>& resv_coeff,
595 const Group::InjectionCMode cmode) const;
596
607 std::size_t getLocalReductionLevel_(const std::vector<std::string>& chain,
608 bool is_production_group,
609 Phase injection_phase) const;
610
611 Scalar getProductionConstraintTarget_(const Group& group,
612 Group::ProductionCMode cmode,
613 const Group::ProductionControls& controls) const;
614
615 Scalar getProductionGroupTargetForMode_(const Group& group, const Group::ProductionCMode cmode) const;
616
617 Scalar getSatelliteRate_(const Group& group,
618 const int phase_pos,
619 const bool res_rates,
620 const bool is_injector) const;
621
625 bool isAutoChokeGroupUnderperforming_(const Group& group) const;
626
628 bool isInGroupChainTopBot_(const std::string& bottom, const std::string& top) const;
629
630 bool isSatelliteGroup_(const Group& group) const;
631
632 Scalar satelliteInjectionRate_(const ScheduleState& sched,
633 const Group& group,
634 const int phase_pos,
635 bool res_rates) const;
636
637 Scalar satelliteProductionRate_(const ScheduleState& sched,
638 const Group& group,
639 const GSatProd::GSatProdGroupProp::Rate rate_comp,
640 bool res_rates) const;
641
642 std::optional<GSatProd::GSatProdGroupProp::Rate> selectRateComponent_(const int phase_pos) const;
643
654 Scalar subtractOtherPhaseResvInjection_(
655 Phase injection_phase,
656 Scalar base_reservoir_rate,
657 const std::vector<Scalar>& group_injection_reservoir_rates) const;
658
659 Scalar sumProductionRateForControlMode_(const Group& group, Group::ProductionCMode cmode) const;
660
661 int updateGroupControlledWellsRecursive_(const std::string& group_name,
662 const bool is_production_group,
663 const Phase injection_phase);
664
665 void updateGroupTargetReductionRecursive_(const Group& group,
666 const bool is_injector,
667 std::vector<Scalar>& group_target_reduction);
668
669 // Validate that GCONINJE top-up phases (RESV/VREP) are consistent
670 // across the group hierarchy. Called once from setCmodeGroup() for FIELD.
671 void validateInjectionTopupPhases_(const Group& group) const;
672
673 // Recursive helper for validateInjectionTopupPhases_()
674 void validateInjectionTopupPhasesRecursive_(const Group& group,
675 std::optional<Phase> inherited_topup_phase,
676 const std::string& inherited_topup_group) const;
677
678 // --- Reservoir coupling private methods ---
679#ifdef RESERVOIR_COUPLING_ENABLED
684 ReservoirCoupling::Phase activePhaseIdxToRescoupPhase_(int phase_pos) const;
685
694 std::unordered_set<std::string> collectMasterGroupHierarchy_() const;
695
714 int getMasterGroupEffectiveGCW_(const std::string& group_name,
715 bool is_production_group) const;
716
730 Scalar getEffectiveProductionLimit_(const std::string& gname,
731 Group::ProductionCMode rate_type,
732 Scalar slave_local_target) const;
733
734 ReservoirCoupling::GrupSlav::FilterFlag getInjectionFilterFlag_(const std::string& group_name,
735 const Phase injection_phase) const;
736
737 ReservoirCoupling::GrupSlav::FilterFlag getProductionFilterFlag_(const std::string& group_name,
738 const Group::ProductionCMode cmode) const;
739
740 Scalar getReservoirCouplingMasterGroupRate_(const Group& group,
741 const int phase_pos,
742 const ReservoirCoupling::RateKind kind) const;
743#endif // RESERVOIR_COUPLING_ENABLED
744
745 const WellState<Scalar, IndexTraits>* well_state_ {nullptr};
746 GroupState<Scalar>* group_state_ {nullptr};
747 const Schedule& schedule_;
748 const SummaryState& summary_state_;
749 const GuideRate& guide_rate_;
750 // NOTE: The deferred logger does not change the object "meaningful" state, so it should be ok to
751 // make it mutable and store a pointer to it here.
752 mutable DeferredLogger* deferred_logger_ {nullptr};
753 // NOTE: The phase usage info seems to be read-only throughout the simulation, so it should be safe
754 // to store a reference to it here.
755 const PhaseUsageInfo<IndexTraits>& phase_usage_info_;
756 const Parallel::Communication& comm_;
757 bool terminal_output_ {false};
758 int report_step_ {0};
759 ReservoirCoupling::Proxy<Scalar> rescoup_{};
760};
761
762// -----------------------------------------------------------------------------
763// Template member function implementations
764// -----------------------------------------------------------------------------
765
766// NOTE: This template member function is defined here in the header because the
767// AverageRegionalPressureType type parameter depends on derived class types that are
768// not available when GroupStateHelper.cpp is compiled. Template functions
769// must be visible at their instantiation point.
770// See GroupStateHelper.cpp for detailed rationale.
771template <typename Scalar, typename IndexTraits>
772template <class AverageRegionalPressureType>
773void
775 const Group& group,
776 const FieldPropsManager& fp,
777 std::map<std::string, std::unique_ptr<AverageRegionalPressureType>>& regional_average_pressure_calculator)
778 const
779{
780 for (const std::string& groupName : group.groups()) {
781 this->setRegionAveragePressureCalculator(this->schedule_.getGroup(groupName, this->report_step_),
782 fp,
783 regional_average_pressure_calculator);
784 }
785 const auto& gpm = group.gpmaint();
786 if (!gpm)
787 return;
788
789 const auto& reg = gpm->region();
790 if (!reg)
791 return;
792
793 if (regional_average_pressure_calculator.count(reg->first) == 0) {
794 const std::string name = (reg->first.rfind("FIP", 0) == 0) ? reg->first : "FIP" + reg->first;
795 const auto& fipnum = fp.get_int(name);
796 regional_average_pressure_calculator[reg->first]
797 = std::make_unique<AverageRegionalPressureType>(fipnum);
798 }
799}
800
801// NOTE: This template member function is defined in the header because the
802// RegionalValues type parameter depends on derived class types that are
803// not available when GroupStateHelper.cpp is compiled. Template functions
804// must be visible at their instantiation point.
805// See GroupStateHelper.cpp for detailed rationale.
806template <typename Scalar, typename IndexTraits>
807template <class RegionalValues>
808void
810 const RegionalValues& regional_values,
811 const double dt)
812{
813 OPM_TIMEFUNCTION();
814 for (const std::string& group_name : group.groups()) {
815 const Group& group_tmp = this->schedule_.getGroup(group_name, this->report_step_);
816 this->updateGpMaintTargetForGroups(group_tmp, regional_values, dt);
817 }
818 const auto& gpm = group.gpmaint();
819 if (!gpm)
820 return;
821
822 const auto& region = gpm->region();
823 if (!region)
824 return;
825 if (this->isReservoirCouplingMasterGroup(group)) {
826 // GPMAINT is not supported for reservoir coupling master groups since master groups do not have
827 // subordinate wells in the master reservoir, so the slaves cannot influence the master reservoir's
828 // average pressure.
829 // Even specifying GPMAINT on a group superior to the master group might not make sense, since if the
830 // superior target is distributed down to the master group with guide rate fractions, adjusting
831 // the master group's target (that is sent to the slave) could only indirectly influence the master
832 // reservoir's average pressure by affecting the guide rate fractions distributed to actual wells
833 // in the master reservoir.
835 std::runtime_error,
836 "GPMAINT is not supported for reservoir coupling master groups.",
837 this->deferredLogger()
838 );
839 return;
840 }
841 else if (this->isReservoirCouplingSlaveGroup(group)) {
842 // GPMAINT is not supported for reservoir coupling slave groups since their targets will be overridden
843 // by the corresponding master group's target anyway.
845 std::runtime_error,
846 "GPMAINT is not supported for reservoir coupling slave groups.",
847 this->deferredLogger()
848 );
849 return;
850 }
851 const auto [name, number] = *region;
852 const Scalar error = gpm->pressure_target() - regional_values.at(name)->pressure(number);
853 Scalar current_rate = 0.0;
854 const auto& pu = this->phase_usage_info_;
855 bool injection = true;
856 Scalar sign = 1.0;
857 switch (gpm->flow_target()) {
858 case GPMaint::FlowTarget::RESV_PROD: {
859 current_rate = -this->groupState().injection_vrep_rate(group.name());
860 injection = false;
861 sign = -1.0;
862 break;
863 }
864 case GPMaint::FlowTarget::RESV_OINJ: {
865 if (pu.phaseIsActive(IndexTraits::oilPhaseIdx)) {
866 const auto io = pu.canonicalToActivePhaseIdx(IndexTraits::oilPhaseIdx);
867 current_rate = this->groupState().injection_reservoir_rates(group.name())[io];
868 }
869 break;
870 }
871 case GPMaint::FlowTarget::RESV_WINJ: {
872 if (pu.phaseIsActive(IndexTraits::waterPhaseIdx)) {
873 const auto iw = pu.canonicalToActivePhaseIdx(IndexTraits::waterPhaseIdx);
874 current_rate = this->groupState().injection_reservoir_rates(group.name())[iw];
875 }
876 break;
877 }
878 case GPMaint::FlowTarget::RESV_GINJ: {
879 if (pu.phaseIsActive(IndexTraits::gasPhaseIdx)) {
880 const auto ig = pu.canonicalToActivePhaseIdx(IndexTraits::gasPhaseIdx);
881 current_rate = this->groupState().injection_reservoir_rates(group.name())[ig];
882 }
883 break;
884 }
885 case GPMaint::FlowTarget::SURF_OINJ: {
886 if (pu.phaseIsActive(IndexTraits::oilPhaseIdx)) {
887 const auto io = pu.canonicalToActivePhaseIdx(IndexTraits::oilPhaseIdx);
888 current_rate = this->groupState().injection_surface_rates(group.name())[io];
889 }
890 break;
891 }
892 case GPMaint::FlowTarget::SURF_WINJ: {
893 if (pu.phaseIsActive(IndexTraits::waterPhaseIdx)) {
894 const auto iw = pu.canonicalToActivePhaseIdx(IndexTraits::waterPhaseIdx);
895 current_rate = this->groupState().injection_surface_rates(group.name())[iw];
896 }
897 break;
898 }
899 case GPMaint::FlowTarget::SURF_GINJ: {
900 if (pu.phaseIsActive(IndexTraits::gasPhaseIdx)) {
901 const auto ig = pu.canonicalToActivePhaseIdx(IndexTraits::gasPhaseIdx);
902 current_rate = this->groupState().injection_surface_rates(group.name())[ig];
903 }
904 break;
905 }
906 default:
907 throw std::invalid_argument("Invalid Flow target type in GPMAINT");
908 }
909 auto& gpmaint_state = this->groupState().gpmaint(group.name());
910 // we only activate gpmaint if pressure is lower than the target regional pressure for injectors
911 // (i.e. error > 0) and higher for producers.
912 bool activate = (injection && error > 0) || (!injection && error < 0);
913 Scalar rate = 0.0;
914 if (activate) {
915 rate = gpm->rate(gpmaint_state, current_rate, error, dt);
916 } else {
917 gpm->resetState(gpmaint_state);
918 }
919 this->groupState().update_gpmaint_target(group.name(), std::max(Scalar {0.0}, sign * rate));
920}
921
922} // namespace Opm
923
924#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
GuideRateModel::Target getProductionGuideTargetModeFromControlMode(const Group::ProductionCMode cmode) 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:774
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
ReservoirCouplingMaster< Scalar > & reservoirCouplingMaster()
Definition: GroupStateHelper.hpp:391
GroupState< Scalar > & groupState() const
Definition: GroupStateHelper.hpp:301
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:341
const ReservoirCouplingMaster< Scalar > & reservoirCouplingMaster() const
Definition: GroupStateHelper.hpp:395
const SummaryState & summaryState() const
Definition: GroupStateHelper.hpp:429
const GuideRate & guideRate() const
Definition: GroupStateHelper.hpp:306
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:383
void updateReservoirRatesInjectionGroups(const Group &group)
void updatePreviousGroupProductionRates(const Group &group)
bool isReservoirCouplingSlaveGroup(const Group &group) const
Definition: GroupStateHelper.hpp:333
const WellState< Scalar, IndexTraits > & wellState() const
Definition: GroupStateHelper.hpp:510
constexpr int numPhases() const
Definition: GroupStateHelper.hpp:335
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:363
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:322
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:378
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:446
typename SingleWellState< Scalar, IndexTraits >::GroupTarget GroupTarget
Definition: GroupStateHelper.hpp:196
bool isMasterGroupEligibleForGuideRate(const std::string &group_name) const
Whether a reservoir-coupling master group is eligible to participate in its parent's guide-rate distr...
void updateGpMaintTargetForGroups(const Group &group, const RegionalValues &regional_values, const double dt)
Definition: GroupStateHelper.hpp:809
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
void updateNetworkLeafNodeRates()
int reportStepIdx() const
Definition: GroupStateHelper.hpp:373
void setReportStep(int report_step)
Definition: GroupStateHelper.hpp:424
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:329
const ReservoirCoupling::Proxy< Scalar > & rescoup() const
Definition: GroupStateHelper.hpp:387
bool isReservoirCouplingMaster() const
Definition: GroupStateHelper.hpp:327
const ReservoirCouplingSlave< Scalar > & reservoirCouplingSlave() const
Definition: GroupStateHelper.hpp:402
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:368
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:345
ReservoirCouplingSlave< Scalar > & reservoirCouplingSlave()
Definition: GroupStateHelper.hpp:399
Scalar sumWellPhaseRates(bool res_rates, const Group &group, const int phase_pos, const bool injector, const bool network=false) 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, std::optional< Group::ProductionCMode > target_cmode=std::nullopt) const
void updateREINForGroups(const Group &group, bool sum_rank)
bool isReservoirCouplingSlave() const
Definition: GroupStateHelper.hpp:331
GuideRateModel::Target getProductionGuideTargetMode(const Group &group) const
Get the guide rate target mode for a production group.
void updateNONEProductionGroups()
Set production control to NONE for groups not targeting any well.
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:171
Phase
Phase indices for reservoir coupling, we currently only support black-oil phases (oil,...
Definition: ReservoirCoupling.hpp:159
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:124