RegularizedBrooksCorey.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) 2008-2013 by Andreas Lauser
5  Copyright (C) 2010 by Philipp Nuske
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 REGULARIZED_BROOKS_COREY_HPP
27 #define REGULARIZED_BROOKS_COREY_HPP
28 
29 #include "BrooksCorey.hpp"
31 
33 
34 namespace Opm {
63 template <class TraitsT, class ParamsT = RegularizedBrooksCoreyParams<TraitsT> >
64 class RegularizedBrooksCorey : public TraitsT
65 {
67 
68 public:
69  typedef TraitsT Traits;
70  typedef ParamsT Params;
71  typedef typename Traits::Scalar Scalar;
72 
74  static const int numPhases = Traits::numPhases;
75  static_assert(numPhases == 2,
76  "The regularized Brooks-Corey capillary pressure law only "
77  "applies to the case of two fluid phases");
78 
81  static const bool implementsTwoPhaseApi = true;
82 
85  static const bool implementsTwoPhaseSatApi = true;
86 
89  static const bool isSaturationDependent = true;
90 
93  static const bool isPressureDependent = false;
94 
97  static const bool isTemperatureDependent = false;
98 
101  static const bool isCompositionDependent = false;
102 
103  static_assert(Traits::numPhases == 2,
104  "The number of fluid phases must be two if you want to use "
105  "this material law!");
106 
117  template <class Container, class FluidState>
118  static void capillaryPressures(Container &values, const Params &params, const FluidState &fs)
119  {
120  typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
121 
122  values[Traits::wettingPhaseIdx] = 0.0; // reference phase
123  values[Traits::nonWettingPhaseIdx] = pcnw<FluidState, Evaluation>(params, fs);
124  }
125 
130  template <class Container, class FluidState>
131  static void saturations(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] = Sw<FluidState, Evaluation>(params, fs);
136  values[Traits::nonWettingPhaseIdx] = 1 - values[Traits::wettingPhaseIdx];
137  }
138 
149  template <class Container, class FluidState>
150  static void relativePermeabilities(Container &values, const Params &params, const FluidState &fs)
151  {
152  typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
153 
154  values[Traits::wettingPhaseIdx] = krw<FluidState, Evaluation>(params, fs);
155  values[Traits::nonWettingPhaseIdx] = krn<FluidState, Evaluation>(params, fs);
156  }
157 
172  template <class FluidState, class Evaluation = typename FluidState::Scalar>
173  static Evaluation pcnw(const Params &params, const FluidState &fs)
174  {
175  typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
176 
177  const auto& Sw = FsToolbox::template toLhs<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
178  return twoPhaseSatPcnw(params, Sw);
179  }
180 
181  template <class Evaluation>
182  static Evaluation twoPhaseSatPcnw(const Params &params, const Evaluation& Sw)
183  {
184  const Scalar Sthres = params.pcnwLowSw();
185 
186  if (Sw <= Sthres) {
187  Scalar m = params.pcnwSlopeLow();
188  Scalar pcnw_SwLow = params.pcnwLow();
189  return pcnw_SwLow + m*(Sw - Sthres);
190  }
191  else if (Sw >= 1.0) {
192  Scalar m = params.pcnwSlopeHigh();
193  Scalar pcnw_SwHigh = params.pcnwHigh();
194  return pcnw_SwHigh + m*(Sw - 1.0);
195  }
196 
197  // if the effective saturation is in an 'reasonable'
198  // range, we use the real Brooks-Corey law...
199  return BrooksCorey::twoPhaseSatPcnw(params, Sw);
200  }
201 
208  template <class FluidState, class Evaluation = typename FluidState::Scalar>
209  static Evaluation Sw(const Params &params, const FluidState &fs)
210  {
211  typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
212 
213  const Evaluation& pC =
214  FsToolbox::template toLhs<Evaluation>(fs.pressure(Traits::nonWettingPhaseIdx))
215  - FsToolbox::template toLhs<Evaluation>(fs.pressure(Traits::wettingPhaseIdx));
216  return twoPhaseSatSw(params, pC);
217  }
218 
219  template <class Evaluation>
220  static Evaluation twoPhaseSatSw(const Params &params, const Evaluation& pcnw)
221  {
222  const Scalar Sthres = params.pcnwLowSw();
223 
224  // calculate the saturation which corrosponds to the
225  // saturation in the non-regularized version of the
226  // Brooks-Corey law. If the input capillary pressure is
227  // smaller than the entry pressure, make sure that we will
228  // regularize.
229  Evaluation Sw = 1.5;
230  if (pcnw >= params.entryPressure())
231  Sw = BrooksCorey::twoPhaseSatSw(params, pcnw);
232 
233  // make sure that the capilary pressure observes a
234  // derivative != 0 for 'illegal' saturations. This is
235  // required for example by newton solvers (if the
236  // derivative calculated numerically) in order to get the
237  // saturation moving to the right direction if it
238  // temporarily is in an 'illegal' range.
239  if (Sw <= Sthres) {
240  // invert the low saturation regularization of pcnw()
241  Scalar m = params.pcnwSlopeLow();
242  Scalar pcnw_SwLow = params.pcnwLow();
243  return Sthres + (pcnw - pcnw_SwLow)/m;
244  }
245  else if (Sw > 1.0) {
246  Scalar m = params.pcnwSlopeHigh();
247  Scalar pcnw_SwHigh = params.pcnwHigh();
248  return 1.0 + (pcnw - pcnw_SwHigh)/m;;
249  }
250 
251  return BrooksCorey::twoPhaseSatSw(params, pcnw);
252  }
253 
258  template <class FluidState, class Evaluation = typename FluidState::Scalar>
259  static Evaluation Sn(const Params &params, const FluidState &fs)
260  { return 1 - Sw<FluidState, Evaluation>(params, fs); }
261 
262  template <class Evaluation>
263  static Evaluation twoPhaseSatSn(const Params &params, const Evaluation& pcnw)
264  { return 1 - twoPhaseSatSw(params, pcnw); }
265 
280  template <class FluidState, class Evaluation = typename FluidState::Scalar>
281  static Evaluation krw(const Params &params, const FluidState &fs)
282  {
283  typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
284 
285  const auto& Sw = FsToolbox::template toLhs<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
286  return twoPhaseSatKrw(params, Sw);
287  }
288 
289  template <class Evaluation>
290  static Evaluation twoPhaseSatKrw(const Params &params, const Evaluation& Sw)
291  {
292  typedef MathToolbox<Evaluation> Toolbox;
293 
294  if (Sw <= 0.0)
295  return Toolbox::createConstant(0.0);
296  else if (Sw >= 1.0)
297  return Toolbox::createConstant(1.0);
298 
299  return BrooksCorey::twoPhaseSatKrw(params, Sw);
300  }
301 
302  template <class Evaluation>
303  static Evaluation twoPhaseSatKrwInv(const Params &params, const Evaluation& krw)
304  {
305  typedef MathToolbox<Evaluation> Toolbox;
306 
307  if (krw <= 0.0)
308  return Toolbox::createConstant(0.0);
309  else if (krw >= 1.0)
310  return Toolbox::createConstant(1.0);
311 
312  return BrooksCorey::twoPhaseSatKrwInv(params, krw);
313  }
314 
329  template <class FluidState, class Evaluation = typename FluidState::Scalar>
330  static Evaluation krn(const Params &params, const FluidState &fs)
331  {
332  typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
333 
334  const Evaluation& Sw =
335  1.0 - FsToolbox::template toLhs<Evaluation>(fs.saturation(Traits::nonWettingPhaseIdx));
336  return twoPhaseSatKrn(params, Sw);
337  }
338 
339  template <class Evaluation>
340  static Evaluation twoPhaseSatKrn(const Params &params, const Evaluation& Sw)
341  {
342  typedef MathToolbox<Evaluation> Toolbox;
343 
344  if (Sw >= 1.0)
345  return Toolbox::createConstant(0.0);
346  else if (Sw <= 0.0)
347  return Toolbox::createConstant(1.0);
348 
349  return BrooksCorey::twoPhaseSatKrn(params, Sw);
350  }
351 
352  template <class Evaluation>
353  static Evaluation twoPhaseSatKrnInv(const Params &params, const Evaluation& krn)
354  {
355  typedef MathToolbox<Evaluation> Toolbox;
356 
357  if (krn <= 0.0)
358  return Toolbox::createConstant(1.0);
359  else if (krn >= 1.0)
360  return Toolbox::createConstant(0.0);
361 
362  return BrooksCorey::twoPhaseSatKrnInv(params, krn);
363  }
364 };
365 } // namespace Opm
366 
367 #endif
Implementation of the Brooks-Corey capillary pressure <-> saturation relation.
Definition: Air_Mesitylene.hpp:31
ParamsT Params
Definition: RegularizedBrooksCorey.hpp:70
Class implementing cubic splines.
Implementation of the Brooks-Corey capillary pressure <-> saturation relation.
Definition: BrooksCorey.hpp:54
Parameters that are necessary for the regularization of the Brooks-Corey capillary pressure model...
Implementation of the regularized Brooks-Corey capillary pressure / relative permeability <-> saturat...
Definition: RegularizedBrooksCorey.hpp:64
static const int numPhases
The number of fluid phases.
Definition: RegularizedBrooksCorey.hpp:74
Traits::Scalar Scalar
Definition: RegularizedBrooksCorey.hpp:71
TraitsT Traits
Definition: RegularizedBrooksCorey.hpp:69