changed: avoid use of boost::multi_array in VFPInjTable

This commit is contained in:
Arne Morten Kvarving
2020-02-12 15:04:21 +01:00
parent de39bb5a7c
commit 143c8ec3bf
4 changed files with 39 additions and 53 deletions

View File

@@ -21,20 +21,18 @@
#ifndef OPM_PARSER_ECLIPSE_ECLIPSESTATE_TABLES_VFPINJTABLE_HPP_
#define OPM_PARSER_ECLIPSE_ECLIPSESTATE_TABLES_VFPINJTABLE_HPP_
#include <opm/parser/eclipse/Deck/Deck.hpp>
#include <boost/multi_array.hpp>
#include <array>
#include <vector>
namespace Opm {
class DeckKeyword;
class UnitSystem;
class VFPInjTable {
public:
typedef boost::multi_array<double, 2> array_type;
typedef boost::array<array_type::index, 2> extents;
typedef std::vector<double> array_type;
enum FLO_TYPE {
FLO_OIL=1,
@@ -75,12 +73,10 @@ public:
}
/**
* Returns the data of the table itself. The data is ordered so that
* Returns the data of the table itself. For ordered access
* use operator()(thp_idx, flo_idx)
*
* table = getTable();
* bhp = table[thp_idx][flo_idx];
*
* gives the bottom hole pressure value in the table for the coordinate
* This gives the bottom hole pressure value in the table for the coordinate
* given by
* flo_axis = getFloAxis();
* thp_axis = getTHPAxis();
@@ -93,10 +89,12 @@ public:
}
bool operator==(const VFPInjTable& data) const;
VFPInjTable& operator=(const VFPInjTable& data);
std::array<size_t,2> shape() const;
double operator()(size_t thp_idx, size_t flo_idx) const;
private:
int m_table_num;
double m_datum_depth;
FLO_TYPE m_flo_type;
@@ -108,6 +106,8 @@ private:
array_type m_data;
void check();
double& operator()(size_t thp_idx, size_t flo_idx);
static FLO_TYPE getFloType(std::string flo_string);
static void scaleValues(std::vector<double>& values,

View File

@@ -17,6 +17,8 @@
along with OPM. If not, see <http://www.gnu.org/licenses/>.
*/
#include <cassert>
#include <cmath>
#include <iostream>
#include <opm/parser/eclipse/Deck/DeckItem.hpp>
@@ -69,10 +71,6 @@ VFPInjTable::VFPInjTable(int table_num,
m_flo_data = flo_data;
m_thp_data = thp_data;
extents shape;
shape[0] = data.shape()[0];
shape[1] = data.shape()[1];
m_data.resize(shape);
m_data = data;
check();
@@ -159,11 +157,8 @@ VFPInjTable::VFPInjTable( const DeckKeyword& table, const UnitSystem& deck_unit_
//Finally, read the actual table itself.
size_t nt = m_thp_data.size();
size_t nf = m_flo_data.size();
extents shape;
shape[0] = nt;
shape[1] = nf;
m_data.resize(shape);
std::fill_n(m_data.data(), m_data.num_elements(), std::nan("0"));
m_data.resize(nt*nf);
std::fill_n(m_data.data(), m_data.size(), std::nan("0"));
//Check that size of table matches size of axis:
if (table.size() != nt + 3) {
@@ -190,7 +185,7 @@ VFPInjTable::VFPInjTable( const DeckKeyword& table, const UnitSystem& deck_unit_
std::cerr << "Too large value encountered in VFPINJ in ["
<< t << "," << f << "]=" << value << std::endl;
}
m_data[t][f] = table_scaling_factor*value;
(*this)(t,f) = table_scaling_factor*value;
}
}
@@ -221,15 +216,13 @@ void VFPInjTable::check() {
assert(is_sorted(m_thp_data.begin(), m_thp_data.end()));
//Check data size matches axes
assert(m_data.num_dimensions() == 2);
assert(m_data.shape()[0] == m_thp_data.size());
assert(m_data.shape()[1] == m_flo_data.size());
assert(m_data.size() == m_thp_data.size() * m_flo_data.size());
//Finally, check that all data is within reasonable ranges, defined to be up-to 1.0e10...
typedef array_type::size_type size_type;
for (size_type t=0; t<m_data.shape()[0]; ++t) {
for (size_type f=0; f<m_data.shape()[1]; ++f) {
if (std::isnan(m_data[t][f])) {
for (size_type t = 0; t < m_thp_data.size(); ++t) {
for (size_type f = 0; f < m_flo_data.size(); ++f) {
if (std::isnan((*this)(t,f))) {
//TODO: Replace with proper log message
std::cerr << "VFPINJ element [" << t << "," << f << "] not set!" << std::endl;
throw std::invalid_argument("Missing VFPINJ value");
@@ -328,24 +321,19 @@ bool VFPInjTable::operator==(const VFPInjTable& data) const {
}
VFPInjTable& VFPInjTable::operator=(const VFPInjTable& data) {
m_table_num = data.m_table_num;
m_datum_depth = data.m_datum_depth;
m_flo_type = data.m_flo_type;
m_flo_data = data.m_flo_data;
m_thp_data = data.m_thp_data;
extents shape;
shape[0] = data.m_data.shape()[0];
shape[1] = data.m_data.shape()[1];
m_data.resize(shape);
for (size_t i = 0; i < data.m_data.num_elements(); ++i)
*(m_data.data() + i) = *(data.m_data.data() + i);
return *this;
double VFPInjTable::operator()(size_t thp_idx, size_t flo_idx) const {
return m_data[thp_idx*m_flo_data.size() + flo_idx];
}
double& VFPInjTable::operator()(size_t thp_idx, size_t flo_idx) {
return m_data[thp_idx*m_flo_data.size() + flo_idx];
}
std::array<size_t,2> VFPInjTable::shape() const {
return {m_thp_data.size(), m_flo_data.size()};
}
} //Namespace

View File

@@ -3106,8 +3106,7 @@ VFPINJ \n \
//The data itself
{
typedef Opm::VFPInjTable::array_type::size_type size_type;
const Opm::VFPInjTable::array_type& data = vfpinjTable.getTable();
const size_type* size = data.shape();
const auto size = vfpinjTable.shape();
BOOST_CHECK_EQUAL(size[0], 2);
BOOST_CHECK_EQUAL(size[1], 3);
@@ -3119,7 +3118,7 @@ VFPINJ \n \
for (size_type t = 0; t < size[0]; ++t) {
for (size_type f = 0; f < size[1]; ++f) {
index += 1.0;
BOOST_CHECK_EQUAL(data[t][f], index*conversion_factor);
BOOST_CHECK_EQUAL(vfpinjTable(t,f), index*conversion_factor);
}
}
}

View File

@@ -1140,8 +1140,7 @@ VFPINJ \n\
//The data itself
{
typedef Opm::VFPInjTable::array_type::size_type size_type;
const Opm::VFPInjTable::array_type& data = vfpinjTable.getTable();
const size_type* size = data.shape();
const auto size = vfpinjTable.shape();
BOOST_CHECK_EQUAL(size[0], 2);
BOOST_CHECK_EQUAL(size[1], 3);
@@ -1153,7 +1152,7 @@ VFPINJ \n\
for (size_type t = 0; t < size[0]; ++t) {
for (size_type f = 0; f < size[1]; ++f) {
index += 1.0;
BOOST_CHECK_EQUAL(data[t][f], index*conversion_factor);
BOOST_CHECK_EQUAL(const_cast<const VFPInjTable&>(vfpinjTable)(t,f), index*conversion_factor);
}
}
}