opm-common
Loading...
Searching...
No Matches
GuideRate.hpp
1/*
2 Copyright 2019, 2020 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 GUIDE_RATE_HPP
21#define GUIDE_RATE_HPP
22
23#include <opm/input/eclipse/Schedule/Group/Group.hpp>
24#include <opm/input/eclipse/Schedule/Group/GuideRateModel.hpp>
25
26#include <cstddef>
27#include <ctime>
28#include <limits>
29#include <memory>
30#include <string>
31#include <unordered_map>
32#include <unordered_set>
33#include <utility>
34
35namespace Opm {
36
37class Schedule;
38enum class WellGuideRateTarget : std::uint8_t;
39
40} // namespace Opm
41
42namespace Opm {
43
44class GuideRate
45{
46public:
47 // used for potentials and well rates
48 struct RateVector {
49 RateVector() = default;
50 RateVector(const double orat, const double grat, const double wrat)
51 : oil_rat(orat)
52 , gas_rat(grat)
53 , wat_rat(wrat)
54 {}
55
56 static RateVector serializationTestObject()
57 {
58 return RateVector{1.0, 2.0, 3.0};
59 }
60
61 double eval(const WellGuideRateTarget target) const;
62 double eval(const Group::GuideRateProdTarget target) const;
63 double eval(const GuideRateModel::Target target) const;
64
65 template<class Serializer>
66 void serializeOp(Serializer& serializer)
67 {
68 serializer(oil_rat);
69 serializer(gas_rat);
70 serializer(wat_rat);
71 }
72
73 double oil_rat{0.0};
74 double gas_rat{0.0};
75 double wat_rat{0.0};
76 };
77
78 struct GuideRateValue {
79 GuideRateValue() = default;
80 GuideRateValue(const double t, const double v, const GuideRateModel::Target tg)
81 : sim_time(t)
82 , value (v)
83 , target (tg)
84 {}
85
86 static GuideRateValue serializationTestObject()
87 {
88 return GuideRateValue{1.0, 2.0, GuideRateModel::Target::LIQ};
89 }
90
91 bool operator==(const GuideRateValue& other) const
92 {
93 return (this->sim_time == other.sim_time)
94 && (this->value == other.value);
95 }
96
97 bool operator!=(const GuideRateValue& other) const
98 {
99 return !(*this == other);
100 }
101
102 template<class Serializer>
103 void serializeOp(Serializer& serializer)
104 {
105 serializer(sim_time);
106 serializer(value);
107 serializer(target);
108 }
109
110 double sim_time { std::numeric_limits<double>::lowest() };
111 double value { std::numeric_limits<double>::lowest() };
112 GuideRateModel::Target target { GuideRateModel::Target::NONE };
113 };
114
115 explicit GuideRate(const Schedule& schedule);
116
117 void setSerializationTestData();
118
119 void compute(const std::string& wgname,
120 const std::size_t report_step,
121 const double sim_time,
122 const double oil_pot,
123 const double gas_pot,
124 const double wat_pot);
125
126 void compute(const std::string& wgname,
127 const Phase& phase,
128 const std::size_t report_step,
129 const std::optional<double> guide_rate);
130
131 bool has(const std::string& name) const;
132 bool hasPotentials(const std::string& name) const;
133 bool has(const std::string& name, const Phase& phase) const;
134
141 void erase(const std::string& name);
142
143 double get(const std::string& well, const WellGuideRateTarget target, const RateVector& rates) const;
144 double get(const std::string& group, const Group::GuideRateProdTarget target, const RateVector& rates) const;
145 double get(const std::string& name, const GuideRateModel::Target model_target, const RateVector& rates) const;
146 double get(const std::string& group, const Phase& phase) const;
147
148 double getSI(const std::string& well, const WellGuideRateTarget target, const RateVector& rates) const;
149 double getSI(const std::string& group, const Group::GuideRateProdTarget target, const RateVector& rates) const;
150 double getSI(const std::string& wgname, const GuideRateModel::Target target, const RateVector& rates) const;
151 double getSI(const std::string& group, const Phase& phase) const;
152
153 void init_grvalue(const std::size_t report_step, const std::string& wgname, GuideRateValue value);
154 void init_grvalue_SI(const std::size_t report_step, const std::string& wgname, GuideRateValue value);
155
156 void updateGuideRateExpiration(const double sim_time,
157 const std::size_t report_step);
158
159 template<class Serializer>
160 void serializeOp(Serializer& serializer)
161 {
162 serializer(values);
163 serializer(injection_group_values);
164 serializer(potentials);
165 serializer(potn_groups);
166 serializer(guide_rates_expired);
167 }
168
169private:
170 struct GRValState
171 {
172 GuideRateValue curr{};
173 GuideRateValue prev{};
174
175 static GRValState serializationTestObject()
176 {
177 return GRValState{GuideRateValue::serializationTestObject(),
178 GuideRateValue::serializationTestObject()};
179 }
180
181 template<class Serializer>
182 void serializeOp(Serializer& serializer)
183 {
184 serializer(curr);
185 serializer(prev);
186 }
187 };
188
189 struct pair_hash
190 {
191 template <class T1, class T2>
192 std::size_t operator()(const std::pair<T1, T2>& pair) const
193 {
194 return std::hash<T1>()(pair.first) ^ std::hash<T2>()(pair.second);
195 }
196 };
197
198 using GRValPtr = std::unique_ptr<GRValState>;
199 using pair = std::pair<Phase, std::string>;
200
201 void well_compute(const std::string& wgname,
202 const std::size_t report_step,
203 const double sim_time,
204 const double oil_pot,
205 const double gas_pot,
206 const double wat_pot);
207
208 void group_compute(const std::string& wgname,
209 const std::size_t report_step,
210 const double sim_time,
211 const double oil_pot,
212 const double gas_pot,
213 const double wat_pot);
214
215 double eval_form(const GuideRateModel& model,
216 const std::string& wgId,
217 const double oil_pot,
218 const double gas_pot,
219 const double wat_pot) const;
220
221 void assign_grvalue(const std::string& wgname,
222 const GuideRateModel& model,
223 GuideRateValue&& value);
224 double get_grvalue_result(const GRValState& gr) const;
225
226 const Schedule& schedule;
227
228 std::unordered_map<std::string, GRValPtr> values{};
229 std::unordered_map<pair, double, pair_hash> injection_group_values{};
230 std::unordered_map<std::string, RateVector> potentials{};
231 std::unordered_set<std::string> potn_groups{};
232 bool guide_rates_expired {false};
233};
234
235} // namespace Opm
236
237#endif
void erase(const std::string &name)
Erase the guide rate for a well or group.
Definition GuideRate.cpp:219
Definition Schedule.hpp:101
Class for (de-)serializing.
Definition Serializer.hpp:94
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30
Definition GuideRate.hpp:78
Definition GuideRate.hpp:48