My Project
PiecewiseLinearTwoPhaseMaterial.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_HPP
28#define OPM_PIECEWISE_LINEAR_TWO_PHASE_MATERIAL_HPP
29
31
33
34#include <stdexcept>
35#include <type_traits>
36
37namespace Opm {
48template <class TraitsT, class ParamsT = PiecewiseLinearTwoPhaseMaterialParams<TraitsT> >
49class PiecewiseLinearTwoPhaseMaterial : public TraitsT
50{
51 using ValueVector = typename ParamsT::ValueVector;
52
53public:
55 using Traits = TraitsT;
56
58 using Params = ParamsT;
59
61 using Scalar = typename Traits::Scalar;
62
64 static constexpr int numPhases = Traits::numPhases;
65 static_assert(numPhases == 2,
66 "The piecewise linear two-phase capillary pressure law only"
67 "applies to the case of two fluid phases");
68
71 static constexpr bool implementsTwoPhaseApi = true;
72
75 static constexpr bool implementsTwoPhaseSatApi = true;
76
79 static constexpr bool isSaturationDependent = true;
80
83 static constexpr bool isPressureDependent = false;
84
87 static constexpr bool isTemperatureDependent = false;
88
91 static constexpr bool isCompositionDependent = false;
92
96 template <class Container, class FluidState>
97 static void capillaryPressures(Container& values, const Params& params, const FluidState& fs)
98 {
99 using Evaluation = typename std::remove_reference<decltype(values[0])>::type;
100
101 values[Traits::wettingPhaseIdx] = 0.0; // reference phase
102 values[Traits::nonWettingPhaseIdx] = pcnw<FluidState, Evaluation>(params, fs);
103 }
104
109 template <class Container, class FluidState>
110 static void saturations(Container& /* values */, const Params& /* params */, const FluidState& /* fs */)
111 { throw std::logic_error("Not implemented: saturations()"); }
112
116 template <class Container, class FluidState>
117 static void relativePermeabilities(Container& values, const Params& params, const FluidState& fs)
118 {
119 using Evaluation = typename std::remove_reference<decltype(values[0])>::type;
120
121 values[Traits::wettingPhaseIdx] = krw<FluidState, Evaluation>(params, fs);
122 values[Traits::nonWettingPhaseIdx] = krn<FluidState, Evaluation>(params, fs);
123 }
124
128 template <class FluidState, class Evaluation = typename FluidState::Scalar>
129 static Evaluation pcnw(const Params& params, const FluidState& fs)
130 {
131 const auto& Sw =
132 decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
133
134 return twoPhaseSatPcnw(params, Sw);
135 }
136
140 template <class Evaluation>
141 static Evaluation twoPhaseSatPcnw(const Params& params, const Evaluation& Sw)
142 { return eval_(params.SwPcwnSamples(), params.pcnwSamples(), Sw); }
143
144 template <class Evaluation>
145 static Evaluation twoPhaseSatPcnwInv(const Params& params, const Evaluation& pcnw)
146 { return eval_(params.pcnwSamples(), params.SwPcwnSamples(), pcnw); }
147
151 template <class FluidState, class Evaluation = typename FluidState::Scalar>
152 static Evaluation Sw(const Params& /* params */, const FluidState& /* fs */)
153 { throw std::logic_error("Not implemented: Sw()"); }
154
155 template <class Evaluation>
156 static Evaluation twoPhaseSatSw(const Params& /* params */, const Evaluation& /* pC */)
157 { throw std::logic_error("Not implemented: twoPhaseSatSw()"); }
158
163 template <class FluidState, class Evaluation = typename FluidState::Scalar>
164 static Evaluation Sn(const Params& params, const FluidState& fs)
165 { return 1 - Sw<FluidState, Scalar>(params, fs); }
166
167 template <class Evaluation>
168 static Evaluation twoPhaseSatSn(const Params& params, const Evaluation& pC)
169 { return 1 - twoPhaseSatSw(params, pC); }
170
175 template <class FluidState, class Evaluation = typename FluidState::Scalar>
176 static Evaluation krw(const Params& params, const FluidState& fs)
177 {
178 const auto& Sw =
179 decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
180
181 return twoPhaseSatKrw(params, Sw);
182 }
183
184 template <class Evaluation>
185 static Evaluation twoPhaseSatKrw(const Params& params, const Evaluation& Sw)
186 { return eval_(params.SwKrwSamples(), params.krwSamples(), Sw); }
187
188 template <class Evaluation>
189 static Evaluation twoPhaseSatKrwInv(const Params& params, const Evaluation& krw)
190 { return eval_(params.krwSamples(), params.SwKrwSamples(), krw); }
191
196 template <class FluidState, class Evaluation = typename FluidState::Scalar>
197 static Evaluation krn(const Params& params, const FluidState& fs)
198 {
199 const auto& Sw =
200 decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
201
202 return twoPhaseSatKrn(params, Sw);
203 }
204
205 template <class Evaluation>
206 static Evaluation twoPhaseSatKrn(const Params& params, const Evaluation& Sw)
207 { return eval_(params.SwKrnSamples(), params.krnSamples(), Sw); }
208
209 template <class Evaluation>
210 static Evaluation twoPhaseSatKrnInv(const Params& params, const Evaluation& krn)
211 { return eval_(params.krnSamples(), params.SwKrnSamples(), krn); }
212
213 template <class Evaluation>
214 static size_t findSegmentIndex(const ValueVector& xValues, const Evaluation& x){
215 return findSegmentIndex_(xValues, scalarValue(x));
216 }
217
218 template <class Evaluation>
219 static size_t findSegmentIndexDescending(const ValueVector& xValues, const Evaluation& x){
220 return findSegmentIndexDescending_(xValues, scalarValue(x));
221 }
222
223 template <class Evaluation>
224 static Evaluation eval(const ValueVector& xValues, const ValueVector& yValues, const Evaluation& x, unsigned segIdx){
225 Scalar x0 = xValues[segIdx];
226 Scalar x1 = xValues[segIdx + 1];
227
228 Scalar y0 = yValues[segIdx];
229 Scalar y1 = yValues[segIdx + 1];
230
231 Scalar m = (y1 - y0)/(x1 - x0);
232
233 return y0 + (x - x0)*m;
234 }
235
236private:
237 template <class Evaluation>
238 static Evaluation eval_(const ValueVector& xValues,
239 const ValueVector& yValues,
240 const Evaluation& x)
241 {
242 if (xValues.front() < xValues.back())
243 return evalAscending_(xValues, yValues, x);
244 return evalDescending_(xValues, yValues, x);
245 }
246
247 template <class Evaluation>
248 static Evaluation evalAscending_(const ValueVector& xValues,
249 const ValueVector& yValues,
250 const Evaluation& x)
251 {
252 if (x <= xValues.front())
253 return yValues.front();
254 if (x >= xValues.back())
255 return yValues.back();
256
257 size_t segIdx = findSegmentIndex_(xValues, scalarValue(x));
258
259 return eval(xValues, yValues, x, segIdx);
260 }
261
262 template <class Evaluation>
263 static Evaluation evalDescending_(const ValueVector& xValues,
264 const ValueVector& yValues,
265 const Evaluation& x)
266 {
267 if (x >= xValues.front())
268 return yValues.front();
269 if (x <= xValues.back())
270 return yValues.back();
271
272 size_t segIdx = findSegmentIndexDescending_(xValues, scalarValue(x));
273
274 return eval(xValues, yValues, x, segIdx);
275 }
276
277 template <class Evaluation>
278 static Evaluation evalDeriv_(const ValueVector& xValues,
279 const ValueVector& yValues,
280 const Evaluation& x)
281 {
282 if (x <= xValues.front())
283 return 0.0;
284 if (x >= xValues.back())
285 return 0.0;
286
287 size_t segIdx = findSegmentIndex_(xValues, scalarValue(x));
288
289 Scalar x0 = xValues[segIdx];
290 Scalar x1 = xValues[segIdx + 1];
291
292 Scalar y0 = yValues[segIdx];
293 Scalar y1 = yValues[segIdx + 1];
294
295 return (y1 - y0)/(x1 - x0);
296 }
297 template<class ScalarT>
298 static size_t findSegmentIndex_(const ValueVector& xValues, const ScalarT& x)
299 {
300 assert(xValues.size() > 1); // we need at least two sampling points!
301 size_t n = xValues.size() - 1;
302 if (xValues.back() <= x)
303 return n - 1;
304 else if (x <= xValues.front())
305 return 0;
306
307 // bisection
308 size_t lowIdx = 0, highIdx = n;
309 while (lowIdx + 1 < highIdx) {
310 size_t curIdx = (lowIdx + highIdx)/2;
311 if (xValues[curIdx] < x)
312 lowIdx = curIdx;
313 else
314 highIdx = curIdx;
315 }
316
317 return lowIdx;
318 }
319
320 static size_t findSegmentIndexDescending_(const ValueVector& xValues, Scalar x)
321 {
322 assert(xValues.size() > 1); // we need at least two sampling points!
323 size_t n = xValues.size() - 1;
324 if (x <= xValues.back())
325 return n;
326 else if (xValues.front() <= x)
327 return 0;
328
329 // bisection
330 size_t lowIdx = 0, highIdx = n;
331 while (lowIdx + 1 < highIdx) {
332 size_t curIdx = (lowIdx + highIdx)/2;
333 if (xValues[curIdx] >= x)
334 lowIdx = curIdx;
335 else
336 highIdx = curIdx;
337 }
338
339 return lowIdx;
340 }
341};
342
343} // namespace Opm
344
345#endif
A traits class which provides basic mathematical functions for arbitrary scalar floating point values...
Implementation of a tabulated, piecewise linear capillary pressure law.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:50
TraitsT Traits
The traits class for this material law.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:55
ParamsT Params
The type of the parameter objects for this law.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:58
static constexpr int numPhases
The number of fluid phases.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:64
static Evaluation krn(const Params &params, const FluidState &fs)
The relative permeability for the non-wetting phase of the porous medium.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:197
static Evaluation twoPhaseSatPcnw(const Params &params, const Evaluation &Sw)
The saturation-capillary pressure curve.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:141
static Evaluation krw(const Params &params, const FluidState &fs)
The relative permeability for the wetting phase of the porous medium.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:176
static constexpr bool implementsTwoPhaseSatApi
Specify whether this material law implements the two-phase convenience API which only depends on the ...
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:75
static void relativePermeabilities(Container &values, const Params &params, const FluidState &fs)
The relative permeabilities.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:117
typename Traits::Scalar Scalar
The type of the scalar values for this law.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:61
static Evaluation Sn(const Params &params, const FluidState &fs)
Calculate the non-wetting phase saturations depending on the phase pressures.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:164
static constexpr bool isCompositionDependent
Specify whether the quantities defined by this material law are dependent on the phase composition.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:91
static constexpr bool isSaturationDependent
Specify whether the quantities defined by this material law are saturation dependent.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:79
static void saturations(Container &, const Params &, const FluidState &)
The saturations of the fluid phases starting from their pressure differences.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:110
static constexpr bool isTemperatureDependent
Specify whether the quantities defined by this material law are temperature dependent.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:87
static constexpr bool isPressureDependent
Specify whether the quantities defined by this material law are dependent on the absolute pressure.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:83
static constexpr bool implementsTwoPhaseApi
Specify whether this material law implements the two-phase convenience API.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:71
static Evaluation pcnw(const Params &params, const FluidState &fs)
The capillary pressure-saturation curve.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:129
static void capillaryPressures(Container &values, const Params &params, const FluidState &fs)
The capillary pressure-saturation curve.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:97
static Evaluation Sw(const Params &, const FluidState &)
The saturation-capillary pressure curve.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:152
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition: Exceptions.hpp:30