TableManager.hpp
Go to the documentation of this file.
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 <set>
26
29
33
60
61namespace Opm {
62
64 public:
65 explicit TableManager( const Deck& deck );
66 TableManager() = default;
67
69
71
72 const TableContainer& getTables( const std::string& tableName ) const;
73 const TableContainer& operator[](const std::string& tableName) const;
74 bool hasTables( const std::string& tableName ) const;
75
76 const Tabdims& getTabdims() const;
77 const Eqldims& getEqldims() const;
78 const Aqudims& getAqudims() const;
79 const Regdims& getRegdims() const;
80
81 /*
82 WIll return max{ Tabdims::NTFIP , Regdims::NTFIP }.
83 */
84 size_t numFIPRegions() const;
85
123
130
131 const JFunc& getJFunc() const;
132
133 const std::vector<PvtgTable>& getPvtgTables() const;
134 const std::vector<PvtoTable>& getPvtoTables() const;
135 const std::vector<Rock2dTable>& getRock2dTables() const;
136 const std::vector<Rock2dtrTable>& getRock2dtrTables() const;
139
140 const DenT& WatDenT() const;
141 const DenT& GasDenT() const;
142 const DenT& OilDenT() const;
143 const StandardCond& stCond() const;
144 std::size_t gas_comp_index() const;
145 const PvtwTable& getPvtwTable() const;
146 const std::vector<PvtwsaltTable>& getPvtwSaltTables() const;
147 const std::vector<BrineDensityTable>& getBrineDensityTables() const;
148 const std::vector<SolventDensityTable>& getSolventDensityTables() const;
149
150 const PvcdoTable& getPvcdoTable() const;
153 const RockTable& getRockTable() const;
160 const std::map<int, PlymwinjTable>& getPlymwinjTables() const;
161 const std::map<int, SkprwatTable>& getSkprwatTables() const;
162 const std::map<int, SkprpolyTable>& getSkprpolyTables() const;
163 const std::map<std::string, TableContainer>& getSimpleTables() const;
164
166 bool useImptvd() const;
167
169 bool useEnptvd() const;
170
172 bool useEqlnum() const;
173
175 bool useShrate() const;
176
178 bool useJFunc() const;
179
180 double rtemp() const;
181
182 bool operator==(const TableManager& data) const;
183
184 template<class Serializer>
185 void serializeOp(Serializer& serializer)
186 {
187 auto simpleTables = m_simpleTables;
188 auto split = splitSimpleTable(simpleTables);
189 serializer.map(simpleTables);
190 serializer(split.plyshMax);
191 serializer.map(split.plyshMap);
192 serializer(split.rockMax);
193 serializer.map(split.rockMap);
194 serializer.vector(m_pvtgTables);
195 serializer.vector(m_pvtoTables);
196 serializer.vector(m_rock2dTables);
197 serializer.vector(m_rock2dtrTables);
198 m_pvtwTable.serializeOp(serializer);
199 m_pvcdoTable.serializeOp(serializer);
200 m_densityTable.serializeOp(serializer);
201 m_plyvmhTable.serializeOp(serializer);
202 m_rockTable.serializeOp(serializer);
203 m_plmixparTable.serializeOp(serializer);
204 m_shrateTable.serializeOp(serializer);
205 m_stone1exTable.serializeOp(serializer);
206 m_tlmixparTable.serializeOp(serializer);
207 m_viscrefTable.serializeOp(serializer);
208 m_watdentTable.serializeOp(serializer);
209 serializer.vector(m_pvtwsaltTables);
210 serializer.vector(m_bdensityTables);
211 serializer.vector(m_sdensityTables);
212 serializer.map(m_plymwinjTables);
213 serializer.map(m_skprwatTables);
214 serializer.map(m_skprpolyTables);
215 m_tabdims.serializeOp(serializer);
216 m_regdims.serializeOp(serializer);
217 m_eqldims.serializeOp(serializer);
218 m_aqudims.serializeOp(serializer);
219 serializer(hasImptvd);
220 serializer(hasEnptvd);
221 serializer(hasEqlnum);
222 serializer(hasShrate);
223 serializer(jfunc);
224 oilDenT.serializeOp(serializer);
225 gasDenT.serializeOp(serializer);
226 watDenT.serializeOp(serializer);
227 stcond.serializeOp(serializer);
228 serializer(m_gas_comp_index);
229 serializer(m_rtemp);
230 if (!serializer.isSerializing()) {
231 m_simpleTables = simpleTables;
232 if (split.plyshMax > 0) {
233 TableContainer container(split.plyshMax);
234 for (const auto& it : split.plyshMap) {
235 container.addTable(it.first, it.second);
236 }
237 m_simpleTables.insert(std::make_pair("PLYSHLOG", container));
238 }
239 if (split.rockMax > 0) {
240 TableContainer container(split.rockMax);
241 for (const auto& it : split.rockMap) {
242 container.addTable(it.first, it.second);
243 }
244 m_simpleTables.insert(std::make_pair("ROCKTAB", container));
245 }
246 }
247 }
248
249 private:
250 TableContainer& forceGetTables( const std::string& tableName , size_t numTables);
251
252 void complainAboutAmbiguousKeyword(const Deck& deck, const std::string& keywordName);
253
254 void addTables( const std::string& tableName , size_t numTables);
255 void initSimpleTables(const Deck& deck);
256 void initRTempTables(const Deck& deck);
257 void initDims(const Deck& deck);
258 void initRocktabTables(const Deck& deck);
259 void initGasvisctTables(const Deck& deck);
260
261 void initPlymaxTables(const Deck& deck);
262 void initPlyrockTables(const Deck& deck);
263 void initPlyshlogTables(const Deck& deck);
264
265 void initPlymwinjTables(const Deck& deck);
266 void initSkprwatTables(const Deck& deck);
267 void initSkprpolyTables(const Deck& deck);
268
269 //void initRockTables(const Deck& deck, const std::string& keywordName);
270
271 template <class TableType>
272 void initRockTables(const Deck& deck, const std::string& keywordName, std::vector<TableType>& rocktable ) {
273
274 if (!deck.hasKeyword(keywordName))
275 return;
276
277 if (!deck.hasKeyword("ROCKCOMP")) {
278 OpmLog::error("ROCKCOMP must be present if ROCK2DTR is used");
279 }
280
281 if (!deck.hasKeyword("ROCKWNOD")) {
282 OpmLog::error("ROCKWNOD must be present if ROCK2DTR is used");
283 }
284
285 const auto& rockcompKeyword = deck.getKeyword("ROCKCOMP");
286 const auto& record = rockcompKeyword.getRecord( 0 );
287 size_t numTables = record.getItem("NTROCC").get< int >(0);
288 rocktable.resize(numTables);
289
290 const auto& keyword = deck.getKeyword(keywordName);
291 size_t numEntries = keyword.size();
292 size_t regionIdx = 0;
293 size_t tableIdx = 0;
294 for (unsigned lineIdx = 0; lineIdx < numEntries; ++lineIdx) {
295 if (keyword.getRecord(lineIdx).getItem("PRESSURE").hasValue(0)) {
296 rocktable[regionIdx].init(keyword.getRecord(lineIdx), tableIdx);
297 tableIdx++;
298 } else { // next region
299 tableIdx = 0;
300 regionIdx++;
301 }
302 }
303 assert(regionIdx == numTables - 1 );
304 }
305
306
307 template <class TableType>
308 void initPvtwsaltTables(const Deck& deck, std::vector<TableType>& pvtwtables ) {
309
310 size_t numTables = m_tabdims.getNumPVTTables();
311 pvtwtables.resize(numTables);
312
313 const auto& keyword = deck.getKeyword("PVTWSALT");
314 size_t numEntries = keyword.size();
315 size_t regionIdx = 0;
316 for (unsigned lineIdx = 0; lineIdx < numEntries; lineIdx += 2) {
317 pvtwtables[regionIdx].init(keyword.getRecord(lineIdx), keyword.getRecord(lineIdx+1));
318 ++regionIdx;
319 }
320 assert(regionIdx == numTables);
321 }
322
323 template <class TableType>
324 void initBrineTables(const Deck& deck, std::vector<TableType>& brinetables ) {
325
326 size_t numTables = m_tabdims.getNumPVTTables();
327 brinetables.resize(numTables);
328
329 const auto& keyword = deck.getKeyword("BDENSITY");
330 size_t numEntries = keyword.size();
331 assert(numEntries == numTables);
332 for (unsigned lineIdx = 0; lineIdx < numEntries; ++lineIdx) {
333 brinetables[lineIdx].init(keyword.getRecord(lineIdx));
334 }
335 }
336
337
338 void initSolventTables(const Deck& deck, std::vector<SolventDensityTable>& solventtables);
339
343 template <class TableType>
344 void initSimpleTableContainerWithJFunc(const Deck& deck,
345 const std::string& keywordName,
346 const std::string& tableName,
347 size_t numTables) {
348 if (!deck.hasKeyword(keywordName))
349 return; // the table is not featured by the deck...
350
351 auto& container = forceGetTables(tableName , numTables);
352
353 if (deck.count(keywordName) > 1) {
354 complainAboutAmbiguousKeyword(deck, keywordName);
355 return;
356 }
357
358 const auto& tableKeyword = deck.getKeyword(keywordName);
359 for (size_t tableIdx = 0; tableIdx < tableKeyword.size(); ++tableIdx) {
360 const auto& dataItem = tableKeyword.getRecord( tableIdx ).getItem("DATA");
361 if (dataItem.data_size() > 0) {
362 std::shared_ptr<TableType> table = std::make_shared<TableType>( dataItem, useJFunc() );
363 container.addTable( tableIdx , table );
364 }
365 }
366 }
367
368
369
370 template <class TableType>
371 void initSimpleTableContainer(const Deck& deck,
372 const std::string& keywordName,
373 const std::string& tableName,
374 size_t numTables) {
375 if (!deck.hasKeyword(keywordName))
376 return; // the table is not featured by the deck...
377
378 auto& container = forceGetTables(tableName , numTables);
379
380 if (deck.count(keywordName) > 1) {
381 complainAboutAmbiguousKeyword(deck, keywordName);
382 return;
383 }
384
385 const auto& tableKeyword = deck.getKeyword(keywordName);
386 for (size_t tableIdx = 0; tableIdx < tableKeyword.size(); ++tableIdx) {
387 const auto& dataItem = tableKeyword.getRecord( tableIdx ).getItem("DATA");
388 if (dataItem.data_size() > 0) {
389 std::shared_ptr<TableType> table = std::make_shared<TableType>( dataItem );
390 container.addTable( tableIdx , table );
391 }
392 }
393 }
394
395 template <class TableType>
396 void initSimpleTableContainer(const Deck& deck,
397 const std::string& keywordName,
398 size_t numTables) {
399 initSimpleTableContainer<TableType>(deck , keywordName , keywordName , numTables);
400 }
401
402
403 template <class TableType>
404 void initSimpleTableContainerWithJFunc(const Deck& deck,
405 const std::string& keywordName,
406 size_t numTables) {
407 initSimpleTableContainerWithJFunc<TableType>(deck , keywordName , keywordName , numTables);
408 }
409
410
411
412 template <class TableType>
413 void initSimpleTable(const Deck& deck,
414 const std::string& keywordName,
415 std::vector<TableType>& tableVector) {
416 if (!deck.hasKeyword(keywordName))
417 return; // the table is not featured by the deck...
418
419 if (deck.count(keywordName) > 1) {
420 complainAboutAmbiguousKeyword(deck, keywordName);
421 return;
422 }
423
424 const auto& tableKeyword = deck.getKeyword(keywordName);
425 for (size_t tableIdx = 0; tableIdx < tableKeyword.size(); ++tableIdx) {
426 const auto& dataItem = tableKeyword.getRecord( tableIdx ).getItem("DATA");
427 if (dataItem.data_size() == 0) {
428 // for simple tables, an empty record indicates that the previous table
429 // should be copied...
430 if (tableIdx == 0) {
431 std::string msg = "The first table for keyword "+keywordName+" must be explicitly defined! Ignoring keyword";
432 const auto& location = tableKeyword.location();
433 OpmLog::warning(Log::fileMessage(location, msg));
434 return;
435 }
436 tableVector.push_back(tableVector.back());
437 continue;
438 }
439
440 tableVector.push_back(TableType());
441 tableVector[tableIdx].init(dataItem);
442 }
443 }
444
445
446 template <class TableType>
447 void initFullTables(const Deck& deck,
448 const std::string& keywordName,
449 std::vector<TableType>& tableVector) {
450 if (!deck.hasKeyword(keywordName))
451 return; // the table is not featured by the deck...
452
453 if (deck.count(keywordName) > 1) {
454 complainAboutAmbiguousKeyword(deck, keywordName);
455 return;
456 }
457
458 const auto& tableKeyword = deck.getKeyword(keywordName);
459
460 int numTables = TableType::numTables( tableKeyword );
461 for (int tableIdx = 0; tableIdx < numTables; ++tableIdx)
462 tableVector.emplace_back( tableKeyword , tableIdx );
463 }
464
465 std::map<std::string , TableContainer> m_simpleTables;
466 std::vector<PvtgTable> m_pvtgTables;
467 std::vector<PvtoTable> m_pvtoTables;
468 std::vector<Rock2dTable> m_rock2dTables;
469 std::vector<Rock2dtrTable> m_rock2dtrTables;
470 PvtwTable m_pvtwTable;
471 PvcdoTable m_pvcdoTable;
472 DensityTable m_densityTable;
473 PlyvmhTable m_plyvmhTable;
474 RockTable m_rockTable;
475 PlmixparTable m_plmixparTable;
476 ShrateTable m_shrateTable;
477 Stone1exTable m_stone1exTable;
478 TlmixparTable m_tlmixparTable;
479 ViscrefTable m_viscrefTable;
480 WatdentTable m_watdentTable;
481 std::vector<PvtwsaltTable> m_pvtwsaltTables;
482 std::vector<BrineDensityTable> m_bdensityTables;
483 std::vector<SolventDensityTable> m_sdensityTables;
484 std::map<int, PlymwinjTable> m_plymwinjTables;
485 std::map<int, SkprwatTable> m_skprwatTables;
486 std::map<int, SkprpolyTable> m_skprpolyTables;
487
488 Tabdims m_tabdims;
489 Regdims m_regdims;
490 Eqldims m_eqldims;
491 Aqudims m_aqudims;
492
493 bool hasImptvd = false;// if deck has keyword IMPTVD
494 bool hasEnptvd = false;// if deck has keyword ENPTVD
495 bool hasEqlnum = false;// if deck has keyword EQLNUM
496 bool hasShrate = false;// if deck has keyword SHRATE
497 std::shared_ptr<JFunc> jfunc;
498
499 DenT oilDenT;
500 DenT gasDenT;
501 DenT watDenT;
502 StandardCond stcond;
503 std::size_t m_gas_comp_index;
504 double m_rtemp;
505
506 struct SplitSimpleTables {
507 size_t plyshMax = 0;
508 size_t rockMax = 0;
509 std::map<size_t, std::shared_ptr<PlyshlogTable>> plyshMap;
510 std::map<size_t, std::shared_ptr<RocktabTable>> rockMap;
511 };
512
513 SplitSimpleTables splitSimpleTable(std::map<std::string,TableContainer>& simpleTables);
514 };
515}
516
517
518#endif
519
const char *const string
Definition: cJSON.h:170
Definition: Aqudims.hpp:34
void serializeOp(Serializer &serializer)
Definition: Aqudims.hpp:133
T get(size_t index) const
const DeckRecord & getRecord(size_t index) const
Definition: Deck.hpp:115
bool hasKeyword(const DeckKeyword &keyword) const
DeckKeyword & getKeyword(size_t)
DeckItem & getItem(size_t index)
Definition: DenT.hpp:30
void serializeOp(Serializer &serializer)
Definition: DenT.hpp:62
Definition: Eqldims.hpp:31
void serializeOp(Serializer &serializer)
Definition: Eqldims.hpp:92
Definition: JFunc.hpp:29
static void error(const std::string &message)
static void warning(const std::string &message)
Definition: Regdims.hpp:31
void serializeOp(Serializer &serializer)
Definition: Regdims.hpp:90
Definition: Serializer.hpp:38
Definition: parser/eclipse/EclipseState/Tables/tabdims.hpp:34
size_t getNumPVTTables() const
Definition: parser/eclipse/EclipseState/Tables/tabdims.hpp:83
void serializeOp(Serializer &serializer)
Definition: parser/eclipse/EclipseState/Tables/tabdims.hpp:113
Definition: TableContainer.hpp:34
void addTable(size_t tableNumber, std::shared_ptr< SimpleTable > table)
Definition: TableManager.hpp:63
const TableContainer & getRsvdTables() const
const TableContainer & getTlpmixpaTables() const
const TableContainer & getRocktabTables() const
const Regdims & getRegdims() const
const std::vector< Rock2dTable > & getRock2dTables() const
TableManager()=default
const TableContainer & getEnptvdTables() const
const PlmixparTable & getPlmixparTable() const
const TableContainer & getMsfnTables() const
const TableContainer & getWatvisctTables() const
const TableContainer & getPlyrockTables() const
const std::vector< Rock2dtrTable > & getRock2dtrTables() const
const TableContainer & getMiscTables() const
const TableContainer & getRockwnodTables() const
const Aqudims & getAqudims() const
bool operator==(const TableManager &data) const
const TableContainer & getImkrvdTables() const
const TableContainer & getSlgofTables() const
const std::map< int, PlymwinjTable > & getPlymwinjTables() const
const TableContainer & getTables(const std::string &tableName) const
bool useEqlnum() const
deck has keyword "EQLNUM" — Equilibriation region numbers
const TableContainer & getSgwfnTables() const
const TableContainer & getGasvisctTables() const
const Eqldims & getEqldims() const
size_t numFIPRegions() const
const std::vector< SolventDensityTable > & getSolventDensityTables() const
void serializeOp(Serializer &serializer)
Definition: TableManager.hpp:185
const TableContainer & getSpecrockTables() const
const ShrateTable & getShrateTable() const
const std::map< int, SkprpolyTable > & getSkprpolyTables() const
const std::vector< PvtoTable > & getPvtoTables() const
const TableContainer & getSwfnTables() const
const TableContainer & getSwofTables() const
const TableContainer & getPmiscTables() const
const TableContainer & operator[](const std::string &tableName) const
const TableContainer & getSgofTables() const
const TableContainer & getFoammobTables() const
bool useJFunc() const
deck has keyword "JFUNC" — Use Leverett's J Function for capillary pressure
const DenT & WatDenT() const
const PlyvmhTable & getPlyvmhTable() const
const WatdentTable & getWatdentTable() const
TableManager & operator=(const TableManager &data)
const TableContainer & getPvdgTables() const
const ViscrefTable & getViscrefTable() const
const TableContainer & getSorwmisTables() const
bool useEnptvd() const
deck has keyword "ENPTVD" — Saturation end-point versus depth tables
const TableContainer & getPlydhflfTables() const
const TableContainer & getEnkrvdTables() const
const std::map< std::string, TableContainer > & getSimpleTables() const
bool hasTables(const std::string &tableName) const
const TableContainer & getPdvdTables() const
TableManager(const Deck &deck)
const TableContainer & getFoamadsTables() const
std::size_t gas_comp_index() const
const TableContainer & getPlyshlogTables() const
const TableContainer & getSgcwmisTables() const
bool useImptvd() const
deck has keyword "IMPTVD" — Imbition end-point versus depth tables
const TlmixparTable & getTlmixparTable() const
const DenT & OilDenT() const
const PvcdoTable & getPvcdoTable() const
const Tabdims & getTabdims() const
bool useShrate() const
deck has keyword "SHRATE"
const TableContainer & getPvdsTables() const
const JFunc & getJFunc() const
const TableContainer & getImptvdTables() const
const TableContainer & getPbvdTables() const
const TableContainer & getSaltvdTables() const
const TableContainer & getPlyviscTables() const
const RockTable & getRockTable() const
double rtemp() const
const TableContainer & getPlyadsTables() const
const StandardCond & stCond() const
const TableContainer & getSof3Tables() const
const TableContainer & getAqutabTables() const
const TableContainer & getRvvdTables() const
const Stone1exTable & getStone1exTable() const
const PvtwTable & getPvtwTable() const
const std::map< int, SkprwatTable > & getSkprwatTables() const
const TableContainer & getPvdoTables() const
const std::vector< BrineDensityTable > & getBrineDensityTables() const
const std::vector< PvtwsaltTable > & getPvtwSaltTables() const
const std::vector< PvtgTable > & getPvtgTables() const
const TableContainer & getOilvisctTables() const
const TableContainer & getOverburdTables() const
const TableContainer & getRtempvdTables() const
const TableContainer & getSgfnTables() const
const TableContainer & getSof2Tables() const
static TableManager serializeObject()
const DenT & GasDenT() const
const TableContainer & getSsfnTables() const
const DensityTable & getDensityTable() const
const TableContainer & getSpecheatTables() const
const TableContainer & getPlymaxTables() const
std::string fileMessage(const Location &location, const std::string &msg)
UDAKeyword keyword(UDAControl control)
Definition: A.hpp:4
static std::string data()
Definition: exprtk.hpp:40022
Definition: FlatTable.hpp:45
void serializeOp(Serializer &serializer)
Definition: FlatTable.hpp:17
Definition: FlatTable.hpp:172
Definition: FlatTable.hpp:206
Definition: FlatTable.hpp:147
Definition: FlatTable.hpp:82
Definition: FlatTable.hpp:110
Definition: FlatTable.hpp:231
Definition: StandardCond.hpp:24
void serializeOp(Serializer &serializer)
Definition: StandardCond.hpp:35
Definition: FlatTable.hpp:256
Definition: FlatTable.hpp:284
Definition: FlatTable.hpp:312
Definition: FlatTable.hpp:343