EffToAbsLaw.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 
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_EFF_TO_ABS_LAW_HPP
26 #define OPM_EFF_TO_ABS_LAW_HPP
27 
28 #include "EffToAbsLawParams.hpp"
29 
31 
32 namespace Opm {
67 template <class EffLawT, class ParamsT = EffToAbsLawParams<typename EffLawT::Params, EffLawT::numPhases> >
68 class EffToAbsLaw : public EffLawT::Traits
69 {
70  typedef EffLawT EffLaw;
71 
72 public:
73  typedef typename EffLaw::Traits Traits;
74  typedef ParamsT Params;
75  typedef typename EffLaw::Scalar Scalar;
76 
78  static const int numPhases = EffLaw::numPhases;
79 
82  static const bool implementsTwoPhaseApi = EffLaw::implementsTwoPhaseApi;
83 
86  static const bool implementsTwoPhaseSatApi = EffLaw::implementsTwoPhaseSatApi;
87 
90  static const bool isSaturationDependent = EffLaw::isSaturationDependent;
91 
94  static const bool isPressureDependent = EffLaw::isPressureDependent;
95 
98  static const bool isTemperatureDependent = EffLaw::isTemperatureDependent;
99 
102  static const bool isCompositionDependent = EffLaw::isCompositionDependent;
103 
114  template <class Container, class FluidState>
115  static void capillaryPressures(Container &values, const Params &params, const FluidState &fs)
116  {
117  typedef Opm::SaturationOverlayFluidState<FluidState> OverlayFluidState;
118 
119  OverlayFluidState overlayFs(fs);
120  for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
121  overlayFs.setSaturation(phaseIdx,
122  effectiveSaturation(params,
123  fs.saturation(phaseIdx),
124  phaseIdx));
125  }
126 
127  EffLaw::template capillaryPressures<Container, OverlayFluidState>(values, params, overlayFs);
128  }
129 
140  template <class Container, class FluidState>
141  static void relativePermeabilities(Container &values, const Params &params, const FluidState &fs)
142  {
143  typedef Opm::SaturationOverlayFluidState<FluidState> OverlayFluidState;
144 
145  OverlayFluidState overlayFs(fs);
146  for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
147  overlayFs.setSaturation(phaseIdx,
148  effectiveSaturation(params,
149  fs.saturation(phaseIdx),
150  phaseIdx));
151  }
152 
153  EffLaw::template relativePermeabilities<Container, OverlayFluidState>(values, params, overlayFs);
154  }
155 
167  template <class FluidState, class Evaluation = typename FluidState::Scalar>
168  static Evaluation pcnw(const Params &params, const FluidState &fs)
169  {
170  typedef Opm::SaturationOverlayFluidState<FluidState> OverlayFluidState;
171 
172  static_assert(FluidState::numPhases == numPhases,
173  "The fluid state and the material law must exhibit the same "
174  "number of phases!");
175 
176  OverlayFluidState overlayFs(fs);
177  for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
178  overlayFs.setSaturation(phaseIdx,
179  effectiveSaturation(params,
180  fs.saturation(phaseIdx),
181  phaseIdx));
182  }
183 
184  return EffLaw::template pcnw<OverlayFluidState, Evaluation>(params, overlayFs);
185  }
186 
187  template <class Evaluation>
188  static typename std::enable_if<implementsTwoPhaseSatApi, Evaluation>::type
189  twoPhaseSatPcnw(const Params &params, const Evaluation& SwAbs)
190  {
191  const Evaluation& SwEff = effectiveSaturation(params, SwAbs, Traits::wettingPhaseIdx);
192 
193  return EffLaw::twoPhaseSatPcnw(params, SwEff);
194  }
195 
199  template <class Container, class FluidState>
200  static void saturations(Container &values, const Params &params, const FluidState &fs)
201  {
202  EffLaw::template saturations<Container, FluidState>(values, params, fs);
203  for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
204  values[phaseIdx] = absoluteSaturation(params, values[phaseIdx], phaseIdx);
205  }
206  }
207 
212  template <class FluidState, class Evaluation = typename FluidState::Scalar>
213  static Evaluation Sw(const Params &params, const FluidState &fs)
214  {
215  return absoluteSaturation(params,
216  EffLaw::template Sw<FluidState, Evaluation>(params, fs),
217  Traits::wettingPhaseIdx);
218  }
219 
220  template <class Evaluation>
221  static typename std::enable_if<implementsTwoPhaseSatApi, Evaluation>::type
222  twoPhaseSatSw(const Params &params, const Evaluation& Sw)
223  { return absoluteSaturation(params,
224  EffLaw::twoPhaseSatSw(params, Sw),
225  Traits::wettingPhaseIdx); }
226 
231  template <class FluidState, class Evaluation = typename FluidState::Scalar>
232  static Evaluation Sn(const Params &params, const FluidState &fs)
233  {
234  return absoluteSaturation(params,
235  EffLaw::template Sn<FluidState, Evaluation>(params, fs),
236  Traits::nonWettingPhaseIdx);
237  }
238 
239  template <class Evaluation>
240  static typename std::enable_if<implementsTwoPhaseSatApi, Evaluation>::type
241  twoPhaseSatSn(const Params &params, const Evaluation& Sw)
242  {
243  return absoluteSaturation(params,
244  EffLaw::twoPhaseSatSn(params, Sw),
245  Traits::nonWettingPhaseIdx);
246  }
247 
254  template <class FluidState, class Evaluation = typename FluidState::Scalar>
255  static typename std::enable_if< (Traits::numPhases > 2), Evaluation>::type
256  Sg(const Params &params, const FluidState &fs)
257  {
258  return absoluteSaturation(params,
259  EffLaw::template Sg<FluidState, Evaluation>(params, fs),
260  Traits::gasPhaseIdx);
261  }
262 
272  template <class FluidState, class Evaluation = typename FluidState::Scalar>
273  static Evaluation krw(const Params &params, const FluidState &fs)
274  {
275  typedef Opm::SaturationOverlayFluidState<FluidState> OverlayFluidState;
276 
277  static_assert(FluidState::numPhases == numPhases,
278  "The fluid state and the material law must exhibit the same "
279  "number of phases!");
280 
281  OverlayFluidState overlayFs(fs);
282  for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
283  overlayFs.setSaturation(phaseIdx,
284  effectiveSaturation(params,
285  fs.saturation(phaseIdx),
286  phaseIdx));
287  }
288 
289  return EffLaw::template krw<OverlayFluidState, Evaluation>(params, overlayFs);
290  }
291 
292  template <class Evaluation>
293  static typename std::enable_if<implementsTwoPhaseSatApi, Evaluation>::type
294  twoPhaseSatKrw(const Params &params, const Evaluation& Sw)
295  { return EffLaw::twoPhaseSatKrw(params, effectiveSaturation(params, Sw, Traits::nonWettingPhaseIdx)); }
296 
300  template <class FluidState, class Evaluation = typename FluidState::Scalar>
301  static Evaluation krn(const Params &params, const FluidState &fs)
302  {
303  typedef Opm::SaturationOverlayFluidState<FluidState> OverlayFluidState;
304 
305  static_assert(FluidState::numPhases == numPhases,
306  "The fluid state and the material law must exhibit the same "
307  "number of phases!");
308 
309  OverlayFluidState overlayFs(fs);
310  for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
311  overlayFs.setSaturation(phaseIdx,
312  effectiveSaturation(params,
313  fs.saturation(phaseIdx),
314  phaseIdx));
315  }
316 
317  return EffLaw::template krn<OverlayFluidState, Evaluation>(params, overlayFs);
318  }
319 
320  template <class Evaluation>
321  static typename std::enable_if<implementsTwoPhaseSatApi, Evaluation>::type
322  twoPhaseSatKrn(const Params &params, const Evaluation& Sw)
323  { return EffLaw::twoPhaseSatKrn(params, effectiveSaturation(params, Sw, Traits::nonWettingPhaseIdx)); }
324 
330  template <class FluidState, class Evaluation = typename FluidState::Scalar>
331  static typename std::enable_if< (Traits::numPhases > 2), Evaluation>::type
332  krg(const Params &params, const FluidState &fs)
333  {
334  typedef Opm::SaturationOverlayFluidState<FluidState> OverlayFluidState;
335 
336  static_assert(FluidState::numPhases == numPhases,
337  "The fluid state and the material law must exhibit the same "
338  "number of phases!");
339 
340  OverlayFluidState overlayFs(fs);
341  for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
342  overlayFs.setSaturation(phaseIdx,
343  effectiveSaturation(params,
344  fs.saturation(phaseIdx),
345  phaseIdx));
346  }
347 
348  return EffLaw::template krg<OverlayFluidState, Evaluation>(params, overlayFs);
349  }
350 
354  template <class Evaluation>
355  static Evaluation effectiveSaturation(const Params &params, const Evaluation& S, unsigned phaseIdx)
356  { return (S - params.residualSaturation(phaseIdx))/(1.0 - params.sumResidualSaturations()); }
357 
361  template <class Evaluation>
362  static Evaluation absoluteSaturation(const Params &params, const Evaluation& S, unsigned phaseIdx)
363  { return S*(1.0 - params.sumResidualSaturations()) + params.residualSaturation(phaseIdx); }
364 
365 private:
374  static Scalar dSeff_dSabs_(const Params &params, int /*phaseIdx*/)
375  { return 1.0/(1 - params.sumResidualSaturations()); }
376 
385  static Scalar dSabs_dSeff_(const Params &params, int /*phaseIdx*/)
386  { return 1 - params.sumResidualSaturations(); }
387 };
388 } // namespace Opm
389 
390 #endif
static void saturations(Container &values, const Params &params, const FluidState &fs)
The saturation-capillary pressure curves.
Definition: EffToAbsLaw.hpp:200
static Evaluation Sw(const Params &params, const FluidState &fs)
Calculate wetting liquid phase saturation given that the rest of the fluid state has been initialized...
Definition: EffToAbsLaw.hpp:213
static const bool isSaturationDependent
Definition: EffToAbsLaw.hpp:90
ParamsT Params
Definition: EffToAbsLaw.hpp:74
static std::enable_if< implementsTwoPhaseSatApi, Evaluation >::type twoPhaseSatKrw(const Params &params, const Evaluation &Sw)
Definition: EffToAbsLaw.hpp:294
static std::enable_if< (Traits::numPhases > 2), Evaluation >::type Sg(const Params &params, const FluidState &fs)
Calculate gas phase saturation given that the rest of the fluid state has been initialized.
Definition: EffToAbsLaw.hpp:256
static Evaluation pcnw(const Params &params, const FluidState &fs)
The capillary pressure-saturation curve.
Definition: EffToAbsLaw.hpp:168
Definition: Air_Mesitylene.hpp:31
static Evaluation krw(const Params &params, const FluidState &fs)
The relative permeability for the wetting phase.
Definition: EffToAbsLaw.hpp:273
static const int numPhases
The number of fluid phases.
Definition: EffToAbsLaw.hpp:78
static void relativePermeabilities(Container &values, const Params &params, const FluidState &fs)
The relative permeability-saturation curves depending on absolute saturations.
Definition: EffToAbsLaw.hpp:141
static std::enable_if< implementsTwoPhaseSatApi, Evaluation >::type twoPhaseSatSw(const Params &params, const Evaluation &Sw)
Definition: EffToAbsLaw.hpp:222
static std::enable_if< (Traits::numPhases > 2), Evaluation >::type krg(const Params &params, const FluidState &fs)
The relative permability of the gas phase.
Definition: EffToAbsLaw.hpp:332
static std::enable_if< implementsTwoPhaseSatApi, Evaluation >::type twoPhaseSatSn(const Params &params, const Evaluation &Sw)
Definition: EffToAbsLaw.hpp:241
static Evaluation absoluteSaturation(const Params &params, const Evaluation &S, unsigned phaseIdx)
Convert an effective saturation to an absolute one.
Definition: EffToAbsLaw.hpp:362
static const bool implementsTwoPhaseApi
Definition: EffToAbsLaw.hpp:82
static const bool isCompositionDependent
Definition: EffToAbsLaw.hpp:102
EffLaw::Scalar Scalar
Definition: EffToAbsLaw.hpp:75
static Evaluation krn(const Params &params, const FluidState &fs)
The relative permeability of the non-wetting phase.
Definition: EffToAbsLaw.hpp:301
static const bool isTemperatureDependent
Definition: EffToAbsLaw.hpp:98
A default implementation of the parameters for the adapter class to convert material laws from effect...
static const bool implementsTwoPhaseSatApi
Definition: EffToAbsLaw.hpp:86
This is a fluid state which allows to set the fluid saturations and takes all other quantities from a...
Definition: SaturationOverlayFluidState.hpp:40
EffLaw::Traits Traits
Definition: EffToAbsLaw.hpp:73
static Evaluation effectiveSaturation(const Params &params, const Evaluation &S, unsigned phaseIdx)
Convert an absolute saturation to an effective one.
Definition: EffToAbsLaw.hpp:355
static const bool isPressureDependent
Definition: EffToAbsLaw.hpp:94
static std::enable_if< implementsTwoPhaseSatApi, Evaluation >::type twoPhaseSatKrn(const Params &params, const Evaluation &Sw)
Definition: EffToAbsLaw.hpp:322
This is a fluid state which allows to set the fluid saturations and takes all other quantities from a...
static Evaluation Sn(const Params &params, const FluidState &fs)
Calculate non-wetting liquid phase saturation given that the rest of the fluid state has been initial...
Definition: EffToAbsLaw.hpp:232
This material law takes a material law defined for effective saturations and converts it to a materia...
Definition: EffToAbsLaw.hpp:68
static std::enable_if< implementsTwoPhaseSatApi, Evaluation >::type twoPhaseSatPcnw(const Params &params, const Evaluation &SwAbs)
Definition: EffToAbsLaw.hpp:189
static void capillaryPressures(Container &values, const Params &params, const FluidState &fs)
The capillary pressure-saturation curves depending on absolute saturations.
Definition: EffToAbsLaw.hpp:115