Made FaultCollection reference in State thereby fixing compiler warn

This commit is contained in:
Pål Grønås Drange
2016-05-02 17:50:45 +02:00
parent 63e452167c
commit 5a7dcc5427
8 changed files with 74 additions and 78 deletions

View File

@@ -145,7 +145,7 @@ namespace Opm {
return m_simulationConfig;
}
std::shared_ptr<const FaultCollection> EclipseState::getFaults() const {
const FaultCollection& EclipseState::getFaults() const {
return m_faults;
}
@@ -208,11 +208,11 @@ namespace Opm {
m_transMult->applyMULT(p.getDoubleGridProperty("MULTZ-"), FaceDir::ZMinus);
}
void EclipseState::initFaults(DeckConstPtr deck) {
void EclipseState::initFaults(const Deck& deck) {
EclipseGridConstPtr grid = getInputGrid();
std::shared_ptr<GRIDSection> gridSection = std::make_shared<GRIDSection>( deck );
m_faults = std::make_shared<FaultCollection>(gridSection , grid);
m_faults = FaultCollection(gridSection, grid);
setMULTFLT(gridSection);
if (Section::hasEDIT(deck)) {
@@ -279,17 +279,20 @@ namespace Opm {
}
}
void EclipseState::applyModifierDeck(std::shared_ptr<const Deck> deckptr) {
applyModifierDeck(*deckptr);
}
void EclipseState::applyModifierDeck( std::shared_ptr<const Deck> deck) {
void EclipseState::applyModifierDeck(const Deck& deck) {
using namespace ParserKeywords;
for (const auto& keyword : *deck) {
for (const auto& keyword : deck) {
if (keyword.isKeyword<MULTFLT>()) {
for (const auto& record : keyword) {
const std::string& faultName = record.getItem<MULTFLT::fault>().get< std::string >(0);
auto fault = m_faults->getFault( faultName );
auto& fault = m_faults.getFault( faultName );
double tmpMultFlt = record.getItem<MULTFLT::factor>().get< double >(0);
double oldMultFlt = fault->getTransMult( );
double oldMultFlt = fault.getTransMult( );
double newMultFlt = oldMultFlt * tmpMultFlt;
/*
@@ -312,9 +315,9 @@ namespace Opm {
MULTFLT value in the fault instance.
*/
fault->setTransMult( tmpMultFlt );
fault.setTransMult( tmpMultFlt );
m_transMult->applyMULTFLT( fault );
fault->setTransMult( newMultFlt );
fault.setTransMult( newMultFlt );
}
}
}

View File

@@ -26,6 +26,10 @@
#include <vector>
#include <opm/parser/eclipse/EclipseState/Eclipse3DProperties.hpp>
#include <opm/parser/eclipse/EclipseState/Grid/FaultCollection.hpp>
#include <opm/parser/eclipse/EclipseState/Grid/FaultFace.hpp>
#include <opm/parser/eclipse/EclipseState/Grid/Fault.hpp>
#include <opm/parser/eclipse/EclipseState/Grid/NNC.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/ScheduleEnums.hpp>
#include <opm/parser/eclipse/EclipseState/Tables/TableManager.hpp>
#include <opm/parser/eclipse/Parser/MessageContainer.hpp>
@@ -43,8 +47,6 @@ namespace Opm {
class DeckRecord;
class EclipseGrid;
class Eclipse3DProperties;
class Fault;
class FaultCollection;
class InitConfig;
class IOConfig;
class NNC;
@@ -83,9 +85,9 @@ namespace Opm {
MessageContainer& getMessageContainer();
std::string getTitle() const;
std::shared_ptr<const FaultCollection> getFaults() const;
const FaultCollection& getFaults() const;
std::shared_ptr<const TransMult> getTransMult() const;
std::shared_ptr<const NNC> getNNC() const;
const NNC& getNNC() const;
bool hasNNC() const;
const Eclipse3DProperties& get3DProperties() const;
@@ -97,18 +99,20 @@ namespace Opm {
// the unit system used by the deck. note that it is rarely needed to convert
// units because internally to opm-parser everything is represented by SI
// units...
const UnitSystem& getDeckUnitSystem() const;
void applyModifierDeck( const Deck& deck);
const UnitSystem& getDeckUnitSystem() const;
/// [deprecated]
void applyModifierDeck(std::shared_ptr<const Deck>);
void applyModifierDeck(const Deck& deck);
private:
void initIOConfig(const Deck& deck);
void initIOConfigPostSchedule(const Deck& deck);
void initTransMult();
void initFaults(std::shared_ptr< const Deck > deck);
void initFaults(const Deck& deck);
void setMULTFLT(std::shared_ptr<const Opm::Section> section) const;
void initMULTREGT(std::shared_ptr< const Deck > deck);
void setMULTFLT(std::shared_ptr<const Opm::Section> section);
void initMULTREGT(const Deck& deck);
void complainAboutAmbiguousKeyword(const Deck& deck,
const std::string& keywordName);
@@ -120,8 +124,8 @@ namespace Opm {
std::string m_title;
std::shared_ptr<TransMult> m_transMult;
std::shared_ptr<FaultCollection> m_faults;
std::shared_ptr<NNC> m_nnc;
FaultCollection m_faults;
NNC m_nnc;
const UnitSystem& m_deckUnitSystem;

View File

@@ -34,9 +34,7 @@
namespace Opm {
FaultCollection::FaultCollection()
{
}
{}
FaultCollection::FaultCollection( std::shared_ptr<const GRIDSection> deck, std::shared_ptr<const EclipseGrid> grid) {
@@ -60,58 +58,48 @@ namespace Opm {
static_cast<size_t>(K1) , static_cast<size_t>(K2) ,
faceDir);
if (!hasFault(faultName)) {
std::shared_ptr<Fault> fault = std::make_shared<Fault>( faultName );
addFault( fault );
addFault( faultName );
}
{
std::shared_ptr<Fault> fault = getFault( faultName );
fault->addFace( face );
Fault& fault = getFault( faultName );
fault.addFace( face );
}
}
}
}
size_t FaultCollection::size() const {
return m_faults.size();
}
bool FaultCollection::hasFault(const std::string& faultName) const {
return m_faults.hasKey( faultName );
}
std::shared_ptr<const Fault> FaultCollection::getFault(const std::string& faultName) const {
const Fault& FaultCollection::getFault(const std::string& faultName) const {
return m_faults.get( faultName );
}
std::shared_ptr<Fault> FaultCollection::getFault(const std::string& faultName) {
Fault& FaultCollection::getFault(const std::string& faultName) {
return m_faults.get( faultName );
}
std::shared_ptr<Fault> FaultCollection::getFault(size_t faultIndex) {
Fault& FaultCollection::getFault(size_t faultIndex) {
return m_faults.get( faultIndex );
}
std::shared_ptr<const Fault> FaultCollection::getFault(size_t faultIndex) const {
const Fault& FaultCollection::getFault(size_t faultIndex) const {
return m_faults.get( faultIndex );
}
void FaultCollection::addFault(std::shared_ptr<Fault> fault) {
m_faults.insert(fault->getName() , fault);
void FaultCollection::addFault(const std::string& faultName) {
Fault fault(faultName);
m_faults.insert(fault.getName() , fault);
}
void FaultCollection::setTransMult(const std::string& faultName , double transMult) {
std::shared_ptr<Fault> fault = getFault( faultName );
fault->setTransMult( transMult );
Fault& fault = getFault( faultName );
fault.setTransMult( transMult );
}
}

View File

@@ -24,13 +24,14 @@
#include <memory>
#include <opm/parser/eclipse/EclipseState/Util/OrderedMap.hpp>
#include <opm/parser/eclipse/EclipseState/Grid/Fault.hpp>
namespace Opm {
class Deck;
class GRIDSection;
class EclipseGrid;
class Fault;
class GRIDSection;
@@ -41,15 +42,17 @@ public:
size_t size() const;
bool hasFault(const std::string& faultName) const;
std::shared_ptr<Fault> getFault(const std::string& faultName);
std::shared_ptr<const Fault> getFault(const std::string& faultName) const;
std::shared_ptr<Fault> getFault(size_t faultIndex);
std::shared_ptr<const Fault> getFault(size_t faultIndex) const;
void addFault(std::shared_ptr<Fault> fault);
Fault& getFault(const std::string& faultName);
const Fault& getFault(const std::string& faultName) const;
Fault& getFault(size_t faultIndex);
const Fault& getFault(size_t faultIndex) const;
/// we construct the fault based on faultName. To get the fault: getFault
void addFault(const std::string& faultName);
void setTransMult(const std::string& faultName , double transMult);
private:
OrderedMap<std::shared_ptr<Fault > > m_faults;
OrderedMap<Fault> m_faults;
};
}

View File

@@ -109,10 +109,10 @@ namespace Opm {
}
void TransMult::applyMULTFLT( std::shared_ptr<const Fault> fault) {
double transMult = fault->getTransMult();
void TransMult::applyMULTFLT(const Fault& fault) {
double transMult = fault.getTransMult();
for (auto face_iter = fault->begin(); face_iter != fault->end(); ++face_iter) {
for (auto face_iter = fault.begin(); face_iter != fault.end(); ++face_iter) {
std::shared_ptr<const FaultFace> face = *face_iter;
FaceDir::DirEnum faceDir = face->getDir();
auto& multProperty = getDirectionProperty(faceDir);
@@ -125,10 +125,10 @@ namespace Opm {
}
void TransMult::applyMULTFLT( std::shared_ptr<const FaultCollection> faults) {
for (size_t faultIndex = 0; faultIndex < faults->size(); faultIndex++) {
std::shared_ptr<const Fault> fault = faults->getFault( faultIndex );
applyMULTFLT( fault );
void TransMult::applyMULTFLT(const FaultCollection& faults) {
for (size_t faultIndex = 0; faultIndex < faults.size(); faultIndex++) {
auto& fault = faults.getFault(faultIndex);
applyMULTFLT(fault);
}
}

View File

@@ -51,8 +51,8 @@ namespace Opm {
double getMultiplier(size_t i , size_t j , size_t k, FaceDir::DirEnum faceDir) const;
double getRegionMultiplier( size_t globalCellIndex1, size_t globalCellIndex2, FaceDir::DirEnum faceDir) const;
void applyMULT(const GridProperty<double>& srcMultProp, FaceDir::DirEnum faceDir);
void applyMULTFLT( std::shared_ptr<const FaultCollection> faults);
void applyMULTFLT( std::shared_ptr<const Fault> fault);
void applyMULTFLT(const FaultCollection& faults);
void applyMULTFLT(const Fault& fault);
void setMultregtScanner(std::shared_ptr<const MULTREGTScanner> multregtScanner);
private:

View File

@@ -117,20 +117,18 @@ BOOST_AUTO_TEST_CASE(CreateFaultCollection) {
BOOST_AUTO_TEST_CASE(AddFaultsToCollection) {
Opm::FaultCollection faults;
std::shared_ptr<Opm::Fault> fault = std::make_shared<Opm::Fault>("FAULT");
faults.addFault(fault);
faults.addFault("FAULT");
BOOST_CHECK_EQUAL( faults.size() , 1 );
BOOST_CHECK(faults.hasFault("FAULT"));
std::shared_ptr<Opm::Fault> fault2 = faults.getFault("FAULT");
std::shared_ptr<Opm::Fault> fault0 = faults.getFault(0);
BOOST_CHECK_EQUAL( fault , fault2 );
BOOST_CHECK_EQUAL( fault , fault0 );
const auto& fault1 = faults.getFault("FAULT");
const auto& fault2 = faults.getFault(0);
BOOST_CHECK_EQUAL(fault1.getName(), fault2.getName());
std::shared_ptr<Opm::Fault> faultx = std::make_shared<Opm::Fault>("FAULTX");
faults.addFault(faultx);
faults.addFault("FAULTX");
const auto& faultx = faults.getFault("FAULTX");
BOOST_CHECK_EQUAL( faults.size() , 2 );
BOOST_CHECK(faults.hasFault("FAULTX"));
BOOST_CHECK_EQUAL( faultx , faults.getFault(1));
BOOST_CHECK_EQUAL( faultx.getName() , faults.getFault(1).getName());
}

View File

@@ -284,16 +284,16 @@ BOOST_AUTO_TEST_CASE(GetTransMult) {
BOOST_AUTO_TEST_CASE(GetFaults) {
DeckPtr deck = createDeck();
EclipseState state( deck, ParseContext() );
std::shared_ptr<const FaultCollection> faults = state.getFaults();
const auto& faults = state.getFaults();
BOOST_CHECK( faults->hasFault( "F1" ) );
BOOST_CHECK( faults->hasFault( "F2" ) );
BOOST_CHECK( faults.hasFault( "F1" ) );
BOOST_CHECK( faults.hasFault( "F2" ) );
std::shared_ptr<const Fault> F1 = faults->getFault( "F1" );
std::shared_ptr<const Fault> F2 = faults->getFault( "F2" );
const auto& F1 = faults.getFault( "F1" );
const auto& F2 = faults.getFault( "F2" );
BOOST_CHECK_EQUAL( 0.50, F1->getTransMult() );
BOOST_CHECK_EQUAL( 0.25, F2->getTransMult() );
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 );