My Project
Loading...
Searching...
No Matches
UDQSet.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 UDQSET_HPP
21#define UDQSET_HPP
22
23#include <opm/input/eclipse/Schedule/UDQ/UDQEnums.hpp>
24
25#include <cstddef>
26#include <optional>
27#include <stdexcept>
28#include <string>
29#include <unordered_map>
30#include <vector>
31
32namespace Opm {
33 class RegionSetMatchResult;
34 class SegmentSet;
35} // namespace Opm
36
37namespace Opm {
38
40{
41public:
42 UDQScalar() = default;
43
51 explicit UDQScalar(const double value, const std::size_t num = 0);
52
62 explicit UDQScalar(const std::string& wgname, const std::size_t num = 0);
63
70 void operator+=(const UDQScalar& rhs);
71
78 void operator+=(double rhs);
79
86 void operator*=(const UDQScalar& rhs);
87
94 void operator*=(double rhs);
95
102 void operator/=(const UDQScalar& rhs);
103
110 void operator/=(double rhs);
111
118 void operator-=(const UDQScalar& rhs);
119
126 void operator-=(double rhs);
127
131 operator bool() const;
132
137 void assign(const std::optional<double>& value);
138
143 void assign(double value);
144
146 bool defined() const;
147
151 double get() const;
152
156 const std::optional<double>& value() const { return this->m_value; }
157
159 const std::string& wgname() const { return this->m_wgname; }
160
165 std::size_t number() const { return this->m_num; }
166
171 bool operator==(const UDQScalar& other) const;
172
173public:
175 std::optional<double> m_value{};
176
178 std::string m_wgname{};
179
182 std::size_t m_num = 0;
183};
184
185
187{
188public:
189 // Connections and segments.
191 {
192 std::string name{};
193 std::vector<std::size_t> numbers{};
194
195 bool operator==(const EnumeratedItems& rhs) const;
196 static EnumeratedItems serializationTestObject();
197
198 template <class Serializer>
199 void serializeOp(Serializer& serializer)
200 {
201 serializer(this->name);
202 serializer(this->numbers);
203 }
204 };
205
206 static std::vector<EnumeratedItems>
207 enumerateItems(const RegionSetMatchResult& regionSet);
208
209 static std::vector<EnumeratedItems>
210 enumerateItems(const SegmentSet& segmentSet);
211
217 UDQSet(const std::string& name, UDQVarType var_type);
218
227 UDQSet(const std::string& name, UDQVarType var_type,
228 const std::vector<std::string>& wgnames);
229
240 UDQSet(const std::string& name, UDQVarType var_type,
241 const std::vector<EnumeratedItems>& items);
242
247 UDQSet(const std::string& name, UDQVarType var_type, std::size_t size);
248
253 UDQSet(const std::string& name, std::size_t size);
254
262 static UDQSet scalar(const std::string& name,
263 const std::optional<double>& scalar_value);
264
271 static UDQSet scalar(const std::string& name, double value);
272
276 static UDQSet empty(const std::string& name);
277
283 static UDQSet wells(const std::string& name,
284 const std::vector<std::string>& wells);
285
295 static UDQSet wells(const std::string& name,
296 const std::vector<std::string>& wells,
297 double scalar_value);
298
304 static UDQSet groups(const std::string& name,
305 const std::vector<std::string>& groups);
306
316 static UDQSet groups(const std::string& name,
317 const std::vector<std::string>& groups,
318 double scalar_value);
319
327 static UDQSet field(const std::string& name, double scalar_value);
328
329 static UDQSet segments(const std::string& name,
330 const std::vector<EnumeratedItems>& segments);
331 static UDQSet segments(const std::string& name,
332 const std::vector<EnumeratedItems>& segments,
333 const double scalar_value);
334
335 static UDQSet regions(const std::string& name,
336 const std::vector<EnumeratedItems>& regSetColl);
337 static UDQSet regions(const std::string& name,
338 const std::vector<EnumeratedItems>& regSetColl,
339 const double scalar_value);
340
345 void assign(const std::optional<double>& value);
346
354 void assign(std::size_t index, const std::optional<double>& value);
355
362 void assign(const std::string& wgname, const std::optional<double>& value);
363
373 void assign(const std::string& wgname, std::size_t number, const std::optional<double>& value);
374
379 void assign(double value);
380
388 void assign(std::size_t index, double value);
389
396 void assign(const std::string& wgname, double value);
397
399 bool has(const std::string& name) const;
400
403 std::size_t size() const;
404
411 void operator+=(const UDQSet& rhs);
412
416 void operator+=(double rhs);
417
424 void operator-=(const UDQSet& rhs);
425
429 void operator-=(double rhs);
430
437 void operator*=(const UDQSet& rhs);
438
442 void operator*=(double rhs);
443
450 void operator/=(const UDQSet& rhs);
451
455 void operator/=(double rhs);
456
462 const UDQScalar& operator[](std::size_t index) const;
463
469 const UDQScalar& operator[](const std::string& wgname) const;
470
481 const UDQScalar& operator()(const std::string& well, const std::size_t item) const;
482
484 std::vector<UDQScalar>::const_iterator begin() const;
485
487 std::vector<UDQScalar>::const_iterator end() const;
488
490 std::vector<std::string> wgnames() const;
491
493 std::vector<double> defined_values() const;
494
496 std::size_t defined_size() const;
497
499 const std::string& name() const;
500
502 void name(const std::string& name);
503
506 UDQVarType var_type() const;
507
511 bool operator==(const UDQSet& other) const;
512
513private:
515 std::string m_name;
516
518 UDQVarType m_var_type = UDQVarType::NONE;
519
521 std::vector<UDQScalar> values;
522
524 UDQSet() = default;
525};
526
527
528UDQScalar operator+(const UDQScalar& lhs, const UDQScalar& rhs);
529UDQScalar operator+(const UDQScalar& lhs, double rhs);
530UDQScalar operator+(double lhs, const UDQScalar& rhs);
531
532UDQScalar operator-(const UDQScalar& lhs, const UDQScalar& rhs);
533UDQScalar operator-(const UDQScalar& lhs, double rhs);
534UDQScalar operator-(double lhs, const UDQScalar& rhs);
535
536UDQScalar operator*(const UDQScalar& lhs, const UDQScalar& rhs);
537UDQScalar operator*(const UDQScalar& lhs, double rhs);
538UDQScalar operator*(double lhs, const UDQScalar& rhs);
539
540UDQScalar operator/(const UDQScalar& lhs, const UDQScalar& rhs);
541UDQScalar operator/(const UDQScalar& lhs, double rhs);
542UDQScalar operator/(double lhs, const UDQScalar& rhs);
543
544UDQSet operator+(const UDQSet& lhs, const UDQSet& rhs);
545UDQSet operator+(const UDQSet& lhs, double rhs);
546UDQSet operator+(double lhs, const UDQSet& rhs);
547
548UDQSet operator-(const UDQSet& lhs, const UDQSet& rhs);
549UDQSet operator-(const UDQSet& lhs, double rhs);
550UDQSet operator-(double lhs, const UDQSet& rhs);
551
552UDQSet operator*(const UDQSet& lhs, const UDQSet& rhs);
553UDQSet operator*(const UDQSet& lhs, double rhs);
554UDQSet operator*(double lhs, const UDQSet& rhs);
555
556UDQSet operator/(const UDQSet& lhs, const UDQSet& rhs);
557UDQSet operator/(const UDQSet& lhs, double rhs);
558UDQSet operator/(double lhs, const UDQSet&rhs);
559
560} // namespace Opm
561
562#endif // UDQSET_HPP
Result Set From RegionSetMatcher's Matching Process.
Definition RegionSetMatcher.hpp:44
Class for (de-)serializing.
Definition Serializer.hpp:84
Definition UDQSet.hpp:40
void operator/=(const UDQScalar &rhs)
Divide this UDQ scalar by other.
Definition UDQSet.cpp:100
void operator-=(const UDQScalar &rhs)
Subtract other UDQ scalar from this.
Definition UDQSet.cpp:88
const std::optional< double > & value() const
Retrive contained numeric value.
Definition UDQSet.hpp:156
void assign(const std::optional< double > &value)
Assign numeric value to this UDQ scalar.
Definition UDQSet.cpp:68
std::size_t number() const
Retrive numbered item, typically segment or connection, to which this scalar is associated.
Definition UDQSet.hpp:165
void operator+=(const UDQScalar &rhs)
Add other UDQ scalar to this.
Definition UDQSet.cpp:112
std::optional< double > m_value
Scalar value.
Definition UDQSet.hpp:175
const std::string & wgname() const
Retrive named well/group to which this scalar is associated.
Definition UDQSet.hpp:159
void operator*=(const UDQScalar &rhs)
Multiply UDQ scalar into this.
Definition UDQSet.cpp:124
bool operator==(const UDQScalar &other) const
Equality predicate.
Definition UDQSet.cpp:141
std::string m_wgname
Associated well/group name.
Definition UDQSet.hpp:178
bool defined() const
Predicate for whether or not this UDQ scalar has a defined value.
Definition UDQSet.cpp:51
std::size_t m_num
Numbered item.
Definition UDQSet.hpp:182
double get() const
Retrive contained numeric value.
Definition UDQSet.cpp:56
Definition UDQSet.hpp:187
std::vector< UDQScalar >::const_iterator begin() const
Range-for traversal support (beginning of range)
Definition UDQSet.cpp:539
static UDQSet empty(const std::string &name)
Form an empty UDQ set.
Definition UDQSet.cpp:234
const UDQScalar & operator()(const std::string &well, const std::size_t item) const
Access individual UDQ scalar assiociated to particular named well and numbered sub-entity of that nam...
Definition UDQSet.cpp:520
void operator-=(const UDQSet &rhs)
Subtract UDQ set from this.
Definition UDQSet.cpp:433
static UDQSet field(const std::string &name, double scalar_value)
Form a UDQ set at the field level.
Definition UDQSet.cpp:239
void operator+=(const UDQSet &rhs)
Add other UDQ set into this.
Definition UDQSet.cpp:415
const UDQScalar & operator[](std::size_t index) const
Access individual UDQ scalar at particular index in UDQ set.
Definition UDQSet.cpp:495
const std::string & name() const
Retrive the name of this UDQ set.
Definition UDQSet.cpp:164
std::vector< double > defined_values() const
Retrive the UDQ set's defined values only.
Definition UDQSet.cpp:473
static UDQSet wells(const std::string &name, const std::vector< std::string > &wells)
Form a UDQ set pertaining to a set of named wells.
Definition UDQSet.cpp:246
std::vector< std::string > wgnames() const
Retrive names of entities associate to this UDQ set.
Definition UDQSet.cpp:402
std::vector< UDQScalar >::const_iterator end() const
Range-for traversal support (one past end of range)
Definition UDQSet.cpp:544
bool operator==(const UDQSet &other) const
Equality comparison operator.
Definition UDQSet.cpp:784
void operator*=(const UDQSet &rhs)
Multiply other UDQ set into this.
Definition UDQSet.cpp:437
static UDQSet groups(const std::string &name, const std::vector< std::string > &groups)
Form a UDQ set pertaining to a set of named groups.
Definition UDQSet.cpp:261
UDQVarType var_type() const
Retrive the variable type of this UDQ set (e.g., well, group, field, segment &c).
Definition UDQSet.cpp:397
bool has(const std::string &name) const
Predicate for whether or not named UDQ element exists.
Definition UDQSet.cpp:306
std::size_t defined_size() const
Retrive the UDQ set's number of defined values.
Definition UDQSet.cpp:486
std::size_t size() const
Number of elements in UDQ set.
Definition UDQSet.cpp:315
void assign(const std::optional< double > &value)
Assign value to every element of the UDQ set.
Definition UDQSet.cpp:386
static UDQSet scalar(const std::string &name, const std::optional< double > &scalar_value)
Form a UDQ set pertaining to a single scalar value.
Definition UDQSet.cpp:226
void operator/=(const UDQSet &rhs)
Divide this UDQ set by other UDQ set.
Definition UDQSet.cpp:455
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30
Definition UDQSet.hpp:191