27#ifndef OPM_ECL_STONE1_MATERIAL_HPP
28#define OPM_ECL_STONE1_MATERIAL_HPP
55template <
class TraitsT,
56 class GasOilMaterialLawT,
57 class OilWaterMaterialLawT,
58 class ParamsT = EclStone1MaterialParams<TraitsT, GasOilMaterialLawT, OilWaterMaterialLawT> >
62 using GasOilMaterialLaw = GasOilMaterialLawT;
63 using OilWaterMaterialLaw = OilWaterMaterialLawT;
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.");
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!");
87 using Traits = TraitsT;
88 using Params = ParamsT;
89 using Scalar =
typename Traits::Scalar;
91 static constexpr int numPhases = 3;
92 static constexpr int waterPhaseIdx = Traits::wettingPhaseIdx;
93 static constexpr int oilPhaseIdx = Traits::nonWettingPhaseIdx;
94 static constexpr int gasPhaseIdx = Traits::gasPhaseIdx;
134 template <
class ContainerT,
class Flu
idState>
136 const Params& params,
137 const FluidState& state)
139 using Evaluation =
typename std::remove_reference<
decltype(values[0])>::type;
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]);
155 static void oilWaterHysteresisParams(Scalar& soMax,
158 const Params& params)
160 soMax = 1.0 - params.oilWaterParams().krnSwMdc();
161 swMax = params.oilWaterParams().krwSwMdc();
162 swMin = params.oilWaterParams().pcSwMdc();
163 Valgrind::CheckDefined(soMax);
164 Valgrind::CheckDefined(swMax);
165 Valgrind::CheckDefined(swMin);
175 static void setOilWaterHysteresisParams(
const Scalar& soMax,
180 params.oilWaterParams().update(swMin, swMax, 1.0 - soMax);
190 static void gasOilHysteresisParams(Scalar& sgMax,
193 const Params& params)
195 const auto Swco = params.Swl();
196 sgMax = 1.0 - params.gasOilParams().krnSwMdc() - Swco;
197 shMax = params.gasOilParams().krwSwMdc();
198 soMin = params.gasOilParams().pcSwMdc();
200 Valgrind::CheckDefined(sgMax);
201 Valgrind::CheckDefined(shMax);
202 Valgrind::CheckDefined(soMin);
212 static void setGasOilHysteresisParams(
const Scalar& sgMax,
217 const auto Swco = params.Swl();
218 params.gasOilParams().update(soMin, shMax, 1.0 - sgMax - Swco);
221 static Scalar trappedGasSaturation(
const Params& params,
bool maximumTrapping)
223 const auto Swco = params.Swl();
224 return params.gasOilParams().SnTrapped(maximumTrapping) - Swco;
227 static Scalar trappedOilSaturation(
const Params& params,
bool maximumTrapping)
229 return params.oilWaterParams().SnTrapped(maximumTrapping) + params.gasOilParams().SwTrapped();
232 static Scalar trappedWaterSaturation(
const Params& params)
234 return params.oilWaterParams().SwTrapped();
237 static Scalar strandedGasSaturation(
const Params& params, Scalar
Sg, Scalar Kg)
239 const auto Swco = params.Swl();
240 return params.gasOilParams().SnStranded(
Sg, Kg) - Swco;
251 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
252 static Evaluation
pcgn(
const Params& params,
253 const FluidState& fs)
256 const auto Sw = 1.0 - params.Swl() - decay<Evaluation>(fs.saturation(gasPhaseIdx));
257 return GasOilMaterialLaw::twoPhaseSatPcnw(params.gasOilParams(),
Sw);
269 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
270 static Evaluation
pcnw(
const Params& params,
271 const FluidState& fs)
273 const auto Sw = decay<Evaluation>(fs.saturation(waterPhaseIdx));
274 Valgrind::CheckDefined(
Sw);
276 const auto result = OilWaterMaterialLaw::twoPhaseSatPcnw(params.oilWaterParams(),
Sw);
277 Valgrind::CheckDefined(result);
285 template <
class ContainerT,
class Flu
idState>
290 throw std::logic_error(
"Not implemented: saturations()");
296 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
297 static Evaluation
Sg(
const Params& ,
300 throw std::logic_error(
"Not implemented: Sg()");
306 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
307 static Evaluation
Sn(
const Params& ,
310 throw std::logic_error(
"Not implemented: Sn()");
316 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
317 static Evaluation
Sw(
const Params& ,
320 throw std::logic_error(
"Not implemented: Sw()");
338 template <
class ContainerT,
class Flu
idState>
340 const Params& params,
341 const FluidState& fluidState)
343 using Evaluation =
typename std::remove_reference<
decltype(values[0])>::type;
345 values[waterPhaseIdx] = krw<FluidState, Evaluation>(params, fluidState);
346 values[oilPhaseIdx] = krn<FluidState, Evaluation>(params, fluidState);
347 values[gasPhaseIdx] = krg<FluidState, Evaluation>(params, fluidState);
353 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
354 static Evaluation
krg(
const Params& params,
355 const FluidState& fluidState)
358 const Evaluation
Sw = 1 - params.Swl() - decay<Evaluation>(fluidState.saturation(gasPhaseIdx));
359 return GasOilMaterialLaw::twoPhaseSatKrn(params.gasOilParams(),
Sw);
365 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
366 static Evaluation
krw(
const Params& params,
367 const FluidState& fluidState)
369 const Evaluation
Sw = decay<Evaluation>(fluidState.saturation(waterPhaseIdx));
370 return OilWaterMaterialLaw::twoPhaseSatKrw(params.oilWaterParams(),
Sw);
376 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
377 static Evaluation
krn(
const Params& params,
378 const FluidState& fluidState)
383 const Scalar Swco = params.Swl();
386 const Scalar krocw = params.krocw();
388 const Evaluation
Sw = decay<Evaluation>(fluidState.saturation(waterPhaseIdx));
389 const Evaluation
Sg = decay<Evaluation>(fluidState.saturation(gasPhaseIdx));
391 const Evaluation kro_ow = relpermOilInOilWaterSystem<Evaluation>(params, fluidState);
392 const Evaluation kro_go = relpermOilInOilGasSystem<Evaluation>(params, fluidState);
401 const Evaluation SSw = (
Sw - Swco)/(1.0 - Swco);
402 const Evaluation SSg =
Sg/(1.0 - Swco);
403 const Evaluation SSo = 1.0 - SSw - SSg;
405 if (SSw >= 1.0 || SSg >= 1.0)
408 beta = pow( SSo/((1 - SSw)*(1 - SSg)), params.eta());
411 return max(0.0, min(1.0, beta*kro_ow*kro_go/krocw));
417 template <
class Evaluation,
class Flu
idState>
419 const FluidState& fluidState)
421 const Evaluation
Sg = decay<Evaluation>(fluidState.saturation(gasPhaseIdx));
423 return GasOilMaterialLaw::twoPhaseSatKrw(params.gasOilParams(), 1 -
Sg - params.Swl());
429 template <
class Evaluation,
class Flu
idState>
431 const FluidState& fluidState)
433 const Evaluation
Sw = decay<Evaluation>(fluidState.saturation(waterPhaseIdx));
435 return OilWaterMaterialLaw::twoPhaseSatKrn(params.oilWaterParams(),
Sw);
445 template <
class Flu
idState>
448 const Scalar Swco = params.Swl();
449 const Scalar
Sw = clampSaturation(fluidState, waterPhaseIdx);
450 const Scalar So = clampSaturation(fluidState, oilPhaseIdx);
451 const Scalar
Sg = clampSaturation(fluidState, gasPhaseIdx);
452 bool owChanged = params.oilWaterParams().update(
Sw,
Sw, 1 - So);
453 bool gochanged = params.gasOilParams().update( So,
456 return owChanged || gochanged;
459 template <
class Flu
idState>
460 static Scalar clampSaturation(
const FluidState& fluidState,
const int phaseIndex)
462 OPM_TIMEFUNCTION_LOCAL();
463 const auto sat = scalarValue(fluidState.saturation(phaseIndex));
464 return std::clamp(sat, Scalar{0.0}, Scalar{1.0});
Default implementation for the parameters required by the three-phase capillary pressure/relperm Ston...
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
static constexpr bool isSaturationDependent
Specify whether the quantities defined by this material law are saturation dependent.
Definition EclStone1Material.hpp:106
static constexpr bool implementsTwoPhaseSatApi
Specify whether this material law implements the two-phase convenience API which only depends on the ...
Definition EclStone1Material.hpp:102
static constexpr bool isCompositionDependent
Specify whether the quantities defined by this material law are dependent on the phase composition.
Definition EclStone1Material.hpp:118
static Evaluation krw(const Params ¶ms, const FluidState &fluidState)
The relative permeability of the wetting phase.
Definition EclStone1Material.hpp:366
static Evaluation krn(const Params ¶ms, const FluidState &fluidState)
The relative permeability of the non-wetting (i.e., oil) phase.
Definition EclStone1Material.hpp:377
static void relativePermeabilities(ContainerT &values, const Params ¶ms, const FluidState &fluidState)
The relative permeability of all phases.
Definition EclStone1Material.hpp:339
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 EclStone1Material.hpp:135
static void saturations(ContainerT &, const Params &, const FluidState &)
The inverse of the capillary pressure.
Definition EclStone1Material.hpp:286
static constexpr bool isTemperatureDependent
Specify whether the quantities defined by this material law are temperature dependent.
Definition EclStone1Material.hpp:114
static Evaluation Sg(const Params &, const FluidState &)
The saturation of the gas phase.
Definition EclStone1Material.hpp:297
static Evaluation krg(const Params ¶ms, const FluidState &fluidState)
The relative permeability of the gas phase.
Definition EclStone1Material.hpp:354
static constexpr bool implementsTwoPhaseApi
Specify whether this material law implements the two-phase convenience API.
Definition EclStone1Material.hpp:98
static Evaluation Sw(const Params &, const FluidState &)
The saturation of the wetting (i.e., water) phase.
Definition EclStone1Material.hpp:317
static Evaluation pcgn(const Params ¶ms, const FluidState &fs)
Capillary pressure between the gas and the non-wetting liquid (i.e., oil) phase.
Definition EclStone1Material.hpp:252
static constexpr bool isPressureDependent
Specify whether the quantities defined by this material law are dependent on the absolute pressure.
Definition EclStone1Material.hpp:110
static Evaluation relpermOilInOilWaterSystem(const Params ¶ms, const FluidState &fluidState)
The relative permeability of oil in oil/water system.
Definition EclStone1Material.hpp:430
static bool updateHysteresis(Params ¶ms, const FluidState &fluidState)
Update the hysteresis parameters after a time step.
Definition EclStone1Material.hpp:446
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 EclStone1Material.hpp:270
static Evaluation relpermOilInOilGasSystem(const Params ¶ms, const FluidState &fluidState)
The relative permeability of oil in oil/gas system.
Definition EclStone1Material.hpp:418
static Evaluation Sn(const Params &, const FluidState &)
The saturation of the non-wetting (i.e., oil) phase.
Definition EclStone1Material.hpp:307
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30