EclMultiplexerMaterial.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) 2015 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_ECL_MULTIPLEXER_MATERIAL_HPP
26 #define OPM_ECL_MULTIPLEXER_MATERIAL_HPP
27 
29 #include "EclDefaultMaterial.hpp"
30 #include "EclStone1Material.hpp"
31 #include "EclStone2Material.hpp"
32 #include "EclTwoPhaseMaterial.hpp"
33 
36 #include <opm/common/Exceptions.hpp>
37 #include <opm/common/ErrorMacros.hpp>
38 
39 #include <algorithm>
40 
41 namespace Opm {
42 
49 template <class TraitsT,
50  class GasOilMaterialLawT,
51  class OilWaterMaterialLawT,
52  class ParamsT = EclMultiplexerMaterialParams<TraitsT,
53  GasOilMaterialLawT,
54  OilWaterMaterialLawT> >
55 class EclMultiplexerMaterial : public TraitsT
56 {
57 public:
58  typedef GasOilMaterialLawT GasOilMaterialLaw;
59  typedef OilWaterMaterialLawT OilWaterMaterialLaw;
60 
65 
66  // some safety checks
67  static_assert(TraitsT::numPhases == 3,
68  "The number of phases considered by this capillary pressure "
69  "law is always three!");
70  static_assert(GasOilMaterialLaw::numPhases == 2,
71  "The number of phases considered by the gas-oil capillary "
72  "pressure law must be two!");
73  static_assert(OilWaterMaterialLaw::numPhases == 2,
74  "The number of phases considered by the oil-water capillary "
75  "pressure law must be two!");
76  static_assert(std::is_same<typename GasOilMaterialLaw::Scalar,
77  typename OilWaterMaterialLaw::Scalar>::value,
78  "The two two-phase capillary pressure laws must use the same "
79  "type of floating point values.");
80 
81  typedef TraitsT Traits;
82  typedef ParamsT Params;
83  typedef typename Traits::Scalar Scalar;
84 
85  static const int numPhases = 3;
86  static const int waterPhaseIdx = Traits::wettingPhaseIdx;
87  static const int oilPhaseIdx = Traits::nonWettingPhaseIdx;
88  static const int gasPhaseIdx = Traits::gasPhaseIdx;
89 
92  static const bool implementsTwoPhaseApi = false;
93 
96  static const bool implementsTwoPhaseSatApi = false;
97 
100  static const bool isSaturationDependent = true;
101 
104  static const bool isPressureDependent = false;
105 
108  static const bool isTemperatureDependent = false;
109 
112  static const bool isCompositionDependent = false;
113 
128  template <class ContainerT, class FluidState>
129  static void capillaryPressures(ContainerT &values,
130  const Params &params,
131  const FluidState &fluidState)
132  {
133  switch (params.approach()) {
134  case EclStone1Approach:
135  Stone1Material::capillaryPressures(values,
136  params.template getRealParams<EclStone1Approach>(),
137  fluidState);
138  break;
139 
140  case EclStone2Approach:
141  Stone2Material::capillaryPressures(values,
142  params.template getRealParams<EclStone2Approach>(),
143  fluidState);
144  break;
145 
146  case EclDefaultApproach:
147  DefaultMaterial::capillaryPressures(values,
148  params.template getRealParams<EclDefaultApproach>(),
149  fluidState);
150  break;
151 
152  case EclTwoPhaseApproach:
153  TwoPhaseMaterial::capillaryPressures(values,
154  params.template getRealParams<EclTwoPhaseApproach>(),
155  fluidState);
156  break;
157  }
158  }
159 
169  template <class FluidState, class Evaluation = typename FluidState::Scalar>
170  static Evaluation pcgn(const Params& /* params */,
171  const FluidState& /* fs */)
172  {
173  OPM_THROW(std::logic_error, "Not implemented: pcgn()");
174  }
175 
185  template <class FluidState, class Evaluation = typename FluidState::Scalar>
186  static Evaluation pcnw(const Params& /* params */,
187  const FluidState& /* fs */)
188  {
189  OPM_THROW(std::logic_error, "Not implemented: pcnw()");
190  }
191 
195  template <class ContainerT, class FluidState>
196  static void saturations(ContainerT& /* values */,
197  const Params& /* params */,
198  const FluidState& /* fs */)
199  {
200  OPM_THROW(std::logic_error, "Not implemented: saturations()");
201  }
202 
206  template <class FluidState, class Evaluation = typename FluidState::Scalar>
207  static Evaluation Sg(const Params& /* params */,
208  const FluidState& /* fluidState */)
209  {
210  OPM_THROW(std::logic_error, "Not implemented: Sg()");
211  }
212 
216  template <class FluidState, class Evaluation = typename FluidState::Scalar>
217  static Evaluation Sn(const Params& /* params */,
218  const FluidState& /* fluidState */)
219  {
220  OPM_THROW(std::logic_error, "Not implemented: Sn()");
221  }
222 
226  template <class FluidState, class Evaluation = typename FluidState::Scalar>
227  static Evaluation Sw(const Params& /* params */,
228  const FluidState& /* fluidState */)
229  {
230  OPM_THROW(std::logic_error, "Not implemented: Sw()");
231  }
232 
248  template <class ContainerT, class FluidState>
249  static void relativePermeabilities(ContainerT &values,
250  const Params &params,
251  const FluidState &fluidState)
252  {
253  switch (params.approach()) {
254  case EclStone1Approach:
255  Stone1Material::relativePermeabilities(values,
256  params.template getRealParams<EclStone1Approach>(),
257  fluidState);
258  break;
259 
260  case EclStone2Approach:
261  Stone2Material::relativePermeabilities(values,
262  params.template getRealParams<EclStone2Approach>(),
263  fluidState);
264  break;
265 
266  case EclDefaultApproach:
267  DefaultMaterial::relativePermeabilities(values,
268  params.template getRealParams<EclDefaultApproach>(),
269  fluidState);
270  break;
271 
272  case EclTwoPhaseApproach:
273  TwoPhaseMaterial::relativePermeabilities(values,
274  params.template getRealParams<EclTwoPhaseApproach>(),
275  fluidState);
276  break;
277  }
278  }
279 
283  template <class FluidState, class Evaluation = typename FluidState::Scalar>
284  static Evaluation krg(const Params& /* params */,
285  const FluidState& /* fluidState */)
286  {
287  OPM_THROW(std::logic_error, "Not implemented: krg()");
288  }
289 
293  template <class FluidState, class Evaluation = typename FluidState::Scalar>
294  static Evaluation krw(const Params& /* params */,
295  const FluidState& /* fluidState */)
296  {
297  OPM_THROW(std::logic_error, "Not implemented: krw()");
298  }
299 
303  template <class FluidState, class Evaluation = typename FluidState::Scalar>
304  static Evaluation krn(const Params& /* params */,
305  const FluidState& /* fluidState */)
306  {
307  OPM_THROW(std::logic_error, "Not implemented: krn()");
308  }
309 
310 
318  template <class FluidState>
319  static void updateHysteresis(Params &params, const FluidState &fluidState)
320  {
321  switch (params.approach()) {
322  case EclStone1Approach:
323  Stone1Material::updateHysteresis(params.template getRealParams<EclStone1Approach>(),
324  fluidState);
325  break;
326 
327  case EclStone2Approach:
328  Stone2Material::updateHysteresis(params.template getRealParams<EclStone2Approach>(),
329  fluidState);
330  break;
331 
332  case EclDefaultApproach:
333  DefaultMaterial::updateHysteresis(params.template getRealParams<EclDefaultApproach>(),
334  fluidState);
335  break;
336 
337  case EclTwoPhaseApproach:
338  TwoPhaseMaterial::updateHysteresis(params.template getRealParams<EclTwoPhaseApproach>(),
339  fluidState);
340  break;
341  }
342  }
343 };
344 } // namespace Opm
345 
346 #endif
Implements the second phase capillary pressure/relperm law suggested by Stone as used by the ECLipse ...
Definition: Air_Mesitylene.hpp:31
Definition: EclMultiplexerMaterialParams.hpp:42
Some templates to wrap the valgrind client request macros.
Implements the default three phase capillary pressure law used by the ECLipse simulator.
Definition: EclDefaultMaterial.hpp:59
EclTwoPhaseApproach
Definition: EclTwoPhaseMaterialParams.hpp:33
Implements a multiplexer class that provides all three phase capillary pressure laws used by the ECLi...
Definition: EclMultiplexerMaterial.hpp:55
Opm::EclTwoPhaseMaterial< TraitsT, GasOilMaterialLaw, OilWaterMaterialLaw > TwoPhaseMaterial
Definition: EclMultiplexerMaterial.hpp:64
Multiplexer implementation for the parameters required by the multiplexed three-phase material law...
Implements a multiplexer class that provides ECL saturation functions for twophase simulations...
Definition: EclTwoPhaseMaterial.hpp:55
Definition: EclMultiplexerMaterialParams.hpp:41
Implements the second phase capillary pressure/relperm law suggested by Stone as used by the ECLipse ...
Implements the second phase capillary pressure/relperm law suggested by Stone as used by the ECLipse ...
Definition: EclStone1Material.hpp:60
Implements the second phase capillary pressure/relperm law suggested by Stone as used by the ECLipse ...
Definition: EclStone2Material.hpp:59
Implements a multiplexer class that provides ECL saturation functions for twophase simulations...
OilWaterMaterialLawT OilWaterMaterialLaw
Definition: EclMultiplexerMaterial.hpp:59
Opm::EclDefaultMaterial< TraitsT, GasOilMaterialLaw, OilWaterMaterialLaw > DefaultMaterial
Definition: EclMultiplexerMaterial.hpp:63
Opm::EclStone2Material< TraitsT, GasOilMaterialLaw, OilWaterMaterialLaw > Stone2Material
Definition: EclMultiplexerMaterial.hpp:62
Opm::EclStone1Material< TraitsT, GasOilMaterialLaw, OilWaterMaterialLaw > Stone1Material
Definition: EclMultiplexerMaterial.hpp:61
Definition: EclMultiplexerMaterialParams.hpp:40
GasOilMaterialLawT GasOilMaterialLaw
Definition: EclMultiplexerMaterial.hpp:58
Implements the default three phase capillary pressure law used by the ECLipse simulator.
A traits class which provides basic mathematical functions for arbitrary scalar floating point values...