2016-05-13 09:25:18 -05:00
|
|
|
|
|
|
|
/*
|
|
|
|
Copyright 2014 Statoil IT
|
|
|
|
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 EclipseWriter
|
|
|
|
#include <boost/test/unit_test.hpp>
|
|
|
|
|
|
|
|
#include <opm/output/eclipse/EclipseWriter.hpp>
|
|
|
|
#include <opm/output/eclipse/EclipseReader.hpp>
|
2016-10-03 07:52:39 -05:00
|
|
|
#include <opm/output/data/Cells.hpp>
|
|
|
|
#include <opm/output/data/Wells.hpp>
|
2016-05-13 09:25:18 -05:00
|
|
|
#include <opm/parser/eclipse/Deck/Deck.hpp>
|
|
|
|
#include <opm/parser/eclipse/EclipseState/Grid/EclipseGrid.hpp>
|
|
|
|
#include <opm/parser/eclipse/EclipseState/IOConfig/IOConfig.hpp>
|
|
|
|
#include <opm/parser/eclipse/EclipseState/Schedule/Schedule.hpp>
|
|
|
|
#include <opm/parser/eclipse/EclipseState/Schedule/Well.hpp>
|
|
|
|
#include <opm/parser/eclipse/Parser/ParseContext.hpp>
|
|
|
|
#include <opm/parser/eclipse/Parser/Parser.hpp>
|
|
|
|
#include <opm/parser/eclipse/Utility/Functional.hpp>
|
|
|
|
|
|
|
|
// ERT stuff
|
|
|
|
#include <ert/ecl/ecl_kw.h>
|
|
|
|
#include <ert/ecl/ecl_file.h>
|
2016-09-19 12:28:06 -05:00
|
|
|
#include <ert/ecl/ecl_util.h>
|
2016-05-13 09:25:18 -05:00
|
|
|
#include <ert/ecl/ecl_kw_magic.h>
|
|
|
|
#include <ert/ecl_well/well_info.h>
|
|
|
|
#include <ert/ecl_well/well_state.h>
|
|
|
|
#include <ert/util/test_work_area.h>
|
|
|
|
|
|
|
|
using namespace Opm;
|
|
|
|
|
2016-09-25 11:31:34 -05:00
|
|
|
inline std::string input( const std::string& rst_name = "FIRST_SIM" ) {
|
|
|
|
return std::string(
|
2016-05-13 09:25:18 -05:00
|
|
|
"RUNSPEC\n"
|
|
|
|
"OIL\n"
|
|
|
|
"GAS\n"
|
|
|
|
"WATER\n"
|
|
|
|
"DISGAS\n"
|
|
|
|
"VAPOIL\n"
|
|
|
|
"UNIFOUT\n"
|
|
|
|
"UNIFIN\n"
|
|
|
|
"DIMENS\n"
|
|
|
|
" 10 10 10 /\n"
|
|
|
|
|
|
|
|
"GRID\n"
|
|
|
|
"DXV\n"
|
|
|
|
"10*0.25 /\n"
|
|
|
|
"DYV\n"
|
|
|
|
"10*0.25 /\n"
|
|
|
|
"DZV\n"
|
|
|
|
"10*0.25 /\n"
|
|
|
|
"TOPS\n"
|
|
|
|
"100*0.25 /\n"
|
|
|
|
"\n"
|
|
|
|
|
|
|
|
"SOLUTION\n"
|
|
|
|
"RESTART\n"
|
2016-09-25 11:31:34 -05:00
|
|
|
) + rst_name + std::string(
|
|
|
|
" 1/\n"
|
2016-05-13 09:25:18 -05:00
|
|
|
"\n"
|
|
|
|
|
|
|
|
"START -- 0 \n"
|
|
|
|
"1 NOV 1979 / \n"
|
|
|
|
|
|
|
|
"SCHEDULE\n"
|
|
|
|
"SKIPREST\n"
|
|
|
|
"RPTRST\n"
|
|
|
|
"BASIC=1\n"
|
|
|
|
"/\n"
|
|
|
|
"DATES -- 1\n"
|
|
|
|
" 10 OKT 2008 / \n"
|
|
|
|
"/\n"
|
|
|
|
"WELSPECS\n"
|
|
|
|
" 'OP_1' 'OP' 9 9 1* 'OIL' 1* 1* 1* 1* 1* 1* 1* / \n"
|
|
|
|
" 'OP_2' 'OP' 9 9 1* 'OIL' 1* 1* 1* 1* 1* 1* 1* / \n"
|
|
|
|
"/\n"
|
|
|
|
"COMPDAT\n"
|
|
|
|
" 'OP_1' 9 9 1 1 'OPEN' 1* 32.948 0.311 3047.839 1* 1* 'X' 22.100 / \n"
|
|
|
|
" 'OP_2' 9 9 2 2 'OPEN' 1* 46.825 0.311 4332.346 1* 1* 'X' 22.123 / \n"
|
|
|
|
" 'OP_1' 9 9 3 9 'OPEN' 1* 32.948 0.311 3047.839 1* 1* 'X' 22.100 / \n"
|
|
|
|
"/\n"
|
|
|
|
"WCONPROD\n"
|
|
|
|
"'OP_1' 'OPEN' 'ORAT' 20000 4* 1000 /\n"
|
|
|
|
"/\n"
|
|
|
|
"WCONINJE\n"
|
|
|
|
"'OP_2' 'GAS' 'OPEN' 'RATE' 100 200 400 /\n"
|
|
|
|
"/\n"
|
|
|
|
|
|
|
|
"DATES -- 2\n"
|
|
|
|
" 20 JAN 2011 / \n"
|
|
|
|
"/\n"
|
|
|
|
"WELSPECS\n"
|
|
|
|
" 'OP_3' 'OP' 9 9 1* 'OIL' 1* 1* 1* 1* 1* 1* 1* / \n"
|
|
|
|
"/\n"
|
|
|
|
"COMPDAT\n"
|
|
|
|
" 'OP_3' 9 9 1 1 'OPEN' 1* 32.948 0.311 3047.839 1* 1* 'X' 22.100 / \n"
|
|
|
|
"/\n"
|
|
|
|
"WCONPROD\n"
|
|
|
|
"'OP_3' 'OPEN' 'ORAT' 20000 4* 1000 /\n"
|
|
|
|
"/\n"
|
|
|
|
|
|
|
|
"DATES -- 3\n"
|
|
|
|
" 15 JUN 2013 / \n"
|
|
|
|
"/\n"
|
|
|
|
"COMPDAT\n"
|
|
|
|
" 'OP_2' 9 9 3 9 'OPEN' 1* 32.948 0.311 3047.839 1* 1* 'X' 22.100 / \n"
|
|
|
|
" 'OP_1' 9 9 7 7 'SHUT' 1* 32.948 0.311 3047.839 1* 1* 'X' 22.100 / \n"
|
|
|
|
"/\n"
|
|
|
|
|
|
|
|
"DATES -- 4\n"
|
|
|
|
" 22 APR 2014 / \n"
|
|
|
|
"/\n"
|
|
|
|
"WELSPECS\n"
|
|
|
|
" 'OP_4' 'OP' 9 9 1* 'OIL' 1* 1* 1* 1* 1* 1* 1* / \n"
|
|
|
|
"/\n"
|
|
|
|
"COMPDAT\n"
|
|
|
|
" 'OP_4' 9 9 3 9 'OPEN' 1* 32.948 0.311 3047.839 1* 1* 'X' 22.100 / \n"
|
|
|
|
" 'OP_3' 9 9 3 9 'OPEN' 1* 32.948 0.311 3047.839 1* 1* 'X' 22.100 / \n"
|
|
|
|
"/\n"
|
|
|
|
"WCONPROD\n"
|
|
|
|
"'OP_4' 'OPEN' 'ORAT' 20000 4* 1000 /\n"
|
|
|
|
"/\n"
|
|
|
|
|
|
|
|
"DATES -- 5\n"
|
|
|
|
" 30 AUG 2014 / \n"
|
|
|
|
"/\n"
|
|
|
|
"WELSPECS\n"
|
|
|
|
" 'OP_5' 'OP' 9 9 1* 'OIL' 1* 1* 1* 1* 1* 1* 1* / \n"
|
|
|
|
"/\n"
|
|
|
|
"COMPDAT\n"
|
|
|
|
" 'OP_5' 9 9 3 9 'OPEN' 1* 32.948 0.311 3047.839 1* 1* 'X' 22.100 / \n"
|
|
|
|
"/\n"
|
|
|
|
"WCONPROD\n"
|
|
|
|
"'OP_5' 'OPEN' 'ORAT' 20000 4* 1000 /\n"
|
|
|
|
"/\n"
|
|
|
|
|
|
|
|
"DATES -- 6\n"
|
|
|
|
" 15 SEP 2014 / \n"
|
|
|
|
"/\n"
|
|
|
|
"WCONPROD\n"
|
|
|
|
"'OP_3' 'SHUT' 'ORAT' 20000 4* 1000 /\n"
|
|
|
|
"/\n"
|
|
|
|
|
|
|
|
"DATES -- 7\n"
|
|
|
|
" 9 OCT 2014 / \n"
|
|
|
|
"/\n"
|
|
|
|
"WELSPECS\n"
|
|
|
|
" 'OP_6' 'OP' 9 9 1* 'OIL' 1* 1* 1* 1* 1* 1* 1* / \n"
|
|
|
|
"/\n"
|
|
|
|
"COMPDAT\n"
|
|
|
|
" 'OP_6' 9 9 3 9 'OPEN' 1* 32.948 0.311 3047.839 1* 1* 'X' 22.100 / \n"
|
|
|
|
"/\n"
|
|
|
|
"WCONPROD\n"
|
|
|
|
"'OP_6' 'OPEN' 'ORAT' 20000 4* 1000 /\n"
|
|
|
|
"/\n"
|
|
|
|
"TSTEP -- 8\n"
|
|
|
|
"10 /"
|
2016-09-25 11:31:34 -05:00
|
|
|
"/\n"
|
|
|
|
);
|
|
|
|
}
|
2016-05-13 09:25:18 -05:00
|
|
|
|
2016-09-25 11:31:34 -05:00
|
|
|
namespace Opm {
|
|
|
|
namespace data {
|
2016-05-13 09:25:18 -05:00
|
|
|
|
2016-09-25 11:31:34 -05:00
|
|
|
/*
|
|
|
|
* Some test specific equivalence definitions and pretty-printing. Not fit as a
|
|
|
|
* general purpose implementation, but does its job for testing and
|
|
|
|
* pretty-pringing for debugging purposes.
|
|
|
|
*/
|
|
|
|
|
|
|
|
std::ostream& operator<<( std::ostream& stream, const Rates& r ) {
|
|
|
|
return stream << "{ "
|
|
|
|
<< "wat: " << r.get( Rates::opt::wat, 0.0 ) << ", "
|
|
|
|
<< "oil: " << r.get( Rates::opt::oil, 0.0 ) << ", "
|
|
|
|
<< "gas: " << r.get( Rates::opt::gas, 0.0 ) << " "
|
|
|
|
<< "}";
|
|
|
|
}
|
|
|
|
|
|
|
|
std::ostream& operator<<( std::ostream& stream, const Completion& c ) {
|
|
|
|
return stream << "{ index: "
|
|
|
|
<< c.index << ", "
|
|
|
|
<< c.rates << ", "
|
|
|
|
<< c.pressure << " }";
|
|
|
|
}
|
|
|
|
|
|
|
|
std::ostream& operator<<( std::ostream& stream,
|
|
|
|
const std::map< std::string, Well >& m ) {
|
|
|
|
stream << "\n";
|
|
|
|
|
|
|
|
for( const auto& p : m ) {
|
|
|
|
stream << p.first << ": \n"
|
|
|
|
<< "\t" << "bhp: " << p.second.bhp << "\n"
|
|
|
|
<< "\t" << "temp: " << p.second.temperature << "\n"
|
|
|
|
<< "\t" << "rates: " << p.second.rates << "\n"
|
|
|
|
<< "\t" << "completions: [\n";
|
|
|
|
|
|
|
|
for( const auto& c : p.second.completions )
|
|
|
|
stream << c.second << " ";
|
|
|
|
|
|
|
|
stream << "]\n";
|
|
|
|
}
|
|
|
|
|
|
|
|
return stream;
|
|
|
|
}
|
|
|
|
|
|
|
|
bool operator==( const Rates& lhs, const Rates& rhs ) {
|
|
|
|
using rt = Rates::opt;
|
|
|
|
|
|
|
|
BOOST_CHECK_EQUAL( lhs.has( rt::wat ), rhs.has( rt::wat ) );
|
|
|
|
BOOST_CHECK_EQUAL( lhs.has( rt::oil ), rhs.has( rt::oil ) );
|
|
|
|
BOOST_CHECK_EQUAL( lhs.has( rt::gas ), rhs.has( rt::gas ) );
|
|
|
|
BOOST_CHECK_EQUAL( lhs.has( rt::polymer ), rhs.has( rt::polymer ) );
|
|
|
|
BOOST_CHECK_EQUAL( lhs.get( rt::wat, 0.0 ), rhs.get( rt::wat, 0.0 ) );
|
|
|
|
BOOST_CHECK_EQUAL( lhs.get( rt::oil, 0.0 ), rhs.get( rt::oil, 0.0 ) );
|
|
|
|
BOOST_CHECK_EQUAL( lhs.get( rt::gas, 0.0 ), rhs.get( rt::gas, 0.0 ) );
|
|
|
|
BOOST_CHECK_EQUAL( lhs.get( rt::polymer, 0.0 ), rhs.get( rt::polymer, 0.0 ) );
|
|
|
|
|
|
|
|
return true;
|
|
|
|
}
|
|
|
|
|
|
|
|
bool operator==( const Completion& lhs, const Completion& rhs ) {
|
|
|
|
BOOST_CHECK_EQUAL( lhs.index, rhs.index );
|
|
|
|
BOOST_CHECK_EQUAL( lhs.rates, rhs.rates );
|
|
|
|
BOOST_CHECK_EQUAL( lhs.pressure, rhs.pressure );
|
|
|
|
BOOST_CHECK_EQUAL( lhs.reservoir_rate, rhs.reservoir_rate );
|
2016-05-13 09:25:18 -05:00
|
|
|
|
2016-09-25 11:31:34 -05:00
|
|
|
return true;
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
bool operator==( const Well& lhs, const Well& rhs ) {
|
|
|
|
BOOST_CHECK_EQUAL( lhs.rates, rhs.rates );
|
|
|
|
BOOST_CHECK_EQUAL( lhs.bhp, rhs.bhp );
|
|
|
|
BOOST_CHECK_EQUAL( lhs.temperature, rhs.temperature );
|
2016-09-29 04:58:10 -05:00
|
|
|
BOOST_CHECK_EQUAL( lhs.control, rhs.control );
|
2016-09-25 11:31:34 -05:00
|
|
|
|
|
|
|
for( const auto& p : lhs.completions )
|
|
|
|
BOOST_CHECK_EQUAL( p.second, rhs.completions.at( p.first ) );
|
|
|
|
|
|
|
|
return true;
|
|
|
|
}
|
|
|
|
|
|
|
|
}
|
|
|
|
|
|
|
|
/*
|
|
|
|
* forward declarations of internal functions that we want to expose to tests
|
|
|
|
* but not to users.
|
|
|
|
*/
|
2016-09-29 04:58:10 -05:00
|
|
|
std::vector< double > serialize_XWEL( const data::Wells& wells,
|
|
|
|
int report_step,
|
|
|
|
const std::vector< const Well* > sched_wells,
|
|
|
|
const TableManager& tm,
|
|
|
|
const EclipseGrid& );
|
|
|
|
|
|
|
|
std::vector< int > serialize_IWEL( const data::Wells& wells,
|
|
|
|
const std::vector< const Well* > sched_wells );
|
|
|
|
|
|
|
|
data::Wells restore_wells( const double* xwel_data,
|
|
|
|
size_t xwel_data_size,
|
|
|
|
const int* iwel_data,
|
|
|
|
size_t iwel_data_size,
|
|
|
|
int restart_step,
|
|
|
|
const std::vector< const Well* > sched_wells,
|
|
|
|
const std::vector< data::Rates::opt >& phases,
|
|
|
|
const EclipseGrid& grid );
|
2016-09-25 11:31:34 -05:00
|
|
|
}
|
|
|
|
|
|
|
|
data::Wells mkWells() {
|
|
|
|
data::Rates r1, r2, rc1, rc2, rc3;
|
|
|
|
r1.set( data::Rates::opt::wat, 5.67 );
|
|
|
|
r1.set( data::Rates::opt::oil, 6.78 );
|
|
|
|
r1.set( data::Rates::opt::gas, 7.89 );
|
|
|
|
|
|
|
|
r2.set( data::Rates::opt::wat, 8.90 );
|
|
|
|
r2.set( data::Rates::opt::oil, 9.01 );
|
|
|
|
r2.set( data::Rates::opt::gas, 10.12 );
|
|
|
|
|
|
|
|
rc1.set( data::Rates::opt::wat, 20.41 );
|
|
|
|
rc1.set( data::Rates::opt::oil, 21.19 );
|
|
|
|
rc1.set( data::Rates::opt::gas, 22.41 );
|
|
|
|
|
|
|
|
rc2.set( data::Rates::opt::wat, 23.19 );
|
|
|
|
rc2.set( data::Rates::opt::oil, 24.41 );
|
|
|
|
rc2.set( data::Rates::opt::gas, 25.19 );
|
|
|
|
|
|
|
|
rc3.set( data::Rates::opt::wat, 26.41 );
|
|
|
|
rc3.set( data::Rates::opt::oil, 27.19 );
|
|
|
|
rc3.set( data::Rates::opt::gas, 28.41 );
|
|
|
|
|
|
|
|
data::Well w1, w2;
|
|
|
|
w1.rates = r1;
|
|
|
|
w1.bhp = 1.23;
|
|
|
|
w1.temperature = 3.45;
|
2016-09-29 04:58:10 -05:00
|
|
|
w1.control = 1;
|
2016-09-25 11:31:34 -05:00
|
|
|
|
|
|
|
/*
|
|
|
|
* the completion keys (active indices) and well names correspond to the
|
|
|
|
* input deck. All other entries in the well structures are arbitrary.
|
|
|
|
*/
|
|
|
|
w1.completions[ 88 ] = { 88, rc1, 30.45 };
|
|
|
|
w1.completions[ 288 ] = { 288, rc2, 33.19 };
|
|
|
|
|
|
|
|
w2.rates = r2;
|
|
|
|
w2.bhp = 2.34;
|
|
|
|
w2.temperature = 4.56;
|
2016-09-29 04:58:10 -05:00
|
|
|
w2.control = 2;
|
2016-09-25 11:31:34 -05:00
|
|
|
w2.completions[ 188 ] = { 188, rc3, 36.22 };
|
|
|
|
|
|
|
|
return { { "OP_1", w1 },
|
|
|
|
{ "OP_2", w2 } };
|
2016-05-13 09:25:18 -05:00
|
|
|
}
|
|
|
|
|
|
|
|
data::Solution mkSolution( int numCells ) {
|
2016-10-06 09:31:41 -05:00
|
|
|
|
|
|
|
using measure = UnitSystem::measure;
|
|
|
|
using namespace data;
|
|
|
|
|
|
|
|
data::Solution sol = {
|
|
|
|
{ "PRESSURE", { measure::pressure, std::vector<double>( numCells ), TargetType::RESTART_SOLUTION } },
|
|
|
|
{ "TEMP", { measure::temperature, std::vector<double>( numCells ), TargetType::RESTART_SOLUTION } },
|
|
|
|
{ "SWAT", { measure::identity, std::vector<double>( numCells ), TargetType::RESTART_SOLUTION } },
|
|
|
|
{ "SGAS", { measure::identity, std::vector<double>( numCells ), TargetType::RESTART_SOLUTION } }
|
|
|
|
};
|
2016-05-13 09:25:18 -05:00
|
|
|
|
|
|
|
|
2016-10-06 01:16:33 -05:00
|
|
|
sol.data("PRESSURE").assign( numCells, 6.0 );
|
|
|
|
sol.data("TEMP").assign( numCells, 7.0 );
|
|
|
|
sol.data("SWAT").assign( numCells, 8.0 );
|
|
|
|
sol.data("SGAS").assign( numCells, 9.0 );
|
2016-05-13 09:25:18 -05:00
|
|
|
|
|
|
|
fun::iota rsi( 300, 300 + numCells );
|
|
|
|
fun::iota rvi( 400, 400 + numCells );
|
|
|
|
|
2016-10-06 09:31:41 -05:00
|
|
|
sol.insert( "RS", measure::identity, { rsi.begin(), rsi.end() } , TargetType::RESTART_SOLUTION );
|
|
|
|
sol.insert( "RV", measure::identity, { rvi.begin(), rvi.end() } , TargetType::RESTART_SOLUTION );
|
2016-05-13 09:25:18 -05:00
|
|
|
|
|
|
|
return sol;
|
|
|
|
}
|
|
|
|
|
|
|
|
std::pair< data::Solution, data::Wells >
|
|
|
|
first_sim(test_work_area_type * test_area) {
|
|
|
|
|
|
|
|
std::string eclipse_data_filename = "FIRST_SIM.DATA";
|
|
|
|
test_work_area_copy_file(test_area, eclipse_data_filename.c_str());
|
|
|
|
|
2016-10-13 07:15:58 -05:00
|
|
|
auto eclipseState = Parser::parse( eclipse_data_filename );
|
2016-05-13 09:25:18 -05:00
|
|
|
|
2016-10-13 07:15:58 -05:00
|
|
|
const auto& grid = eclipseState.getInputGrid();
|
2016-05-13 09:25:18 -05:00
|
|
|
auto num_cells = grid.getNX() * grid.getNY() * grid.getNZ();
|
|
|
|
|
2016-09-07 11:08:05 -05:00
|
|
|
EclipseWriter eclWriter( eclipseState, grid);
|
2016-09-19 12:28:06 -05:00
|
|
|
auto start_time = ecl_util_make_date( 1, 11, 1979 );
|
|
|
|
auto first_step = ecl_util_make_date( 10, 10, 2008 );
|
2016-05-13 09:25:18 -05:00
|
|
|
|
|
|
|
auto sol = mkSolution( num_cells );
|
|
|
|
auto wells = mkWells();
|
|
|
|
|
|
|
|
eclWriter.writeTimeStep( 1,
|
2016-08-20 11:50:43 -05:00
|
|
|
false,
|
2016-05-13 09:25:18 -05:00
|
|
|
first_step - start_time,
|
2016-08-20 11:50:43 -05:00
|
|
|
sol, wells);
|
2016-05-13 09:25:18 -05:00
|
|
|
|
|
|
|
return { sol, wells };
|
|
|
|
}
|
|
|
|
|
|
|
|
std::pair< data::Solution, data::Wells > second_sim() {
|
2016-09-25 11:31:34 -05:00
|
|
|
auto eclipseState = Parser::parseData( input() );
|
2016-05-13 09:25:18 -05:00
|
|
|
|
2016-10-13 07:15:58 -05:00
|
|
|
const auto& grid = eclipseState.getInputGrid();
|
2016-05-13 09:25:18 -05:00
|
|
|
auto num_cells = grid.getNX() * grid.getNY() * grid.getNZ();
|
|
|
|
|
|
|
|
return init_from_restart_file( eclipseState, num_cells );
|
|
|
|
}
|
|
|
|
|
|
|
|
void compare( std::pair< data::Solution, data::Wells > fst,
|
|
|
|
std::pair< data::Solution, data::Wells > snd ) {
|
|
|
|
|
2016-10-06 01:16:33 -05:00
|
|
|
for( auto key : { "PRESSURE", "TEMP", "SWAT", "SGAS",
|
|
|
|
"RS", "RV" } ) {
|
2016-05-13 09:25:18 -05:00
|
|
|
|
2016-10-06 01:16:33 -05:00
|
|
|
auto first = fst.first.data( key ).begin();
|
|
|
|
auto second = snd.first.data( key ).begin();
|
2016-05-13 09:25:18 -05:00
|
|
|
|
2016-10-06 01:16:33 -05:00
|
|
|
for( ; first != fst.first.data( key ).end(); ++first, ++second )
|
2016-05-13 09:25:18 -05:00
|
|
|
BOOST_CHECK_CLOSE( *first, *second, 0.00001 );
|
|
|
|
}
|
|
|
|
|
2016-09-25 11:31:34 -05:00
|
|
|
BOOST_CHECK_EQUAL( fst.second, snd.second );
|
2016-05-13 09:25:18 -05:00
|
|
|
}
|
|
|
|
|
|
|
|
BOOST_AUTO_TEST_CASE(EclipseReadWriteWellStateData) {
|
|
|
|
test_work_area_type * test_area = test_work_area_alloc("test_Restart");
|
|
|
|
|
|
|
|
auto state1 = first_sim(test_area);
|
|
|
|
auto state2 = second_sim();
|
|
|
|
compare(state1, state2);
|
|
|
|
|
|
|
|
test_work_area_free(test_area);
|
|
|
|
}
|
2016-09-25 11:31:34 -05:00
|
|
|
|
|
|
|
BOOST_AUTO_TEST_CASE(OPM_XWEL) {
|
|
|
|
auto es = Parser::parseData( input( "XWEL" ) );
|
|
|
|
const auto& sched = es.getSchedule();
|
|
|
|
const auto& grid = es.getInputGrid();
|
|
|
|
const auto& tm = es.getTableManager();
|
|
|
|
|
|
|
|
std::vector< data::Rates::opt > phases {
|
|
|
|
data::Rates::opt::wat,
|
|
|
|
data::Rates::opt::oil,
|
|
|
|
data::Rates::opt::gas,
|
|
|
|
};
|
|
|
|
|
|
|
|
const auto wells = mkWells();
|
|
|
|
const auto& sched_wells = sched.getWells( 1 );
|
2016-09-29 04:58:10 -05:00
|
|
|
const auto xwel = serialize_XWEL( wells, 1, sched_wells, tm, grid );
|
|
|
|
const auto iwel = serialize_IWEL( wells, sched_wells );
|
|
|
|
|
|
|
|
const auto restored_wells = restore_wells( xwel.data(), xwel.size(),
|
|
|
|
iwel.data(), iwel.size(),
|
|
|
|
1,
|
|
|
|
sched.getWells( 1 ),
|
|
|
|
phases,
|
|
|
|
grid );
|
2016-09-25 11:31:34 -05:00
|
|
|
|
|
|
|
BOOST_CHECK_EQUAL( wells, restored_wells );
|
|
|
|
}
|