EclStone1Material.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-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_STONE1_MATERIAL_HPP
26 #define OPM_ECL_STONE1_MATERIAL_HPP
27 
29 
32 
33 #include <opm/common/Exceptions.hpp>
34 #include <opm/common/ErrorMacros.hpp>
35 
36 #include <algorithm>
37 #include <cmath>
38 
39 namespace Opm {
40 
54 template <class TraitsT,
55  class GasOilMaterialLawT,
56  class OilWaterMaterialLawT,
57  class ParamsT = EclStone1MaterialParams<TraitsT,
58  typename GasOilMaterialLawT::Params,
59  typename OilWaterMaterialLawT::Params> >
60 class EclStone1Material : public TraitsT
61 {
62 public:
63  typedef GasOilMaterialLawT GasOilMaterialLaw;
64  typedef OilWaterMaterialLawT OilWaterMaterialLaw;
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  static_assert(GasOilMaterialLaw::implementsTwoPhaseSatApi,
82  "The gas-oil material law must implement the two-phase saturation "
83  "only API to for the default Ecl capillary pressure law!");
84  static_assert(OilWaterMaterialLaw::implementsTwoPhaseSatApi,
85  "The oil-water material law must implement the two-phase saturation "
86  "only API to for the default Ecl capillary pressure law!");
87 
88  typedef TraitsT Traits;
89  typedef ParamsT Params;
90  typedef typename Traits::Scalar Scalar;
91 
92  static const int numPhases = 3;
93  static const int waterPhaseIdx = Traits::wettingPhaseIdx;
94  static const int oilPhaseIdx = Traits::nonWettingPhaseIdx;
95  static const int gasPhaseIdx = Traits::gasPhaseIdx;
96 
99  static const bool implementsTwoPhaseApi = false;
100 
103  static const bool implementsTwoPhaseSatApi = false;
104 
107  static const bool isSaturationDependent = true;
108 
111  static const bool isPressureDependent = false;
112 
115  static const bool isTemperatureDependent = false;
116 
119  static const bool isCompositionDependent = false;
120 
135  template <class ContainerT, class FluidState>
136  static void capillaryPressures(ContainerT &values,
137  const Params &params,
138  const FluidState &state)
139  {
140  typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
141  values[gasPhaseIdx] = pcgn<FluidState, Evaluation>(params, state);
142  values[oilPhaseIdx] = 0;
143  values[waterPhaseIdx] = - pcnw<FluidState, Evaluation>(params, state);
144  Valgrind::CheckDefined(values[gasPhaseIdx]);
145  Valgrind::CheckDefined(values[oilPhaseIdx]);
146  Valgrind::CheckDefined(values[waterPhaseIdx]);
147  }
148 
158  template <class FluidState, class Evaluation = typename FluidState::Scalar>
159  static Evaluation pcgn(const Params &params,
160  const FluidState &fs)
161  {
163 
164  const auto& Sw = 1.0 - FsToolbox::template toLhs<Evaluation>(fs.saturation(gasPhaseIdx));
165  return GasOilMaterialLaw::twoPhaseSatPcnw(params.gasOilParams(), Sw);
166  }
167 
177  template <class FluidState, class Evaluation = typename FluidState::Scalar>
178  static Evaluation pcnw(const Params &params,
179  const FluidState &fs)
180  {
181  typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
182 
183  const auto& Sw = FsToolbox::template toLhs<Evaluation>(fs.saturation(waterPhaseIdx));
185  const auto& result = OilWaterMaterialLaw::twoPhaseSatPcnw(params.oilWaterParams(), Sw);
186  Valgrind::CheckDefined(result);
187  return result;
188  }
189 
193  template <class ContainerT, class FluidState>
194  static void saturations(ContainerT& /* values */,
195  const Params& /* params */,
196  const FluidState& /* fluidState */)
197  {
198  OPM_THROW(std::logic_error, "Not implemented: saturations()");
199  }
200 
204  template <class FluidState, class Evaluation = typename FluidState::Scalar>
205  static Evaluation Sg(const Params& /* params */,
206  const FluidState& /* fluidState */)
207  {
208  OPM_THROW(std::logic_error, "Not implemented: Sg()");
209  }
210 
214  template <class FluidState, class Evaluation = typename FluidState::Scalar>
215  static Evaluation Sn(const Params& /* params */,
216  const FluidState& /* fluidState */)
217  {
218  OPM_THROW(std::logic_error, "Not implemented: Sn()");
219  }
220 
224  template <class FluidState, class Evaluation = typename FluidState::Scalar>
225  static Evaluation Sw(const Params& /* params */,
226  const FluidState& /* fluidState */)
227  {
228  OPM_THROW(std::logic_error, "Not implemented: Sw()");
229  }
230 
246  template <class ContainerT, class FluidState>
247  static void relativePermeabilities(ContainerT &values,
248  const Params &params,
249  const FluidState &fluidState)
250  {
251  typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
252 
253  values[waterPhaseIdx] = krw<FluidState, Evaluation>(params, fluidState);
254  values[oilPhaseIdx] = krn<FluidState, Evaluation>(params, fluidState);
255  values[gasPhaseIdx] = krg<FluidState, Evaluation>(params, fluidState);
256  }
257 
261  template <class FluidState, class Evaluation = typename FluidState::Scalar>
262  static Evaluation krg(const Params &params,
263  const FluidState &fluidState)
264  {
265  typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
266 
267  const Evaluation& Sw = 1 - FsToolbox::template toLhs<Evaluation>(fluidState.saturation(gasPhaseIdx));
268  return GasOilMaterialLaw::twoPhaseSatKrn(params.gasOilParams(), Sw);
269  }
270 
274  template <class FluidState, class Evaluation = typename FluidState::Scalar>
275  static Evaluation krw(const Params &params,
276  const FluidState &fluidState)
277  {
278  typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
279 
280  const Evaluation& Sw = FsToolbox::template toLhs<Evaluation>(fluidState.saturation(waterPhaseIdx));
281  return OilWaterMaterialLaw::twoPhaseSatKrw(params.oilWaterParams(), Sw);
282  }
283 
287  template <class FluidState, class Evaluation = typename FluidState::Scalar>
288  static Evaluation krn(const Params &params,
289  const FluidState &fluidState)
290  {
291  typedef MathToolbox<Evaluation> Toolbox;
292  typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
293 
294  // the Eclipse docu is inconsistent here: In some places the connate water
295  // saturation is represented by "Swl", in others "Swco" is used.
296  Scalar Swco = params.Swl();
297 
298  Scalar Sowcr = params.Sowcr();
299  Scalar Sogcr = params.Sogcr();
300  Scalar Som = std::min(Sowcr, Sogcr); // minimum residual oil saturation
301 
302  Scalar eta = params.eta(); // exponent of the beta term
303 
304  const Evaluation& Sw = FsToolbox::template toLhs<Evaluation>(fluidState.saturation(waterPhaseIdx));
305  const Evaluation& So = FsToolbox::template toLhs<Evaluation>(fluidState.saturation(oilPhaseIdx));
306  const Evaluation& Sg = FsToolbox::template toLhs<Evaluation>(fluidState.saturation(gasPhaseIdx));
307 
308  Evaluation SSw;
309  if (Sw > Swco)
310  SSw = (Sw - Swco)/(1 - Swco - Som);
311  else
312  SSw = 0.0;
313 
314  Evaluation SSo;
315  if (So > Som)
316  SSo = (So - Som)/(1 - Swco - Som);
317  else
318  SSo = 0.0;
319 
320  Evaluation SSg = Sg/(1 - Swco - Som);
321 
322  Scalar krocw = OilWaterMaterialLaw::twoPhaseSatKrn(params.oilWaterParams(), Swco);
323  Evaluation krow = OilWaterMaterialLaw::twoPhaseSatKrn(params.oilWaterParams(), Sw);
324  Evaluation krog = GasOilMaterialLaw::twoPhaseSatKrw(params.gasOilParams(), 1 - Sg);
325 
326  Evaluation beta = Toolbox::pow(SSo/((1 - SSw)*(1 - SSg)), eta);
327  return beta*krow*krog/krocw;
328  }
329 
337  template <class FluidState>
338  static void updateHysteresis(Params &params, const FluidState &fluidState)
339  {
340  typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
341 
342  Scalar Sw = FsToolbox::value(fluidState.saturation(waterPhaseIdx));
343  Scalar Sg = FsToolbox::value(fluidState.saturation(gasPhaseIdx));
344 
345  params.oilWaterParams().update(/*pcSw=*/Sw, /*krwSw=*/Sw, /*krnSw=*/Sw);
346  params.gasOilParams().update(/*pcSw=*/1 - Sg, /*krwSw=*/1 - Sg, /*krnSw=*/1 - Sg);
347  }
348 };
349 } // namespace Opm
350 
351 #endif
bool CheckDefined(const T &value OPM_UNUSED)
Make valgrind complain if any of the memory occupied by an object is undefined.
Definition: Valgrind.hpp:74
Definition: MathToolbox.hpp:39
Definition: Air_Mesitylene.hpp:31
GasOilMaterialLawT GasOilMaterialLaw
Definition: EclStone1Material.hpp:63
OilWaterMaterialLawT OilWaterMaterialLaw
Definition: EclStone1Material.hpp:64
Some templates to wrap the valgrind client request macros.
Implements the second phase capillary pressure/relperm law suggested by Stone as used by the ECLipse ...
Definition: EclStone1Material.hpp:60
Evaluation< Scalar, VarSetTag, numVars > min(const Evaluation< Scalar, VarSetTag, numVars > &x1, const Evaluation< Scalar, VarSetTag, numVars > &x2)
Definition: Math.hpp:61
Evaluation< Scalar, VarSetTag, numVars > pow(const Evaluation< Scalar, VarSetTag, numVars > &base, Scalar exp)
Definition: Math.hpp:312
Default implementation for the parameters required by the three-phase capillary pressure/relperm Ston...
A traits class which provides basic mathematical functions for arbitrary scalar floating point values...