WellInterfaceFluidSystem.hpp
Go to the documentation of this file.
1/*
2 Copyright 2017 SINTEF Digital, Mathematics and Cybernetics.
3 Copyright 2017 Statoil ASA.
4 Copyright 2017 IRIS
5 Copyright 2019 Norce
6
7 This file is part of the Open Porous Media project (OPM).
8
9 OPM is free software: you can redistribute it and/or modify
10 it under the terms of the GNU General Public License as published by
11 the Free Software Foundation, either version 3 of the License, or
12 (at your option) any later version.
13
14 OPM is distributed in the hope that it will be useful,
15 but WITHOUT ANY WARRANTY; without even the implied warranty of
16 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17 GNU General Public License for more details.
18
19 You should have received a copy of the GNU General Public License
20 along with OPM. If not, see <http://www.gnu.org/licenses/>.
21*/
22
23
24#ifndef OPM_WELLINTERFACE_FLUID_SYSTEM_HEADER_INCLUDED
25#define OPM_WELLINTERFACE_FLUID_SYSTEM_HEADER_INCLUDED
26
28
30
31#include <limits>
32#include <optional>
33#include <vector>
34
35namespace Opm
36{
37namespace RateConverter
38{
39 template <class FluidSystem, class Region> class SurfaceToReservoirVoidage;
40}
41
42class Group;
43template<class Scalar> class GroupState;
44class Schedule;
45struct RatioLimitCheckReport;
46template<class Scalar> class SingleWellState;
47template<class Scalar> class WellState;
48
49template<class FluidSystem>
50class WellInterfaceFluidSystem : public WellInterfaceGeneric<typename FluidSystem::Scalar>
51{
52protected:
53 using RateConverterType = RateConverter::
54 SurfaceToReservoirVoidage<FluidSystem, std::vector<int>>;
55 // to indicate a invalid completion
56 static constexpr int INVALIDCOMPLETION = std::numeric_limits<int>::max();
57
58public:
59 using Scalar = typename FluidSystem::Scalar;
60
61 int flowPhaseToModelPhaseIdx(const int phaseIdx) const;
62
63 static constexpr int Water = BlackoilPhases::Aqua;
64 static constexpr int Oil = BlackoilPhases::Liquid;
65 static constexpr int Gas = BlackoilPhases::Vapour;
66
68 {
69 return rateConverter_;
70 }
71
72protected:
73 WellInterfaceFluidSystem(const Well& well,
74 const ParallelWellInfo<Scalar>& parallel_well_info,
75 const int time_step,
76 const RateConverterType& rate_converter,
77 const int pvtRegionIdx,
78 const int num_components,
79 const int num_phases,
80 const int index_of_well,
81 const std::vector<PerforationData<Scalar>>& perf_data);
82
83 // updating the voidage rates in well_state when requested
84 void calculateReservoirRates(const bool co2store, SingleWellState<Scalar>& ws) const;
85
87 const SummaryState& summaryState,
88 DeferredLogger& deferred_logger,
89 const std::optional<Well::InjectionControls>& inj_controls = std::nullopt,
90 const std::optional<Well::ProductionControls>& prod_controls = std::nullopt) const;
91
93 const GroupState<Scalar>& group_state,
94 const Schedule& schedule,
95 const SummaryState& summaryState,
96 DeferredLogger& deferred_logger) const;
97
99 const GroupState<Scalar>& group_state,
100 const Schedule& schedule,
101 const SummaryState& summaryState,
102 DeferredLogger& deferred_logger) const;
103
104 std::optional<Scalar>
105 getGroupInjectionTargetRate(const Group& group,
106 const WellState<Scalar>& well_state,
107 const GroupState<Scalar>& group_state,
108 const Schedule& schedule,
109 const SummaryState& summaryState,
110 const InjectorType& injectorType,
111 Scalar efficiencyFactor,
112 DeferredLogger& deferred_logger) const;
113
114 Scalar
116 const WellState<Scalar>& well_state,
117 const GroupState<Scalar>& group_state,
118 const Schedule& schedule,
119 const SummaryState& summaryState,
120 Scalar efficiencyFactor,
121 DeferredLogger& deferred_logger) const;
122
123 bool zeroGroupRateTarget(const SummaryState& summary_state,
124 const Schedule& schedule,
125 const WellState<Scalar>& well_state,
126 const GroupState<Scalar>& group_state,
127 DeferredLogger& deferredLogger) const;
128
129 // For the conversion between the surface volume rate and reservoir voidage rate
131};
132
133}
134
135#endif // OPM_WELLINTERFACE_FLUID_SYSTEM_HEADER_INCLUDED
@ Liquid
Definition: BlackoilPhases.hpp:42
@ Aqua
Definition: BlackoilPhases.hpp:42
@ Vapour
Definition: BlackoilPhases.hpp:42
Definition: DeferredLogger.hpp:57
Definition: GroupState.hpp:38
Class encapsulating some information about parallel wells.
Definition: ParallelWellInfo.hpp:186
Definition: RateConverter.hpp:71
Definition: SingleWellState.hpp:42
Definition: WellInterfaceFluidSystem.hpp:51
bool checkConstraints(WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, const Schedule &schedule, const SummaryState &summaryState, DeferredLogger &deferred_logger) const
std::optional< Scalar > getGroupInjectionTargetRate(const Group &group, const WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, const Schedule &schedule, const SummaryState &summaryState, const InjectorType &injectorType, Scalar efficiencyFactor, DeferredLogger &deferred_logger) const
WellInterfaceFluidSystem(const Well &well, const ParallelWellInfo< Scalar > &parallel_well_info, const int time_step, const RateConverterType &rate_converter, const int pvtRegionIdx, const int num_components, const int num_phases, const int index_of_well, const std::vector< PerforationData< Scalar > > &perf_data)
const RateConverterType & rateConverter() const
Definition: WellInterfaceFluidSystem.hpp:67
bool checkGroupConstraints(WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, const Schedule &schedule, const SummaryState &summaryState, DeferredLogger &deferred_logger) const
bool checkIndividualConstraints(SingleWellState< Scalar > &ws, const SummaryState &summaryState, DeferredLogger &deferred_logger, const std::optional< Well::InjectionControls > &inj_controls=std::nullopt, const std::optional< Well::ProductionControls > &prod_controls=std::nullopt) const
bool zeroGroupRateTarget(const SummaryState &summary_state, const Schedule &schedule, const WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferredLogger) const
static constexpr int Oil
Definition: WellInterfaceFluidSystem.hpp:64
void calculateReservoirRates(const bool co2store, SingleWellState< Scalar > &ws) const
int flowPhaseToModelPhaseIdx(const int phaseIdx) const
Scalar getGroupProductionTargetRate(const Group &group, const WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, const Schedule &schedule, const SummaryState &summaryState, Scalar efficiencyFactor, DeferredLogger &deferred_logger) const
static constexpr int Water
Definition: WellInterfaceFluidSystem.hpp:63
typename FluidSystem::Scalar Scalar
Definition: WellInterfaceFluidSystem.hpp:59
static constexpr int Gas
Definition: WellInterfaceFluidSystem.hpp:65
const RateConverterType & rateConverter_
Definition: WellInterfaceFluidSystem.hpp:130
static constexpr int INVALIDCOMPLETION
Definition: WellInterfaceFluidSystem.hpp:56
Definition: WellInterfaceGeneric.hpp:51
int pvtRegionIdx() const
Definition: WellInterfaceGeneric.hpp:119
Definition: WellState.hpp:62
Definition: blackoilboundaryratevector.hh:37
Static data associated with a well perforation.
Definition: PerforationData.hpp:30