My Project
Loading...
Searching...
No Matches
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#include <opm/common/TimingMacros.hpp>
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 OPM_TIMEFUNCTION_LOCAL();
100 using Evaluation = typename std::remove_reference<decltype(values[0])>::type;
101
102 values[Traits::wettingPhaseIdx] = 0.0; // reference phase
103 values[Traits::nonWettingPhaseIdx] = pcnw<FluidState, Evaluation>(params, fs);
104 }
105
110 template <class Container, class FluidState>
111 static void saturations(Container& /* values */, const Params& /* params */, const FluidState& /* fs */)
112 { throw std::logic_error("Not implemented: saturations()"); }
113
117 template <class Container, class FluidState>
118 static void relativePermeabilities(Container& values, const Params& params, const FluidState& fs)
119 {
120 OPM_TIMEFUNCTION_LOCAL();
121 using Evaluation = typename std::remove_reference<decltype(values[0])>::type;
122
123 values[Traits::wettingPhaseIdx] = krw<FluidState, Evaluation>(params, fs);
124 values[Traits::nonWettingPhaseIdx] = krn<FluidState, Evaluation>(params, fs);
125 }
126
130 template <class FluidState, class Evaluation = typename FluidState::Scalar>
131 static Evaluation pcnw(const Params& params, const FluidState& fs)
132 {
133 OPM_TIMEFUNCTION_LOCAL();
134 const auto& Sw =
135 decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
136
137 return twoPhaseSatPcnw(params, Sw);
138 }
139
143 template <class Evaluation>
144 static Evaluation twoPhaseSatPcnw(const Params& params, const Evaluation& Sw)
145 {
146 OPM_TIMEFUNCTION_LOCAL();
147 return eval_(params.SwPcwnSamples(), params.pcnwSamples(), Sw);
148 }
149
150 template <class Evaluation>
151 static Evaluation twoPhaseSatPcnwInv(const Params& params, const Evaluation& pcnw)
152 { return eval_(params.pcnwSamples(), params.SwPcwnSamples(), pcnw); }
153
157 template <class FluidState, class Evaluation = typename FluidState::Scalar>
158 static Evaluation Sw(const Params& /* params */, const FluidState& /* fs */)
159 { throw std::logic_error("Not implemented: Sw()"); }
160
161 template <class Evaluation>
162 static Evaluation twoPhaseSatSw(const Params& /* params */, const Evaluation& /* pC */)
163 { throw std::logic_error("Not implemented: twoPhaseSatSw()"); }
164
169 template <class FluidState, class Evaluation = typename FluidState::Scalar>
170 static Evaluation Sn(const Params& params, const FluidState& fs)
171 { return 1 - Sw<FluidState, Scalar>(params, fs); }
172
173 template <class Evaluation>
174 static Evaluation twoPhaseSatSn(const Params& params, const Evaluation& pC)
175 { return 1 - twoPhaseSatSw(params, pC); }
176
181 template <class FluidState, class Evaluation = typename FluidState::Scalar>
182 static Evaluation krw(const Params& params, const FluidState& fs)
183 {
184 OPM_TIMEFUNCTION_LOCAL();
185 const auto& Sw =
186 decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
187
188 return twoPhaseSatKrw(params, Sw);
189 }
190
191 template <class Evaluation>
192 static Evaluation twoPhaseSatKrw(const Params& params, const Evaluation& Sw)
193 {
194 OPM_TIMEFUNCTION_LOCAL();
195 return eval_(params.SwKrwSamples(), params.krwSamples(), Sw);
196 }
197
198 template <class Evaluation>
199 static Evaluation twoPhaseSatKrwInv(const Params& params, const Evaluation& krw)
200 {
201 OPM_TIMEFUNCTION_LOCAL();
202 return eval_(params.krwSamples(), params.SwKrwSamples(), krw);
203 }
204
209 template <class FluidState, class Evaluation = typename FluidState::Scalar>
210 static Evaluation krn(const Params& params, const FluidState& fs)
211 {
212 OPM_TIMEFUNCTION_LOCAL();
213 const auto& Sw =
214 decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
215
216 return twoPhaseSatKrn(params, Sw);
217 }
218
219 template <class Evaluation>
220 static Evaluation twoPhaseSatKrn(const Params& params, const Evaluation& Sw)
221 {
222 OPM_TIMEFUNCTION_LOCAL();
223 return eval_(params.SwKrnSamples(), params.krnSamples(), Sw);
224 }
225
226 template <class Evaluation>
227 static Evaluation twoPhaseSatKrnInv(const Params& params, const Evaluation& krn)
228 { return eval_(params.krnSamples(), params.SwKrnSamples(), krn); }
229
230 template <class Evaluation>
231 static size_t findSegmentIndex(const ValueVector& xValues, const Evaluation& x){
232 return findSegmentIndex_(xValues, scalarValue(x));
233 }
234
235 template <class Evaluation>
236 static size_t findSegmentIndexDescending(const ValueVector& xValues, const Evaluation& x){
237 return findSegmentIndexDescending_(xValues, scalarValue(x));
238 }
239
240 template <class Evaluation>
241 static Evaluation eval(const ValueVector& xValues, const ValueVector& yValues, const Evaluation& x, unsigned segIdx){
242 Scalar x0 = xValues[segIdx];
243 Scalar x1 = xValues[segIdx + 1];
244
245 Scalar y0 = yValues[segIdx];
246 Scalar y1 = yValues[segIdx + 1];
247
248 Scalar m = (y1 - y0)/(x1 - x0);
249
250 return y0 + (x - x0)*m;
251 }
252
253private:
254 template <class Evaluation>
255 static Evaluation eval_(const ValueVector& xValues,
256 const ValueVector& yValues,
257 const Evaluation& x)
258 {
259 OPM_TIMEFUNCTION_LOCAL();
260 if (xValues.front() < xValues.back())
261 return evalAscending_(xValues, yValues, x);
262 return evalDescending_(xValues, yValues, x);
263 }
264
265 template <class Evaluation>
266 static Evaluation evalAscending_(const ValueVector& xValues,
267 const ValueVector& yValues,
268 const Evaluation& x)
269 {
270 OPM_TIMEFUNCTION_LOCAL();
271 if (x <= xValues.front())
272 return yValues.front();
273 if (x >= xValues.back())
274 return yValues.back();
275
276 size_t segIdx = findSegmentIndex_(xValues, scalarValue(x));
277
278 return eval(xValues, yValues, x, segIdx);
279 }
280
281 template <class Evaluation>
282 static Evaluation evalDescending_(const ValueVector& xValues,
283 const ValueVector& yValues,
284 const Evaluation& x)
285 {
286 OPM_TIMEFUNCTION_LOCAL();
287 if (x >= xValues.front())
288 return yValues.front();
289 if (x <= xValues.back())
290 return yValues.back();
291
292 size_t segIdx = findSegmentIndexDescending_(xValues, scalarValue(x));
293
294 return eval(xValues, yValues, x, segIdx);
295 }
296
297 template <class Evaluation>
298 static Evaluation evalDeriv_(const ValueVector& xValues,
299 const ValueVector& yValues,
300 const Evaluation& x)
301 {
302 OPM_TIMEFUNCTION_LOCAL();
303 if (x <= xValues.front())
304 return 0.0;
305 if (x >= xValues.back())
306 return 0.0;
307
308 size_t segIdx = findSegmentIndex_(xValues, scalarValue(x));
309
310 Scalar x0 = xValues[segIdx];
311 Scalar x1 = xValues[segIdx + 1];
312
313 Scalar y0 = yValues[segIdx];
314 Scalar y1 = yValues[segIdx + 1];
315
316 return (y1 - y0)/(x1 - x0);
317 }
318 template<class ScalarT>
319 static size_t findSegmentIndex_(const ValueVector& xValues, const ScalarT& x)
320 {
321 OPM_TIMEFUNCTION_LOCAL();
322 assert(xValues.size() > 1); // we need at least two sampling points!
323 size_t n = xValues.size() - 1;
324 if (xValues.back() <= x)
325 return n - 1;
326 else if (x <= xValues.front())
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 static size_t findSegmentIndexDescending_(const ValueVector& xValues, Scalar x)
343 {
344 OPM_TIMEFUNCTION_LOCAL();
345 assert(xValues.size() > 1); // we need at least two sampling points!
346 size_t n = xValues.size() - 1;
347 if (x <= xValues.back())
348 return n;
349 else if (xValues.front() <= x)
350 return 0;
351
352 // bisection
353 size_t lowIdx = 0, highIdx = n;
354 while (lowIdx + 1 < highIdx) {
355 size_t curIdx = (lowIdx + highIdx)/2;
356 if (xValues[curIdx] >= x)
357 lowIdx = curIdx;
358 else
359 highIdx = curIdx;
360 }
361
362 return lowIdx;
363 }
364};
365
366} // namespace Opm
367
368#endif
A traits class which provides basic mathematical functions for arbitrary scalar floating point values...
Specification of the material parameters for a two-phase material law which uses a table and piecewis...
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:210
static Evaluation twoPhaseSatPcnw(const Params &params, const Evaluation &Sw)
The saturation-capillary pressure curve.
Definition PiecewiseLinearTwoPhaseMaterial.hpp:144
static Evaluation krw(const Params &params, const FluidState &fs)
The relative permeability for the wetting phase of the porous medium.
Definition PiecewiseLinearTwoPhaseMaterial.hpp:182
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:118
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:170
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:111
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:131
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:158
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30