Serialization: Add Initial Support for Segment Information

This commit extends Opm::data::Wells to include a set of output
vectors for well segment information.  At present we define output
structures for segment rates and segment pressures.  The immediate
use case is properly assigning restart vector items RSEG[8 .. 11],
but these same values are also usable for outputting the summary
vectors SPR, SOFR, SGFR, and SWFR.  Future expansion is likely.
This commit is contained in:
Bård Skaflestad 2018-10-04 15:18:35 +02:00
parent ae0cb3e5ad
commit 80154e8f5f
2 changed files with 92 additions and 4 deletions

View File

@ -20,12 +20,14 @@
#ifndef OPM_OUTPUT_WELLS_HPP
#define OPM_OUTPUT_WELLS_HPP
#include <algorithm>
#include <cstddef>
#include <initializer_list>
#include <map>
#include <algorithm>
#include <stdexcept>
#include <string>
#include <type_traits>
#include <unordered_map>
#include <vector>
namespace Opm {
@ -120,6 +122,18 @@ namespace Opm {
void read(MessageBufferType& buffer);
};
struct Segment {
Rates rates;
double pressure;
std::size_t segIndex;
template <class MessageBufferType>
void write(MessageBufferType& buffer) const;
template <class MessageBufferType>
void read(MessageBufferType& buffer);
};
struct Well {
Rates rates;
double bhp;
@ -127,6 +141,7 @@ namespace Opm {
double temperature;
int control;
std::vector< Connection > connections;
std::unordered_map<std::size_t, Segment> segments;
inline bool flowing() const noexcept;
template <class MessageBufferType>
@ -300,6 +315,13 @@ namespace Opm {
buffer.write(this->effective_Kh);
}
template <class MessageBufferType>
void Segment::write(MessageBufferType& buffer) const {
buffer.write(this->segIndex);
this->rates.write(buffer);
buffer.write(this->pressure);
}
template <class MessageBufferType>
void Well::write(MessageBufferType& buffer) const {
this->rates.write(buffer);
@ -311,6 +333,16 @@ namespace Opm {
buffer.write(size);
for (const Connection& comp : this->connections)
comp.write(buffer);
{
const auto nSeg =
static_cast<unsigned int>(this->segments.size());
buffer.write(nSeg);
for (const auto& seg : this->segments) {
seg.second.write(buffer);
}
}
}
template <class MessageBufferType>
@ -341,6 +373,13 @@ namespace Opm {
buffer.read(this->effective_Kh);
}
template <class MessageBufferType>
void Segment::read(MessageBufferType& buffer) {
buffer.read(this->segIndex);
this->rates.read(buffer);
buffer.read(this->pressure);
}
template <class MessageBufferType>
void Well::read(MessageBufferType& buffer) {
this->rates.read(buffer);
@ -348,6 +387,8 @@ namespace Opm {
buffer.read(this->thp);
buffer.read(this->temperature);
buffer.read(this->control);
// Connection information
unsigned int size = 0.0; //this->connections.size();
buffer.read(size);
this->connections.resize(size);
@ -356,9 +397,27 @@ namespace Opm {
auto& comp = this->connections[ i ];
comp.read(buffer);
}
// Segment information (if applicable)
const auto nSeg = [&buffer]() -> unsigned int
{
auto n = 0u;
buffer.read(n);
return n;
}();
for (auto segID = 0*nSeg; segID < nSeg; ++segID) {
auto seg = Segment{};
seg.read(buffer);
const auto segIndex = seg.segIndex;
this->segments.emplace(segIndex, std::move(seg));
}
}
}
}
}} // Opm::data
#endif //OPM_OUTPUT_WELLS_HPP

View File

@ -23,6 +23,7 @@
#include <boost/test/unit_test.hpp>
#include <boost/date_time/posix_time/posix_time.hpp>
#include <exception>
#include <stdexcept>
#include <ert/ecl/ecl_sum.h>
@ -43,6 +44,8 @@
#include <opm/parser/eclipse/Parser/ParseContext.hpp>
#include <opm/parser/eclipse/Parser/Parser.hpp>
#include <opm/parser/eclipse/Units/Units.hpp>
using namespace Opm;
using rt = data::Rates::opt;
@ -141,6 +144,15 @@ static data::Wells result_wells() {
crates3.set( rt::reservoir_oil, 300.7 / day );
crates3.set( rt::reservoir_gas, 300.8 / day );
// Segment vectors
const auto sm3_pr_day = unit::cubic(unit::meter) / day;
auto segment = ::Opm::data::Segment{};
segment.rates.set(rt::wat, 123.45*sm3_pr_day);
segment.rates.set(rt::oil, 543.21*sm3_pr_day);
segment.rates.set(rt::gas, 1729.496*sm3_pr_day);
segment.pressure = 314.159*unit::barsa;
segment.segIndex = 1;
/*
The global index assigned to the completion must be manually
syncronized with the global index in the COMPDAT keyword in the
@ -151,10 +163,15 @@ static data::Wells result_wells() {
data::Connection well2_comp2 { 101, crates3, 1.11 , 123.4, 150.6, 0.001, 0.89, 100.0};
data::Connection well3_comp1 { 2 , crates3, 1.11 , 123.4, 456.78, 0.0, 0.15, 432.1};
/*
The completions
*/
data::Well well1 { rates1, 0.1 * ps, 0.2 * ps, 0.3 * ps, 1, { {well1_comp1} } };
data::Well well1 {
rates1, 0.1 * ps, 0.2 * ps, 0.3 * ps, 1,
{ {well1_comp1} },
{ { segment.segIndex, segment } },
};
data::Well well2 { rates2, 1.1 * ps, 1.2 * ps, 1.3 * ps, 2, { {well2_comp1 , well2_comp2} } };
data::Well well3 { rates3, 2.1 * ps, 2.2 * ps, 2.3 * ps, 3, { {well3_comp1} } };
@ -1178,6 +1195,18 @@ BOOST_AUTO_TEST_CASE(READ_WRITE_WELLDATA) {
BOOST_CHECK_CLOSE( wellRatesCopy.get( "W_1" , rt::wat) , wellRates.get( "W_1" , rt::wat), 1e-16);
BOOST_CHECK_CLOSE( wellRatesCopy.get( "W_2" , 101 , rt::wat) , wellRates.get( "W_2" , 101 , rt::wat), 1e-16);
const auto sm3_pr_day = unit::cubic(unit::meter) / day;
const auto& seg = wellRatesCopy.at("W_1").segments.at(1);
BOOST_CHECK_CLOSE(seg.rates.get(rt::wat), 123.45*sm3_pr_day, 1.0e-10);
BOOST_CHECK_CLOSE(seg.rates.get(rt::oil), 543.21*sm3_pr_day, 1.0e-10);
BOOST_CHECK_CLOSE(seg.rates.get(rt::gas), 1729.496*sm3_pr_day, 1.0e-10);
BOOST_CHECK_CLOSE(seg.pressure, 314.159*unit::barsa, 1.0e-10);
BOOST_CHECK_EQUAL(seg.segIndex, 1);
// No data for segment 10 of well W_2 (or no such segment).
const auto& W2 = wellRatesCopy.at("W_2");
BOOST_CHECK_THROW(W2.segments.at(10), std::out_of_range);
}
BOOST_AUTO_TEST_CASE(efficiency_factor) {