SingleWellState.hpp
Go to the documentation of this file.
1/*
2 Copyright 2021 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
20#ifndef OPM_SINGLE_WELL_STATE_HEADER_INCLUDED
21#define OPM_SINGLE_WELL_STATE_HEADER_INCLUDED
22
23#include <functional>
24#include <vector>
25
26#include <opm/input/eclipse/Schedule/Well/WellEnums.hpp>
27#include <opm/input/eclipse/Schedule/Events.hpp>
28
29#include <opm/material/fluidsystems/PhaseUsageInfo.hpp>
30
35
36namespace Opm {
37
38template<class Scalar> struct PerforationData;
39class SummaryState;
40class Well;
41
42template<typename Scalar, typename IndexTraits>
44public:
48
49 SingleWellState(const std::string& name,
50 const ParallelWellInfo<Scalar>& pinfo,
52 bool is_producer,
54 const std::vector<PerforationData<Scalar>>& perf_input,
55 Scalar temp);
56
58
59 template<class Serializer>
60 void serializeOp(Serializer& serializer)
61 {
62 serializer(name);
63 serializer(status);
64 serializer(producer);
65 serializer(bhp);
66 serializer(thp);
67 serializer(pressure_first_connection);
68 serializer(temperature);
69 serializer(energy_rate);
70 serializer(efficiency_scaling_factor);
71 serializer(phase_mixing_rates);
72 serializer(well_potentials);
73 serializer(productivity_index);
74 serializer(implicit_ipr_a);
75 serializer(implicit_ipr_b);
76 serializer(surface_rates);
77 serializer(reservoir_rates);
78 serializer(prev_surface_rates);
79 serializer(trivial_group_target);
80 serializer(segments);
81 serializer(events);
82 serializer(injection_cmode);
83 serializer(production_cmode);
84 serializer(filtrate_conc);
85 serializer(perf_data);
86 serializer(primaryvar);
87 serializer(alq_state);
88 serializer(group_target);
90 }
91
92 bool operator==(const SingleWellState&) const;
93
94 std::string name;
95 std::reference_wrapper<const ParallelWellInfo<Scalar>> parallel_info;
96
97 WellStatus status{WellStatus::OPEN};
100 Scalar bhp{0};
101 Scalar thp{0};
103
104 // thermal related
105 Scalar temperature{0};
106 Scalar energy_rate{0.};
107
109
110 // filtration injection concentration
111 Scalar filtrate_conc{0};
112
113 std::array<Scalar,4> phase_mixing_rates{};
119 };
120
121 struct GroupTarget {
122 std::string group_name;
124
125 bool operator==(const GroupTarget& other) const {
126 return group_name == other.group_name && target_value == other.target_value;
127 }
128
129 template<class Serializer>
130 void serializeOp(Serializer& serializer) {
131 serializer(group_name);
132 serializer(target_value);
133 }
134 };
135
136 std::vector<Scalar> well_potentials;
137 std::vector<Scalar> productivity_index;
138 std::vector<Scalar> implicit_ipr_a;
139 std::vector<Scalar> implicit_ipr_b;
140 std::vector<Scalar> surface_rates;
141 std::vector<Scalar> reservoir_rates;
142 std::vector<Scalar> prev_surface_rates;
145 std::optional<GroupTarget> group_target;
147 Events events;
148 WellInjectorCMode injection_cmode{WellInjectorCMode::CMODE_UNDEFINED};
149 WellProducerCMode production_cmode{WellProducerCMode::CMODE_UNDEFINED};
150 std::vector<Scalar> primaryvar;
152 // This is used to indicate whether the well was shut before applying an action
153 // if it was SHUT, even the action set the well to OPEN, the data in the well state
154 // is not well-defined. We do not use it to overwrite the current well state.
156
163 void reset_connection_factors(const std::vector<PerforationData<Scalar>>& new_perf_data);
164 void update_producer_targets(const Well& ecl_well, const SummaryState& st);
165 void update_injector_targets(const Well& ecl_well, const SummaryState& st);
171 void update_type_and_targets(const Well& ecl_well, const SummaryState& st);
172 bool updateStatus(WellStatus status);
173 void init_timestep(const SingleWellState& other);
174 void shut();
175 void stop();
176 void open();
177
178 // The sum_xxx_rates() functions sum over all connection rates of pertinent
179 // types. In the case of distributed wells this involves an MPI
180 // communication.
181 Scalar sum_solvent_rates() const;
182 Scalar sum_polymer_rates() const;
183 Scalar sum_brine_rates() const;
184 Scalar sum_microbial_rates() const;
185 Scalar sum_oxygen_rates() const;
186 Scalar sum_urea_rates() const;
187 Scalar sum_wat_mass_rates() const;
188
189 Scalar sum_filtrate_rate() const;
190 Scalar sum_filtrate_total() const;
191
192private:
193 Scalar sum_connection_rates(const std::vector<Scalar>& connection_rates) const;
194};
195
196}
197
198#endif
Definition: ALQState.hpp:32
Class encapsulating some information about parallel wells.
Definition: ParallelWellInfo.hpp:198
Definition: PerfData.hpp:34
Definition: GasLiftGroupInfo.hpp:37
Definition: SegmentState.hpp:34
Definition: SingleWellState.hpp:43
WellProducerCMode production_cmode
Definition: SingleWellState.hpp:149
Scalar thp
Definition: SingleWellState.hpp:101
std::vector< Scalar > prev_surface_rates
Definition: SingleWellState.hpp:142
void reset_connection_factors(const std::vector< PerforationData< Scalar > > &new_perf_data)
static SingleWellState serializationTestObject(const ParallelWellInfo< Scalar > &pinfo)
std::array< Scalar, 4 > phase_mixing_rates
Definition: SingleWellState.hpp:113
std::vector< Scalar > surface_rates
Definition: SingleWellState.hpp:140
SingleWellState(const std::string &name, const ParallelWellInfo< Scalar > &pinfo, const PhaseUsageInfo< IndexTraits > &pu, bool is_producer, Scalar pressure_first_connection, const std::vector< PerforationData< Scalar > > &perf_input, Scalar temp)
std::optional< GroupTarget > group_target
Definition: SingleWellState.hpp:145
WellStatus status
Definition: SingleWellState.hpp:97
void init_timestep(const SingleWellState &other)
Scalar sum_solvent_rates() const
Scalar sum_oxygen_rates() const
std::vector< Scalar > productivity_index
Definition: SingleWellState.hpp:137
Scalar pressure_first_connection
Definition: SingleWellState.hpp:102
Scalar bhp
Definition: SingleWellState.hpp:100
std::reference_wrapper< const ParallelWellInfo< Scalar > > parallel_info
Definition: SingleWellState.hpp:95
void update_injector_targets(const Well &ecl_well, const SummaryState &st)
ALQState< Scalar > alq_state
Definition: SingleWellState.hpp:151
Scalar sum_filtrate_total() const
Scalar sum_microbial_rates() const
std::vector< Scalar > implicit_ipr_a
Definition: SingleWellState.hpp:138
void serializeOp(Serializer &serializer)
Definition: SingleWellState.hpp:60
Scalar sum_filtrate_rate() const
Scalar sum_wat_mass_rates() const
std::vector< Scalar > implicit_ipr_b
Definition: SingleWellState.hpp:139
RateIndices
Definition: SingleWellState.hpp:114
@ vaporized_water
Definition: SingleWellState.hpp:118
@ dissolved_gas_in_water
Definition: SingleWellState.hpp:116
@ dissolved_gas
Definition: SingleWellState.hpp:115
@ vaporized_oil
Definition: SingleWellState.hpp:117
static const int gasPhaseIdx
Definition: SingleWellState.hpp:47
bool operator==(const SingleWellState &) const
std::vector< Scalar > primaryvar
Definition: SingleWellState.hpp:150
Scalar sum_brine_rates() const
bool was_shut_before_action_applied
Definition: SingleWellState.hpp:155
SegmentState< Scalar > segments
Definition: SingleWellState.hpp:146
void update_type_and_targets(const Well &ecl_well, const SummaryState &st)
update the type of the well and the targets.
Scalar filtrate_conc
Definition: SingleWellState.hpp:111
std::vector< Scalar > well_potentials
Definition: SingleWellState.hpp:136
Scalar energy_rate
Definition: SingleWellState.hpp:106
bool updateStatus(WellStatus status)
static const int waterPhaseIdx
Definition: SingleWellState.hpp:45
std::string name
Definition: SingleWellState.hpp:94
Scalar sum_polymer_rates() const
Events events
Definition: SingleWellState.hpp:147
Scalar temperature
Definition: SingleWellState.hpp:105
bool trivial_group_target
Definition: SingleWellState.hpp:144
WellInjectorCMode injection_cmode
Definition: SingleWellState.hpp:148
PhaseUsageInfo< IndexTraits > pu
Definition: SingleWellState.hpp:99
bool producer
Definition: SingleWellState.hpp:98
static const int oilPhaseIdx
Definition: SingleWellState.hpp:46
void update_producer_targets(const Well &ecl_well, const SummaryState &st)
Scalar sum_urea_rates() const
Scalar efficiency_scaling_factor
Definition: SingleWellState.hpp:108
PerfData< Scalar > perf_data
Definition: SingleWellState.hpp:143
std::vector< Scalar > reservoir_rates
Definition: SingleWellState.hpp:141
Definition: blackoilbioeffectsmodules.hh:43
Static data associated with a well perforation.
Definition: PerforationData.hpp:30
Definition: SingleWellState.hpp:121
bool operator==(const GroupTarget &other) const
Definition: SingleWellState.hpp:125
Scalar target_value
Definition: SingleWellState.hpp:123
void serializeOp(Serializer &serializer)
Definition: SingleWellState.hpp:130
std::string group_name
Definition: SingleWellState.hpp:122