SpinParser  1.0
Classes | Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
HMP::LoadManager Class Referenceabstract

Common load manager interface for both, the MPI server rank and slave ranks. More...

#include <LoadManager.hpp>

Inheritance diagram for HMP::LoadManager:
Inheritance graph
[legend]

Classes

struct  Chunk
 Definition of a small chunk of work. Used to communicate workload between different LoadManager instances. More...
 
struct  DataStack
 Concrete derivation of DataStackBase, managing an array of fundamental data types (or tuples thereof, if DataStackBase::typeMultiplicity is greater than one). More...
 
struct  DataStackBase
 Abstract base class for DataStack types. More...
 

Public Member Functions

virtual ~LoadManager ()
 Destroy the LoadManager object.
 
template<class StackT >
StackIdentifier addMasterStackExplicit (StackT *const data, const int size, const std::function< StackT(int)> &calculator, const int recommendedChunkSizeMultiple=1, const int recommendedChunksPerRank=10, const bool autoBroadcast=false)
 Create an explicit data stack and attach it to the LoadManager. More...
 
template<class StackT >
StackIdentifier addMasterStackImplicit (StackT *const data, const int size, const std::function< void(int)> &calculator, const int typeMultiplicity=1, const int recommendedChunkSizeMultiple=1, const int recommendedChunksPerRank=10, const bool autoBroadcast=false)
 Create an implicit data stack and attach it to the LoadManager. More...
 
template<class StackT >
StackIdentifier addSlaveStack (StackT *const data, const int size, const StackIdentifier master, const int typeMultiplicity=1)
 Create a slave data stack and attach it to the LoadManager. More...
 
template<class StackT >
StackIdentifier addPassiveStack (StackT *const data, const int size)
 Create a passive data stack and attach it to the LoadManager. More...
 
virtual void calculate (const StackIdentifier *stackIds, const int size)=0
 Calculate a list of stacks, where the stack identifiers are provided in list form. More...
 
void calculate (const std::initializer_list< StackIdentifier > &stackIds)
 Calculate a list of stacks, where the stack identifiers are provided in initializer list form. More...
 
void calculate (const StackIdentifier stackId)
 Calculate a single stack. More...
 
void calculateAll ()
 Calculate all stacks.
 
void broadcast (const StackIdentifier *stackIds, const int size)
 Broadcast a list of stacks, where the stack identifiers are provided in list form. More...
 
void broadcast (const std::initializer_list< StackIdentifier > &stackIds)
 Broadcast a list of stacks, where the stack identifiers are provided in initializer list form. More...
 
void broadcast (const StackIdentifier stackId)
 Broadcast a single stack. More...
 
void broadcastAll ()
 Broadcast all stacks.
 
virtual void printRuntimeStatistics () const
 Virtual implementation to print runtime statistics, including information on the efficiency of LoadManager instances running on different MPI ranks.
 

Protected Member Functions

 LoadManager (const int serverRank, const MPI_Comm communicator)
 Construct a new Load Manager object. More...
 
virtual StackIdentifier _registerStack (DataStackBase *stack)
 Register a DataStackBase with the LoadManager. By registering the stack, the LoadManager assumes responsibility to synchronize data between MPI ranks as necessary. More...
 
void _calculateChunk (const Chunk &chunk)
 Calculate the workload defined by a specific chunk. More...
 

Protected Attributes

std::vector< DataStackBase * > _stacks
 List of all registered stacks.
 
int _serverRank
 Designated MPI master rank.
 
int _rank
 MPI rank of the current LoadManager instance.
 
int _commSize
 MPI communicator size.
 
MPI_Comm _communicator
 MPI communicator to operate on.
 

Detailed Description

Common load manager interface for both, the MPI server rank and slave ranks.

This class provides an implementation of an automatic load balancing scheme for hybrid shared/distributed memory parallelization. Intra-node parallelization is implemented via OpenMP. Inter-node parallelization is implemented via MPI. A LoadManager instance is created on every MPI rank by calling the LoadManager::newLoadManager routine; This will create a LoadManagerMaster instance on the designated MPI master rank, and LoadManagerSlave instances on all remaining ranks. The work load is registered with the LoadManager in the form of DataStacks. Each data stack is associated with an array of numbers and with a functional description (referred to as calculators) on how to perform computations on the individual array entries. Different types of data stack exist with different calculator properties, see DataStackBase::StackType for more information. The LoadManager interface provides functions to create and register the different stack types with the LoadManager.

The role of the LoadManagerMaster is to break down the collective computation of all array entries into smaller chunks of work, to distribute them either locally or across multiple MPI ranks, and to retrieve the result on the master rank. In order to execute the calculators on a specific stack, the LoadManager provides the LoadManager::calculate interface. In order to make the results available also on other MPI ranks, the LoadManager::broadcast interface is provided. Various different types of data stacks exist, see the appropriate member functions for adding new stacks.

See also
LoadManager::newLoadManager
LoadManager::calculate
LoadManager::broadcast

Constructor & Destructor Documentation

◆ LoadManager()

HMP::LoadManager::LoadManager ( const int  serverRank,
const MPI_Comm  communicator 
)
inlineprotected

Construct a new Load Manager object.

Parameters
serverRankThe MPI rank to take on the master role.
communicatorThe MPI communicator to operate on.

Member Function Documentation

◆ _calculateChunk()

void HMP::LoadManager::_calculateChunk ( const Chunk chunk)
inlineprotected

Calculate the workload defined by a specific chunk.

Parameters
chunkProvides the workload definition.

◆ _registerStack()

virtual StackIdentifier HMP::LoadManager::_registerStack ( DataStackBase stack)
inlineprotectedvirtual

