VFPProdTable.hpp
Go to the documentation of this file.
1/*
2 Copyright 2015 SINTEF ICT, Applied Mathematics.
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_PARSER_ECLIPSE_ECLIPSESTATE_TABLES_VFPPRODTABLE_HPP_
21#define OPM_PARSER_ECLIPSE_ECLIPSESTATE_TABLES_VFPPRODTABLE_HPP_
22
23
24#include <array>
25#include <vector>
26
27namespace Opm {
28
29 class DeckItem;
30 class DeckKeyword;
31 class UnitSystem;
32
37public:
38 typedef std::vector<double> array_type;
39
40 enum FLO_TYPE {
45 };
46
47 enum WFR_TYPE {
52 };
53
54 enum GFR_TYPE {
59 };
60
61 enum ALQ_TYPE {
70 };
71
73 VFPProdTable( const DeckKeyword& table, const UnitSystem& deck_unit_system);
74 VFPProdTable(int table_num,
75 double datum_depth,
76 FLO_TYPE flo_type,
77 WFR_TYPE wfr_type,
78 GFR_TYPE gfr_type,
79 ALQ_TYPE alq_type,
80 const std::vector<double>& flo_data,
81 const std::vector<double>& thp_data,
82 const std::vector<double>& wfr_data,
83 const std::vector<double>& gfr_data,
84 const std::vector<double>& alq_data,
85 const array_type& data);
86
88
89 inline int getTableNum() const {
90 return m_table_num;
91 }
92
93 inline double getDatumDepth() const {
94 return m_datum_depth;
95 }
96
97 inline FLO_TYPE getFloType() const {
98 return m_flo_type;
99 }
100
101 inline WFR_TYPE getWFRType() const {
102 return m_wfr_type;
103 }
104
105 inline GFR_TYPE getGFRType() const {
106 return m_gfr_type;
107 }
108
109 inline ALQ_TYPE getALQType() const {
110 return m_alq_type;
111 }
112
113 inline const std::vector<double>& getFloAxis() const {
114 return m_flo_data;
115 }
116
117 inline const std::vector<double>& getTHPAxis() const {
118 return m_thp_data;
119 }
120
121 inline const std::vector<double>& getWFRAxis() const {
122 return m_wfr_data;
123 }
124
125 inline const std::vector<double>& getGFRAxis() const {
126 return m_gfr_data;
127 }
128
129 inline const std::vector<double>& getALQAxis() const {
130 return m_alq_data;
131 }
132
147 inline const array_type& getTable() const {
148 return m_data;
149 }
150
151 bool operator==(const VFPProdTable& data) const;
152
153 std::array<size_t,5> shape() const;
154
155 double operator()(size_t thp_idx, size_t wfr_idx, size_t gfr_idx, size_t alq_idx, size_t flo_idx) const;
156
157 template<class Serializer>
158 void serializeOp(Serializer& serializer)
159 {
160 serializer(m_table_num);
161 serializer(m_datum_depth);
162 serializer(m_flo_type);
163 serializer(m_wfr_type);
164 serializer(m_gfr_type);
165 serializer(m_alq_type);
166 serializer(m_flo_data);
167 serializer(m_thp_data);
168 serializer(m_wfr_data);
169 serializer(m_gfr_data);
170 serializer(m_alq_data);
171 serializer(m_data);
172 }
173
174private:
175 int m_table_num;
176 double m_datum_depth;
177 FLO_TYPE m_flo_type;
178 WFR_TYPE m_wfr_type;
179 GFR_TYPE m_gfr_type;
180 ALQ_TYPE m_alq_type;
181
182 std::vector<double> m_flo_data;
183 std::vector<double> m_thp_data;
184 std::vector<double> m_wfr_data;
185 std::vector<double> m_gfr_data;
186 std::vector<double> m_alq_data;
187
188 array_type m_data;
189
190 void check(const DeckKeyword& table, const double factor);
191
192 double& operator()(size_t thp_idx, size_t wfr_idx, size_t gfr_idx, size_t alq_idx, size_t flo_idx);
193
194 static void scaleValues(std::vector<double>& values,
195 const double& scaling_factor);
196
197 static void convertFloToSI(const FLO_TYPE& type,
198 std::vector<double>& values,
199 const UnitSystem& unit_system);
200 static void convertTHPToSI(std::vector<double>& values,
201 const UnitSystem& unit_system);
202 static void convertWFRToSI(const WFR_TYPE& type,
203 std::vector<double>& values,
204 const UnitSystem& unit_system);
205 static void convertGFRToSI(const GFR_TYPE& type,
206 std::vector<double>& values,
207 const UnitSystem& unit_system);
208};
209
210
211
212}
213
214
215#endif
Definition: DeckKeyword.hpp:38
Definition: Serializer.hpp:38
Definition: UnitSystem.hpp:32
Definition: VFPProdTable.hpp:36
ALQ_TYPE getALQType() const
Definition: VFPProdTable.hpp:109
const std::vector< double > & getFloAxis() const
Definition: VFPProdTable.hpp:113
std::vector< double > array_type
Definition: VFPProdTable.hpp:38
FLO_TYPE getFloType() const
Definition: VFPProdTable.hpp:97
void serializeOp(Serializer &serializer)
Definition: VFPProdTable.hpp:158
VFPProdTable(int table_num, double datum_depth, FLO_TYPE flo_type, WFR_TYPE wfr_type, GFR_TYPE gfr_type, ALQ_TYPE alq_type, const std::vector< double > &flo_data, const std::vector< double > &thp_data, const std::vector< double > &wfr_data, const std::vector< double > &gfr_data, const std::vector< double > &alq_data, const array_type &data)
ALQ_TYPE
Definition: VFPProdTable.hpp:61
@ ALQ_UNDEF
Definition: VFPProdTable.hpp:68
@ ALQ_BEAN
Definition: VFPProdTable.hpp:67
@ ALQ_GRAT
Definition: VFPProdTable.hpp:62
@ ALQ_IGLR
Definition: VFPProdTable.hpp:63
@ ALQ_COMP
Definition: VFPProdTable.hpp:66
@ ALQ_PUMP
Definition: VFPProdTable.hpp:65
@ ALQ_TGLR
Definition: VFPProdTable.hpp:64
@ ALQ_INVALID
Definition: VFPProdTable.hpp:69
WFR_TYPE
Definition: VFPProdTable.hpp:47
@ WFR_WOR
Definition: VFPProdTable.hpp:48
@ WFR_INVALID
Definition: VFPProdTable.hpp:51
@ WFR_WCT
Definition: VFPProdTable.hpp:49
@ WFR_WGR
Definition: VFPProdTable.hpp:50
const std::vector< double > & getALQAxis() const
Definition: VFPProdTable.hpp:129
GFR_TYPE
Definition: VFPProdTable.hpp:54
@ GFR_OGR
Definition: VFPProdTable.hpp:57
@ GFR_INVALID
Definition: VFPProdTable.hpp:58
@ GFR_GLR
Definition: VFPProdTable.hpp:56
@ GFR_GOR
Definition: VFPProdTable.hpp:55
const std::vector< double > & getTHPAxis() const
Definition: VFPProdTable.hpp:117
std::array< size_t, 5 > shape() const
static VFPProdTable serializeObject()
bool operator==(const VFPProdTable &data) const
WFR_TYPE getWFRType() const
Definition: VFPProdTable.hpp:101
const std::vector< double > & getWFRAxis() const
Definition: VFPProdTable.hpp:121
const std::vector< double > & getGFRAxis() const
Definition: VFPProdTable.hpp:125
double getDatumDepth() const
Definition: VFPProdTable.hpp:93
const array_type & getTable() const
Definition: VFPProdTable.hpp:147
GFR_TYPE getGFRType() const
Definition: VFPProdTable.hpp:105
double operator()(size_t thp_idx, size_t wfr_idx, size_t gfr_idx, size_t alq_idx, size_t flo_idx) const
FLO_TYPE
Definition: VFPProdTable.hpp:40
@ FLO_OIL
Definition: VFPProdTable.hpp:41
@ FLO_GAS
Definition: VFPProdTable.hpp:43
@ FLO_INVALID
Definition: VFPProdTable.hpp:44
@ FLO_LIQ
Definition: VFPProdTable.hpp:42
VFPProdTable(const DeckKeyword &table, const UnitSystem &deck_unit_system)
int getTableNum() const
Definition: VFPProdTable.hpp:89
Definition: A.hpp:4
static std::string data()
Definition: exprtk.hpp:40022