EclHysteresisTwoPhaseLaw.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_ECL_HYSTERESIS_TWO_PHASE_LAW_HPP
28 #define OPM_ECL_HYSTERESIS_TWO_PHASE_LAW_HPP
29 
31 
32 namespace Opm {
38 template <class EffectiveLawT,
39  class ParamsT = EclHysteresisTwoPhaseLawParams<EffectiveLawT> >
40 class EclHysteresisTwoPhaseLaw : public EffectiveLawT::Traits
41 {
42 public:
43  typedef EffectiveLawT EffectiveLaw;
44  typedef typename EffectiveLaw::Params EffectiveLawParams;
45 
46  typedef typename EffectiveLaw::Traits Traits;
47  typedef ParamsT Params;
48  typedef typename EffectiveLaw::Scalar Scalar;
49 
50  enum { wettingPhaseIdx = Traits::wettingPhaseIdx };
51  enum { nonWettingPhaseIdx = Traits::nonWettingPhaseIdx };
52 
54  static const int numPhases = EffectiveLaw::numPhases;
55  static_assert(numPhases == 2,
56  "The endpoint scaling applies to the nested twophase laws, not to "
57  "the threephase one!");
58 
61  static const bool implementsTwoPhaseApi = true;
62 
63  static_assert(EffectiveLaw::implementsTwoPhaseApi,
64  "The material laws put into EclEpsTwoPhaseLaw must implement the "
65  "two-phase material law API!");
66 
69  static const bool implementsTwoPhaseSatApi = true;
70 
71  static_assert(EffectiveLaw::implementsTwoPhaseSatApi,
72  "The material laws put into EclEpsTwoPhaseLaw must implement the "
73  "two-phase material law saturation API!");
74 
77  static const bool isSaturationDependent = true;
78 
81  static const bool isPressureDependent = false;
82 
85  static const bool isTemperatureDependent = false;
86 
89  static const bool isCompositionDependent = false;
90 
101  template <class Container, class FluidState>
102  static void capillaryPressures(Container& /* values */,
103  const Params& /* params */,
104  const FluidState& /* fs */)
105  {
106  OPM_THROW(NotImplemented,
107  "The capillaryPressures(fs) method is not yet implemented");
108  }
109 
120  template <class Container, class FluidState>
121  static void relativePermeabilities(Container& /* values */,
122  const Params& /* params */,
123  const FluidState& /* fs */)
124  {
125  OPM_THROW(NotImplemented,
126  "The pcnw(fs) method is not yet implemented");
127  }
128 
140  template <class FluidState, class Evaluation = typename FluidState::Scalar>
141  static Evaluation pcnw(const Params& /* params */,
142  const FluidState& /* fs */)
143  {
144  OPM_THROW(NotImplemented,
145  "The pcnw(fs) method is not yet implemented");
146  }
147 
148  template <class Evaluation>
149  static Evaluation twoPhaseSatPcnw(const Params& params, const Evaluation& Sw)
150  {
151  // TODO: capillary pressure hysteresis
152  return EffectiveLaw::twoPhaseSatPcnw(params.drainageParams(), Sw);
153 /*
154  if (!params.config().enableHysteresis() || params.config().pcHysteresisModel() < 0)
155  return EffectiveLaw::twoPhaseSatPcnw(params.drainageParams(), Sw);
156 
157  if (Sw < params.SwMdc())
158  return EffectiveLaw::twoPhaseSatPcnw(params.drainageParams(), Sw);
159 
160  const Evaluation& SwEff = Sw;
161 
162  //return EffectiveLaw::twoPhaseSatPcnw(params.imbibitionParams(), SwEff);
163  return EffectiveLaw::twoPhaseSatPcnw(params.drainageParams(), SwEff);
164 */
165  }
166 
170  template <class Container, class FluidState>
171  static void saturations(Container& /* values */,
172  const Params& /* params */,
173  const FluidState& /* fs */)
174  {
175  OPM_THROW(NotImplemented,
176  "The saturations(fs) method is not yet implemented");
177  }
178 
183  template <class FluidState, class Evaluation = typename FluidState::Scalar>
184  static Evaluation Sw(const Params& /* params */,
185  const FluidState& /* fs */)
186  {
187  OPM_THROW(NotImplemented,
188  "The Sw(fs) method is not yet implemented");
189  }
190 
191  template <class Evaluation>
192  static Evaluation twoPhaseSatSw(const Params& /* params */,
193  const Evaluation& /* pc */)
194  {
195  OPM_THROW(NotImplemented,
196  "The twoPhaseSatSw(pc) method is not yet implemented");
197  }
198 
203  template <class FluidState, class Evaluation = typename FluidState::Scalar>
204  static Evaluation Sn(const Params& /* params */,
205  const FluidState& /* fs */)
206  {
207  OPM_THROW(NotImplemented,
208  "The Sn(pc) method is not yet implemented");
209  }
210 
211  template <class Evaluation>
212  static Evaluation twoPhaseSatSn(const Params& /* params */,
213  const Evaluation& /* pc */)
214  {
215  OPM_THROW(NotImplemented,
216  "The twoPhaseSatSn(pc) method is not yet implemented");
217  }
218 
228  template <class FluidState, class Evaluation = typename FluidState::Scalar>
229  static Evaluation krw(const Params& /* params */,
230  const FluidState& /* fs */)
231  {
232  OPM_THROW(NotImplemented,
233  "The krw(fs) method is not yet implemented");
234  }
235 
236  template <class Evaluation>
237  static Evaluation twoPhaseSatKrw(const Params& params, const Evaluation& Sw)
238  {
239 
240  // if no relperm hysteresis is enabled, use the drainage curve
241  if (!params.config().enableHysteresis() || params.config().krHysteresisModel() < 0)
242  return EffectiveLaw::twoPhaseSatKrw(params.drainageParams(), Sw);
243 
244  // if it is enabled, use either the drainage or the imbibition curve. if the
245  // imbibition curve is used, the saturation must be shifted.
246  if (Sw <= params.krwSwMdc())
247  return EffectiveLaw::twoPhaseSatKrw(params.drainageParams(), Sw);
248 
249  return EffectiveLaw::twoPhaseSatKrw(params.imbibitionParams(),
250  Sw + params.deltaSwImbKrw());
251  }
252 
256  template <class FluidState, class Evaluation = typename FluidState::Scalar>
257  static Evaluation krn(const Params& /* params */,
258  const FluidState& /* fs */)
259  {
260  OPM_THROW(NotImplemented,
261  "The krn(fs) method is not yet implemented");
262  }
263 
264  template <class Evaluation>
265  static Evaluation twoPhaseSatKrn(const Params& params, const Evaluation& Sw)
266  {
267  // if no relperm hysteresis is enabled, use the drainage curve
268  if (!params.config().enableHysteresis() || params.config().krHysteresisModel() < 0)
269  return EffectiveLaw::twoPhaseSatKrn(params.drainageParams(), Sw);
270 
271  // if it is enabled, use either the drainage or the imbibition curve. if the
272  // imbibition curve is used, the saturation must be shifted.
273  if (Sw <= params.krnSwMdc())
274  return EffectiveLaw::twoPhaseSatKrn(params.drainageParams(), Sw);
275 
276  return EffectiveLaw::twoPhaseSatKrn(params.imbibitionParams(),
277  Sw + params.deltaSwImbKrn());
278  }
279 };
280 } // namespace Opm
281 
282 #endif
static Evaluation krw(const Params &, const FluidState &)
The relative permeability for the wetting phase.
Definition: EclHysteresisTwoPhaseLaw.hpp:229
This material law implements the hysteresis model of the ECL file format.
Definition: EclHysteresisTwoPhaseLaw.hpp:40
static const bool isCompositionDependent
Specify whether the quantities defined by this material law are dependent on the phase composition...
Definition: EclHysteresisTwoPhaseLaw.hpp:89
static Evaluation Sw(const Params &, const FluidState &)
Calculate wetting liquid phase saturation given that the rest of the fluid state has been initialized...
Definition: EclHysteresisTwoPhaseLaw.hpp:184
static void capillaryPressures(Container &, const Params &, const FluidState &)
The capillary pressure-saturation curves depending on absolute saturations.
Definition: EclHysteresisTwoPhaseLaw.hpp:102
static const bool implementsTwoPhaseSatApi
Specify whether this material law implements the two-phase convenience API which only depends on the ...
Definition: EclHysteresisTwoPhaseLaw.hpp:69
static void saturations(Container &, const Params &, const FluidState &)
The saturation-capillary pressure curves.
Definition: EclHysteresisTwoPhaseLaw.hpp:171
Definition: Air_Mesitylene.hpp:33
static const bool isTemperatureDependent
Specify whether the quantities defined by this material law are temperature dependent.
Definition: EclHysteresisTwoPhaseLaw.hpp:85
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: EclHysteresisTwoPhaseLaw.hpp:204
static const int numPhases
The number of fluid phases.
Definition: EclHysteresisTwoPhaseLaw.hpp:54
static Evaluation krn(const Params &, const FluidState &)
The relative permeability of the non-wetting phase.
Definition: EclHysteresisTwoPhaseLaw.hpp:257
static const bool isSaturationDependent
Specify whether the quantities defined by this material law are saturation dependent.
Definition: EclHysteresisTwoPhaseLaw.hpp:77
static void relativePermeabilities(Container &, const Params &, const FluidState &)
The relative permeability-saturation curves depending on absolute saturations.
Definition: EclHysteresisTwoPhaseLaw.hpp:121
static Evaluation pcnw(const Params &, const FluidState &)
The capillary pressure-saturation curve.
Definition: EclHysteresisTwoPhaseLaw.hpp:141
static const bool implementsTwoPhaseApi
Specify whether this material law implements the two-phase convenience API.
Definition: EclHysteresisTwoPhaseLaw.hpp:61
static const bool isPressureDependent
Specify whether the quantities defined by this material law are dependent on the absolute pressure...
Definition: EclHysteresisTwoPhaseLaw.hpp:81
A default implementation of the parameters for the material law which implements the ECL relative per...