opm-common
Loading...
Searching...
No Matches
TableManager.hpp
1/*
2 Copyright 2015 Statoil ASA.
3 Copyright 2018 IRIS
4
5 This file is part of the Open Porous Media project (OPM).
6
7 OPM is free software: you can redistribute it and/or modify
8 it under the terms of the GNU General Public License as published by
9 the Free Software Foundation, either version 3 of the License, or
10 (at your option) any later version.
11
12 OPM is distributed in the hope that it will be useful,
13 but WITHOUT ANY WARRANTY; without even the implied warranty of
14 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 GNU General Public License for more details.
16
17 You should have received a copy of the GNU General Public License
18 along with OPM. If not, see <http://www.gnu.org/licenses/>.
19 */
20
21#ifndef OPM_TABLE_MANAGER_HPP
22#define OPM_TABLE_MANAGER_HPP
23
24#include <cassert>
25#include <optional>
26#include <set>
27
28#include <opm/input/eclipse/EclipseState/Tables/DenT.hpp>
29#include <opm/input/eclipse/EclipseState/Tables/JouleThomson.hpp>
30#include <opm/input/eclipse/EclipseState/Tables/PvtgTable.hpp>
31#include <opm/input/eclipse/EclipseState/Tables/PvtgwTable.hpp>
32#include <opm/input/eclipse/EclipseState/Tables/PvtgwoTable.hpp>
33#include <opm/input/eclipse/EclipseState/Tables/PvtoTable.hpp>
34#include <opm/input/eclipse/EclipseState/Tables/PvtsolTable.hpp>
35#include <opm/input/eclipse/EclipseState/Tables/RocktabTable.hpp>
36#include <opm/input/eclipse/EclipseState/Tables/Rock2dTable.hpp>
37#include <opm/input/eclipse/EclipseState/Tables/Rock2dtrTable.hpp>
38#include <opm/input/eclipse/EclipseState/Tables/PlyshlogTable.hpp>
39#include <opm/input/eclipse/EclipseState/Tables/PvtwsaltTable.hpp>
40#include <opm/input/eclipse/EclipseState/Tables/RwgsaltTable.hpp>
41#include <opm/input/eclipse/EclipseState/Tables/BrineDensityTable.hpp>
42#include <opm/input/eclipse/EclipseState/Tables/SolventDensityTable.hpp>
43#include <opm/input/eclipse/EclipseState/Tables/StandardCond.hpp>
44#include <opm/input/eclipse/EclipseState/Tables/FlatTable.hpp>
45#include <opm/input/eclipse/EclipseState/Tables/SorwmisTable.hpp>
46#include <opm/input/eclipse/EclipseState/Tables/SgcwmisTable.hpp>
47#include <opm/input/eclipse/EclipseState/Tables/MiscTable.hpp>
48#include <opm/input/eclipse/EclipseState/Tables/PmiscTable.hpp>
49#include <opm/input/eclipse/EclipseState/Tables/MsfnTable.hpp>
50#include <opm/input/eclipse/EclipseState/Tables/JFunc.hpp>
51#include <opm/input/eclipse/EclipseState/Tables/Tabdims.hpp>
52#include <opm/input/eclipse/EclipseState/Tables/TableContainer.hpp>
53#include <opm/input/eclipse/EclipseState/Tables/Aqudims.hpp>
54#include <opm/input/eclipse/EclipseState/Tables/PlymwinjTable.hpp>
55#include <opm/input/eclipse/EclipseState/Tables/SkprwatTable.hpp>
56#include <opm/input/eclipse/EclipseState/Tables/SkprpolyTable.hpp>
57#include <opm/input/eclipse/EclipseState/Tables/Eqldims.hpp>
58#include <opm/input/eclipse/EclipseState/Tables/Regdims.hpp>
59#include <opm/input/eclipse/EclipseState/Tables/TLMixpar.hpp>
60#include <opm/input/eclipse/EclipseState/Tables/Ppcwmax.hpp>
61
62namespace Opm {
63
64 class Deck;
65
66 class TableManager {
67 public:
68 explicit TableManager( const Deck& deck );
69 TableManager() = default;
70
71 static TableManager serializationTestObject();
72
73 const TableContainer& getTables( const std::string& tableName ) const;
74 const TableContainer& operator[](const std::string& tableName) const;
75 bool hasTables( const std::string& tableName ) const;
76
77 const Tabdims& getTabdims() const;
78 const Eqldims& getEqldims() const;
79 const Aqudims& getAqudims() const;
80 const Regdims& getRegdims() const;
81 const TLMixpar& getTLMixpar() const;
82 const Ppcwmax& getPpcwmax() const;
83 /*
84 WIll return max{ Tabdims::NTFIP , Regdims::NTFIP }.
85 */
86 size_t numFIPRegions() const;
87
88 const TableContainer& getSwofTables() const;
89 const TableContainer& getSgwfnTables() const;
90 const TableContainer& getSof2Tables() const;
91 const TableContainer& getSof3Tables() const;
92 const TableContainer& getSgofTables() const;
93 const TableContainer& getSlgofTables() const;
94 const TableContainer& getSwfnTables() const;
95 const TableContainer& getSgfnTables() const;
96 const TableContainer& getSsfnTables() const;
97 const TableContainer& getRsvdTables() const;
98 const TableContainer& getRvvdTables() const;
99 const TableContainer& getRvwvdTables() const;
100 const TableContainer& getPbvdTables() const;
101 const TableContainer& getPdvdTables() const;
102 const TableContainer& getSaltvdTables() const;
103 const TableContainer& getSaltpvdTables() const;
104 const TableContainer& getSaltsolTables() const;
105 const TableContainer& getPcfactTables() const;
106 const TableContainer& getPermfactTables() const;
107 const TableContainer& getBiofilmTables() const;
108 const TableContainer& getDiffMICPTables() const;
109 const TableContainer& getEnkrvdTables() const;
110 const TableContainer& getEnptvdTables() const;
111 const TableContainer& getImkrvdTables() const;
112 const TableContainer& getImptvdTables() const;
113 const TableContainer& getPvdgTables() const;
114 const TableContainer& getPvdoTables() const;
115 const TableContainer& getRsconstTables() const;
116 const TableContainer& getPvdsTables() const;
117 const TableContainer& getSpecheatTables() const;
118 const TableContainer& getSpecrockTables() const;
119 const TableContainer& getWatvisctTables() const;
120 const TableContainer& getOilvisctTables() const;
121 const TableContainer& getGasvisctTables() const;
122 const TableContainer& getRtempvdTables() const;
123 const TableContainer& getZmfvdTables() const;
124 const TableContainer& getRocktabTables() const;
125 const TableContainer& getPlyadsTables() const;
126 const TableContainer& getPlyviscTables() const;
127 const TableContainer& getPlydhflfTables() const;
128 const TableContainer& getPlymaxTables() const;
129 const TableContainer& getPlyrockTables() const;
130 const TableContainer& getPlyshlogTables() const;
131 const TableContainer& getAqutabTables() const;
132 const TableContainer& getFoamadsTables() const;
133 const TableContainer& getFoammobTables() const;
134
135 const TableContainer& getSorwmisTables() const;
136 const TableContainer& getSgcwmisTables() const;
137 const TableContainer& getMiscTables() const;
138 const TableContainer& getPmiscTables() const;
139 const TableContainer& getMsfnTables() const;
140 const TableContainer& getTlpmixpaTables() const;
141
142 const TableContainer& getWsfTables() const;
143 const TableContainer& getGsfTables() const;
144
145 const JFunc& getJFunc() const;
146
147 const std::vector<PvtgTable>& getPvtgTables() const;
148 const std::vector<PvtgwTable>& getPvtgwTables() const;
149 const std::vector<PvtgwoTable>& getPvtgwoTables() const;
150 const std::vector<PvtoTable>& getPvtoTables() const;
151 const std::vector<PvtsolTable>& getPvtsolTables() const;
152 const std::vector<Rock2dTable>& getRock2dTables() const;
153 const std::vector<Rock2dtrTable>& getRock2dtrTables() const;
154 const TableContainer& getRockwnodTables() const;
155 const TableContainer& getOverburdTables() const;
156
157 const DenT& WatDenT() const;
158 const DenT& GasDenT() const;
159 const DenT& OilDenT() const;
160 const JouleThomson& WatJT() const;
161 const JouleThomson& GasJT() const;
162 const JouleThomson& OilJT() const;
163 const StandardCond& stCond() const;
164 std::size_t gas_comp_index() const;
165 const PvtwTable& getPvtwTable() const;
166 const std::vector<PvtwsaltTable>& getPvtwSaltTables() const;
167 const std::vector<RwgsaltTable>& getRwgSaltTables() const;
168 const std::vector<BrineDensityTable>& getBrineDensityTables() const;
169 const std::vector<SolventDensityTable>& getSolventDensityTables() const;
170
171 const PvcdoTable& getPvcdoTable() const;
172 const DensityTable& getDensityTable() const;
173 const DiffCoeffTable& getDiffusionCoefficientTable() const;
174 const DiffCoeffWatTable& getDiffusionCoefficientWaterTable() const;
175 const DiffCoeffGasTable& getDiffusionCoefficientGasTable() const;
176 const PlyvmhTable& getPlyvmhTable() const;
177 const RockTable& getRockTable() const;
178 const ViscrefTable& getViscrefTable() const;
179 const PlmixparTable& getPlmixparTable() const;
180 const ShrateTable& getShrateTable() const;
181 const Stone1exTable& getStone1exTable() const;
182 const WatdentTable& getWatdentTable() const;
183 const SgofletTable& getSgofletTable() const;
184 const SwofletTable& getSwofletTable() const;
185 const std::map<int, PlymwinjTable>& getPlymwinjTables() const;
186 const std::map<int, SkprwatTable>& getSkprwatTables() const;
187 const std::map<int, SkprpolyTable>& getSkprpolyTables() const;
188 const std::map<std::string, TableContainer>& getSimpleTables() const;
189
191 bool useImptvd() const;
192
194 bool useEnptvd() const;
195
197 bool useEqlnum() const;
198
200 bool useShrate() const;
201
203 bool useJFunc() const;
204
205 double rtemp() const;
206
207 double salinity() const;
208
209 bool diffMoleFraction() const;
210
211 bool operator==(const TableManager& data) const;
212
213 template<class Serializer>
214 void serializeOp(Serializer& serializer)
215 {
216 auto simpleTables = m_simpleTables;
217 auto split = splitSimpleTable(simpleTables);
218 serializer(simpleTables);
219 serializer(split.plyshMax);
220 serializer(split.plyshMap);
221 serializer(split.rockMax);
222 serializer(split.rockMap);
223 serializer(m_pvtgTables);
224 serializer(m_pvtgwTables);
225 serializer(m_pvtgwoTables);
226 serializer(m_pvtoTables);
227 serializer(m_pvtsolTables);
228 serializer(m_rock2dTables);
229 serializer(m_rock2dtrTables);
230 serializer(m_pvtwTable);
231 serializer(m_pvcdoTable);
232 serializer(m_densityTable);
233 serializer(m_diffCoeffTable);
234 serializer(m_diffCoeffWatTable);
235 serializer(m_diffCoeffGasTable);
236 serializer(m_plyvmhTable);
237 serializer(m_rockTable);
238 serializer(m_plmixparTable);
239 serializer(m_shrateTable);
240 serializer(m_stone1exTable);
241 serializer(m_viscrefTable);
242 serializer(m_watdentTable);
243 serializer(m_sgofletTable);
244 serializer(m_swofletTable);
245 serializer(m_pvtwsaltTables);
246 serializer(m_rwgsaltTables);
247 serializer(m_bdensityTables);
248 serializer(m_sdensityTables);
249 serializer(m_plymwinjTables);
250 serializer(m_skprwatTables);
251 serializer(m_skprpolyTables);
252 serializer(m_tabdims);
253 serializer(m_regdims);
254 serializer(m_eqldims);
255 serializer(m_aqudims);
256 serializer(hasImptvd);
257 serializer(hasEnptvd);
258 serializer(hasEqlnum);
259 serializer(hasShrate);
260 serializer(jfunc);
261 serializer(oilDenT);
262 serializer(gasDenT);
263 serializer(watDenT);
264 serializer(oilJT);
265 serializer(gasJT);
266 serializer(watJT);
267 serializer(stcond);
268 serializer(m_gas_comp_index);
269 serializer(m_rtemp);
270 serializer(m_salinity);
271 serializer(m_diff_mole_fraction);
272 serializer(m_tlmixpar);
273 serializer(m_ppcwmax);
274 if (!serializer.isSerializing()) {
275 m_simpleTables = simpleTables;
276 if (split.plyshMax > 0) {
277 TableContainer container(split.plyshMax);
278 for (const auto& it : split.plyshMap) {
279 container.addTable(it.first, it.second);
280 }
281 m_simpleTables.insert(std::make_pair("PLYSHLOG", container));
282 }
283 if (split.rockMax > 0) {
284 TableContainer container(split.rockMax);
285 for (const auto& it : split.rockMap) {
286 container.addTable(it.first, it.second);
287 }
288 m_simpleTables.insert(std::make_pair("ROCKTAB", container));
289 }
290 }
291 }
292
293 private:
294 TableContainer& forceGetTables( const std::string& tableName , size_t numTables);
295
296 void complainAboutAmbiguousKeyword(const Deck& deck, const std::string& keywordName);
297
298 void addTables( const std::string& tableName , size_t numTables);
299 void initSimpleTables(const Deck& deck);
300 void initRTempTables(const Deck& deck);
301 void initZmfvdTables(const Deck& deck);
302 void initDims(const Deck& deck);
303 void initRocktabTables(const Deck& deck);
304
305 void initPlymaxTables(const Deck& deck);
306 void initRsconstTables(const Deck& deck);
307 void initPlyrockTables(const Deck& deck);
308 void initPlyshlogTables(const Deck& deck);
309
310 void initPlymwinjTables(const Deck& deck);
311 void initSkprwatTables(const Deck& deck);
312 void initSkprpolyTables(const Deck& deck);
313
314 //void initRockTables(const Deck& deck, const std::string& keywordName);
315
316 template <class TableType>
317 void initRockTables(const Deck& deck, const std::string& keywordName, std::vector<TableType>& rocktable );
318
319 template <class TableType>
320 void initPvtwsaltTables(const Deck& deck, std::vector<TableType>& pvtwtables );
321
322 template <class TableType>
323 void initRwgsaltTables(const Deck& deck, std::vector<TableType>& rwgtables );
324
325 template <class TableType>
326 void initBrineTables(const Deck& deck, std::vector<TableType>& brinetables );
327
328 void initSolventTables(const Deck& deck, std::vector<SolventDensityTable>& solventtables);
329
333 template <class TableType>
334 void initSimpleTableContainerWithJFunc(const Deck& deck,
335 const std::string& keywordName,
336 const std::string& tableName,
337 size_t numTables);
338
339 template <class TableType>
340 void initSimpleTableContainer(const Deck& deck,
341 const std::string& keywordName,
342 const std::string& tableName,
343 size_t numTables);
344
345 template <class TableType>
346 void initSimpleTableContainer(const Deck& deck,
347 const std::string& keywordName,
348 size_t numTables);
349
350 template <class TableType>
351 void initSimpleTableContainerWithJFunc(const Deck& deck,
352 const std::string& keywordName,
353 size_t numTables);
354
355 template <class TableType>
356 void initSimpleTable(const Deck& deck,
357 const std::string& keywordName,
358 std::vector<TableType>& tableVector);
359
360 template <class TableType>
361 void initFullTables(const Deck& deck,
362 const std::string& keywordName,
363 std::vector<TableType>& tableVector);
364
365 void checkPVTOMonotonicity(const Deck& deck) const;
366
367 void logPVTOMonotonicityFailure(const Deck& deck,
368 const std::size_t tableID,
369 const std::vector<PvtoTable::FlippedFVF>& flipped_Bo) const;
370
371 std::map<std::string , TableContainer> m_simpleTables;
372 std::vector<PvtgTable> m_pvtgTables;
373 std::vector<PvtgwTable> m_pvtgwTables;
374 std::vector<PvtgwoTable> m_pvtgwoTables;
375 std::vector<PvtoTable> m_pvtoTables;
376 std::vector<PvtsolTable> m_pvtsolTables;
377 std::vector<Rock2dTable> m_rock2dTables;
378 std::vector<Rock2dtrTable> m_rock2dtrTables;
379 PvtwTable m_pvtwTable;
380 PvcdoTable m_pvcdoTable;
381 DensityTable m_densityTable;
382 DiffCoeffTable m_diffCoeffTable;
383 DiffCoeffWatTable m_diffCoeffWatTable;
384 DiffCoeffGasTable m_diffCoeffGasTable;
385 PlyvmhTable m_plyvmhTable;
386 RockTable m_rockTable;
387 PlmixparTable m_plmixparTable;
388 ShrateTable m_shrateTable;
389 Stone1exTable m_stone1exTable;
390 ViscrefTable m_viscrefTable;
391 WatdentTable m_watdentTable;
392 SgofletTable m_sgofletTable;
393 SwofletTable m_swofletTable;
394 std::vector<PvtwsaltTable> m_pvtwsaltTables;
395 std::vector<RwgsaltTable> m_rwgsaltTables;
396 std::vector<BrineDensityTable> m_bdensityTables;
397 std::vector<SolventDensityTable> m_sdensityTables;
398 std::map<int, PlymwinjTable> m_plymwinjTables;
399 std::map<int, SkprwatTable> m_skprwatTables;
400 std::map<int, SkprpolyTable> m_skprpolyTables;
401
402 Tabdims m_tabdims;
403 Regdims m_regdims;
404 Eqldims m_eqldims;
405 Aqudims m_aqudims;
406 TLMixpar m_tlmixpar;
407 Ppcwmax m_ppcwmax;
408
409 bool hasImptvd = false;// if deck has keyword IMPTVD
410 bool hasEnptvd = false;// if deck has keyword ENPTVD
411 bool hasEqlnum = false;// if deck has keyword EQLNUM
412 bool hasShrate = false;// if deck has keyword SHRATE
413 std::optional<JFunc> jfunc;
414
415 DenT oilDenT;
416 DenT gasDenT;
417 DenT watDenT;
418 JouleThomson oilJT;
419 JouleThomson gasJT;
420 JouleThomson watJT;
421 StandardCond stcond;
422 std::size_t m_gas_comp_index = 77;
423 double m_rtemp {288.7056}; // 60 Fahrenheit in Kelvin
424 double m_salinity {0.0};
425 bool m_diff_mole_fraction {true};
426
427 struct SplitSimpleTables {
428 size_t plyshMax = 0;
429 size_t rockMax = 0;
430 std::map<size_t, std::shared_ptr<PlyshlogTable>> plyshMap;
431 std::map<size_t, std::shared_ptr<RocktabTable>> rockMap;
432 };
433
434 SplitSimpleTables splitSimpleTable(std::map<std::string,TableContainer>& simpleTables);
435 };
436}
437
438
439#endif
Definition Aqudims.hpp:34
Definition Deck.hpp:46
Definition DenT.hpp:31
Definition Eqldims.hpp:32
Definition JFunc.hpp:28
Definition JouleThomson.hpp:30
Definition Ppcwmax.hpp:51
Definition Regdims.hpp:36
Class for (de-)serializing.
Definition Serializer.hpp:94
bool isSerializing() const
Returns true if we are currently doing a serialization operation.
Definition Serializer.hpp:207
Definition TLMixpar.hpp:54
Definition Tabdims.hpp:36
The TableContainer class implements a simple map:
Definition TableContainer.hpp:61
bool useEqlnum() const
deck has keyword "EQLNUM" — Equilibriation region numbers
Definition TableManager.cpp:1319
bool useJFunc() const
deck has keyword "JFUNC" — Use Leverett's J Function for capillary pressure
Definition TableManager.cpp:1327
bool useEnptvd() const
deck has keyword "ENPTVD" — Saturation end-point versus depth tables
Definition TableManager.cpp:1315
bool useImptvd() const
deck has keyword "IMPTVD" — Imbition end-point versus depth tables
Definition TableManager.cpp:1311
bool useShrate() const
deck has keyword "SHRATE"
Definition TableManager.cpp:1323
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30
Scalar Brine< Scalar, H2O >::salinity
Default value for the salinity of the brine (dimensionless).
Definition Brine.hpp:391
Definition FlatTable.hpp:124
Definition FlatTable.hpp:230
Definition FlatTable.hpp:173
Definition FlatTable.hpp:202
Definition FlatTable.hpp:363
Definition FlatTable.hpp:397
Definition FlatTable.hpp:338
Definition FlatTable.hpp:268
Definition FlatTable.hpp:299
Definition FlatTable.hpp:617
Definition FlatTable.hpp:422
Definition StandardCond.hpp:24
Definition FlatTable.hpp:447
Definition FlatTable.hpp:607
Definition FlatTable.hpp:503
Definition FlatTable.hpp:534