InterRegFlows.hpp
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*
4 Copyright 2022 Equinor AS
5
6 This file is part of the Open Porous Media project (OPM).
7
8 OPM is free software: you can redistribute it and/or modify
9 it under the terms of the GNU General Public License as published by
10 the Free Software Foundation, either version 3 of the License, or
11 (at your option) any later version.
12
13 OPM is distributed in the hope that it will be useful,
14 but WITHOUT ANY WARRANTY; without even the implied warranty of
15 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16 GNU General Public License for more details.
17
18 You should have received a copy of the GNU General Public License
19 along with OPM. If not, see <http://www.gnu.org/licenses/>.
20*/
21
22#ifndef OPM_INTERREG_FLOWS_MODULE_HPP
23#define OPM_INTERREG_FLOWS_MODULE_HPP
24
25#include <opm/output/data/InterRegFlowMap.hpp>
26
27#include <algorithm>
28#include <cstddef>
29#include <string>
30#include <utility>
31#include <vector>
32
38
39namespace Opm {
40
41 class InterRegFlowMap;
42
47 {
48 public:
50 struct Cell {
52 int activeIndex{-1};
53
56
58 bool isInterior{true};
59 };
60
61 friend class InterRegFlowMap;
62
66 explicit InterRegFlowMapSingleFIP(const std::vector<int>& region);
67
82 void addConnection(const Cell& source,
83 const Cell& destination,
84 const data::InterRegFlowMap::FlowRates& rates);
85
91 void compress();
92
94 void clear();
95
99 const data::InterRegFlowMap& getInterRegFlows() const;
100
102 std::size_t getLocalMaxRegionID() const;
103
112 bool assignGlobalMaxRegionID(const std::size_t regID);
113
122 template <class MessageBufferType>
123 void write(MessageBufferType& buffer) const
124 {
125 // Note: this->region_ is *intentionally* omitted here.
126 buffer.write(this->maxGlobalRegionID_);
127 this->iregFlow_.write(buffer);
128 }
129
144 template <class MessageBufferType>
146 {
147 auto otherMaxRegionID = 0 * this->maxGlobalRegionID_;
148 buffer.read(otherMaxRegionID);
149 this->maxGlobalRegionID_ = std::max(this->maxGlobalRegionID_, otherMaxRegionID);
150
151 this->iregFlow_.read(buffer);
152 this->isReadFromStream_ = true;
153 }
154
155 private:
157 std::vector<int> region_{};
158
160 std::size_t maxLocalRegionID_{0};
161
164 std::size_t maxGlobalRegionID_{0};
165
167 data::InterRegFlowMap iregFlow_{};
168
171 bool isReadFromStream_{false};
172
174 InterRegFlowMapSingleFIP() = default;
175 };
176
179 {
180 public:
186 std::string name;
187
189 std::reference_wrapper<const std::vector<int>> definition;
190 };
191
194
196 InterRegFlowMap() = default;
197
204 static InterRegFlowMap
205 createMapFromNames(std::vector<std::string> names);
206
217 explicit InterRegFlowMap(const std::size_t numCells,
218 const std::vector<SingleRegion>& regions,
219 const std::size_t declaredMaxRegID = 0);
220
221 InterRegFlowMap(const InterRegFlowMap& rhs) = default;
222 InterRegFlowMap(InterRegFlowMap&& rhs) noexcept = default;
223
225 InterRegFlowMap& operator=(InterRegFlowMap&& rhs) noexcept = default;
226
242 void addConnection(const Cell& source,
243 const Cell& destination,
244 const data::InterRegFlowMap::FlowRates& rates);
245
251 void compress();
252
254 void clear();
255
259 const std::vector<std::string>& names() const;
260
264 std::vector<data::InterRegFlowMap> getInterRegFlows() const;
265
267 std::vector<std::size_t> getLocalMaxRegionID() const;
268
277 bool assignGlobalMaxRegionID(const std::vector<std::size_t>& regID);
278
280 bool readIsConsistent() const;
281
283 {
284 return !this->regionMaps_.empty();
285 }
286
295 template <class MessageBufferType>
296 void write(MessageBufferType& buffer) const
297 {
298 buffer.write(this->names_.size());
299 for (const auto& name : this->names_) {
300 buffer.write(name);
301 }
302
303 for (const auto& map : this->regionMaps_) {
304 map.write(buffer);
305 }
306 }
307
322 template <class MessageBufferType>
324 {
325 const auto names = this->readNames(buffer);
326
327 if (names == this->names_) {
328 // Input stream holds the maps in expected order (common
329 // case). Read the maps and merge with internal values.
330 for (auto& map : this->regionMaps_) {
331 map.read(buffer);
332 }
333 }
334 else {
335 // Input stream does not hold the maps in expected order (or
336 // different number of maps). Unexpected. Read the values
337 // from the input stream, but do not merge with internal
338 // values.
339 auto map = InterRegFlowMapSingleFIP {
340 std::vector<int>(this->numCells_, 1)
341 };
342
343 const auto numMaps = names.size();
344 for (auto n = 0*numMaps; n < numMaps; ++n) {
345 map.read(buffer);
346 }
347
348 this->readIsConsistent_ = false;
349 }
350 }
351
352 private:
355 std::vector<InterRegFlowMapSingleFIP> regionMaps_{};
356
359 std::vector<std::string> names_;
360
363 std::size_t numCells_{0};
364
367 bool readIsConsistent_{true};
368
372 template <class MessageBufferType>
373 std::vector<std::string> readNames(MessageBufferType& buffer) const
374 {
375 auto numNames = 0 * this->names_.size();
376 buffer.read(numNames);
377
378 auto names = std::vector<std::string>(numNames);
379 for (auto name = 0*numNames; name < numNames; ++name) {
380 buffer.read(names[name]);
381 }
382
383 return names;
384 }
385 };
386} // namespace Opm
387
388#endif // OPM_INTERREG_FLOWS_MODULE_HPP
Inter-region flow accumulation maps for all region definition arrays.
Definition: InterRegFlows.hpp:179
void read(MessageBufferType &buffer)
Definition: InterRegFlows.hpp:323
InterRegFlowMap(const InterRegFlowMap &rhs)=default
InterRegFlowMap & operator=(const InterRegFlowMap &rhs)=default
void write(MessageBufferType &buffer) const
Definition: InterRegFlows.hpp:296
InterRegFlowMap(InterRegFlowMap &&rhs) noexcept=default
std::vector< data::InterRegFlowMap > getInterRegFlows() const
InterRegFlowMap()=default
Default constructor.
bool wantInterRegflowSummary() const
Definition: InterRegFlows.hpp:282
const std::vector< std::string > & names() const
std::vector< std::size_t > getLocalMaxRegionID() const
Retrieve maximum FIP region ID on local MPI rank.
void addConnection(const Cell &source, const Cell &destination, const data::InterRegFlowMap::FlowRates &rates)
void clear()
Clear all internal buffers, but preserve allocated capacity.
InterRegFlowMap & operator=(InterRegFlowMap &&rhs) noexcept=default
bool readIsConsistent() const
Whether or not previous read() operation succeeded.
static InterRegFlowMap createMapFromNames(std::vector< std::string > names)
InterRegFlowMap(const std::size_t numCells, const std::vector< SingleRegion > &regions, const std::size_t declaredMaxRegID=0)
bool assignGlobalMaxRegionID(const std::vector< std::size_t > &regID)
Definition: InterRegFlows.hpp:47
std::size_t getLocalMaxRegionID() const
Retrieve maximum FIP region ID on local MPI rank.
void read(MessageBufferType &buffer)
Definition: InterRegFlows.hpp:145
const data::InterRegFlowMap & getInterRegFlows() const
InterRegFlowMapSingleFIP(const std::vector< int > &region)
void clear()
Clear all internal buffers, but preserve allocated capacity.
void write(MessageBufferType &buffer) const
Definition: InterRegFlows.hpp:123
bool assignGlobalMaxRegionID(const std::size_t regID)
void addConnection(const Cell &source, const Cell &destination, const data::InterRegFlowMap::FlowRates &rates)
Definition: BlackoilPhases.hpp:27
typename P2PCommunicatorType::MessageBufferType MessageBufferType
Definition: CollectDataOnIORank_impl.hpp:83
Definition: InterRegFlows.hpp:184
std::reference_wrapper< const std::vector< int > > definition
Region definition array.
Definition: InterRegFlows.hpp:189
std::string name
Region definition array name.
Definition: InterRegFlows.hpp:186
Minimal characteristics of a cell from a simulation grid.
Definition: InterRegFlows.hpp:50
int cartesianIndex
Cell's global cell ID.
Definition: InterRegFlows.hpp:55
int activeIndex
Cell's active index on local rank.
Definition: InterRegFlows.hpp:52
bool isInterior
Whether or not cell is interior to local rank.
Definition: InterRegFlows.hpp:58