All Classes Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
SplineTwoPhaseMaterial.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_SPLINE_TWO_PHASE_MATERIAL_HPP
28 #define OPM_SPLINE_TWO_PHASE_MATERIAL_HPP
29 
31 
32 #include <opm/common/ErrorMacros.hpp>
33 #include <opm/common/Exceptions.hpp>
34 
35 #include <algorithm>
36 #include <cmath>
37 #include <cassert>
38 
39 namespace Opm {
46 template <class TraitsT, class ParamsT = SplineTwoPhaseMaterialParams<TraitsT> >
47 class SplineTwoPhaseMaterial : public TraitsT
48 {
49  typedef typename ParamsT::SamplePoints SamplePoints;
50 
51 public:
53  typedef TraitsT Traits;
54 
56  typedef ParamsT Params;
57 
59  typedef typename Traits::Scalar Scalar;
60 
62  static const int numPhases = Traits::numPhases;
63  static_assert(numPhases == 2,
64  "The piecewise linear two-phase capillary pressure law only"
65  "applies to the case of two fluid phases");
66 
69  static const bool implementsTwoPhaseApi = true;
70 
73  static const bool implementsTwoPhaseSatApi = true;
74 
77  static const bool isSaturationDependent = true;
78 
81  static const bool isPressureDependent = false;
82 
85  static const bool isTemperatureDependent = false;
86 
89  static const bool isCompositionDependent = false;
90 
94  template <class Container, class FluidState>
95  static void capillaryPressures(Container& values, const Params& params, const FluidState& fluidState)
96  {
97  typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
98 
99  values[Traits::wettingPhaseIdx] = 0.0; // reference phase
100  values[Traits::nonWettingPhaseIdx] = pcnw<FluidState, Evaluation>(params, fluidState);
101  }
102 
107  template <class Container, class FluidState>
108  static void saturations(Container& /*values*/, const Params& /*params*/, const FluidState& /*fluidState*/)
109  { OPM_THROW(std::logic_error, "Not implemented: saturations()"); }
110 
114  template <class Container, class FluidState>
115  static void relativePermeabilities(Container& values, const Params& params, const FluidState& fluidState)
116  {
117  typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
118 
119  values[Traits::wettingPhaseIdx] = krw<FluidState, Evaluation>(params, fluidState);
120  values[Traits::nonWettingPhaseIdx] = krn<FluidState, Evaluation>(params, fluidState);
121  }
122 
126  template <class FluidState, class Evaluation = typename FluidState::Scalar>
127  static Evaluation pcnw(const Params& params, const FluidState& fluidState)
128  {
129  const Evaluation& Sw =
130  Opm::decay<Evaluation>(fluidState.saturation(Traits::wettingPhaseIdx));
131 
132  return twoPhaseSatPcnw(params, Sw);
133  }
134 
138  template <class Evaluation>
139  static Evaluation twoPhaseSatPcnw(const Params& params, const Evaluation& Sw)
140  {
141  // this assumes that the capillary pressure is monotonically decreasing
142  const auto& pcnwSpline = params.pcnwSpline();
143  if (Sw <= pcnwSpline.xAt(0))
144  return Evaluation(pcnwSpline.valueAt(0));
145  if (Sw >= pcnwSpline.xAt(pcnwSpline.numSamples() - 1))
146  return Evaluation(pcnwSpline.valueAt(pcnwSpline.numSamples() - 1));
147 
148  return pcnwSpline.eval(Sw);
149  }
150 
151  template <class Evaluation>
152  static Evaluation twoPhaseSatPcnwInv(const Params& params, const Evaluation& pcnw)
153  {
154  static const Evaluation nil(0.0);
155 
156  // this assumes that the capillary pressure is monotonically decreasing
157  const auto& pcnwSpline = params.pcnwSpline();
158  if (pcnw >= pcnwSpline.valueAt(0))
159  return Evaluation(pcnwSpline.xAt(0));
160  if (pcnw <= pcnwSpline.y(pcnwSpline.numSamples() - 1))
161  return Evaluation(pcnwSpline.xAt(pcnwSpline.numSamples() - 1));
162 
163  // the intersect() method of splines is a bit slow, but this code path is not too
164  // time critical...
165  return pcnwSpline.intersect(/*a=*/nil, /*b=*/nil, /*c=*/nil, /*d=*/pcnw);
166  }
167 
171  template <class FluidState, class Evaluation = typename FluidState::Scalar>
172  static Evaluation Sw(const Params& /*params*/, const FluidState& /*fluidState*/)
173  { OPM_THROW(std::logic_error, "Not implemented: Sw()"); }
174 
175  template <class Evaluation>
176  static Evaluation twoPhaseSatSw(const Params& /*params*/, const Evaluation& /*pC*/)
177  { OPM_THROW(std::logic_error, "Not implemented: twoPhaseSatSw()"); }
178 
183  template <class FluidState, class Evaluation = typename FluidState::Scalar>
184  static Evaluation Sn(const Params& params, const FluidState& fluidState)
185  { return 1 - Sw<FluidState, Evaluation>(params, fluidState); }
186 
187  template <class Evaluation>
188  static Evaluation twoPhaseSatSn(const Params& params, const Evaluation& pC)
189  { return 1 - twoPhaseSatSw(params, pC); }
190 
195  template <class FluidState, class Evaluation = typename FluidState::Scalar>
196  static Evaluation krw(const Params& params, const FluidState& fluidState)
197  {
198  const Evaluation& Sw =
199  Opm::decay<Evaluation>(fluidState.saturation(Traits::wettingPhaseIdx));
200 
201  return twoPhaseSatKrw(params, Sw);
202  }
203 
204  template <class Evaluation>
205  static Evaluation twoPhaseSatKrw(const Params& params, const Evaluation& Sw)
206  {
207  const auto& krwSpline = params.krwSpline();
208  if (Sw <= krwSpline.xAt(0))
209  return Evaluation(krwSpline.valueAt(0));
210  if (Sw >= krwSpline.xAt(krwSpline.numSamples() - 1))
211  return Evaluation(krwSpline.valueAt(krwSpline.numSamples() - 1));
212 
213  return krwSpline.eval(Sw);
214  }
215 
216  template <class Evaluation>
217  static Evaluation twoPhaseSatKrwInv(const Params& params, const Evaluation& krw)
218  {
219  static const Evaluation nil(0.0);
220 
221  const auto& krwSpline = params.krwSpline();
222  if (krw <= krwSpline.valueAt(0))
223  return Evaluation(krwSpline.xAt(0));
224  if (krw >= krwSpline.valueAt(krwSpline.numSamples() - 1))
225  return Evaluation(krwSpline.xAt(krwSpline.numSamples() - 1));
226 
227  return krwSpline.intersect(/*a=*/nil, /*b=*/nil, /*c=*/nil, /*d=*/krw);
228  }
229 
234  template <class FluidState, class Evaluation = typename FluidState::Scalar>
235  static Evaluation krn(const Params& params, const FluidState& fluidState)
236  {
237  const Evaluation& Sn =
238  Opm::decay<Evaluation>(fluidState.saturation(Traits::nonWettingPhaseIdx));
239 
240  return twoPhaseSatKrn(params, 1.0 - Sn);
241  }
242 
243  template <class Evaluation>
244  static Evaluation twoPhaseSatKrn(const Params& params, const Evaluation& Sw)
245  {
246  const auto& krnSpline = params.krnSpline();
247  if (Sw <= krnSpline.xAt(0))
248  return Evaluation(krnSpline.valueAt(0));
249  if (Sw >= krnSpline.xAt(krnSpline.numSamples() - 1))
250  return Evaluation(krnSpline.valueAt(krnSpline.numSamples() - 1));
251 
252  return krnSpline.eval(Sw);
253  }
254 
255  template <class Evaluation>
256  static Evaluation twoPhaseSatKrnInv(const Params& params, const Evaluation& krn)
257  {
258  static const Evaluation nil(0.0);
259 
260  const auto& krnSpline = params.krnSpline();
261  if (krn >= krnSpline.valueAt(0))
262  return Evaluation(krnSpline.xAt(0));
263  if (krn <= krnSpline.valueAt(krnSpline.numSamples() - 1))
264  return Evaluation(krnSpline.xAt(krnSpline.numSamples() - 1));
265 
266  return krnSpline.intersect(/*a=*/nil, /*b=*/nil, /*c=*/nil, /*d=*/krn);
267  }
268 };
269 } // namespace Opm
270 
271 #endif
static Evaluation Sn(const Params &params, const FluidState &fluidState)
Calculate the non-wetting phase saturations depending on the phase pressures.
Definition: SplineTwoPhaseMaterial.hpp:184
static void relativePermeabilities(Container &values, const Params &params, const FluidState &fluidState)
The relative permeabilities.
Definition: SplineTwoPhaseMaterial.hpp:115
static Evaluation twoPhaseSatPcnw(const Params &params, const Evaluation &Sw)
The saturation-capillary pressure curve.
Definition: SplineTwoPhaseMaterial.hpp:139
static const bool isCompositionDependent
Specify whether the quantities defined by this material law are dependent on the phase composition...
Definition: SplineTwoPhaseMaterial.hpp:89
static const bool isSaturationDependent
Specify whether the quantities defined by this material law are saturation dependent.
Definition: SplineTwoPhaseMaterial.hpp:77
TraitsT Traits
The traits class for this material law.
Definition: SplineTwoPhaseMaterial.hpp:53
static Evaluation pcnw(const Params &params, const FluidState &fluidState)
The capillary pressure-saturation curve.
Definition: SplineTwoPhaseMaterial.hpp:127
static void capillaryPressures(Container &values, const Params &params, const FluidState &fluidState)
The capillary pressure-saturation curve.
Definition: SplineTwoPhaseMaterial.hpp:95
Implementation of a tabulated capillary pressure and relperm law which uses spline curves as interpol...
Definition: SplineTwoPhaseMaterial.hpp:47
static Evaluation krn(const Params &params, const FluidState &fluidState)
The relative permeability for the non-wetting phase of the porous medium.
Definition: SplineTwoPhaseMaterial.hpp:235
static const int numPhases
The number of fluid phases.
Definition: SplineTwoPhaseMaterial.hpp:62
ParamsT Params
The type of the parameter objects for this law.
Definition: SplineTwoPhaseMaterial.hpp:56
Traits::Scalar Scalar
The type of the scalar values for this law.
Definition: SplineTwoPhaseMaterial.hpp:59
Specification of the material parameters for a two-phase material law which uses a table and spline-b...
static const bool isTemperatureDependent
Specify whether the quantities defined by this material law are temperature dependent.
Definition: SplineTwoPhaseMaterial.hpp:85
static const bool isPressureDependent
Specify whether the quantities defined by this material law are dependent on the absolute pressure...
Definition: SplineTwoPhaseMaterial.hpp:81
static void saturations(Container &, const Params &, const FluidState &)
The saturations of the fluid phases starting from their pressure differences.
Definition: SplineTwoPhaseMaterial.hpp:108
static const bool implementsTwoPhaseApi
Specify whether this material law implements the two-phase convenience API.
Definition: SplineTwoPhaseMaterial.hpp:69
static const bool implementsTwoPhaseSatApi
Specify whether this material law implements the two-phase convenience API which only depends on the ...
Definition: SplineTwoPhaseMaterial.hpp:73
static Evaluation krw(const Params &params, const FluidState &fluidState)
The relative permeability for the wetting phase of the porous medium.
Definition: SplineTwoPhaseMaterial.hpp:196
static Evaluation Sw(const Params &, const FluidState &)
The saturation-capillary pressure curve.
Definition: SplineTwoPhaseMaterial.hpp:172