WellGroupHelpers.hpp
Go to the documentation of this file.
1/*
2 Copyright 2019 Norce.
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
20
21#ifndef OPM_WELLGROUPHELPERS_HEADER_INCLUDED
22#define OPM_WELLGROUPHELPERS_HEADER_INCLUDED
23
24#include "opm/material/fluidsystems/PhaseUsageInfo.hpp"
25
26#include <opm/input/eclipse/Schedule/Group/GuideRate.hpp>
27#include <opm/input/eclipse/Schedule/Group/GSatProd.hpp>
28#include <opm/input/eclipse/EclipseState/Grid/FieldPropsManager.hpp>
30#include <opm/input/eclipse/Schedule/ScheduleState.hpp>
31
32#include <map>
33#include <string>
34#include <vector>
35
36namespace Opm {
37
38class DeferredLogger;
39class Group;
40template<class Scalar> class GroupState;
41namespace Network { class ExtNetwork; }
42class Schedule;
43template<class Scalar> class VFPProdProperties;
44template<typename Scalar, typename IndexTraits> class WellState;
45class FieldPropsManager;
46
47namespace Network { class ExtNetwork; }
48
49template<typename Scalar, typename IndexTraits>
51{
52public:
53
55
56 static Scalar sumWellPhaseRates(bool res_rates,
57 const Opm::Group& group,
58 const Opm::Schedule& schedule,
59 const WellStateType& wellState,
60 const SummaryState& summaryState,
61 const int reportStepIdx,
62 const int phasePos,
63 const bool injector,
64 const bool network = false);
65
66 static Scalar satelliteInjectionRate(const ScheduleState& sched,
67 const Group& group,
69 const int phase_pos,
70 bool res_rates);
71
72 static Scalar satelliteProductionRate(const SummaryState& summaryState,
73 const ScheduleState& sched,
74 const Group& group,
75 const GSatProd::GSatProdGroupProp::Rate rateComp,
76 bool res_rates);
77
78 static std::optional<GSatProd::GSatProdGroupProp::Rate>
79 selectRateComponent(const PhaseUsageInfo<IndexTraits>& pu, const int phasePos);
80
81 static void setCmodeGroup(const Group& group,
82 const Schedule& schedule,
83 const SummaryState& summaryState,
84 const int reportStepIdx,
85 GroupState<Scalar>& group_state);
86
87 static void accumulateGroupEfficiencyFactor(const Group& group,
88 const Schedule& schedule,
89 const int reportStepIdx,
90 Scalar& factor);
91
92 static Scalar sumWellSurfaceRates(const Group& group,
93 const Schedule& schedule,
94 const WellStateType& wellState,
95 const int reportStepIdx,
96 const int phasePos,
97 const bool injector,
98 const SummaryState& summaryState);
99
101 static std::pair<std::optional<std::string>, Scalar>
102 worstOffendingWell(const Group& group,
103 const Schedule& schedule,
104 const int reportStepIdx,
105 const Group::ProductionCMode& offendedControl,
106 const Parallel::Communication& comm,
107 const WellStateType& wellState,
108 DeferredLogger& deferred_logger);
109
110 static Scalar sumWellResRates(const Group& group,
111 const Schedule& schedule,
112 const WellStateType& wellState,
113 const int reportStepIdx,
114 const int phasePos,
115 const bool injector,
116 const SummaryState& summaryState);
117
118 static Scalar sumSolventRates(const Group& group,
119 const Schedule& schedule,
120 const WellStateType& wellState,
121 const int reportStepIdx,
122 const bool injector);
123
124 static void updateGroupTargetReduction(const Group& group,
125 const Schedule& schedule,
126 const int reportStepIdx,
127 const bool isInjector,
128 const GuideRate& guide_rate,
129 const WellStateType& wellState,
130 const SummaryState& summaryState,
131 GroupState<Scalar>& group_state,
132 std::vector<Scalar>& groupTargetReduction);
133
134 static void updateVREPForGroups(const Group& group,
135 const Schedule& schedule,
136 const int reportStepIdx,
137 const WellStateType& wellState,
138 GroupState<Scalar>& group_state,
139 const SummaryState& summaryState);
140
141 template <class RegionalValues>
142 static void updateGpMaintTargetForGroups(const Group& group,
143 const Schedule& schedule,
144 const RegionalValues& regional_values,
145 const int reportStepIdx,
146 const double dt,
147 const WellStateType& well_state,
148 GroupState<Scalar>& group_state);
149
150 static void updateReservoirRatesInjectionGroups(const Group& group,
151 const Schedule& schedule,
152 const int reportStepIdx,
153 const WellStateType& wellState,
154 GroupState<Scalar>& group_state,
155 const SummaryState& summaryState);
156
157 static void updateSurfaceRatesInjectionGroups(const Group& group,
158 const Schedule& schedule,
159 const int reportStepIdx,
160 const WellStateType& wellState,
161 GroupState<Scalar>& group_state,
162 const SummaryState& summaryState);
163
164 static void updateWellRates(const Group& group,
165 const Schedule& schedule,
166 const int reportStepIdx,
167 const WellStateType& wellStateNupcol,
168 WellStateType& wellState);
169
170 static void updateGroupProductionRates(const Group& group,
171 const Schedule& schedule,
172 const int reportStepIdx,
173 const WellStateType& wellState,
174 GroupState<Scalar>& group_state,
175 const SummaryState& summaryState);
176
177 static void updateNetworkLeafNodeProductionRates(const Schedule& schedule,
178 const int reportStepIdx,
179 const WellStateType& wellState,
180 GroupState<Scalar>& group_state,
181 const SummaryState& summaryState);
182
183
184 static void updateWellRatesFromGroupTargetScale(const Scalar scale,
185 const Group& group,
186 const Schedule& schedule,
187 const int reportStepIdx,
188 bool isInjector,
189 const GroupState<Scalar>& group_state,
190 WellStateType& wellState);
191
192 static void updateREINForGroups(const Group& group,
193 const Schedule& schedule,
194 const int reportStepIdx,
195 const SummaryState& st,
196 const WellStateType& wellState,
197 GroupState<Scalar>& group_state,
198 bool sum_rank);
199
200
201 static std::map<std::string, Scalar>
202 computeNetworkPressures(const Network::ExtNetwork& network,
203 const WellStateType& well_state,
204 const GroupState<Scalar>& group_state,
205 const VFPProdProperties<Scalar>& vfp_prod_props,
206 const Schedule& schedule,
207 const Parallel::Communication& comm,
208 const int report_time_step);
209
210 static GuideRate::RateVector
212 const std::string& name);
213
214 static GuideRate::RateVector
217 const std::string& group_name);
218
219 static Scalar getGuideRate(const std::string& name,
220 const Schedule& schedule,
221 const WellStateType& wellState,
222 const GroupState<Scalar>& group_state,
223 const int reportStepIdx,
224 const GuideRate* guideRate,
225 const GuideRateModel::Target target);
226
227 static Scalar getGuideRateInj(const std::string& name,
228 const Schedule& schedule,
229 const WellStateType& wellState,
230 const GroupState<Scalar>& group_state,
231 const int reportStepIdx,
232 const GuideRate* guideRate,
233 const GuideRateModel::Target target,
234 const Phase& injectionPhase);
235
238 static int updateGroupControlledWells(const Schedule& schedule,
239 const WellStateType& well_state,
240 GroupState<Scalar>& group_state,
241 const SummaryState& summary_state,
242 const GuideRate* guideRate,
243 const int report_step,
244 const std::string& group_name,
245 const bool is_production_group,
246 const Phase injection_phase);
247
249 static int groupControlledWells(const Schedule& schedule,
250 const WellStateType& well_state,
251 const GroupState<Scalar>& group_state,
252 const int report_step,
253 const std::string& group_name,
254 const std::string& always_included_child,
255 const bool is_production_group,
256 const Phase injection_phase);
257
258 static std::pair<bool, Scalar>
259 checkGroupConstraintsInj(const std::string& name,
260 const std::string& parent,
261 const Group& group,
262 const WellStateType& wellState,
263 const GroupState<Scalar>& group_state,
264 const int reportStepIdx,
265 const GuideRate* guideRate,
266 const Scalar* rates,
267 Phase injectionPhase,
268 const Scalar efficiencyFactor,
269 const Schedule& schedule,
270 const SummaryState& summaryState,
271 const std::vector<Scalar>& resv_coeff,
272 const bool check_guide_rate,
273 DeferredLogger& deferred_logger);
274
275 static Scalar
276 getWellGroupTargetInjector(const std::string& name,
277 const std::string& parent,
278 const Group& group,
279 const WellStateType& wellState,
280 const GroupState<Scalar>& group_state,
281 const int reportStepIdx,
282 const GuideRate* guideRate,
283 const Scalar* rates,
284 Phase injectionPhase,
285 const Scalar efficiencyFactor,
286 const Schedule& schedule,
287 const SummaryState& summaryState,
288 const std::vector<Scalar>& resv_coeff,
289 DeferredLogger& deferred_logger);
290
291 static std::vector<std::string>
292 groupChainTopBot(const std::string& bottom,
293 const std::string& top,
294 const Schedule& schedule,
295 const int report_step);
296
297
298 // check if well/group bottom is a sub well/group of the group top
299 static bool
300 isInGroupChainTopBot(const std::string& bottom,
301 const std::string& top,
302 const Schedule& schedule,
303 const int report_step);
304
305 static std::string
306 control_group(const Group& group,
307 const GroupState<Scalar>& group_state,
308 const int reportStepIdx,
309 const Schedule& schedule);
310
311 static std::pair<bool, Scalar>
312 checkGroupConstraintsProd(const std::string& name,
313 const std::string& parent,
314 const Group& group,
315 const WellStateType& wellState,
316 const GroupState<Scalar>& group_state,
317 const int reportStepIdx,
318 const GuideRate* guideRate,
319 const Scalar* rates,
320 const Scalar efficiencyFactor,
321 const Schedule& schedule,
322 const SummaryState& summaryState,
323 const std::vector<Scalar>& resv_coeff,
324 const bool check_guide_rate,
325 DeferredLogger& deferred_logger);
326 static Scalar
327 getWellGroupTargetProducer(const std::string& name,
328 const std::string& parent,
329 const Group& group,
330 const WellStateType& wellState,
331 const GroupState<Scalar>& group_state,
332 const int reportStepIdx,
333 const GuideRate* guideRate,
334 const Scalar* rates,
335 const Scalar efficiencyFactor,
336 const Schedule& schedule,
337 const SummaryState& summaryState,
338 const std::vector<Scalar>& resv_coeff,
339 DeferredLogger& deferred_logger);
340
341 template <class AverageRegionalPressureType>
342 static void setRegionAveragePressureCalculator(const Group& group,
343 const Schedule& schedule,
344 const int reportStepIdx,
345 const FieldPropsManager& fp,
346 std::map<std::string, std::unique_ptr<AverageRegionalPressureType>>& regionalAveragePressureCalculator);
347};
348
349} // namespace Opm
350
351#endif
Definition: DeferredLogger.hpp:57
Definition: GroupState.hpp:41
Definition: GasLiftGroupInfo.hpp:37
Definition: VFPProdProperties.hpp:38
Definition: WellGroupHelpers.hpp:51
static void updateNetworkLeafNodeProductionRates(const Schedule &schedule, const int reportStepIdx, const WellStateType &wellState, GroupState< Scalar > &group_state, const SummaryState &summaryState)
static void updateGroupProductionRates(const Group &group, const Schedule &schedule, const int reportStepIdx, const WellStateType &wellState, GroupState< Scalar > &group_state, const SummaryState &summaryState)
static void setCmodeGroup(const Group &group, const Schedule &schedule, const SummaryState &summaryState, const int reportStepIdx, GroupState< Scalar > &group_state)
static void updateWellRates(const Group &group, const Schedule &schedule, const int reportStepIdx, const WellStateType &wellStateNupcol, WellStateType &wellState)
static bool isInGroupChainTopBot(const std::string &bottom, const std::string &top, const Schedule &schedule, const int report_step)
static void updateSurfaceRatesInjectionGroups(const Group &group, const Schedule &schedule, const int reportStepIdx, const WellStateType &wellState, GroupState< Scalar > &group_state, const SummaryState &summaryState)
static std::string control_group(const Group &group, const GroupState< Scalar > &group_state, const int reportStepIdx, const Schedule &schedule)
static std::pair< bool, Scalar > checkGroupConstraintsProd(const std::string &name, const std::string &parent, const Group &group, const WellStateType &wellState, const GroupState< Scalar > &group_state, const int reportStepIdx, const GuideRate *guideRate, const Scalar *rates, const Scalar efficiencyFactor, const Schedule &schedule, const SummaryState &summaryState, const std::vector< Scalar > &resv_coeff, const bool check_guide_rate, DeferredLogger &deferred_logger)
static Scalar sumWellResRates(const Group &group, const Schedule &schedule, const WellStateType &wellState, const int reportStepIdx, const int phasePos, const bool injector, const SummaryState &summaryState)
static Scalar getWellGroupTargetProducer(const std::string &name, const std::string &parent, const Group &group, const WellStateType &wellState, const GroupState< Scalar > &group_state, const int reportStepIdx, const GuideRate *guideRate, const Scalar *rates, const Scalar efficiencyFactor, const Schedule &schedule, const SummaryState &summaryState, const std::vector< Scalar > &resv_coeff, DeferredLogger &deferred_logger)
static Scalar sumWellPhaseRates(bool res_rates, const Opm::Group &group, const Opm::Schedule &schedule, const WellStateType &wellState, const SummaryState &summaryState, const int reportStepIdx, const int phasePos, const bool injector, const bool network=false)
static std::pair< bool, Scalar > checkGroupConstraintsInj(const std::string &name, const std::string &parent, const Group &group, const WellStateType &wellState, const GroupState< Scalar > &group_state, const int reportStepIdx, const GuideRate *guideRate, const Scalar *rates, Phase injectionPhase, const Scalar efficiencyFactor, const Schedule &schedule, const SummaryState &summaryState, const std::vector< Scalar > &resv_coeff, const bool check_guide_rate, DeferredLogger &deferred_logger)
static std::map< std::string, Scalar > computeNetworkPressures(const Network::ExtNetwork &network, const WellStateType &well_state, const GroupState< Scalar > &group_state, const VFPProdProperties< Scalar > &vfp_prod_props, const Schedule &schedule, const Parallel::Communication &comm, const int report_time_step)
static void updateVREPForGroups(const Group &group, const Schedule &schedule, const int reportStepIdx, const WellStateType &wellState, GroupState< Scalar > &group_state, const SummaryState &summaryState)
static void updateGpMaintTargetForGroups(const Group &group, const Schedule &schedule, const RegionalValues &regional_values, const int reportStepIdx, const double dt, const WellStateType &well_state, GroupState< Scalar > &group_state)
static Scalar satelliteInjectionRate(const ScheduleState &sched, const Group &group, const PhaseUsageInfo< IndexTraits > &pu, const int phase_pos, bool res_rates)
static void setRegionAveragePressureCalculator(const Group &group, const Schedule &schedule, const int reportStepIdx, const FieldPropsManager &fp, std::map< std::string, std::unique_ptr< AverageRegionalPressureType > > &regionalAveragePressureCalculator)
static void updateREINForGroups(const Group &group, const Schedule &schedule, const int reportStepIdx, const SummaryState &st, const WellStateType &wellState, GroupState< Scalar > &group_state, bool sum_rank)
static Scalar getWellGroupTargetInjector(const std::string &name, const std::string &parent, const Group &group, const WellStateType &wellState, const GroupState< Scalar > &group_state, const int reportStepIdx, const GuideRate *guideRate, const Scalar *rates, Phase injectionPhase, const Scalar efficiencyFactor, const Schedule &schedule, const SummaryState &summaryState, const std::vector< Scalar > &resv_coeff, DeferredLogger &deferred_logger)
static std::pair< std::optional< std::string >, Scalar > worstOffendingWell(const Group &group, const Schedule &schedule, const int reportStepIdx, const Group::ProductionCMode &offendedControl, const Parallel::Communication &comm, const WellStateType &wellState, DeferredLogger &deferred_logger)
Returns the name of the worst offending well and its fraction (i.e. violated_phase / preferred_phase)
static void accumulateGroupEfficiencyFactor(const Group &group, const Schedule &schedule, const int reportStepIdx, Scalar &factor)
static Scalar sumSolventRates(const Group &group, const Schedule &schedule, const WellStateType &wellState, const int reportStepIdx, const bool injector)
static Scalar getGuideRateInj(const std::string &name, const Schedule &schedule, const WellStateType &wellState, const GroupState< Scalar > &group_state, const int reportStepIdx, const GuideRate *guideRate, const GuideRateModel::Target target, const Phase &injectionPhase)
static int updateGroupControlledWells(const Schedule &schedule, const WellStateType &well_state, GroupState< Scalar > &group_state, const SummaryState &summary_state, const GuideRate *guideRate, const int report_step, const std::string &group_name, const bool is_production_group, const Phase injection_phase)
static Scalar sumWellSurfaceRates(const Group &group, const Schedule &schedule, const WellStateType &wellState, const int reportStepIdx, const int phasePos, const bool injector, const SummaryState &summaryState)
static Scalar satelliteProductionRate(const SummaryState &summaryState, const ScheduleState &sched, const Group &group, const GSatProd::GSatProdGroupProp::Rate rateComp, bool res_rates)
static std::vector< std::string > groupChainTopBot(const std::string &bottom, const std::string &top, const Schedule &schedule, const int report_step)
static Scalar getGuideRate(const std::string &name, const Schedule &schedule, const WellStateType &wellState, const GroupState< Scalar > &group_state, const int reportStepIdx, const GuideRate *guideRate, const GuideRateModel::Target target)
static void updateGroupTargetReduction(const Group &group, const Schedule &schedule, const int reportStepIdx, const bool isInjector, const GuideRate &guide_rate, const WellStateType &wellState, const SummaryState &summaryState, GroupState< Scalar > &group_state, std::vector< Scalar > &groupTargetReduction)
static void updateWellRatesFromGroupTargetScale(const Scalar scale, const Group &group, const Schedule &schedule, const int reportStepIdx, bool isInjector, const GroupState< Scalar > &group_state, WellStateType &wellState)
static GuideRate::RateVector getWellRateVector(const WellStateType &well_state, const std::string &name)
static GuideRate::RateVector getProductionGroupRateVector(const GroupState< Scalar > &group_state, const PhaseUsageInfo< IndexTraits > &pu, const std::string &group_name)
static void updateReservoirRatesInjectionGroups(const Group &group, const Schedule &schedule, const int reportStepIdx, const WellStateType &wellState, GroupState< Scalar > &group_state, const SummaryState &summaryState)
static int groupControlledWells(const Schedule &schedule, const WellStateType &well_state, const GroupState< Scalar > &group_state, const int report_step, const std::string &group_name, const std::string &always_included_child, const bool is_production_group, const Phase injection_phase)
returns the number of wells that are actively under group control for a given group with name given b...
static std::optional< GSatProd::GSatProdGroupProp::Rate > selectRateComponent(const PhaseUsageInfo< IndexTraits > &pu, const int phasePos)
Definition: WellState.hpp:66
Dune::Communication< MPIComm > Communication
Definition: ParallelCommunication.hpp:30
Definition: blackoilboundaryratevector.hh:39