VanGenuchten.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  Copyright (C) 2010 by Bernd Flemisch
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 OPM_VAN_GENUCHTEN_HPP
27 #define OPM_VAN_GENUCHTEN_HPP
28 
29 #include "VanGenuchtenParams.hpp"
30 
32 
33 #include <algorithm>
34 #include <cmath>
35 #include <cassert>
36 
37 namespace Opm {
53 template <class TraitsT, class ParamsT = VanGenuchtenParams<TraitsT> >
54 class VanGenuchten : public TraitsT
55 {
56 public:
58  typedef TraitsT Traits;
59 
61  typedef ParamsT Params;
62 
64  typedef typename Traits::Scalar Scalar;
65 
67  static const int numPhases = Traits::numPhases;
68  static_assert(numPhases == 2,
69  "The van Genuchten capillary pressure law only "
70  "applies to the case of two fluid phases");
71 
74  static const bool implementsTwoPhaseApi = true;
75 
78  static const bool implementsTwoPhaseSatApi = true;
79 
82  static const bool isSaturationDependent = true;
83 
86  static const bool isPressureDependent = false;
87 
90  static const bool isTemperatureDependent = false;
91 
94  static const bool isCompositionDependent = false;
95 
112  template <class Container, class FluidState>
113  static void capillaryPressures(Container &values, const Params &params, const FluidState &fs)
114  {
115  typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
116 
117  values[Traits::wettingPhaseIdx] = 0.0; // reference phase
118  values[Traits::nonWettingPhaseIdx] = pcnw<FluidState, Evaluation>(params, fs);
119  }
120 
125  template <class Container, class FluidState>
126  static void saturations(Container &values, const Params &params, const FluidState &fs)
127  {
128  typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
129 
130  values[Traits::wettingPhaseIdx] = Sw<FluidState, Evaluation>(params, fs);
131  values[Traits::nonWettingPhaseIdx] = 1 - values[Traits::wettingPhaseIdx];
132  }
133 
144  template <class Container, class FluidState>
145  static void relativePermeabilities(Container &values, const Params &params, const FluidState &fs)
146  {
147  typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
148 
149  values[Traits::wettingPhaseIdx] = krw<FluidState, Evaluation>(params, fs);
150  values[Traits::nonWettingPhaseIdx] = krn<FluidState, Evaluation>(params, fs);
151  }
152 
167  template <class FluidState, class Evaluation = typename FluidState::Scalar>
168  static Evaluation pcnw(const Params &params, const FluidState &fs)
169  {
170  typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
171 
172  const Evaluation& Sw =
173  FsToolbox::template toLhs<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
174 
175  assert(0 <= Sw && Sw <= 1);
176 
177  return twoPhaseSatPcnw(params, Sw);
178  }
179 
194  template <class Evaluation>
195  static Evaluation twoPhaseSatPcnw(const Params &params, const Evaluation& Sw)
196  {
197  typedef MathToolbox<Evaluation> Toolbox;
198 
199  return Toolbox::pow(Toolbox::pow(Sw, -1.0/params.vgM()) - 1, 1.0/params.vgN())/params.vgAlpha();
200  }
201 
214  template <class FluidState, class Evaluation = typename FluidState::Scalar>
215  static Evaluation Sw(const Params &params, const FluidState &fs)
216  {
217  typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
218 
219  Evaluation pC =
220  FsToolbox::template toLhs<Evaluation>(fs.pressure(Traits::nonWettingPhaseIdx))
221  - FsToolbox::template toLhs<Evaluation>(fs.pressure(Traits::wettingPhaseIdx));
222  return twoPhaseSatSw(params, pC);
223  }
224 
225  template <class Evaluation>
226  static Evaluation twoPhaseSatSw(const Params &params, const Evaluation& pC)
227  {
228  typedef MathToolbox<Evaluation> Toolbox;
229 
230  assert(pC >= 0);
231 
232  return Toolbox::pow(Toolbox::pow(params.vgAlpha()*pC, params.vgN()) + 1, -params.vgM());
233  }
234 
239  template <class FluidState, class Evaluation = typename FluidState::Scalar>
240  static Evaluation Sn(const Params &params, const FluidState &fs)
241  { return 1 - Sw<FluidState, Evaluation>(params, fs); }
242 
243  template <class Evaluation>
244  static Evaluation twoPhaseSatSn(const Params &params, const Evaluation& pC)
245  { return 1 - twoPhaseSatSw(params, pC); }
246 
257  template <class FluidState, class Evaluation = typename FluidState::Scalar>
258  static Evaluation krw(const Params &params, const FluidState &fs)
259  {
260  typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
261 
262  const Evaluation& Sw =
263  FsToolbox::template toLhs<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
264 
265  return twoPhaseSatKrw(params, Sw);
266  }
267 
268  template <class Evaluation>
269  static Evaluation twoPhaseSatKrw(const Params &params, const Evaluation& Sw)
270  {
271  typedef MathToolbox<Evaluation> Toolbox;
272 
273  assert(0.0 <= Sw && Sw <= 1.0);
274 
275  Evaluation r = 1.0 - Toolbox::pow(1.0 - Toolbox::pow(Sw, 1/params.vgM()), params.vgM());
276  return Toolbox::sqrt(Sw)*r*r;
277  }
278 
288  template <class FluidState, class Evaluation = typename FluidState::Scalar>
289  static Evaluation krn(const Params &params, const FluidState &fs)
290  {
291  typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
292 
293  const Evaluation& Sw =
294  1.0 - FsToolbox::template toLhs<Evaluation>(fs.saturation(Traits::nonWettingPhaseIdx));
295 
296  return twoPhaseSatKrn(params, Sw);
297  }
298 
299  template <class Evaluation>
300  static Evaluation twoPhaseSatKrn(const Params &params, Evaluation Sw)
301  {
302  typedef MathToolbox<Evaluation> Toolbox;
303 
304  assert(0 <= Sw && Sw <= 1);
305 
306  return
307  Toolbox::pow(1 - Sw, 1.0/3) *
308  Toolbox::pow(1 - Toolbox::pow(Sw, 1/params.vgM()), 2*params.vgM());
309  }
310 };
311 } // namespace Opm
312 
313 #endif
Definition: Air_Mesitylene.hpp:31
Implementation of the van Genuchten capillary pressure - saturation relation.
Definition: VanGenuchten.hpp:54
Evaluation< Scalar, VarSetTag, numVars > sqrt(const Evaluation< Scalar, VarSetTag, numVars > &x)
Definition: Math.hpp:278
ParamsT Params
The type of the parameter objects for this law.
Definition: VanGenuchten.hpp:61
static const int numPhases
The number of fluid phases.
Definition: VanGenuchten.hpp:67
Evaluation< Scalar, VarSetTag, numVars > pow(const Evaluation< Scalar, VarSetTag, numVars > &base, Scalar exp)
Definition: Math.hpp:312
Traits::Scalar Scalar
The type of the scalar values for this law.
Definition: VanGenuchten.hpp:64
Specification of the material parameters for the van Genuchten constitutive relations.
TraitsT Traits
The traits class for this material law.
Definition: VanGenuchten.hpp:58
A traits class which provides basic mathematical functions for arbitrary scalar floating point values...