Merge pull request #446 from bska/eval-2p-krox
Expose Two-Phase Oil Relperm As Individual Operation
This commit is contained in:
commit
819d8379ed
@ -355,26 +355,60 @@ public:
|
||||
const Evaluation Sg = Opm::decay<Evaluation>(fluidState.saturation(gasPhaseIdx));
|
||||
|
||||
const Evaluation Sw_ow = Sg + Sw;
|
||||
const Evaluation So_go = 1.0 - Sw_ow;
|
||||
const Evaluation kro_ow = OilWaterMaterialLaw::twoPhaseSatKrn(params.oilWaterParams(), Sw_ow);
|
||||
const Evaluation kro_go = GasOilMaterialLaw::twoPhaseSatKrw(params.gasOilParams(), So_go);
|
||||
const Evaluation kro_ow = relpermOilInOilWaterSystem<Evaluation>(params, fluidState);
|
||||
const Evaluation kro_go = relpermOilInOilGasSystem<Evaluation>(params, fluidState);
|
||||
|
||||
// avoid the division by zero: chose a regularized kro which is used if Sw - Swco
|
||||
// < epsilon/2 and interpolate between the oridinary and the regularized kro between
|
||||
// epsilon and epsilon/2
|
||||
const Scalar epsilon = 1e-5;
|
||||
if (Opm::scalarValue(Sw_ow) - Swco < epsilon) {
|
||||
const Evaluation kro2 = (kro_ow + kro_go)/2;;
|
||||
const Evaluation kro2 = (kro_ow + kro_go)/2;
|
||||
if (Opm::scalarValue(Sw_ow) - Swco > epsilon/2) {
|
||||
const Evaluation kro1 = (Sg*kro_go + (Sw - Swco)*kro_ow)/(Sw_ow - Swco);
|
||||
const Evaluation alpha = (epsilon - (Sw_ow - Swco))/(epsilon/2);
|
||||
|
||||
return kro2*alpha + kro1*(1 - alpha);
|
||||
}
|
||||
|
||||
return kro2;
|
||||
}
|
||||
else
|
||||
return (Sg*kro_go + (Sw - Swco)*kro_ow)/(Sw_ow - Swco);
|
||||
|
||||
return (Sg*kro_go + (Sw - Swco)*kro_ow) / (Sw_ow - Swco);
|
||||
}
|
||||
|
||||
/*!
|
||||
* \brief The relative permeability of oil in oil/gas system.
|
||||
*/
|
||||
template <class Evaluation, class FluidState>
|
||||
static Evaluation relpermOilInOilGasSystem(const Params& params,
|
||||
const FluidState& fluidState)
|
||||
{
|
||||
const Evaluation Sw =
|
||||
Opm::max(Evaluation{ params.Swl() },
|
||||
Opm::decay<Evaluation>(fluidState.saturation(waterPhaseIdx)));
|
||||
|
||||
const Evaluation Sg = Opm::decay<Evaluation>(fluidState.saturation(gasPhaseIdx));
|
||||
const Evaluation So_go = 1.0 - (Sg + Sw);
|
||||
|
||||
return GasOilMaterialLaw::twoPhaseSatKrw(params.gasOilParams(), So_go);
|
||||
}
|
||||
|
||||
/*!
|
||||
* \brief The relative permeability of oil in oil/water system.
|
||||
*/
|
||||
template <class Evaluation, class FluidState>
|
||||
static Evaluation relpermOilInOilWaterSystem(const Params& params,
|
||||
const FluidState& fluidState)
|
||||
{
|
||||
const Evaluation Sw =
|
||||
Opm::max(Evaluation{ params.Swl() },
|
||||
Opm::decay<Evaluation>(fluidState.saturation(waterPhaseIdx)));
|
||||
|
||||
const Evaluation Sg = Opm::decay<Evaluation>(fluidState.saturation(gasPhaseIdx));
|
||||
const Evaluation Sw_ow = Sg + Sw;
|
||||
|
||||
return OilWaterMaterialLaw::twoPhaseSatKrn(params.oilWaterParams(), Sw_ow);
|
||||
}
|
||||
|
||||
/*!
|
||||
|
@ -443,6 +443,66 @@ public:
|
||||
}
|
||||
}
|
||||
|
||||
/*!
|
||||
* \brief The relative permeability of oil in oil/gas system.
|
||||
*/
|
||||
template <class Evaluation, class FluidState>
|
||||
static Evaluation relpermOilInOilGasSystem(const Params& params,
|
||||
const FluidState& fluidState)
|
||||
{
|
||||
switch (params.approach()) {
|
||||
case EclMultiplexerApproach::EclStone1Approach:
|
||||
return Stone1Material::template relpermOilInOilGasSystem<Evaluation>
|
||||
(params.template getRealParams<EclMultiplexerApproach::EclStone1Approach>(),
|
||||
fluidState);
|
||||
|
||||
case EclMultiplexerApproach::EclStone2Approach:
|
||||
return Stone2Material::template relpermOilInOilGasSystem<Evaluation>
|
||||
(params.template getRealParams<EclMultiplexerApproach::EclStone2Approach>(),
|
||||
fluidState);
|
||||
|
||||
case EclMultiplexerApproach::EclDefaultApproach:
|
||||
return DefaultMaterial::template relpermOilInOilGasSystem<Evaluation>
|
||||
(params.template getRealParams<EclMultiplexerApproach::EclDefaultApproach>(),
|
||||
fluidState);
|
||||
|
||||
default:
|
||||
throw std::logic_error {
|
||||
"relpermOilInOilGasSystem() is specific to three phases"
|
||||
};
|
||||
}
|
||||
}
|
||||
|
||||
/*!
|
||||
* \brief The relative permeability of oil in oil/water system.
|
||||
*/
|
||||
template <class Evaluation, class FluidState>
|
||||
static Evaluation relpermOilInOilWaterSystem(const Params& params,
|
||||
const FluidState& fluidState)
|
||||
{
|
||||
switch (params.approach()) {
|
||||
case EclMultiplexerApproach::EclStone1Approach:
|
||||
return Stone1Material::template relpermOilInOilWaterSystem<Evaluation>
|
||||
(params.template getRealParams<EclMultiplexerApproach::EclStone1Approach>(),
|
||||
fluidState);
|
||||
|
||||
case EclMultiplexerApproach::EclStone2Approach:
|
||||
return Stone2Material::template relpermOilInOilWaterSystem<Evaluation>
|
||||
(params.template getRealParams<EclMultiplexerApproach::EclStone2Approach>(),
|
||||
fluidState);
|
||||
|
||||
case EclMultiplexerApproach::EclDefaultApproach:
|
||||
return DefaultMaterial::template relpermOilInOilWaterSystem<Evaluation>
|
||||
(params.template getRealParams<EclMultiplexerApproach::EclDefaultApproach>(),
|
||||
fluidState);
|
||||
|
||||
default:
|
||||
throw std::logic_error {
|
||||
"relpermOilInOilWaterSystem() is specific to three phases"
|
||||
};
|
||||
}
|
||||
}
|
||||
|
||||
/*!
|
||||
* \brief The relative permeability of the gas phase.
|
||||
*/
|
||||
|
@ -359,11 +359,11 @@ public:
|
||||
// oil relperm at connate water saturations (with Sg=0)
|
||||
const Scalar krocw = params.krocw();
|
||||
|
||||
const Evaluation& Sw = Opm::decay<Evaluation>(fluidState.saturation(waterPhaseIdx));
|
||||
const Evaluation& Sg = Opm::decay<Evaluation>(fluidState.saturation(gasPhaseIdx));
|
||||
const Evaluation Sw = Opm::decay<Evaluation>(fluidState.saturation(waterPhaseIdx));
|
||||
const Evaluation Sg = Opm::decay<Evaluation>(fluidState.saturation(gasPhaseIdx));
|
||||
|
||||
const Evaluation kro_ow = OilWaterMaterialLaw::twoPhaseSatKrn(params.oilWaterParams(), Sw);
|
||||
const Evaluation kro_go = GasOilMaterialLaw::twoPhaseSatKrw(params.gasOilParams(), 1 - Sg - Swco);
|
||||
const Evaluation kro_ow = relpermOilInOilWaterSystem<Evaluation>(params, fluidState);
|
||||
const Evaluation kro_go = relpermOilInOilGasSystem<Evaluation>(params, fluidState);
|
||||
|
||||
Evaluation beta;
|
||||
if (Sw <= Swco)
|
||||
@ -385,6 +385,30 @@ public:
|
||||
return Opm::max(0.0, Opm::min(1.0, beta*kro_ow*kro_go/krocw));
|
||||
}
|
||||
|
||||
/*!
|
||||
* \brief The relative permeability of oil in oil/gas system.
|
||||
*/
|
||||
template <class Evaluation, class FluidState>
|
||||
static Evaluation relpermOilInOilGasSystem(const Params& params,
|
||||
const FluidState& fluidState)
|
||||
{
|
||||
const Evaluation Sg = Opm::decay<Evaluation>(fluidState.saturation(gasPhaseIdx));
|
||||
|
||||
return GasOilMaterialLaw::twoPhaseSatKrw(params.gasOilParams(), 1 - Sg - params.Swl());
|
||||
}
|
||||
|
||||
/*!
|
||||
* \brief The relative permeability of oil in oil/water system.
|
||||
*/
|
||||
template <class Evaluation, class FluidState>
|
||||
static Evaluation relpermOilInOilWaterSystem(const Params& params,
|
||||
const FluidState& fluidState)
|
||||
{
|
||||
const Evaluation Sw = Opm::decay<Evaluation>(fluidState.saturation(waterPhaseIdx));
|
||||
|
||||
return OilWaterMaterialLaw::twoPhaseSatKrn(params.oilWaterParams(), Sw);
|
||||
}
|
||||
|
||||
/*!
|
||||
* \brief Update the hysteresis parameters after a time step.
|
||||
*
|
||||
|
@ -353,18 +353,46 @@ public:
|
||||
const FluidState& fluidState)
|
||||
{
|
||||
const Scalar Swco = params.Swl();
|
||||
|
||||
const Evaluation Sw = Opm::decay<Evaluation>(fluidState.saturation(waterPhaseIdx));
|
||||
const Evaluation Sg = Opm::decay<Evaluation>(fluidState.saturation(gasPhaseIdx));
|
||||
|
||||
const Scalar krocw = OilWaterMaterialLaw::twoPhaseSatKrn(params.oilWaterParams(), Swco);
|
||||
const Evaluation krow = OilWaterMaterialLaw::twoPhaseSatKrn(params.oilWaterParams(), Sw);
|
||||
const Evaluation krow = relpermOilInOilWaterSystem<Evaluation>(params, fluidState);
|
||||
const Evaluation krw = OilWaterMaterialLaw::twoPhaseSatKrw(params.oilWaterParams(), Sw);
|
||||
const Evaluation krg = GasOilMaterialLaw::twoPhaseSatKrn(params.gasOilParams(), 1 - Swco - Sg);
|
||||
const Evaluation krog = GasOilMaterialLaw::twoPhaseSatKrw(params.gasOilParams(), 1 - Swco - Sg);
|
||||
const Evaluation krog = relpermOilInOilGasSystem<Evaluation>(params, fluidState);
|
||||
|
||||
return Opm::max(krocw * ((krow/krocw + krw) * (krog/krocw + krg) - krw - krg), Evaluation{0});
|
||||
}
|
||||
|
||||
|
||||
/*!
|
||||
* \brief The relative permeability of oil in oil/gas system.
|
||||
*/
|
||||
template <class Evaluation, class FluidState>
|
||||
static Evaluation relpermOilInOilGasSystem(const Params& params,
|
||||
const FluidState& fluidState)
|
||||
{
|
||||
const Scalar Swco = params.Swl();
|
||||
const Evaluation Sg = Opm::decay<Evaluation>(fluidState.saturation(gasPhaseIdx));
|
||||
|
||||
return GasOilMaterialLaw::twoPhaseSatKrw(params.gasOilParams(), 1 - Swco - Sg);
|
||||
}
|
||||
|
||||
|
||||
/*!
|
||||
* \brief The relative permeability of oil in oil/water system.
|
||||
*/
|
||||
template <class Evaluation, class FluidState>
|
||||
static Evaluation relpermOilInOilWaterSystem(const Params& params,
|
||||
const FluidState& fluidState)
|
||||
{
|
||||
const Evaluation Sw = Opm::decay<Evaluation>(fluidState.saturation(waterPhaseIdx));
|
||||
|
||||
return OilWaterMaterialLaw::twoPhaseSatKrn(params.oilWaterParams(), Sw);
|
||||
}
|
||||
|
||||
/*!
|
||||
* \brief Update the hysteresis parameters after a time step.
|
||||
*
|
||||
|
Loading…
Reference in New Issue
Block a user