Register a DataStackBase with the LoadManager. By registering the stack, the LoadManager assumes responsibility to synchronize data between MPI ranks as necessary.

Parameters
stackThe stack to be registered.
Returns
StackIdentifier Unique identifier of the newly registered stack; Used e.g. for LoadManager::calculate calls.

Reimplemented in HMP::LoadManagerSlave, and HMP::LoadManagerMaster.

◆ addMasterStackExplicit()

template<class StackT >
StackIdentifier HMP::LoadManager::addMasterStackExplicit ( StackT *const  data,
const int  size,
const std::function< StackT(int)> &  calculator,
const int  recommendedChunkSizeMultiple = 1,
const int  recommendedChunksPerRank = 10,
const bool  autoBroadcast = false 
)
inline

Create an explicit data stack and attach it to the LoadManager.

Template Parameters
StackTThe fundamental data type of the DataStack to be created.
Parameters
dataData array on which the stack operates. The allocated size should be at least size * typeof(StackT).
sizeNumber of elements in the stack.
calculatorCalculator object associated with the stack.
recommendedChunkSizeMultipleSuggested quantization of elements per work chunk.
recommendedChunksPerRankSuggested number of work chunks per MPI rank.
autoBroadcastEnable or disable auto broadcast.
Returns
StackIdentifier Id of the newly generated stack as registered with the LoadManager.
See also
DataStackBase::StackType::Explicit
DataStack

◆ addMasterStackImplicit()

template<class StackT >
StackIdentifier HMP::LoadManager::addMasterStackImplicit ( StackT *const  data,
const int  size,
const std::function< void(int)> &  calculator,
const int  typeMultiplicity = 1,
const int  recommendedChunkSizeMultiple = 1,
const int  recommendedChunksPerRank = 10,
const bool  autoBroadcast = false 
)
inline

Create an implicit data stack and attach it to the LoadManager.

Template Parameters
StackTThe fundamental data type of the DataStack to be created.
Parameters
dataData array on which the stack operates. The allocated size should be at least size * typeMultiplicity * typeof(StackT).
sizeNumber of elements (or element tuples) in the stack.
calculatorCalculator object associated with the stack.
typeMultiplicityElement tuple size.
recommendedChunkSizeMultipleSuggested quantization of elements per work chunk.
recommendedChunksPerRankSuggested number of work chunks per MPI rank.
autoBroadcastEnable or disable auto broadcast.
Returns
StackIdentifier Id of the newly generated stack as registered with the LoadManager.
See also
DataStackBase::StackType::Implicit
DataStack

◆ addPassiveStack()

template<class StackT >
StackIdentifier HMP::LoadManager::addPassiveStack ( StackT *const  data,
const int  size 
)
inline

Create a passive data stack and attach it to the LoadManager.

Template Parameters
StackTThe fundamental data type of the DataStack to be created.
Parameters
dataData array on which the stack operates. The allocated size should be at least size * typeMultiplicity * typeof(StackT).
sizeNumber of elements (or element tuples) in the stack.
Returns
StackIdentifier Id of the newly generated stack as registered with the LoadManager.
See also
DataStackBase::StackType::Passive
DataStack

◆ addSlaveStack()

template<class StackT >
StackIdentifier HMP::LoadManager::addSlaveStack ( StackT *const  data,
const int  size,
const StackIdentifier  master,
const int  typeMultiplicity = 1 
)
inline

Create a slave data stack and attach it to the LoadManager.

Template Parameters
StackTThe fundamental data type of the DataStack to be created.
Parameters
dataData array on which the stack operates. The allocated size should be at least size * typeMultiplicity * typeof(StackT).
sizeNumber of elements (or element tuples) in the stack.
masterStack id of the associated implicit master stack.
typeMultiplicityElement tuple size.
Returns
StackIdentifier Id of the newly generated stack as registered with the LoadManager.
See also
DataStackBase::StackType::Slave
DataStack

◆ broadcast() [1/3]

void HMP::LoadManager::broadcast ( const StackIdentifier stackIds,
const int  size 
)
inline

Broadcast a list of stacks, where the stack identifiers are provided in list form.

Parameters
stackIdsPointer to the first StackIdentifier.
sizeNumber of stacks.

◆ broadcast() [2/3]

void HMP::LoadManager::broadcast ( const StackIdentifier  stackId)
inline

Broadcast a single stack.

Parameters
stackIdStackIdentifier of the stack to be broadcasted.

◆ broadcast() [3/3]

void HMP::LoadManager::broadcast ( const std::initializer_list< StackIdentifier > &  stackIds)
inline

Broadcast a list of stacks, where the stack identifiers are provided in initializer list form.

Parameters
stackIdsInitializer list of StackIdentifiers.

◆ calculate() [1/3]

virtual void HMP::LoadManager::calculate ( const StackIdentifier stackIds,
const int  size 
)
pure virtual

Calculate a list of stacks, where the stack identifiers are provided in list form.

Parameters
stackIdsPointer to the first StackIdentifier.
sizeNumber of stacks.

Implemented in HMP::LoadManagerSlave, and HMP::LoadManagerMaster.

◆ calculate() [2/3]

void HMP::LoadManager::calculate ( const StackIdentifier  stackId)
inline

Calculate a single stack.

Parameters
stackIdStackIdentifier of the stack to be calculated.

◆ calculate() [3/3]

void HMP::LoadManager::calculate ( const std::initializer_list< StackIdentifier > &  stackIds)
inline

Calculate a list of stacks, where the stack identifiers are provided in initializer list form.

Parameters
stackIdsInitializer list of StackIdentifiers.

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