Inter-region flow accumulation maps for all region definition arrays. More...

#include <InterRegFlows.hpp>

Classes

struct  SingleRegion
 

Public Types

using Cell = InterRegFlowMapSingleFIP::Cell
 Characteristics of a cell from a simulation grid. More...
 

Public Member Functions

 InterRegFlowMap ()=default
 Default constructor. More...
 
 InterRegFlowMap (const std::size_t numCells, const std::vector< SingleRegion > &regions, const std::size_t declaredMaxRegID=0)
 
 InterRegFlowMap (const InterRegFlowMap &rhs)=default
 
 InterRegFlowMap (InterRegFlowMap &&rhs) noexcept=default
 
InterRegFlowMapoperator= (const InterRegFlowMap &rhs)=default
 
InterRegFlowMapoperator= (InterRegFlowMap &&rhs) noexcept=default
 
void addConnection (const Cell &source, const Cell &destination, const data::InterRegFlowMap::FlowRates &rates)
 
void compress ()
 
void clear ()
 Clear all internal buffers, but preserve allocated capacity. More...
 
const std::vector< std::string > & names () const
 
std::vector< data::InterRegFlowMap > getInterRegFlows () const
 
std::vector< std::size_t > getLocalMaxRegionID () const
 Retrieve maximum FIP region ID on local MPI rank. More...
 
bool assignGlobalMaxRegionID (const std::vector< std::size_t > &regID)
 
bool readIsConsistent () const
 Whether or not previous read() operation succeeded. More...
 
bool wantInterRegflowSummary () const
 
template<class MessageBufferType >
void write (MessageBufferType &buffer) const
 
template<class MessageBufferType >
void read (MessageBufferType &buffer)
 

Static Public Member Functions

static InterRegFlowMap createMapFromNames (std::vector< std::string > names)
 

Detailed Description

Inter-region flow accumulation maps for all region definition arrays.

Member Typedef Documentation

◆ Cell

Characteristics of a cell from a simulation grid.

Constructor & Destructor Documentation

◆ InterRegFlowMap() [1/4]

Opm::InterRegFlowMap::InterRegFlowMap ( )
default

Default constructor.

◆ InterRegFlowMap() [2/4]

Opm::InterRegFlowMap::InterRegFlowMap ( const std::size_t  numCells,
const std::vector< SingleRegion > &  regions,
const std::size_t  declaredMaxRegID = 0 
)
explicit

Constructor.

Parameters
[in]numCellsNumber of cells on local MPI rank, including overlap cells if applicable.
[in]regionsAll applicable region definition arrays.
[in]declaredMaxRegIDDeclared maximum region ID in the run-typically from the TABDIMS and/or REGDIMS keywords. Used for sizing internal data structures if greater than zero.

◆ InterRegFlowMap() [3/4]

Opm::InterRegFlowMap::InterRegFlowMap ( const InterRegFlowMap rhs)
default

◆ InterRegFlowMap() [4/4]

Opm::InterRegFlowMap::InterRegFlowMap ( InterRegFlowMap &&  rhs)
defaultnoexcept

Member Function Documentation

◆ addConnection()

void Opm::InterRegFlowMap::addConnection ( const Cell source,
const Cell destination,
const data::InterRegFlowMap::FlowRates &  rates 
)

Add flow rate connection between regions for all region definitions.

Parameters
[in]sourceCell from which the flow nominally originates.
[in]destinationCell into which flow nominally goes.
[in]ratesFlow rates associated to single connection.

If both cells are in the same region, or if neither cell is interior to this MPI rank, then this function does nothing. If one cell is interior to this MPI rank and the other isn't, then this function will include the flow rate contribution if and only if the cell with the smallest associate region ID is interior to this MPI rank.

Referenced by Opm::OutputBlackOilModule< TypeTag >::processFluxes().

◆ assignGlobalMaxRegionID()

bool Opm::InterRegFlowMap::assignGlobalMaxRegionID ( const std::vector< std::size_t > &  regID)

Assign maximum FIP region ID across all MPI ranks.

Fails if global maximum is smaller than local maximum region ID.

Parameters
[in]regIDGlobal maximum FIP region ID for this FIP definition array across all MPI ranks.
Returns
Whether or not assignment succeeded.

◆ clear()

void Opm::InterRegFlowMap::clear ( )

Clear all internal buffers, but preserve allocated capacity.

Referenced by Opm::OutputBlackOilModule< TypeTag >::initializeFluxData().

◆ compress()

void Opm::InterRegFlowMap::compress ( )

Form CSR adjacency matrix representation of input graph from connections established in previous calls to addConnection().

Number of rows in the CSR representation is the maximum FIP region ID.

Referenced by Opm::EclWriter< TypeTag >::evalSummaryState(), and Opm::OutputBlackOilModule< TypeTag >::finalizeFluxData().

◆ createMapFromNames()

static InterRegFlowMap Opm::InterRegFlowMap::createMapFromNames ( std::vector< std::string >  names)
static

Special purpose constructor for global object being collected on the I/O rank.

Only knows about the FIP region set names.

Parameters
[in]namesSorted sequence of FIP region names.

◆ getInterRegFlows()

std::vector< data::InterRegFlowMap > Opm::InterRegFlowMap::getInterRegFlows ( ) const

Get read-only access to the underlying CSR representation.

Mostly intended for summary output purposes.

◆ getLocalMaxRegionID()

std::vector< std::size_t > Opm::InterRegFlowMap::getLocalMaxRegionID ( ) const

Retrieve maximum FIP region ID on local MPI rank.

◆ names()

const std::vector< std::string > & Opm::InterRegFlowMap::names ( ) const

Names of all applicable region definition arrays.

Mostly intended for summary output purposes.

Referenced by read().

◆ operator=() [1/2]

InterRegFlowMap & Opm::InterRegFlowMap::operator= ( const InterRegFlowMap rhs)
default

◆ operator=() [2/2]

InterRegFlowMap & Opm::InterRegFlowMap::operator= ( InterRegFlowMap &&  rhs)
defaultnoexcept

◆ read()

template<class MessageBufferType >
void Opm::InterRegFlowMap::read ( MessageBufferType buffer)
inline

Reconstitute internal object representation from MPI message buffer

This object (

*this

) is not usable in subsequent calls to

void addConnection(const Cell &source, const Cell &destination, const data::InterRegFlowMap::FlowRates &rates)

following a call to member function

read()
void read(MessageBufferType &buffer)
Definition: InterRegFlows.hpp:323

.

Template Parameters
MessageBufferTypeLinear MPI message buffer. API should be similar to Dune::MessageBufferIF
Parameters
[in,out]bufferLinear MPI message buffer instance. Function reads a partially linearised representation of
*this
from the buffer contents and advances the buffer's read position.

References names().

Referenced by Opm::PackUnpackInterRegFlows::unpack().

◆ readIsConsistent()

bool Opm::InterRegFlowMap::readIsConsistent ( ) const

Whether or not previous read() operation succeeded.

◆ wantInterRegflowSummary()

bool Opm::InterRegFlowMap::wantInterRegflowSummary ( ) const
inline

◆ write()

template<class MessageBufferType >
void Opm::InterRegFlowMap::write ( MessageBufferType buffer) const
inline

Serialise internal representation to MPI message buffer

Template Parameters
MessageBufferTypeLinear MPI message buffer. API should be similar to Dune::MessageBufferIF
Parameters
[in,out]bufferLinear MPI message buffer instance. Function appends a partially linearised representation of
*this
to the buffer contents.

Referenced by Opm::PackUnpackInterRegFlows::pack().


The documentation for this class was generated from the following file: