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
311 bool isRank0() const
312 {
313 return this->well_state_->isRank0();
314 }
315
316 bool isReservoirCouplingMaster() const { return rescoup_.isMaster(); }
317
318 bool isReservoirCouplingMasterGroup(const Group& group) const { return rescoup_.isMasterGroup(group.name()); }
319
320 bool isReservoirCouplingSlave() const { return rescoup_.isSlave(); }
321
322 bool isReservoirCouplingSlaveGroup(const Group& group) const { return rescoup_.isSlaveGroup(group.name()); }
323
324 constexpr int numPhases() const {
325 return this->wellState().numPhases();
326 }
327
328 int phaseToActivePhaseIdx(const Phase phase) const;
329
331 return this->phase_usage_info_;
332 }
333
335 {
336 return GroupStateGuard(*this, group_state);
337 }
338
352 ScopedLoggerGuard pushLogger(bool do_mpi_gather = true) const
353 {
354 return ScopedLoggerGuard(*this, do_mpi_gather);
355 }
356
358 {
359 return WellStateGuard(*this, well_state);
360 }
361
362 int reportStepIdx() const
363 {
364 return report_step_;
365 }
366
367 const Schedule& schedule() const
368 {
369 return this->schedule_;
370 }
371
373 return rescoup_;
374 }
375
377 return rescoup_;
378 }
379
381 return rescoup_.master();
382 }
383
385 return rescoup_.master();
386 }
387
389 return rescoup_.slave();
390 }
392 return rescoup_.slave();
393 }
394
395#ifdef RESERVOIR_COUPLING_ENABLED
396 void setReservoirCouplingMaster(ReservoirCouplingMaster<Scalar>* master) {
397 rescoup_.setMaster(master);
398 }
399 void setReservoirCouplingSlave(ReservoirCouplingSlave<Scalar>* slave) {
400 rescoup_.setSlave(slave);
401 }
402#endif
403
404 void setCmodeGroup(const Group& group);
405
406 template <class AverageRegionalPressureType>
407 void
409 const FieldPropsManager& fp,
410 std::map<std::string, std::unique_ptr<AverageRegionalPressureType>>&
411 regional_average_pressure_calculator) const;
412
413 void setReportStep(int report_step)
414 {
415 report_step_ = report_step;
416 }
417
418 const SummaryState& summaryState() const
419 {
420 return this->summary_state_;
421 }
422
423 Scalar sumSolventRates(const Group& group, const bool is_injector) const;
424
425 Scalar sumWellResRates(const Group& group, const int phase_pos, const bool injector) const;
426
427 Scalar sumWellSurfaceRates(const Group& group, const int phase_pos, const bool injector) const;
428
429 Scalar sumWellPhaseRates(bool res_rates,
430 const Group& group,
431 const int phase_pos,
432 const bool injector,
433 const bool network = false) const;
434
435 bool terminalOutput() const
436 {
437 return this->terminal_output_;
438 }
439
440 template <class RegionalValues>
441 void updateGpMaintTargetForGroups(const Group& group,
442 const RegionalValues& regional_values,
443 const double dt);
444
447 int updateGroupControlledWells(const bool is_production_group,
448 const Phase injection_phase);
449
450 void updateGroupProductionRates(const Group& group);
451
452 void updatePreviousGroupProductionRates(const Group& group);
453
454 void updateGroupTargetReduction(const Group& group,
455 const bool is_injector);
456
458
471
472 void updateREINForGroups(const Group& group, bool sum_rank);
473
474 void updateReservoirRatesInjectionGroups(const Group& group);
475
476#ifdef RESERVOIR_COUPLING_ENABLED
486 void updateSlaveGroupCmodesFromMaster();
487#endif
488
490
491 void updateSurfaceRatesInjectionGroups(const Group& group);
492
493 void updateVREPForGroups(const Group& group);
494
495 void updateWellRates(const Group& group,
496 const WellState<Scalar, IndexTraits>& well_state_nupcol,
497 WellState<Scalar, IndexTraits>& well_state) const;
498
500 {
501 return *this->well_state_;
502 }
503
504 void updateWellRatesFromGroupTargetScale(const Scalar scale,
505 const Group& group,
506 bool is_injector,
507 WellState<Scalar, IndexTraits>& well_state) const;
508
510 std::pair<std::optional<std::string>, Scalar>
511 worstOffendingWell(const Group& group,
512 const Group::ProductionCMode& offended_control) const;
513
514private:
515
537 template<typename ReductionLambda, typename FractionLambda>
538 Scalar applyReductionsAndFractions_(const std::vector<std::string>& chain,
539 Scalar orig_target,
540 Scalar current_rate_available,
541 std::size_t local_reduction_level,
542 bool is_production_group,
543 Phase injection_phase,
544 ReductionLambda&& local_reduction_lambda,
545 FractionLambda&& local_fraction_lambda,
546 bool do_addback) const;
547
548 std::pair<Group::ProductionCMode, Scalar>
549 checkProductionRateConstraint_(const Group& group,
550 Group::ProductionCMode cmode,
551 Group::ProductionCMode currentControl,
552 Scalar target,
553 Scalar current_rate) const;
554
562 std::unordered_set<std::string> collectTargetedProductionGroups_() const;
563
574 Scalar computeAddbackEfficiency_(const std::vector<std::string>& chain,
575 std::size_t local_reduction_level) const;
576
577 std::string controlGroup_(const Group& group) const;
578
579 GuideRate::RateVector getGuideRateVector_(const std::vector<Scalar>& rates) const;
580
581 Scalar getInjectionGroupTargetForMode_(const Group& group,
582 const Phase& injection_phase,
583 const std::vector<Scalar>& resv_coeff,
584 const Group::InjectionCMode cmode) const;
585
596 std::size_t getLocalReductionLevel_(const std::vector<std::string>& chain,
597 bool is_production_group,
598 Phase injection_phase) const;
599
600 Scalar getProductionConstraintTarget_(const Group& group,
601 Group::ProductionCMode cmode,
602 const Group::ProductionControls& controls) const;
603
604 Scalar getProductionGroupTargetForMode_(const Group& group, const Group::ProductionCMode cmode) const;
605
606 Scalar getSatelliteRate_(const Group& group,
607 const int phase_pos,
608 const bool res_rates,
609 const bool is_injector) const;
610
614 bool isAutoChokeGroupUnderperforming_(const Group& group) const;
615
617 bool isInGroupChainTopBot_(const std::string& bottom, const std::string& top) const;
618
619 bool isSatelliteGroup_(const Group& group) const;
620
621 Scalar satelliteInjectionRate_(const ScheduleState& sched,
622 const Group& group,
623 const int phase_pos,
624 bool res_rates) const;
625
626 Scalar satelliteProductionRate_(const ScheduleState& sched,
627 const Group& group,
628 const GSatProd::GSatProdGroupProp::Rate rate_comp,
629 bool res_rates) const;
630
631 std::optional<GSatProd::GSatProdGroupProp::Rate> selectRateComponent_(const int phase_pos) const;
632
643 Scalar subtractOtherPhaseResvInjection_(
644 Phase injection_phase,
645 Scalar base_reservoir_rate,
646 const std::vector<Scalar>& group_injection_reservoir_rates) const;
647
648 Scalar sumProductionRateForControlMode_(const Group& group, Group::ProductionCMode cmode) const;
649
650 int updateGroupControlledWellsRecursive_(const std::string& group_name,
651 const bool is_production_group,
652 const Phase injection_phase);
653
654 void updateGroupTargetReductionRecursive_(const Group& group,
655 const bool is_injector,
656 std::vector<Scalar>& group_target_reduction);
657
658 // Validate that GCONINJE top-up phases (RESV/VREP) are consistent
659 // across the group hierarchy. Called once from setCmodeGroup() for FIELD.
660 void validateInjectionTopupPhases_(const Group& group) const;
661
662 // Recursive helper for validateInjectionTopupPhases_()
663 void validateInjectionTopupPhasesRecursive_(const Group& group,
664 std::optional<Phase> inherited_topup_phase,
665 const std::string& inherited_topup_group) const;
666
667 // --- Reservoir coupling private methods ---
668#ifdef RESERVOIR_COUPLING_ENABLED
673 ReservoirCoupling::Phase activePhaseIdxToRescoupPhase_(int phase_pos) const;
674
683 std::unordered_set<std::string> collectMasterGroupHierarchy_() const;
684
698 Scalar getEffectiveProductionLimit_(const std::string& gname,
699 Group::ProductionCMode rate_type,
700 Scalar slave_local_target) const;
701
702 ReservoirCoupling::GrupSlav::FilterFlag getInjectionFilterFlag_(const std::string& group_name,
703 const Phase injection_phase) const;
704
705 ReservoirCoupling::GrupSlav::FilterFlag getProductionFilterFlag_(const std::string& group_name,
706 const Group::ProductionCMode cmode) const;
707
708 Scalar getReservoirCouplingMasterGroupRate_(const Group& group,
709 const int phase_pos,
710 const ReservoirCoupling::RateKind kind) const;
711#endif // RESERVOIR_COUPLING_ENABLED
712
713 const WellState<Scalar, IndexTraits>* well_state_ {nullptr};
714 GroupState<Scalar>* group_state_ {nullptr};
715 const Schedule& schedule_;
716 const SummaryState& summary_state_;
717 const GuideRate& guide_rate_;
718 // NOTE: The deferred logger does not change the object "meaningful" state, so it should be ok to
719 // make it mutable and store a pointer to it here.
720 mutable DeferredLogger* deferred_logger_ {nullptr};
721 // NOTE: The phase usage info seems to be read-only throughout the simulation, so it should be safe
722 // to store a reference to it here.
723 const PhaseUsageInfo<IndexTraits>& phase_usage_info_;
724 const Parallel::Communication& comm_;
725 bool terminal_output_ {false};
726 int report_step_ {0};
727 ReservoirCoupling::Proxy<Scalar> rescoup_{};
728};
729
730// -----------------------------------------------------------------------------
731// Template member function implementations
732// -----------------------------------------------------------------------------
733
734// NOTE: This template member function is defined here in the header because the
735// AverageRegionalPressureType type parameter depends on derived class types that are
736// not available when GroupStateHelper.cpp is compiled. Template functions
737// must be visible at their instantiation point.
738// See GroupStateHelper.cpp for detailed rationale.
739template <typename Scalar, typename IndexTraits>
740template <class AverageRegionalPressureType>
741void
743 const Group& group,
744 const FieldPropsManager& fp,
745 std::map<std::string, std::unique_ptr<AverageRegionalPressureType>>& regional_average_pressure_calculator)
746 const
747{
748 for (const std::string& groupName : group.groups()) {
749 this->setRegionAveragePressureCalculator(this->schedule_.getGroup(groupName, this->report_step_),
750 fp,
751 regional_average_pressure_calculator);
752 }
753 const auto& gpm = group.gpmaint();
754 if (!gpm)
755 return;
756
757 const auto& reg = gpm->region();
758 if (!reg)
759 return;
760
761 if (regional_average_pressure_calculator.count(reg->first) == 0) {
762 const std::string name = (reg->first.rfind("FIP", 0) == 0) ? reg->first : "FIP" + reg->first;
763 const auto& fipnum = fp.get_int(name);
764 regional_average_pressure_calculator[reg->first]
765 = std::make_unique<AverageRegionalPressureType>(fipnum);
766 }
767}
768
769// NOTE: This template member function is defined in the header because the
770// RegionalValues type parameter depends on derived class types that are
771// not available when GroupStateHelper.cpp is compiled. Template functions
772// must be visible at their instantiation point.
773// See GroupStateHelper.cpp for detailed rationale.
774template <typename Scalar, typename IndexTraits>
775template <class RegionalValues>
776void
778 const RegionalValues& regional_values,
779 const double dt)
780{
781 OPM_TIMEFUNCTION();
782 for (const std::string& group_name : group.groups()) {
783 const Group& group_tmp = this->schedule_.getGroup(group_name, this->report_step_);
784 this->updateGpMaintTargetForGroups(group_tmp, regional_values, dt);
785 }
786 const auto& gpm = group.gpmaint();
787 if (!gpm)
788 return;
789
790 const auto& region = gpm->region();
791 if (!region)
792 return;
793 if (this->isReservoirCouplingMasterGroup(group)) {
794 // GPMAINT is not supported for reservoir coupling master groups since master groups do not have
795 // subordinate wells in the master reservoir, so the slaves cannot influence the master reservoir's
796 // average pressure.
797 // Even specifying GPMAINT on a group superior to the master group might not make sense, since if the
798 // superior target is distributed down to the master group with guide rate fractions, adjusting
799 // the master group's target (that is sent to the slave) could only indirectly influence the master
800 // reservoir's average pressure by affecting the guide rate fractions distributed to actual wells
801 // in the master reservoir.
803 std::runtime_error,
804 "GPMAINT is not supported for reservoir coupling master groups.",
805 this->deferredLogger()
806 );
807 return;
808 }
809 else if (this->isReservoirCouplingSlaveGroup(group)) {
810 // GPMAINT is not supported for reservoir coupling slave groups since their targets will be overridden
811 // by the corresponding master group's target anyway.
813 std::runtime_error,
814 "GPMAINT is not supported for reservoir coupling slave groups.",
815 this->deferredLogger()
816 );
817 return;
818 }
819 const auto [name, number] = *region;
820 const Scalar error = gpm->pressure_target() - regional_values.at(name)->pressure(number);
821 Scalar current_rate = 0.0;
822 const auto& pu = this->phase_usage_info_;
823 bool injection = true;
824 Scalar sign = 1.0;
825 switch (gpm->flow_target()) {
826 case GPMaint::FlowTarget::RESV_PROD: {
827 current_rate = -this->groupState().injection_vrep_rate(group.name());
828 injection = false;
829 sign = -1.0;
830 break;
831 }
832 case GPMaint::FlowTarget::RESV_OINJ: {
833 if (pu.phaseIsActive(IndexTraits::oilPhaseIdx)) {
834 const auto io = pu.canonicalToActivePhaseIdx(IndexTraits::oilPhaseIdx);
835 current_rate = this->groupState().injection_reservoir_rates(group.name())[io];
836 }
837 break;
838 }
839 case GPMaint::FlowTarget::RESV_WINJ: {
840 if (pu.phaseIsActive(IndexTraits::waterPhaseIdx)) {
841 const auto iw = pu.canonicalToActivePhaseIdx(IndexTraits::waterPhaseIdx);
842 current_rate = this->groupState().injection_reservoir_rates(group.name())[iw];
843 }
844 break;
845 }
846 case GPMaint::FlowTarget::RESV_GINJ: {
847 if (pu.phaseIsActive(IndexTraits::gasPhaseIdx)) {
848 const auto ig = pu.canonicalToActivePhaseIdx(IndexTraits::gasPhaseIdx);
849 current_rate = this->groupState().injection_reservoir_rates(group.name())[ig];
850 }
851 break;
852 }
853 case GPMaint::FlowTarget::SURF_OINJ: {
854 if (pu.phaseIsActive(IndexTraits::oilPhaseIdx)) {
855 const auto io = pu.canonicalToActivePhaseIdx(IndexTraits::oilPhaseIdx);
856 current_rate = this->groupState().injection_surface_rates(group.name())[io];
857 }
858 break;
859 }
860 case GPMaint::FlowTarget::SURF_WINJ: {
861 if (pu.phaseIsActive(IndexTraits::waterPhaseIdx)) {
862 const auto iw = pu.canonicalToActivePhaseIdx(IndexTraits::waterPhaseIdx);
863 current_rate = this->groupState().injection_surface_rates(group.name())[iw];
864 }
865 break;
866 }
867 case GPMaint::FlowTarget::SURF_GINJ: {
868 if (pu.phaseIsActive(IndexTraits::gasPhaseIdx)) {
869 const auto ig = pu.canonicalToActivePhaseIdx(IndexTraits::gasPhaseIdx);
870 current_rate = this->groupState().injection_surface_rates(group.name())[ig];
871 }
872 break;
873 }
874 default:
875 throw std::invalid_argument("Invalid Flow target type in GPMAINT");
876 }
877 auto& gpmaint_state = this->groupState().gpmaint(group.name());
878 // we only activate gpmaint if pressure is lower than the target regional pressure for injectors
879 // (i.e. error > 0) and higher for producers.
880 bool activate = (injection && error > 0) || (!injection && error < 0);
881 Scalar rate = 0.0;
882 if (activate) {
883 rate = gpm->rate(gpmaint_state, current_rate, error, dt);
884 } else {
885 gpm->resetState(gpmaint_state);
886 }
887 this->groupState().update_gpmaint_target(group.name(), std::max(Scalar {0.0}, sign * rate));
888}
889
890} // namespace Opm
891
892#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:742
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:380
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:330
const ReservoirCouplingMaster< Scalar > & reservoirCouplingMaster() const
Definition: GroupStateHelper.hpp:384
const SummaryState & summaryState() const
Definition: GroupStateHelper.hpp:418
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:372
void updateReservoirRatesInjectionGroups(const Group &group)
void updatePreviousGroupProductionRates(const Group &group)
bool isReservoirCouplingSlaveGroup(const Group &group) const
Definition: GroupStateHelper.hpp:322
const WellState< Scalar, IndexTraits > & wellState() const
Definition: GroupStateHelper.hpp:499
constexpr int numPhases() const
Definition: GroupStateHelper.hpp:324
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:352
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:311
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:367
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:435
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:777
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:362
void setReportStep(int report_step)
Definition: GroupStateHelper.hpp:413
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:318
const ReservoirCoupling::Proxy< Scalar > & rescoup() const
Definition: GroupStateHelper.hpp:376
bool isReservoirCouplingMaster() const
Definition: GroupStateHelper.hpp:316
const ReservoirCouplingSlave< Scalar > & reservoirCouplingSlave() const
Definition: GroupStateHelper.hpp:391
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:357
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:334
ReservoirCouplingSlave< Scalar > & reservoirCouplingSlave()
Definition: GroupStateHelper.hpp:388
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:320
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: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:124