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(dissolved_gas_rate);
62 serializer(dissolved_gas_rate_in_water);
63 serializer(vaporized_oil_rate);
64 serializer(vaporized_wat_rate);
65 serializer(well_potentials);
66 serializer(productivity_index);
67 serializer(surface_rates);
68 serializer(reservoir_rates);
69 serializer(trivial_target);
70 serializer(segments);
71 serializer(events);
72 serializer(injection_cmode);
73 serializer(production_cmode);
74 serializer(perf_data);
75 }
76
77 bool operator==(const SingleWellState&) const;
78
79 std::string name;
80 std::reference_wrapper<const ParallelWellInfo> parallel_info;
81
82 WellStatus status{WellStatus::OPEN};
83 bool producer;
84 PhaseUsage pu;
85 double bhp{0};
86 double thp{0};
87 double temperature{0};
88 double dissolved_gas_rate{0};
89 double dissolved_gas_rate_in_water{0};
90 double vaporized_oil_rate{0};
91 double vaporized_wat_rate{0};
92 std::vector<double> well_potentials;
93 std::vector<double> productivity_index;
94 std::vector<double> surface_rates;
95 std::vector<double> reservoir_rates;
96 PerfData perf_data;
97 bool trivial_target;
98 SegmentState segments;
99 Events events;
100 WellInjectorCMode injection_cmode{WellInjectorCMode::CMODE_UNDEFINED};
101 WellProducerCMode production_cmode{WellProducerCMode::CMODE_UNDEFINED};
102
103
110 void reset_connection_factors(const std::vector<PerforationData>& new_perf_data);
111 void update_producer_targets(const Well& ecl_well, const SummaryState& st);
112 void update_injector_targets(const Well& ecl_well, const SummaryState& st);
113 void update_targets(const Well& ecl_well, const SummaryState& st);
114 void updateStatus(WellStatus status);
115 void init_timestep(const SingleWellState& other);
116 void shut();
117 void stop();
118 void open();
119
120 // The sum_xxx_rates() functions sum over all connection rates of pertinent
121 // types. In the case of distributed wells this involves an MPI
122 // communication.
123 double sum_solvent_rates() const;
124 double sum_polymer_rates() const;
125 double sum_brine_rates() const;
126
127private:
128 double sum_connection_rates(const std::vector<double>& connection_rates) const;
129};
130
131
132}
133
134
135
136#endif
Definition AquiferInterface.hpp:35
Class encapsulating some information about parallel wells.
Definition ParallelWellInfo.hpp:184
Definition PerfData.hpp:30
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:119
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition BlackoilPhases.hpp:27
Definition BlackoilPhases.hpp:46