My Project
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 
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;
103 
104  for (unsigned phaseIdx = 0; phaseIdx < Traits::numPhases; ++phaseIdx) {
105  const Evaluation& S =
106  decay<Evaluation>(state.saturation(phaseIdx));
107  Valgrind::CheckDefined(S);
108 
109  values[phaseIdx] =
110  S*params.pcMaxSat(phaseIdx) +
111  (1.0 - S)*params.pcMinSat(phaseIdx);
112  }
113  }
114 
118  template <class ContainerT, class FluidState>
119  static void saturations(ContainerT& /*values*/,
120  const Params& /*params*/,
121  const FluidState& /*state*/)
122  {
123  throw std::runtime_error("Not implemented: LinearMaterial::saturations()");
124  }
125 
129  template <class ContainerT, class FluidState>
130  static void relativePermeabilities(ContainerT& values,
131  const Params& /*params*/,
132  const FluidState& state)
133  {
134  typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
135 
136  for (unsigned phaseIdx = 0; phaseIdx < Traits::numPhases; ++phaseIdx) {
137  const Evaluation& S =
138  decay<Evaluation>(state.saturation(phaseIdx));
139  Valgrind::CheckDefined(S);
140 
141  values[phaseIdx] = max(min(S,1.0),0.0);
142  }
143  }
144 
148  template <class FluidState, class Evaluation = typename FluidState::Scalar>
149  static Evaluation pcnw(const Params& params, const FluidState& fs)
150  {
151  const Evaluation& Sw =
152  decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
153  Valgrind::CheckDefined(Sw);
154 
155  const Evaluation& wPhasePressure =
156  Sw*params.pcMaxSat(Traits::wettingPhaseIdx) +
157  (1.0 - Sw)*params.pcMinSat(Traits::wettingPhaseIdx);
158 
159  const Evaluation& Sn =
160  decay<Evaluation>(fs.saturation(Traits::nonWettingPhaseIdx));
161  Valgrind::CheckDefined(Sn);
162 
163  const Evaluation& nPhasePressure =
164  Sn*params.pcMaxSat(Traits::nonWettingPhaseIdx) +
165  (1.0 - Sn)*params.pcMinSat(Traits::nonWettingPhaseIdx);
166 
167  return nPhasePressure - wPhasePressure;
168  }
169 
170 
171  template <class Evaluation = Scalar>
172  static typename std::enable_if<Traits::numPhases == 2, Evaluation>::type
173  twoPhaseSatPcnw(const Params& params, const Evaluation& Sw)
174  {
175  const Evaluation& wPhasePressure =
176  Sw*params.pcMaxSat(Traits::wettingPhaseIdx) +
177  (1.0 - Sw)*params.pcMinSat(Traits::wettingPhaseIdx);
178 
179  const Evaluation& nPhasePressure =
180  (1.0 - Sw)*params.pcMaxSat(Traits::nonWettingPhaseIdx) +
181  Sw*params.pcMinSat(Traits::nonWettingPhaseIdx);
182 
183  return nPhasePressure - wPhasePressure;
184  }
185 
190  template <class FluidState, class Evaluation = typename FluidState::Scalar>
191  static Evaluation Sw(const Params& /*params*/, const FluidState& /*fs*/)
192  { throw std::runtime_error("Not implemented: Sw()"); }
193 
194  template <class Evaluation = Scalar>
195  static typename std::enable_if<Traits::numPhases == 2, Evaluation>::type
196  twoPhaseSatSw(const Params& /*params*/, const Evaluation& /*Sw*/)
197  { throw std::runtime_error("Not implemented: twoPhaseSatSw()"); }
198 
203  template <class FluidState, class Evaluation = typename FluidState::Scalar>
204  static Evaluation Sn(const Params& /*params*/, const FluidState& /*fs*/)
205  { throw std::runtime_error("Not implemented: Sn()"); }
206 
207  template <class Evaluation = Scalar>
208  static typename std::enable_if<Traits::numPhases == 2, Evaluation>::type
209  twoPhaseSatSn(const Params& /*params*/, const Evaluation& /*Sw*/)
210  { throw std::runtime_error("Not implemented: twoPhaseSatSn()"); }
211 
218  template <class FluidState, class Evaluation = typename FluidState::Scalar>
219  static typename std::enable_if<Traits::numPhases == 3, Evaluation>::type
220  Sg(const Params& /*params*/, const FluidState& /*fs*/)
221  { throw std::runtime_error("Not implemented: Sg()"); }
222 
226  template <class FluidState, class Evaluation = typename FluidState::Scalar>
227  static Evaluation krw(const Params& /*params*/, const FluidState& fs)
228  {
229  const Evaluation& Sw =
230  decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
231  return max(0.0, min(1.0, Sw));
232  }
233 
234  template <class Evaluation = Scalar>
235  static typename std::enable_if<Traits::numPhases == 2, Evaluation>::type
236  twoPhaseSatKrw(const Params& /*params*/, const Evaluation& Sw)
237  { return max(0.0, min(1.0, Sw)); }
238 
242  template <class FluidState, class Evaluation = typename FluidState::Scalar>
243  static Evaluation krn(const Params& /*params*/, const FluidState& fs)
244  {
245  const Evaluation& Sn =
246  decay<Evaluation>(fs.saturation(Traits::nonWettingPhaseIdx));
247  return max(0.0, min(1.0, Sn));
248  }
249 
250  template <class Evaluation = Scalar>
251  static typename std::enable_if<Traits::numPhases == 2, Evaluation>::type
252  twoPhaseSatKrn(const Params& /*params*/, const Evaluation& Sw)
253  {
254  return max(0.0, min(1.0, Sw));
255  }
256 
262  template <class FluidState, class Evaluation = typename FluidState::Scalar>
263  static typename std::enable_if<Traits::numPhases == 3, Evaluation>::type
264  krg(const Params& /*params*/, const FluidState& fs)
265  {
266  const Evaluation& Sg =
267  decay<Evaluation>(fs.saturation(Traits::gasPhaseIdx));
268  return max(0.0, min(1.0, Sg));
269  }
270 
276  template <class FluidState, class Evaluation = Scalar>
277  static typename std::enable_if<Traits::numPhases == 3, Evaluation>::type
278  pcgn(const Params& params, const FluidState& fs)
279  {
280  const Evaluation& Sn =
281  decay<Evaluation>(fs.saturation(Traits::nonWettingPhaseIdx));
282  Valgrind::CheckDefined(Sn);
283 
284  const Evaluation& nPhasePressure =
285  Sn*params.pcMaxSat(Traits::nonWettingPhaseIdx) +
286  (1.0 - Sn)*params.pcMinSat(Traits::nonWettingPhaseIdx);
287 
288  const Evaluation& Sg =
289  decay<Evaluation>(fs.saturation(Traits::gasPhaseIdx));
290  Valgrind::CheckDefined(Sg);
291 
292  const Evaluation& gPhasePressure =
293  Sg*params.pcMaxSat(Traits::gasPhaseIdx) +
294  (1.0 - Sg)*params.pcMinSat(Traits::gasPhaseIdx);
295 
296  return gPhasePressure - nPhasePressure;
297  }
298 };
299 } // namespace Opm
300 
301 #endif
Reference implementation of params for the linear M-phase material material.
A traits class which provides basic mathematical functions for arbitrary scalar floating point values...
Some templates to wrap the valgrind client request macros.
Implements a linear saturation-capillary pressure relation.
Definition: LinearMaterial.hpp:52
static Evaluation pcnw(const Params &params, const FluidState &fs)
The difference between the pressures of the non-wetting and wetting phase.
Definition: LinearMaterial.hpp:149
static Evaluation krn(const Params &, const FluidState &fs)
The relative permability of the liquid non-wetting phase.
Definition: LinearMaterial.hpp:243
static void relativePermeabilities(ContainerT &values, const Params &, const FluidState &state)
The relative permeability of all phases.
Definition: LinearMaterial.hpp:130
static const bool implementsTwoPhaseSatApi
Specify whether this material law implements the two-phase convenience API which only depends on the ...
Definition: LinearMaterial.hpp:67
static const bool isPressureDependent
Specify whether the quantities defined by this material law are dependent on the absolute pressure.
Definition: LinearMaterial.hpp:75
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:278
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:220
static const bool isSaturationDependent
Specify whether the quantities defined by this material law are saturation dependent.
Definition: LinearMaterial.hpp:71
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:264
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:191
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:204
static const bool isTemperatureDependent
Specify whether the quantities defined by this material law are temperature dependent.
Definition: LinearMaterial.hpp:79
static void saturations(ContainerT &, const Params &, const FluidState &)
The inverse of the capillary pressure.
Definition: LinearMaterial.hpp:119
static const bool isCompositionDependent
Specify whether the quantities defined by this material law are dependent on the phase composition.
Definition: LinearMaterial.hpp:83
static const bool implementsTwoPhaseApi
Specify whether this material law implements the two-phase convenience API.
Definition: LinearMaterial.hpp:63
static void capillaryPressures(ContainerT &values, const Params &params, const FluidState &state)
The linear capillary pressure-saturation curve.
Definition: LinearMaterial.hpp:98
static const int numPhases
The number of fluid phases.
Definition: LinearMaterial.hpp:59
static Evaluation krw(const Params &, const FluidState &fs)
The relative permability of the wetting phase.
Definition: LinearMaterial.hpp:227