opm-simulators
Loading...
Searching...
No Matches
WellInterfaceFluidSystem.hpp
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
27#include <opm/simulators/wells/WellInterfaceGeneric.hpp>
28
29#include <limits>
30#include <optional>
31#include <vector>
32
33namespace Opm
34{
35namespace RateConverter
36{
37 template <class FluidSystem, class Region> class SurfaceToReservoirVoidage;
38}
39
40class Group;
41template<class Scalar> class GroupState;
42class Schedule;
43struct RatioLimitCheckReport;
44template<typename Scalar, typename IndexTraits> class SingleWellState;
45template<typename Scalar, typename IndexTraits> class WellState;
46
47template<class FluidSystem>
48class WellInterfaceFluidSystem : public WellInterfaceGeneric<typename FluidSystem::Scalar, typename FluidSystem::IndexTraitsType>
49{
50protected:
51 using RateConverterType = RateConverter::
52 SurfaceToReservoirVoidage<FluidSystem, std::vector<int>>;
53 // to indicate a invalid completion
54 static constexpr int INVALIDCOMPLETION = std::numeric_limits<int>::max();
55
56public:
57 using Scalar = typename FluidSystem::Scalar;
58 using IndexTraits = typename FluidSystem::IndexTraitsType;
59 using ModelParameters = typename WellInterfaceGeneric<Scalar, IndexTraits>::ModelParameters;
60
61 static constexpr int Water = IndexTraits::waterPhaseIdx;
62 static constexpr int Oil = IndexTraits::oilPhaseIdx;
63 static constexpr int Gas = IndexTraits::gasPhaseIdx;
64
65 const RateConverterType& rateConverter() const
66 {
67 return rateConverter_;
68 }
69
70protected:
71 WellInterfaceFluidSystem(const Well& well,
72 const ParallelWellInfo<Scalar>& parallel_well_info,
73 const int time_step,
74 const ModelParameters& param,
75 const RateConverterType& rate_converter,
76 const int pvtRegionIdx,
77 const int num_conservation_quantities,
78 const int num_phases,
79 const int index_of_well,
80 const std::vector<PerforationData<Scalar>>& perf_data);
81
82 // updating the voidage rates in well_state when requested
83 void calculateReservoirRates(const bool co2store, SingleWellState<Scalar, IndexTraits>& ws) const;
84
85 bool checkIndividualConstraints(SingleWellState<Scalar, IndexTraits>& ws,
86 const SummaryState& summaryState,
87 DeferredLogger& deferred_logger,
88 const std::optional<Well::InjectionControls>& inj_controls = std::nullopt,
89 const std::optional<Well::ProductionControls>& prod_controls = std::nullopt) const;
90
91 bool checkGroupConstraints(WellState<Scalar, IndexTraits>& well_state,
92 const GroupState<Scalar>& group_state,
93 const Schedule& schedule,
94 const SummaryState& summaryState,
95 const bool check_guide_rate,
96 DeferredLogger& deferred_logger) const;
97
98 bool checkConstraints(WellState<Scalar, IndexTraits>& well_state,
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, IndexTraits>& 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
115 getGroupProductionTargetRate(const Group& group,
116 const WellState<Scalar, IndexTraits>& 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, IndexTraits>& 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
130 const RateConverterType& rateConverter_;
131};
132
133}
134
135#endif // OPM_WELLINTERFACE_FLUID_SYSTEM_HEADER_INCLUDED
Definition DeferredLogger.hpp:57
Definition GroupState.hpp:41
Class encapsulating some information about parallel wells.
Definition ParallelWellInfo.hpp:198
Convert component rates at surface conditions to phase (voidage) rates at reservoir conditions.
Definition RateConverter.hpp:70
Definition SingleWellState.hpp:43
The state of a set of wells, tailored for use by the fully implicit blackoil simulator.
Definition WellState.hpp:66
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:43
Static data associated with a well perforation.
Definition PerforationData.hpp:30