BlackoilWellModelGasLift.hpp
Go to the documentation of this file.
1/*
2 Copyright 2016 SINTEF ICT, Applied Mathematics.
3 Copyright 2016 - 2017 Statoil ASA.
4 Copyright 2017 Dr. Blatt - HPC-Simulation-Software & Services
5 Copyright 2016 - 2018 IRIS AS
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#ifndef OPM_BLACKOILWELLMODEL_GASLIFT_HEADER_INCLUDED
23#define OPM_BLACKOILWELLMODEL_GASLIFT_HEADER_INCLUDED
24
28
29#include <memory>
30#include <map>
31#include <string>
32
33namespace Opm {
34
35class DeferredLogger;
36template<class Scalar> class GroupState;
37template<typename Scalar, typename IndexTraits> class WellState;
38template<class TypeTag> class WellInterface;
39
40template<typename Scalar, typename IndexTraits>
42{
43public:
44 using GLiftOptWells = std::map<std::string, std::unique_ptr<GasLiftSingleWellGeneric<Scalar, IndexTraits>>>;
45 using GLiftProdWells = std::map<std::string, const WellInterfaceGeneric<Scalar, IndexTraits>*>;
46 using GLiftWellStateMap = std::map<std::string, std::unique_ptr<GasLiftWellState<Scalar>>>;
49
50 explicit BlackoilWellModelGasLiftGeneric(bool terminal_output)
51 : terminal_output_(terminal_output)
52 {}
53
54 static constexpr bool glift_debug = false;
55
56 void gliftDebug(const std::string& msg,
57 DeferredLogger& deferred_logger) const;
58
59 bool terminalOutput() const { return terminal_output_; }
60
61 template<class Serializer>
62 void serializeOp(Serializer& serializer)
63 {
64 serializer(last_glift_opt_time_);
65 }
66
68 { return this->last_glift_opt_time_ == that.last_glift_opt_time_; }
69
70protected:
71 void gliftDebugShowALQ(const std::vector<WellInterfaceGeneric<Scalar, IndexTraits>*>& well_container,
72 const WellState<Scalar, IndexTraits>& wellState,
73 DeferredLogger& deferred_logger);
74
76 const Schedule& schedule,
77 const SummaryState& summaryState,
79 GroupState<Scalar>& groupState,
80 GLiftProdWells& prod_wells,
81 GLiftOptWells& glift_wells,
84 const int episodeIndex,
85 DeferredLogger& deferred_logger);
86
88 double last_glift_opt_time_ = -1.0;
89};
90
92template<typename TypeTag>
94 public BlackoilWellModelGasLiftGeneric<GetPropType<TypeTag, Properties::Scalar>,
95 typename GetPropType<TypeTag, Properties::FluidSystem>::IndexTraitsType>
96{
97public:
100 using IndexTraits = typename FluidSystem::IndexTraitsType;
102 using Base::glift_debug;
109 using WellInterfacePtr = std::shared_ptr<WellInterface<TypeTag>>;
111
112 explicit BlackoilWellModelGasLift(bool terminal_output)
113 : Base(terminal_output)
114 {}
115
116 static void initGliftEclWellMap(const std::vector<WellInterfacePtr>& well_container,
117 GLiftEclWells& ecl_well_map);
118
119 bool maybeDoGasLiftOptimize(const Simulator& simulator,
120 const std::vector<WellInterfacePtr>& well_container,
121 const std::map<std::string, Scalar>& node_pressures,
122 const bool updatePotentials,
123 WellStateType& wellState,
124 GroupState<Scalar>& groupState,
125 DeferredLogger& deferred_logger);
126
127private:
128 void gasLiftOptimizationStage1(const Simulator& simulator,
129 const std::vector<WellInterfacePtr>& well_container,
130 WellStateType& wellState,
131 GroupState<Scalar>& groupState,
132 GLiftProdWells& prod_wells,
133 GLiftOptWells& glift_wells,
135 GLiftWellStateMap& state_map,
136 DeferredLogger& deferred_logger);
137
138 // cannot be const since it accesses the non-const WellState
139 void gasLiftOptimizationStage1SingleWell(WellInterface<TypeTag>* well,
140 const Simulator& simulator,
141 WellStateType& wellState,
142 GroupState<Scalar>& groupState,
143 GLiftProdWells& prod_wells,
144 GLiftOptWells& glift_wells,
146 GLiftWellStateMap& state_map,
147 GLiftSyncGroups& groups_to_sync,
148 DeferredLogger& deferred_logger);
149
150 void updateWellPotentials(const Simulator& simulator,
151 const std::vector<WellInterfacePtr>& well_container,
152 const std::map<std::string, Scalar>& node_pressures,
154 DeferredLogger& deferred_logger);
155
156};
157
158} // namespace Opm
159
161
162#endif
Defines a type tags and some fundamental properties all models.
Definition: BlackoilWellModelGasLift.hpp:42
std::map< std::string, std::unique_ptr< GasLiftWellState< Scalar > > > GLiftWellStateMap
Definition: BlackoilWellModelGasLift.hpp:46
void gliftDebug(const std::string &msg, DeferredLogger &deferred_logger) const
std::map< std::string, const WellInterfaceGeneric< Scalar, IndexTraits > * > GLiftProdWells
Definition: BlackoilWellModelGasLift.hpp:45
typename GasLiftSingleWellGeneric< Scalar, IndexTraits >::GLiftSyncGroups GLiftSyncGroups
Definition: BlackoilWellModelGasLift.hpp:48
BlackoilWellModelGasLiftGeneric(bool terminal_output)
Definition: BlackoilWellModelGasLift.hpp:50
bool terminalOutput() const
Definition: BlackoilWellModelGasLift.hpp:59
static constexpr bool glift_debug
Definition: BlackoilWellModelGasLift.hpp:54
void gliftDebugShowALQ(const std::vector< WellInterfaceGeneric< Scalar, IndexTraits > * > &well_container, const WellState< Scalar, IndexTraits > &wellState, DeferredLogger &deferred_logger)
typename GasLiftGroupInfo< Scalar, IndexTraits >::GLiftEclWells GLiftEclWells
Definition: BlackoilWellModelGasLift.hpp:47
void serializeOp(Serializer &serializer)
Definition: BlackoilWellModelGasLift.hpp:62
bool terminal_output_
Definition: BlackoilWellModelGasLift.hpp:87
bool operator==(const BlackoilWellModelGasLiftGeneric &that) const
Definition: BlackoilWellModelGasLift.hpp:67
void gasLiftOptimizationStage2(const Parallel::Communication &comm, const Schedule &schedule, const SummaryState &summaryState, WellState< Scalar, IndexTraits > &wellState, GroupState< Scalar > &groupState, GLiftProdWells &prod_wells, GLiftOptWells &glift_wells, GasLiftGroupInfo< Scalar, IndexTraits > &group_info, GLiftWellStateMap &map, const int episodeIndex, DeferredLogger &deferred_logger)
double last_glift_opt_time_
Definition: BlackoilWellModelGasLift.hpp:88
std::map< std::string, std::unique_ptr< GasLiftSingleWellGeneric< Scalar, IndexTraits > > > GLiftOptWells
Definition: BlackoilWellModelGasLift.hpp:44
Class for handling the gaslift in the blackoil well model.
Definition: BlackoilWellModelGasLift.hpp:96
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: BlackoilWellModelGasLift.hpp:108
bool maybeDoGasLiftOptimize(const Simulator &simulator, const std::vector< WellInterfacePtr > &well_container, const std::map< std::string, Scalar > &node_pressures, const bool updatePotentials, WellStateType &wellState, GroupState< Scalar > &groupState, DeferredLogger &deferred_logger)
Definition: BlackoilWellModelGasLift_impl.hpp:43
static void initGliftEclWellMap(const std::vector< WellInterfacePtr > &well_container, GLiftEclWells &ecl_well_map)
Definition: BlackoilWellModelGasLift_impl.hpp:306
typename Base::GLiftWellStateMap GLiftWellStateMap
Definition: BlackoilWellModelGasLift.hpp:107
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: BlackoilWellModelGasLift.hpp:99
typename FluidSystem::IndexTraitsType IndexTraits
Definition: BlackoilWellModelGasLift.hpp:100
typename GasLiftSingleWellGeneric< Scalar, IndexTraits >::GLiftSyncGroups GLiftSyncGroups
Definition: BlackoilWellModelGasLift.hpp:106
typename Base::GLiftOptWells GLiftOptWells
Definition: BlackoilWellModelGasLift.hpp:104
std::shared_ptr< WellInterface< TypeTag > > WellInterfacePtr
Definition: BlackoilWellModelGasLift.hpp:109
BlackoilWellModelGasLift(bool terminal_output)
Definition: BlackoilWellModelGasLift.hpp:112
typename GasLiftGroupInfo< Scalar, IndexTraits >::GLiftEclWells GLiftEclWells
Definition: BlackoilWellModelGasLift.hpp:103
typename Base::GLiftProdWells GLiftProdWells
Definition: BlackoilWellModelGasLift.hpp:105
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: BlackoilWellModelGasLift.hpp:98
WellState< Scalar, IndexTraits > WellStateType
Definition: BlackoilWellModelGasLift.hpp:110
Definition: DeferredLogger.hpp:57
Definition: GasLiftGroupInfo.hpp:46
std::map< std::string, std::pair< const Well *, int > > GLiftEclWells
Definition: GasLiftGroupInfo.hpp:64
std::set< int > GLiftSyncGroups
Definition: GasLiftSingleWellGeneric.hpp:57
Definition: GroupState.hpp:41
Manages the initializing and running of time dependent problems.
Definition: simulator.hh:84
Definition: WellInterfaceGeneric.hpp:53
Definition: WellInterface.hpp:76
Definition: WellState.hpp:66
Dune::Communication< MPIComm > Communication
Definition: ParallelCommunication.hpp:30
Definition: blackoilboundaryratevector.hh:39
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:233