Treat All-Defaulted PVTW Record as Copy of Previous

A simulation model may choose to give PVTW data as

    PVTW
      1.0 1.0 1.0e-5 0.2 0.0 /
      /  -- record 2 (copied from record 1)

if, for instance, the oil and/or gas tables are different in regions
1 and 2, but the water is the same.  In this case we must properly
copy record 1 into record 2 and essentially recreate the table.

To this end, decouple the 'PvtwTable' from the 'FlatTable' machinery
and make the former into an independent type containing vector<>
instead of inheriting from vector<>.  Implement the default->copy
behaviour in the new PvtwTable::PvtwTable(const DeckKeyword&)
constructor.
This commit is contained in:
Bård Skaflestad 2022-06-17 17:52:46 +02:00
parent 00ace58e6c
commit 433cc4d649
3 changed files with 154 additions and 3 deletions

View File

@ -1,6 +1,10 @@
#ifndef OPM_FLAT_TABLE_HPP
#define OPM_FLAT_TABLE_HPP
#include <cstddef>
#include <initializer_list>
#include <vector>
namespace Opm {
class DeckKeyword;
@ -156,13 +160,44 @@ struct PVTWRecord {
}
};
struct PvtwTable : public FlatTable< PVTWRecord > {
using FlatTable< PVTWRecord >::FlatTable;
class PvtwTable
{
public:
PvtwTable() = default;
explicit PvtwTable(const DeckKeyword& kw);
explicit PvtwTable(std::initializer_list<PVTWRecord> records);
auto size() const { return this->table_.size(); }
bool empty() const { return this->table_.empty(); }
const PVTWRecord& operator[](const std::size_t tableID) const
{
return this->table_[tableID];
}
const PVTWRecord& at(const std::size_t tableID) const
{
return this->table_.at(tableID);
}
bool operator==(const PvtwTable& other) const
{
return this->table_ == other.table_;
}
static PvtwTable serializeObject()
{
return PvtwTable({{1.0, 2.0, 3.0, 4.0, 5.0}});
}
template <class Serializer>
void serializeOp(Serializer& serializer)
{
serializer.vector(this->table_);
}
private:
std::vector<PVTWRecord> table_{};
};
struct ROCKRecord {

View File

@ -86,6 +86,16 @@
#include <opm/input/eclipse/EclipseState/Tables/WatvisctTable.hpp>
#include <opm/input/eclipse/EclipseState/Tables/AqutabTable.hpp>
#include <opm/common/utility/OpmInputError.hpp>
#include <algorithm>
#include <cstddef>
#include <stdexcept>
#include <stddef.h>
#include <fmt/format.h>
namespace Opm {
PvtgTable::PvtgTable( const DeckKeyword& keyword, size_t tableIdx ) :
@ -1547,8 +1557,58 @@ struct flat_props< PVCDORecord, N > {
}
};
bool all_defaulted(const DeckRecord& record)
{
return std::all_of(record.begin(), record.end(),
[](const DeckItem& item)
{
const auto& vstat = item.getValueStatus();
return std::all_of(vstat.begin(), vstat.end(), &value::defaulted);
});
}
} // Anonymous namespace
// ------------------------------------------------------------------------
PvtwTable::PvtwTable(const DeckKeyword& kw)
{
if (kw.name() != ParserKeywords::PVTW::keywordName) {
throw std::invalid_argument {
fmt::format("Keyword {} cannot be used to "
"initialise {} table structures", kw.name(),
ParserKeywords::PVTW::keywordName)
};
}
this->table_.reserve(kw.size());
for (const auto& record : kw) {
if (all_defaulted(record)) {
// All-defaulted records imply PVTW in region R is equal to PVTW
// in region R-1. PVTW must not be defaulted in region 1 (i.e.,
// when PVTNUM=1).
if (this->table_.empty()) {
throw OpmInputError {
"First record cannot be defaulted",
kw.location()
};
}
this->table_.push_back(this->table_.back());
}
else {
this->table_.push_back(flat_get<PVTWRecord>(record, mkseq<PVTWRecord::size>{}));
}
}
}
PvtwTable::PvtwTable(std::initializer_list<PVTWRecord> records)
: table_(records)
{}
// ------------------------------------------------------------------------
template< typename T >
FlatTable< T >::FlatTable( const DeckKeyword& kw ) :
std::vector< T >( flat_records< T >( kw, mkseq< T::size >{} ) )
@ -1557,7 +1617,6 @@ FlatTable< T >::FlatTable( const DeckKeyword& kw ) :
template FlatTable< DENSITYRecord >::FlatTable( const DeckKeyword& );
template FlatTable< GRAVITYRecord >::FlatTable( const DeckKeyword& );
template FlatTable< DiffCoeffRecord >::FlatTable( const DeckKeyword& );
template FlatTable< PVTWRecord >::FlatTable( const DeckKeyword& );
template FlatTable< PVCDORecord >::FlatTable( const DeckKeyword& );
template FlatTable< ROCKRecord >::FlatTable( const DeckKeyword& );
template FlatTable< PlyvmhRecord >::FlatTable( const DeckKeyword& );

View File

@ -38,6 +38,7 @@
#include <opm/input/eclipse/EclipseState/Tables/Tabdims.hpp>
#include <opm/input/eclipse/EclipseState/Tables/PlyadsTable.hpp>
#include <opm/input/eclipse/EclipseState/Tables/PlymaxTable.hpp>
#include <opm/input/eclipse/EclipseState/Tables/FlatTable.hpp>
#include <opm/input/eclipse/EclipseState/Tables/FoamadsTable.hpp>
#include <opm/input/eclipse/EclipseState/Tables/FoammobTable.hpp>
#include <opm/input/eclipse/EclipseState/Tables/PbvdTable.hpp>
@ -653,6 +654,62 @@ BOOST_AUTO_TEST_CASE(FoammobTable_Tests) {
BOOST_AUTO_TEST_CASE(PvtwTable_Tests) {
// PVT tables from opm-tests/model5/include/pvt_live_oil_dgas.ecl .
const auto deck = Opm::Parser{}.parseString(R"(RUNSPEC
OIL
WATER
TABDIMS
1 2 /
PROPS
DENSITY
924.1 1026.0 1.03446 /
924.1 1026.0 1.03446 /
PVTW
79.0 1.02643 0.37876E-04 0.39831 0.74714E-04 /
/
END
)");
const auto tmgr = Opm::TableManager { deck };
const auto& pvtw = tmgr.getPvtwTable();
BOOST_REQUIRE_EQUAL(pvtw.size(), std::size_t{2});
{
const auto& t1 = pvtw[0];
BOOST_CHECK_CLOSE(t1.reference_pressure, 7.9e6, 1.0e-8);
BOOST_CHECK_CLOSE(t1.volume_factor, 1.02643, 1.0e-8);
BOOST_CHECK_CLOSE(t1.compressibility, 0.37876e-9, 1.0e-8);
BOOST_CHECK_CLOSE(t1.viscosity, 0.39831e-3, 1.0e-8);
BOOST_CHECK_CLOSE(t1.viscosibility, 0.74714e-9, 1.0e-9);
}
{
const auto& t2 = pvtw[1];
BOOST_CHECK_CLOSE(t2.reference_pressure, 7.9e6, 1.0e-8);
BOOST_CHECK_CLOSE(t2.volume_factor, 1.02643, 1.0e-8);
BOOST_CHECK_CLOSE(t2.compressibility, 0.37876e-9, 1.0e-8);
BOOST_CHECK_CLOSE(t2.viscosity, 0.39831e-3, 1.0e-8);
BOOST_CHECK_CLOSE(t2.viscosibility, 0.74714e-9, 1.0e-9);
}
const auto& dens = tmgr.getDensityTable();
BOOST_REQUIRE_EQUAL(dens.size(), std::size_t{2});
{
const auto& t1 = dens[0];
BOOST_CHECK_CLOSE(t1.oil, 924.1, 1.0e-8);
BOOST_CHECK_CLOSE(t1.gas, 1.03446, 1.0e-8);
BOOST_CHECK_CLOSE(t1.water, 1026.0, 1.0e-8);
}
{
const auto& t2 = dens[1];
BOOST_CHECK_CLOSE(t2.oil, 924.1, 1.0e-8);
BOOST_CHECK_CLOSE(t2.gas, 1.03446, 1.0e-8);
BOOST_CHECK_CLOSE(t2.water, 1026.0, 1.0e-8);
}
}
/**
* Tests "happy path" for a VFPPROD table, i.e., when everything goes well