20#ifndef OPM_AGGREGATE_UDQ_DATA_HPP
21#define OPM_AGGREGATE_UDQ_DATA_HPP
25#include <opm/io/eclipse/PaddedOutputString.hpp>
41namespace Opm::RestartIO::Helpers {
48 void captureDeclaredUDQData(
const Schedule& sched,
49 const std::size_t simStep,
51 const std::vector<int>& inteHead);
53 const std::vector<int>& getIUDQ()
const
55 return this->iUDQ_.
data();
59 const std::optional<WindowedArray<int>>&
getIUAD()
const
64 const std::vector<EclIO::PaddedOutputString<8>>& getZUDN()
const
66 return this->zUDN_.data();
69 const std::vector<EclIO::PaddedOutputString<8>>& getZUDL()
const
71 return this->zUDL_.data();
76 const std::optional<WindowedArray<int>>&
getIGPH()
const
82 const std::optional<WindowedArray<int>>&
getIUAP()
const
88 const std::optional<WindowedArray<double>>&
getDUDF()
const
94 const std::optional<WindowedArray<double>>&
getDUDG()
const
100 const std::optional<WindowedArray<double>>&
getDUDW()
const
115 std::optional<WindowedArray<int>> iUAD_;
131 std::optional<WindowedArray<int>> iGPH_;
136 std::optional<WindowedArray<int>> iUAP_;
141 std::optional<WindowedArray<double>> dUDF_{};
147 std::optional<WindowedArray<double>> dUDG_{};
153 std::optional<WindowedArray<double>> dUDW_{};
155 void collectUserDefinedQuantities(
const std::vector<UDQInput>& udqInput,
156 const std::vector<int>& inteHead);
158 void collectUserDefinedArguments(
const Schedule& sched,
159 const std::size_t simStep,
160 const std::vector<int>& inteHead);
162 void collectFieldUDQValues(
const std::vector<UDQInput>& udqInput,
163 const UDQState& udq_state,
164 const int expectNumFieldUDQs);
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);
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);
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