My Project
Loading...
Searching...
No Matches
AggregateUDQData.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_UDQ_DATA_HPP
21#define OPM_AGGREGATE_UDQ_DATA_HPP
22
24
25#include <opm/io/eclipse/PaddedOutputString.hpp>
26
27#include <cstddef>
28#include <optional>
29#include <string>
30#include <vector>
31
32namespace Opm {
33 class Group;
34 class Schedule;
35 class UDQConfig;
36 class UDQDims;
37 class UDQInput;
38 class UDQState;
39} // Opm
40
41namespace Opm::RestartIO::Helpers {
42
44{
45public:
46 explicit AggregateUDQData(const UDQDims& udqDims);
47
48 void captureDeclaredUDQData(const Schedule& sched,
49 const std::size_t simStep,
50 const UDQState& udqState,
51 const std::vector<int>& inteHead);
52
53 const std::vector<int>& getIUDQ() const
54 {
55 return this->iUDQ_.data();
56 }
57
59 const std::optional<WindowedArray<int>>& getIUAD() const
60 {
61 return this->iUAD_;
62 }
63
64 const std::vector<EclIO::PaddedOutputString<8>>& getZUDN() const
65 {
66 return this->zUDN_.data();
67 }
68
69 const std::vector<EclIO::PaddedOutputString<8>>& getZUDL() const
70 {
71 return this->zUDL_.data();
72 }
73
76 const std::optional<WindowedArray<int>>& getIGPH() const
77 {
78 return this->iGPH_;
79 }
80
82 const std::optional<WindowedArray<int>>& getIUAP() const
83 {
84 return this->iUAP_;
85 }
86
88 const std::optional<WindowedArray<double>>& getDUDF() const
89 {
90 return this->dUDF_;
91 }
92
94 const std::optional<WindowedArray<double>>& getDUDG() const
95 {
96 return this->dUDG_;
97 }
98
100 const std::optional<WindowedArray<double>>& getDUDW() const
101 {
102 return this->dUDW_;
103 }
104
105private:
109 WindowedArray<int> iUDQ_;
110
115 std::optional<WindowedArray<int>> iUAD_;
116
121
126
131 std::optional<WindowedArray<int>> iGPH_;
132
136 std::optional<WindowedArray<int>> iUAP_;
137
141 std::optional<WindowedArray<double>> dUDF_{};
142
147 std::optional<WindowedArray<double>> dUDG_{};
148
153 std::optional<WindowedArray<double>> dUDW_{};
154
155 void collectUserDefinedQuantities(const std::vector<UDQInput>& udqInput,
156 const std::vector<int>& inteHead);
157
158 void collectUserDefinedArguments(const Schedule& sched,
159 const std::size_t simStep,
160 const std::vector<int>& inteHead);
161
162 void collectFieldUDQValues(const std::vector<UDQInput>& udqInput,
163 const UDQState& udq_state,
164 const int expectNumFieldUDQs);
165
166 void collectGroupUDQValues(const std::vector<UDQInput>& udqInput,
167 const UDQState& udqState,
168 const std::size_t ngmax,
169 const std::vector<const Group*>& groups,
170 const int expectedNumGroupUDQs);
171
172 void collectWellUDQValues(const std::vector<UDQInput>& udqInput,
173 const UDQState& udqState,
174 const std::size_t nwmax,
175 const std::vector<std::string>& wells,
176 const int expectedNumWellUDQs);
177};
178
179} // Opm::RestartIO::Helpers
180
181#endif // OPM_AGGREGATE_UDQ_DATA_HPP
Provide facilities to simplify constructing restart vectors such as IWEL or RSEG.
Definition AggregateUDQData.hpp:44
const std::optional< WindowedArray< int > > & getIUAP() const
Associate well/group IDs for IUAD. Nullopt if no UDAs in use.
Definition AggregateUDQData.hpp:82
const std::optional< WindowedArray< int > > & getIUAD() const
Retrieve UDA descriptive data. Nullopt if no UDAs in use.
Definition AggregateUDQData.hpp:59
const std::optional< WindowedArray< int > > & getIGPH() const
Retrive group level injection phase UDAs.
Definition AggregateUDQData.hpp:76
const std::optional< WindowedArray< double > > & getDUDW() const
Retrieve values of well level UDQs. Nullopt if no such UDQs exist.
Definition AggregateUDQData.hpp:100
const std::optional< WindowedArray< double > > & getDUDG() const
Retrieve values of group level UDQs. Nullopt if no such UDQs exist.
Definition AggregateUDQData.hpp:94
const std::optional< WindowedArray< double > > & getDUDF() const
Retrieve values of field level UDQs. Nullopt if no such UDQs exist.
Definition AggregateUDQData.hpp:88
Provide read-only and read/write access to constantly sized portions/windows of a linearised buffer w...
Definition WindowedArray.hpp:50
const std::vector< T > & data() const
Get read-only access to full, linearised data items for all windows.
Definition WindowedArray.hpp:131
Definition Schedule.hpp:89
Definition UDQDims.hpp:37
Definition UDQState.hpp:37
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30