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  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_HYSTERESIS_TWO_PHASE_LAW_HPP
26 #define OPM_ECL_HYSTERESIS_TWO_PHASE_LAW_HPP
27 
29 
30 namespace Opm {
36 template <class EffectiveLawT,
37  class ParamsT = EclHysteresisTwoPhaseLawParams<EffectiveLawT> >
38 class EclHysteresisTwoPhaseLaw : public EffectiveLawT::Traits
39 {
40 public:
41  typedef EffectiveLawT EffectiveLaw;
42  typedef typename EffectiveLaw::Params EffectiveLawParams;
43 
44  typedef typename EffectiveLaw::Traits Traits;
45  typedef ParamsT Params;
46  typedef typename EffectiveLaw::Scalar Scalar;
47 
48  enum { wettingPhaseIdx = Traits::wettingPhaseIdx };
49  enum { nonWettingPhaseIdx = Traits::nonWettingPhaseIdx };
50 
52  static const int numPhases = EffectiveLaw::numPhases;
53  static_assert(numPhases == 2,
54  "The endpoint scaling applies to the nested twophase laws, not to "
55  "the threephase one!");
56 
59  static const bool implementsTwoPhaseApi = true;
60 
61  static_assert(EffectiveLaw::implementsTwoPhaseApi,
62  "The material laws put into EclEpsTwoPhaseLaw must implement the "
63  "two-phase material law API!");
64 
67  static const bool implementsTwoPhaseSatApi = true;
68 
69  static_assert(EffectiveLaw::implementsTwoPhaseSatApi,
70  "The material laws put into EclEpsTwoPhaseLaw must implement the "
71  "two-phase material law saturation API!");
72 
75  static const bool isSaturationDependent = true;
76 
79  static const bool isPressureDependent = false;
80 
83  static const bool isTemperatureDependent = false;
84 
87  static const bool isCompositionDependent = false;
88 
99  template <class Container, class FluidState>
100  static void capillaryPressures(Container& /* values */,
101  const Params& /* params */,
102  const FluidState& /* fs */)
103  {
104  OPM_THROW(NotImplemented,
105  "The capillaryPressures(fs) method is not yet implemented");
106  }
107 
118  template <class Container, class FluidState>
119  static void relativePermeabilities(Container& /* values */,
120  const Params& /* params */,
121  const FluidState& /* fs */)
122  {
123  OPM_THROW(NotImplemented,
124  "The pcnw(fs) method is not yet implemented");
125  }
126 
138  template <class FluidState, class Evaluation = typename FluidState::Scalar>
139  static Evaluation pcnw(const Params& /* params */,
140  const FluidState& /* fs */)
141  {
142  OPM_THROW(NotImplemented,
143  "The pcnw(fs) method is not yet implemented");
144  }
145 
146  template <class Evaluation>
147  static Evaluation twoPhaseSatPcnw(const Params &params, const Evaluation& Sw)
148  {
149  // TODO: capillary pressure hysteresis
150  return EffectiveLaw::twoPhaseSatPcnw(params.drainageParams(), Sw);
151 /*
152  if (!params.config().enableHysteresis() || params.config().pcHysteresisModel() < 0)
153  return EffectiveLaw::twoPhaseSatPcnw(params.drainageParams(), Sw);
154 
155  if (Sw < params.SwMdc())
156  return EffectiveLaw::twoPhaseSatPcnw(params.drainageParams(), Sw);
157 
158  const Evaluation& SwEff = Sw;
159 
160  //return EffectiveLaw::twoPhaseSatPcnw(params.imbibitionParams(), SwEff);
161  return EffectiveLaw::twoPhaseSatPcnw(params.drainageParams(), SwEff);
162 */
163  }
164 
168  template <class Container, class FluidState>
169  static void saturations(Container& /* values */,
170  const Params& /* params */,
171  const FluidState& /* fs */)
172  {
173  OPM_THROW(NotImplemented,
174  "The saturations(fs) method is not yet implemented");
175  }
176 
181  template <class FluidState, class Evaluation = typename FluidState::Scalar>
182  static Evaluation Sw(const Params& /* params */,
183  const FluidState& /* fs */)
184  {
185  OPM_THROW(NotImplemented,
186  "The Sw(fs) method is not yet implemented");
187  }
188 
189  template <class Evaluation>
190  static Evaluation twoPhaseSatSw(const Params& /* params */,
191  const Evaluation& /* pc */)
192  {
193  OPM_THROW(NotImplemented,
194  "The twoPhaseSatSw(pc) method is not yet implemented");
195  }
196 
201  template <class FluidState, class Evaluation = typename FluidState::Scalar>
202  static Evaluation Sn(const Params& /* params */,
203  const FluidState& /* fs */)
204  {
205  OPM_THROW(NotImplemented,
206  "The Sn(pc) method is not yet implemented");
207  }
208 
209  template <class Evaluation>
210  static Evaluation twoPhaseSatSn(const Params& /* params */,
211  const Evaluation& /* pc */)
212  {
213  OPM_THROW(NotImplemented,
214  "The twoPhaseSatSn(pc) method is not yet implemented");
215  }
216 
226  template <class FluidState, class Evaluation = typename FluidState::Scalar>
227  static Evaluation krw(const Params& /* params */,
228  const FluidState& /* fs */)
229  {
230  OPM_THROW(NotImplemented,
231  "The krw(fs) method is not yet implemented");
232  }
233 
234  template <class Evaluation>
235  static Evaluation twoPhaseSatKrw(const Params &params, const Evaluation& Sw)
236  {
237 
238  // if no relperm hysteresis is enabled, use the drainage curve
239  if (!params.config().enableHysteresis() || params.config().krHysteresisModel() < 0)
240  return EffectiveLaw::twoPhaseSatKrw(params.drainageParams(), Sw);
241 
242  // if it is enabled, use either the drainage or the imbibition curve. if the
243  // imbibition curve is used, the saturation must be shifted.
244  if (Sw <= params.krwSwMdc())
245  return EffectiveLaw::twoPhaseSatKrw(params.drainageParams(), Sw);
246 
247  return EffectiveLaw::twoPhaseSatKrw(params.imbibitionParams(),
248  Sw + params.deltaSwImbKrw());
249  }
250 
254  template <class FluidState, class Evaluation = typename FluidState::Scalar>
255  static Evaluation krn(const Params& /* params */,
256  const FluidState& /* fs */)
257  {
258  OPM_THROW(NotImplemented,
259  "The krn(fs) method is not yet implemented");
260  }
261 
262  template <class Evaluation>
263  static Evaluation twoPhaseSatKrn(const Params &params, const Evaluation& Sw)
264  {
265  // if no relperm hysteresis is enabled, use the drainage curve
266  if (!params.config().enableHysteresis() || params.config().krHysteresisModel() < 0)
267  return EffectiveLaw::twoPhaseSatKrn(params.drainageParams(), Sw);
268 
269  // if it is enabled, use either the drainage or the imbibition curve. if the
270  // imbibition curve is used, the saturation must be shifted.
271  if (Sw <= params.krnSwMdc())
272  return EffectiveLaw::twoPhaseSatKrn(params.drainageParams(), Sw);
273 
274  return EffectiveLaw::twoPhaseSatKrn(params.imbibitionParams(),
275  Sw + params.deltaSwImbKrn());
276  }
277 };
278 } // namespace Opm
279 
280 #endif
EffectiveLaw::Traits Traits
Definition: EclHysteresisTwoPhaseLaw.hpp:44
Definition: Air_Mesitylene.hpp:31
EffectiveLawT EffectiveLaw
Definition: EclHysteresisTwoPhaseLaw.hpp:41
Definition: EclHysteresisTwoPhaseLaw.hpp:48
EffectiveLaw::Params EffectiveLawParams
Definition: EclHysteresisTwoPhaseLaw.hpp:42
ParamsT Params
Definition: EclHysteresisTwoPhaseLaw.hpp:45
A default implementation of the parameters for the material law which implements the ECL relative per...
EffectiveLaw::Scalar Scalar
Definition: EclHysteresisTwoPhaseLaw.hpp:46
Definition: EclHysteresisTwoPhaseLaw.hpp:49
This material law implements the hysteresis model of the ECL file format.
Definition: EclHysteresisTwoPhaseLaw.hpp:38
static const int numPhases
The number of fluid phases.
Definition: EclHysteresisTwoPhaseLaw.hpp:52