changed: simplify OilVaporizationProperties

there is no vappars per pvt region. do not store it as such either.
This commit is contained in:
Arne Morten Kvarving
2020-02-28 10:54:46 +01:00
parent 8a4240b2fe
commit 710a94709e
4 changed files with 20 additions and 40 deletions

View File

@@ -43,18 +43,16 @@ namespace Opm
OilVaporizationProperties();
explicit OilVaporizationProperties(const size_t numPvtReginIdx);
OilVaporizationProperties(OilVaporization type,
const std::vector<double>& vap1,
const std::vector<double>& vap2,
double vap1,
double vap2,
const std::vector<double>& maxDRSDT,
const std::vector<bool>& maxDRSDT_allCells,
const std::vector<double>& maxDRVDT);
static void updateDRSDT(Opm::OilVaporizationProperties& ovp, const std::vector<double>& maxDRSDT, const std::vector<std::string>& option);
static void updateDRVDT(Opm::OilVaporizationProperties& ovp, const std::vector<double>& maxDRVDT);
static void updateVAPPARS(Opm::OilVaporizationProperties& ovp, const std::vector<double>& vap1, const std::vector<double>& vap2);
static void updateVAPPARS(Opm::OilVaporizationProperties& ovp, double vap1, double vap2);
OilVaporization getType() const;
double getVap1(const size_t pvtRegionIdx) const;
double getVap2(const size_t pvtRegionIdx) const;
double getMaxDRSDT(const size_t pvtRegionIdx) const;
double getMaxDRVDT(const size_t pvtRegionIdx) const;
bool getOption(const size_t pvtRegionIdx) const;
@@ -63,8 +61,8 @@ namespace Opm
bool defined() const;
size_t numPvtRegions() const {return m_maxDRSDT.size();}
const std::vector<double>& vap1() const;
const std::vector<double>& vap2() const;
double vap1() const;
double vap2() const;
const std::vector<double>& maxDRSDT() const;
const std::vector<bool>& maxDRSDT_allCells() const;
const std::vector<double>& maxDRVDT() const;
@@ -78,8 +76,8 @@ namespace Opm
private:
OilVaporization m_type = OilVaporization::UNDEF;
std::vector<double> m_vap1;
std::vector<double> m_vap2;
double m_vap1;
double m_vap2;
std::vector<double> m_maxDRSDT;
std::vector<bool> m_maxDRSDT_allCells;
std::vector<double> m_maxDRVDT;

View File

@@ -23,19 +23,20 @@ namespace Opm {
OilVaporizationProperties::OilVaporizationProperties()
{
m_type = OilVaporization::UNDEF;
m_vap1 = m_vap2 = -1.0;
}
OilVaporizationProperties::OilVaporizationProperties(const size_t numPvtRegionIdx):
m_vap1(numPvtRegionIdx, -1.0),
m_vap2(numPvtRegionIdx, -1.0),
m_vap1(-1.0),
m_vap2(-1.0),
m_maxDRSDT(numPvtRegionIdx, -1.0),
m_maxDRSDT_allCells(numPvtRegionIdx),
m_maxDRVDT(numPvtRegionIdx, -1.0)
{ }
OilVaporizationProperties::OilVaporizationProperties(OilVaporization type,
const std::vector<double>& vap1,
const std::vector<double>& vap2,
double vap1,
double vap2,
const std::vector<double>& maxDRSDT,
const std::vector<bool>& maxDRSDT_allCells,
const std::vector<double>& maxDRVDT):
@@ -75,22 +76,6 @@ namespace Opm {
return m_type;
}
double OilVaporizationProperties::getVap1(const size_t pvtRegionIdx) const{
if (m_type == OilVaporization::VAPPARS){
return m_vap1[pvtRegionIdx];
}else{
throw std::logic_error("Only valid if type is VAPPARS");
}
}
double OilVaporizationProperties::getVap2(const size_t pvtRegionIdx) const{
if (m_type == OilVaporization::VAPPARS){
return m_vap2[pvtRegionIdx];
}else{
throw std::logic_error("Only valid if type is VAPPARS");
}
}
void OilVaporizationProperties::updateDRSDT(OilVaporizationProperties& ovp, const std::vector<double>& maximums, const std::vector<std::string>& options){
ovp.m_type = OilVaporization::DRDT;
ovp.m_maxDRSDT = maximums;
@@ -110,7 +95,7 @@ namespace Opm {
ovp.m_maxDRVDT = maximums;
}
void OilVaporizationProperties::updateVAPPARS(OilVaporizationProperties& ovp, const std::vector<double>& vap1, const std::vector<double>& vap2){
void OilVaporizationProperties::updateVAPPARS(OilVaporizationProperties& ovp, double vap1, double vap2){
ovp.m_type = OilVaporization::VAPPARS;
ovp.m_vap1 = vap1;
ovp.m_vap2 = vap2;
@@ -152,11 +137,11 @@ namespace Opm {
return !(*this == rhs);
}
const std::vector<double>& OilVaporizationProperties::vap1() const {
double OilVaporizationProperties::vap1() const {
return m_vap1;
}
const std::vector<double>& OilVaporizationProperties::vap2() const {
double OilVaporizationProperties::vap2() const {
return m_vap2;
}

View File

@@ -753,14 +753,11 @@ std::pair<std::time_t, std::size_t> restart_info(const RestartIO::RstState * rst
}
void Schedule::handleVAPPARS( const DeckKeyword& keyword, size_t currentStep){
size_t numPvtRegions = m_runspec.tabdims().getNumPVTTables();
std::vector<double> vap(numPvtRegions);
std::vector<double> density(numPvtRegions);
for( const auto& record : keyword ) {
std::fill(vap.begin(), vap.end(), record.getItem("OIL_VAP_PROPENSITY").get< double >(0));
std::fill(density.begin(), density.end(), record.getItem("OIL_DENSITY_PROPENSITY").get< double >(0));
double vap1 = record.getItem("OIL_VAP_PROPENSITY").get< double >(0);
double vap2 = record.getItem("OIL_DENSITY_PROPENSITY").get< double >(0);
OilVaporizationProperties ovp = this->m_oilvaporizationproperties.get(currentStep);
OilVaporizationProperties::updateVAPPARS(ovp, vap, density);
OilVaporizationProperties::updateVAPPARS(ovp, vap1, vap2);
this->m_oilvaporizationproperties.update( currentStep, ovp );
}
}

View File

@@ -1465,9 +1465,9 @@ BOOST_AUTO_TEST_CASE(createDeckWithVAPPARS) {
BOOST_CHECK_EQUAL(schedule.hasOilVaporizationProperties(), true);
const OilVaporizationProperties& ovap = schedule.getOilVaporizationProperties(currentStep);
BOOST_CHECK(ovap.getType() == OilVaporizationProperties::OilVaporization::VAPPARS);
double vap1 = ovap.getVap1(0);
double vap1 = ovap.vap1();
BOOST_CHECK_EQUAL(2, vap1);
double vap2 = ovap.getVap2(0);
double vap2 = ovap.vap2();
BOOST_CHECK_EQUAL(0.100, vap2);
BOOST_CHECK_EQUAL(false, ovap.drsdtActive());
BOOST_CHECK_EQUAL(false, ovap.drvdtActive());