parser/eclipse/EclipseState/Schedule/MSW/segment.hpp
Go to the documentation of this file.
1/*
2 Copyright 2015 SINTEF ICT, Applied Mathematics.
3
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 3 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18*/
19
20#ifndef SEGMENT_HPP_HEADER_INCLUDED
21#define SEGMENT_HPP_HEADER_INCLUDED
22
23#include <memory>
24#include <vector>
25
26namespace Opm {
27 class SpiralICD;
28 class Valve;
29
30 namespace RestartIO {
31 struct RstSegment;
32 }
33}
34
35namespace Opm {
36
37 class Segment {
38 public:
39
40 enum class SegmentType {
41 REGULAR,
42 SICD,
43 AICD, // Not really supported - just included to complete the enum
44 VALVE,
45 };
46
48
49 Segment(const Segment& src, double new_depth, double new_length, double new_volume);
50 Segment(const Segment& src, double new_depth, double new_length);
51 Segment(const Segment& src, double new_volume);
52 Segment(int segment_number_in, int branch_in, int outlet_segment_in, double length_in, double depth_in,
53 double internal_diameter_in, double roughness_in, double cross_area_in, double volume_in, bool data_ready_in, SegmentType segment_type_in);
54 Segment(const RestartIO::RstSegment& rst_segment);
55
57
58 int segmentNumber() const;
59 int branchNumber() const;
60 int outletSegment() const;
61 double totalLength() const;
62 double depth() const;
63 double internalDiameter() const;
64 double roughness() const;
65 double crossArea() const;
66 double volume() const;
67 bool dataReady() const;
68
70 int ecl_type_id() const;
71
72
73 const std::vector<int>& inletSegments() const;
74
75 static double invalidValue();
76 static SegmentType type_from_int(int ecl_id);
77
78 bool operator==( const Segment& ) const;
79 bool operator!=( const Segment& ) const;
80
81 const std::shared_ptr<SpiralICD>& spiralICD() const;
82 const Valve* valve() const;
83
84 void updateSpiralICD(const SpiralICD& spiral_icd);
85 void updateValve(const Valve& valve, const double segment_length);
86 void addInletSegment(const int segment_number);
87
88 template<class Serializer>
89 void serializeOp(Serializer& serializer)
90 {
91 serializer(m_segment_number);
92 serializer(m_branch);
93 serializer(m_outlet_segment);
94 serializer(m_inlet_segments);
95 serializer(m_total_length);
96 serializer(m_depth);
97 serializer(m_internal_diameter);
98 serializer(m_roughness);
99 serializer(m_cross_area);
100 serializer(m_volume);
101 serializer(m_data_ready);
102 serializer(m_segment_type);
103 serializer(m_spiral_icd);
104 serializer(m_valve);
105 }
106
107 private:
108 // segment number
109 // it should work as a ID.
110 int m_segment_number;
111 // branch number
112 // for top segment, it should always be 1
113 int m_branch;
114 // the outlet junction segment
115 // for top segment, it should be -1
116 int m_outlet_segment;
117 // the segments whose outlet segments are the current segment
118 std::vector<int> m_inlet_segments;
119 // length of the segment node to the bhp reference point.
120 // when reading in from deck, with 'INC',
121 // it will be incremental length before processing.
122 // After processing and in the class Well, it always stores the 'ABS' value.
123 // which means the total_length
124 double m_total_length;
125 // depth of the nodes to the bhp reference point
126 // when reading in from deck, with 'INC',
127 // it will be the incremental depth before processing.
128 // in the class Well, it always stores the 'ABS' value.
129 // TODO: to check if it is good to use 'ABS' always.
130 double m_depth;
131 // tubing internal diameter
132 // or the equivalent diameter for annular cross-sections
133 // for top segment, it is UNDEFINED
134 // we use invalid_value for the top segment
135 double m_internal_diameter;
136 // effective roughness of the tubing
137 // used to calculate the Fanning friction factor
138 // for top segment, it is UNDEFINED
139 // we use invalid_value for the top segment
140 double m_roughness;
141 // cross-sectional area for fluid flow
142 // not defined for the top segment,
143 // we use invalid_value for the top segment.
144 double m_cross_area;
145 // valume of the segment;
146 // it is defined for top segment.
147 // TODO: to check if the definition is the same with other segments.
148 double m_volume;
149 // indicate if the data related to 'INC' or 'ABS' is ready
150 // the volume will be updated at a final step.
151 bool m_data_ready;
152
153 // indicate the type of the segment
154 // regular, spiral ICD, or Valve.
155 SegmentType m_segment_type;
156
157 // information related to SpiralICD. It is nullptr for segments are not
158 // spiral ICD type
159 std::shared_ptr<SpiralICD> m_spiral_icd;
160
161 // information related to sub-critical valve. It is nullptr for segments are not
162 // of type of Valve
163 std::shared_ptr<Valve> m_valve;
164
165 // We are not handling the length of segment projected onto the X-axis and Y-axis.
166 // They are not used in the simulations and we are not supporting the plotting.
167 // There are other three properties for segment related to thermal conduction,
168 // while they are not supported by the keyword at the moment.
169 };
170
171 inline bool isRegular(const Segment& segment)
172 {
173 return segment.segmentType() == Segment::SegmentType::REGULAR;
174 }
175
176 inline bool isSpiralICD(const Segment& segment)
177 {
178 return segment.segmentType() == Segment::SegmentType::SICD;
179 }
180
181 inline bool isValve(const Segment& segment)
182 {
183 return segment.segmentType() == Segment::SegmentType::VALVE;
184 }
185}
186
187#endif
Definition: parser/eclipse/EclipseState/Schedule/MSW/segment.hpp:37
Segment(const Segment &src, double new_depth, double new_length, double new_volume)
double crossArea() const
Segment(const Segment &src, double new_depth, double new_length)
double totalLength() const
const std::shared_ptr< SpiralICD > & spiralICD() const
Segment(const RestartIO::RstSegment &rst_segment)
double volume() const
int segmentNumber() const
static Segment serializeObject()
void updateSpiralICD(const SpiralICD &spiral_icd)
void updateValve(const Valve &valve, const double segment_length)
SegmentType
Definition: parser/eclipse/EclipseState/Schedule/MSW/segment.hpp:40
int branchNumber() const
double depth() const
Segment(int segment_number_in, int branch_in, int outlet_segment_in, double length_in, double depth_in, double internal_diameter_in, double roughness_in, double cross_area_in, double volume_in, bool data_ready_in, SegmentType segment_type_in)
Segment(const Segment &src, double new_volume)
int outletSegment() const
bool operator==(const Segment &) const
void addInletSegment(const int segment_number)
bool operator!=(const Segment &) const
const std::vector< int > & inletSegments() const
double internalDiameter() const
int ecl_type_id() const
const Valve * valve() const
double roughness() const
void serializeOp(Serializer &serializer)
Definition: parser/eclipse/EclipseState/Schedule/MSW/segment.hpp:89
static double invalidValue()
static SegmentType type_from_int(int ecl_id)
bool dataReady() const
SegmentType segmentType() const
Definition: Serializer.hpp:38
Definition: SpiralICD.hpp:36
Definition: Valve.hpp:37
Definition: A.hpp:4
bool isValve(const Segment &segment)
Definition: parser/eclipse/EclipseState/Schedule/MSW/segment.hpp:181
bool isRegular(const Segment &segment)
Definition: parser/eclipse/EclipseState/Schedule/MSW/segment.hpp:171
bool isSpiralICD(const Segment &segment)
Definition: parser/eclipse/EclipseState/Schedule/MSW/segment.hpp:176
Definition: io/eclipse/rst/segment.hpp:32