EclStone2Material.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_STONE2_MATERIAL_HPP
26 #define OPM_ECL_STONE2_MATERIAL_HPP
27 
29 
32 
33 #include <opm/common/Exceptions.hpp>
34 #include <opm/common/ErrorMacros.hpp>
35 
36 #include <algorithm>
37 
38 namespace Opm {
39 
53 template <class TraitsT,
54  class GasOilMaterialLawT,
55  class OilWaterMaterialLawT,
56  class ParamsT = EclStone2MaterialParams<TraitsT,
57  typename GasOilMaterialLawT::Params,
58  typename OilWaterMaterialLawT::Params> >
59 class EclStone2Material : public TraitsT
60 {
61 public:
62  typedef GasOilMaterialLawT GasOilMaterialLaw;
63  typedef OilWaterMaterialLawT OilWaterMaterialLaw;
64 
65  // some safety checks
66  static_assert(TraitsT::numPhases == 3,
67  "The number of phases considered by this capillary pressure "
68  "law is always three!");
69  static_assert(GasOilMaterialLaw::numPhases == 2,
70  "The number of phases considered by the gas-oil capillary "
71  "pressure law must be two!");
72  static_assert(OilWaterMaterialLaw::numPhases == 2,
73  "The number of phases considered by the oil-water capillary "
74  "pressure law must be two!");
75  static_assert(std::is_same<typename GasOilMaterialLaw::Scalar,
76  typename OilWaterMaterialLaw::Scalar>::value,
77  "The two two-phase capillary pressure laws must use the same "
78  "type of floating point values.");
79 
80  static_assert(GasOilMaterialLaw::implementsTwoPhaseSatApi,
81  "The gas-oil material law must implement the two-phase saturation "
82  "only API to for the default Ecl capillary pressure law!");
83  static_assert(OilWaterMaterialLaw::implementsTwoPhaseSatApi,
84  "The oil-water material law must implement the two-phase saturation "
85  "only API to for the default Ecl capillary pressure law!");
86 
87  typedef TraitsT Traits;
88  typedef ParamsT Params;
89  typedef typename Traits::Scalar Scalar;
90 
91  static const int numPhases = 3;
92  static const int waterPhaseIdx = Traits::wettingPhaseIdx;
93  static const int oilPhaseIdx = Traits::nonWettingPhaseIdx;
94  static const int gasPhaseIdx = Traits::gasPhaseIdx;
95 
98  static const bool implementsTwoPhaseApi = false;
99 
102  static const bool implementsTwoPhaseSatApi = false;
103 
106  static const bool isSaturationDependent = true;
107 
110  static const bool isPressureDependent = false;
111 
114  static const bool isTemperatureDependent = false;
115 
118  static const bool isCompositionDependent = false;
119 
134  template <class ContainerT, class FluidState>
135  static void capillaryPressures(ContainerT &values,
136  const Params &params,
137  const FluidState &state)
138  {
139  typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
140  values[gasPhaseIdx] = pcgn<FluidState, Evaluation>(params, state);
141  values[oilPhaseIdx] = 0;
142  values[waterPhaseIdx] = - pcnw<FluidState, Evaluation>(params, state);
143  Valgrind::CheckDefined(values[gasPhaseIdx]);
144  Valgrind::CheckDefined(values[oilPhaseIdx]);
145  Valgrind::CheckDefined(values[waterPhaseIdx]);
146  }
147 
157  template <class FluidState, class Evaluation = typename FluidState::Scalar>
158  static Evaluation pcgn(const Params &params,
159  const FluidState &fs)
160  {
162 
163  const auto& Sw = 1.0 - FsToolbox::template toLhs<Evaluation>(fs.saturation(gasPhaseIdx));
164  return GasOilMaterialLaw::twoPhaseSatPcnw(params.gasOilParams(), Sw);
165  }
166 
176  template <class FluidState, class Evaluation = typename FluidState::Scalar>
177  static Evaluation pcnw(const Params &params,
178  const FluidState &fs)
179  {
180  typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
181 
182  const auto& Sw = FsToolbox::template toLhs<Evaluation>(fs.saturation(waterPhaseIdx));
184  const auto& result = OilWaterMaterialLaw::twoPhaseSatPcnw(params.oilWaterParams(), Sw);
185  Valgrind::CheckDefined(result);
186  return result;
187  }
188 
192  template <class ContainerT, class FluidState>
193  static void saturations(ContainerT& /*values */,
194  const Params& /* params */,
195  const FluidState& /* fluidState */)
196  {
197  OPM_THROW(std::logic_error, "Not implemented: saturations()");
198  }
199 
203  template <class FluidState, class Evaluation = typename FluidState::Scalar>
204  static Evaluation Sg(const Params& /* params */,
205  const FluidState& /* fluidState */)
206  {
207  OPM_THROW(std::logic_error, "Not implemented: Sg()");
208  }
209 
213  template <class FluidState, class Evaluation = typename FluidState::Scalar>
214  static Evaluation Sn(const Params& /* params */,
215  const FluidState& /* fluidState */)
216  {
217  OPM_THROW(std::logic_error, "Not implemented: Sn()");
218  }
219 
223  template <class FluidState, class Evaluation = typename FluidState::Scalar>
224  static Evaluation Sw(const Params& /* params */,
225  const FluidState& /* fluidState */)
226  {
227  OPM_THROW(std::logic_error, "Not implemented: Sw()");
228  }
229 
245  template <class ContainerT, class FluidState>
246  static void relativePermeabilities(ContainerT &values,
247  const Params &params,
248  const FluidState &fluidState)
249  {
250  typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
251 
252  values[waterPhaseIdx] = krw<FluidState, Evaluation>(params, fluidState);
253  values[oilPhaseIdx] = krn<FluidState, Evaluation>(params, fluidState);
254  values[gasPhaseIdx] = krg<FluidState, Evaluation>(params, fluidState);
255  }
256 
260  template <class FluidState, class Evaluation = typename FluidState::Scalar>
261  static Evaluation krg(const Params &params,
262  const FluidState &fluidState)
263  {
264  typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
265 
266  const Evaluation& Sw = 1 - FsToolbox::template toLhs<Evaluation>(fluidState.saturation(gasPhaseIdx));
267  return GasOilMaterialLaw::twoPhaseSatKrn(params.gasOilParams(), Sw);
268  }
269 
273  template <class FluidState, class Evaluation = typename FluidState::Scalar>
274  static Evaluation krw(const Params &params,
275  const FluidState &fluidState)
276  {
277  typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
278 
279  const Evaluation& Sw = FsToolbox::template toLhs<Evaluation>(fluidState.saturation(waterPhaseIdx));
280  return OilWaterMaterialLaw::twoPhaseSatKrw(params.oilWaterParams(), Sw);
281  }
282 
286  template <class FluidState, class Evaluation = typename FluidState::Scalar>
287  static Evaluation krn(const Params &params,
288  const FluidState &fluidState)
289  {
290  typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
291 
292  Scalar Swco = params.Swl();
293  const Evaluation& Sw = FsToolbox::template toLhs<Evaluation>(fluidState.saturation(waterPhaseIdx));
294  const Evaluation& Sg = FsToolbox::template toLhs<Evaluation>(fluidState.saturation(gasPhaseIdx));
295 
296  Scalar krocw = OilWaterMaterialLaw::twoPhaseSatKrn(params.oilWaterParams(), Swco);
297  Evaluation krow = OilWaterMaterialLaw::twoPhaseSatKrn(params.oilWaterParams(), Sw);
298  Evaluation krw = OilWaterMaterialLaw::twoPhaseSatKrw(params.oilWaterParams(), Sw);
299  Evaluation krg = GasOilMaterialLaw::twoPhaseSatKrn(params.gasOilParams(), 1 - Sg);
300  Evaluation krog = GasOilMaterialLaw::twoPhaseSatKrw(params.gasOilParams(), 1 - Sg);
301 
302  return krocw*((krow/krocw + krw)*(krog/krocw + krg) - krw - krg);
303  }
304 
312  template <class FluidState>
313  static void updateHysteresis(Params &params, const FluidState &fluidState)
314  {
315  typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
316 
317  Scalar Sw = FsToolbox::value(fluidState.saturation(waterPhaseIdx));
318  Scalar Sg = FsToolbox::value(fluidState.saturation(gasPhaseIdx));
319 
320  params.oilWaterParams().update(/*pcSw=*/Sw, /*krwSw=*/Sw, /*krnSw=*/Sw);
321  params.gasOilParams().update(/*pcSw=*/1 - Sg, /*krwSw=*/1 - Sg, /*krnSw=*/1 - Sg);
322  }
323 };
324 } // namespace Opm
325 
326 #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
Some templates to wrap the valgrind client request macros.
OilWaterMaterialLawT OilWaterMaterialLaw
Definition: EclStone2Material.hpp:63
Implements the second phase capillary pressure/relperm law suggested by Stone as used by the ECLipse ...
Definition: EclStone2Material.hpp:59
GasOilMaterialLawT GasOilMaterialLaw
Definition: EclStone2Material.hpp:62
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...