22#ifndef ECL_INTERREG_FLOWS_MODULE_HPP
23#define ECL_INTERREG_FLOWS_MODULE_HPP
25#include <opm/output/data/InterRegFlowMap.hpp>
41 class EclInterRegFlowMap;
84 const data::InterRegFlowMap::FlowRates& rates);
122 template <
class MessageBufferType>
126 buffer.write(this->maxGlobalRegionID_);
127 this->iregFlow_.write(
buffer);
144 template <
class MessageBufferType>
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};
214 const std::vector<SingleRegion>&
regions);
239 const data::InterRegFlowMap::FlowRates& rates);
254 const std::vector<std::string>&
names()
const;
277 bool wantInterRegflowSummary()
const
279 return !this->regionMaps_.empty();
290 template <
class MessageBufferType>
293 buffer.write(this->names_.size());
294 for (
const auto& name : this->names_) {
298 for (
const auto&
map : this->regionMaps_) {
317 template <
class MessageBufferType>
322 if (
names == this->names_) {
325 for (
auto&
map : this->regionMaps_) {
335 std::vector<int>(this->numCells_, 1)
343 this->readIsConsistent_ =
false;
350 std::vector<EclInterRegFlowMapSingleFIP> regionMaps_{};
354 std::vector<std::string> names_;
358 std::size_t numCells_{0};
362 bool readIsConsistent_{
true};
367 template <
class MessageBufferType>
368 std::vector<std::string> readNames(MessageBufferType& buffer)
const
370 auto numNames = 0 * this->names_.size();
371 buffer.read(numNames);
373 auto names = std::vector<std::string>(numNames);
374 for (
auto name = 0*numNames; name < numNames; ++name) {
375 buffer.read(
names[name]);
Definition AquiferInterface.hpp:35
Form CSR adjacency matrix representation of inter-region flow rate graph provided as a list of connec...
Definition EclInterRegFlows.hpp:47
void clear()
Clear all internal buffers, but preserve allocated capacity.
Definition EclInterRegFlows.cpp:83
void addConnection(const Cell &source, const Cell &destination, const data::InterRegFlowMap::FlowRates &rates)
Add flow rate connection between regions.
Definition EclInterRegFlows.cpp:47
const data::InterRegFlowMap & getInterRegFlows() const
Get read-only access to the underlying CSR representation.
Definition EclInterRegFlows.cpp:90
bool assignGlobalMaxRegionID(const std::size_t regID)
Assign maximum FIP region ID across all MPI ranks.
Definition EclInterRegFlows.cpp:102
std::size_t getLocalMaxRegionID() const
Retrieve maximum FIP region ID on local MPI rank.
Definition EclInterRegFlows.cpp:95
void compress()
Form CSR adjacency matrix representation of input graph from connections established in previous call...
Definition EclInterRegFlows.cpp:78
EclInterRegFlowMapSingleFIP(const std::vector< int > ®ion)
Constructor.
Definition EclInterRegFlows.cpp:32
void write(MessageBufferType &buffer) const
Serialise internal representation to MPI message buffer.
Definition EclInterRegFlows.hpp:123
void read(MessageBufferType &buffer)
Reconstitute internal object representation from MPI message buffer.
Definition EclInterRegFlows.hpp:145
Inter-region flow accumulation maps for all region definition arrays.
Definition EclInterRegFlows.hpp:179
static EclInterRegFlowMap createMapFromNames(std::vector< std::string > names)
Special purpose constructor for global object being collected on the I/O rank.
Definition EclInterRegFlows.cpp:119
bool assignGlobalMaxRegionID(const std::vector< std::size_t > ®ID)
Assign maximum FIP region ID across all MPI ranks.
Definition EclInterRegFlows.cpp:205
std::vector< std::size_t > getLocalMaxRegionID() const
Retrieve maximum FIP region ID on local MPI rank.
Definition EclInterRegFlows.cpp:191
EclInterRegFlowMap()=default
Default constructor.
void read(MessageBufferType &buffer)
Reconstitute internal object representation from MPI message buffer.
Definition EclInterRegFlows.hpp:318
bool readIsConsistent() const
Whether or not previous read() operation succeeded.
Definition EclInterRegFlows.cpp:224
const std::vector< std::string > & names() const
Names of all applicable region definition arrays.
Definition EclInterRegFlows.cpp:172
void compress()
Form CSR adjacency matrix representation of input graph from connections established in previous call...
Definition EclInterRegFlows.cpp:155
std::vector< data::InterRegFlowMap > getInterRegFlows() const
Get read-only access to the underlying CSR representation.
Definition EclInterRegFlows.cpp:178
void addConnection(const Cell &source, const Cell &destination, const data::InterRegFlowMap::FlowRates &rates)
Add flow rate connection between regions for all region definitions.
Definition EclInterRegFlows.cpp:146
void write(MessageBufferType &buffer) const
Serialise internal representation to MPI message buffer.
Definition EclInterRegFlows.hpp:291
void clear()
Clear all internal buffers, but preserve allocated capacity.
Definition EclInterRegFlows.cpp:162
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition BlackoilPhases.hpp:27
Minimal characteristics of a cell from a simulation grid.
Definition EclInterRegFlows.hpp:50
int cartesianIndex
Cell's global cell ID.
Definition EclInterRegFlows.hpp:55
int activeIndex
Cell's active index on local rank.
Definition EclInterRegFlows.hpp:52
bool isInterior
Whether or not cell is interior to local rank.
Definition EclInterRegFlows.hpp:58
Minimal representation of a single named region defintion.
Definition EclInterRegFlows.hpp:184
std::string name
Region definition array name.
Definition EclInterRegFlows.hpp:186
std::reference_wrapper< const std::vector< int > > definition
Region definition array.
Definition EclInterRegFlows.hpp:189