Merge pull request #909 from andlaus/implement_multregp

Implement support for MULTREGP
This commit is contained in:
Joakim Hove
2016-09-05 18:11:18 +02:00
committed by GitHub
2 changed files with 124 additions and 2 deletions

View File

@@ -33,7 +33,6 @@
#include <opm/parser/eclipse/EclipseState/Tables/TableManager.hpp>
#include <opm/parser/eclipse/Utility/String.hpp>
namespace Opm {
namespace {
@@ -50,9 +49,26 @@ namespace Opm {
}
}
/// initPORV uses doubleGridProperties: PORV, PORO, NTG, MULTPV
// this funcion applies a single pore volume multiplier (i.e., it deals with a single
// record of the MULTREGP keyword)
void applyPorosityRegionMultiplier_(std::vector<double>& porev,
const std::vector<int>& regionId,
int multRegionId,
double multValue)
{
for (unsigned i = 0; i < porev.size(); ++i) {
if (regionId[i] == multRegionId)
porev[i] *= multValue;
}
}
/// this function initializes the pore volume of all cells. it uses the raw keyword
/// 'MULTREGP', the integer grid properties 'FLUXNUM', 'MULTNUM' and 'OPERNUM' as
/// well as the double grid properties 'PORV', 'PORO', 'NTG' and 'MULTPV'
void initPORV( std::vector<double>& values,
const Deck* deck,
const EclipseGrid* eclipseGrid,
const GridProperties<int>* intGridProperties,
const GridProperties<double>* doubleGridProperties)
{
const auto iter = std::find_if( values.begin(), values.end(), []( double value ) { return std::isnan(value); });
@@ -81,6 +97,54 @@ namespace Opm {
for (size_t globalIndex = 0; globalIndex < multpvData.size(); globalIndex++)
values[globalIndex] *= multpvData[globalIndex];
}
// deal with the region multiplier for porosity
if (deck->hasKeyword("MULTREGP")) {
const DeckKeyword& multregpKeyword = deck->getKeyword("MULTREGP");
for (unsigned recordIdx = 0; recordIdx < multregpKeyword.size(); ++recordIdx) {
const DeckRecord& multregpRecord = multregpKeyword.getRecord(recordIdx);
int regionId = multregpRecord.getItem("REGION").template get<int>(0);
std::string regionType = multregpRecord.getItem("REGION_TYPE").template get<std::string>(0);
double multValue = multregpRecord.getItem("MULTIPLIER").template get<double>(0);
uppercase(regionType, regionType);
// deal with the "convenience" feature of ECL: if the region index is
// zero or negative, the record is ignored.
if (regionId <= 0)
continue;
// implement a (documented) ECL bug: the region index must be unique
// (i.e., it is impossible to specify multipliers for different
// region types with the same index). Also, only the last occurence
// of each region counts.
unsigned record2Idx = recordIdx + 1;
for (; record2Idx < multregpKeyword.size(); ++record2Idx) {
const DeckRecord& multregpRecord2 = multregpKeyword.getRecord(record2Idx);
int region2Idx = multregpRecord2.getItem("REGION").template get<int>(0);
if (region2Idx == regionId)
break;
}
if (record2Idx < multregpKeyword.size())
// the region was specified twice
continue;
if (regionType == "M") {
const GridProperty<int>& multnumProp = intGridProperties->getKeyword("MULTNUM");
applyPorosityRegionMultiplier_(values, multnumProp.getData(), regionId, multValue);
}
else if (regionType == "F") {
const GridProperty<int>& fluxnumProp = intGridProperties->getKeyword("FLUXNUM");
applyPorosityRegionMultiplier_(values, fluxnumProp.getData(), regionId, multValue);
}
else if (regionType == "O") {
const GridProperty<int>& opernumProp = intGridProperties->getKeyword("OPERNUM");
applyPorosityRegionMultiplier_(values, opernumProp.getData(), regionId, multValue);
}
else
throw std::logic_error("Unknown or illegal region type for MULTREGP keyword: '"+regionType+"'");
}
}
}
@@ -396,7 +460,9 @@ namespace Opm {
{
auto initPORVProcessor = std::bind(&initPORV,
std::placeholders::_1,
&deck,
&eclipseGrid,
&m_intGridProperties,
&m_doubleGridProperties);
m_doubleGridProperties.postAddKeyword( "PORV",

View File

@@ -228,6 +228,39 @@ static Opm::DeckPtr createDeckWithNTG() {
return parser->parseString(deckData, Opm::ParseContext()) ;
}
static Opm::DeckPtr createDeckWithMULTREGP() {
const char *deckData =
"RUNSPEC\n"
"\n"
"DIMENS\n"
" 10 10 10 /\n"
"GRID\n"
"DX\n"
"1000*0.25 /\n"
"DYV\n"
"10*0.25 /\n"
"DZ\n"
"1000*0.25 /\n"
"TOPS\n"
"100*0.25 /\n"
"PORV\n"
"1000*77 /\n"
"MULTNUM\n"
"100*1 700*2 200*3 / \n"
"FLUXNUM\n"
"200*1 700*2 100*3 / \n"
"MULTREGP\n"
"1 10.0 / \n"
"1 2.0 F/ \n"
"0 123.0 F/ \n" // ignored
"2 20.0 M / \n"
"/\n"
"\n";
Opm::ParserPtr parser(new Opm::Parser());
return parser->parseString(deckData, Opm::ParseContext()) ;
}
BOOST_AUTO_TEST_CASE(PORV_cartesianDeck) {
/* Check that an exception is raised if we try to create a PORV field without PORO. */
Opm::DeckPtr deck = createCARTDeck();
@@ -324,6 +357,29 @@ BOOST_AUTO_TEST_CASE(PORV_multpvAndNtg) {
BOOST_CHECK_CLOSE( cell_volume * poro*multpv*NTG , porv.iget(9,9,9) , 0.001);
}
BOOST_AUTO_TEST_CASE(PORV_multregp) {
Opm::DeckPtr deck = createDeckWithMULTREGP();
const auto props = getProps(deck);
const auto& porv = props.getDoubleGridProperty("PORV");
const auto& porvData = porv.getData();
double basePorv = 77.0;
// the cumlative effect of the base pore volume is multiplied 2 for the cells in
// region 1 of FLUXNUM and by 20 for region 2 of MULTNUM. This means that the first
// 100 cels are multiplied by 2, then follow 100 cells multiplied by 2*20, then
// 600 cells multiplied by 20 while the last 200 cells are not modified at all.
for (int i=0; i < 100; ++i)
BOOST_CHECK_CLOSE(porvData[i], basePorv*2, 1e-8);
for (int i=100; i < 200; ++i)
BOOST_CHECK_CLOSE(porvData[i], basePorv*2*20, 1e-8);
for (int i=200; i < 800; ++i)
BOOST_CHECK_CLOSE(porvData[i], basePorv*20, 1e-8);
for (int i=800; i < 1000; ++i)
BOOST_CHECK_CLOSE(porvData[i], basePorv, 1e-8);
}
static Opm::DeckPtr createDeckNakedGRID() {