All Classes Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
BrooksCorey.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_BROOKS_COREY_HPP
28 #define OPM_BROOKS_COREY_HPP
29 
30 #include "BrooksCoreyParams.hpp"
31 
33 #include <opm/common/ErrorMacros.hpp>
34 #include <opm/common/Exceptions.hpp>
35 
36 #include <algorithm>
37 #include <cassert>
38 #include <cmath>
39 
40 namespace Opm {
53 template <class TraitsT, class ParamsT = BrooksCoreyParams<TraitsT> >
54 class BrooksCorey : public TraitsT
55 {
56 public:
57  typedef TraitsT Traits;
58  typedef ParamsT Params;
59  typedef typename Traits::Scalar Scalar;
60 
62  static const int numPhases = Traits::numPhases;
63  static_assert(numPhases == 2,
64  "The Brooks-Corey capillary pressure law only applies "
65  "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 
91  static_assert(Traits::numPhases == 2,
92  "The number of fluid phases must be two if you want to use "
93  "this material law!");
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  {
114  typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
115 
116  values[Traits::wettingPhaseIdx] = Sw<FluidState, Evaluation>(params, fs);
117  values[Traits::nonWettingPhaseIdx] = 1.0 - values[Traits::wettingPhaseIdx];
118  }
119 
130  template <class Container, class FluidState>
131  static void relativePermeabilities(Container& values, const Params& params, const FluidState& fs)
132  {
133  typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
134 
135  values[Traits::wettingPhaseIdx] = krw<FluidState, Evaluation>(params, fs);
136  values[Traits::nonWettingPhaseIdx] = krn<FluidState, Evaluation>(params, fs);
137  }
138 
152  template <class FluidState, class Evaluation = typename FluidState::Scalar>
153  static Evaluation pcnw(const Params& params, const FluidState& fs)
154  {
155  const Evaluation& Sw =
156  Opm::decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
157 
158  assert(0.0 <= Sw && Sw <= 1.0);
159 
160  return twoPhaseSatPcnw(params, Sw);
161  }
162 
163  template <class Evaluation>
164  static Evaluation twoPhaseSatPcnw(const Params& params, const Evaluation& Sw)
165  {
166  assert(0.0 <= Sw && Sw <= 1.0);
167 
168  return params.entryPressure()*Opm::pow(Sw, -1/params.lambda());
169  }
170 
171  template <class Evaluation>
172  static Evaluation twoPhaseSatPcnwInv(const Params& params, const Evaluation& pcnw)
173  {
174  assert(pcnw > 0.0);
175 
176  return Opm::pow(params.entryPressure()/pcnw, -params.lambda());
177  }
178 
191  template <class FluidState, class Evaluation = typename FluidState::Scalar>
192  static Evaluation Sw(const Params& params, const FluidState& fs)
193  {
194  Evaluation pC =
195  Opm::decay<Evaluation>(fs.pressure(Traits::nonWettingPhaseIdx))
196  - Opm::decay<Evaluation>(fs.pressure(Traits::wettingPhaseIdx));
197  return twoPhaseSatSw(params, pC);
198  }
199 
200  template <class Evaluation>
201  static Evaluation twoPhaseSatSw(const Params& params, const Evaluation& pc)
202  {
203  assert(pc > 0.0); // if we don't assume that, std::pow will screw up!
204 
205  return Opm::pow(pc/params.entryPressure(), -params.lambda());
206  }
207 
212  template <class FluidState, class Evaluation = typename FluidState::Scalar>
213  static Evaluation Sn(const Params& params, const FluidState& fs)
214  { return 1.0 - Sw<FluidState, Evaluation>(params, fs); }
215 
216  template <class Evaluation>
217  static Evaluation twoPhaseSatSn(const Params& params, const Evaluation& pc)
218  { return 1.0 - twoPhaseSatSw(params, pc); }
227  template <class FluidState, class Evaluation = typename FluidState::Scalar>
228  static Evaluation krw(const Params& params, const FluidState& fs)
229  {
230  const auto& Sw =
231  Opm::decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
232 
233  return twoPhaseSatKrw(params, Sw);
234  }
235 
236  template <class Evaluation>
237  static Evaluation twoPhaseSatKrw(const Params& params, const Evaluation& Sw)
238  {
239  assert(0.0 <= Sw && Sw <= 1.0);
240 
241  return Opm::pow(Sw, 2.0/params.lambda() + 3.0);
242  }
243 
244  template <class Evaluation>
245  static Evaluation twoPhaseSatKrwInv(const Params& params, const Evaluation& krw)
246  {
247  return Opm::pow(krw, 1.0/(2.0/params.lambda() + 3.0));
248  }
249 
258  template <class FluidState, class Evaluation = typename FluidState::Scalar>
259  static Evaluation krn(const Params& params, const FluidState& fs)
260  {
261  const Evaluation& Sw =
262  1.0 - Opm::decay<Evaluation>(fs.saturation(Traits::nonWettingPhaseIdx));
263 
264  return twoPhaseSatKrn(params, Sw);
265  }
266 
267  template <class Evaluation>
268  static Evaluation twoPhaseSatKrn(const Params& params, const Evaluation& Sw)
269  {
270  assert(0.0 <= Sw && Sw <= 1.0);
271 
272  Scalar exponent = 2.0/params.lambda() + 1.0;
273  const Evaluation Sn = 1.0 - Sw;
274  return Sn*Sn*(1. - Opm::pow(Sw, exponent));
275  }
276 
277  template <class Evaluation>
278  static Evaluation twoPhaseSatKrnInv(const Params& params, const Evaluation& krn)
279  {
280  // since inverting the formula for krn is hard to do analytically, we use the
281  // Newton-Raphson method
282  Evaluation Sw = 0.5;
283  Scalar eps = 1e-10;
284  for (int i = 0; i < 20; ++i) {
285  Evaluation f = krn - twoPhaseSatKrn(params, Sw);
286  Evaluation fStar = krn - twoPhaseSatKrn(params, Sw + eps);
287  Evaluation fPrime = (fStar - f)/eps;
288 
289  Evaluation delta = f/fPrime;
290  Sw -= delta;
291  if (Sw < 0)
292  Sw = 0.0;
293  if (Opm::abs(delta) < 1e-10)
294  return Sw;
295  }
296 
297  OPM_THROW(NumericalProblem,
298  "Couldn't invert the Brooks-Corey non-wetting phase"
299  " relperm within 20 iterations");
300  }
301 };
302 } // namespace Opm
303 
304 #endif // BROOKS_COREY_HPP
static const bool implementsTwoPhaseApi
Specify whether this material law implements the two-phase convenience API.
Definition: BrooksCorey.hpp:69
A traits class which provides basic mathematical functions for arbitrary scalar floating point values...
Specification of the material parameters for the Brooks-Corey constitutive relations.
static Evaluation pcnw(const Params &params, const FluidState &fs)
The capillary pressure-saturation curve according to Brooks and Corey.
Definition: BrooksCorey.hpp:153
static void capillaryPressures(Container &values, const Params &params, const FluidState &fs)
The capillary pressure-saturation curves.
Definition: BrooksCorey.hpp:99
static void relativePermeabilities(Container &values, const Params &params, const FluidState &fs)
The relative permeability-saturation curves.
Definition: BrooksCorey.hpp:131
static Evaluation Sn(const Params &params, const FluidState &fs)
Calculate the non-wetting phase saturations depending on the phase pressures.
Definition: BrooksCorey.hpp:213
static const bool isCompositionDependent
Specify whether the quantities defined by this material law are dependent on the phase composition...
Definition: BrooksCorey.hpp:89
static const bool isTemperatureDependent
Specify whether the quantities defined by this material law are temperature dependent.
Definition: BrooksCorey.hpp:85
static void saturations(Container &values, const Params &params, const FluidState &fs)
Calculate the saturations of the phases starting from their pressure differences. ...
Definition: BrooksCorey.hpp:112
Implementation of the Brooks-Corey capillary pressure &lt;-&gt; saturation relation.
Definition: BrooksCorey.hpp:54
static Evaluation krw(const Params &params, const FluidState &fs)
The relative permeability for the wetting phase of the medium implied by the Brooks-Corey parameteriz...
Definition: BrooksCorey.hpp:228
static const bool isSaturationDependent
Specify whether the quantities defined by this material law are saturation dependent.
Definition: BrooksCorey.hpp:77
static const bool implementsTwoPhaseSatApi
Specify whether this material law implements the two-phase convenience API which only depends on the ...
Definition: BrooksCorey.hpp:73
static const int numPhases
The number of fluid phases to which this material law applies.
Definition: BrooksCorey.hpp:62
static Evaluation Sw(const Params &params, const FluidState &fs)
The saturation-capillary pressure curve according to Brooks &amp; Corey.
Definition: BrooksCorey.hpp:192
static const bool isPressureDependent
Specify whether the quantities defined by this material law are dependent on the absolute pressure...
Definition: BrooksCorey.hpp:81
static Evaluation krn(const Params &params, const FluidState &fs)
The relative permeability for the non-wetting phase of the medium as implied by the Brooks-Corey para...
Definition: BrooksCorey.hpp:259