21 #ifndef OPM_WELLGROUPHELPERS_HEADER_INCLUDED
22 #define OPM_WELLGROUPHELPERS_HEADER_INCLUDED
24 #include <opm/input/eclipse/Schedule/Group/GuideRate.hpp>
36 namespace Network {
class ExtNetwork; }
39 class VFPProdProperties;
42 namespace Network {
class ExtNetwork; }
44 namespace WellGroupHelpers
49 void setCmodeGroup(
const Group& group,
50 const Schedule& schedule,
51 const SummaryState& summaryState,
52 const int reportStepIdx,
54 GroupState& group_state);
56 void accumulateGroupEfficiencyFactor(
const Group& group,
57 const Schedule& schedule,
58 const int reportStepIdx,
61 double sumWellSurfaceRates(
const Group& group,
62 const Schedule& schedule,
63 const WellState& wellState,
64 const int reportStepIdx,
68 double sumWellResRates(
const Group& group,
69 const Schedule& schedule,
70 const WellState& wellState,
71 const int reportStepIdx,
75 double sumSolventRates(
const Group& group,
76 const Schedule& schedule,
77 const WellState& wellState,
78 const int reportStepIdx,
81 void updateGroupTargetReduction(
const Group& group,
82 const Schedule& schedule,
83 const int reportStepIdx,
84 const bool isInjector,
86 const GuideRate& guide_rate,
87 const WellState& wellState,
88 GroupState& group_state,
89 std::vector<double>& groupTargetReduction);
92 void updateGuideRates(
const Group& group,
93 const Schedule& schedule,
94 const SummaryState& summary_state,
98 WellState& well_state,
99 const GroupState& group_state,
101 GuideRate* guide_rate,
102 std::vector<double>& pot,
105 template <
class Comm>
106 void updateGuideRateForProductionGroups(
const Group& group,
107 const Schedule& schedule,
108 const PhaseUsage& pu,
109 const int reportStepIdx,
110 const double& simTime,
111 WellState& wellState,
112 const GroupState& group_state,
114 GuideRate* guideRate,
115 std::vector<double>& pot);
117 template <
class Comm>
118 void updateGuideRatesForWells(
const Schedule& schedule,
119 const PhaseUsage& pu,
120 const int reportStepIdx,
121 const double& simTime,
122 const WellState& wellState,
124 GuideRate* guideRate);
126 void updateGuideRatesForInjectionGroups(
const Group& group,
127 const Schedule& schedule,
128 const SummaryState& summaryState,
130 const int reportStepIdx,
131 const WellState& wellState,
132 const GroupState& group_state,
133 GuideRate* guideRate,
136 void updateVREPForGroups(
const Group& group,
137 const Schedule& schedule,
138 const int reportStepIdx,
139 const WellState& wellState,
140 GroupState& group_state);
142 void updateReservoirRatesInjectionGroups(
const Group& group,
143 const Schedule& schedule,
144 const int reportStepIdx,
145 const WellState& wellState,
146 GroupState& group_state);
148 void updateSurfaceRatesInjectionGroups(
const Group& group,
149 const Schedule& schedule,
150 const int reportStepIdx,
151 const WellState& wellState,
152 GroupState& group_state);
154 void updateWellRates(
const Group& group,
155 const Schedule& schedule,
156 const int reportStepIdx,
157 const WellState& wellStateNupcol,
158 WellState& wellState);
160 void updateGroupProductionRates(
const Group& group,
161 const Schedule& schedule,
162 const int reportStepIdx,
163 const WellState& wellState,
164 GroupState& group_state);
166 void updateWellRatesFromGroupTargetScale(
const double scale,
168 const Schedule& schedule,
169 const int reportStepIdx,
171 const GroupState& group_state,
172 WellState& wellState);
174 void updateREINForGroups(
const Group& group,
175 const Schedule& schedule,
176 const int reportStepIdx,
177 const PhaseUsage& pu,
178 const SummaryState& st,
179 const WellState& wellState,
180 GroupState& group_state);
182 template <
class RegionalValues>
183 void updateGpMaintTargetForGroups(
const Group& group,
184 const Schedule& schedule,
185 const RegionalValues& regional_values,
186 const int reportStepIdx,
188 const WellState& well_state,
189 GroupState& group_state);
191 std::map<std::string, double>
192 computeNetworkPressures(
const Opm::Network::ExtNetwork& network,
193 const WellState& well_state,
194 const GroupState& group_state,
195 const VFPProdProperties& vfp_prod_props,
196 const Schedule& schedule,
197 const int report_time_step);
199 GuideRate::RateVector
200 getWellRateVector(
const WellState& well_state,
const PhaseUsage& pu,
const std::string& name);
202 GuideRate::RateVector
203 getProductionGroupRateVector(
const GroupState& group_state,
const PhaseUsage& pu,
const std::string& group_name);
205 double getGuideRate(
const std::string& name,
206 const Schedule& schedule,
207 const WellState& wellState,
208 const GroupState& group_state,
209 const int reportStepIdx,
210 const GuideRate* guideRate,
211 const GuideRateModel::Target target,
212 const PhaseUsage& pu);
215 double getGuideRateInj(
const std::string& name,
216 const Schedule& schedule,
217 const WellState& wellState,
218 const GroupState& group_state,
219 const int reportStepIdx,
220 const GuideRate* guideRate,
221 const GuideRateModel::Target target,
222 const Phase& injectionPhase,
223 const PhaseUsage& pu);
225 int groupControlledWells(
const Schedule& schedule,
226 const WellState& well_state,
227 const GroupState& group_state,
228 const int report_step,
229 const std::string& group_name,
230 const std::string& always_included_child,
231 const bool is_production_group,
232 const Phase injection_phase);
241 const int report_step,
242 const GuideRate* guide_rate,
243 const GuideRateModel::Target target,
245 const bool is_producer,
246 const Phase injection_phase);
247 double fraction(
const std::string& name,
const std::string& control_group_name,
const bool always_include_this);
248 double localFraction(
const std::string& name,
const std::string& always_included_child);
251 std::string parent(
const std::string& name);
252 double guideRateSum(
const Group& group,
const std::string& always_included_child);
253 double guideRate(
const std::string& name,
const std::string& always_included_child);
254 int groupControlledWells(
const std::string& group_name,
const std::string& always_included_child);
255 GuideRate::RateVector getGroupRateVector(
const std::string& group_name);
256 const Schedule& schedule_;
260 const GuideRate* guide_rate_;
261 GuideRateModel::Target target_;
264 Phase injection_phase_;
268 std::pair<bool, double> checkGroupConstraintsInj(
const std::string& name,
269 const std::string& parent,
273 const int reportStepIdx,
274 const GuideRate* guideRate,
276 Phase injectionPhase,
278 const double efficiencyFactor,
279 const Schedule& schedule,
280 const SummaryState& summaryState,
281 const std::vector<double>& resv_coeff,
289 std::vector<std::string> groupChainTopBot(
const std::string& bottom,
290 const std::string& top,
291 const Schedule& schedule,
292 const int report_step);
297 std::pair<bool, double> checkGroupConstraintsProd(
const std::string& name,
298 const std::string& parent,
302 const int reportStepIdx,
303 const GuideRate* guideRate,
306 const double efficiencyFactor,
307 const Schedule& schedule,
308 const SummaryState& summaryState,
309 const std::vector<double>& resv_coeff,
Definition: DeferredLogger.hpp:57
Definition: GroupState.hpp:34
Definition: WellGroupHelpers.hpp:236
The state of a set of wells, tailored for use by the fully implicit blackoil simulator.
Definition: WellState.hpp:56
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition: BlackoilPhases.hpp:27
Definition: BlackoilPhases.hpp:46