LinearMaterial.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_LINEAR_MATERIAL_HPP
26 #define OPM_LINEAR_MATERIAL_HPP
27 
28 #include "LinearMaterialParams.hpp"
29 
32 #include <opm/common/Exceptions.hpp>
33 #include <opm/common/ErrorMacros.hpp>
34 
35 #include <algorithm>
36 #include <type_traits>
37 
38 namespace Opm {
39 
50 template <class TraitsT, class ParamsT = LinearMaterialParams<TraitsT> >
51 class LinearMaterial : public TraitsT
52 {
53 public:
54  typedef TraitsT Traits;
55  typedef ParamsT Params;
56  typedef typename Traits::Scalar Scalar;
57 
59  static const int numPhases = Traits::numPhases;
60 
63  static const bool implementsTwoPhaseApi = (numPhases == 2);
64 
67  static const bool implementsTwoPhaseSatApi = (numPhases == 2);
68 
71  static const bool isSaturationDependent = true;
72 
75  static const bool isPressureDependent = false;
76 
79  static const bool isTemperatureDependent = false;
80 
83  static const bool isCompositionDependent = false;
84 
97  template <class ContainerT, class FluidState>
98  static void capillaryPressures(ContainerT &values,
99  const Params &params,
100  const FluidState &state)
101  {
102  typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
104 
105  for (unsigned phaseIdx = 0; phaseIdx < Traits::numPhases; ++phaseIdx) {
106  const Evaluation& S =
107  FsToolbox::template toLhs<Evaluation>(state.saturation(phaseIdx));
109 
110  values[phaseIdx] =
111  S*params.pcMaxSat(phaseIdx) +
112  (1.0 - S)*params.pcMinSat(phaseIdx);
113  }
114  }
115 
119  template <class ContainerT, class FluidState>
120  static void saturations(ContainerT &/*values*/,
121  const Params &/*params*/,
122  const FluidState &/*state*/)
123  {
124  OPM_THROW(std::runtime_error, "Not implemented: LinearMaterial::saturations()");
125  }
126 
130  template <class ContainerT, class FluidState>
131  static void relativePermeabilities(ContainerT &values,
132  const Params &/*params*/,
133  const FluidState &state)
134  {
135  typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
136  typedef MathToolbox<Evaluation> Toolbox;
138 
139  for (unsigned phaseIdx = 0; phaseIdx < Traits::numPhases; ++phaseIdx) {
140  const Evaluation& S =
141  FsToolbox::template toLhs<Evaluation>(state.saturation(phaseIdx));
143 
144  values[phaseIdx] = Toolbox::max(Toolbox::min(S,1.0),0.0);
145  }
146  }
147 
151  template <class FluidState, class Evaluation = typename FluidState::Scalar>
152  static Evaluation pcnw(const Params &params, const FluidState &fs)
153  {
155  const Evaluation& Sw =
156  FsToolbox::template toLhs<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
158 
159  const Evaluation& wPhasePressure =
160  Sw*params.pcMaxSat(Traits::wettingPhaseIdx) +
161  (1.0 - Sw)*params.pcMinSat(Traits::wettingPhaseIdx);
162 
163  const Evaluation& Sn =
164  FsToolbox::template toLhs<Evaluation>(fs.saturation(Traits::nonWettingPhaseIdx));
166 
167  const Evaluation& nPhasePressure =
168  Sn*params.pcMaxSat(Traits::nonWettingPhaseIdx) +
169  (1.0 - Sn)*params.pcMinSat(Traits::nonWettingPhaseIdx);
170 
171  return nPhasePressure - wPhasePressure;
172  }
173 
174 
175  template <class Evaluation = Scalar>
176  static typename std::enable_if<Traits::numPhases == 2, Evaluation>::type
177  twoPhaseSatPcnw(const Params &params, const Evaluation& Sw)
178  {
179  const Evaluation& wPhasePressure =
180  Sw*params.pcMaxSat(Traits::wettingPhaseIdx) +
181  (1.0 - Sw)*params.pcMinSat(Traits::wettingPhaseIdx);
182 
183  const Evaluation& nPhasePressure =
184  (1.0 - Sw)*params.pcMaxSat(Traits::nonWettingPhaseIdx) +
185  Sw*params.pcMinSat(Traits::nonWettingPhaseIdx);
186 
187  return nPhasePressure - wPhasePressure;
188  }
189 
194  template <class FluidState, class Evaluation = typename FluidState::Scalar>
195  static Evaluation Sw(const Params &/*params*/, const FluidState &/*fs*/)
196  { OPM_THROW(std::runtime_error, "Not implemented: Sw()"); }
197 
198  template <class Evaluation = Scalar>
199  static typename std::enable_if<Traits::numPhases == 2, Evaluation>::type
200  twoPhaseSatSw(const Params &/*params*/, const Evaluation& /*Sw*/)
201  { OPM_THROW(std::runtime_error, "Not implemented: twoPhaseSatSw()"); }
202 
207  template <class FluidState, class Evaluation = typename FluidState::Scalar>
208  static Evaluation Sn(const Params &/*params*/, const FluidState &/*fs*/)
209  { OPM_THROW(std::runtime_error, "Not implemented: Sn()"); }
210 
211  template <class Evaluation = Scalar>
212  static typename std::enable_if<Traits::numPhases == 2, Evaluation>::type
213  twoPhaseSatSn(const Params &/*params*/, const Evaluation& /*Sw*/)
214  { OPM_THROW(std::runtime_error, "Not implemented: twoPhaseSatSn()"); }
215 
222  template <class FluidState, class Evaluation = typename FluidState::Scalar>
223  static typename std::enable_if<Traits::numPhases == 3, Evaluation>::type
224  Sg(const Params &/*params*/, const FluidState &/*fs*/)
225  { OPM_THROW(std::runtime_error, "Not implemented: Sg()"); }
226 
230  template <class FluidState, class Evaluation = typename FluidState::Scalar>
231  static Evaluation krw(const Params &/*params*/, const FluidState &fs)
232  {
233  typedef MathToolbox<Evaluation> Toolbox;
235 
236  const Evaluation& Sw =
237  FsToolbox::template toLhs<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
238  return Toolbox::max(0.0, Toolbox::min(1.0, Sw));
239  }
240 
241  template <class Evaluation = Scalar>
242  static typename std::enable_if<Traits::numPhases == 2, Evaluation>::type
243  twoPhaseSatKrw(const Params &/*params*/, const Evaluation& Sw)
244  {
245  typedef MathToolbox<Evaluation> Toolbox;
246 
247  return Toolbox::max(0.0, Toolbox::min(1.0, Sw));
248  }
249 
253  template <class FluidState, class Evaluation = typename FluidState::Scalar>
254  static Evaluation krn(const Params &/*params*/, const FluidState &fs)
255  {
256  typedef MathToolbox<Evaluation> Toolbox;
258 
259  const Evaluation& Sn =
260  FsToolbox::template toLhs<Evaluation>(fs.saturation(Traits::nonWettingPhaseIdx));
261  return Toolbox::max(0.0, Toolbox::min(1.0, Sn));
262  }
263 
264  template <class Evaluation = Scalar>
265  static typename std::enable_if<Traits::numPhases == 2, Evaluation>::type
266  twoPhaseSatKrn(const Params &/*params*/, const Evaluation& Sw)
267  {
268  typedef MathToolbox<Evaluation> Toolbox;
269 
270  return Toolbox::max(0.0, Toolbox::min(1.0, Sw));
271  }
272 
278  template <class FluidState, class Evaluation = typename FluidState::Scalar>
279  static typename std::enable_if<Traits::numPhases == 3, Evaluation>::type
280  krg(const Params &/*params*/, const FluidState &fs)
281  {
282  typedef MathToolbox<Evaluation> Toolbox;
284 
285  const Evaluation& Sg =
286  FsToolbox::template toLhs<Evaluation>(fs.saturation(Traits::gasPhaseIdx));
287  return Toolbox::max(0.0, Toolbox::min(1.0, Sg));
288  }
289 
295  template <class FluidState, class Evaluation = Scalar>
296  static typename std::enable_if<Traits::numPhases == 3, Evaluation>::type
297  pcgn(const Params &params, const FluidState &fs)
298  {
300 
301  const Evaluation& Sn =
302  FsToolbox::template toLhs<Evaluation>(fs.saturation(Traits::nonWettingPhaseIdx));
304 
305  const Evaluation& nPhasePressure =
306  Sn*params.pcMaxSat(Traits::nonWettingPhaseIdx) +
307  (1.0 - Sn)*params.pcMinSat(Traits::nonWettingPhaseIdx);
308 
309  const Evaluation& Sg =
310  FsToolbox::template toLhs<Evaluation>(fs.saturation(Traits::gasPhaseIdx));
312 
313  const Evaluation& gPhasePressure =
314  Sg*params.pcMaxSat(Traits::gasPhaseIdx) +
315  (1.0 - Sg)*params.pcMinSat(Traits::gasPhaseIdx);
316 
317  return gPhasePressure - nPhasePressure;
318  }
319 };
320 } // namespace Opm
321 
322 #endif
Reference implementation of params for the linear M-phase material material.
bool CheckDefined(const T &value OPM_UNUSED)
Make valgrind complain if any of the memory occupied by an object is undefined.
Definition: Valgrind.hpp:74
static std::enable_if< Traits::numPhases==2, Evaluation >::type twoPhaseSatSw(const Params &, const Evaluation &)
Definition: LinearMaterial.hpp:200
static std::enable_if< Traits::numPhases==2, Evaluation >::type twoPhaseSatKrw(const Params &, const Evaluation &Sw)
Definition: LinearMaterial.hpp:243
Definition: MathToolbox.hpp:39
static void saturations(ContainerT &, const Params &, const FluidState &)
The inverse of the capillary pressure.
Definition: LinearMaterial.hpp:120
Definition: Air_Mesitylene.hpp:31
static void relativePermeabilities(ContainerT &values, const Params &, const FluidState &state)
The relative permeability of all phases.
Definition: LinearMaterial.hpp:131
static Evaluation Sw(const Params &, const FluidState &)
Calculate wetting phase saturation given that the rest of the fluid state has been initialized...
Definition: LinearMaterial.hpp:195
static void capillaryPressures(ContainerT &values, const Params &params, const FluidState &state)
The linear capillary pressure-saturation curve.
Definition: LinearMaterial.hpp:98
static const bool isSaturationDependent
Definition: LinearMaterial.hpp:71
Some templates to wrap the valgrind client request macros.
Evaluation< Scalar, VarSetTag, numVars > max(const Evaluation< Scalar, VarSetTag, numVars > &x1, const Evaluation< Scalar, VarSetTag, numVars > &x2)
Definition: Math.hpp:114
static std::enable_if< Traits::numPhases==3, Evaluation >::type krg(const Params &, const FluidState &fs)
The relative permability of the gas phase.
Definition: LinearMaterial.hpp:280
static std::enable_if< Traits::numPhases==2, Evaluation >::type twoPhaseSatPcnw(const Params &params, const Evaluation &Sw)
Definition: LinearMaterial.hpp:177
static const bool isPressureDependent
Definition: LinearMaterial.hpp:75
static const bool isTemperatureDependent
Definition: LinearMaterial.hpp:79
Traits::Scalar Scalar
Definition: LinearMaterial.hpp:56
static const bool implementsTwoPhaseSatApi
Definition: LinearMaterial.hpp:67
Implements a linear saturation-capillary pressure relation.
Definition: LinearMaterial.hpp:51
Evaluation< Scalar, VarSetTag, numVars > min(const Evaluation< Scalar, VarSetTag, numVars > &x1, const Evaluation< Scalar, VarSetTag, numVars > &x2)
Definition: Math.hpp:61
ParamsT Params
Definition: LinearMaterial.hpp:55
TraitsT Traits
Definition: LinearMaterial.hpp:54
static const int numPhases
The number of fluid phases.
Definition: LinearMaterial.hpp:59
static const bool implementsTwoPhaseApi
Definition: LinearMaterial.hpp:63
static std::enable_if< Traits::numPhases==3, Evaluation >::type Sg(const Params &, const FluidState &)
Calculate gas phase saturation given that the rest of the fluid state has been initialized.
Definition: LinearMaterial.hpp:224
static std::enable_if< Traits::numPhases==2, Evaluation >::type twoPhaseSatSn(const Params &, const Evaluation &)
Definition: LinearMaterial.hpp:213
static Evaluation krw(const Params &, const FluidState &fs)
The relative permability of the wetting phase.
Definition: LinearMaterial.hpp:231
static Evaluation Sn(const Params &, const FluidState &)
Calculate non-wetting liquid phase saturation given that the rest of the fluid state has been initial...
Definition: LinearMaterial.hpp:208
static Evaluation krn(const Params &, const FluidState &fs)
The relative permability of the liquid non-wetting phase.
Definition: LinearMaterial.hpp:254
static std::enable_if< Traits::numPhases==3, Evaluation >::type pcgn(const Params &params, const FluidState &fs)
The difference between the pressures of the gas and the non-wetting phase.
Definition: LinearMaterial.hpp:297
static const bool isCompositionDependent
Definition: LinearMaterial.hpp:83
static std::enable_if< Traits::numPhases==2, Evaluation >::type twoPhaseSatKrn(const Params &, const Evaluation &Sw)
Definition: LinearMaterial.hpp:266
static Evaluation pcnw(const Params &params, const FluidState &fs)
The difference between the pressures of the non-wetting and wetting phase.
Definition: LinearMaterial.hpp:152
A traits class which provides basic mathematical functions for arbitrary scalar floating point values...