My Project
Loading...
Searching...
No Matches
SingleWellState.hpp
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/simulators/wells/SegmentState.hpp>
30#include <opm/simulators/wells/PerfData.hpp>
31#include <opm/simulators/wells/ParallelWellInfo.hpp>
32#include <opm/core/props/BlackoilPhases.hpp>
33
34namespace Opm {
35
36struct PerforationData;
37class SummaryState;
38class Well;
39
41public:
42 SingleWellState(const std::string& name,
43 const ParallelWellInfo& pinfo,
44 bool is_producer,
46 const std::vector<PerforationData>& perf_input,
47 const PhaseUsage& pu,
48 double temp);
49
50 static SingleWellState serializationTestObject(const ParallelWellInfo& pinfo);
51
52 template<class Serializer>
53 void serializeOp(Serializer& serializer)
54 {
55 serializer(name);
56 serializer(status);
57 serializer(producer);
58 serializer(bhp);
59 serializer(thp);
60 serializer(temperature);
61 serializer(phase_mixing_rates);
62 serializer(well_potentials);
63 serializer(productivity_index);
64 serializer(surface_rates);
65 serializer(reservoir_rates);
66 serializer(prev_surface_rates);
67 serializer(trivial_target);
68 serializer(segments);
69 serializer(events);
70 serializer(injection_cmode);
71 serializer(production_cmode);
72 serializer(filtrate_conc);
73 serializer(perf_data);
74 }
75
76 bool operator==(const SingleWellState&) const;
77
78 std::string name;
79 std::reference_wrapper<const ParallelWellInfo> parallel_info;
80
81 WellStatus status{WellStatus::OPEN};
82 bool producer;
83 PhaseUsage pu;
84 double bhp{0};
85 double thp{0};
86 double temperature{0};
87
88 // filtration injection concentration
89 double filtrate_conc{0};
90
91 std::array<double,4> phase_mixing_rates{};
92 enum RateIndices {
93 dissolved_gas = 0,
94 dissolved_gas_in_water = 1,
95 vaporized_oil = 2,
96 vaporized_water = 3
97 };
98
99 std::vector<double> well_potentials;
100 std::vector<double> productivity_index;
101 std::vector<double> surface_rates;
102 std::vector<double> reservoir_rates;
103 std::vector<double> prev_surface_rates;
104 PerfData perf_data;
105 bool trivial_target;
106 SegmentState segments;
107 Events events;
108 WellInjectorCMode injection_cmode{WellInjectorCMode::CMODE_UNDEFINED};
109 WellProducerCMode production_cmode{WellProducerCMode::CMODE_UNDEFINED};
110
111
118 void reset_connection_factors(const std::vector<PerforationData>& new_perf_data);
119 void update_producer_targets(const Well& ecl_well, const SummaryState& st);
120 void update_injector_targets(const Well& ecl_well, const SummaryState& st);
121 void update_targets(const Well& ecl_well, const SummaryState& st);
122 void updateStatus(WellStatus status);
123 void init_timestep(const SingleWellState& other);
124 void shut();
125 void stop();
126 void open();
127
128 // The sum_xxx_rates() functions sum over all connection rates of pertinent
129 // types. In the case of distributed wells this involves an MPI
130 // communication.
131 double sum_solvent_rates() const;
132 double sum_polymer_rates() const;
133 double sum_brine_rates() const;
134
135 double sum_filtrate_rate() const;
136 double sum_filtrate_total() const;
137
138private:
139 double sum_connection_rates(const std::vector<double>& connection_rates) const;
140};
141
142
143}
144
145
146
147#endif
Definition AquiferInterface.hpp:35
Class encapsulating some information about parallel wells.
Definition ParallelWellInfo.hpp:184
Definition PerfData.hpp:32
Definition SegmentState.hpp:35
Definition SingleWellState.hpp:40
void reset_connection_factors(const std::vector< PerforationData > &new_perf_data)
Special purpose method to support dynamically rescaling a well's CTFs through WELPI.
Definition SingleWellState.cpp:121
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition BlackoilPhases.hpp:27
Definition BlackoilPhases.hpp:46