My Project
Loading...
Searching...
No Matches
Segment.hpp
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 <opm/input/eclipse/Schedule/MSW/AICD.hpp>
24#include <opm/input/eclipse/Schedule/MSW/SICD.hpp>
25#include <opm/input/eclipse/Schedule/MSW/Valve.hpp>
26
27#include <optional>
28#include <variant>
29#include <vector>
30
31namespace Opm { namespace RestartIO {
32 struct RstSegment;
33}} // namespace Opm::RestartIO
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
47 Segment();
48
49 Segment(const Segment& src, double new_depth, double new_length, double new_volume, double new_x, double new_y);
50 Segment(const Segment& src, double new_depth, double new_length, double new_x, double new_y);
51 Segment(const Segment& src, double new_depth, double new_length, double new_volume);
52 Segment(const Segment& src, double new_depth, double new_length);
53 Segment(const Segment& src, double new_volume);
54 Segment(const int segment_number_in,
55 const int branch_in,
56 const int outlet_segment_in,
57 const double length_in,
58 const double depth_in,
59 const double internal_diameter_in,
60 const double roughness_in,
61 const double cross_area_in,
62 const double volume_in,
63 const bool data_ready_in,
64 const double x_in,
65 const double y_in);
66
67 explicit Segment(const RestartIO::RstSegment& rst_segment);
68
69 static Segment serializationTestObject();
70
71 int segmentNumber() const;
72 int branchNumber() const;
73 int outletSegment() const;
74 double perfLength() const;
75 double totalLength() const;
76 double node_X() const;
77 double node_Y() const;
78 double depth() const;
79 double internalDiameter() const;
80 double roughness() const;
81 double crossArea() const;
82 double volume() const;
83 bool dataReady() const;
84
85 SegmentType segmentType() const;
86 int ecl_type_id() const;
87
88 const std::vector<int>& inletSegments() const;
89
90 static double invalidValue();
91
92 bool operator==( const Segment& ) const;
93 bool operator!=( const Segment& ) const;
94
95 const SICD& spiralICD() const;
96 const AutoICD& autoICD() const;
97 const Valve& valve() const;
98
99 void updatePerfLength(double perf_length);
100 void updateSpiralICD(const SICD& spiral_icd);
101 void updateAutoICD(const AutoICD& aicd);
102 void updateValve(const Valve& valve, const double segment_length);
103 void updateValve(const Valve& valve);
104 void addInletSegment(const int segment_number);
105
106 bool isRegular() const
107 {
108 return std::holds_alternative<RegularSegment>(this->m_icd);
109 }
110
111 inline bool isSpiralICD() const
112 {
113 return std::holds_alternative<SICD>(this->m_icd);
114 }
115
116 inline bool isAICD() const
117 {
118 return std::holds_alternative<AutoICD>(this->m_icd);
119 }
120
121 inline bool isValve() const
122 {
123 return std::holds_alternative<Valve>(this->m_icd);
124 }
125
126 template<class Serializer>
127 void serializeOp(Serializer& serializer)
128 {
129 serializer(m_segment_number);
130 serializer(m_branch);
131 serializer(m_outlet_segment);
132 serializer(m_inlet_segments);
133 serializer(m_total_length);
134 serializer(m_depth);
135 serializer(m_internal_diameter);
136 serializer(m_roughness);
137 serializer(m_cross_area);
138 serializer(m_volume);
139 serializer(m_data_ready);
140 serializer(m_x);
141 serializer(m_y);
142 serializer(m_perf_length);
143 serializer(m_icd);
144 }
145
146 private:
147 // The current serialization of std::variant<> requires that all the
148 // types in the variant have a serializeOp() method. We introduce
149 // this RegularSegment to meet this requirement. Ideally, the ICD
150 // variant<> would use std::monostate to represent non-ICD segments.
151 struct RegularSegment
152 {
153 template <class Serializer>
154 void serializeOp(Serializer&) {}
155
156 static RegularSegment serializationTestObject() { return {}; }
157
158 bool operator==(const RegularSegment&) const { return true; }
159 };
160
161 // segment number
162 // it should work as a ID.
163 int m_segment_number;
164
165 // branch number
166 // for top segment, it should always be 1
167 int m_branch;
168
169 // the outlet junction segment
170 // for top segment, it should be -1
171 int m_outlet_segment;
172
173 // the segments whose outlet segments are the current segment
174 std::vector<int> m_inlet_segments;
175
176 // length of the segment node to the bhp reference point.
177 // when reading in from deck, with 'INC',
178 // it will be incremental length before processing.
179 // After processing and in the class Well, it always stores the 'ABS' value.
180 // which means the total_length
181 double m_total_length;
182
183 // depth of the nodes to the bhp reference point
184 // when reading in from deck, with 'INC',
185 // it will be the incremental depth before processing.
186 // in the class Well, it always stores the 'ABS' value.
187 // TODO: to check if it is good to use 'ABS' always.
188 double m_depth;
189
190 // tubing internal diameter
191 // or the equivalent diameter for annular cross-sections
192 // for top segment, it is UNDEFINED
193 // we use invalid_value for the top segment
194 double m_internal_diameter;
195
196 // effective roughness of the tubing
197 // used to calculate the Fanning friction factor
198 // for top segment, it is UNDEFINED
199 // we use invalid_value for the top segment
200 double m_roughness;
201
202 // cross-sectional area for fluid flow
203 // not defined for the top segment,
204 // we use invalid_value for the top segment.
205 double m_cross_area;
206
207 // valume of the segment;
208 // it is defined for top segment.
209 // TODO: to check if the definition is the same with other segments.
210 double m_volume;
211
212 // indicate if the data related to 'INC' or 'ABS' is ready
213 // the volume will be updated at a final step.
214 bool m_data_ready;
215
216 // Length of segment projected onto the X axis. Not used in
217 // simulations, but needed for the SEG option in WRFTPLT.
218 double m_x{};
219
220 // Length of segment projected onto the Y axis. Not used in
221 // simulations, but needed for the SEG option in WRFTPLT.
222 double m_y{};
223
224 std::optional<double> m_perf_length;
225 std::variant<RegularSegment, SICD, AutoICD, Valve> m_icd;
226
227 // There are three other properties for the segment pertaining to
228 // thermal conduction. These are not currently supported.
229
230 void updateValve__(Valve& valve, const double segment_length);
231 };
232
233}
234
235#endif
Definition AICD.hpp:44
Definition SICD.hpp:46
Definition Segment.hpp:37
Class for (de-)serializing.
Definition Serializer.hpp:91
Definition Valve.hpp:63
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30
Definition segment.hpp:34