ERT-700: Included Fault treatment in EclipseState

This commit is contained in:
Joakim Hove
2014-07-05 13:47:22 +02:00
parent 5411f9ec74
commit c4cd20f9f2
3 changed files with 101 additions and 1 deletions

View File

@@ -40,6 +40,7 @@ namespace Opm {
initTitle(deck);
initProperties(deck);
initTransMult(deck);
initFaults(deck);
}
@@ -52,6 +53,10 @@ namespace Opm {
return schedule;
}
std::shared_ptr<const FaultCollection> EclipseState::getFaults() const {
return m_faults;
}
std::shared_ptr<const TransMult> EclipseState::getTransMult() const {
return m_transMult;
}
@@ -69,6 +74,62 @@ namespace Opm {
m_transMult = std::make_shared<TransMult>( grid->getNX() , grid->getNY() , grid->getNZ());
}
void EclipseState::initFaults(DeckConstPtr deck) {
EclipseGridConstPtr grid = getEclipseGrid();
m_faults = std::make_shared<FaultCollection>( grid->getNX() , grid->getNY() , grid->getNZ());
std::shared_ptr<Opm::GRIDSection> gridSection(new Opm::GRIDSection(deck) );
std::shared_ptr<Opm::EDITSection> editSection(new Opm::EDITSection(deck) );
for (size_t index=0; index < gridSection->count("FAULTS"); index++) {
DeckKeywordConstPtr faultsKeyword = gridSection->getKeyword("FAULTS" , index);
for (auto iter = faultsKeyword->begin(); iter != faultsKeyword->end(); ++iter) {
DeckRecordConstPtr faultRecord = *iter;
const std::string& faultName = faultRecord->getItem(0)->getString(0);
int I1 = faultRecord->getItem(1)->getInt(0) - 1;
int I2 = faultRecord->getItem(2)->getInt(0) - 1;
int J1 = faultRecord->getItem(3)->getInt(0) - 1;
int J2 = faultRecord->getItem(4)->getInt(0) - 1;
int K1 = faultRecord->getItem(5)->getInt(0) - 1;
int K2 = faultRecord->getItem(6)->getInt(0) - 1;
FaceDir::DirEnum faceDir = FaceDir::FromString( faultRecord->getItem(7)->getString(0) );
std::shared_ptr<const FaultFace> face = std::make_shared<const FaultFace>(grid->getNX() , grid->getNY() , grid->getNZ(),
static_cast<size_t>(I1) , static_cast<size_t>(I2) ,
static_cast<size_t>(J1) , static_cast<size_t>(J2) ,
static_cast<size_t>(K1) , static_cast<size_t>(K2) ,
faceDir);
if (!m_faults->hasFault(faultName)) {
std::shared_ptr<Fault> fault = std::make_shared<Fault>( faultName );
m_faults->addFault( fault );
}
{
std::shared_ptr<Fault> fault = m_faults->getFault( faultName );
fault->addFace( face );
}
}
}
setMULTFLT( gridSection );
setMULTFLT( editSection );
m_transMult->applyMULTFLT( m_faults );
}
void EclipseState::setMULTFLT(std::shared_ptr<const Section> section) const {
for (size_t index=0; index < section->count("MULTFLT"); index++) {
DeckKeywordConstPtr faultsKeyword = section->getKeyword("MULTFLT" , index);
for (auto iter = faultsKeyword->begin(); iter != faultsKeyword->end(); ++iter) {
DeckRecordConstPtr faultRecord = *iter;
const std::string& faultName = faultRecord->getItem(0)->getString(0);
double multFlt = faultRecord->getItem(1)->getRawDouble(0);
m_faults->setTransMult( faultName , multFlt );
}
}
}
void EclipseState::initEclipseGrid(DeckConstPtr deck) {
std::shared_ptr<Opm::GRIDSection> gridSection(new Opm::GRIDSection(deck) );

View File

@@ -28,6 +28,7 @@
#include <opm/parser/eclipse/EclipseState/Grid/GridProperties.hpp>
#include <opm/parser/eclipse/EclipseState/Grid/TransMult.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/ScheduleEnums.hpp>
#include <opm/parser/eclipse/EclipseState/Grid/FaultCollection.hpp>
#include <set>
#include <memory>
@@ -50,6 +51,7 @@ namespace Opm {
void loadGridPropertyFromDeckKeyword(std::shared_ptr<const Box> inputBox , DeckKeywordConstPtr deckKeyword);
std::shared_ptr<const FaultCollection> getFaults() const;
std::shared_ptr<const TransMult> getTransMult() const;
private:
void initSchedule(DeckConstPtr deck);
@@ -58,6 +60,8 @@ namespace Opm {
void initTitle(DeckConstPtr deck);
void initProperties(DeckConstPtr deck);
void initTransMult(DeckConstPtr deck);
void initFaults(DeckConstPtr deck);
void setMULTFLT(std::shared_ptr<const Section> section) const;
double getSIScaling(const std::string &dimensionString) const;
@@ -82,6 +86,7 @@ namespace Opm {
std::shared_ptr<GridProperties<int> > m_intGridProperties;
std::shared_ptr<GridProperties<double> > m_doubleGridProperties;
std::shared_ptr<TransMult> m_transMult;
std::shared_ptr<FaultCollection> m_faults;
};
typedef std::shared_ptr<EclipseState> EclipseStatePtr;

View File

@@ -51,7 +51,18 @@ static DeckPtr createDeck() {
"1000*0.25 /\n"
"TOPS\n"
"1000*0.25 /\n"
"FAULTS \n"
" 'F1' 1 1 1 4 1 4 'X' / \n"
" 'F2' 5 5 1 4 1 4 'X-' / \n"
"/\n"
"MULTFLT \n"
" 'F1' 0.50 / \n"
" 'F2' 0.50 / \n"
"/\n"
"EDIT\n"
"MULTFLT /\n"
" 'F2' 0.25 / \n"
"/\n"
"OIL\n"
"\n"
"GAS\n"
@@ -145,6 +156,29 @@ BOOST_AUTO_TEST_CASE(GetTransMult) {
std::shared_ptr<const TransMult> transMult = state.getTransMult();
BOOST_CHECK_EQUAL( 1.0 , transMult->getMultiplier(0,0,0,FaceDir::XPlus));
BOOST_CHECK_EQUAL( 1.0 , transMult->getMultiplier(1,0,0,FaceDir::XPlus));
BOOST_CHECK_THROW(transMult->getMultiplier(1000 , FaceDir::XPlus) , std::invalid_argument);
}
BOOST_AUTO_TEST_CASE(GetFaults) {
DeckPtr deck = createDeck();
EclipseState state(deck);
std::shared_ptr<const FaultCollection> faults = state.getFaults();
BOOST_CHECK( faults->hasFault("F1") );
BOOST_CHECK( faults->hasFault("F2") );
std::shared_ptr<Fault> F1 = faults->getFault("F1");
std::shared_ptr<Fault> F2 = faults->getFault("F2");
BOOST_CHECK_EQUAL( 0.50 , F1->getTransMult());
BOOST_CHECK_EQUAL( 0.25 , F2->getTransMult());
std::shared_ptr<const TransMult> transMult = state.getTransMult();
BOOST_CHECK_EQUAL( transMult->getMultiplier(0 , 0 , 0 , FaceDir::XPlus) , 0.50 );
BOOST_CHECK_EQUAL( transMult->getMultiplier(4 , 3 , 0 , FaceDir::XMinus) , 0.25 );
BOOST_CHECK_EQUAL( transMult->getMultiplier(4 , 3 , 0 , FaceDir::ZPlus) , 1.00 );
}