parser/eclipse/EclipseState/Schedule/Well/connection.hpp
Go to the documentation of this file.
1/*
2 Copyright 2013 Statoil ASA.
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
21#ifndef COMPLETION_HPP_
22#define COMPLETION_HPP_
23
24#include <array>
25#include <cstddef>
26#include <map>
27#include <memory>
28#include <string>
29#include <vector>
30
33
34namespace Opm {
35
36namespace RestartIO {
37 struct RstConnection;
38}
39
40 class DeckKeyword;
41 class DeckRecord;
42 class EclipseGrid;
43 class FieldPropsManager;
44
45 class Connection {
46 public:
47 enum class State {
48 OPEN = 1,
49 SHUT = 2,
50 AUTO = 3 // Seems like the AUTO state can not be serialized to restart files.
51 };
52
53 static const std::string State2String( State enumValue );
54 static State StateFromString( const std::string& stringValue );
55
56
57 enum class Direction{
58 X = 1,
59 Y = 2,
60 Z = 3
61 };
62
63 static std::string Direction2String(const Direction enumValue);
64 static Direction DirectionFromString(const std::string& stringValue);
65
66
67 enum class Order {
68 DEPTH,
69 INPUT,
70 TRACK
71 };
72
73 static const std::string Order2String( Order enumValue );
74 static Order OrderFromString(const std::string& comporderStringValue);
75
76 enum class CTFKind {
78 Defaulted,
79 };
80
81
83 Connection(int i, int j , int k ,
84 std::size_t global_index,
85 int complnum,
86 double depth,
88 double CF,
89 double Kh,
90 double rw,
91 double r0,
92 double skin_factor,
93 const int satTableId,
94 const Direction direction,
95 const CTFKind ctf_kind,
96 const std::size_t sort_value,
97 const double segDistStart,
98 const double segDistEnd,
99 const bool defaultSatTabId);
100
101 Connection(const RestartIO::RstConnection& rst_connection, const EclipseGrid& grid, const FieldPropsManager& fp);
102
104
105 bool attachedToSegment() const;
106 bool sameCoordinate(const int i, const int j, const int k) const;
107 int getI() const;
108 int getJ() const;
109 int getK() const;
110 std::size_t global_index() const;
111 State state() const;
112 Direction dir() const;
113 double depth() const;
114 int satTableId() const;
115 int complnum() const;
116 int segment() const;
117 double CF() const;
118 double Kh() const;
119 double rw() const;
120 double r0() const;
121 double skinFactor() const;
122 CTFKind kind() const;
123
125 void setComplnum(int compnum);
126 void scaleWellPi(double wellPi);
127 void updateSegmentRST(int segment_number_arg,
128 double center_depth_arg);
129 void updateSegment(int segment_number_arg,
130 double center_depth_arg,
131 std::size_t compseg_insert_index,
132 double start,
133 double end);
134 std::size_t sort_value() const;
135 const bool& getDefaultSatTabId() const;
136 void setDefaultSatTabId(bool id);
137 const double& getSegDistStart() const;
138 const double& getSegDistEnd() const;
141 {
142 return this->m_ctfkind == CTFKind::DeckValue;
143 }
144
145 bool operator==( const Connection& ) const;
146 bool operator!=( const Connection& ) const;
147
148 template<class Serializer>
149 void serializeOp(Serializer& serializer)
150 {
151 serializer(direction);
152 serializer(center_depth);
153 serializer(open_state);
154 serializer(sat_tableId);
155 serializer(m_complnum);
156 serializer(m_CF);
157 serializer(m_Kh);
158 serializer(m_rw);
159 serializer(m_r0);
160 serializer(m_skin_factor);
161 serializer(ijk);
162 serializer(m_global_index);
163 serializer(m_ctfkind);
164 serializer(m_sort_value);
165 serializer(m_segDistStart);
166 serializer(m_segDistEnd);
167 serializer(m_defaultSatTabId);
168 serializer(segment_number);
169 }
170
171 private:
172 Direction direction;
173 double center_depth;
174 State open_state;
175 int sat_tableId;
176 int m_complnum;
177 double m_CF;
178 double m_Kh;
179 double m_rw;
180 double m_r0;
181 double m_skin_factor;
182
183 std::array<int,3> ijk;
184 CTFKind m_ctfkind;
185 std::size_t m_global_index;
186 /*
187 The sort_value member is a peculiar quantity. The connections are
188 assembled in the WellConnections class. During the lifetime of the
189 connections there are three different sort orders which are all
190 relevant:
191
192 input: This is the ordering implied be the order of the
193 connections in the input deck.
194
195 simulation: This is the ordering the connections have in
196 WellConnections container during the simulation and RFT output.
197
198 restart: This is the ordering the connections have when they are
199 written out to a restart file.
200
201 Exactly what consitutes input, simulation and restart ordering, and
202 how the connections transition between the three during application
203 lifetime is different from MSW and normal wells.
204
205 normal wells: For normal wells the simulation order is given by the
206 COMPORD keyword, and then when the connections are serialized to the
207 restart file they are written in input order; i.e. we have:
208
209 input == restart and simulation given COMPORD
210
211 To recover the input order when creating the restart files the
212 sort_value member corresponds to the insert index for normal wells.
213
214 MSW wells: For MSW wells the wells simulator order[*] is given by
215 COMPSEGS keyword, the COMPORD keyword is ignored. The connections are
216 sorted in WellConnections::order() and then retain that order for all
217 eternity i.e.
218
219 input and simulation == restart
220
221 Now the important point is that the COMPSEGS detail used to perform
222 this sorting is not available when loading from a restart file, but
223 then the connections are already sorted correctly. I.e. *after* a
224 restart we will have:
225
226 input(from restart) == simulation == restart
227
228 The sort_value member is used to sort the connections into restart
229 ordering. In the case of normal wells this corresponds to recovering
230 the input order, whereas for MSW wells this is equivalent to the
231 simulation order.
232
233 [*]: For MSW wells the topology is given by the segments and entered
234 explicitly, so the truth is probably that the storage order
235 during simulation makes no difference?
236 */
237
238 std::size_t m_sort_value;
239 double m_segDistStart;
240 double m_segDistEnd;
241 bool m_defaultSatTabId;
242
243 // related segment number
244 // 0 means the completion is not related to segment
245 int segment_number = 0;
246
247 static std::string CTFKindToString(const CTFKind);
248 };
249}
250
251#endif /* COMPLETION_HPP_ */
252
const char *const string
Definition: cJSON.h:170
Definition: parser/eclipse/EclipseState/Schedule/Well/connection.hpp:45
bool sameCoordinate(const int i, const int j, const int k) const
bool operator!=(const Connection &) const
int complnum() const
static Connection serializeObject()
Direction dir() const
std::size_t sort_value() const
Connection(const RestartIO::RstConnection &rst_connection, const EclipseGrid &grid, const FieldPropsManager &fp)
State state() const
static std::string Direction2String(const Direction enumValue)
int segment() const
double skinFactor() const
void setComplnum(int compnum)
Direction
Definition: parser/eclipse/EclipseState/Schedule/Well/connection.hpp:57
bool attachedToSegment() const
static Order OrderFromString(const std::string &comporderStringValue)
const bool & getDefaultSatTabId() const
void updateSegment(int segment_number_arg, double center_depth_arg, std::size_t compseg_insert_index, double start, double end)
static State StateFromString(const std::string &stringValue)
CTFKind kind() const
std::string str() const
double r0() const
int satTableId() const
double CF() const
CTFKind
Definition: parser/eclipse/EclipseState/Schedule/Well/connection.hpp:76
const double & getSegDistEnd() const
void setState(State state)
double depth() const
Order
Definition: parser/eclipse/EclipseState/Schedule/Well/connection.hpp:67
static Direction DirectionFromString(const std::string &stringValue)
static const std::string Order2String(Order enumValue)
std::size_t global_index() const
static const std::string State2String(State enumValue)
bool ctfAssignedFromInput() const
Definition: parser/eclipse/EclipseState/Schedule/Well/connection.hpp:140
bool operator==(const Connection &) const
int getJ() const
const double & getSegDistStart() const
void setDefaultSatTabId(bool id)
int getK() const
double Kh() const
State
Definition: parser/eclipse/EclipseState/Schedule/Well/connection.hpp:47
void serializeOp(Serializer &serializer)
Definition: parser/eclipse/EclipseState/Schedule/Well/connection.hpp:149
void scaleWellPi(double wellPi)
int getI() const
double rw() const
void updateSegmentRST(int segment_number_arg, double center_depth_arg)
Connection(int i, int j, int k, std::size_t global_index, int complnum, double depth, State state, double CF, double Kh, double rw, double r0, double skin_factor, const int satTableId, const Direction direction, const CTFKind ctf_kind, const std::size_t sort_value, const double segDistStart, const double segDistEnd, const bool defaultSatTabId)
Definition: DeckValue.hpp:29
Definition: EclipseGrid.hpp:54
Definition: FieldPropsManager.hpp:33
Definition: Serializer.hpp:38
Definition: A.hpp:4
@ end
Definition: ActionValue.hpp:20
Definition: io/eclipse/rst/connection.hpp:33