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