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  Copyright (C) 2009-2013 by Andreas Lauser
5  Copyright (C) 2015 by IRIS AS
6 
7  This file is part of the Open Porous Media project (OPM).
8 
9  OPM is free software: you can redistribute it and/or modify
10  it under the terms of the GNU General Public License as published by
11  the Free Software Foundation, either version 2 of the License, or
12  (at your option) any later version.
13 
14  OPM is distributed in the hope that it will be useful,
15  but WITHOUT ANY WARRANTY; without even the implied warranty of
16  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17  GNU General Public License for more details.
18 
19  You should have received a copy of the GNU General Public License
20  along with OPM. If not, see <http://www.gnu.org/licenses/>.
21 */
26 #ifndef OPM_PIECEWISE_LINEAR_TWO_PHASE_MATERIAL_HPP
27 #define OPM_PIECEWISE_LINEAR_TWO_PHASE_MATERIAL_HPP
28 
30 
31 #include <opm/common/ErrorMacros.hpp>
32 #include <opm/common/Exceptions.hpp>
34 
35 #include <algorithm>
36 #include <cmath>
37 #include <cassert>
38 
39 namespace Opm {
50 template <class TraitsT, class ParamsT = PiecewiseLinearTwoPhaseMaterialParams<TraitsT> >
51 class PiecewiseLinearTwoPhaseMaterial : public TraitsT
52 {
53  typedef typename ParamsT::ValueVector ValueVector;
54 
55 public:
57  typedef TraitsT Traits;
58 
60  typedef ParamsT Params;
61 
63  typedef typename Traits::Scalar Scalar;
64 
66  static const int numPhases = Traits::numPhases;
67  static_assert(numPhases == 2,
68  "The piecewise linear two-phase capillary pressure law only"
69  "applies to the case of two fluid phases");
70 
73  static const bool implementsTwoPhaseApi = true;
74 
77  static const bool implementsTwoPhaseSatApi = true;
78 
81  static const bool isSaturationDependent = true;
82 
85  static const bool isPressureDependent = false;
86 
89  static const bool isTemperatureDependent = false;
90 
93  static const bool isCompositionDependent = false;
94 
98  template <class Container, class FluidState>
99  static void capillaryPressures(Container &values, const Params &params, const FluidState &fs)
100  {
101  typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
102 
103  values[Traits::wettingPhaseIdx] = 0.0; // reference phase
104  values[Traits::nonWettingPhaseIdx] = pcnw<FluidState, Evaluation>(params, fs);
105  }
106 
111  template <class Container, class FluidState>
112  static void saturations(Container& /* values */, const Params& /* params */, const FluidState& /* fs */)
113  { OPM_THROW(std::logic_error, "Not implemented: saturations()"); }
114 
118  template <class Container, class FluidState>
119  static void relativePermeabilities(Container &values, const Params &params, const FluidState &fs)
120  {
121  typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
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  typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
134  const auto& Sw =
135  FsToolbox::template toLhs<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  { return eval_(params.SwPcwnSamples(), params.pcnwSamples(), Sw); }
146 
147  template <class Evaluation>
148  static Evaluation twoPhaseSatPcnwInv(const Params &params, const Evaluation& pcnw)
149  { return eval_(params.pcnwSamples(), params.SwPcwnSamples(), pcnw); }
150 
154  template <class FluidState, class Evaluation = typename FluidState::Scalar>
155  static Evaluation Sw(const Params& /* params */, const FluidState& /* fs */)
156  { OPM_THROW(std::logic_error, "Not implemented: Sw()"); }
157 
158  template <class Evaluation>
159  static Evaluation twoPhaseSatSw(const Params& /* params */, const Evaluation& /* pC */)
160  { OPM_THROW(std::logic_error, "Not implemented: twoPhaseSatSw()"); }
161 
166  template <class FluidState, class Evaluation = typename FluidState::Scalar>
167  static Evaluation Sn(const Params &params, const FluidState &fs)
168  { return 1 - Sw<FluidState, Scalar>(params, fs); }
169 
170  template <class Evaluation>
171  static Evaluation twoPhaseSatSn(const Params &params, const Evaluation& pC)
172  { return 1 - twoPhaseSatSw(params, pC); }
173 
178  template <class FluidState, class Evaluation = typename FluidState::Scalar>
179  static Evaluation krw(const Params &params, const FluidState &fs)
180  {
181  typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
182  const auto& Sw =
183  FsToolbox::template toLhs<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
184 
185  return twoPhaseSatKrw(params, Sw);
186  }
187 
188  template <class Evaluation>
189  static Evaluation twoPhaseSatKrw(const Params &params, const Evaluation& Sw)
190  { return eval_(params.SwKrwSamples(), params.krwSamples(), Sw); }
191 
192  template <class Evaluation>
193  static Evaluation twoPhaseSatKrwInv(const Params &params, const Evaluation& krw)
194  { return eval_(params.krwSamples(), params.SwKrwSamples(), krw); }
195 
200  template <class FluidState, class Evaluation = typename FluidState::Scalar>
201  static Evaluation krn(const Params &params, const FluidState &fs)
202  {
203  typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
204  const auto& Sw =
205  FsToolbox::template toLhs<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
206 
207  return twoPhaseSatKrn(params, Sw);
208  }
209 
210  template <class Evaluation>
211  static Evaluation twoPhaseSatKrn(const Params &params, const Evaluation& Sw)
212  { return eval_(params.SwKrnSamples(), params.krnSamples(), Sw); }
213 
214  template <class Evaluation>
215  static Evaluation twoPhaseSatKrnInv(const Params &params, const Evaluation& krn)
216  { return eval_(params.krnSamples(), params.SwKrnSamples(), krn); }
217 
218 private:
219  template <class Evaluation>
220  static Evaluation eval_(const ValueVector &xValues,
221  const ValueVector &yValues,
222  const Evaluation& x)
223  {
224  if (xValues.front() < xValues.back())
225  return evalAscending_(xValues, yValues, x);
226  return evalDescending_(xValues, yValues, x);
227  }
228 
229  template <class Evaluation>
230  static Evaluation evalAscending_(const ValueVector &xValues,
231  const ValueVector &yValues,
232  const Evaluation& x)
233  {
234  typedef MathToolbox<Evaluation> Toolbox;
235 
236  if (x <= xValues.front())
237  return yValues.front();
238  if (x >= xValues.back())
239  return yValues.back();
240 
241  size_t segIdx = findSegmentIndex_(xValues, Toolbox::value(x));
242 
243  Scalar x0 = xValues[segIdx];
244  Scalar x1 = xValues[segIdx + 1];
245 
246  Scalar y0 = yValues[segIdx];
247  Scalar y1 = yValues[segIdx + 1];
248 
249  Scalar m = (y1 - y0)/(x1 - x0);
250 
251  return y0 + (x - x0)*m;
252  }
253 
254  template <class Evaluation>
255  static Evaluation evalDescending_(const ValueVector &xValues,
256  const ValueVector &yValues,
257  const Evaluation& x)
258  {
259  typedef MathToolbox<Evaluation> Toolbox;
260 
261  if (x >= xValues.front())
262  return yValues.front();
263  if (x <= xValues.back())
264  return yValues.back();
265 
266  size_t segIdx = findSegmentIndexDescending_(xValues, Toolbox::value(x));
267 
268  Scalar x0 = xValues[segIdx];
269  Scalar x1 = xValues[segIdx + 1];
270 
271  Scalar y0 = yValues[segIdx];
272  Scalar y1 = yValues[segIdx + 1];
273 
274  Scalar m = (y1 - y0)/(x1 - x0);
275 
276  return y0 + (x - x0)*m;
277  }
278 
279  template <class Evaluation>
280  static Evaluation evalDeriv_(const ValueVector &xValues,
281  const ValueVector &yValues,
282  const Evaluation& x)
283  {
284  typedef MathToolbox<Evaluation> Toolbox;
285 
286  if (Toolbox::value(x) <= xValues.front())
287  return Toolbox::createConstant(0.0);
288  if (Toolbox::value(x) >= xValues.back())
289  return Toolbox::createConstant(0.0);
290 
291  size_t segIdx = findSegmentIndex_(xValues, Toolbox::value(x));
292 
293  Scalar x0 = xValues[segIdx];
294  Scalar x1 = xValues[segIdx + 1];
295 
296  Scalar y0 = yValues[segIdx];
297  Scalar y1 = yValues[segIdx + 1];
298 
299  return Toolbox::createConstant((y1 - y0)/(x1 - x0));
300  }
301 
302  static size_t findSegmentIndex_(const ValueVector &xValues, Scalar x)
303  {
304  assert(xValues.size() > 1); // we need at least two sampling points!
305  size_t n = xValues.size() - 1;
306  if (xValues.back() <= x)
307  return n - 1;
308  else if (x <= xValues.front())
309  return 0;
310 
311  // bisection
312  size_t lowIdx = 0, highIdx = n;
313  while (lowIdx + 1 < highIdx) {
314  size_t curIdx = (lowIdx + highIdx)/2;
315  if (xValues[curIdx] < x)
316  lowIdx = curIdx;
317  else
318  highIdx = curIdx;
319  }
320 
321  return lowIdx;
322  }
323 
324  static size_t findSegmentIndexDescending_(const ValueVector &xValues, Scalar x)
325  {
326  assert(xValues.size() > 1); // we need at least two sampling points!
327  size_t n = xValues.size() - 1;
328  if (x <= xValues.back())
329  return n;
330  else if (xValues.front() <= x)
331  return 0;
332 
333  // bisection
334  size_t lowIdx = 0, highIdx = n;
335  while (lowIdx + 1 < highIdx) {
336  size_t curIdx = (lowIdx + highIdx)/2;
337  if (xValues[curIdx] >= x)
338  lowIdx = curIdx;
339  else
340  highIdx = curIdx;
341  }
342 
343  return lowIdx;
344  }
345 };
346 } // namespace Opm
347 
348 #endif
TraitsT Traits
The traits class for this material law.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:57
ParamsT Params
The type of the parameter objects for this law.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:60
Definition: Air_Mesitylene.hpp:31
static const int numPhases
The number of fluid phases.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:66
Traits::Scalar Scalar
The type of the scalar values for this law.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:63
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:51
A traits class which provides basic mathematical functions for arbitrary scalar floating point values...