opm-common
Loading...
Searching...
No Matches
UDQConfig.hpp
1/*
2 Copyright 2019 Equinor 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 UDQINPUT_HPP_
21#define UDQINPUT_HPP_
22
23#include <opm/input/eclipse/Schedule/UDQ/UDQAssign.hpp>
24#include <opm/input/eclipse/Schedule/UDQ/UDQDefine.hpp>
25#include <opm/input/eclipse/Schedule/UDQ/UDQEnums.hpp>
26#include <opm/input/eclipse/Schedule/UDQ/UDQFunctionTable.hpp>
27#include <opm/input/eclipse/Schedule/UDQ/UDQInput.hpp>
28#include <opm/input/eclipse/Schedule/UDQ/UDQParams.hpp>
29#include <opm/input/eclipse/Schedule/UDQ/UDT.hpp>
30
31#include <opm/input/eclipse/EclipseState/Util/OrderedMap.hpp>
32#include <opm/input/eclipse/EclipseState/Util/IOrderSet.hpp>
33
34#include <array>
35#include <cstddef>
36#include <functional>
37#include <map>
38#include <memory>
39#include <optional>
40#include <string>
41#include <unordered_map>
42#include <unordered_set>
43#include <utility>
44#include <vector>
45
46namespace Opm {
47
48 class DeckRecord;
49 class KeywordLocation;
50 class GroupOrder;
51 class RegionSetMatcher;
52 class Schedule;
53 class SegmentMatcher;
54 class SummaryState;
55 class UDQState;
56 class WellMatcher;
57
58} // namespace Opm
59
60namespace Opm::RestartIO {
61 struct RstState;
62 class RstUDQ;
63} // namespace Opm::RestartIO
64
65namespace Opm {
66
69 {
70 public:
76 {
77 public:
82 {
83 public:
88 using Iter = std::vector<std::string>::const_iterator;
89
93 StringRange() = default;
94
103 StringRange(Iter first, Iter last)
104 : first_ { first }
105 , last_ { last }
106 {}
107
109 auto begin() const { return this->first_; }
110
112 auto end() const { return this->last_; }
113
117 std::vector<std::string> asVector() const
118 {
119 return { this->begin(), this->end() };
120 }
121
122 private:
124 Iter first_{};
125
127 Iter last_{};
128 };
129
133 DynamicSelector() = default;
134
143 DynamicSelector& wells(std::vector<std::string>::const_iterator first,
144 std::vector<std::string>::const_iterator last)
145 {
146 this->wells_.emplace(first, last);
147 return *this;
148 }
149
156 {
157 this->wells_.emplace(srange);
158 return *this;
159 }
160
164 const std::optional<StringRange> wells() const
165 {
166 return this->wells_;
167 }
168
169 private:
173 std::optional<StringRange> wells_{};
174 };
175
177 using RegionSetMatcherFactory = std::function<std::unique_ptr<RegionSetMatcher>()>;
178
180 using SegmentMatcherFactory = std::function<std::unique_ptr<SegmentMatcher>()>;
181
183 UDQConfig() = default;
184
190 explicit UDQConfig(const UDQParams& params);
191
199 UDQConfig(const UDQParams& params, const RestartIO::RstState& rst_state);
200
203
212 const std::string& unit(const std::string& key) const;
213
222 bool has_unit(const std::string& keyword) const;
223
230 bool has_keyword(const std::string& keyword) const;
231
248 std::optional<std::string>
249 add_record(SegmentMatcherFactory create_segment_matcher,
250 const DeckRecord& record,
251 const KeywordLocation& location,
252 std::size_t report_step,
253 const std::optional<DynamicSelector>& dynamic_selector = std::nullopt);
254
263 void add_unit(const std::string& keyword,
264 const std::string& unit);
265
280 void add_update(const std::string& keyword,
281 std::size_t report_step,
282 const KeywordLocation& location,
283 const std::vector<std::string>& data);
284
307 void add_assign(const std::string& quantity,
308 SegmentMatcherFactory create_segment_matcher,
309 const std::vector<std::string>& selector,
310 double value,
311 std::size_t report_step,
312 const std::optional<DynamicSelector>& dynamic_selector = std::nullopt);
313
330 void add_define(const std::string& quantity,
331 const KeywordLocation& location,
332 const std::vector<std::string>& expression,
333 std::size_t report_step);
334
342 void add_table(const std::string& name, UDT udt);
343
355
371
389 void eval_assign(const WellMatcher& wm,
390 const GroupOrder& go,
391 SegmentMatcherFactory create_segment_matcher,
392 SummaryState& st,
393 UDQState& udq_state) const;
394
419 void eval(std::size_t report_step,
420 const WellMatcher& wm,
421 const GroupOrder& go,
422 SegmentMatcherFactory create_segment_matcher,
423 RegionSetMatcherFactory create_region_matcher,
424 SummaryState& st,
425 UDQState& udq_state) const;
426
436 const UDQDefine& define(const std::string& key) const;
437
446 const UDQAssign& assign(const std::string& key) const;
447
450 std::vector<UDQDefine> definitions() const;
451
459 std::vector<UDQDefine> definitions(UDQVarType var_type) const;
460
464 std::vector<UDQInput> input() const;
465
470 void exportTypeCount(std::array<int, static_cast<std::size_t>(UDQVarType::NumTypes)>& count) const;
471
475 std::size_t size() const;
476
485 UDQInput operator[](const std::string& keyword) const;
486
495 UDQInput operator[](std::size_t insert_index) const;
496
498 std::vector<UDQAssign> assignments() const;
499
507 std::vector<UDQAssign> assignments(UDQVarType var_type) const;
508
510 const UDQParams& params() const;
511
513 const UDQFunctionTable& function_table() const;
514
516 const std::unordered_map<std::string, UDT>& tables() const;
517
524 bool operator==(const UDQConfig& config) const;
525
533 void required_summary(std::unordered_set<std::string>& summary_keys) const;
534
540 template<class Serializer>
541 void serializeOp(Serializer& serializer)
542 {
543 serializer(udq_params);
544 serializer(m_definitions);
545 serializer(m_assignments);
546 serializer(m_tables);
547 serializer(units);
548 serializer(define_order);
549 serializer(input_index);
550 serializer(type_count);
551 serializer(pending_assignments_);
552
553 // The UDQFunction table is constant up to udq_params, so we can
554 // just construct a new instance here.
555 if (!serializer.isSerializing()) {
556 udqft = UDQFunctionTable(udq_params);
557 }
558 }
559
560 private:
564 UDQParams udq_params;
565
567 UDQFunctionTable udqft;
568
569 // The following data structures are constrained by compatibility
570 // requirements in our simulation restart files. In particular we
571 // need to control the keyword ordering. In this class the ordering
572 // is maintained mainly by the input_index map which tracks the
573 // insertion order of each keyword and whether the keyword (UDQ
574 // name) is currently DEFINE'ed or ASSIGN'ed.
575
580 std::unordered_map<std::string, UDQDefine> m_definitions{};
581
585 std::unordered_map<std::string, UDQAssign> m_assignments{};
586
590 std::unordered_map<std::string, UDT> m_tables{};
591
598 std::unordered_map<std::string, std::string> units{};
599
603 IOrderSet<std::string> define_order;
604
606 OrderedMap<UDQIndex> input_index;
607
609 std::map<UDQVarType, std::size_t> type_count{};
610
616 mutable std::vector<std::string> pending_assignments_{};
617
626 void add_node(const std::string& quantity, UDQAction action);
627
634 void add_assign(const RestartIO::RstUDQ& udq, const std::size_t report_step);
635
643 void add_define(const RestartIO::RstUDQ& udq, const std::size_t report_step);
644
651 void eval_assign(UDQContext& context) const;
652
664 void eval_define(std::size_t report_step,
665 const UDQState& udq_state,
666 UDQContext& context) const;
667
685 void add_enumerated_assign(const std::string& quantity,
686 SegmentMatcherFactory create_segment_matcher,
687 const std::vector<std::string>& selector,
688 double value,
689 std::size_t report_step);
690
709 void add_assign_wellvar(const std::string& quantity,
710 const std::vector<std::string>& selector,
711 double value,
712 std::size_t report_step,
713 const std::optional<DynamicSelector>& dynamic_selector);
714
724 template <typename... Args>
725 void add_assign_impl(const std::string& quantity,
726 Args&&... args)
727 {
728 // Note: Duplicate 'quantity' arguments is intentional here.
729 // The first is the key in 'm_assignments', while the second is
730 // the first argument to the UDQAssign constructor.
731 const auto& [asgnPos, inserted] =
732 this->m_assignments.try_emplace(quantity, quantity, args...);
733
734 if (! inserted) {
735 // We already have an assignment object for 'quantity'. Add
736 // a new assignment record to that object.
737 asgnPos->second.add_record(std::forward<Args>(args)...);
738 }
739 }
740 };
741
742} // namespace Opm
743
744#endif // UDQINPUT_HPP_
Definition DeckRecord.hpp:32
Collection of group names with built-in ordering.
Definition NameOrder.hpp:69
Definition KeywordLocation.hpp:27
Definition Schedule.hpp:101
Class for (de-)serializing.
Definition Serializer.hpp:94
bool isSerializing() const
Returns true if we are currently doing a serialization operation.
Definition Serializer.hpp:207
Definition SummaryState.hpp:73
Representation of a UDQ ASSIGN statement.
Definition UDQAssign.hpp:41
Random access range of string values.
Definition UDQConfig.hpp:82
StringRange(Iter first, Iter last)
Constructor.
Definition UDQConfig.hpp:103
StringRange()=default
Default constructor.
auto end() const
End of range's elements.
Definition UDQConfig.hpp:112
std::vector< std::string > asVector() const
Convert value range to a std::vector.
Definition UDQConfig.hpp:117
std::vector< std::string >::const_iterator Iter
Random access iterator.
Definition UDQConfig.hpp:88
auto begin() const
Beginning of range's elements.
Definition UDQConfig.hpp:109
DynamicSelector & wells(const StringRange &srange)
Include a range of well names into the dynamic selector.
Definition UDQConfig.hpp:155
const std::optional< StringRange > wells() const
Retrieve current dynamic well set.
Definition UDQConfig.hpp:164
DynamicSelector & wells(std::vector< std::string >::const_iterator first, std::vector< std::string >::const_iterator last)
Include a range of well names into the dynamic selector.
Definition UDQConfig.hpp:143
DynamicSelector()=default
Default constructor.
std::size_t size() const
Total number of active DEFINE and ASSIGN statements.
Definition UDQConfig.cpp:552
const UDQDefine & define(const std::string &key) const
Retrieve defining expression and evaluation object for a single UDQ.
Definition UDQConfig.cpp:485
const UDQParams & params() const
Retrieve run's active UDQ parameters.
Definition UDQConfig.cpp:641
void add_assign(const std::string &quantity, SegmentMatcherFactory create_segment_matcher, const std::vector< std::string > &selector, double value, std::size_t report_step, const std::optional< DynamicSelector > &dynamic_selector=std::nullopt)
Incorporate a UDQ assignment.
Definition UDQConfig.cpp:370
void eval_assign(const WellMatcher &wm, const GroupOrder &go, SegmentMatcherFactory create_segment_matcher, SummaryState &st, UDQState &udq_state) const
Apply all pending assignments.
Definition UDQConfig.cpp:447
static UDQConfig serializationTestObject()
Create a serialisation test object.
Definition UDQConfig.cpp:247
void eval(std::size_t report_step, const WellMatcher &wm, const GroupOrder &go, SegmentMatcherFactory create_segment_matcher, RegionSetMatcherFactory create_region_matcher, SummaryState &st, UDQState &udq_state) const
Compute new values for all UDQs.
Definition UDQConfig.cpp:464
void serializeOp(Serializer &serializer)
Convert between byte array and object representation.
Definition UDQConfig.hpp:541
const std::unordered_map< std::string, UDT > & tables() const
Retrieve run's active user defined tables.
Definition UDQConfig.cpp:651
void add_define(const std::string &quantity, const KeywordLocation &location, const std::vector< std::string > &expression, std::size_t report_step)
Incorporate a UDQ defining expressions.
Definition UDQConfig.cpp:402
UDQConfig()=default
Default constructor.
void exportTypeCount(std::array< int, static_cast< std::size_t >(UDQVarType::NumTypes)> &count) const
Export count of all known UDQ categories in the current run.
Definition UDQConfig.cpp:543
bool clear_update_next_for_new_report_step()
Clear "UPDATE NEXT" flags for all pertinent UDQ definitions.
Definition UDQConfig.cpp:435
void add_unit(const std::string &keyword, const std::string &unit)
Incorporate a unit string for a UDQ.
Definition UDQConfig.cpp:325
bool has_keyword(const std::string &keyword) const
Query whether or not a particular UDQ exists in the collection.
Definition UDQConfig.cpp:278
std::optional< std::string > add_record(SegmentMatcherFactory create_segment_matcher, const DeckRecord &record, const KeywordLocation &location, std::size_t report_step, const std::optional< DynamicSelector > &dynamic_selector=std::nullopt)
Incorporate a single UDQ record into the known collection.
Definition UDQConfig.cpp:285
std::vector< UDQDefine > definitions() const
Retrieve defining expressions and evaluation objects for all known UDQs.
Definition UDQConfig.cpp:495
std::function< std::unique_ptr< RegionSetMatcher >()> RegionSetMatcherFactory
Factory function for constructing region set matchers.
Definition UDQConfig.hpp:177
std::vector< UDQAssign > assignments() const
Retrieve pending assignment objects for all known UDQs.
Definition UDQConfig.cpp:612
bool clear_pending_assignments()
Clear all pending assignments.
Definition UDQConfig.cpp:426
bool has_unit(const std::string &keyword) const
Query whether or not a particular UDQ has an associated unit string.
Definition UDQConfig.cpp:273
void add_update(const std::string &keyword, std::size_t report_step, const KeywordLocation &location, const std::vector< std::string > &data)
Incorporate update status change for a UDQ.
Definition UDQConfig.cpp:340
void add_table(const std::string &name, UDT udt)
Incorporate a user defined table.
Definition UDQConfig.cpp:421
const UDQAssign & assign(const std::string &key) const
Retrieve any pending assignment object for a single UDQ.
Definition UDQConfig.cpp:490
UDQInput operator[](const std::string &keyword) const
Unprocessed input object for named quantity.
Definition UDQConfig.cpp:564
const std::string & unit(const std::string &key) const
Retrieve unit string for a particular UDQ.
Definition UDQConfig.cpp:263
std::function< std::unique_ptr< SegmentMatcher >()> SegmentMatcherFactory
Factory function for constructing segment set matchers.
Definition UDQConfig.hpp:180
std::vector< UDQInput > input() const
Retrieve unprocessed input objects for all UDQs.
Definition UDQConfig.cpp:524
void required_summary(std::unordered_set< std::string > &summary_keys) const
Export all summary vectors needed to compute values for the current collection of user defined quanti...
Definition UDQConfig.cpp:671
bool operator==(const UDQConfig &config) const
Equality predicate.
Definition UDQConfig.cpp:656
const UDQFunctionTable & function_table() const
Retrieve run's active UDQ function table.
Definition UDQConfig.cpp:646
Definition UDQDefine.hpp:51
Definition UDQFunctionTable.hpp:32
Definition UDQInput.hpp:84
Definition UDQParams.hpp:31
Definition UDQState.hpp:40
Definition UDT.hpp:30
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30
Definition state.hpp:57