Merge pull request #368 from qilicun/wpolymer

add WPOLYMER to schedule section.
This commit is contained in:
Joakim Hove 2014-12-02 12:21:40 +01:00
commit 00dc724aba
9 changed files with 203 additions and 0 deletions

View File

@ -70,6 +70,7 @@ EclipseState/Schedule/Schedule.cpp
EclipseState/Schedule/Well.cpp
EclipseState/Schedule/WellProductionProperties.cpp
EclipseState/Schedule/WellInjectionProperties.cpp
EclipseState/Schedule/WellPolymerProperties.cpp
EclipseState/Schedule/WellSet.cpp
EclipseState/Schedule/Group.cpp
EclipseState/Schedule/Completion.cpp
@ -134,6 +135,7 @@ EclipseState/Schedule/Schedule.hpp
EclipseState/Schedule/Well.hpp
EclipseState/Schedule/WellProductionProperties.hpp
EclipseState/Schedule/WellInjectionProperties.hpp
EclipseState/Schedule/WellPolymerProperties.hpp
EclipseState/Schedule/WellSet.hpp
EclipseState/Schedule/Group.hpp
EclipseState/Schedule/DynamicState.hpp

View File

@ -21,6 +21,7 @@
#include <opm/parser/eclipse/EclipseState/Schedule/TimeMap.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/WellProductionProperties.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/WellInjectionProperties.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/WellPolymerProperties.hpp>
#include <opm/parser/eclipse/Deck/Section.hpp>
#include <boost/algorithm/string.hpp>
@ -83,6 +84,9 @@ namespace Opm {
if (keyword->name() == "WCONINJE")
handleWCONINJE(deck, keyword, parserLog, currentStep);
if (keyword->name() == "WPOLYMER")
handleWPOLYMER(deck, keyword, parserLog, currentStep);
if (keyword->name() == "WCONINJH")
handleWCONINJH(deck, keyword, parserLog, currentStep);
@ -321,6 +325,35 @@ namespace Opm {
}
void Schedule::handleWPOLYMER(DeckConstPtr deck, DeckKeywordConstPtr keyword, ParserLogPtr /*parserLog*/, size_t currentStep) {
for (size_t recordNr = 0; recordNr < keyword->size(); recordNr++) {
DeckRecordConstPtr record = keyword->getRecord(recordNr);
const std::string& wellNamePattern = record->getItem("WELL")->getTrimmedString(0);
std::vector<WellPtr> wells = getWells(wellNamePattern);
for (auto wellIter=wells.begin(); wellIter != wells.end(); ++wellIter) {
WellPtr well = *wellIter;
WellPolymerProperties properties(well->getPolymerPropertiesCopy(currentStep));
properties.m_polymerConcentration = record->getItem("POLYMER_CONCENTRATION")->getSIDouble(0);
properties.m_saltConcentration = record->getItem("SALT_CONCENTRATION")->getSIDouble(0);
auto group_polymer_item = record->getItem("GROUP_POLYMER_CONCENTRATION");
auto group_salt_item = record->getItem("GROUP_SALT_CONCENTRATION");
if (!group_polymer_item->defaultApplied(0)) {
throw std::logic_error("Sorry explicit setting of \`GROUP_POLYMER_CONCENTRATION\` is not supported!");
}
if (!group_salt_item->defaultApplied(0)) {
throw std::logic_error("Sorry explicit setting of \`GROUP_SALT_CONCENTRATION\` is not supported!");
}
well->setPolymerProperties(currentStep, properties);
}
}
}
void Schedule::handleWCONINJH(DeckConstPtr deck, DeckKeywordConstPtr keyword, ParserLogPtr /*parserLog*/, size_t currentStep) {
for (size_t recordNr = 0; recordNr < keyword->size(); recordNr++) {
DeckRecordConstPtr record = keyword->getRecord(recordNr);
@ -617,6 +650,7 @@ namespace Opm {
}
}
bool Schedule::convertEclipseStringToBool(const std::string& eclipseString) {
std::string lowerTrimmed = boost::algorithm::to_lower_copy(eclipseString);
boost::algorithm::trim(lowerTrimmed);

View File

@ -79,6 +79,7 @@ namespace Opm
void handleWGRUPCON(DeckKeywordConstPtr keyword, ParserLogPtr parserLog, size_t currentStep);
void handleCOMPDAT(DeckKeywordConstPtr keyword, ParserLogPtr parserLog, size_t currentStep);
void handleWCONINJE(DeckConstPtr deck, DeckKeywordConstPtr keyword, ParserLogPtr parserLog, size_t currentStep);
void handleWPOLYMER(DeckConstPtr deck, DeckKeywordConstPtr keyword, ParserLogPtr parserLog, size_t currentStep);
void handleWCONINJH(DeckConstPtr deck, DeckKeywordConstPtr keyword, ParserLogPtr parserLog, size_t currentStep);
void handleWELOPEN(DeckKeywordConstPtr keyword, ParserLogPtr parserLog, size_t currentStep);
void handleGCONINJE(DeckConstPtr deck, DeckKeywordConstPtr keyword, ParserLogPtr parserLog, size_t currentStep);

View File

@ -38,6 +38,7 @@ namespace Opm {
m_completions( new DynamicState<CompletionSetConstPtr>( timeMap , CompletionSetConstPtr( new CompletionSet()) )),
m_productionProperties( new DynamicState<WellProductionProperties>(timeMap, WellProductionProperties() )),
m_injectionProperties( new DynamicState<WellInjectionProperties>(timeMap, WellInjectionProperties() )),
m_polymerProperties( new DynamicState<WellPolymerProperties>(timeMap, WellPolymerProperties() )),
m_groupName( new DynamicState<std::string>( timeMap , "" )),
m_headI(headI),
m_headJ(headJ),
@ -60,6 +61,7 @@ namespace Opm {
m_completions( new DynamicState<CompletionSetConstPtr>( timeMap , CompletionSetConstPtr( new CompletionSet()) )),
m_productionProperties( new DynamicState<WellProductionProperties>(timeMap, WellProductionProperties() )),
m_injectionProperties( new DynamicState<WellInjectionProperties>(timeMap, WellInjectionProperties() )),
m_polymerProperties( new DynamicState<WellPolymerProperties>(timeMap, WellPolymerProperties() )),
m_groupName( new DynamicState<std::string>( timeMap , "" )),
m_headI(headI),
m_headJ(headJ),
@ -102,6 +104,19 @@ namespace Opm {
return m_injectionProperties->at(timeStep);
}
void Well::setPolymerProperties(size_t timeStep , const WellPolymerProperties newProperties) {
m_isProducer->add(timeStep , false);
m_polymerProperties->add(timeStep, newProperties);
}
WellPolymerProperties Well::getPolymerPropertiesCopy(size_t timeStep) const {
return m_polymerProperties->get(timeStep);
}
const WellPolymerProperties& Well::getPolymerProperties(size_t timeStep) const {
return m_polymerProperties->at(timeStep);
}
bool Well::hasBeenDefined(size_t timeStep) const {
if (timeStep < m_creationTimeStep)
return false;

View File

@ -27,6 +27,7 @@
#include <opm/parser/eclipse/EclipseState/Schedule/Completion.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/WellProductionProperties.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/WellInjectionProperties.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/WellPolymerProperties.hpp>
#include <boost/optional.hpp>
@ -83,6 +84,10 @@ namespace Opm {
WellInjectionProperties getInjectionPropertiesCopy(size_t timeStep) const;
const WellInjectionProperties& getInjectionProperties(size_t timeStep) const;
void setPolymerProperties(size_t timeStep , const WellPolymerProperties properties);
WellPolymerProperties getPolymerPropertiesCopy(size_t timeStep) const;
const WellPolymerProperties& getPolymerProperties(size_t timeStep) const;
private:
size_t m_creationTimeStep;
std::string m_name;
@ -98,6 +103,7 @@ namespace Opm {
std::shared_ptr<DynamicState<CompletionSetConstPtr> > m_completions;
std::shared_ptr<DynamicState<WellProductionProperties> > m_productionProperties;
std::shared_ptr<DynamicState<WellInjectionProperties> > m_injectionProperties;
std::shared_ptr<DynamicState<WellPolymerProperties> > m_polymerProperties;
std::shared_ptr<DynamicState<std::string> > m_groupName;
// WELSPECS data - assumes this is not dynamic

View File

@ -0,0 +1,12 @@
#include <opm/parser/eclipse/EclipseState/Schedule/WellPolymerProperties.hpp>
#include <string>
#include <vector>
namespace Opm {
WellPolymerProperties::WellPolymerProperties() {
m_polymerConcentration = 0.0;
m_saltConcentration = 0.0;
}
}

View File

@ -0,0 +1,34 @@
/*
Copyright 2014 Statoil ASA.
This file is part of the Open Porous Media project (OPM).
OPM is free software: you can redistribute it and/or modify
it under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OPM is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
GNU General Public License for more details.
You should have received a copy of the GNU General Public License
along with OPM. If not, see <http://www.gnu.org/licenses/>.
*/
#ifndef WELLPOLYMERPROPERTIES_HPP_HEADER_INCLUDED
#define WELLPOLYMERPROPERTIES_HPP_HEADER_INCLUDED
namespace Opm {
struct WellPolymerProperties {
double m_polymerConcentration;
double m_saltConcentration;
WellPolymerProperties();
};
}
#endif

View File

@ -592,3 +592,46 @@ BOOST_AUTO_TEST_CASE(WELLS_SHUT) {
BOOST_CHECK_EQUAL( WellCommon::StatusEnum::SHUT , well3->getStatus(2));
}
BOOST_AUTO_TEST_CASE(WellTestWPOLYMER) {
ParserPtr parser(new Parser());
boost::filesystem::path scheduleFile("testdata/integration_tests/SCHEDULE/SCHEDULE_POLYMER");
DeckPtr deck = parser->parseFile(scheduleFile.string());
ScheduleConstPtr sched(new Schedule(deck));
BOOST_CHECK_EQUAL(4U, sched->numWells());
BOOST_CHECK(sched->hasWell("INJE01"));
BOOST_CHECK(sched->hasWell("PROD01"));
WellConstPtr well1 = sched->getWell("INJE01");
BOOST_CHECK( well1->isInjector(0));
{
const WellPolymerProperties& props_well10 = well1->getPolymerProperties(0);
BOOST_CHECK_CLOSE(1.5*Metric::PolymerDensity, props_well10.m_polymerConcentration, 0.0001);
const WellPolymerProperties& props_well11 = well1->getPolymerProperties(1);
BOOST_CHECK_CLOSE(1.0*Metric::PolymerDensity, props_well11.m_polymerConcentration, 0.0001);
const WellPolymerProperties& props_well12 = well1->getPolymerProperties(2);
BOOST_CHECK_CLOSE(0.1*Metric::PolymerDensity, props_well12.m_polymerConcentration, 0.0001);
}
WellConstPtr well2 = sched->getWell("INJE02");
BOOST_CHECK( well2->isInjector(0));
{
const WellPolymerProperties& props_well20 = well2->getPolymerProperties(0);
BOOST_CHECK_CLOSE(2.0*Metric::PolymerDensity, props_well20.m_polymerConcentration, 0.0001);
const WellPolymerProperties& props_well21 = well2->getPolymerProperties(1);
BOOST_CHECK_CLOSE(1.5*Metric::PolymerDensity, props_well21.m_polymerConcentration, 0.0001);
const WellPolymerProperties& props_well22 = well2->getPolymerProperties(2);
BOOST_CHECK_CLOSE(0.2*Metric::PolymerDensity, props_well22.m_polymerConcentration, 0.0001);
}
WellConstPtr well3 = sched->getWell("INJE03");
BOOST_CHECK( well3->isInjector(0));
{
const WellPolymerProperties& props_well30 = well3->getPolymerProperties(0);
BOOST_CHECK_CLOSE(2.5*Metric::PolymerDensity, props_well30.m_polymerConcentration, 0.0001);
const WellPolymerProperties& props_well31 = well3->getPolymerProperties(1);
BOOST_CHECK_CLOSE(2.0*Metric::PolymerDensity, props_well31.m_polymerConcentration, 0.0001);
const WellPolymerProperties& props_well32 = well3->getPolymerProperties(2);
BOOST_CHECK_CLOSE(0.3*Metric::PolymerDensity, props_well32.m_polymerConcentration, 0.0001);
}
}

View File

@ -0,0 +1,56 @@
START
10 MAI 2007 /
SCHEDULE
WELSPECS
'INJE01' 'I' 1 1 1* 'WATER' /
'INJE02' 'I' 5 1 1* 'WATER' /
'INJE03' 'I' 10 1 1* 'WATER' /
'PROD01' 'P' 20 1 1* 'OIL' 7* /
/
COMPDAT
'INJE01' 1 1 2 2 'OPEN' 1* 200. 0.5 /
'INJE02' 5 1 2 2 'OPEN' 1* 200. 0.5 /
'INJE03' 10 1 2 2 'OPEN' 1* 200. 0.5 /
'PROD01' 20 1 1 1 'OPEN' 1* 200. 0.5 /
/
WCONINJE
'INJE01' 'WATER' 'OPEN' 'RATE' 800.00 1* 1000 /
'INJE02' 'WATER' 'OPEN' 'RATE' 800.00 1* 1000 /
'INJE03' 'WATER' 'OPEN' 'RATE' 800.00 1* 1000 /
/
WCONPROD
'PROD01' 'OPEN' 'BHP' 5* 200 /
/
WPOLYMER
'INJE01' 1.5 0.0 /
'INJE02' 2.0 0.0 /
'INJE03' 2.5 0.0 /
/
TSTEP
10
/
WPOLYMER
'INJE01' 1.0 0.0 /
'INJE02' 1.5 0.0 /
'INJE03' 2.0 0.0 /
/
TSTEP
10
/
WPOLYMER
'INJE01' 0.1 0.0 /
'INJE02' 0.2 0.0 /
'INJE03' 0.3 0.0 /
/
TSTEP
10
/