All Classes Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
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  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_LINEAR_MATERIAL_HPP
28 #define OPM_LINEAR_MATERIAL_HPP
29 
30 #include "LinearMaterialParams.hpp"
31 
33 #include <opm/common/Valgrind.hpp>
34 #include <opm/common/Exceptions.hpp>
35 #include <opm/common/ErrorMacros.hpp>
36 
37 #include <algorithm>
38 #include <type_traits>
39 
40 namespace Opm {
41 
52 template <class TraitsT, class ParamsT = LinearMaterialParams<TraitsT> >
53 class LinearMaterial : 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 
65  static const bool implementsTwoPhaseApi = (numPhases == 2);
66 
69  static const bool implementsTwoPhaseSatApi = (numPhases == 2);
70 
73  static const bool isSaturationDependent = true;
74 
77  static const bool isPressureDependent = false;
78 
81  static const bool isTemperatureDependent = false;
82 
85  static const bool isCompositionDependent = false;
86 
99  template <class ContainerT, class FluidState>
100  static void capillaryPressures(ContainerT& values,
101  const Params& params,
102  const FluidState& state)
103  {
104  typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
105 
106  for (unsigned phaseIdx = 0; phaseIdx < Traits::numPhases; ++phaseIdx) {
107  const Evaluation& S =
108  Opm::decay<Evaluation>(state.saturation(phaseIdx));
109  Valgrind::CheckDefined(S);
110 
111  values[phaseIdx] =
112  S*params.pcMaxSat(phaseIdx) +
113  (1.0 - S)*params.pcMinSat(phaseIdx);
114  }
115  }
116 
120  template <class ContainerT, class FluidState>
121  static void saturations(ContainerT& /*values*/,
122  const Params& /*params*/,
123  const FluidState& /*state*/)
124  {
125  OPM_THROW(std::runtime_error, "Not implemented: LinearMaterial::saturations()");
126  }
127 
131  template <class ContainerT, class FluidState>
132  static void relativePermeabilities(ContainerT& values,
133  const Params& /*params*/,
134  const FluidState& state)
135  {
136  typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
137 
138  for (unsigned phaseIdx = 0; phaseIdx < Traits::numPhases; ++phaseIdx) {
139  const Evaluation& S =
140  Opm::decay<Evaluation>(state.saturation(phaseIdx));
141  Valgrind::CheckDefined(S);
142 
143  values[phaseIdx] = Opm::max(Opm::min(S,1.0),0.0);
144  }
145  }
146 
150  template <class FluidState, class Evaluation = typename FluidState::Scalar>
151  static Evaluation pcnw(const Params& params, const FluidState& fs)
152  {
153  const Evaluation& Sw =
154  Opm::decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
155  Valgrind::CheckDefined(Sw);
156 
157  const Evaluation& wPhasePressure =
158  Sw*params.pcMaxSat(Traits::wettingPhaseIdx) +
159  (1.0 - Sw)*params.pcMinSat(Traits::wettingPhaseIdx);
160 
161  const Evaluation& Sn =
162  Opm::decay<Evaluation>(fs.saturation(Traits::nonWettingPhaseIdx));
163  Valgrind::CheckDefined(Sn);
164 
165  const Evaluation& nPhasePressure =
166  Sn*params.pcMaxSat(Traits::nonWettingPhaseIdx) +
167  (1.0 - Sn)*params.pcMinSat(Traits::nonWettingPhaseIdx);
168 
169  return nPhasePressure - wPhasePressure;
170  }
171 
172 
173  template <class Evaluation = Scalar>
174  static typename std::enable_if<Traits::numPhases == 2, Evaluation>::type
175  twoPhaseSatPcnw(const Params& params, const Evaluation& Sw)
176  {
177  const Evaluation& wPhasePressure =
178  Sw*params.pcMaxSat(Traits::wettingPhaseIdx) +
179  (1.0 - Sw)*params.pcMinSat(Traits::wettingPhaseIdx);
180 
181  const Evaluation& nPhasePressure =
182  (1.0 - Sw)*params.pcMaxSat(Traits::nonWettingPhaseIdx) +
183  Sw*params.pcMinSat(Traits::nonWettingPhaseIdx);
184 
185  return nPhasePressure - wPhasePressure;
186  }
187 
192  template <class FluidState, class Evaluation = typename FluidState::Scalar>
193  static Evaluation Sw(const Params& /*params*/, const FluidState& /*fs*/)
194  { OPM_THROW(std::runtime_error, "Not implemented: Sw()"); }
195 
196  template <class Evaluation = Scalar>
197  static typename std::enable_if<Traits::numPhases == 2, Evaluation>::type
198  twoPhaseSatSw(const Params& /*params*/, const Evaluation& /*Sw*/)
199  { OPM_THROW(std::runtime_error, "Not implemented: twoPhaseSatSw()"); }
200 
205  template <class FluidState, class Evaluation = typename FluidState::Scalar>
206  static Evaluation Sn(const Params& /*params*/, const FluidState& /*fs*/)
207  { OPM_THROW(std::runtime_error, "Not implemented: Sn()"); }
208 
209  template <class Evaluation = Scalar>
210  static typename std::enable_if<Traits::numPhases == 2, Evaluation>::type
211  twoPhaseSatSn(const Params& /*params*/, const Evaluation& /*Sw*/)
212  { OPM_THROW(std::runtime_error, "Not implemented: twoPhaseSatSn()"); }
213 
220  template <class FluidState, class Evaluation = typename FluidState::Scalar>
221  static typename std::enable_if<Traits::numPhases == 3, Evaluation>::type
222  Sg(const Params& /*params*/, const FluidState& /*fs*/)
223  { OPM_THROW(std::runtime_error, "Not implemented: Sg()"); }
224 
228  template <class FluidState, class Evaluation = typename FluidState::Scalar>
229  static Evaluation krw(const Params& /*params*/, const FluidState& fs)
230  {
231  const Evaluation& Sw =
232  Opm::decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
233  return Opm::max(0.0, Opm::min(1.0, Sw));
234  }
235 
236  template <class Evaluation = Scalar>
237  static typename std::enable_if<Traits::numPhases == 2, Evaluation>::type
238  twoPhaseSatKrw(const Params& /*params*/, const Evaluation& Sw)
239  { return Opm::max(0.0, Opm::min(1.0, Sw)); }
240 
244  template <class FluidState, class Evaluation = typename FluidState::Scalar>
245  static Evaluation krn(const Params& /*params*/, const FluidState& fs)
246  {
247  const Evaluation& Sn =
248  Opm::decay<Evaluation>(fs.saturation(Traits::nonWettingPhaseIdx));
249  return Opm::max(0.0, Opm::min(1.0, Sn));
250  }
251 
252  template <class Evaluation = Scalar>
253  static typename std::enable_if<Traits::numPhases == 2, Evaluation>::type
254  twoPhaseSatKrn(const Params& /*params*/, const Evaluation& Sw)
255  {
256  return Opm::max(0.0, Opm::min(1.0, Sw));
257  }
258 
264  template <class FluidState, class Evaluation = typename FluidState::Scalar>
265  static typename std::enable_if<Traits::numPhases == 3, Evaluation>::type
266  krg(const Params& /*params*/, const FluidState& fs)
267  {
268  const Evaluation& Sg =
269  Opm::decay<Evaluation>(fs.saturation(Traits::gasPhaseIdx));
270  return Opm::max(0.0, Opm::min(1.0, Sg));
271  }
272 
278  template <class FluidState, class Evaluation = Scalar>
279  static typename std::enable_if<Traits::numPhases == 3, Evaluation>::type
280  pcgn(const Params& params, const FluidState& fs)
281  {
282  const Evaluation& Sn =
283  Opm::decay<Evaluation>(fs.saturation(Traits::nonWettingPhaseIdx));
284  Valgrind::CheckDefined(Sn);
285 
286  const Evaluation& nPhasePressure =
287  Sn*params.pcMaxSat(Traits::nonWettingPhaseIdx) +
288  (1.0 - Sn)*params.pcMinSat(Traits::nonWettingPhaseIdx);
289 
290  const Evaluation& Sg =
291  Opm::decay<Evaluation>(fs.saturation(Traits::gasPhaseIdx));
292  Valgrind::CheckDefined(Sg);
293 
294  const Evaluation& gPhasePressure =
295  Sg*params.pcMaxSat(Traits::gasPhaseIdx) +
296  (1.0 - Sg)*params.pcMinSat(Traits::gasPhaseIdx);
297 
298  return gPhasePressure - nPhasePressure;
299  }
300 };
301 } // namespace Opm
302 
303 #endif
static const bool isCompositionDependent
Specify whether the quantities defined by this material law are dependent on the phase composition...
Definition: LinearMaterial.hpp:85
A traits class which provides basic mathematical functions for arbitrary scalar floating point values...
static const bool isSaturationDependent
Specify whether the quantities defined by this material law are saturation dependent.
Definition: LinearMaterial.hpp:73
static const bool implementsTwoPhaseSatApi
Specify whether this material law implements the two-phase convenience API which only depends on the ...
Definition: LinearMaterial.hpp:69
static const int numPhases
The number of fluid phases.
Definition: LinearMaterial.hpp:61
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: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:151
static Evaluation krw(const Params &, const FluidState &fs)
The relative permability of the wetting phase.
Definition: LinearMaterial.hpp:229
static void relativePermeabilities(ContainerT &values, const Params &, const FluidState &state)
The relative permeability of all phases.
Definition: LinearMaterial.hpp:132
static void saturations(ContainerT &, const Params &, const FluidState &)
The inverse of the capillary pressure.
Definition: LinearMaterial.hpp:121
Reference implementation of params for the linear M-phase material material.
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:193
static Evaluation krn(const Params &, const FluidState &fs)
The relative permability of the liquid non-wetting phase.
Definition: LinearMaterial.hpp:245
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:280
static void capillaryPressures(ContainerT &values, const Params &params, const FluidState &state)
The linear capillary pressure-saturation curve.
Definition: LinearMaterial.hpp:100
static const bool isTemperatureDependent
Specify whether the quantities defined by this material law are temperature dependent.
Definition: LinearMaterial.hpp:81
Implements a linear saturation-capillary pressure relation.
Definition: LinearMaterial.hpp:53
static const bool isPressureDependent
Specify whether the quantities defined by this material law are dependent on the absolute pressure...
Definition: LinearMaterial.hpp:77
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:222
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:206
static const bool implementsTwoPhaseApi
Specify whether this material law implements the two-phase convenience API.
Definition: LinearMaterial.hpp:65