SatFunc: Add Means of Retrieving Select Function Values

This commit introduces a new helper structure, RawFunctionValues,
which collects unscaled saturation function values that are needed
for vertical scaling of saturation functions using keywords such as

    KRO, KRORG, KRGR, PCW

and their hysteretic and directional counterparts.  We also
introduce a new helper function, getRawFunctionValues, which
extracts those values from the function tables in TableManager.

Add a set of unit tests to exercise the new feature.
This commit is contained in:
Bård Skaflestad 2020-09-16 00:13:58 +02:00
parent 07780a793f
commit df1f03e4a8
3 changed files with 285 additions and 0 deletions

View File

@ -30,31 +30,128 @@ namespace Opm {
namespace Opm { namespace satfunc {
/// Collection of unscaled/raw saturation range endpoints extracted
/// directly from tables of tabulated saturation functions.
struct RawTableEndPoints
{
/// Connate saturation endpoints
struct {
/// Connate gas saturation. One value for each saturation region.
/// All zero values unless gas is an active phase.
std::vector<double> gas;
/// Connate water saturation. One value for each saturation region.
/// All zero values unless water is an active phase.
std::vector<double> water;
} connate;
/// Critical saturation endpoints
struct {
/// Critical saturation of oil in oil/gas two-phase system. One
/// value for each saturation region. All zero values unless oil
/// is an active phase.
std::vector<double> oil_in_gas;
/// Critical saturation of oil in oil/water two-phase system. One
/// value for each saturation region. All zero values unless oil
/// is an active phase.
std::vector<double> oil_in_water;
/// Critical saturation of gas. One value for each saturation
/// region. All zero values unless oil is an active phase.
std::vector<double> gas;
/// Critical saturation of water. One value for each saturation
/// region. All zero values unless oil is an active phase.
std::vector<double> water;
} critical;
/// Maximum saturation endpoints
struct {
/// Maximum gas saturation value. One value for each saturation
/// region All zero values unless gas is an active phase.
std::vector<double> gas;
/// Maximum water saturation value. One value for each saturation
/// region All zero values unless gas is an active phase.
std::vector<double> water;
} maximum;
};
/// Collection of unscaled/raw saturation function value range endpoints
/// extracted directly from tables of tabulated saturation functions.
struct RawFunctionValues
{
/// Function values for relative permeability of oil.
struct {
/// Maximum relative permeability value of oil in both oil/gas and
/// oil/water two-phase systems. One value for each saturation
/// region. All zero values unless oil is an active phase.
std::vector<double> max;
/// Relative permeability of oil at critical saturation of
/// displacing phase in oil/gas two-phase system. One value for
/// each saturation region. All zero values unless oil is
/// an active phase.
std::vector<double> rg;
/// Relative permeability of oil at critical saturation of
/// displacing phase in oil/water two-phase system. One value
/// for each saturation region. All zero values unless oil is
/// an active phase.
std::vector<double> rw;
} kro;
/// Function values for relative permeability of gas.
struct {
/// Maximum relative permeability value of gas. One value for
/// each saturation region. All zero values unless gas is
/// an active phase.
std::vector<double> max;
/// Relative permeability of gas at critical saturation of
/// displacing phase. One value for each saturation region.
/// All zero values unless gas is an active phase.
std::vector<double> r;
} krg;
/// Function values for relative permeability of gas.
struct {
/// Maximum relative permeability value of water. One value for
/// each saturation region. All zero values unless water is
/// an active phase.
std::vector<double> max;
/// Relative permeability of water at critical saturation of
/// displacing phase. One value for each saturation region.
/// All zero values unless water is an active phase.
std::vector<double> r;
} krw;
/// Maximum capillary function values.
struct {
/// Maximum gas/oil capillary pressure value (Pg - Po). One
/// value for each saturation region. All zero values unless
/// both gas and oil are active phase.
std::vector<double> g;
/// Maximum oil/eater capillary pressure value (Po - Pw). One
/// value for each saturation region. All zero values unless
/// both oil and water are active phase.
std::vector<double> w;
} pc;
};
std::shared_ptr<RawTableEndPoints>
getRawTableEndpoints(const Opm::TableManager& tm,
const Opm::Phases& phases,
const double tolcrit);
std::shared_ptr<RawFunctionValues>
getRawFunctionValues(const Opm::TableManager& tm,
const Opm::Phases& phases,
const RawTableEndPoints& ep);
std::vector<double> init(const std::string& kewyord,
const TableManager& tables,
const Phases& phases,

View File

@ -1,4 +1,7 @@
/*
Copyright 2019-2020 Equinor ASA
Copyright 2014 Andreas Lauser
This file is part of the Open Porous Media project (OPM).
OPM is free software: you can redistribute it and/or modify it under the terms
@ -1595,6 +1598,29 @@ Opm::satfunc::getRawTableEndpoints(const Opm::TableManager& tm,
return ep;
}
std::shared_ptr<Opm::satfunc::RawFunctionValues>
Opm::satfunc::getRawFunctionValues(const Opm::TableManager& tm,
const Opm::Phases& phases,
const RawTableEndPoints& ep)
{
auto fval = std::make_shared<RawFunctionValues>();
fval->kro.max = findMaxKro(tm, phases);
fval->kro.rg = findKrorg(tm, phases, ep);
fval->kro.rw = findKrorw(tm, phases, ep);
fval->krg.max = findMaxKrg(tm, phases);
fval->krg.r = findKrgr(tm, phases, ep);
fval->krw.max = findMaxKrw(tm, phases);
fval->krw.r = findKrwr(tm, phases, ep);
fval->pc.g = findMaxPcog(tm, phases);
fval->pc.w = findMaxPcow(tm, phases);
return fval;
}
std::vector<double>
Opm::satfunc::init(const std::string& keyword,
const TableManager& tables,

View File

@ -1218,6 +1218,166 @@ BOOST_AUTO_TEST_CASE(RawTableEndPoints_Family_II_TolCrit_Large) {
// =====================================================================
BOOST_AUTO_TEST_CASE(RawFunctionValues_Family_I_Tolcrit_Zero) {
const auto es = ::Opm::EclipseState {
::Opm::Parser{}.parseString(satfunc_model_setup() + satfunc_family_I() + end())
};
const auto& tm = es.getTableManager();
const auto& ph = es.runspec().phases();
const auto tolcrit = 0.0;
const auto rtepPtr = satfunc::getRawTableEndpoints(tm, ph, tolcrit);
const auto rfuncPtr = satfunc::getRawFunctionValues(tm, ph, *rtepPtr);
// Oil values (TOLCRIT = 0.0)
BOOST_CHECK_CLOSE(rfuncPtr->kro.rw [0], 1.0, 1.0e-10); // Krow(So=1-Swcr-Sgl) = Krow(So=0.928996)
BOOST_CHECK_CLOSE(rfuncPtr->kro.rg [0], 0.896942, 1.0e-10); // Krog(So=1-Sgcr-Swl) = Krog(So=0.898996)
BOOST_CHECK_CLOSE(rfuncPtr->kro.max[0], 1.0, 1.0e-10); // Krow(So=Somax) = Krog(So=Somax)
// Gas values
BOOST_CHECK_CLOSE(rfuncPtr->krg.r [0], 0.866135, 1.0e-10); // Krg(Sg=1-Sogcr-Swl) = Krg(Sg=0.858996)
BOOST_CHECK_CLOSE(rfuncPtr->krg.max[0], 1.0 , 1.0e-10); // Krg(Sgmax) = Krg(Sg=0.928996)
// Water values
BOOST_CHECK_CLOSE(rfuncPtr->krw.r [0], 0.911429, 1.0e-10); // Krw(Sw=1-Sowcr-Sgl) = Krw(Sw=0.791004)
BOOST_CHECK_CLOSE(rfuncPtr->krw.max[0], 1.0 , 1.0e-10); // Krw(Swmax) = Krw(Sw=1)
}
BOOST_AUTO_TEST_CASE(RawFunctionValues_Family_II_Tolcrit_Zero) {
const auto es = ::Opm::EclipseState {
::Opm::Parser{}.parseString(satfunc_model_setup() + satfunc_family_II() + end())
};
const auto& tm = es.getTableManager();
const auto& ph = es.runspec().phases();
const auto tolcrit = 0.0;
const auto rtepPtr = satfunc::getRawTableEndpoints(tm, ph, tolcrit);
const auto rfuncPtr = satfunc::getRawFunctionValues(tm, ph, *rtepPtr);
// Oil values
BOOST_CHECK_CLOSE(rfuncPtr->kro.rw [0], 1.0, 1.0e-10); // Krow(So=1-Swcr-Sgl) = Krow(So=0.928996)
BOOST_CHECK_CLOSE(rfuncPtr->kro.rg [0], 0.896942, 1.0e-10); // Krog(So=1-Sgcr-Swl) = Krog(So=0.898996)
BOOST_CHECK_CLOSE(rfuncPtr->kro.max[0], 1.0, 1.0e-10); // Krow(So=Somax) = Krog(So=Somax)
// Gas values
BOOST_CHECK_CLOSE(rfuncPtr->krg.r [0], 0.866135, 1.0e-10); // Krg(Sg=1-Sogcr-Swl) = Krg(Sg=0.858996)
BOOST_CHECK_CLOSE(rfuncPtr->krg.max[0], 1.0 , 1.0e-10); // Krg(Sgmax) = Krg(Sg=0.928996)
// Water values
BOOST_CHECK_CLOSE(rfuncPtr->krw.r [0], 0.911429, 1.0e-10); // Krw(Sw=1-Sowcr-Sgl) = Krw(Sw=0.791004)
BOOST_CHECK_CLOSE(rfuncPtr->krw.max[0], 1.0 , 1.0e-10); // Krw(Swmax) = Krw(Sw=1)
}
BOOST_AUTO_TEST_CASE(RawFunctionValues_Family_I_Tolcrit_Default) {
const auto es = ::Opm::EclipseState {
::Opm::Parser{}.parseString(satfunc_model_setup() + satfunc_family_I() + end())
};
const auto& tm = es.getTableManager();
const auto& ph = es.runspec().phases();
const auto tolcrit = es.runspec().saturationFunctionControls()
.minimumRelpermMobilityThreshold(); // 1.0e-6.
const auto rtepPtr = satfunc::getRawTableEndpoints(tm, ph, tolcrit);
const auto rfuncPtr = satfunc::getRawFunctionValues(tm, ph, *rtepPtr);
// Oil values
BOOST_CHECK_CLOSE(rfuncPtr->kro.rw [0], 0.882459, 1.0e-10); // Krow(So=1-Swcr-Sgl) = Krow(So=0.908996)
BOOST_CHECK_CLOSE(rfuncPtr->kro.rg [0], 0.896942, 1.0e-10); // Krog(So=1-Sgcr-Swl) = Krog(So=0.898996)
BOOST_CHECK_CLOSE(rfuncPtr->kro.max[0], 1.0, 1.0e-10); // Krow(So=Somax) = Krog(So=Somax)
// Gas values
BOOST_CHECK_CLOSE(rfuncPtr->krg.r [0], 0.866135, 1.0e-10); // Krg(Sg=1-Sogcr-Swl) = Krg(Sg=0.858996)
BOOST_CHECK_CLOSE(rfuncPtr->krg.max[0], 1.0 , 1.0e-10); // Krg(Sgmax) = Krg(Sg=0.928996)
// Water values
BOOST_CHECK_CLOSE(rfuncPtr->krw.r [0], 0.835916, 1.0e-10); // Krw(Sw=1-Sowcr-Sgl) = Krw(Sw=0.771004)
BOOST_CHECK_CLOSE(rfuncPtr->krw.max[0], 1.0 , 1.0e-10); // Krw(Swmax) = Krw(Sw=1)
}
BOOST_AUTO_TEST_CASE(RawFunctionValues_Family_II_Tolcrit_Default) {
const auto es = ::Opm::EclipseState {
::Opm::Parser{}.parseString(satfunc_model_setup() + satfunc_family_II() + end())
};
const auto& tm = es.getTableManager();
const auto& ph = es.runspec().phases();
const auto tolcrit = es.runspec().saturationFunctionControls()
.minimumRelpermMobilityThreshold(); // 1.0e-6.
const auto rtepPtr = satfunc::getRawTableEndpoints(tm, ph, tolcrit);
const auto rfuncPtr = satfunc::getRawFunctionValues(tm, ph, *rtepPtr);
// Oil values
BOOST_CHECK_CLOSE(rfuncPtr->kro.rw [0], 0.882459, 1.0e-10); // Krow(So=1-Swcr-Sgl) = Krow(So=0.908996)
BOOST_CHECK_CLOSE(rfuncPtr->kro.rg [0], 0.896942, 1.0e-10); // Krog(So=1-Sgcr-Swl) = Krog(So=0.898996)
BOOST_CHECK_CLOSE(rfuncPtr->kro.max[0], 1.0 , 1.0e-10); // Krow(So=Somax) = Krog(So=Somax)
// Gas values
BOOST_CHECK_CLOSE(rfuncPtr->krg.r [0], 0.866135, 1.0e-10); // Krg(Sg=1-Sogcr-Swl) = Krg(Sg=0.858996)
BOOST_CHECK_CLOSE(rfuncPtr->krg.max[0], 1.0 , 1.0e-10); // Krg(Sgmax) = Krg(Sg=0.928996)
// Water values
BOOST_CHECK_CLOSE(rfuncPtr->krw.r [0], 0.835916, 1.0e-10); // Krw(Sw=1-Sowcr-Sgl) = Krw(Sw=0.771004)
BOOST_CHECK_CLOSE(rfuncPtr->krw.max[0], 1.0 , 1.0e-10); // Krw(Swmax) = Krw(Sw=1)
}
BOOST_AUTO_TEST_CASE(RawFunctionValues_Family_I_Tolcrit_Large) {
const auto es = ::Opm::EclipseState {
::Opm::Parser{}.parseString(satfunc_model_setup() + satfunc_family_I() + end())
};
const auto& tm = es.getTableManager();
const auto& ph = es.runspec().phases();
const auto tolcrit = 0.01;
const auto rtepPtr = satfunc::getRawTableEndpoints(tm, ph, tolcrit);
const auto rfuncPtr = satfunc::getRawFunctionValues(tm, ph, *rtepPtr);
// Oil values
BOOST_CHECK_CLOSE(rfuncPtr->kro.rw [0], 0.328347, 1.0e-10); // Krow(So=1-Swcr-Sgl) = Krow(So=0.768996)
BOOST_CHECK_CLOSE(rfuncPtr->kro.rg [0], 0.712749, 1.0e-10); // Krog(So=1-Sgcr-Swl) = Krog(So=0.838996)
BOOST_CHECK_CLOSE(rfuncPtr->kro.max[0], 1.0, 1.0e-10); // Krow(So=Somax) = Krog(So=Somax)
// Gas values
BOOST_CHECK_CLOSE(rfuncPtr->krg.r [0], 0.578171, 1.0e-10); // Krg(Sg=1-Sogcr-Swl) = Krg(Sg=0.690000)
BOOST_CHECK_CLOSE(rfuncPtr->krg.max[0], 1.0 , 1.0e-10); // Krg(Sgmax) = Krg(Sg=0.928996)
// Water values
BOOST_CHECK_CLOSE(rfuncPtr->krw.r [0], 0.261115, 1.0e-10); // Krw(Sw=1-Sowcr-Sgl) = Krw(Sw=0.551004)
BOOST_CHECK_CLOSE(rfuncPtr->krw.max[0], 1.0 , 1.0e-10); // Krw(Swmax) = Krw(Sw=1)
}
BOOST_AUTO_TEST_CASE(RawFunctionValues_Family_II_Tolcrit_Large) {
const auto es = ::Opm::EclipseState {
::Opm::Parser{}.parseString(satfunc_model_setup() + satfunc_family_II() + end())
};
const auto& tm = es.getTableManager();
const auto& ph = es.runspec().phases();
const auto tolcrit = 0.01;
const auto rtepPtr = satfunc::getRawTableEndpoints(tm, ph, tolcrit);
const auto rfuncPtr = satfunc::getRawFunctionValues(tm, ph, *rtepPtr);
// Oil values
BOOST_CHECK_CLOSE(rfuncPtr->kro.rw [0], 0.328347, 1.0e-10); // Krow(So=1-Swcr-Sgl) = Krow(So=0.768996)
BOOST_CHECK_CLOSE(rfuncPtr->kro.rg [0], 0.712749, 1.0e-10); // Krog(So=1-Sgcr-Swl) = Krog(So=0.838996)
BOOST_CHECK_CLOSE(rfuncPtr->kro.max[0], 1.0, 1.0e-10); // Krow(So=Somax) = Krog(So=Somax)
// Gas values
BOOST_CHECK_CLOSE(rfuncPtr->krg.r [0], 0.562914, 1.0e-10); // Krg(Sg=1-Sogcr-Swl) = Krg(Sg=0.680000)
BOOST_CHECK_CLOSE(rfuncPtr->krg.max[0], 1.0 , 1.0e-10); // Krg(Sgmax) = Krg(Sg=0.928996)
// Water values
BOOST_CHECK_CLOSE(rfuncPtr->krw.r [0], 0.261115, 1.0e-10); // Krw(Sw=1-Sowcr-Sgl) = Krw(Sw=0.551004)
BOOST_CHECK_CLOSE(rfuncPtr->krw.max[0], 1.0 , 1.0e-10); // Krw(Swmax) = Krw(Sw=1)
}
// =====================================================================
BOOST_AUTO_TEST_CASE(SatFunc_EndPts_Family_I_TolCrit_Zero) {
const auto es = ::Opm::EclipseState {
::Opm::Parser{}.parseString(satfunc_model_setup() + tolCrit(0.0) + satfunc_family_I() + end())
@ -1576,6 +1736,8 @@ BOOST_AUTO_TEST_CASE(SatFunc_EndPts_Family_II_TolCrit_Large) {
}
}
// =====================================================================
BOOST_AUTO_TEST_CASE(GET_FIPXYZ) {
std::string deck_string = R"(
GRID