WellGroupHelper.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_WELLGROUPHELPER_HEADER_INCLUDED
20#define OPM_WELLGROUPHELPER_HEADER_INCLUDED
21#if HAVE_MPI
22#define RESERVOIR_COUPLING_ENABLED
23#endif
24#ifdef RESERVOIR_COUPLING_ENABLED
28#endif
29
30#include <opm/common/TimingMacros.hpp>
31#include <opm/input/eclipse/EclipseState/Grid/FieldPropsManager.hpp>
32#include <opm/input/eclipse/Schedule/Group/GPMaint.hpp>
33#include <opm/input/eclipse/Schedule/Group/GSatProd.hpp>
34#include <opm/input/eclipse/Schedule/Group/GuideRate.hpp>
35#include <opm/input/eclipse/Schedule/Schedule.hpp>
36#include <opm/input/eclipse/Schedule/ScheduleState.hpp>
37#include <opm/input/eclipse/Schedule/SummaryState.hpp>
38#include <opm/material/fluidsystems/PhaseUsageInfo.hpp>
43
45
46#include <algorithm>
47#include <map>
48#include <memory>
49#include <optional>
50#include <stdexcept>
51#include <string>
52#include <vector>
53
54namespace Opm
55{
56
57template <typename Scalar, typename IndexTraits>
59{
60public:
61 // RAII guard for temporarily setting wellstate pointer
63 {
64 public:
66 : wgHelper_ {wgHelper}
67 , previous_state_ptr_ {wgHelper_.well_state_}
68 {
69 // Set the new state directly
70 wgHelper_.well_state_ = &well_state;
71 }
72
74 {
75 // Restore the previous state
76 wgHelper_.well_state_ = previous_state_ptr_;
77 }
78
79 // Delete copy and move operations
84
85 private:
86 WellGroupHelper& wgHelper_;
87 const WellState<Scalar, IndexTraits>* previous_state_ptr_;
88 };
89
90 // RAII guard for temporarily setting groupstate pointer
92 {
93 public:
95 : wgHelper_ {wgHelper}
96 , previous_state_ptr_ {wgHelper_.group_state_}
97 {
98 // Set the new state directly
99 wgHelper_.group_state_ = &group_state;
100 }
101
103 {
104 // Restore the previous state
105 wgHelper_.group_state_ = previous_state_ptr_;
106 }
107
108 // Delete copy and move operations
113
114 private:
115 WellGroupHelper& wgHelper_;
116 const GroupState<Scalar>* previous_state_ptr_;
117 };
118
120 GroupState<Scalar>& group_state,
121 const Schedule& schedule,
122 const SummaryState& summary_state,
123 const GuideRate& guide_rate,
124 const PhaseUsageInfo<IndexTraits>& phase_usage_info);
125
126 void accumulateGroupEfficiencyFactor(const Group& group, Scalar& factor) const;
127
128 std::pair<bool, Scalar> checkGroupConstraintsInj(const std::string& name,
129 const std::string& parent,
130 const Group& group,
131 const Scalar* rates,
132 const Phase injection_phase,
133 const Scalar efficiency_factor,
134 const std::vector<Scalar>& resv_coeff,
135 const bool check_guide_rate,
136 DeferredLogger& deferred_logger) const;
137
138 std::pair<bool, Scalar> checkGroupConstraintsProd(const std::string& name,
139 const std::string& parent,
140 const Group& group,
141 const Scalar* rates,
142 const Scalar efficiency_factor,
143 const std::vector<Scalar>& resv_coeff,
144 const bool check_guide_rate,
145 DeferredLogger& deferred_logger) const;
146
147 std::map<std::string, Scalar> computeNetworkPressures(const Network::ExtNetwork& network,
148 const VFPProdProperties<Scalar>& vfp_prod_props,
149 const Parallel::Communication& comm) const;
150
151 Scalar getGuideRate(const std::string& name, const GuideRateModel::Target target) const;
152
153 GuideRate::RateVector getProductionGroupRateVector(const std::string& group_name) const;
154
155 std::optional<Scalar> getWellGroupTargetInjector(const std::string& name,
156 const std::string& parent,
157 const Group& group,
158 const Scalar* rates,
159 const Phase injection_phase,
160 const Scalar efficiency_factor,
161 const std::vector<Scalar>& resv_coeff,
162 DeferredLogger& deferred_logger) const;
163
164 std::optional<Scalar> getWellGroupTargetProducer(const std::string& name,
165 const std::string& parent,
166 const Group& group,
167 const Scalar* rates,
168 const Scalar efficiency_factor,
169 const std::vector<Scalar>& resv_coeff,
170 DeferredLogger& deferred_logger) const;
171
172 GuideRate::RateVector getWellRateVector(const std::string& name) const;
173
174 std::vector<std::string> groupChainTopBot(const std::string& bottom, const std::string& top) const;
175
178 int groupControlledWells(const std::string& group_name,
179 const std::string& always_included_child,
180 const bool is_production_group,
181 const Phase injection_phase) const;
182
184 {
185 return *this->group_state_;
186 }
187
189 {
190 return this->phase_usage_info_;
191 }
192
194 {
195 return GroupStateGuard(*this, group_state);
196 }
197
199 {
200 return WellStateGuard(*this, well_state);
201 }
202
203 int reportStepIdx() const
204 {
205 return report_step_;
206 }
207
208 const Schedule& schedule() const
209 {
210 return this->schedule_;
211 }
212
213 const GuideRate& guideRate() const
214 {
215 return this->guide_rate_;
216 }
217
219 return this->phase_usage_info_;
220 }
221
222#ifdef RESERVOIR_COUPLING_ENABLED
224 {
225 return *this->reservoir_coupling_master_;
226 }
227
229 {
230 return *this->reservoir_coupling_slave_;
231 }
232#endif
233
234 void setCmodeGroup(const Group& group, GroupState<Scalar>& group_state) const;
235
236 template <class AverageRegionalPressureType>
237 void
239 const FieldPropsManager& fp,
240 std::map<std::string, std::unique_ptr<AverageRegionalPressureType>>&
241 regional_average_pressure_calculator) const;
242
243 void setReportStep(int report_step)
244 {
245 report_step_ = report_step;
246 }
247
248#ifdef RESERVOIR_COUPLING_ENABLED
250 {
251 reservoir_coupling_master_ = reservoir_coupling_master;
252 }
253
255 {
256 reservoir_coupling_slave_ = reservoir_coupling_slave;
257 }
258#endif
259
260 const SummaryState& summaryState() const
261 {
262 return this->summary_state_;
263 }
264
265 Scalar sumSolventRates(const Group& group, const bool is_injector) const;
266
267 Scalar sumWellResRates(const Group& group, const int phase_pos, const bool injector) const;
268
269 Scalar sumWellSurfaceRates(const Group& group, const int phase_pos, const bool injector) const;
270
271 Scalar sumWellPhaseRates(bool res_rates,
272 const Group& group,
273 const int phase_pos,
274 const bool injector,
275 const bool network = false) const;
276
277 template <class RegionalValues>
278 void updateGpMaintTargetForGroups(const Group& group,
279 const RegionalValues& regional_values,
280 const double dt,
281 GroupState<Scalar>& group_state) const;
282
285 int updateGroupControlledWells(const bool is_production_group,
286 const Phase injection_phase,
287 GroupState<Scalar>& group_state,
288 DeferredLogger& deferred_logger) const;
289
290 void updateGroupProductionRates(const Group& group, GroupState<Scalar>& group_state) const;
291
292 void updateGroupTargetReduction(const Group& group,
293 const bool is_injector,
294 GroupState<Scalar>& group_state) const;
295
297
298 void updateREINForGroups(const Group& group, bool sum_rank, GroupState<Scalar>& group_state) const;
299
300 void updateReservoirRatesInjectionGroups(const Group& group, GroupState<Scalar>& group_state) const;
301
302 void updateVREPForGroups(const Group& group, GroupState<Scalar>& group_state) const;
303
305
306 void updateSurfaceRatesInjectionGroups(const Group& group, GroupState<Scalar>& group_state) const;
307
308 void updateWellRates(const Group& group,
309 const WellState<Scalar, IndexTraits>& well_state_nupcol,
310 WellState<Scalar, IndexTraits>& well_state) const;
311
313 {
314 return *this->well_state_;
315 }
316
317 void updateWellRatesFromGroupTargetScale(const Scalar scale,
318 const Group& group,
319 bool is_injector,
320 WellState<Scalar, IndexTraits>& well_state) const;
321
323 std::pair<std::optional<std::string>, Scalar>
324 worstOffendingWell(const Group& group,
325 const Group::ProductionCMode& offended_control,
326 const Parallel::Communication& comm,
327 DeferredLogger& deferred_logger) const;
328
329private:
330 std::string controlGroup_(const Group& group) const;
331
332 GuideRate::RateVector getGuideRateVector_(const std::vector<Scalar>& rates) const;
333
335 bool isInGroupChainTopBot_(const std::string& bottom, const std::string& top) const;
336
337 int phaseToActivePhaseIdx_(const Phase phase) const;
338
339 Scalar satelliteInjectionRate_(const ScheduleState& sched,
340 const Group& group,
341 const int phase_pos,
342 bool res_rates) const;
343
344 Scalar satelliteProductionRate_(const ScheduleState& sched,
345 const Group& group,
346 const GSatProd::GSatProdGroupProp::Rate rate_comp,
347 bool res_rates) const;
348
349 std::optional<GSatProd::GSatProdGroupProp::Rate> selectRateComponent_(const int phase_pos) const;
350
351 int updateGroupControlledWellsRecursive_(const std::string& group_name,
352 const bool is_production_group,
353 const Phase injection_phase,
354 GroupState<Scalar>& group_state,
355 DeferredLogger& deferred_logger) const;
356
357 void updateGroupTargetReductionRecursive_(const Group& group,
358 const bool is_injector,
359 std::vector<Scalar>& group_target_reduction,
360 GroupState<Scalar>& group_state) const;
361
362 const WellState<Scalar, IndexTraits>* well_state_ {nullptr};
363 const GroupState<Scalar>* group_state_ {nullptr};
364 const Schedule& schedule_;
365 const SummaryState& summary_state_;
366 const GuideRate& guide_rate_;
367 int report_step_ {0};
368 // NOTE: The phase usage info seems to be read-only throughout the simulation, so it should be safe
369 // to store a reference to it here.
370 const PhaseUsageInfo<IndexTraits>& phase_usage_info_;
371#ifdef RESERVOIR_COUPLING_ENABLED
372 ReservoirCouplingMaster<Scalar>* reservoir_coupling_master_{nullptr};
373 ReservoirCouplingSlave<Scalar>* reservoir_coupling_slave_{nullptr};
374#endif
375};
376
377// -----------------------------------------------------------------------------
378// Template member function implementations
379// -----------------------------------------------------------------------------
380
381// NOTE: This template member function is defined here in the header because the
382// AverageRegionalPressureType type parameter depends on derived class types that are
383// not available when WellGroupHelper.cpp is compiled. Template functions
384// must be visible at their instantiation point.
385// See WellGroupHelper.cpp for detailed rationale.
386template <typename Scalar, typename IndexTraits>
387template <class AverageRegionalPressureType>
388void
390 const Group& group,
391 const FieldPropsManager& fp,
392 std::map<std::string, std::unique_ptr<AverageRegionalPressureType>>& regional_average_pressure_calculator)
393 const
394{
395 for (const std::string& groupName : group.groups()) {
396 this->setRegionAveragePressureCalculator(this->schedule_.getGroup(groupName, this->report_step_),
397 fp,
398 regional_average_pressure_calculator);
399 }
400 const auto& gpm = group.gpmaint();
401 if (!gpm)
402 return;
403
404 const auto& reg = gpm->region();
405 if (!reg)
406 return;
407
408 if (regional_average_pressure_calculator.count(reg->first) == 0) {
409 const std::string name = (reg->first.rfind("FIP", 0) == 0) ? reg->first : "FIP" + reg->first;
410 const auto& fipnum = fp.get_int(name);
411 regional_average_pressure_calculator[reg->first]
412 = std::make_unique<AverageRegionalPressureType>(fipnum);
413 }
414}
415
416// NOTE: This template member function is defined in the header because the
417// RegionalValues type parameter depends on derived class types that are
418// not available when WellGroupHelper.cpp is compiled. Template functions
419// must be visible at their instantiation point.
420// See WellGroupHelper.cpp for detailed rationale.
421template <typename Scalar, typename IndexTraits>
422template <class RegionalValues>
423void
425 const RegionalValues& regional_values,
426 const double dt,
427 GroupState<Scalar>& group_state) const
428{
429 OPM_TIMEFUNCTION();
430 for (const std::string& group_name : group.groups()) {
431 const Group& group_tmp = this->schedule_.getGroup(group_name, this->report_step_);
432 this->updateGpMaintTargetForGroups(group_tmp, regional_values, dt, group_state);
433 }
434 const auto& gpm = group.gpmaint();
435 if (!gpm)
436 return;
437
438 const auto& region = gpm->region();
439 if (!region)
440 return;
441
442 const auto [name, number] = *region;
443 const Scalar error = gpm->pressure_target() - regional_values.at(name)->pressure(number);
444 Scalar current_rate = 0.0;
445 const auto& pu = this->phase_usage_info_;
446 bool injection = true;
447 Scalar sign = 1.0;
448 switch (gpm->flow_target()) {
449 case GPMaint::FlowTarget::RESV_PROD: {
450 current_rate = -this->groupState().injection_vrep_rate(group.name());
451 injection = false;
452 sign = -1.0;
453 break;
454 }
455 case GPMaint::FlowTarget::RESV_OINJ: {
456 if (pu.phaseIsActive(IndexTraits::oilPhaseIdx)) {
457 const auto io = pu.canonicalToActivePhaseIdx(IndexTraits::oilPhaseIdx);
458 current_rate = group_state.injection_reservoir_rates(group.name())[io];
459 }
460 break;
461 }
462 case GPMaint::FlowTarget::RESV_WINJ: {
463 if (pu.phaseIsActive(IndexTraits::waterPhaseIdx)) {
464 const auto iw = pu.canonicalToActivePhaseIdx(IndexTraits::waterPhaseIdx);
465 current_rate = group_state.injection_reservoir_rates(group.name())[iw];
466 }
467 break;
468 }
469 case GPMaint::FlowTarget::RESV_GINJ: {
470 if (pu.phaseIsActive(IndexTraits::gasPhaseIdx)) {
471 const auto ig = pu.canonicalToActivePhaseIdx(IndexTraits::gasPhaseIdx);
472 current_rate = group_state.injection_reservoir_rates(group.name())[ig];
473 }
474 break;
475 }
476 case GPMaint::FlowTarget::SURF_OINJ: {
477 if (pu.phaseIsActive(IndexTraits::oilPhaseIdx)) {
478 const auto io = pu.canonicalToActivePhaseIdx(IndexTraits::oilPhaseIdx);
479 current_rate = group_state.injection_surface_rates(group.name())[io];
480 }
481 break;
482 }
483 case GPMaint::FlowTarget::SURF_WINJ: {
484 if (pu.phaseIsActive(IndexTraits::waterPhaseIdx)) {
485 const auto iw = pu.canonicalToActivePhaseIdx(IndexTraits::waterPhaseIdx);
486 current_rate = group_state.injection_surface_rates(group.name())[iw];
487 }
488 break;
489 }
490 case GPMaint::FlowTarget::SURF_GINJ: {
491 if (pu.phaseIsActive(IndexTraits::gasPhaseIdx)) {
492 const auto ig = pu.canonicalToActivePhaseIdx(IndexTraits::gasPhaseIdx);
493 current_rate = group_state.injection_surface_rates(group.name())[ig];
494 }
495 break;
496 }
497 default:
498 throw std::invalid_argument("Invalid Flow target type in GPMAINT");
499 }
500 auto& gpmaint_state = group_state.gpmaint(group.name());
501 // we only activate gpmaint if pressure is lower than the target regional pressure for injectors
502 // (i.e. error > 0) and higher for producers.
503 bool activate = (injection && error > 0) || (!injection && error < 0);
504 Scalar rate = 0.0;
505 if (activate) {
506 rate = gpm->rate(gpmaint_state, current_rate, error, dt);
507 } else {
508 gpm->resetState(gpmaint_state);
509 }
510 group_state.update_gpmaint_target(group.name(), std::max(Scalar {0.0}, sign * rate));
511}
512
513} // namespace Opm
514
515#endif
Definition: DeferredLogger.hpp:57
Definition: GroupState.hpp:41
const std::vector< Scalar > & injection_reservoir_rates(const std::string &gname) const
const std::vector< Scalar > & injection_surface_rates(const std::string &gname) const
GPMaint::State & gpmaint(const std::string &gname)
void update_gpmaint_target(const std::string &gname, Scalar target)
Definition: GasLiftGroupInfo.hpp:37
Definition: ReservoirCouplingMaster.hpp:38
Definition: ReservoirCouplingSlave.hpp:40
Definition: VFPProdProperties.hpp:38
Definition: WellGroupHelper.hpp:92
GroupStateGuard & operator=(GroupStateGuard &&)=delete
GroupStateGuard(WellGroupHelper &wgHelper, GroupState< Scalar > &group_state)
Definition: WellGroupHelper.hpp:94
GroupStateGuard & operator=(const GroupStateGuard &)=delete
GroupStateGuard(const GroupStateGuard &)=delete
GroupStateGuard(GroupStateGuard &&)=delete
~GroupStateGuard()
Definition: WellGroupHelper.hpp:102
Definition: WellGroupHelper.hpp:63
WellStateGuard(WellGroupHelper &wgHelper, WellState< Scalar, IndexTraits > &well_state)
Definition: WellGroupHelper.hpp:65
WellStateGuard(WellStateGuard &&)=delete
~WellStateGuard()
Definition: WellGroupHelper.hpp:73
WellStateGuard(const WellStateGuard &)=delete
WellStateGuard & operator=(const WellStateGuard &)=delete
WellStateGuard & operator=(WellStateGuard &&)=delete
Definition: WellGroupHelper.hpp:59
void setReservoirCouplingSlave(ReservoirCouplingSlave< Scalar > *reservoir_coupling_slave)
Definition: WellGroupHelper.hpp:254
const GroupState< Scalar > & groupState() const
Definition: WellGroupHelper.hpp:183
WellGroupHelper(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)
Scalar sumWellPhaseRates(bool res_rates, const Group &group, const int phase_pos, const bool injector, const bool network=false) const
void updateVREPForGroups(const Group &group, GroupState< Scalar > &group_state) const
const Schedule & schedule() const
Definition: WellGroupHelper.hpp:208
std::vector< std::string > groupChainTopBot(const std::string &bottom, const std::string &top) const
Scalar sumSolventRates(const Group &group, const bool is_injector) const
const PhaseUsageInfo< IndexTraits > & phaseUsage() const
Definition: WellGroupHelper.hpp:218
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, DeferredLogger &deferred_logger) const
Scalar sumWellSurfaceRates(const Group &group, const int phase_pos, const bool injector) const
void updateGroupTargetReduction(const Group &group, const bool is_injector, GroupState< Scalar > &group_state) const
void updateREINForGroups(const Group &group, bool sum_rank, GroupState< Scalar > &group_state) const
GuideRate::RateVector getProductionGroupRateVector(const std::string &group_name) const
ReservoirCouplingMaster< Scalar > & reservoirCouplingMaster()
Definition: WellGroupHelper.hpp:223
int groupControlledWells(const std::string &group_name, const std::string &always_included_child, const bool is_production_group, const Phase injection_phase) const
void setReportStep(int report_step)
Definition: WellGroupHelper.hpp:243
void updateWellRatesFromGroupTargetScale(const Scalar scale, const Group &group, bool is_injector, WellState< Scalar, IndexTraits > &well_state) const
int updateGroupControlledWells(const bool is_production_group, const Phase injection_phase, GroupState< Scalar > &group_state, DeferredLogger &deferred_logger) const
const PhaseUsageInfo< IndexTraits > & phaseUsageInfo() const
Definition: WellGroupHelper.hpp:188
void setRegionAveragePressureCalculator(const Group &group, const FieldPropsManager &fp, std::map< std::string, std::unique_ptr< AverageRegionalPressureType > > &regional_average_pressure_calculator) const
Definition: WellGroupHelper.hpp:389
void updateGroupProductionRates(const Group &group, GroupState< Scalar > &group_state) const
Scalar sumWellResRates(const Group &group, const int phase_pos, const bool injector) const
ReservoirCouplingSlave< Scalar > & reservoirCouplingSlave()
Definition: WellGroupHelper.hpp:228
std::optional< Scalar > 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, DeferredLogger &deferred_logger) const
const SummaryState & summaryState() const
Definition: WellGroupHelper.hpp:260
void updateWellRates(const Group &group, const WellState< Scalar, IndexTraits > &well_state_nupcol, WellState< Scalar, IndexTraits > &well_state) const
void setCmodeGroup(const Group &group, GroupState< Scalar > &group_state) const
int reportStepIdx() const
Definition: WellGroupHelper.hpp:203
GuideRate::RateVector getWellRateVector(const std::string &name) const
void accumulateGroupEfficiencyFactor(const Group &group, Scalar &factor) const
GroupStateGuard pushGroupState(GroupState< Scalar > &group_state)
Definition: WellGroupHelper.hpp:193
const GuideRate & guideRate() const
Definition: WellGroupHelper.hpp:213
void updateNetworkLeafNodeProductionRates(GroupState< Scalar > &group_state) const
std::map< std::string, Scalar > computeNetworkPressures(const Network::ExtNetwork &network, const VFPProdProperties< Scalar > &vfp_prod_props, const Parallel::Communication &comm) 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, DeferredLogger &deferred_logger) const
const WellState< Scalar, IndexTraits > & wellState() const
Definition: WellGroupHelper.hpp:312
Scalar getGuideRate(const std::string &name, const GuideRateModel::Target target) const
std::optional< Scalar > 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, DeferredLogger &deferred_logger) const
std::pair< std::optional< std::string >, Scalar > worstOffendingWell(const Group &group, const Group::ProductionCMode &offended_control, const Parallel::Communication &comm, DeferredLogger &deferred_logger) const
Returns the name of the worst offending well and its fraction (i.e. violated_phase / preferred_phase)
void setReservoirCouplingMaster(ReservoirCouplingMaster< Scalar > *reservoir_coupling_master)
Definition: WellGroupHelper.hpp:249
void updateState(WellState< Scalar, IndexTraits > &well_state, GroupState< Scalar > &group_state)
WellStateGuard pushWellState(WellState< Scalar, IndexTraits > &well_state)
Definition: WellGroupHelper.hpp:198
void updateGpMaintTargetForGroups(const Group &group, const RegionalValues &regional_values, const double dt, GroupState< Scalar > &group_state) const
Definition: WellGroupHelper.hpp:424
void updateReservoirRatesInjectionGroups(const Group &group, GroupState< Scalar > &group_state) const
void updateSurfaceRatesInjectionGroups(const Group &group, GroupState< Scalar > &group_state) const
Definition: WellState.hpp:66
Dune::Communication< MPIComm > Communication
Definition: ParallelCommunication.hpp:30
Phase
Phase indices for reservoir coupling, we currently only support black-oil phases (oil,...
Definition: ReservoirCoupling.hpp:141
Definition: blackoilbioeffectsmodules.hh:43