21 #ifndef UDQASSIGN_HPP_
22 #define UDQASSIGN_HPP_
25 #include <unordered_set>
28 #include <opm/input/eclipse/Schedule/UDQ/UDQSet.hpp>
29 #include <opm/input/eclipse/Schedule/UDQ/UDQEnums.hpp>
44 std::vector<std::string> input_selector;
45 std::unordered_set<std::string> rst_selector;
47 std::size_t report_step;
51 AssignRecord(
const std::vector<std::string>& selector,
double value_arg, std::size_t report_step_arg)
52 : input_selector(selector)
54 , report_step(report_step_arg)
57 AssignRecord(
const std::unordered_set<std::string>& selector,
double value_arg, std::size_t report_step_arg)
58 : rst_selector(selector)
60 , report_step(report_step_arg)
63 void eval(
UDQSet& values)
const {
64 if (this->input_selector.empty() && this->rst_selector.empty())
65 values.assign( this->value );
67 if (this->rst_selector.empty())
68 values.assign(this->input_selector[0], this->value);
70 for (
const auto& wgname : this->rst_selector)
71 values.assign(wgname, this->value);
77 return input_selector == data.input_selector &&
78 rst_selector == data.rst_selector &&
79 report_step == data.report_step &&
83 template<
class Serializer>
86 serializer(input_selector);
87 serializer(rst_selector);
89 serializer(report_step);
94 UDQAssign(
const std::string& keyword,
const std::vector<std::string>& selector,
double value, std::size_t report_step);
95 UDQAssign(
const std::string& keyword,
const std::unordered_set<std::string>& selector,
double value, std::size_t report_step);
97 static UDQAssign serializationTestObject();
99 const std::string& keyword()
const;
100 UDQVarType var_type()
const;
101 void add_record(
const std::vector<std::string>& selector,
double value, std::size_t report_step);
102 void add_record(
const std::unordered_set<std::string>& rst_selector,
double value, std::size_t report_step);
103 UDQSet eval(
const std::vector<std::string>& wells)
const;
105 std::size_t report_step()
const;
107 bool operator==(
const UDQAssign& data)
const;
109 template<
class Serializer>
112 serializer(m_keyword);
113 serializer(m_var_type);
118 std::string m_keyword;
119 UDQVarType m_var_type;
120 std::vector<AssignRecord> records;
Class for (de-)serializing.
Definition: Serializer.hpp:75
Definition: UDQAssign.hpp:34
Definition: UDQSet.hpp:63
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition: Exceptions.hpp:29
Definition: UDQAssign.hpp:43