27 #ifndef OPM_ECL_DEFAULT_MATERIAL_HPP
28 #define OPM_ECL_DEFAULT_MATERIAL_HPP
37 #include <type_traits>
54 template <
class TraitsT,
55 class GasOilMaterialLawT,
56 class OilWaterMaterialLawT,
57 class ParamsT = EclDefaultMaterialParams<TraitsT,
58 typename GasOilMaterialLawT::Params,
59 typename OilWaterMaterialLawT::Params> >
63 using GasOilMaterialLaw = GasOilMaterialLawT;
64 using OilWaterMaterialLaw = OilWaterMaterialLawT;
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.");
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!");
88 using Traits = TraitsT;
89 using Params = ParamsT;
90 using Scalar =
typename Traits::Scalar;
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;
135 template <
class ContainerT,
class Flu
idState>
137 const Params& params,
138 const FluidState& state)
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);
145 Valgrind::CheckDefined(values[gasPhaseIdx]);
146 Valgrind::CheckDefined(values[oilPhaseIdx]);
147 Valgrind::CheckDefined(values[waterPhaseIdx]);
156 static void oilWaterHysteresisParams(Scalar& pcSwMdc,
158 const Params& params)
160 pcSwMdc = params.oilWaterParams().pcSwMdc();
161 krnSwMdc = params.oilWaterParams().krnSwMdc();
163 Valgrind::CheckDefined(pcSwMdc);
164 Valgrind::CheckDefined(krnSwMdc);
173 static void setOilWaterHysteresisParams(
const Scalar& pcSwMdc,
174 const Scalar& krnSwMdc,
177 constexpr
const double krwSw = 2.0;
178 params.oilWaterParams().update(pcSwMdc, krwSw, krnSwMdc);
187 static void gasOilHysteresisParams(Scalar& pcSwMdc,
189 const Params& params)
191 const auto Swco = params.Swl();
195 pcSwMdc = std::min(params.gasOilParams().pcSwMdc() + Swco, Scalar{2.0});
196 krnSwMdc = std::min(params.gasOilParams().krnSwMdc() + Swco, Scalar{2.0});
198 Valgrind::CheckDefined(pcSwMdc);
199 Valgrind::CheckDefined(krnSwMdc);
208 static void setGasOilHysteresisParams(
const Scalar& pcSwMdc,
209 const Scalar& krnSwMdc,
213 const auto Swco = params.Swl();
214 constexpr
const double krwSw = 2.0;
215 params.gasOilParams().update(pcSwMdc - Swco, krwSw, krnSwMdc - Swco);
227 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
228 static Evaluation
pcgn(
const Params& params,
229 const FluidState& fs)
232 const auto Sw = 1.0 - params.Swl() - decay<Evaluation>(fs.saturation(gasPhaseIdx));
233 return GasOilMaterialLaw::twoPhaseSatPcnw(params.gasOilParams(),
Sw);
245 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
246 static Evaluation
pcnw(
const Params& params,
247 const FluidState& fs)
249 const auto Sw = decay<Evaluation>(fs.saturation(waterPhaseIdx));
250 return OilWaterMaterialLaw::twoPhaseSatPcnw(params.oilWaterParams(),
Sw);
256 template <
class ContainerT,
class Flu
idState>
261 throw std::logic_error(
"Not implemented: saturations()");
267 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
268 static Evaluation
Sg(
const Params& ,
271 throw std::logic_error(
"Not implemented: Sg()");
277 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
278 static Evaluation
Sn(
const Params& ,
281 throw std::logic_error(
"Not implemented: Sn()");
287 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
288 static Evaluation
Sw(
const Params& ,
291 throw std::logic_error(
"Not implemented: Sw()");
309 template <
class ContainerT,
class Flu
idState>
311 const Params& params,
312 const FluidState& fluidState)
314 using Evaluation =
typename std::remove_reference<decltype(values[0])>::type;
316 values[waterPhaseIdx] = krw<FluidState, Evaluation>(params, fluidState);
317 values[oilPhaseIdx] = krn<FluidState, Evaluation>(params, fluidState);
318 values[gasPhaseIdx] = krg<FluidState, Evaluation>(params, fluidState);
324 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
325 static Evaluation
krg(
const Params& params,
326 const FluidState& fluidState)
329 const Evaluation
Sw = 1.0 - params.Swl() - decay<Evaluation>(fluidState.saturation(gasPhaseIdx));
330 return GasOilMaterialLaw::twoPhaseSatKrn(params.gasOilParams(),
Sw);
336 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
337 static Evaluation
krw(
const Params& params,
338 const FluidState& fluidState)
340 const Evaluation
Sw = decay<Evaluation>(fluidState.saturation(waterPhaseIdx));
341 return OilWaterMaterialLaw::twoPhaseSatKrw(params.oilWaterParams(),
Sw);
347 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
348 static Evaluation
krn(
const Params& params,
349 const FluidState& fluidState)
351 const Scalar Swco = params.Swl();
353 const Evaluation
Sw =
354 max(Evaluation(Swco),
355 decay<Evaluation>(fluidState.saturation(waterPhaseIdx)));
357 const Evaluation
Sg = decay<Evaluation>(fluidState.saturation(gasPhaseIdx));
359 const Evaluation Sw_ow =
Sg +
Sw;
360 const Evaluation kro_ow = relpermOilInOilWaterSystem<Evaluation>(params, fluidState);
361 const Evaluation kro_go = relpermOilInOilGasSystem<Evaluation>(params, fluidState);
366 constexpr
const Scalar epsilon = 1e-5;
367 if (scalarValue(Sw_ow) - Swco < epsilon) {
368 const Evaluation kro2 = (kro_ow + kro_go)/2;
369 if (scalarValue(Sw_ow) - Swco > epsilon/2) {
370 const Evaluation kro1 = (
Sg*kro_go + (
Sw - Swco)*kro_ow)/(Sw_ow - Swco);
371 const Evaluation alpha = (epsilon - (Sw_ow - Swco))/(epsilon/2);
373 return kro2*alpha + kro1*(1 - alpha);
379 return (
Sg*kro_go + (
Sw - Swco)*kro_ow) / (Sw_ow - Swco);
385 template <
class Evaluation,
class Flu
idState>
387 const FluidState& fluidState)
389 const Evaluation
Sw =
390 max(Evaluation{ params.Swl() },
391 decay<Evaluation>(fluidState.saturation(waterPhaseIdx)));
393 const Evaluation
Sg = decay<Evaluation>(fluidState.saturation(gasPhaseIdx));
394 const Evaluation So_go = 1.0 - (
Sg +
Sw);
396 return GasOilMaterialLaw::twoPhaseSatKrw(params.gasOilParams(), So_go);
402 template <
class Evaluation,
class Flu
idState>
404 const FluidState& fluidState)
406 const Evaluation
Sw =
407 max(Evaluation{ params.Swl() },
408 decay<Evaluation>(fluidState.saturation(waterPhaseIdx)));
410 const Evaluation
Sg = decay<Evaluation>(fluidState.saturation(gasPhaseIdx));
411 const Evaluation Sw_ow =
Sg +
Sw;
413 return OilWaterMaterialLaw::twoPhaseSatKrn(params.oilWaterParams(), Sw_ow);
423 template <
class Flu
idState>
426 const Scalar Swco = params.Swl();
428 const Scalar
Sw = clampSaturation(fluidState, waterPhaseIdx);
429 const Scalar So = clampSaturation(fluidState, oilPhaseIdx);
430 const Scalar
Sg = clampSaturation(fluidState, gasPhaseIdx);
432 if (params.inconsistentHysteresisUpdate()) {
442 params.oilWaterParams().update(
Sw,
446 params.gasOilParams().update( 1.0 - Swco -
Sg,
451 const Scalar Sw_ow =
Sg + std::max(Swco,
Sw);
452 const Scalar So_go = 1.0 - Sw_ow;
454 params.oilWaterParams().update(
Sw,
458 params.gasOilParams().update( 1.0 - Swco -
Sg,
464 template <
class Flu
idState>
465 static Scalar clampSaturation(
const FluidState& fluidState,
const int phaseIndex)
467 const auto sat = scalarValue(fluidState.saturation(phaseIndex));
468 return std::clamp(sat, Scalar{0.0}, Scalar{1.0});
Default implementation for the parameters required by the default three-phase capillary pressure mode...
Some templates to wrap the valgrind client request macros.
Implements the default three phase capillary pressure law used by the ECLipse simulator.
Definition: EclDefaultMaterial.hpp:61
static Evaluation relpermOilInOilWaterSystem(const Params ¶ms, const FluidState &fluidState)
The relative permeability of oil in oil/water system.
Definition: EclDefaultMaterial.hpp:403
static Evaluation krg(const Params ¶ms, const FluidState &fluidState)
The relative permeability of the gas phase.
Definition: EclDefaultMaterial.hpp:325
static Evaluation pcgn(const Params ¶ms, const FluidState &fs)
Capillary pressure between the gas and the non-wetting liquid (i.e., oil) phase.
Definition: EclDefaultMaterial.hpp:228
static constexpr bool isCompositionDependent
Specify whether the quantities defined by this material law are dependent on the phase composition.
Definition: EclDefaultMaterial.hpp:119
static constexpr bool isTemperatureDependent
Specify whether the quantities defined by this material law are temperature dependent.
Definition: EclDefaultMaterial.hpp:115
static Evaluation Sg(const Params &, const FluidState &)
The saturation of the gas phase.
Definition: EclDefaultMaterial.hpp:268
static Evaluation relpermOilInOilGasSystem(const Params ¶ms, const FluidState &fluidState)
The relative permeability of oil in oil/gas system.
Definition: EclDefaultMaterial.hpp:386
static constexpr bool implementsTwoPhaseApi
Specify whether this material law implements the two-phase convenience API.
Definition: EclDefaultMaterial.hpp:99
static void capillaryPressures(ContainerT &values, const Params ¶ms, const FluidState &state)
Implements the default three phase capillary pressure law used by the ECLipse simulator.
Definition: EclDefaultMaterial.hpp:136
static void saturations(ContainerT &, const Params &, const FluidState &)
The inverse of the capillary pressure.
Definition: EclDefaultMaterial.hpp:257
static Evaluation Sn(const Params &, const FluidState &)
The saturation of the non-wetting (i.e., oil) phase.
Definition: EclDefaultMaterial.hpp:278
static Evaluation krn(const Params ¶ms, const FluidState &fluidState)
The relative permeability of the non-wetting (i.e., oil) phase.
Definition: EclDefaultMaterial.hpp:348
static Evaluation pcnw(const Params ¶ms, const FluidState &fs)
Capillary pressure between the non-wetting liquid (i.e., oil) and the wetting liquid (i....
Definition: EclDefaultMaterial.hpp:246
static constexpr bool isSaturationDependent
Specify whether the quantities defined by this material law are saturation dependent.
Definition: EclDefaultMaterial.hpp:107
static void relativePermeabilities(ContainerT &values, const Params ¶ms, const FluidState &fluidState)
The relative permeability of all phases.
Definition: EclDefaultMaterial.hpp:310
static Evaluation krw(const Params ¶ms, const FluidState &fluidState)
The relative permeability of the wetting phase.
Definition: EclDefaultMaterial.hpp:337
static Evaluation Sw(const Params &, const FluidState &)
The saturation of the wetting (i.e., water) phase.
Definition: EclDefaultMaterial.hpp:288
static constexpr bool isPressureDependent
Specify whether the quantities defined by this material law are dependent on the absolute pressure.
Definition: EclDefaultMaterial.hpp:111
static constexpr bool implementsTwoPhaseSatApi
Specify whether this material law implements the two-phase convenience API which only depends on the ...
Definition: EclDefaultMaterial.hpp:103
static void updateHysteresis(Params ¶ms, const FluidState &fluidState)
Update the hysteresis parameters after a time step.
Definition: EclDefaultMaterial.hpp:424