22#ifndef OPM_INTERREG_FLOWS_MODULE_HPP
23#define OPM_INTERREG_FLOWS_MODULE_HPP
25#include <opm/output/data/InterRegFlowMap.hpp>
41 class InterRegFlowMap;
83 const Cell& destination,
84 const data::InterRegFlowMap::FlowRates& rates);
122 template <
class MessageBufferType>
123 void write(MessageBufferType& buffer)
const
126 buffer.write(this->maxGlobalRegionID_);
127 this->iregFlow_.write(buffer);
144 template <
class MessageBufferType>
145 void read(MessageBufferType& buffer)
147 auto otherMaxRegionID = 0 * this->maxGlobalRegionID_;
148 buffer.read(otherMaxRegionID);
149 this->maxGlobalRegionID_ = std::max(this->maxGlobalRegionID_, otherMaxRegionID);
151 this->iregFlow_.read(buffer);
152 this->isReadFromStream_ =
true;
157 std::vector<int> region_{};
160 std::size_t maxLocalRegionID_{0};
164 std::size_t maxGlobalRegionID_{0};
167 data::InterRegFlowMap iregFlow_{};
171 bool isReadFromStream_{
false};
218 const std::vector<SingleRegion>& regions,
219 const std::size_t declaredMaxRegID = 0);
243 const Cell& destination,
244 const data::InterRegFlowMap::FlowRates& rates);
259 const std::vector<std::string>&
names()
const;
282 bool wantInterRegflowSummary()
const
284 return !this->regionMaps_.empty();
295 template <
class MessageBufferType>
296 void write(MessageBufferType& buffer)
const
298 buffer.write(this->names_.size());
299 for (
const auto& name : this->names_) {
303 for (
const auto& map : this->regionMaps_) {
322 template <
class MessageBufferType>
323 void read(MessageBufferType& buffer)
325 const auto names = this->readNames(buffer);
327 if (
names == this->names_) {
330 for (
auto& map : this->regionMaps_) {
340 std::vector<int>(this->numCells_, 1)
343 const auto numMaps =
names.size();
344 for (
auto n = 0*numMaps; n < numMaps; ++n) {
348 this->readIsConsistent_ =
false;
355 std::vector<InterRegFlowMapSingleFIP> regionMaps_{};
359 std::vector<std::string> names_;
363 std::size_t numCells_{0};
367 bool readIsConsistent_{
true};
372 template <
class MessageBufferType>
373 std::vector<std::string> readNames(MessageBufferType& buffer)
const
375 auto numNames = 0 * this->names_.size();
376 buffer.read(numNames);
378 auto names = std::vector<std::string>(numNames);
379 for (
auto name = 0*numNames; name < numNames; ++name) {
380 buffer.read(
names[name]);
Form CSR adjacency matrix representation of inter-region flow rate graph provided as a list of connec...
Definition InterRegFlows.hpp:47
std::size_t getLocalMaxRegionID() const
Retrieve maximum FIP region ID on local MPI rank.
Definition InterRegFlows.cpp:95
void read(MessageBufferType &buffer)
Reconstitute internal object representation from MPI message buffer.
Definition InterRegFlows.hpp:145
void compress()
Form CSR adjacency matrix representation of input graph from connections established in previous call...
Definition InterRegFlows.cpp:78
InterRegFlowMapSingleFIP(const std::vector< int > ®ion)
Constructor.
Definition InterRegFlows.cpp:32
void clear()
Clear all internal buffers, but preserve allocated capacity.
Definition InterRegFlows.cpp:83
const data::InterRegFlowMap & getInterRegFlows() const
Get read-only access to the underlying CSR representation.
Definition InterRegFlows.cpp:90
void write(MessageBufferType &buffer) const
Serialise internal representation to MPI message buffer.
Definition InterRegFlows.hpp:123
bool assignGlobalMaxRegionID(const std::size_t regID)
Assign maximum FIP region ID across all MPI ranks.
Definition InterRegFlows.cpp:102
void addConnection(const Cell &source, const Cell &destination, const data::InterRegFlowMap::FlowRates &rates)
Add flow rate connection between regions.
Definition InterRegFlows.cpp:47
Inter-region flow accumulation maps for all region definition arrays.
Definition InterRegFlows.hpp:179
void read(MessageBufferType &buffer)
Reconstitute internal object representation from MPI message buffer.
Definition InterRegFlows.hpp:323
static InterRegFlowMap createMapFromNames(std::vector< std::string > names)
Special purpose constructor for global object being collected on the I/O rank.
Definition InterRegFlows.cpp:122
void write(MessageBufferType &buffer) const
Serialise internal representation to MPI message buffer.
Definition InterRegFlows.hpp:296
void compress()
Form CSR adjacency matrix representation of input graph from connections established in previous call...
Definition InterRegFlows.cpp:165
InterRegFlowMap()=default
Default constructor.
const std::vector< std::string > & names() const
Names of all applicable region definition arrays.
Definition InterRegFlows.cpp:182
std::vector< std::size_t > getLocalMaxRegionID() const
Retrieve maximum FIP region ID on local MPI rank.
Definition InterRegFlows.cpp:201
void addConnection(const Cell &source, const Cell &destination, const data::InterRegFlowMap::FlowRates &rates)
Add flow rate connection between regions for all region definitions.
Definition InterRegFlows.cpp:156
void clear()
Clear all internal buffers, but preserve allocated capacity.
Definition InterRegFlows.cpp:172
bool readIsConsistent() const
Whether or not previous read() operation succeeded.
Definition InterRegFlows.cpp:234
std::vector< data::InterRegFlowMap > getInterRegFlows() const
Get read-only access to the underlying CSR representation.
Definition InterRegFlows.cpp:188
bool assignGlobalMaxRegionID(const std::vector< std::size_t > ®ID)
Assign maximum FIP region ID across all MPI ranks.
Definition InterRegFlows.cpp:215
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilboundaryratevector.hh:37
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
Minimal representation of a single named region defintion.
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