Framework for adding tables to INIT file - PVTO.

This commit is contained in:
Joakim Hove 2016-11-16 16:12:59 +01:00
parent 82f1277fe8
commit 180057c062
5 changed files with 374 additions and 0 deletions

View File

@ -0,0 +1,51 @@
/*
Copyright 2016 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 OUTPUT_TABLES_HPP
#define OUTPUT_TABLES_HPP
#include <vector>
#include <ert/ecl/FortIO.hpp>
#include <ert/ecl/EclKW.hpp>
#include <opm/parser/eclipse/EclipseState/Tables/PvtoTable.hpp>
namespace Opm {
class UnitSystem;
class Tables {
public:
Tables( const UnitSystem& units_);
void fwrite( ERT::FortIO& fortio ) const;
void addPVTO( const std::vector<PvtoTable>& pvtoTables);
private:
void addData( size_t offset_index , const std::vector<double>& new_data);
const UnitSystem& units;
ERT::EclKW<int> tabdims;
std::vector<double> data;
};
}
#endif

View File

@ -39,6 +39,7 @@
#include <opm/parser/eclipse/EclipseState/Schedule/Well.hpp>
#include <opm/parser/eclipse/Utility/Functional.hpp>
#include <opm/output/eclipse/Summary.hpp>
#include <opm/output/eclipse/Tables.hpp>
#include <opm/output/eclipse/RegionCache.hpp>
#include <cstdlib>
@ -507,6 +508,12 @@ void EclipseWriter::Impl::writeINITFile( const data::Solution& simProps, const N
}
}
// Write tables
{
Tables tables( this->es.getUnits() );
tables.addPVTO( this->es.getTableManager().getPvtoTables() );
tables.fwrite( fortio );
}
// Write all integer field properties from the input deck.
{

View File

@ -0,0 +1,133 @@
/*
Copyright 2016 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/>.
*/
#include <ert/ecl/FortIO.hpp>
#include <ert/ecl/EclKW.hpp>
#include <ert/ecl/ecl_kw_magic.h>
#include <opm/parser/eclipse/Units/UnitSystem.hpp>
#include <opm/output/eclipse/Tables.hpp>
namespace Opm {
Tables::Tables( const UnitSystem& units_ ) :
units( units_ ),
tabdims( "TABDIMS" , TABDIMS_SIZE )
{
}
void Tables::addData( size_t offset_index, const std::vector<double>& new_data) {
this->tabdims[ offset_index ] = this->data.size();
this->data.insert( this->data.end() , new_data.begin() , new_data.end());
}
namespace {
struct PvtxDims {
size_t num_tables;
size_t outer_size;
size_t inner_size;
size_t num_columns = 3;
size_t data_size;
};
template <class TableType>
PvtxDims tableDims( const std::vector<TableType>& pvtxTables) {
PvtxDims dims;
dims.num_tables = pvtxTables.size();
dims.inner_size = 0;
dims.outer_size = 0;
for (const auto& table : pvtxTables) {
dims.outer_size = std::max( dims.outer_size, table.size());
for (const auto& underSatTable : table)
dims.inner_size = std::max(dims.inner_size, underSatTable.numRows() );
}
dims.data_size = dims.num_tables * dims.outer_size * dims.inner_size * dims.num_columns;
return dims;
}
}
void Tables::addPVTO( const std::vector<PvtoTable>& pvtoTables)
{
const double default_value = 2e20;
PvtxDims dims = tableDims( pvtoTables );
this->tabdims[ TABDIMS_NTPVTO_ITEM ] = dims.num_tables;
this->tabdims[ TABDIMS_NRPVTO_ITEM ] = dims.outer_size;
this->tabdims[ TABDIMS_NPPVTO_ITEM ] = dims.inner_size;
{
std::vector<double> pvtoData( dims.data_size , default_value );
std::vector<double> rs_values( dims.num_tables * dims.outer_size , default_value );
size_t composition_stride = dims.inner_size;
size_t table_stride = dims.outer_size * composition_stride;
size_t column_stride = table_stride * pvtoTables.size();
size_t table_index = 0;
for (const auto& table : pvtoTables) {
size_t composition_index = 0;
for (const auto& underSatTable : table) {
const auto& p = underSatTable.getColumn("P");
const auto& bo = underSatTable.getColumn("BO");
const auto& mu = underSatTable.getColumn("MU");
for (size_t row = 0; row < p.size(); row++) {
size_t data_index = row + composition_stride * composition_index + table_stride * table_index;
pvtoData[ data_index ] = units.from_si( UnitSystem::measure::pressure, p[row]);
pvtoData[ data_index + column_stride ] = 1.0 / bo[row];
pvtoData[ data_index + 2*column_stride] = units.from_si( UnitSystem::measure::viscosity , mu[row]) / bo[row];
}
composition_index++;
}
/*
The RS values which apply for one inner table each
are added as a separate data vector to the TABS
array.
*/
{
const auto& sat_table = table.getSaturatedTable();
const auto& rs = sat_table.getColumn("RS");
for (size_t index = 0; index < rs.size(); index++)
rs_values[index + table_index * dims.outer_size ] = rs[index];
}
table_index++;
}
addData( TABDIMS_IBPVTO_OFFSET_ITEM , pvtoData );
addData( TABDIMS_JBPVTO_OFFSET_ITEM , rs_values );
}
}
void Tables::fwrite( ERT::FortIO& fortio ) const
{
tabdims.fwrite( fortio );
{
ERT::EclKW<double> tab( "TAB" , this->data );
tab.fwrite( fortio );
}
}
}

