24#include <unordered_map>
27#include <opm/input/eclipse/Schedule/UDQ/UDQSet.hpp>
41 bool has(
const std::string& key)
const;
44 bool has_well_var(
const std::string& well,
const std::string& key)
const;
45 bool has_group_var(
const std::string& group,
const std::string& key)
const;
47 double get(
const std::string& key)
const;
48 double get_group_var(
const std::string& well,
const std::string& var)
const;
49 double get_well_var(
const std::string& well,
const std::string& var)
const;
50 void add_define(std::size_t report_step,
const std::string& udq_key,
const UDQSet& result);
51 void add_assign(std::size_t report_step,
const std::string& udq_key,
const UDQSet& result);
52 bool assign(std::size_t report_step,
const std::string& udq_key)
const;
53 bool define(
const std::string& udq_key, std::pair<UDQUpdate, std::size_t> update_status)
const;
54 double undefined_value()
const;
56 bool operator==(
const UDQState& other)
const;
58 static UDQState serializationTestObject() {
61 st.scalar_values = {{
"FU1", 100}, {
"FU2", 200}};
62 st.assignments = {{
"GU1", 99}, {
"GU2", 199}};
63 st.defines = {{
"DU1", 299}, {
"DU2", 399}};
65 st.well_values.emplace(
"W1", std::unordered_map<std::string, double>{{
"U1", 100}, {
"U2", 200}});
66 st.well_values.emplace(
"W2", std::unordered_map<std::string, double>{{
"U1", 700}, {
"32", 600}});
68 st.group_values.emplace(
"G1", std::unordered_map<std::string, double>{{
"U1", 100}, {
"U2", 200}});
69 st.group_values.emplace(
"G2", std::unordered_map<std::string, double>{{
"U1", 700}, {
"32", 600}});
74 template<
class Serializer>
75 void pack_unpack_wgmap(
Serializer& serializer, std::unordered_map<std::string, std::unordered_map<std::string, double>>& wgmap) {
76 std::size_t map_size = wgmap.size();
79 for (
auto& [udq_key, values] : wgmap) {
84 for (std::size_t index=0; index < map_size; index++) {
86 std::unordered_map<std::string, double> inner_map;
88 serializer(inner_map);
90 wgmap.emplace(udq_key, inner_map);
95 template<
class Serializer>
98 serializer(this->undef_value);
99 serializer(this->scalar_values);
100 serializer(this->assignments);
101 serializer(this->defines);
103 pack_unpack_wgmap(serializer, this->well_values);
104 pack_unpack_wgmap(serializer, this->group_values);
109 void add(
const std::string& udq_key,
const UDQSet& result);
110 double get_wg_var(
const std::string& well,
const std::string& key, UDQVarType var_type)
const;
112 std::unordered_map<std::string, double> scalar_values;
113 std::unordered_map<std::string, std::unordered_map<std::string, double>> well_values;
114 std::unordered_map<std::string, std::unordered_map<std::string, double>> group_values;
115 std::unordered_map<std::string, std::size_t> assignments;
116 std::unordered_map<std::string, std::size_t> defines;
Class for (de-)serializing.
Definition: Serializer.hpp:75
bool isSerializing() const
Returns true if we are currently doing a serialization operation.
Definition: Serializer.hpp:174
Definition: UDQSet.hpp:63
Definition: UDQState.hpp:36
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition: Exceptions.hpp:29