opm-common
Loading...
Searching...
No Matches
AggregateWellData.hpp
1/*
2 Copyright (c) 2018 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#ifndef OPM_AGGREGATE_WELL_DATA_HPP
21#define OPM_AGGREGATE_WELL_DATA_HPP
22
24
25#include <opm/io/eclipse/PaddedOutputString.hpp>
26
27#include <cstddef>
28#include <string>
29#include <vector>
30
31namespace Opm {
32 class Schedule;
33 class EclipseGrid;
34 class SummaryState;
35 class UnitSystem;
36 class WellTestState;
37 class TracerConfig;
38 namespace Action {
39 class State;
40 }
41} // Opm
42
43namespace Opm { namespace data {
44 class Wells;
45}} // Opm::data
46
47namespace Opm { namespace RestartIO { namespace Helpers {
48
49 class AggregateWellData
50 {
51 public:
52 explicit AggregateWellData(const std::vector<int>& inteHead);
53
54 void captureDeclaredWellData(const Schedule& sched,
55 const TracerConfig& tracer,
56 const std::size_t sim_step,
57 const Opm::Action::State& action_state,
58 const Opm::WellTestState& wtest_state,
59 const Opm::SummaryState& smry,
60 const std::vector<int>& inteHead);
61
62 void captureDeclaredWellData(const Schedule& sched,
63 const EclipseGrid& grid,
64 const TracerConfig& tracer,
65 const std::size_t sim_step,
66 const Opm::Action::State& action_state,
67 const Opm::WellTestState& wtest_state,
68 const Opm::SummaryState& smry,
69 const std::vector<int>& inteHead);
70
71 void captureDeclaredWellDataLGR(const Schedule& sched,
72 const EclipseGrid& grid,
73 const TracerConfig& tracer,
74 const std::size_t sim_step,
75 const Opm::Action::State& action_state,
76 const Opm::WellTestState& wtest_state,
77 const Opm::SummaryState& smry,
78 const std::vector<int>& inteHead,
79 const std::string& lgr_tag);
80
81 void captureDynamicWellData(const Opm::Schedule& sched,
82 const TracerConfig& tracer,
83 const std::size_t sim_step,
84 const Opm::data::Wells& xw,
85 const Opm::SummaryState& smry);
86
87 void captureDynamicWellDataLGR(const Opm::Schedule& sched,
88 const TracerConfig& tracer,
89 const std::size_t sim_step,
90 const Opm::data::Wells& xw,
91 const Opm::SummaryState& smry,
92 const std::string& lgr_tag);
93
94
96 const std::vector<int>& getIWell() const
97 {
98 return this->iWell_.data();
99 }
100
102 const std::vector<float>& getSWell() const
103 {
104 return this->sWell_.data();
105 }
106
108 const std::vector<double>& getXWell() const
109 {
110 return this->xWell_.data();
111 }
112
114 const std::vector<EclIO::PaddedOutputString<8>>& getZWell() const
115 {
116 return this->zWell_.data();
117 }
118
120 const std::vector<int>& getLGWell() const
121 {
122 return this->lgWell_.data();
123 }
124
125
126
127 private:
129 WindowedArray<int> iWell_;
130
133
136
139
141 WindowedArray<int> lgWell_;
142
144 int nWGMax_;
145 };
146
147}}} // Opm::RestartIO::Helpers
148
149#endif // OPM_AGGREGATE_WELL_DATA_HPP
Provide facilities to simplify constructing restart vectors such as IWEL or RSEG.
Management information about the current run's ACTION system, especially concerning the number of tim...
Definition State.hpp:51
About cell information and dimension: The actual grid information is held in a pointer to an ERT ecl_...
Definition EclipseGrid.hpp:61
const std::vector< int > & getIWell() const
Retrieve Integer Well Data Array.
Definition AggregateWellData.hpp:96
const std::vector< double > & getXWell() const
Retrieve Floating-Point (Double Precision) Well Data Array.
Definition AggregateWellData.hpp:108
const std::vector< float > & getSWell() const
Retrieve Floating-Point (Real) Well Data Array.
Definition AggregateWellData.hpp:102
const std::vector< int > & getLGWell() const
Retrieve Interger LGWell Data Array.
Definition AggregateWellData.hpp:120
const std::vector< EclIO::PaddedOutputString< 8 > > & getZWell() const
Retrieve Character Well Data Array.
Definition AggregateWellData.hpp:114
Provide read-only and read/write access to constantly sized portions/windows of a linearised buffer w...
Definition WindowedArray.hpp:50
Definition Schedule.hpp:101
Definition SummaryState.hpp:72
Definition TracerConfig.hpp:33
Definition UnitSystem.hpp:34
Definition WellTestState.hpp:65
Definition Wells.hpp:1199
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30