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