My Project
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  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_EFF_TO_ABS_LAW_HPP
28 #define OPM_EFF_TO_ABS_LAW_HPP
29 
30 #include "EffToAbsLawParams.hpp"
31 
33 
34 namespace Opm {
69 template <class EffLawT, class ParamsT = EffToAbsLawParams<typename EffLawT::Params, EffLawT::numPhases> >
70 class EffToAbsLaw : public EffLawT::Traits
71 {
72  typedef EffLawT EffLaw;
73 
74 public:
75  typedef typename EffLaw::Traits Traits;
76  typedef ParamsT Params;
77  typedef typename EffLaw::Scalar Scalar;
78 
80  static const int numPhases = EffLaw::numPhases;
81 
84  static const bool implementsTwoPhaseApi = EffLaw::implementsTwoPhaseApi;
85 
88  static const bool implementsTwoPhaseSatApi = EffLaw::implementsTwoPhaseSatApi;
89 
92  static const bool isSaturationDependent = EffLaw::isSaturationDependent;
93 
96  static const bool isPressureDependent = EffLaw::isPressureDependent;
97 
100  static const bool isTemperatureDependent = EffLaw::isTemperatureDependent;
101 
104  static const bool isCompositionDependent = EffLaw::isCompositionDependent;
105 
116  template <class Container, class FluidState>
117  static void capillaryPressures(Container& values, const Params& params, const FluidState& fs)
118  {
119  typedef SaturationOverlayFluidState<FluidState> OverlayFluidState;
120 
121  OverlayFluidState overlayFs(fs);
122  for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
123  overlayFs.setSaturation(phaseIdx,
124  effectiveSaturation(params,
125  fs.saturation(phaseIdx),
126  phaseIdx));
127  }
128 
129  EffLaw::template capillaryPressures<Container, OverlayFluidState>(values, params, overlayFs);
130  }
131 
142  template <class Container, class FluidState>
143  static void relativePermeabilities(Container& values, const Params& params, const FluidState& fs)
144  {
145  typedef SaturationOverlayFluidState<FluidState> OverlayFluidState;
146 
147  OverlayFluidState overlayFs(fs);
148  for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
149  overlayFs.setSaturation(phaseIdx,
150  effectiveSaturation(params,
151  fs.saturation(phaseIdx),
152  phaseIdx));
153  }
154 
155  EffLaw::template relativePermeabilities<Container, OverlayFluidState>(values, params, overlayFs);
156  }
157 
169  template <class FluidState, class Evaluation = typename FluidState::Scalar>
170  static Evaluation pcnw(const Params& params, const FluidState& fs)
171  {
172  typedef SaturationOverlayFluidState<FluidState> OverlayFluidState;
173 
174  static_assert(FluidState::numPhases == numPhases,
175  "The fluid state and the material law must exhibit the same "
176  "number of phases!");
177 
178  OverlayFluidState overlayFs(fs);
179  for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
180  overlayFs.setSaturation(phaseIdx,
181  effectiveSaturation(params,
182  fs.saturation(phaseIdx),
183  phaseIdx));
184  }
185 
186  return EffLaw::template pcnw<OverlayFluidState, Evaluation>(params, overlayFs);
187  }
188 
189  template <class Evaluation>
190  static typename std::enable_if<implementsTwoPhaseSatApi, Evaluation>::type
191  twoPhaseSatPcnw(const Params& params, const Evaluation& SwAbs)
192  {
193  const Evaluation& SwEff = effectiveSaturation(params, SwAbs, Traits::wettingPhaseIdx);
194 
195  return EffLaw::twoPhaseSatPcnw(params, SwEff);
196  }
197 
201  template <class Container, class FluidState>
202  static void saturations(Container& values, const Params& params, const FluidState& fs)
203  {
204  EffLaw::template saturations<Container, FluidState>(values, params, fs);
205  for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
206  values[phaseIdx] = absoluteSaturation(params, values[phaseIdx], phaseIdx);
207  }
208  }
209 
214  template <class FluidState, class Evaluation = typename FluidState::Scalar>
215  static Evaluation Sw(const Params& params, const FluidState& fs)
216  {
217  return absoluteSaturation(params,
218  EffLaw::template Sw<FluidState, Evaluation>(params, fs),
219  Traits::wettingPhaseIdx);
220  }
221 
222  template <class Evaluation>
223  static typename std::enable_if<implementsTwoPhaseSatApi, Evaluation>::type
224  twoPhaseSatSw(const Params& params, const Evaluation& Sw)
225  { return absoluteSaturation(params,
226  EffLaw::twoPhaseSatSw(params, Sw),
227  Traits::wettingPhaseIdx); }
228 
233  template <class FluidState, class Evaluation = typename FluidState::Scalar>
234  static Evaluation Sn(const Params& params, const FluidState& fs)
235  {
236  return absoluteSaturation(params,
237  EffLaw::template Sn<FluidState, Evaluation>(params, fs),
238  Traits::nonWettingPhaseIdx);
239  }
240 
241  template <class Evaluation>
242  static typename std::enable_if<implementsTwoPhaseSatApi, Evaluation>::type
243  twoPhaseSatSn(const Params& params, const Evaluation& Sw)
244  {
245  return absoluteSaturation(params,
246  EffLaw::twoPhaseSatSn(params, Sw),
247  Traits::nonWettingPhaseIdx);
248  }
249 
256  template <class FluidState, class Evaluation = typename FluidState::Scalar>
257  static typename std::enable_if< (Traits::numPhases > 2), Evaluation>::type
258  Sg(const Params& params, const FluidState& fs)
259  {
260  return absoluteSaturation(params,
261  EffLaw::template Sg<FluidState, Evaluation>(params, fs),
262  Traits::gasPhaseIdx);
263  }
264 
274  template <class FluidState, class Evaluation = typename FluidState::Scalar>
275  static Evaluation krw(const Params& params, const FluidState& fs)
276  {
277  typedef SaturationOverlayFluidState<FluidState> OverlayFluidState;
278 
279  static_assert(FluidState::numPhases == numPhases,
280  "The fluid state and the material law must exhibit the same "
281  "number of phases!");
282 
283  OverlayFluidState overlayFs(fs);
284  for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
285  overlayFs.setSaturation(phaseIdx,
286  effectiveSaturation(params,
287  fs.saturation(phaseIdx),
288  phaseIdx));
289  }
290 
291  return EffLaw::template krw<OverlayFluidState, Evaluation>(params, overlayFs);
292  }
293 
294  template <class Evaluation>
295  static typename std::enable_if<implementsTwoPhaseSatApi, Evaluation>::type
296  twoPhaseSatKrw(const Params& params, const Evaluation& Sw)
297  { return EffLaw::twoPhaseSatKrw(params, effectiveSaturation(params, Sw, Traits::wettingPhaseIdx)); }
298 
302  template <class FluidState, class Evaluation = typename FluidState::Scalar>
303  static Evaluation krn(const Params& params, const FluidState& fs)
304  {
305  typedef SaturationOverlayFluidState<FluidState> OverlayFluidState;
306 
307  static_assert(FluidState::numPhases == numPhases,
308  "The fluid state and the material law must exhibit the same "
309  "number of phases!");
310 
311  OverlayFluidState overlayFs(fs);
312  for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
313  overlayFs.setSaturation(phaseIdx,
314  effectiveSaturation(params,
315  fs.saturation(phaseIdx),
316  phaseIdx));
317  }
318 
319  return EffLaw::template krn<OverlayFluidState, Evaluation>(params, overlayFs);
320  }
321 
322  template <class Evaluation>
323  static typename std::enable_if<implementsTwoPhaseSatApi, Evaluation>::type
324  twoPhaseSatKrn(const Params& params, const Evaluation& Sw)
325  { return EffLaw::twoPhaseSatKrn(params, effectiveSaturation(params, Sw, Traits::wettingPhaseIdx)); }
326 
332  template <class FluidState, class Evaluation = typename FluidState::Scalar>
333  static typename std::enable_if< (Traits::numPhases > 2), Evaluation>::type
334  krg(const Params& params, const FluidState& fs)
335  {
336  typedef SaturationOverlayFluidState<FluidState> OverlayFluidState;
337 
338  static_assert(FluidState::numPhases == numPhases,
339  "The fluid state and the material law must exhibit the same "
340  "number of phases!");
341 
342  OverlayFluidState overlayFs(fs);
343  for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
344  overlayFs.setSaturation(phaseIdx,
345  effectiveSaturation(params,
346  fs.saturation(phaseIdx),
347  phaseIdx));
348  }
349 
350  return EffLaw::template krg<OverlayFluidState, Evaluation>(params, overlayFs);
351  }
352 
356  template <class Evaluation>
357  static Evaluation effectiveSaturation(const Params& params, const Evaluation& S, unsigned phaseIdx)
358  { return (S - params.residualSaturation(phaseIdx))/(1.0 - params.sumResidualSaturations()); }
359 
363  template <class Evaluation>
364  static Evaluation absoluteSaturation(const Params& params, const Evaluation& S, unsigned phaseIdx)
365  { return S*(1.0 - params.sumResidualSaturations()) + params.residualSaturation(phaseIdx); }
366 
367 private:
376  static Scalar dSeff_dSabs_(const Params& params, int /*phaseIdx*/)
377  { return 1.0/(1 - params.sumResidualSaturations()); }
378 
387  static Scalar dSabs_dSeff_(const Params& params, int /*phaseIdx*/)
388  { return 1 - params.sumResidualSaturations(); }
389 };
390 } // namespace Opm
391 
392 #endif
A default implementation of the parameters for the adapter class to convert material laws from effect...
This is a fluid state which allows to set the fluid saturations and takes all other quantities from a...
This material law takes a material law defined for effective saturations and converts it to a materia...
Definition: EffToAbsLaw.hpp:71
static const bool implementsTwoPhaseApi
Specify whether this material law implements the two-phase convenience API.
Definition: EffToAbsLaw.hpp:84
static const bool isPressureDependent
Specify whether the quantities defined by this material law are dependent on the absolute pressure.
Definition: EffToAbsLaw.hpp:96
static const bool isTemperatureDependent
Specify whether the quantities defined by this material law are temperature dependent.
Definition: EffToAbsLaw.hpp:100
static void capillaryPressures(Container &values, const Params &params, const FluidState &fs)
The capillary pressure-saturation curves depending on absolute saturations.
Definition: EffToAbsLaw.hpp:117
static const bool isCompositionDependent
Specify whether the quantities defined by this material law are dependent on the phase composition.
Definition: EffToAbsLaw.hpp:104
static Evaluation pcnw(const Params &params, const FluidState &fs)
The capillary pressure-saturation curve.
Definition: EffToAbsLaw.hpp:170
static Evaluation krn(const Params &params, const FluidState &fs)
The relative permeability of the non-wetting phase.
Definition: EffToAbsLaw.hpp:303
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:234
static const bool isSaturationDependent
Specify whether the quantities defined by this material law are saturation dependent.
Definition: EffToAbsLaw.hpp:92
static Evaluation krw(const Params &params, const FluidState &fs)
The relative permeability for the wetting phase.
Definition: EffToAbsLaw.hpp:275
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:258
static const bool implementsTwoPhaseSatApi
Specify whether this material law implements the two-phase convenience API which only depends on the ...
Definition: EffToAbsLaw.hpp:88
static Evaluation absoluteSaturation(const Params &params, const Evaluation &S, unsigned phaseIdx)
Convert an effective saturation to an absolute one.
Definition: EffToAbsLaw.hpp:364
static const int numPhases
The number of fluid phases.
Definition: EffToAbsLaw.hpp:80
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:334
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:215
static Evaluation effectiveSaturation(const Params &params, const Evaluation &S, unsigned phaseIdx)
Convert an absolute saturation to an effective one.
Definition: EffToAbsLaw.hpp:357
static void relativePermeabilities(Container &values, const Params &params, const FluidState &fs)
The relative permeability-saturation curves depending on absolute saturations.
Definition: EffToAbsLaw.hpp:143
static void saturations(Container &values, const Params &params, const FluidState &fs)
The saturation-capillary pressure curves.
Definition: EffToAbsLaw.hpp:202
This is a fluid state which allows to set the fluid saturations and takes all other quantities from a...
Definition: SaturationOverlayFluidState.hpp:44