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
32
34
35namespace Opm {
36
37template<class Scalar> struct PerforationData;
38class SummaryState;
39class Well;
40
41template<class Scalar>
43public:
44 SingleWellState(const std::string& name,
45 const ParallelWellInfo<Scalar>& pinfo,
46 bool is_producer,
47 Scalar presssure_first_connection,
48 const std::vector<PerforationData<Scalar>>& perf_input,
49 const PhaseUsage& pu,
50 Scalar temp);
51
53
54 template<class Serializer>
55 void serializeOp(Serializer& serializer)
56 {
57 serializer(name);
58 serializer(status);
59 serializer(producer);
60 serializer(bhp);
61 serializer(thp);
62 serializer(temperature);
63 serializer(phase_mixing_rates);
64 serializer(well_potentials);
65 serializer(productivity_index);
66 serializer(implicit_ipr_a);
67 serializer(implicit_ipr_b);
68 serializer(surface_rates);
69 serializer(reservoir_rates);
70 serializer(prev_surface_rates);
71 serializer(trivial_target);
72 serializer(segments);
73 serializer(events);
74 serializer(injection_cmode);
75 serializer(production_cmode);
76 serializer(filtrate_conc);
77 serializer(perf_data);
78 }
79
80 bool operator==(const SingleWellState&) const;
81
82 std::string name;
83 std::reference_wrapper<const ParallelWellInfo<Scalar>> parallel_info;
84
85 WellStatus status{WellStatus::OPEN};
88 Scalar bhp{0};
89 Scalar thp{0};
90 Scalar temperature{0};
91
92 // filtration injection concentration
93 Scalar filtrate_conc{0};
94
95 std::array<Scalar,4> phase_mixing_rates{};
101 };
102
103 std::vector<Scalar> well_potentials;
104 std::vector<Scalar> productivity_index;
105 std::vector<Scalar> implicit_ipr_a;
106 std::vector<Scalar> implicit_ipr_b;
107 std::vector<Scalar> surface_rates;
108 std::vector<Scalar> reservoir_rates;
109 std::vector<Scalar> prev_surface_rates;
113 Events events;
114 WellInjectorCMode injection_cmode{WellInjectorCMode::CMODE_UNDEFINED};
115 WellProducerCMode production_cmode{WellProducerCMode::CMODE_UNDEFINED};
116
117
124 void reset_connection_factors(const std::vector<PerforationData<Scalar>>& new_perf_data);
125 void update_producer_targets(const Well& ecl_well, const SummaryState& st);
126 void update_injector_targets(const Well& ecl_well, const SummaryState& st);
127 void update_targets(const Well& ecl_well, const SummaryState& st);
128 void updateStatus(WellStatus status);
129 void init_timestep(const SingleWellState& other);
130 void shut();
131 void stop();
132 void open();
133
134 // The sum_xxx_rates() functions sum over all connection rates of pertinent
135 // types. In the case of distributed wells this involves an MPI
136 // communication.
137 Scalar sum_solvent_rates() const;
138 Scalar sum_polymer_rates() const;
139 Scalar sum_brine_rates() const;
140
141 Scalar sum_filtrate_rate() const;
142 Scalar sum_filtrate_total() const;
143
144private:
145 Scalar sum_connection_rates(const std::vector<Scalar>& connection_rates) const;
146};
147
148}
149
150#endif
Class encapsulating some information about parallel wells.
Definition: ParallelWellInfo.hpp:186
Definition: PerfData.hpp:33
Definition: SegmentState.hpp:34
Definition: SingleWellState.hpp:42
void updateStatus(WellStatus status)
Scalar sum_polymer_rates() const
std::vector< Scalar > surface_rates
Definition: SingleWellState.hpp:107
PhaseUsage pu
Definition: SingleWellState.hpp:87
Scalar sum_filtrate_total() const
std::vector< Scalar > implicit_ipr_b
Definition: SingleWellState.hpp:106
Events events
Definition: SingleWellState.hpp:113
std::vector< Scalar > prev_surface_rates
Definition: SingleWellState.hpp:109
Scalar sum_brine_rates() const
WellInjectorCMode injection_cmode
Definition: SingleWellState.hpp:114
SegmentState< Scalar > segments
Definition: SingleWellState.hpp:112
std::string name
Definition: SingleWellState.hpp:82
void init_timestep(const SingleWellState &other)
WellProducerCMode production_cmode
Definition: SingleWellState.hpp:115
std::vector< Scalar > well_potentials
Definition: SingleWellState.hpp:103
bool trivial_target
Definition: SingleWellState.hpp:111
std::vector< Scalar > implicit_ipr_a
Definition: SingleWellState.hpp:105
Scalar thp
Definition: SingleWellState.hpp:89
bool producer
Definition: SingleWellState.hpp:86
void update_targets(const Well &ecl_well, const SummaryState &st)
void update_injector_targets(const Well &ecl_well, const SummaryState &st)
Scalar filtrate_conc
Definition: SingleWellState.hpp:93
PerfData< Scalar > perf_data
Definition: SingleWellState.hpp:110
Scalar sum_solvent_rates() const
RateIndices
Definition: SingleWellState.hpp:96
@ dissolved_gas_in_water
Definition: SingleWellState.hpp:98
@ vaporized_oil
Definition: SingleWellState.hpp:99
@ dissolved_gas
Definition: SingleWellState.hpp:97
@ vaporized_water
Definition: SingleWellState.hpp:100
bool operator==(const SingleWellState &) const
Scalar bhp
Definition: SingleWellState.hpp:88
SingleWellState(const std::string &name, const ParallelWellInfo< Scalar > &pinfo, bool is_producer, Scalar presssure_first_connection, const std::vector< PerforationData< Scalar > > &perf_input, const PhaseUsage &pu, Scalar temp)
void reset_connection_factors(const std::vector< PerforationData< Scalar > > &new_perf_data)
static SingleWellState serializationTestObject(const ParallelWellInfo< Scalar > &pinfo)
Scalar temperature
Definition: SingleWellState.hpp:90
void update_producer_targets(const Well &ecl_well, const SummaryState &st)
Scalar sum_filtrate_rate() const
WellStatus status
Definition: SingleWellState.hpp:85
std::reference_wrapper< const ParallelWellInfo< Scalar > > parallel_info
Definition: SingleWellState.hpp:83
std::vector< Scalar > reservoir_rates
Definition: SingleWellState.hpp:108
std::vector< Scalar > productivity_index
Definition: SingleWellState.hpp:104
void serializeOp(Serializer &serializer)
Definition: SingleWellState.hpp:55
std::array< Scalar, 4 > phase_mixing_rates
Definition: SingleWellState.hpp:95
Definition: blackoilboundaryratevector.hh:37
Static data associated with a well perforation.
Definition: PerforationData.hpp:30
Definition: BlackoilPhases.hpp:46