Add method WellSegments::segmentDepthChange()
This commit is contained in:
@@ -437,6 +437,7 @@ if(ENABLE_ECL_OUTPUT)
|
||||
tests/SPE1CASE2.X0060
|
||||
tests/PYACTION.DATA
|
||||
tests/act1.py
|
||||
tests/MSW.DATA
|
||||
tests/EXIT_TEST.DATA
|
||||
tests/action_syntax_error.py
|
||||
tests/action_missing_run.py
|
||||
|
||||
@@ -89,6 +89,7 @@ namespace Opm {
|
||||
bool operator!=( const WellSegments& ) const;
|
||||
|
||||
double segmentLength(const int segment_number) const;
|
||||
double segmentDepthChange(const int segment_number) const;
|
||||
|
||||
// it returns true if there is no error encountered during the update
|
||||
bool updateWSEGSICD(const std::vector<std::pair<int, SpiralICD> >& sicd_pairs);
|
||||
|
||||
@@ -442,19 +442,13 @@ namespace Opm {
|
||||
}
|
||||
|
||||
double WellSegments::segmentLength(const int segment_number) const {
|
||||
double segment_length;
|
||||
|
||||
const Segment& segment = getFromSegmentNumber(segment_number);
|
||||
if (segment_number == 1) {// top segment
|
||||
segment_length = segment.totalLength();
|
||||
} else {
|
||||
// other segments
|
||||
const int outlet_segment_number = segment.outletSegment();
|
||||
const Segment &outlet_segment = getFromSegmentNumber(outlet_segment_number);
|
||||
|
||||
segment_length = segment.totalLength() - outlet_segment.totalLength();
|
||||
}
|
||||
const Segment& segment = this->getFromSegmentNumber(segment_number);
|
||||
if (segment_number == 1) // top segment
|
||||
return segment.totalLength();
|
||||
|
||||
// other segments
|
||||
const Segment &outlet_segment = getFromSegmentNumber(segment.outletSegment());
|
||||
const double segment_length = segment.totalLength() - outlet_segment.totalLength();
|
||||
if (segment_length <= 0.)
|
||||
throw std::runtime_error(" non positive segment length is obtained for segment "
|
||||
+ std::to_string(segment_number));
|
||||
@@ -462,6 +456,18 @@ namespace Opm {
|
||||
return segment_length;
|
||||
}
|
||||
|
||||
|
||||
double WellSegments::segmentDepthChange(const int segment_number) const {
|
||||
const Segment& segment = getFromSegmentNumber(segment_number);
|
||||
if (segment_number == 1) // top segment
|
||||
return segment.depth();
|
||||
|
||||
// other segments
|
||||
const Segment &outlet_segment = this->getFromSegmentNumber(segment.outletSegment());
|
||||
return segment.depth() - outlet_segment.depth();
|
||||
}
|
||||
|
||||
|
||||
bool WellSegments::updateWSEGSICD(const std::vector<std::pair<int, SpiralICD> >& sicd_pairs) {
|
||||
if (m_comp_pressure_drop == CompPressureDrop::H__) {
|
||||
const std::string msg = "to use spiral ICD segment you have to activate the frictional pressure drop calculation";
|
||||
|
||||
@@ -33,6 +33,7 @@
|
||||
#include <opm/parser/eclipse/Deck/DeckRecord.hpp>
|
||||
#include <opm/parser/eclipse/Deck/DeckKeyword.hpp>
|
||||
|
||||
#include <opm/parser/eclipse/EclipseState/Schedule/Schedule.hpp>
|
||||
#include <opm/parser/eclipse/EclipseState/Schedule/MSW/SpiralICD.hpp>
|
||||
#include <opm/parser/eclipse/EclipseState/Schedule/MSW/Valve.hpp>
|
||||
#include <opm/parser/eclipse/EclipseState/Schedule/Well/Connection.hpp>
|
||||
@@ -439,3 +440,24 @@ BOOST_AUTO_TEST_CASE(testwsegvalv) {
|
||||
BOOST_CHECK_EQUAL(segment2.roughness(), valv2->pipeRoughness());
|
||||
BOOST_CHECK_EQUAL(segment2.crossArea(), valv2->pipeCrossArea());
|
||||
}
|
||||
|
||||
|
||||
BOOST_AUTO_TEST_CASE(MSW_SEGMENT_LENGTH) {
|
||||
Opm::Parser parser;
|
||||
Opm::Deck deck = parser.parseFile("MSW.DATA");
|
||||
Opm::EclipseState st(deck);
|
||||
Opm::Schedule sched(deck, st);
|
||||
|
||||
|
||||
const auto& well = sched.getWell("PROD01", 0);
|
||||
const auto& segments = well.getSegments();
|
||||
BOOST_CHECK_CLOSE( segments.segmentLength(1), 2512.50, 1e-5);
|
||||
BOOST_CHECK_CLOSE( segments.segmentLength(2), 25, 1e-5);
|
||||
BOOST_CHECK_CLOSE( segments.segmentLength(6), 25, 1e-5);
|
||||
BOOST_CHECK_CLOSE( segments.segmentLength(7), 200, 1e-5);
|
||||
|
||||
BOOST_CHECK_CLOSE( segments.segmentDepthChange(1), 2512.50, 1e-5);
|
||||
BOOST_CHECK_CLOSE( segments.segmentDepthChange(2), 22, 1e-5);
|
||||
BOOST_CHECK_CLOSE( segments.segmentDepthChange(6), 21, 1e-5);
|
||||
BOOST_CHECK_CLOSE( segments.segmentDepthChange(7), 4, 1e-5);
|
||||
}
|
||||
|
||||
Reference in New Issue
Block a user