20#ifndef OPM_AGGREGATE_UDQ_DATA_HPP
21#define OPM_AGGREGATE_UDQ_DATA_HPP
25#include <opm/io/eclipse/PaddedOutputString.hpp>
42namespace Opm::RestartIO::Helpers {
49 void captureDeclaredUDQData(
const Schedule& sched,
50 const std::size_t simStep,
52 const std::vector<int>& inteHead);
54 const std::vector<int>& getIUDQ()
const
56 return this->iUDQ_.
data();
60 const std::optional<WindowedArray<int>>&
getIUAD()
const
65 const std::vector<EclIO::PaddedOutputString<8>>& getZUDN()
const
67 return this->zUDN_.data();
70 const std::vector<EclIO::PaddedOutputString<8>>& getZUDL()
const
72 return this->zUDL_.data();
77 const std::optional<WindowedArray<int>>&
getIGPH()
const
83 const std::optional<WindowedArray<int>>&
getIUAP()
const
89 const std::optional<WindowedArray<double>>&
getDUDF()
const
95 const std::optional<WindowedArray<double>>&
getDUDG()
const
101 const std::optional<WindowedMatrix<double>>&
getDUDS()
const
107 const std::optional<WindowedArray<double>>&
getDUDW()
const
122 std::optional<WindowedArray<int>> iUAD_{};
127 WindowedArray<EclIO::PaddedOutputString<8>> zUDN_;
132 WindowedArray<EclIO::PaddedOutputString<8>> zUDL_;
138 std::optional<WindowedArray<int>> iGPH_{};
143 std::optional<WindowedArray<int>> iUAP_{};
148 std::optional<WindowedArray<double>> dUDF_{};
154 std::optional<WindowedArray<double>> dUDG_{};
161 std::optional<WindowedMatrix<double>> dUDS_{};
167 std::optional<WindowedArray<double>> dUDW_{};
169 void collectUserDefinedQuantities(
const std::vector<UDQInput>& udqInput,
170 const std::vector<int>& inteHead);
172 void collectUserDefinedArguments(
const Schedule& sched,
173 const std::size_t simStep,
174 const std::vector<int>& inteHead);
176 void collectFieldUDQValues(
const std::vector<UDQInput>& udqInput,
177 const UDQState& udq_state,
178 const int expectNumFieldUDQs);
180 void collectGroupUDQValues(
const std::vector<UDQInput>& udqInput,
181 const UDQState& udqState,
182 const std::size_t ngmax,
183 const std::vector<const Group*>& groups,
184 const int expectedNumGroupUDQs);
186 void collectSegmentUDQValues(
const std::vector<UDQInput>& udqInput,
187 const UDQState& udqState,
188 const std::vector<std::string>& msWells);
190 void collectWellUDQValues(
const std::vector<UDQInput>& udqInput,
191 const UDQState& udqState,
192 const std::size_t nwmax,
193 const std::vector<std::string>& wells,
194 const int expectedNumWellUDQs);
202 void collectIUAD(
const UDQActive& udqActive,
203 const std::size_t expectNumIUAD);
211 void collectIUAP(
const std::vector<int>& wgIndex,
212 const std::size_t expectNumIUAP);
221 void collectIGPH(
const std::vector<int>& phase_vector,
222 const std::size_t expectNumIGPH);
Provide facilities to simplify constructing restart vectors such as IWEL or RSEG.
Definition AggregateUDQData.hpp:45
const std::optional< WindowedArray< int > > & getIUAP() const
Associate well/group IDs for IUAD. Nullopt if no UDAs in use.
Definition AggregateUDQData.hpp:83
const std::optional< WindowedArray< int > > & getIUAD() const
Retrieve UDA descriptive data. Nullopt if no UDAs in use.
Definition AggregateUDQData.hpp:60
const std::optional< WindowedArray< int > > & getIGPH() const
Retrive group level injection phase UDAs.
Definition AggregateUDQData.hpp:77
const std::optional< WindowedArray< double > > & getDUDW() const
Retrieve values of well level UDQs. Nullopt if no such UDQs exist.
Definition AggregateUDQData.hpp:107
const std::optional< WindowedArray< double > > & getDUDG() const
Retrieve values of group level UDQs. Nullopt if no such UDQs exist.
Definition AggregateUDQData.hpp:95
const std::optional< WindowedMatrix< double > > & getDUDS() const
Retrieve values of segment level UDQs. Nullopt if no such UDQs exist.
Definition AggregateUDQData.hpp:101
const std::optional< WindowedArray< double > > & getDUDF() const
Retrieve values of field level UDQs. Nullopt if no such UDQs exist.
Definition AggregateUDQData.hpp:89
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:136
Definition Schedule.hpp:101
Collection of UDQ and UDA related dimension queries.
Definition UDQDims.hpp:40
Definition UDQState.hpp:40
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30