Add implementation for autonomous ICD: AutoICD

This commit is contained in:
Joakim Hove
2020-06-11 09:32:54 +02:00
parent 73a244ceda
commit 224911bd10
6 changed files with 328 additions and 25 deletions

View File

@@ -46,6 +46,148 @@
#include <opm/parser/eclipse/Parser/ParseContext.hpp>
BOOST_AUTO_TEST_CASE(AICDWellTest) {
auto dir = Opm::Connection::Direction::Z;
const auto kind = Opm::Connection::CTFKind::DeckValue;
Opm::WellConnections connection_set(Opm::Connection::Order::TRACK, 10,10);
Opm::EclipseGrid grid(20,20,20);
connection_set.add(Opm::Connection( 19, 0, 0,grid.getGlobalIndex(19,0,0), 1, 0.0, Opm::Connection::State::OPEN , 200, 17.29, 0.25, 0.0, 0.0, 0, dir, kind, 0, true) );
connection_set.add(Opm::Connection( 19, 0, 1,grid.getGlobalIndex(19,0,1), 1, 0.0, Opm::Connection::State::OPEN , 200, 17.29, 0.25, 0.0, 0.0, 0, dir, kind, 0, true) );
connection_set.add(Opm::Connection( 19, 0, 2,grid.getGlobalIndex(19,0,2), 1, 0.0, Opm::Connection::State::OPEN , 200, 17.29, 0.25, 0.0, 0.0, 0, dir, kind, 0, true) );
connection_set.add(Opm::Connection( 18, 0, 1,grid.getGlobalIndex(18,0,1), 1, 0.0, Opm::Connection::State::OPEN , 200, 17.29, 0.25, 0.0, 0.0, 0, Opm::Connection::Direction::X, kind, 0, true) );
connection_set.add(Opm::Connection( 17, 0, 1,grid.getGlobalIndex(17,0,1), 1, 0.0, Opm::Connection::State::OPEN , 200, 17.29, 0.25, 0.0, 0.0, 0, Opm::Connection::Direction::X, kind, 0, true) );
connection_set.add(Opm::Connection( 16, 0, 1,grid.getGlobalIndex(16,0,1), 1, 0.0, Opm::Connection::State::OPEN , 200, 17.29, 0.25, 0.0, 0.0, 0, Opm::Connection::Direction::X, kind, 0, true) );
connection_set.add(Opm::Connection( 15, 0, 1,grid.getGlobalIndex(15,0,1), 1, 0.0, Opm::Connection::State::OPEN , 200, 17.29, 0.25, 0.0, 0.0, 0, Opm::Connection::Direction::X, kind, 0, true) );
BOOST_CHECK_EQUAL( 7U , connection_set.size() );
const std::string compsegs_string = R"(
WELSEGS
'PROD01' 2512.5 2512.5 1.0e-5 'ABS' 'HF-' 'HO' /
2 2 1 1 2537.5 2537.5 0.3 0.00010 /
3 3 1 2 2562.5 2562.5 0.2 0.00010 /
4 4 2 2 2737.5 2537.5 0.2 0.00010 /
6 6 2 4 3037.5 2539.5 0.2 0.00010 /
7 7 2 6 3337.5 2534.5 0.2 0.00010 /
8 8 3 7 3337.6 2534.5 0.2 0.00015 /
/
COMPSEGS
PROD01 /
20 1 1 1 2512.5 2525.0 /
20 1 2 1 2525.0 2550.0 /
20 1 3 1 2550.0 2575.0 /
19 1 2 2 2637.5 2837.5 /
18 1 2 2 2837.5 3037.5 /
17 1 2 2 3037.5 3237.5 /
16 1 2 3 3237.5 3437.5 /
/
WSEGAICD
'PROD01' 8 8 0.002 -0.7 1* 1* 0.6 1* 1* 2* 1.0 1.0 'SHUT' /
/
)";
Opm::Parser parser;
Opm::Deck deck = parser.parseString(compsegs_string);
const Opm::DeckKeyword compsegs = deck.getKeyword("COMPSEGS");
BOOST_CHECK_EQUAL( 8U, compsegs.size() );
const Opm::DeckKeyword welsegs = deck.getKeyword("WELSEGS");
Opm::WellSegments segment_set(welsegs);
BOOST_CHECK_EQUAL(7U, segment_set.size());
Opm::ErrorGuard errorGuard;
Opm::ParseContext parseContext;
parseContext.update(Opm::ParseContext::SCHEDULE_COMPSEGS_INVALID, Opm::InputError::THROW_EXCEPTION);
parseContext.update(Opm::ParseContext::SCHEDULE_COMPSEGS_NOT_SUPPORTED, Opm::InputError::THROW_EXCEPTION);
const auto& [new_connection_set, new_segment_set] = Opm::Compsegs::processCOMPSEGS(compsegs, connection_set, segment_set, grid, parseContext, errorGuard);
// checking the ICD segment
const Opm::DeckKeyword wsegaicd = deck.getKeyword("WSEGAICD");
const auto aicd_map = Opm::AutoICD::fromWSEGAICD(wsegaicd);
BOOST_CHECK_EQUAL(1U, aicd_map.size());
const auto it = aicd_map.begin();
const std::string& well_name = it->first;
BOOST_CHECK_EQUAL(well_name, "PROD01");
const auto& aicd_vector = it->second;
BOOST_CHECK_EQUAL(1U, aicd_vector.size());
const int segment_number = aicd_vector[0].first;
const Opm::AutoICD& aicd0 = aicd_vector[0].second;
BOOST_CHECK_EQUAL(8, segment_number);
Opm::Segment segment = segment_set.getFromSegmentNumber(segment_number);
segment.updateAutoICD(aicd0);
BOOST_CHECK(Opm::Segment::SegmentType::AICD == segment.segmentType());
auto aicd = segment.autoICD();
BOOST_CHECK_GT(aicd.maxAbsoluteRate(), 1.e99);
BOOST_CHECK(aicd.status()==Opm::ICDStatus::SHUT);
// 0.002 bars*day*day/Volume^2
BOOST_CHECK_EQUAL(aicd.strength(), 0.002*1.e5*86400.*86400.);
BOOST_CHECK_EQUAL(aicd.length(), -0.7);
BOOST_CHECK_EQUAL(aicd.densityCalibration(), 1000.25);
// 1.45 cp
BOOST_CHECK_EQUAL(aicd.viscosityCalibration(), 1.45 * 0.001);
BOOST_CHECK_EQUAL(aicd.criticalValue(), 0.6);
BOOST_CHECK_EQUAL(aicd.widthTransitionRegion(), 0.05);
BOOST_CHECK_EQUAL(aicd.maxViscosityRatio(), 5.0);
BOOST_CHECK_EQUAL(aicd.methodFlowScaling(), -1);
// the scaling factor has not been updated properly, so it will throw
BOOST_CHECK_THROW(aicd.scalingFactor(), std::runtime_error);
const int outlet_segment_number = segment.outletSegment();
const double outlet_segment_length = segment_set.segmentLength(outlet_segment_number);
// only one connection attached to the outlet segment in this case
const Opm::Connection& connection = new_connection_set.getFromIJK(15, 0, 1);
const auto& perf_range = connection.perf_range();
const auto connection_length = perf_range->second - perf_range->first;
aicd.updateScalingFactor(outlet_segment_length, connection_length);
// updated, so it should not throw
BOOST_CHECK_NO_THROW(aicd.scalingFactor());
BOOST_CHECK_EQUAL(0.7, aicd.scalingFactor());
BOOST_CHECK_EQUAL(7U, new_connection_set.size());
const Opm::Connection& connection1 = new_connection_set.get(0);
const int segment_number_connection1 = connection1.segment();
const double center_depth_connection1 = connection1.depth();
BOOST_CHECK_EQUAL(segment_number_connection1, 1);
BOOST_CHECK_EQUAL(center_depth_connection1, 2512.5);
const Opm::Connection& connection3 = new_connection_set.get(2);
const int segment_number_connection3 = connection3.segment();
const double center_depth_connection3 = connection3.depth();
BOOST_CHECK_EQUAL(segment_number_connection3, 3);
BOOST_CHECK_EQUAL(center_depth_connection3, 2562.5);
const Opm::Connection& connection5 = new_connection_set.get(4);
const int segment_number_connection5 = connection5.segment();
const double center_depth_connection5 = connection5.depth();
BOOST_CHECK_EQUAL(segment_number_connection5, 6);
BOOST_CHECK_CLOSE(center_depth_connection5, 2538.83, 0.001);
const Opm::Connection& connection6 = new_connection_set.get(5);
const int segment_number_connection6 = connection6.segment();
const double center_depth_connection6 = connection6.depth();
BOOST_CHECK_EQUAL(segment_number_connection6, 6);
BOOST_CHECK_CLOSE(center_depth_connection6, 2537.83, 0.001);
const Opm::Connection& connection7 = new_connection_set.get(6);
const int segment_number_connection7 = connection7.segment();
const double center_depth_connection7 = connection7.depth();
BOOST_CHECK_EQUAL(segment_number_connection7, 8);
BOOST_CHECK_EQUAL(center_depth_connection7, 2534.5);
}
BOOST_AUTO_TEST_CASE(MultisegmentWellTest) {
@@ -64,30 +206,31 @@ BOOST_AUTO_TEST_CASE(MultisegmentWellTest) {
BOOST_CHECK_EQUAL( 7U , connection_set.size() );
const std::string compsegs_string =
"WELSEGS \n"
"'PROD01' 2512.5 2512.5 1.0e-5 'ABS' 'HF-' 'HO' /\n"
"2 2 1 1 2537.5 2537.5 0.3 0.00010 /\n"
"3 3 1 2 2562.5 2562.5 0.2 0.00010 /\n"
"4 4 2 2 2737.5 2537.5 0.2 0.00010 /\n"
"6 6 2 4 3037.5 2539.5 0.2 0.00010 /\n"
"7 7 2 6 3337.5 2534.5 0.2 0.00010 /\n"
"8 8 3 7 3337.6 2534.5 0.2 0.00015 /\n"
"/\n"
"\n"
"COMPSEGS\n"
"PROD01 / \n"
"20 1 1 1 2512.5 2525.0 /\n"
"20 1 2 1 2525.0 2550.0 /\n"
"20 1 3 1 2550.0 2575.0 /\n"
"19 1 2 2 2637.5 2837.5 /\n"
"18 1 2 2 2837.5 3037.5 /\n"
"17 1 2 2 3037.5 3237.5 /\n"
"16 1 2 3 3237.5 3437.5 /\n"
"/\n"
"WSEGSICD\n"
"'PROD01' 8 8 0.002 -0.7 1* 1* 0.6 1* 1* 2* 'SHUT' /\n"
"/\n";
const std::string compsegs_string = R"(
WELSEGS
'PROD01' 2512.5 2512.5 1.0e-5 'ABS' 'HF-' 'HO' /
2 2 1 1 2537.5 2537.5 0.3 0.00010 /
3 3 1 2 2562.5 2562.5 0.2 0.00010 /
4 4 2 2 2737.5 2537.5 0.2 0.00010 /
6 6 2 4 3037.5 2539.5 0.2 0.00010 /
7 7 2 6 3337.5 2534.5 0.2 0.00010 /
8 8 3 7 3337.6 2534.5 0.2 0.00015 /
/
COMPSEGS
PROD01 /
20 1 1 1 2512.5 2525.0 /
20 1 2 1 2525.0 2550.0 /
20 1 3 1 2550.0 2575.0 /
19 1 2 2 2637.5 2837.5 /
18 1 2 2 2837.5 3037.5 /
17 1 2 2 3037.5 3237.5 /
16 1 2 3 3237.5 3437.5 /
/
WSEGSICD
'PROD01' 8 8 0.002 -0.7 1* 1* 0.6 1* 1* 2* 'SHUT' /
/
)";
Opm::Parser parser;
Opm::Deck deck = parser.parseString(compsegs_string);