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 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_STONE2_MATERIAL_HPP
28#define OPM_ECL_STONE2_MATERIAL_HPP
29
31
34
35#include <algorithm>
36#include <stdexcept>
37#include <type_traits>
38
39namespace Opm {
40
54template <class TraitsT,
55 class GasOilMaterialLawT,
56 class OilWaterMaterialLawT,
57 class ParamsT = EclStone2MaterialParams<TraitsT,
58 typename GasOilMaterialLawT::Params,
59 typename OilWaterMaterialLawT::Params> >
60class EclStone2Material : public TraitsT
61{
62public:
63 using GasOilMaterialLaw = GasOilMaterialLawT;
64 using OilWaterMaterialLaw = OilWaterMaterialLawT;
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 using Traits = TraitsT;
89 using Params = ParamsT;
90 using Scalar = typename Traits::Scalar;
91
92 static constexpr int numPhases = 3;
93 static constexpr int waterPhaseIdx = Traits::wettingPhaseIdx;
94 static constexpr int oilPhaseIdx = Traits::nonWettingPhaseIdx;
95 static constexpr int gasPhaseIdx = Traits::gasPhaseIdx;
96
99 static constexpr bool implementsTwoPhaseApi = false;
100
103 static constexpr bool implementsTwoPhaseSatApi = false;
104
107 static constexpr bool isSaturationDependent = true;
108
111 static constexpr bool isPressureDependent = false;
112
115 static constexpr bool isTemperatureDependent = false;
116
119 static constexpr 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 using Evaluation = typename std::remove_reference<decltype(values[0])>::type;
141 values[gasPhaseIdx] = pcgn<FluidState, Evaluation>(params, state);
142 values[oilPhaseIdx] = 0;
143 values[waterPhaseIdx] = - pcnw<FluidState, Evaluation>(params, state);
147 }
148
149 /*
150 * Hysteresis parameters for oil-water
151 * @see EclHysteresisTwoPhaseLawParams::pcSwMdc(...)
152 * @see EclHysteresisTwoPhaseLawParams::krnSwMdc(...)
153 * \param params Parameters
154 */
155 static void oilWaterHysteresisParams(Scalar& pcSwMdc,
156 Scalar& krnSwMdc,
157 const Params& params)
158 {
159 pcSwMdc = params.oilWaterParams().pcSwMdc();
160 krnSwMdc = params.oilWaterParams().krnSwMdc();
161
162 Valgrind::CheckDefined(pcSwMdc);
163 Valgrind::CheckDefined(krnSwMdc);
164 }
165
166 /*
167 * Hysteresis parameters for oil-water
168 * @see EclHysteresisTwoPhaseLawParams::pcSwMdc(...)
169 * @see EclHysteresisTwoPhaseLawParams::krnSwMdc(...)
170 * \param params Parameters
171 */
172 static void setOilWaterHysteresisParams(const Scalar& pcSwMdc,
173 const Scalar& krnSwMdc,
174 Params& params)
175 {
176 constexpr const double krwSw = 2.0; //Should not be used
177 params.oilWaterParams().update(pcSwMdc, krwSw, krnSwMdc);
178 }
179
180 /*
181 * Hysteresis parameters for gas-oil
182 * @see EclHysteresisTwoPhaseLawParams::pcSwMdc(...)
183 * @see EclHysteresisTwoPhaseLawParams::krnSwMdc(...)
184 * \param params Parameters
185 */
186 static void gasOilHysteresisParams(Scalar& pcSwMdc,
187 Scalar& krnSwMdc,
188 const Params& params)
189 {
190 const auto Swco = params.Swl();
191
192 // Pretend oil saturation is 'Swco' larger than it really is in
193 // order to infer correct maximum Sg values in output layer.
194 pcSwMdc = std::min(params.gasOilParams().pcSwMdc() + Swco, Scalar{2.0});
195 krnSwMdc = std::min(params.gasOilParams().krnSwMdc() + Swco, Scalar{2.0});
196
197 Valgrind::CheckDefined(pcSwMdc);
198 Valgrind::CheckDefined(krnSwMdc);
199 }
200
201 /*
202 * Hysteresis parameters for gas-oil
203 * @see EclHysteresisTwoPhaseLawParams::pcSwMdc(...)
204 * @see EclHysteresisTwoPhaseLawParams::krnSwMdc(...)
205 * \param params Parameters
206 */
207 static void setGasOilHysteresisParams(const Scalar& pcSwMdc,
208 const Scalar& krnSwMdc,
209 Params& params)
210 {
211 // Maximum attainable oil saturation is 1-SWL
212 const auto Swco = params.Swl();
213 constexpr const double krwSw = 2.0; //Should not be used
214 params.gasOilParams().update(pcSwMdc - Swco, krwSw, krnSwMdc - Swco);
215 }
216
226 template <class FluidState, class Evaluation = typename FluidState::Scalar>
227 static Evaluation pcgn(const Params& params,
228 const FluidState& fs)
229 {
230 // Maximum attainable oil saturation is 1-SWL.
231 const auto Sw = 1.0 - params.Swl() - decay<Evaluation>(fs.saturation(gasPhaseIdx));
232 return GasOilMaterialLaw::twoPhaseSatPcnw(params.gasOilParams(), Sw);
233 }
234
244 template <class FluidState, class Evaluation = typename FluidState::Scalar>
245 static Evaluation pcnw(const Params& params,
246 const FluidState& fs)
247 {
248 const auto Sw = decay<Evaluation>(fs.saturation(waterPhaseIdx));
250
251 const auto result = OilWaterMaterialLaw::twoPhaseSatPcnw(params.oilWaterParams(), Sw);
253
254 return result;
255 }
256
260 template <class ContainerT, class FluidState>
261 static void saturations(ContainerT& /*values */,
262 const Params& /* params */,
263 const FluidState& /* fluidState */)
264 {
265 throw std::logic_error("Not implemented: saturations()");
266 }
267
271 template <class FluidState, class Evaluation = typename FluidState::Scalar>
272 static Evaluation Sg(const Params& /* params */,
273 const FluidState& /* fluidState */)
274 {
275 throw std::logic_error("Not implemented: Sg()");
276 }
277
281 template <class FluidState, class Evaluation = typename FluidState::Scalar>
282 static Evaluation Sn(const Params& /* params */,
283 const FluidState& /* fluidState */)
284 {
285 throw std::logic_error("Not implemented: Sn()");
286 }
287
291 template <class FluidState, class Evaluation = typename FluidState::Scalar>
292 static Evaluation Sw(const Params& /* params */,
293 const FluidState& /* fluidState */)
294 {
295 throw std::logic_error("Not implemented: Sw()");
296 }
297
313 template <class ContainerT, class FluidState>
314 static void relativePermeabilities(ContainerT& values,
315 const Params& params,
316 const FluidState& fluidState)
317 {
318 using Evaluation = typename std::remove_reference<decltype(values[0])>::type;
319
320 values[waterPhaseIdx] = krw<FluidState, Evaluation>(params, fluidState);
321 values[oilPhaseIdx] = krn<FluidState, Evaluation>(params, fluidState);
322 values[gasPhaseIdx] = krg<FluidState, Evaluation>(params, fluidState);
323 }
324
328 template <class FluidState, class Evaluation = typename FluidState::Scalar>
329 static Evaluation krg(const Params& params,
330 const FluidState& fluidState)
331 {
332 // Maximum attainable oil saturation is 1-SWL.
333 const Evaluation Sw = 1 - params.Swl() - decay<Evaluation>(fluidState.saturation(gasPhaseIdx));
334 return GasOilMaterialLaw::twoPhaseSatKrn(params.gasOilParams(), Sw);
335 }
336
340 template <class FluidState, class Evaluation = typename FluidState::Scalar>
341 static Evaluation krw(const Params& params,
342 const FluidState& fluidState)
343 {
344 const Evaluation Sw = decay<Evaluation>(fluidState.saturation(waterPhaseIdx));
345 return OilWaterMaterialLaw::twoPhaseSatKrw(params.oilWaterParams(), Sw);
346 }
347
351 template <class FluidState, class Evaluation = typename FluidState::Scalar>
352 static Evaluation krn(const Params& params,
353 const FluidState& fluidState)
354 {
355 const Scalar Swco = params.Swl();
356
357 const Evaluation Sw = decay<Evaluation>(fluidState.saturation(waterPhaseIdx));
358 const Evaluation Sg = decay<Evaluation>(fluidState.saturation(gasPhaseIdx));
359
360 const Scalar krocw = OilWaterMaterialLaw::twoPhaseSatKrn(params.oilWaterParams(), Swco);
361 const Evaluation krow = relpermOilInOilWaterSystem<Evaluation>(params, fluidState);
362 const Evaluation krw = OilWaterMaterialLaw::twoPhaseSatKrw(params.oilWaterParams(), Sw);
363 const Evaluation krg = GasOilMaterialLaw::twoPhaseSatKrn(params.gasOilParams(), 1 - Swco - Sg);
364 const Evaluation krog = relpermOilInOilGasSystem<Evaluation>(params, fluidState);
365
366 return max(krocw * ((krow/krocw + krw) * (krog/krocw + krg) - krw - krg), Evaluation{0});
367 }
368
369
373 template <class Evaluation, class FluidState>
374 static Evaluation relpermOilInOilGasSystem(const Params& params,
375 const FluidState& fluidState)
376 {
377 const Scalar Swco = params.Swl();
378 const Evaluation Sg = decay<Evaluation>(fluidState.saturation(gasPhaseIdx));
379
380 return GasOilMaterialLaw::twoPhaseSatKrw(params.gasOilParams(), 1 - Swco - Sg);
381 }
382
383
387 template <class Evaluation, class FluidState>
388 static Evaluation relpermOilInOilWaterSystem(const Params& params,
389 const FluidState& fluidState)
390 {
391 const Evaluation Sw = decay<Evaluation>(fluidState.saturation(waterPhaseIdx));
392
393 return OilWaterMaterialLaw::twoPhaseSatKrn(params.oilWaterParams(), Sw);
394 }
395
403 template <class FluidState>
404 static void updateHysteresis(Params& params, const FluidState& fluidState)
405 {
406 const Scalar Swco = params.Swl();
407 const Scalar Sw = scalarValue(fluidState.saturation(waterPhaseIdx));
408 const Scalar Sg = scalarValue(fluidState.saturation(gasPhaseIdx));
409
410 params.oilWaterParams().update(/*pcSw=*/Sw, /*krwSw=*/Sw, /*krnSw=*/Sw);
411 params.gasOilParams().update(/*pcSw=*/ 1.0 - Swco - Sg,
412 /*krwSw=*/ 1.0 - Swco - Sg,
413 /*krnSw=*/ 1.0 - Swco - Sg);
414 }
415};
416
417} // namespace Opm
418
419#endif
A traits class which provides basic mathematical functions for arbitrary scalar floating point values...
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: EclStone2Material.hpp:61
static constexpr bool isPressureDependent
Definition: EclStone2Material.hpp:111
static Evaluation Sw(const Params &, const FluidState &)
The saturation of the wetting (i.e., water) phase.
Definition: EclStone2Material.hpp:292
TraitsT Traits
Definition: EclStone2Material.hpp:88
static Evaluation pcgn(const Params &params, const FluidState &fs)
Capillary pressure between the gas and the non-wetting liquid (i.e., oil) phase.
Definition: EclStone2Material.hpp:227
static Evaluation pcnw(const Params &params, const FluidState &fs)
Capillary pressure between the non-wetting liquid (i.e., oil) and the wetting liquid (i....
Definition: EclStone2Material.hpp:245
static void updateHysteresis(Params &params, const FluidState &fluidState)
Update the hysteresis parameters after a time step.
Definition: EclStone2Material.hpp:404
static void capillaryPressures(ContainerT &values, const Params &params, const FluidState &state)
Implements the default three phase capillary pressure law used by the ECLipse simulator.
Definition: EclStone2Material.hpp:136
static void setOilWaterHysteresisParams(const Scalar &pcSwMdc, const Scalar &krnSwMdc, Params &params)
Definition: EclStone2Material.hpp:172
static Evaluation relpermOilInOilWaterSystem(const Params &params, const FluidState &fluidState)
The relative permeability of oil in oil/water system.
Definition: EclStone2Material.hpp:388
static constexpr bool implementsTwoPhaseSatApi
Definition: EclStone2Material.hpp:103
static constexpr bool isTemperatureDependent
Definition: EclStone2Material.hpp:115
GasOilMaterialLawT GasOilMaterialLaw
Definition: EclStone2Material.hpp:63
typename Traits::Scalar Scalar
Definition: EclStone2Material.hpp:90
static void setGasOilHysteresisParams(const Scalar &pcSwMdc, const Scalar &krnSwMdc, Params &params)
Definition: EclStone2Material.hpp:207
static constexpr bool isCompositionDependent
Definition: EclStone2Material.hpp:119
static Evaluation Sg(const Params &, const FluidState &)
The saturation of the gas phase.
Definition: EclStone2Material.hpp:272
static void saturations(ContainerT &, const Params &, const FluidState &)
The inverse of the capillary pressure.
Definition: EclStone2Material.hpp:261
static constexpr bool implementsTwoPhaseApi
Definition: EclStone2Material.hpp:99
static Evaluation relpermOilInOilGasSystem(const Params &params, const FluidState &fluidState)
The relative permeability of oil in oil/gas system.
Definition: EclStone2Material.hpp:374
static constexpr int waterPhaseIdx
Definition: EclStone2Material.hpp:93
static constexpr int oilPhaseIdx
Definition: EclStone2Material.hpp:94
static Evaluation Sn(const Params &, const FluidState &)
The saturation of the non-wetting (i.e., oil) phase.
Definition: EclStone2Material.hpp:282
static Evaluation krn(const Params &params, const FluidState &fluidState)
The relative permeability of the non-wetting (i.e., oil) phase.
Definition: EclStone2Material.hpp:352
OilWaterMaterialLawT OilWaterMaterialLaw
Definition: EclStone2Material.hpp:64
static void oilWaterHysteresisParams(Scalar &pcSwMdc, Scalar &krnSwMdc, const Params &params)
Definition: EclStone2Material.hpp:155
ParamsT Params
Definition: EclStone2Material.hpp:89
static void gasOilHysteresisParams(Scalar &pcSwMdc, Scalar &krnSwMdc, const Params &params)
Definition: EclStone2Material.hpp:186
static constexpr int gasPhaseIdx
Definition: EclStone2Material.hpp:95
static constexpr bool isSaturationDependent
Definition: EclStone2Material.hpp:107
static Evaluation krw(const Params &params, const FluidState &fluidState)
The relative permeability of the wetting phase.
Definition: EclStone2Material.hpp:341
static Evaluation krg(const Params &params, const FluidState &fluidState)
The relative permeability of the gas phase.
Definition: EclStone2Material.hpp:329
static void relativePermeabilities(ContainerT &values, const Params &params, const FluidState &fluidState)
The relative permeability of all phases.
Definition: EclStone2Material.hpp:314
static constexpr int numPhases
Definition: EclStone2Material.hpp:92
bool CheckDefined(const T &value)
Make valgrind complain if any of the memory occupied by an object is undefined.
Definition: Valgrind.hpp:74
Definition: Air_Mesitylene.hpp:34
ReturnEval_< Evaluation1, Evaluation2 >::type min(const Evaluation1 &arg1, const Evaluation2 &arg2)
Definition: MathToolbox.hpp:346
ReturnEval_< Evaluation1, Evaluation2 >::type max(const Evaluation1 &arg1, const Evaluation2 &arg2)
Definition: MathToolbox.hpp:341
auto scalarValue(const Evaluation &val) -> decltype(MathToolbox< Evaluation >::scalarValue(val))
Definition: MathToolbox.hpp:335