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/EclipseState/Grid/FieldPropsManager.hpp>
27
28#include <map>
29#include <string>
30#include <vector>
31
32namespace Opm {
33
34class DeferredLogger;
35class Group;
36template<class Scalar> class GroupState;
37namespace Network { class ExtNetwork; }
38struct PhaseUsage;
39class Schedule;
41template<class Scalar> class WellState;
42class FieldPropsManager;
43
44namespace Network { class ExtNetwork; }
45
46template<class Scalar>
48{
49public:
50 static void setCmodeGroup(const Group& group,
51 const Schedule& schedule,
52 const SummaryState& summaryState,
53 const int reportStepIdx,
54 GroupState<Scalar>& group_state);
55
56 static void accumulateGroupEfficiencyFactor(const Group& group,
57 const Schedule& schedule,
58 const int reportStepIdx,
59 Scalar& factor);
60
61 static Scalar sumWellSurfaceRates(const Group& group,
62 const Schedule& schedule,
63 const WellState<Scalar>& wellState,
64 const int reportStepIdx,
65 const int phasePos,
66 const bool injector);
67
69 static std::pair<std::optional<std::string>, Scalar>
70 worstOffendingWell(const Group& group,
71 const Schedule& schedule,
72 const int reportStepIdx,
73 const Group::ProductionCMode& offendedControl,
74 const PhaseUsage& pu,
75 const Parallel::Communication& comm,
76 const WellState<Scalar>& wellState,
77 DeferredLogger& deferred_logger);
78
79 static Scalar sumWellResRates(const Group& group,
80 const Schedule& schedule,
81 const WellState<Scalar>& wellState,
82 const int reportStepIdx,
83 const int phasePos,
84 const bool injector);
85
86 static Scalar sumSolventRates(const Group& group,
87 const Schedule& schedule,
88 const WellState<Scalar>& wellState,
89 const int reportStepIdx,
90 const bool injector);
91
92 static void updateGroupTargetReduction(const Group& group,
93 const Schedule& schedule,
94 const int reportStepIdx,
95 const bool isInjector,
96 const PhaseUsage& pu,
97 const GuideRate& guide_rate,
98 const WellState<Scalar>& wellState,
99 GroupState<Scalar>& group_state,
100 std::vector<Scalar>& groupTargetReduction);
101
102 static void updateGuideRates(const Group& group,
103 const Schedule& schedule,
104 const SummaryState& summary_state,
105 const PhaseUsage& pu,
106 int report_step,
107 double sim_time,
108 WellState<Scalar>& well_state,
109 const GroupState<Scalar>& group_state,
110 const Parallel::Communication& comm,
111 GuideRate* guide_rate,
112 std::vector<Scalar>& pot,
113 DeferredLogger& deferred_logger);
114
115 static void updateGuideRateForProductionGroups(const Group& group,
116 const Schedule& schedule,
117 const PhaseUsage& pu,
118 const int reportStepIdx,
119 const double& simTime,
120 WellState<Scalar>& wellState,
121 const GroupState<Scalar>& group_state,
122 const Parallel::Communication& comm,
123 GuideRate* guideRate,
124 std::vector<Scalar>& pot);
125
126 static void updateGuideRatesForWells(const Schedule& schedule,
127 const PhaseUsage& pu,
128 const int reportStepIdx,
129 const double& simTime,
130 const WellState<Scalar>& wellState,
131 const Parallel::Communication& comm,
132 GuideRate* guideRate);
133
134 static void updateGuideRatesForInjectionGroups(const Group& group,
135 const Schedule& schedule,
136 const SummaryState& summaryState,
137 const PhaseUsage& pu,
138 const int reportStepIdx,
139 const WellState<Scalar>& wellState,
140 const GroupState<Scalar>& group_state,
141 GuideRate* guideRate,
142 DeferredLogger& deferred_logger);
143
144 static void updateVREPForGroups(const Group& group,
145 const Schedule& schedule,
146 const int reportStepIdx,
147 const WellState<Scalar>& wellState,
148 GroupState<Scalar>& group_state);
149
150 template <class RegionalValues>
151 static void updateGpMaintTargetForGroups(const Group& group,
152 const Schedule& schedule,
153 const RegionalValues& regional_values,
154 const int reportStepIdx,
155 const double dt,
156 const WellState<Scalar>& well_state,
157 GroupState<Scalar>& group_state);
158
159 static void updateReservoirRatesInjectionGroups(const Group& group,
160 const Schedule& schedule,
161 const int reportStepIdx,
162 const WellState<Scalar>& wellState,
163 GroupState<Scalar>& group_state);
164
165 static void updateSurfaceRatesInjectionGroups(const Group& group,
166 const Schedule& schedule,
167 const int reportStepIdx,
168 const WellState<Scalar>& wellState,
169 GroupState<Scalar>& group_state);
170
171 static void updateWellRates(const Group& group,
172 const Schedule& schedule,
173 const int reportStepIdx,
174 const WellState<Scalar>& wellStateNupcol,
175 WellState<Scalar>& wellState);
176
177 static void updateGroupProductionRates(const Group& group,
178 const Schedule& schedule,
179 const int reportStepIdx,
180 const WellState<Scalar>& wellState,
181 GroupState<Scalar>& group_state);
182
183 static void updateWellRatesFromGroupTargetScale(const Scalar scale,
184 const Group& group,
185 const Schedule& schedule,
186 const int reportStepIdx,
187 bool isInjector,
188 const GroupState<Scalar>& group_state,
189 WellState<Scalar>& wellState);
190
191 static void updateREINForGroups(const Group& group,
192 const Schedule& schedule,
193 const int reportStepIdx,
194 const PhaseUsage& pu,
195 const SummaryState& st,
196 const WellState<Scalar>& 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 WellState<Scalar>& well_state,
204 const GroupState<Scalar>& group_state,
205 const VFPProdProperties& vfp_prod_props,
206 const Schedule& schedule,
207 const int report_time_step);
208
209 static GuideRate::RateVector
211 const PhaseUsage& pu,
212 const std::string& name);
213
214 static GuideRate::RateVector
216 const PhaseUsage& pu,
217 const std::string& group_name);
218
219 static Scalar getGuideRate(const std::string& name,
220 const Schedule& schedule,
221 const WellState<Scalar>& wellState,
222 const GroupState<Scalar>& group_state,
223 const int reportStepIdx,
224 const GuideRate* guideRate,
225 const GuideRateModel::Target target,
226 const PhaseUsage& pu);
227
228 static Scalar getGuideRateInj(const std::string& name,
229 const Schedule& schedule,
230 const WellState<Scalar>& wellState,
231 const GroupState<Scalar>& group_state,
232 const int reportStepIdx,
233 const GuideRate* guideRate,
234 const GuideRateModel::Target target,
235 const Phase& injectionPhase,
236 const PhaseUsage& pu);
237
238 static int groupControlledWells(const Schedule& schedule,
239 const WellState<Scalar>& well_state,
240 const GroupState<Scalar>& group_state,
241 const int report_step,
242 const std::string& group_name,
243 const std::string& always_included_child,
244 const bool is_production_group,
245 const Phase injection_phase);
246
247 static std::pair<bool, Scalar>
248 checkGroupConstraintsInj(const std::string& name,
249 const std::string& parent,
250 const Group& group,
251 const WellState<Scalar>& wellState,
252 const GroupState<Scalar>& group_state,
253 const int reportStepIdx,
254 const GuideRate* guideRate,
255 const Scalar* rates,
256 Phase injectionPhase,
257 const PhaseUsage& pu,
258 const Scalar efficiencyFactor,
259 const Schedule& schedule,
260 const SummaryState& summaryState,
261 const std::vector<Scalar>& resv_coeff,
262 DeferredLogger& deferred_logger);
263
264 static std::vector<std::string>
265 groupChainTopBot(const std::string& bottom,
266 const std::string& top,
267 const Schedule& schedule,
268 const int report_step);
269
270 static std::pair<bool, Scalar>
271 checkGroupConstraintsProd(const std::string& name,
272 const std::string& parent,
273 const Group& group,
274 const WellState<Scalar>& wellState,
275 const GroupState<Scalar>& group_state,
276 const int reportStepIdx,
277 const GuideRate* guideRate,
278 const Scalar* rates,
279 const PhaseUsage& pu,
280 const Scalar efficiencyFactor,
281 const Schedule& schedule,
282 const SummaryState& summaryState,
283 const std::vector<Scalar>& resv_coeff,
284 DeferredLogger& deferred_logger);
285
286 template <class AverageRegionalPressureType>
287 static void setRegionAveragePressureCalculator(const Group& group,
288 const Schedule& schedule,
289 const int reportStepIdx,
290 const FieldPropsManager& fp,
291 const PhaseUsage& pu,
292 std::map<std::string, std::unique_ptr<AverageRegionalPressureType>>& regionalAveragePressureCalculator);
293};
294
295} // namespace Opm
296
297#endif
Definition: DeferredLogger.hpp:57
Definition: GroupState.hpp:35
Definition: VFPProdProperties.hpp:37
Definition: WellGroupHelpers.hpp:48
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, GroupState< Scalar > &group_state, std::vector< Scalar > &groupTargetReduction)
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::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, DeferredLogger &deferred_logger)
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 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, DeferredLogger &deferred_logger)
static void updateSurfaceRatesInjectionGroups(const Group &group, const Schedule &schedule, const int reportStepIdx, const WellState< Scalar > &wellState, GroupState< Scalar > &group_state)
static std::map< std::string, Scalar > computeNetworkPressures(const Network::ExtNetwork &network, const WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, const VFPProdProperties &vfp_prod_props, const Schedule &schedule, const int report_time_step)
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 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 int groupControlledWells(const Schedule &schedule, const WellState< Scalar > &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)
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 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:62
Dune::Communication< MPIComm > Communication
Definition: ParallelCommunication.hpp:30
Definition: BlackoilPhases.hpp:27
Definition: BlackoilPhases.hpp:46