72
tests/table_deck.DATA Normal file
View File

@ -0,0 +1,72 @@
RUNSPEC
TITLE
SIMPLE TEST
DIMENS
10 10 10 /
OIL
WATER
GAS
DISGAS
VAPOIL
METRIC
EQLDIMS
1 100 2 1 1 /
TABDIMS
1 1 33 60 16 60 /
GRID
DX
1000*1 /
DY
1000*1 /
DZ
1000*1 /
TOPS
100*1 /
PROPS ===============================================================
PVTO
-- RSO PRESSURE B-OIL VISCOSITY
-- (BAR) (CP)
20.59 50.00 1.10615 1.180
75.00 1.10164 1.247
100.00 1.09744 1.315
125.00 1.09351 1.384
150.00 1.08984 1.453 /
28.19 70.00 1.12522 1.066
95.00 1.12047 1.124
120.00 1.11604 1.182
170.00 1.10804 1.300 /
/
PVTG
--
20.00 0.00002448 0.061895 0.01299
0.00001224 0.061810 0.01300
0.00000000 0.061725 0.01300 /
40.00 0.00000628 0.030252 0.01383
0.00000314 0.030249 0.01383
0.00000000 0.030245 0.01383 /
/

111
tests/test_Tables.cpp Normal file
View File

