My Project
Loading...
Searching...
No Matches
UDQAssign.hpp
1/*
2 Copyright 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 UDQASSIGN_HPP_
21#define UDQASSIGN_HPP_
22
23#include <opm/input/eclipse/Schedule/UDQ/UDQEnums.hpp>
24#include <opm/input/eclipse/Schedule/UDQ/UDQSet.hpp>
25
26#include <cstddef>
27#include <functional>
28#include <string>
29#include <unordered_set>
30#include <utility>
31#include <vector>
32
33namespace Opm::RestartIO {
34 class RstUDQ;
35}
36
37namespace Opm {
38
41{
42public:
44 using VString = std::vector<std::string>;
45
48 using VEnumItems = std::vector<UDQSet::EnumeratedItems>;
49
52 using WGNameMatcher = std::function<VString(const VString&)>;
53
56 using ItemMatcher = std::function<VEnumItems(const VString&)>;
57
59 UDQAssign() = default;
60
76 UDQAssign(const std::string& keyword,
77 const VString& selector,
78 double value,
79 std::size_t report_step);
80
95 UDQAssign(const std::string& keyword,
96 const VEnumItems& selector,
97 double value,
98 std::size_t report_step);
99
116 UDQAssign(const std::string& keyword,
117 VEnumItems&& selector,
118 double value,
119 std::size_t report_step);
120
133 UDQAssign(const std::string& keyword,
134 const RestartIO::RstUDQ& assignRst,
135 const std::size_t report_step);
136
139
141 const std::string& keyword() const;
142
144 UDQVarType var_type() const;
145
159 void add_record(const VString& selector,
160 double value,
161 std::size_t report_step);
162
175 void add_record(const VEnumItems& selector,
176 double value,
177 std::size_t report_step);
178
193 void add_record(VEnumItems&& selector,
194 double value,
195 std::size_t report_step);
196
208 void add_record(const RestartIO::RstUDQ& assignRst,
209 const std::size_t report_step);
210
222 UDQSet eval(const VEnumItems& items) const;
223
235 UDQSet eval(const VString& wells) const;
236
244 UDQSet eval() const;
245
260 UDQSet eval(const VString& wgNames, WGNameMatcher matcher) const;
261
276 UDQSet eval(const VEnumItems& items, ItemMatcher matcher) const;
277
284 std::size_t report_step() const;
285
292 bool operator==(const UDQAssign& data) const;
293
299 template<class Serializer>
300 void serializeOp(Serializer& serializer)
301 {
302 serializer(m_keyword);
303 serializer(m_var_type);
304 serializer(records);
305 }
306
307private:
308 // If the same keyword is assigned several times the different
309 // assignment records are assembled in one UDQAssign instance. This is
310 // an attempt to support restart in a situation where a full UDQ ASSIGN
311 // statement can be swapped with a UDQ DEFINE statement.
312 struct AssignRecord
313 {
321 VString input_selector{};
322
330 VEnumItems numbered_selector{};
331
333 double value{};
334
339 std::size_t report_step{};
340
342 AssignRecord() = default;
343
357 AssignRecord(const VString& selector,
358 const double value_arg,
359 const std::size_t report_step_arg)
360 : input_selector(selector)
361 , value (value_arg)
362 , report_step (report_step_arg)
363 {}
364
377 AssignRecord(const VEnumItems& selector,
378 const double value_arg,
379 const std::size_t report_step_arg)
380 : numbered_selector(selector)
381 , value (value_arg)
382 , report_step (report_step_arg)
383 {}
384
399 AssignRecord(VEnumItems&& selector,
400 const double value_arg,
401 const std::size_t report_step_arg)
402 : numbered_selector(std::move(selector))
403 , value (value_arg)
404 , report_step (report_step_arg)
405 {}
406
413 void eval(UDQSet& values) const;
414
422 void eval(WGNameMatcher matcher, UDQSet& values) const;
423
430 void eval(ItemMatcher matcher, UDQSet& values) const;
431
438 bool operator==(const AssignRecord& data) const;
439
445 template<class Serializer>
446 void serializeOp(Serializer& serializer)
447 {
448 serializer(this->input_selector);
449 serializer(this->numbered_selector);
450 serializer(this->value);
451 serializer(this->report_step);
452 }
453
454 private:
463 void assignEnumeration(const VEnumItems& items, UDQSet& values) const;
464 };
465
467 std::string m_keyword{};
468
470 UDQVarType m_var_type{UDQVarType::NONE};
471
473 std::vector<AssignRecord> records{};
474
484 void add_well_or_group_records(const RestartIO::RstUDQ& assignRst,
485 const std::size_t report_step);
486
495 void add_segment_records(const RestartIO::RstUDQ& assignRst,
496 const std::size_t report_step);
497};
498
499} // namespace Opm
500
501#endif // UDQASSIGN_HPP_
Class for (de-)serializing.
Definition Serializer.hpp:94
Representation of a UDQ ASSIGN statement.
Definition UDQAssign.hpp:41
std::size_t report_step() const
Time at which this assignment happens.
Definition UDQAssign.cpp:214
std::function< VString(const VString &)> WGNameMatcher
Call-back function type for a well/group name matcher.
Definition UDQAssign.hpp:52
std::vector< std::string > VString
Type alias for a vector of strings. Simplifies function signatures.
Definition UDQAssign.hpp:44
UDQSet eval() const
Construct scalar UDQ set for a scalar UDQ assignment.
Definition UDQAssign.cpp:257
const std::string & keyword() const
Name of UDQ to which this assignment applies.
Definition UDQAssign.cpp:204
static UDQAssign serializationTestObject()
Create a serialisation test object.
Definition UDQAssign.cpp:138
void add_record(const VString &selector, double value, std::size_t report_step)
Add new record to existing UDQ assignment.
Definition UDQAssign.cpp:151
std::function< VEnumItems(const VString &)> ItemMatcher
Call-back function type for a matcher of enumerated items.
Definition UDQAssign.hpp:56
UDQSet eval(const VString &wells) const
Apply current assignment to a selection of named items.
void serializeOp(Serializer &serializer)
Convert between byte array and object representation.
Definition UDQAssign.hpp:300
void add_record(const VEnumItems &selector, double value, std::size_t report_step)
Add new record to existing UDQ assignment.
std::vector< UDQSet::EnumeratedItems > VEnumItems
Type alias for a vector of enumerated items.
Definition UDQAssign.hpp:48
UDQVarType var_type() const
Kind of UDQ to which this assignment applies.
Definition UDQAssign.cpp:209
UDQAssign()=default
Default constructor.
bool operator==(const UDQAssign &data) const
Equality predicate.
Definition UDQAssign.cpp:312
UDQAssign(const std::string &keyword, const VEnumItems &selector, double value, std::size_t report_step)
Constructor.
Definition UDQSet.hpp:187
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30