Add support for SOF2 for 2p runs.

Make it possible to use SOF2 in combination with either
SGFN og SWFN for gas-oil and water-oil case respectivly.
This commit is contained in:
Tor Harald Sandve 2017-11-23 10:46:55 +01:00
parent cb75a7ad6f
commit f6146db122
3 changed files with 264 additions and 21 deletions

View File

@ -37,6 +37,7 @@
#include <opm/parser/eclipse/EclipseState/Tables/SgfnTable.hpp>
#include <opm/parser/eclipse/EclipseState/Tables/SgofTable.hpp>
#include <opm/parser/eclipse/EclipseState/Tables/SlgofTable.hpp>
#include <opm/parser/eclipse/EclipseState/Tables/Sof2Table.hpp>
#include <opm/parser/eclipse/EclipseState/Tables/Sof3Table.hpp>
#include <opm/parser/eclipse/EclipseState/Tables/SwfnTable.hpp>
#include <opm/parser/eclipse/EclipseState/Tables/SwofTable.hpp>
@ -221,7 +222,6 @@ struct EclEpsScalingPointsInfo
const Opm::EclipseState& eclState,
unsigned satRegionIdx)
{
// TODO: support for the SOF2/SOF3 keyword family
const auto& tables = eclState.getTableManager();
const TableContainer& swofTables = tables.getSwofTables();
const TableContainer& sgofTables = tables.getSgofTables();
@ -229,6 +229,8 @@ struct EclEpsScalingPointsInfo
const TableContainer& swfnTables = tables.getSwfnTables();
const TableContainer& sgfnTables = tables.getSgfnTables();
const TableContainer& sof3Tables = tables.getSof3Tables();
const TableContainer& sof2Tables = tables.getSof2Tables();
bool hasWater = deck.hasKeyword("WATER");
bool hasGas = deck.hasKeyword("GAS");
@ -238,31 +240,50 @@ struct EclEpsScalingPointsInfo
Swl = 0.0;
Swu = 0.0;
Swcr = 0.0;
if (!sgofTables.empty())
extractUnscaledSgof_(sgofTables.getTable<SgofTable>(satRegionIdx));
bool family1 = (!sgofTables.empty() || !slgofTables.empty());
bool family2 = !sgfnTables.empty() && !sof2Tables.empty();
if (family1) {
if (!sgofTables.empty())
extractUnscaledSgof_(sgofTables.getTable<SgofTable>(satRegionIdx));
else {
assert(!slgofTables.empty());
extractUnscaledSlgof_(slgofTables.getTable<SlgofTable>(satRegionIdx));
}
} else if (family2) {
extractUnscaledSgfn_(sgfnTables.getTable<SgfnTable>(satRegionIdx));
extractUnscaledSof2_(sof2Tables.getTable<Sof2Table>(satRegionIdx));
}
else {
assert(!slgofTables.empty());
extractUnscaledSlgof_(slgofTables.getTable<SlgofTable>(satRegionIdx));
throw std::domain_error("No valid saturation keyword family specified");
}
return;
}
else if (!hasGas) {
assert(!swofTables.empty());
Sgl = 0.0;
Sgu = 0.0;
Sgcr = 0.0;
extractUnscaledSwof_(swofTables.getTable<SwofTable>(satRegionIdx));
bool family1 = !swofTables.empty();
bool family2 = !swfnTables.empty() && !sof2Tables.empty();
if (family1) {
extractUnscaledSwof_(swofTables.getTable<SwofTable>(satRegionIdx));
} else if (family2) {
extractUnscaledSwfn_(swfnTables.getTable<SwfnTable>(satRegionIdx));
extractUnscaledSof2_(sof2Tables.getTable<Sof2Table>(satRegionIdx));
}
else {
throw std::domain_error("No valid saturation keyword family specified");
}
return;
}
bool family1 = (!sgofTables.empty() || !slgofTables.empty()) && !swofTables.empty();
bool family2 = !swfnTables.empty() && !sgfnTables.empty() && !sof3Tables.empty();
// so far, only water-oil and oil-gas simulations are supported, i.e.,
// there's no gas-water yet.
if (!hasWater || !hasGas || !hasOil)
throw std::domain_error("The specified phase configuration is not suppored");
bool family1 = (!sgofTables.empty() || !slgofTables.empty()) && !swofTables.empty();
bool family2 = !swfnTables.empty() && !sgfnTables.empty() && !sof3Tables.empty();
if (family1) {
extractUnscaledSwof_(swofTables.getTable<SwofTable>(satRegionIdx));
@ -587,6 +608,32 @@ private:
maxKrow = sof3Table.getKrowColumn().back();
maxKrog = sof3Table.getKrogColumn().back();
}
void extractUnscaledSof2_(const Opm::Sof2Table& sof2Table)
{
// connate oil saturations
Sowl = sof2Table.getSoColumn().front() + Sgl;
Sogl = sof2Table.getSoColumn().front() + Swl;
// maximum oil saturations
Sowu = sof2Table.getSoColumn().back();
// critical oil saturation of oil-water system or
// critical oil saturation of gas-oil system
Sowcr = 0.0;
for (size_t rowIdx = 0 ; rowIdx < sof2Table.numRows(); ++ rowIdx) {
if (sof2Table.getKroColumn()[rowIdx] > 0) {
break;
};
Sowcr = sof2Table.getSoColumn()[rowIdx];
}
Sogcr = Sowcr;
// maximum relative oil permeabilities
maxKrow = sof2Table.getKroColumn().back();
maxKrog = maxKrow;
}
#endif // HAVE_OPM_PARSER
void extractGridPropertyValue_(Scalar& targetValue,

View File

@ -673,6 +673,8 @@ private:
const TableContainer& swfnTables = tableManager.getSwfnTables();
const TableContainer& sgfnTables = tableManager.getSgfnTables();
const TableContainer& sof3Tables = tableManager.getSof3Tables();
const TableContainer& sof2Tables = tableManager.getSof2Tables();
bool hasGas = deck.hasKeyword("GAS");
bool hasOil = deck.hasKeyword("OIL");
@ -683,16 +685,12 @@ private:
if (!hasGas) {
// oil-water case
family1 = !swofTables.empty();
if (!family1)
throw std::runtime_error("only SWOF is supporeted for oil-water two-phase simulations");
//family2 = !swfnTables.empty() && !sof2Tables.empty();
family2 = !swfnTables.empty() && !sof2Tables.empty();
}
else if (!hasWater) {
// oil-gas case
family1 = !sgofTables.empty();
if (!family1)
throw std::runtime_error("only SGOF is supporeted for oil-gas two-phase simulations");
//family2 = !sgfnTables.empty() && !sof2Tables.empty();
family2 = !sgfnTables.empty() && !sof2Tables.empty();
}
else if (!hasOil) {
// water-gas case
@ -762,12 +760,22 @@ private:
case FamilyII:
{
const Sof3Table& sof3Table = tableManager.getSof3Tables().getTable<Sof3Table>( satRegionIdx );
const SgfnTable& sgfnTable = tableManager.getSgfnTables().getTable<SgfnTable>( satRegionIdx );
readGasOilEffectiveParametersFamily2_(effParams,
Swco,
sof3Table,
sgfnTable);
bool hasWater = deck.hasKeyword("WATER");
if (!hasWater) {
// oil and gas case
const Sof2Table& sof2Table = tableManager.getSof2Tables().getTable<Sof2Table>( satRegionIdx );
readGasOilEffectiveParametersFamily2_(effParams,
Swco,
sof2Table,
sgfnTable);
} else {
const Sof3Table& sof3Table = tableManager.getSof3Tables().getTable<Sof3Table>( satRegionIdx );
readGasOilEffectiveParametersFamily2_(effParams,
Swco,
sof3Table,
sgfnTable);
}
break;
}
@ -832,6 +840,24 @@ private:
effParams.finalize();
}
void readGasOilEffectiveParametersFamily2_(GasOilEffectiveTwoPhaseParams& effParams,
Scalar /* Swco */,
const Opm::Sof2Table& sof2Table,
const Opm::SgfnTable& sgfnTable)
{
// convert the saturations of the SGFN keyword from gas to oil saturations
std::vector<double> SoSamples(sgfnTable.numRows());
std::vector<double> SoColumn = sof2Table.getColumn("SO").vectorCopy();
for (size_t sampleIdx = 0; sampleIdx < sgfnTable.numRows(); ++ sampleIdx) {
SoSamples[sampleIdx] = 1 - sgfnTable.get("SG", sampleIdx);
}
effParams.setKrwSamples(SoColumn, sof2Table.getColumn("KRO").vectorCopy());
effParams.setKrnSamples(SoSamples, sgfnTable.getColumn("KRG").vectorCopy());
effParams.setPcnwSamples(SoSamples, sgfnTable.getColumn("PCOG").vectorCopy());
effParams.finalize();
}
template <class Container>
void readOilWaterEffectiveParameters_(Container& dest,
const Opm::Deck& deck,

View File

@ -287,7 +287,122 @@ static const char* hysterDeckString =
"0.85 0.98 0.000 0\n"
"0.88 0.984 0.000 0 /\n";
static const char* fam1DeckStringGasOil =
"RUNSPEC\n"
"\n"
"DIMENS\n"
" 10 10 3 /\n"
"\n"
"TABDIMS\n"
"/\n"
"\n"
"OIL\n"
"GAS\n"
"\n"
"DISGAS\n"
"\n"
"FIELD\n"
"\n"
"GRID\n"
"\n"
"DX\n"
" 300*1000 /\n"
"DY\n"
" 300*1000 /\n"
"DZ\n"
" 100*20 100*30 100*50 /\n"
"\n"
"TOPS\n"
" 100*8325 /\n"
"\n"
"\n"
"PROPS\n"
"\n"
"\n"
"SGOF\n"
"0 0 1 0\n"
"0.001 0 1 0\n"
"0.02 0 0.997 0\n"
"0.05 0.005 0.980 0\n"
"0.12 0.025 0.700 0\n"
"0.2 0.075 0.350 0\n"
"0.25 0.125 0.200 0\n"
"0.3 0.190 0.090 0\n"
"0.4 0.410 0.021 0\n"
"0.45 0.60 0.010 0\n"
"0.5 0.72 0.001 0\n"
"0.6 0.87 0.0001 0\n"
"0.7 0.94 0.000 0\n"
"0.85 0.98 0.000 0\n"
"0.88 0.984 0.000 0 /\n";
static const char* fam2DeckStringGasOil =
"RUNSPEC\n"
"\n"
"DIMENS\n"
" 10 10 3 /\n"
"\n"
"TABDIMS\n"
"/\n"
"\n"
"OIL\n"
"GAS\n"
"\n"
"DISGAS\n"
"\n"
"FIELD\n"
"\n"
"GRID\n"
"\n"
"DX\n"
" 300*1000 /\n"
"DY\n"
" 300*1000 /\n"
"DZ\n"
" 100*20 100*30 100*50 /\n"
"\n"
"TOPS\n"
" 100*8325 /\n"
"\n"
"\n"
"PROPS\n"
"\n"
"PVTW\n"
" 4017.55 1.038 3.22E-6 0.318 0.0 /\n"
"\n"
"\n"
"SGFN\n"
"0 0 0\n"
"0.001 0 0\n"
"0.02 0 0\n"
"0.05 0.005 0\n"
"0.12 0.025 0\n"
"0.2 0.075 0\n"
"0.25 0.125 0\n"
"0.3 0.190 0\n"
"0.4 0.410 0\n"
"0.45 0.60 0\n"
"0.5 0.72 0\n"
"0.6 0.87 0\n"
"0.7 0.94 0\n"
"0.85 0.98 0\n"
"0.88 0.984 0 /\n"
"\n"
"SOF2\n"
"0.12 0.000 \n"
"0.15 0.000 \n"
"0.3 0.000 \n"
"0.4 0.0001 \n"
"0.5 0.001 \n"
"0.55 0.010 \n"
"0.6 0.021 \n"
"0.7 0.090 \n"
"0.8 0.350 \n"
"0.88 0.700 \n"
"0.95 0.980 \n"
"0.98 0.997 \n"
"0.999 1 \n"
"1.0 1 \n /\n";
template <class Scalar>
inline void testAll()
@ -450,6 +565,61 @@ inline void testAll()
}
}
}
// Gas oil
{
const auto fam1Deck = parser.parseString(fam1DeckStringGasOil, parseContext);
const Opm::EclipseState fam1EclState(fam1Deck, parseContext);
MaterialLawManager fam1materialLawManager;
fam1materialLawManager.initFromDeck(fam1Deck, fam1EclState, compressedToCartesianIdx);
const auto fam2Deck = parser.parseString(fam2DeckStringGasOil, parseContext);
const Opm::EclipseState fam2EclState(fam2Deck, parseContext);
Opm::EclMaterialLawManager<MaterialTraits> fam2MaterialLawManager;
fam2MaterialLawManager.initFromDeck(fam2Deck, fam2EclState, compressedToCartesianIdx);
for (unsigned elemIdx = 0; elemIdx < n; ++ elemIdx) {
for (int i = 0; i < 100; ++ i) {
Scalar Sw = 0;
Scalar So = Scalar(i)/100;
Scalar Sg = 1 - Sw - So;
FluidState fs;
fs.setSaturation(waterPhaseIdx, Sw);
fs.setSaturation(oilPhaseIdx, So);
fs.setSaturation(gasPhaseIdx, Sg);
Scalar pcFam1[numPhases] = { 0.0, 0.0 };
Scalar pcFam2[numPhases] = { 0.0, 0.0 };
MaterialLaw::capillaryPressures(pcFam1,
fam1materialLawManager.materialLawParams(elemIdx),
fs);
MaterialLaw::capillaryPressures(pcFam2,
fam2MaterialLawManager.materialLawParams(elemIdx),
fs);
Scalar krFam1[numPhases] = { 0.0, 0.0 };
Scalar krFam2[numPhases] = { 0.0, 0.0 };
MaterialLaw::relativePermeabilities(krFam1,
fam1materialLawManager.materialLawParams(elemIdx),
fs);
MaterialLaw::relativePermeabilities(krFam2,
fam2MaterialLawManager.materialLawParams(elemIdx),
fs);
for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
if (std::abs(pcFam1[phaseIdx] - pcFam2[phaseIdx]) > 1e-5)
OPM_THROW(std::logic_error,
"Discrepancy between capillary pressure of family 1 and family 2 keywords");
if (std::abs(krFam1[phaseIdx] - krFam2[phaseIdx]) > 1e-1)
OPM_THROW(std::logic_error,
"Discrepancy between relative permeabilities of family 1 and family 2 keywords");
}
}
}
}
}
}