Opm::ParallelNLDDPartitioningZoltan Class Reference

#include <ParallelNLDDPartitioningZoltan.hpp>

Public Types

using GlobalCellID = std::function< int(int)>
 
using ZoltanParamMap = std::map< std::string, std::string >
 Collection of parameters to control the partitioning procedure. More...
 

Public Member Functions

 ParallelNLDDPartitioningZoltan (const Parallel::Communication comm, const std::size_t numElements, const GlobalCellID &globalCell)
 
void registerConnection (std::size_t c1, std::size_t c2)
 
void forceSameDomain (std::vector< int > &&cells)
 
std::vector< int > partitionElements (const ZoltanParamMap &params) const
 

Detailed Description

Partition connectivity graph into non-overlapping domains using the Zoltan graph partitioning software package. Primarily intended for use in NLDD-based non-linear solves.

Member Typedef Documentation

◆ GlobalCellID

using Opm::ParallelNLDDPartitioningZoltan::GlobalCellID = std::function<int(int)>

Callback type for mapping a vertex/cell ID to a globally unique ID.

◆ ZoltanParamMap

using Opm::ParallelNLDDPartitioningZoltan::ZoltanParamMap = std::map<std::string, std::string>

Collection of parameters to control the partitioning procedure.

Constructor & Destructor Documentation

◆ ParallelNLDDPartitioningZoltan()

Opm::ParallelNLDDPartitioningZoltan::ParallelNLDDPartitioningZoltan ( const Parallel::Communication  comm,
const std::size_t  numElements,
const GlobalCellID globalCell 
)
inlineexplicit

Constructor.

Parameters
[in]commMPI communication object. Needed by Zoltan.
[in]numElementsNumber of potential vertices in connectivity graph. Typically the total number of cells on the current rank, i.e., both owned cells and overlap cells.
[in]globalCellCallback for mapping (local) vertex IDs to globally unique vertex IDs.

Member Function Documentation

◆ forceSameDomain()

void Opm::ParallelNLDDPartitioningZoltan::forceSameDomain ( std::vector< int > &&  cells)
inline

Force collection of cells to be in same result domain.

Mostly as a means to ensuring wells do not intersect multiple domains/blocks.

Parameters
[in]cellsCell collection. Typically those cells which are intersected by a single well.

◆ partitionElements()

std::vector< int > Opm::ParallelNLDDPartitioningZoltan::partitionElements ( const ZoltanParamMap params) const

Partition connectivity graph using Zoltan graph partitioning package.

Honours any prescribed requirement that certain cells be placed in a single domain/block.

Parameters
[in]paramsParameters for Zoltan. Override default settings.
Returns
Partition vector of size numElements. Reachable vertices/cells are partitioned into N blocks numbered 0..N-1. Unreachable vertices get a block ID of -1.

◆ registerConnection()

void Opm::ParallelNLDDPartitioningZoltan::registerConnection ( std::size_t  c1,
std::size_t  c2 
)
inline

Insert directed graph edge between two vertices.

Parameters
[in]c1Source vertex.
]in]c2 Sink/destination vertex.

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