From 7ce91cf5044ee5ba1c51161951249f48c4a60255 Mon Sep 17 00:00:00 2001 From: Joakim Hove Date: Tue, 28 Jan 2014 14:31:51 +0100 Subject: [PATCH] Added well->Injectortype --- opm/parser/eclipse/EclipseState/Schedule/Well.cpp | 9 +++++++++ opm/parser/eclipse/EclipseState/Schedule/Well.hpp | 3 +++ .../eclipse/EclipseState/Schedule/tests/WellTests.cpp | 11 +++++++++++ 3 files changed, 23 insertions(+) diff --git a/opm/parser/eclipse/EclipseState/Schedule/Well.cpp b/opm/parser/eclipse/EclipseState/Schedule/Well.cpp index b30c4bc9c..070ecaed6 100644 --- a/opm/parser/eclipse/EclipseState/Schedule/Well.cpp +++ b/opm/parser/eclipse/EclipseState/Schedule/Well.cpp @@ -35,6 +35,7 @@ namespace Opm { m_reservoirInjectionRate(new DynamicState(timeMap, 0.0)), m_BHPLimit(new DynamicState(timeMap , 0.0)), m_THPLimit(new DynamicState(timeMap , 0.0)), + m_injectorType(new DynamicState(timeMap, InjectorType::WATER)), m_inPredictionMode(new DynamicState(timeMap, true)), m_isProducer(new DynamicState(timeMap, true)) , m_completions( new DynamicState( timeMap , CompletionSetConstPtr( new CompletionSet()) )), @@ -76,6 +77,14 @@ namespace Opm { m_THPLimit->add(timeStep, THPLimit); } + InjectorType::InjectorEnum Well::getInjectorType(size_t timeStep) const { + return m_injectorType->get(timeStep); + } + + void Well::setInjectorType(size_t timeStep, InjectorType::InjectorEnum injectorType) { + m_injectorType->add(timeStep , injectorType); + } + double Well::getOilRate(size_t timeStep) const { return m_oilRate->get(timeStep); diff --git a/opm/parser/eclipse/EclipseState/Schedule/Well.hpp b/opm/parser/eclipse/EclipseState/Schedule/Well.hpp index 70a010127..a57861cab 100644 --- a/opm/parser/eclipse/EclipseState/Schedule/Well.hpp +++ b/opm/parser/eclipse/EclipseState/Schedule/Well.hpp @@ -54,6 +54,8 @@ namespace Opm { void setBHPLimit(size_t timeStep, double BHPLimit); double getTHPLimit(size_t timeStep) const; void setTHPLimit(size_t timeStep, double THPLimit); + InjectorType::InjectorEnum getInjectorType(size_t timeStep) const; + void setInjectorType(size_t timeStep, InjectorType::InjectorEnum injectorType); int getHeadI() const; @@ -81,6 +83,7 @@ namespace Opm { std::shared_ptr > m_reservoirInjectionRate; std::shared_ptr > m_BHPLimit; std::shared_ptr > m_THPLimit; + std::shared_ptr > m_injectorType; std::shared_ptr > m_inPredictionMode; std::shared_ptr > m_isProducer; diff --git a/opm/parser/eclipse/EclipseState/Schedule/tests/WellTests.cpp b/opm/parser/eclipse/EclipseState/Schedule/tests/WellTests.cpp index bddbdcc02..f06a9c5ce 100644 --- a/opm/parser/eclipse/EclipseState/Schedule/tests/WellTests.cpp +++ b/opm/parser/eclipse/EclipseState/Schedule/tests/WellTests.cpp @@ -29,6 +29,7 @@ #include #include +#include #include #include @@ -253,3 +254,13 @@ BOOST_AUTO_TEST_CASE(XHPLimitDefault) { BOOST_CHECK_EQUAL( 200 , well.getTHPLimit( 5 )); } + + +BOOST_AUTO_TEST_CASE(InjectorType) { + Opm::TimeMapPtr timeMap = createXDaysTimeMap(10); + Opm::Well well("WELL1", 1, 2, 2334.32, timeMap, 0); + + well.setInjectorType( 1 , Opm::InjectorType::WATER ); + BOOST_CHECK_EQUAL( Opm::InjectorType::WATER , well.getInjectorType( 5 )); +} +