PiecewiseLinearTwoPhaseMaterialParams.hpp
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
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 2 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 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
27#ifndef OPM_PIECEWISE_LINEAR_TWO_PHASE_MATERIAL_PARAMS_HPP
28#define OPM_PIECEWISE_LINEAR_TWO_PHASE_MATERIAL_PARAMS_HPP
29
30#include <algorithm>
31#include <cassert>
32#include <vector>
33
35
36namespace Opm {
43template<class TraitsT>
45{
46 using Scalar = typename TraitsT::Scalar;
47
48public:
49 using ValueVector = std::vector<Scalar>;
50
51 using Traits = TraitsT;
52
54 {
55 }
56
61 void finalize()
62 {
64
65 // revert the order of the sampling points if they were given
66 // in reverse direction.
67 if (SwPcwnSamples_.front() > SwPcwnSamples_.back())
68 swapOrder_(SwPcwnSamples_, pcwnSamples_);
69
70 if (SwKrwSamples_.front() > SwKrwSamples_.back())
71 swapOrder_(SwKrwSamples_, krwSamples_);
72
73
74 if (SwKrnSamples_.front() > SwKrnSamples_.back())
75 swapOrder_(SwKrnSamples_, krnSamples_);
76
77 }
78
83 { EnsureFinalized::check(); return SwKrwSamples_; }
84
89 { EnsureFinalized::check(); return SwKrnSamples_; }
90
95 { EnsureFinalized::check(); return SwPcwnSamples_; }
96
103 { EnsureFinalized::check(); return pcwnSamples_; }
104
110 template <class Container>
111 void setPcnwSamples(const Container& SwValues, const Container& values)
112 {
113 assert(SwValues.size() == values.size());
114
115 size_t n = SwValues.size();
116 SwPcwnSamples_.resize(n);
117 pcwnSamples_.resize(n);
118
119 std::copy(SwValues.begin(), SwValues.end(), SwPcwnSamples_.begin());
120 std::copy(values.begin(), values.end(), pcwnSamples_.begin());
121 }
122
129 const ValueVector& krwSamples() const
130 { EnsureFinalized::check(); return krwSamples_; }
131
138 template <class Container>
139 void setKrwSamples(const Container& SwValues, const Container& values)
140 {
141 assert(SwValues.size() == values.size());
142
143 size_t n = SwValues.size();
144 SwKrwSamples_.resize(n);
145 krwSamples_.resize(n);
146
147 std::copy(SwValues.begin(), SwValues.end(), SwKrwSamples_.begin());
148 std::copy(values.begin(), values.end(), krwSamples_.begin());
149 }
150
157 const ValueVector& krnSamples() const
158 { EnsureFinalized::check(); return krnSamples_; }
159
166 template <class Container>
167 void setKrnSamples(const Container& SwValues, const Container& values)
168 {
169 assert(SwValues.size() == values.size());
170
171 size_t n = SwValues.size();
172 SwKrnSamples_.resize(n);
173 krnSamples_.resize(n);
174
175 std::copy(SwValues.begin(), SwValues.end(), SwKrnSamples_.begin());
176 std::copy(values.begin(), values.end(), krnSamples_.begin());
177 }
178
179private:
180 void swapOrder_(ValueVector& swValues, ValueVector& values) const
181 {
182 if (swValues.front() > values.back()) {
183 for (unsigned origSampleIdx = 0;
184 origSampleIdx < swValues.size() / 2;
185 ++ origSampleIdx)
186 {
187 size_t newSampleIdx = swValues.size() - origSampleIdx - 1;
188
189 std::swap(swValues[origSampleIdx], swValues[newSampleIdx]);
190 std::swap(values[origSampleIdx], values[newSampleIdx]);
191 }
192 }
193 }
194
195 ValueVector SwPcwnSamples_;
196 ValueVector SwKrwSamples_;
197 ValueVector SwKrnSamples_;
198 ValueVector pcwnSamples_;
199 ValueVector krwSamples_;
200 ValueVector krnSamples_;
201};
202} // namespace Opm
203
204#endif
Default implementation for asserting finalization of parameter objects.
Definition: EnsureFinalized.hpp:47
void finalize()
Mark the object as finalized.
Definition: EnsureFinalized.hpp:75
void check() const
Definition: EnsureFinalized.hpp:63
Specification of the material parameters for a two-phase material law which uses a table and piecewis...
Definition: PiecewiseLinearTwoPhaseMaterialParams.hpp:45
const ValueVector & SwKrnSamples() const
Return the wetting-phase saturation values of all sampling points.
Definition: PiecewiseLinearTwoPhaseMaterialParams.hpp:88
const ValueVector & SwKrwSamples() const
Return the wetting-phase saturation values of all sampling points.
Definition: PiecewiseLinearTwoPhaseMaterialParams.hpp:82
void finalize()
Calculate all dependent quantities once the independent quantities of the parameter object have been ...
Definition: PiecewiseLinearTwoPhaseMaterialParams.hpp:61
const ValueVector & krwSamples() const
Return the sampling points for the relative permeability curve of the wetting phase.
Definition: PiecewiseLinearTwoPhaseMaterialParams.hpp:129
void setKrwSamples(const Container &SwValues, const Container &values)
Set the sampling points for the relative permeability curve of the wetting phase.
Definition: PiecewiseLinearTwoPhaseMaterialParams.hpp:139
const ValueVector & pcnwSamples() const
Return the sampling points for the capillary pressure curve.
Definition: PiecewiseLinearTwoPhaseMaterialParams.hpp:102
void setKrnSamples(const Container &SwValues, const Container &values)
Set the sampling points for the relative permeability curve of the non-wetting phase.
Definition: PiecewiseLinearTwoPhaseMaterialParams.hpp:167
std::vector< Scalar > ValueVector
Definition: PiecewiseLinearTwoPhaseMaterialParams.hpp:49
PiecewiseLinearTwoPhaseMaterialParams()
Definition: PiecewiseLinearTwoPhaseMaterialParams.hpp:53
void setPcnwSamples(const Container &SwValues, const Container &values)
Set the sampling points for the capillary pressure curve.
Definition: PiecewiseLinearTwoPhaseMaterialParams.hpp:111
const ValueVector & krnSamples() const
Return the sampling points for the relative permeability curve of the non-wetting phase.
Definition: PiecewiseLinearTwoPhaseMaterialParams.hpp:157
const ValueVector & SwPcwnSamples() const
Return the wetting-phase saturation values of all sampling points.
Definition: PiecewiseLinearTwoPhaseMaterialParams.hpp:94
TraitsT Traits
Definition: PiecewiseLinearTwoPhaseMaterialParams.hpp:51
Definition: Air_Mesitylene.hpp:34