SwfnTable.hpp
Go to the documentation of this file.
1 /*
2  Copyright (C) 2014 IRIS AS
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 #ifndef OPM_PARSER_SWFN_TABLE_HPP
20 #define OPM_PARSER_SWFN_TABLE_HPP
21 
22 #include "SimpleTable.hpp"
23 
24 namespace Opm {
25  // forward declaration
26  class TableManager;
27 
28  class SwfnTable : public SimpleTable {
29 
30 
31  friend class TableManager;
32 
37  void init(Opm::DeckItemConstPtr item)
38  {
39  SimpleTable::init(item,
40  std::vector<std::string>{"SW", "KRW", "PCOW"});
41 
43  SimpleTable::checkMonotonic("SW", /*isAscending=*/true);
46  SimpleTable::checkMonotonic("KRW", /*isAscending=*/true, /*strict=*/false);
47  SimpleTable::checkMonotonic("PCOW", /*isAscending=*/false, /*strict=*/false);
48  }
49 
50  public:
51  SwfnTable() = default;
52 
53 #ifdef BOOST_TEST_MODULE
54  // DO NOT TRY TO CALL THIS METHOD! it is only for the unit tests!
55  void initFORUNITTESTONLY(Opm::DeckItemConstPtr item)
56  { init(item); }
57 #endif
58 
63 
64  const std::vector<double> &getSwColumn() const
65  { return SimpleTable::getColumn(0); }
66 
67  const std::vector<double> &getKrwColumn() const
68  { return SimpleTable::getColumn(1); }
69 
70  // this column is p_o - p_w (non-wetting phase pressure minus
71  // wetting phase pressure for a given water saturation)
72  const std::vector<double> &getPcowColumn() const
73  { return SimpleTable::getColumn(2); }
74  };
75 }
76 
77 #endif
78 
size_t numColumns() const
void init(Opm::DeckItemConstPtr deckItem, const std::vector< std::string > &columnNames)
Read simple tables from keywords like SWOF.
const std::vector< double > & getPcowColumn() const
Definition: SwfnTable.hpp:72
Definition: Deck.hpp:29
const std::vector< double > & getColumn(const std::string &name) const
void checkNonDefaultable(const std::string &columnName)
const std::vector< double > & getKrwColumn() const
Definition: SwfnTable.hpp:67
std::shared_ptr< const DeckItem > DeckItemConstPtr
Definition: DeckItem.hpp:127
Definition: TableManager.hpp:66
const std::vector< double > & getSwColumn() const
Definition: SwfnTable.hpp:64
size_t numRows() const
static size_t numTables(Opm::DeckKeywordConstPtr keyword)
Returns the number of tables in a keyword.
double evaluate(const std::string &columnName, double xPos) const
Evaluate a column of the table at a given position.
void checkMonotonic(const std::string &columnName, bool isAscending, bool isStrictlyMonotonic=true)
Definition: SimpleTable.hpp:32
SwfnTable()=default
void applyDefaultsLinear(const std::string &columnName)
Definition: SwfnTable.hpp:28