My Project
VFPInjTable.hpp
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 
21 #ifndef OPM_PARSER_ECLIPSE_ECLIPSESTATE_TABLES_VFPINJTABLE_HPP_
22 #define OPM_PARSER_ECLIPSE_ECLIPSESTATE_TABLES_VFPINJTABLE_HPP_
23 
24 
25 #include <array>
26 #include <vector>
27 #include <string>
28 
29 #include <opm/common/OpmLog/KeywordLocation.hpp>
30 
31 namespace Opm {
32 
33 class DeckKeyword;
34 class UnitSystem;
35 
36 class VFPInjTable {
37 public:
38 
39  enum class FLO_TYPE {
40  FLO_OIL=1,
41  FLO_WAT,
42  FLO_GAS,
43  };
44 
45 
46  VFPInjTable();
47  VFPInjTable(const DeckKeyword& table, const UnitSystem& deck_unit_system);
48 
49  static VFPInjTable serializeObject();
50 
51  inline const KeywordLocation& location() const {
52  return this->m_location;
53  }
54 
55  inline int getTableNum() const {
56  return m_table_num;
57  }
58 
59  // The name() method is added to simplify serialization.
60  inline int name() const {
61  return m_table_num;
62  }
63 
64  inline double getDatumDepth() const {
65  return m_datum_depth;
66  }
67 
68  inline FLO_TYPE getFloType() const {
69  return m_flo_type;
70  }
71 
72  inline const std::vector<double>& getFloAxis() const {
73  return m_flo_data;
74  }
75 
76  inline const std::vector<double>& getTHPAxis() const {
77  return m_thp_data;
78  }
79 
92  inline const std::vector<double>& getTable() const {
93  return m_data;
94  }
95 
96  bool operator==(const VFPInjTable& data) const;
97 
98  std::array<size_t,2> shape() const;
99 
100  double operator()(size_t thp_idx, size_t flo_idx) const;
101 
102  template<class Serializer>
103  void serializeOp(Serializer& serializer)
104  {
105  serializer(m_table_num);
106  serializer(m_datum_depth);
107  serializer(m_flo_type);
108  serializer(m_flo_data);
109  serializer(m_thp_data);
110  serializer(m_data);
111  m_location.serializeOp(serializer);
112  }
113 
114 private:
115  int m_table_num;
116  double m_datum_depth;
117  FLO_TYPE m_flo_type;
118 
119  std::vector<double> m_flo_data;
120  std::vector<double> m_thp_data;
121 
122 
123  std::vector<double> m_data;
124  KeywordLocation m_location;
125 
126  void check();
127 
128  double& operator()(size_t thp_idx, size_t flo_idx);
129 
130  static FLO_TYPE getFloType(std::string flo_string);
131 
132  static void scaleValues(std::vector<double>& values,
133  const double& scaling_factor);
134 
135  static void convertFloToSI(const FLO_TYPE& type,
136  std::vector<double>& values,
137  const UnitSystem& unit_system);
138  static void convertTHPToSI(std::vector<double>& values,
139  const UnitSystem& unit_system);
140 };
141 
142 
143 }
144 
145 
146 #endif
Definition: DeckKeyword.hpp:36
Definition: KeywordLocation.hpp:27
Definition: Serializer.hpp:38
Definition: UnitSystem.hpp:34
Definition: VFPInjTable.hpp:36
const std::vector< double > & getTable() const
Returns the data of the table itself.
Definition: VFPInjTable.hpp:92
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition: Exceptions.hpp:29