@ -0,0 +1,111 @@
/*
Copyright 2016 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/>.
*/
#include "config.h"
#if HAVE_DYNAMIC_BOOST_TEST
#define BOOST_TEST_DYN_LINK
#endif
#define BOOST_TEST_MODULE Wells
#include <boost/test/unit_test.hpp>
#include <boost/date_time/posix_time/posix_time.hpp>
#include <stdexcept>
#include <fstream>
#include <ert/ecl/ecl_kw_magic.h>
#include <ert/ecl/ecl_sum.h>
#include <ert/ecl/ecl_file.h>
#include <ert/util/util.h>
#include <ert/util/TestArea.hpp>
#include <opm/output/data/Wells.hpp>
#include <opm/output/data/Cells.hpp>
#include <opm/output/eclipse/Tables.hpp>
#include <opm/parser/eclipse/Deck/Deck.hpp>
#include <opm/parser/eclipse/Units/UnitSystem.hpp>
#include <opm/parser/eclipse/EclipseState/Grid/EclipseGrid.hpp>
#include <opm/parser/eclipse/EclipseState/EclipseState.hpp>
#include <opm/parser/eclipse/EclipseState/SummaryConfig/SummaryConfig.hpp>
#include <opm/parser/eclipse/Parser/ParseContext.hpp>
#include <opm/parser/eclipse/Parser/Parser.hpp>
using namespace Opm;
struct setup {
Deck deck;
EclipseState es;
ERT::TestArea ta;
setup( const std::string& name, const std::string& path , const ParseContext& parseContext = ParseContext( )) :
deck( Parser().parseFile( path, parseContext ) ),
es( deck, ParseContext() ),
ta( ERT::TestArea("test_tables") )
{
}
};
BOOST_AUTO_TEST_CASE(Test_PVTO) {
setup cfg( "PVTO" , "table_deck.DATA");
Tables tables( cfg.es.getUnits() );
tables.addPVTO( cfg.es.getTableManager().getPvtoTables() );
{
ERT::FortIO f("TEST.INIT" , std::fstream::out);
tables.fwrite( f );
}
{
ecl_file_type * f = ecl_file_open("TEST.INIT" , 0 );
BOOST_CHECK( ecl_file_has_kw( f , "TABDIMS" ));
BOOST_CHECK( ecl_file_has_kw( f , "TAB" ));
{
const ecl_kw_type * tabdims = ecl_file_iget_named_kw( f , "TABDIMS" , 0 );
const ecl_kw_type * tab = ecl_file_iget_named_kw( f , "TAB" , 0 );
int offset = ecl_kw_iget_int( tabdims , TABDIMS_IBPVTO_OFFSET_ITEM );
int rs_offset = ecl_kw_iget_int( tabdims , TABDIMS_JBPVTO_OFFSET_ITEM );
int column_stride = ecl_kw_iget_int( tabdims , TABDIMS_NRPVTO_ITEM ) * ecl_kw_iget_int( tabdims , TABDIMS_NPPVTO_ITEM ) * ecl_kw_iget_int( tabdims , TABDIMS_NTPVTO_ITEM );
BOOST_CHECK_EQUAL( 2, ecl_kw_iget_int( tabdims , TABDIMS_NRPVTO_ITEM ) );
BOOST_CHECK_EQUAL( 5, ecl_kw_iget_int( tabdims , TABDIMS_NPPVTO_ITEM ) );
BOOST_CHECK_EQUAL( 1, ecl_kw_iget_int( tabdims , TABDIMS_NTPVTO_ITEM ) );
BOOST_CHECK_CLOSE(50.0 , ecl_kw_iget_double( tab , offset ), 1e-6 );
BOOST_CHECK_CLOSE(1.0 / 1.10615 , ecl_kw_iget_double( tab , offset + column_stride), 1e-6 );
BOOST_CHECK_CLOSE(1.18 / 1.10615 , ecl_kw_iget_double( tab , offset + 2*column_stride ), 1e-6 );
BOOST_CHECK_CLOSE(150.0 , ecl_kw_iget_double( tab , 4 + offset ), 1e-6 );
BOOST_CHECK_CLOSE(1.0 / 1.08984 , ecl_kw_iget_double( tab , 4 + offset + column_stride), 1e-6 );
BOOST_CHECK_CLOSE(1.453 / 1.08984 , ecl_kw_iget_double( tab , 4 + offset + 2*column_stride ), 1e-6 );
BOOST_CHECK_CLOSE(20.59 , ecl_kw_iget_double( tab , rs_offset ), 1e-6 );
BOOST_CHECK_CLOSE(28.19 , ecl_kw_iget_double( tab , rs_offset + 1), 1e-6 );
}
ecl_file_close( f );
}
}