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  Copyright (C) 2009-2013 by Andreas Lauser
5 
6  This file is part of the Open Porous Media project (OPM).
7 
8  OPM is free software: you can redistribute it and/or modify
9  it under the terms of the GNU General Public License as published by
10  the Free Software Foundation, either version 2 of the License, or
11  (at your option) any later version.
12 
13  OPM is distributed in the hope that it will be useful,
14  but WITHOUT ANY WARRANTY; without even the implied warranty of
15  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16  GNU General Public License for more details.
17 
18  You should have received a copy of the GNU General Public License
19  along with OPM. If not, see <http://www.gnu.org/licenses/>.
20 */
25 #ifndef OPM_SPLINE_TWO_PHASE_MATERIAL_HPP
26 #define OPM_SPLINE_TWO_PHASE_MATERIAL_HPP
27 
29 
30 #include <opm/common/ErrorMacros.hpp>
31 #include <opm/common/Exceptions.hpp>
32 
33 #include <algorithm>
34 #include <cmath>
35 #include <cassert>
36 
37 namespace Opm {
44 template <class TraitsT, class ParamsT = SplineTwoPhaseMaterialParams<TraitsT> >
45 class SplineTwoPhaseMaterial : public TraitsT
46 {
47  typedef typename ParamsT::SamplePoints SamplePoints;
48 
49 public:
51  typedef TraitsT Traits;
52 
54  typedef ParamsT Params;
55 
57  typedef typename Traits::Scalar Scalar;
58 
60  static const int numPhases = Traits::numPhases;
61  static_assert(numPhases == 2,
62  "The piecewise linear two-phase capillary pressure law only"
63  "applies to the case of two fluid phases");
64 
67  static const bool implementsTwoPhaseApi = true;
68 
71  static const bool implementsTwoPhaseSatApi = true;
72 
75  static const bool isSaturationDependent = true;
76 
79  static const bool isPressureDependent = false;
80 
83  static const bool isTemperatureDependent = false;
84 
87  static const bool isCompositionDependent = false;
88 
92  template <class Container, class FluidState>
93  static void capillaryPressures(Container &values, const Params &params, const FluidState &fluidState)
94  {
95  typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
96 
97  values[Traits::wettingPhaseIdx] = 0.0; // reference phase
98  values[Traits::nonWettingPhaseIdx] = pcnw<FluidState, Evaluation>(params, fluidState);
99  }
100 
105  template <class Container, class FluidState>
106  static void saturations(Container &/*values*/, const Params &/*params*/, const FluidState &/*fluidState*/)
107  { OPM_THROW(std::logic_error, "Not implemented: saturations()"); }
108 
112  template <class Container, class FluidState>
113  static void relativePermeabilities(Container &values, const Params &params, const FluidState &fluidState)
114  {
115  typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
116 
117  values[Traits::wettingPhaseIdx] = krw<FluidState, Evaluation>(params, fluidState);
118  values[Traits::nonWettingPhaseIdx] = krn<FluidState, Evaluation>(params, fluidState);
119  }
120 
124  template <class FluidState, class Evaluation = typename FluidState::Scalar>
125  static Evaluation pcnw(const Params &params, const FluidState &fluidState)
126  {
127  typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
128 
129  const Evaluation& Sw =
130  FsToolbox::template toLhs<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  if (Sw <= params.pcnwSpline().xMin())
143  return Evaluation(params.pcnwSpline().yFirst());
144  if (Sw >= params.pcnwSpline().xMax())
145  return Evaluation(params.pcnwSpline().yLast());
146 
147  return params.pcnwSpline().eval(Sw);
148  }
149 
150  template <class Evaluation>
151  static Evaluation twoPhaseSatPcnwInv(const Params &params, const Evaluation& pcnw)
152  {
153  static const Evaluation nil(0.0);
154 
155  // this assumes that the capillary pressure is monotonically decreasing
156  if (pcnw >= params.pcnwSpline().yFirst())
157  return Evaluation(params.pcnwSpline().xMin());
158  if (pcnw <= params.pcnwSpline().yLast())
159  return Evaluation(params.pcnwSpline().xMax());
160 
161  // the intersect() method of splines is a bit slow, but this code path is not too
162  // time critical...
163  return params.pcnwSpline().intersect(/*a=*/nil, /*b=*/nil, /*c=*/nil, /*d=*/pcnw);
164  }
165 
169  template <class FluidState, class Evaluation = typename FluidState::Scalar>
170  static Evaluation Sw(const Params &/*params*/, const FluidState &/*fluidState*/)
171  { OPM_THROW(std::logic_error, "Not implemented: Sw()"); }
172 
173  template <class Evaluation>
174  static Evaluation twoPhaseSatSw(const Params &/*params*/, const Evaluation& /*pC*/)
175  { OPM_THROW(std::logic_error, "Not implemented: twoPhaseSatSw()"); }
176 
181  template <class FluidState, class Evaluation = typename FluidState::Scalar>
182  static Evaluation Sn(const Params &params, const FluidState &fluidState)
183  { return 1 - Sw<FluidState, Evaluation>(params, fluidState); }
184 
185  template <class Evaluation>
186  static Evaluation twoPhaseSatSn(const Params &params, const Evaluation& pC)
187  { return 1 - twoPhaseSatSw(params, pC); }
188 
193  template <class FluidState, class Evaluation = typename FluidState::Scalar>
194  static Evaluation krw(const Params &params, const FluidState &fluidState)
195  {
196  typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
197 
198  const Evaluation& Sw =
199  FsToolbox::template toLhs<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  if (Sw <= params.krnSpline().xMin())
208  return Evaluation(params.krwSpline().yFirst());
209  if (Sw >= params.krnSpline().xMax())
210  return Evaluation(params.krwSpline().yLast());
211 
212  return params.krwSpline().eval(Sw);
213  }
214 
215  template <class Evaluation>
216  static Evaluation twoPhaseSatKrwInv(const Params &params, const Evaluation& krw)
217  {
218  static const Evaluation nil(0.0);
219 
220  if (krw <= params.krwSpline().yFirst())
221  return Evaluation(params.krwSpline().xMin());
222  if (krw >= params.krwSpline().yLast())
223  return Evaluation(params.krwSpline().xMax());
224 
225  return params.krwSpline().intersect(/*a=*/nil, /*b=*/nil, /*c=*/nil, /*d=*/krw);
226  }
227 
232  template <class FluidState, class Evaluation = typename FluidState::Scalar>
233  static Evaluation krn(const Params &params, const FluidState &fluidState)
234  {
235  typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
236 
237  const Evaluation& Sn =
238  FsToolbox::template toLhs<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  if (Sw <= params.krnSpline().xMin())
247  return Evaluation(params.krnSpline().yFirst());
248  if (Sw >= params.krnSpline().xMax())
249  return Evaluation(params.krnSpline().yLast());
250 
251  return params.krnSpline().eval(Sw);
252  }
253 
254  template <class Evaluation>
255  static Evaluation twoPhaseSatKrnInv(const Params &params, const Evaluation& krn)
256  {
257  static const Evaluation nil(0.0);
258 
259  if (krn >= params.krnSpline().yFirst())
260  return Evaluation(params.krnSpline().xMin());
261  if (krn <= params.krnSpline().yLast())
262  return Evaluation(params.krnSpline().xMax());
263 
264  return params.krnSpline().intersect(/*a=*/nil, /*b=*/nil, /*c=*/nil, /*d=*/krn);
265  }
266 };
267 } // namespace Opm
268 
269 #endif
Definition: Air_Mesitylene.hpp:31
TraitsT Traits
The traits class for this material law.
Definition: SplineTwoPhaseMaterial.hpp:51
static const int numPhases
The number of fluid phases.
Definition: SplineTwoPhaseMaterial.hpp:60
Specification of the material parameters for a two-phase material law which uses a table and spline-b...
Implementation of a tabulated capillary pressure and relperm law which uses spline curves as interpol...
Definition: SplineTwoPhaseMaterial.hpp:45
ParamsT Params
The type of the parameter objects for this law.
Definition: SplineTwoPhaseMaterial.hpp:54
Traits::Scalar Scalar
The type of the scalar values for this law.
Definition: SplineTwoPhaseMaterial.hpp:57