Dune Core Modules (2.6.0)

Parallel Computing based on Indexsets

Provides classes for syncing distributed indexed data structures. More...

Files

file  communicator.hh
 Provides utility classes for syncing distributed data via MPI communication.
 
file  indexset.hh
 Provides a map between global and local indices.
 
file  indicessyncer.hh
 Class for adding missing indices of a distributed index set in a local communication.
 
file  interface.hh
 Provides classes for building the communication interface between remote indices.
 
file  localindex.hh
 Provides classes for use as the local index in ParallelIndexSet.
 
file  plocalindex.hh
 Provides classes for use as the local index in ParallelIndexSet for distributed computing.
 
file  remoteindices.hh
 Classes describing a distributed indexset.
 
file  selection.hh
 Provides classes for selecting indices based on attribute flags.
 
file  variablesizecommunicator.hh
 A communicator that only needs to know the number of elements per index at the sender side.
 

Classes

struct  Dune::SizeOne
 Flag for marking indexed data structures where data at each index is of the same size. More...
 
struct  Dune::VariableSize
 Flag for marking indexed data structures where the data at each index may be a variable multiple of another type. More...
 
struct  Dune::CommPolicy< V >
 Default policy used for communicating an indexed type. More...
 
class  Dune::BufferedCommunicator
 A communicator that uses buffers to gather and scatter the data to be send or received. More...
 
class  Dune::IndexPair< TG, TL >
 A pair consisting of a global and local index. More...
 
struct  Dune::MPITraits< T >
 A traits class describing the mapping of types onto MPI_Datatypes. More...
 
class  Dune::InvalidIndexSetState
 Exception indicating that the index set is not in the expected state. More...
 
class  Dune::GlobalLookupIndexSet< I >
 Decorates an index set with the possibility to find a global index that is mapped to a specific local. More...
 
class  Dune::ParallelIndexSet< TG, TL, N >
 Manager class for the mapping between local indices and globally unique indices. More...
 
class  Dune::ParallelIndexSet< TG, TL, N >::iterator
 The iterator over the pairs. More...
 
class  Dune::IndicesSyncer< T >
 Class for recomputing missing indices of a distributed index set. More...
 
class  Dune::Interface
 Communication interface between remote and local indices. More...
 
class  Dune::LocalIndex
 An index present on the local process. More...
 
class  Dune::ParallelLocalIndex< T >
 An index present on the local process with an additional attribute flag. More...
 
class  Dune::MPITraits< ParallelLocalIndex< T > >
 
class  Dune::MPITraits< IndexPair< TG, ParallelLocalIndex< TA > > >
 
class  Dune::RemoteIndices< T, A >
 The indices present on remote processes. More...
 
class  Dune::RemoteIndex< T1, T2 >
 Information about an index residing on another processor. More...
 
class  Dune::RemoteIndexListModifier< T, A, mode >
 Modifier for adding and/or deleting remote indices from the remote index list. More...
 
class  Dune::CollectiveIterator< T, A >
 A collective iterator for moving over the remote indices for all processes collectively. More...
 
class  Dune::OwnerOverlapCopyCommunication< GlobalIdType, LocalIdType >
 A class setting up standard communication for a two-valued attribute set with owner/overlap/copy semantics. More...
 
class  Dune::SelectionIterator< TS, TG, TL, N >
 A const iterator over an uncached selection. More...
 
class  Dune::UncachedSelection< TS, TG, TL, N >
 An uncached selection of indices. More...
 
class  Dune::Selection< TS, TG, TL, N >
 A cached selection of indices. More...
 
class  Dune::VariableSizeCommunicator< Allocator >
 A buffered communicator where the amount of data sent does not have to be known a priori. More...
 

Typedefs

typedef TG Dune::IndexPair< TG, TL >::GlobalIndex
 the type of the global index. More...
 
typedef TL Dune::IndexPair< TG, TL >::LocalIndex
 the type of the local index. More...
 
typedef TG Dune::ParallelIndexSet< TG, TL, N >::GlobalIndex
 the type of the global index. This type has to provide at least a operator< for sorting.
 
typedef TL Dune::ParallelIndexSet< TG, TL, N >::LocalIndex
 The type of the local index, e.g. ParallelLocalIndex. More...
 
typedef Dune::IndexPair< GlobalIndex, LocalIndexDune::ParallelIndexSet< TG, TL, N >::IndexPair
 The type of the pair stored.
 
typedef ArrayList< IndexPair, N >::const_iterator Dune::ParallelIndexSet< TG, TL, N >::const_iterator
 The constant iterator over the pairs.
 
typedef I Dune::GlobalLookupIndexSet< I >::ParallelIndexSet
 The type of the index set.
 
typedef ParallelIndexSet::LocalIndex Dune::GlobalLookupIndexSet< I >::LocalIndex
 The type of the local index.
 
typedef ParallelIndexSet::GlobalIndex Dune::GlobalLookupIndexSet< I >::GlobalIndex
 The type of the global index.
 
typedef ParallelIndexSet::const_iterator Dune::GlobalLookupIndexSet< I >::const_iterator
 The iterator over the index pairs.
 
typedef T Dune::IndicesSyncer< T >::ParallelIndexSet
 The type of the index set.
 
typedef ParallelIndexSet::IndexPair Dune::IndicesSyncer< T >::IndexPair
 The type of the index pair.
 
typedef ParallelIndexSet::GlobalIndex Dune::IndicesSyncer< T >::GlobalIndex
 Type of the global index used in the index set.
 
typedef ParallelIndexSet::LocalIndex::Attribute Dune::IndicesSyncer< T >::Attribute
 Type of the attribute used in the index set.
 
typedef Dune::RemoteIndices< ParallelIndexSetDune::IndicesSyncer< T >::RemoteIndices
 Type of the remote indices.
 

Enumerations

enum  Dune::ParallelIndexSetState { Dune::GROUND , Dune::RESIZE }
 The states the index set can be in. More...
 
enum  { Dune::ParallelIndexSet< TG, TL, N >::arraySize = (N>0) ? N : 1 }
 
enum  Dune::LocalIndexState
 The states avaiable for the local indices. More...
 

Functions

template<class TG , class TL >
std::ostream & Dune::operator<< (std::ostream &os, const IndexPair< TG, TL > &pair)
 Print an index pair. More...
 
 Dune::IndexPair< TG, TL >::IndexPair (const GlobalIndex &global, const LocalIndex &local)
 Constructs a new Pair. More...
 
 Dune::IndexPair< TG, TL >::IndexPair ()
 Construct a new Pair.
 
 Dune::IndexPair< TG, TL >::IndexPair (const GlobalIndex &global)
 Constructs a new Pair. More...
 
const GlobalIndexDune::IndexPair< TG, TL >::global () const
 Get the global index. More...
 
LocalIndexDune::IndexPair< TG, TL >::local ()
 Get the local index. More...
 
const LocalIndexDune::IndexPair< TG, TL >::local () const
 Get the local index. More...
 
void Dune::IndexPair< TG, TL >::setLocal (int index)
 Set the local index. More...
 
 Dune::ParallelIndexSet< TG, TL, N >::ParallelIndexSet ()
 Constructor.
 
const ParallelIndexSetStateDune::ParallelIndexSet< TG, TL, N >::state ()
 Get the state the index set is in. More...
 
void Dune::ParallelIndexSet< TG, TL, N >::beginResize ()
 Indicate that the index set is to be resized. More...
 
void Dune::ParallelIndexSet< TG, TL, N >::add (const GlobalIndex &global)
 Add an new index to the set. More...
 
void Dune::ParallelIndexSet< TG, TL, N >::add (const GlobalIndex &global, const LocalIndex &local)
 Add an new index to the set. More...
 
void Dune::ParallelIndexSet< TG, TL, N >::markAsDeleted (const iterator &position)
 Mark an index as deleted. More...
 
void Dune::ParallelIndexSet< TG, TL, N >::endResize ()
 Indicate that the resizing finishes. More...
 
IndexPairDune::ParallelIndexSet< TG, TL, N >::operator[] (const GlobalIndex &global)
 Find the index pair with a specific global id. More...
 
IndexPairDune::ParallelIndexSet< TG, TL, N >::at (const GlobalIndex &global)
 Find the index pair with a specific global id. More...
 
const IndexPairDune::ParallelIndexSet< TG, TL, N >::operator[] (const GlobalIndex &global) const
 Find the index pair with a specific global id. More...
 
const IndexPairDune::ParallelIndexSet< TG, TL, N >::at (const GlobalIndex &global) const
 Find the index pair with a specific global id. More...
 
iterator Dune::ParallelIndexSet< TG, TL, N >::begin ()
 Get an iterator over the indices positioned at the first index. More...
 
iterator Dune::ParallelIndexSet< TG, TL, N >::end ()
 Get an iterator over the indices positioned after the last index. More...
 
const_iterator Dune::ParallelIndexSet< TG, TL, N >::begin () const
 Get an iterator over the indices positioned at the first index. More...
 
const_iterator Dune::ParallelIndexSet< TG, TL, N >::end () const
 Get an iterator over the indices positioned after the last index. More...
 
void Dune::ParallelIndexSet< TG, TL, N >::renumberLocal ()
 Renumbers the local index numbers. More...
 
int Dune::ParallelIndexSet< TG, TL, N >::seqNo () const
 Get the internal sequence number. More...
 
size_t Dune::ParallelIndexSet< TG, TL, N >::size () const
 Get the total number (public and nonpublic) indices. More...
 
template<class TG , class TL , int N>
std::ostream & Dune::operator<< (std::ostream &os, const ParallelIndexSet< TG, TL, N > &indexSet)
 Print an index set. More...
 
 Dune::GlobalLookupIndexSet< I >::GlobalLookupIndexSet (const ParallelIndexSet &indexset, std::size_t size)
 Constructor. More...
 
 Dune::GlobalLookupIndexSet< I >::GlobalLookupIndexSet (const ParallelIndexSet &indexset)
 Constructor. More...
 
 Dune::GlobalLookupIndexSet< I >::~GlobalLookupIndexSet ()
 Destructor.
 
const IndexPairDune::GlobalLookupIndexSet< I >::operator[] (const GlobalIndex &global) const
 Find the index pair with a specific global id. More...
 
const IndexPairDune::GlobalLookupIndexSet< I >::pair (const std::size_t &local) const
 Get the index pair corresponding to a local index.
 
const_iterator Dune::GlobalLookupIndexSet< I >::begin () const
 Get an iterator over the indices positioned at the first index. More...
 
const_iterator Dune::GlobalLookupIndexSet< I >::end () const
 Get an iterator over the indices positioned after the last index. More...
 
int Dune::GlobalLookupIndexSet< I >::seqNo () const
 Get the internal sequence number. More...
 
size_t Dune::GlobalLookupIndexSet< I >::size () const
 Get the total number (public and nonpublic) indices. More...
 
 Dune::IndicesSyncer< T >::IndicesSyncer (ParallelIndexSet &indexSet, RemoteIndices &remoteIndices)
 Constructor. More...
 
void Dune::IndicesSyncer< T >::sync ()
 Sync the index set. More...
 
template<typename T1 >
void Dune::IndicesSyncer< T >::sync (T1 &numberer)
 Synce the index set and assign local numbers to new indices. More...
 
std::size_t Dune::IndicesSyncer< T >::DefaultNumberer::operator() (const GlobalIndex &global)
 Provide the lcoal index, always std::numeric_limits<size_t>::max() More...
 
 Dune::IndicesSyncer< T >::Iterators::Iterators (RemoteIndexList &remoteIndices, GlobalIndexList &globalIndices, BoolList &booleans)
 Constructor. More...
 
 Dune::IndicesSyncer< T >::Iterators::Iterators ()
 Default constructor.
 
Iterators & Dune::IndicesSyncer< T >::Iterators::operator++ ()
 Increment all iteraors.
 
void Dune::IndicesSyncer< T >::Iterators::insert (const RemoteIndex &index, const std::pair< GlobalIndex, Attribute > &global)
 Insert a new remote index to the underlying remote index list. More...
 
RemoteIndexDune::IndicesSyncer< T >::Iterators::remoteIndex () const
 Get the remote index at current position. More...
 
std::pair< GlobalIndex, Attribute > & Dune::IndicesSyncer< T >::Iterators::globalIndexPair () const
 Get the global index of the remote index at current position. More...
 
bool Dune::IndicesSyncer< T >::Iterators::isOld () const
 Was this entry already in the remote index list before the sync process? More...
 
void Dune::IndicesSyncer< T >::Iterators::reset (RemoteIndexList &remoteIndices, GlobalIndexList &globalIndices, BoolList &booleans)
 Reset all the underlying iterators. More...
 
bool Dune::IndicesSyncer< T >::Iterators::isNotAtEnd () const
 Are we not at the end of the list? More...
 
bool Dune::IndicesSyncer< T >::Iterators::isAtEnd () const
 Are we at the end of the list? More...
 
template<typename T , typename A , typename A1 >
void Dune::storeGlobalIndicesOfRemoteIndices (std::map< int, SLList< std::pair< typename T::GlobalIndex, typename T::LocalIndex::Attribute >, A > > &globalMap, const RemoteIndices< T, A1 > &remoteIndices)
 Stores the corresponding global indices of the remote index information. More...
 
template<typename T , typename A , typename A1 >
void Dune::repairLocalIndexPointers (std::map< int, SLList< std::pair< typename T::GlobalIndex, typename T::LocalIndex::Attribute >, A > > &globalMap, RemoteIndices< T, A1 > &remoteIndices, const T &indexSet)
 Repair the pointers to the local indices in the remote indices. More...
 
template<class T >
std::ostream & Dune::operator<< (std::ostream &os, const ParallelLocalIndex< T > &index)
 Print the local index to a stream. More...
 
template<class R , class T1 , class T2 , class Op , bool send>
void Dune::InterfaceBuilder::buildInterface (const R &remoteIndices, const T1 &sourceFlags, const T2 &destFlags, Op &functor) const
 Builds the interface between remote processes. More...
 
MPI_Comm Dune::Interface::communicator () const
 Get the MPI Communicator.
 
const InformationMapDune::Interface::interfaces () const
 Get information about the interfaces. More...
 
InformationMapDune::Interface::interfaces ()
 Get information about the interfaces. More...
 
void Dune::Interface::print () const
 Print the interface to std::out for debugging.
 
template<typename R , typename T1 , typename T2 >
void Dune::Interface::build (const R &remoteIndices, const T1 &sourceFlags, const T2 &destFlags)
 Builds the interface. More...
 
void Dune::Interface::free ()
 Frees memory allocated during the build.
 
virtual Dune::Interface::~Interface ()
 Destructor.
 
const std::size_t & Dune::LocalIndex::local () const
 get the local index. More...
 
 Dune::LocalIndex::operator std::size_t () const
 Convert to the local index represented by an int.
 
LocalIndexDune::LocalIndex::operator= (std::size_t index)
 Assign a new local index. More...
 
LocalIndexState Dune::LocalIndex::state () const
 Get the state. More...
 
void Dune::LocalIndex::setState (LocalIndexState state)
 Set the state. More...
 
 Dune::ParallelLocalIndex< T >::ParallelLocalIndex (const Attribute &attribute, bool isPublic)
 Constructor. More...
 
 Dune::ParallelLocalIndex< T >::ParallelLocalIndex (size_t localIndex, const Attribute &attribute, bool isPublic=true)
 Constructor. More...
 
 Dune::ParallelLocalIndex< T >::ParallelLocalIndex ()
 Parameterless constructor. More...
 
const Attribute Dune::ParallelLocalIndex< T >::attribute () const
 Get the attribute of the index. More...
 
void Dune::ParallelLocalIndex< T >::setAttribute (const Attribute &attribute)
 Set the attribute of the index. More...
 
size_t Dune::ParallelLocalIndex< T >::local () const
 get the local index. More...
 
 Dune::ParallelLocalIndex< T >::operator size_t () const
 Convert to the local index represented by an int.
 
ParallelLocalIndex< Attribute > & Dune::ParallelLocalIndex< T >::operator= (size_t index)
 Assign a new local index. More...
 
bool Dune::ParallelLocalIndex< T >::isPublic () const
 Check whether the index might also be known other processes. More...
 
LocalIndexState Dune::ParallelLocalIndex< T >::state () const
 Get the state. More...
 
void Dune::ParallelLocalIndex< T >::setState (const LocalIndexState &state)
 Set the state. More...
 
void Dune::Selection< TS, TG, TL, N >::setIndexSet (const ParallelIndexSet &indexset)
 Set the index set of the selection. More...
 
const_iterator Dune::Selection< TS, TG, TL, N >::begin () const
 Get the index set we are a selection for. More...
 
const_iterator Dune::Selection< TS, TG, TL, N >::end () const
 Get an iterator over the selected indices. More...
 
void Dune::Selection< TS, TG, TL, N >::free ()
 Free allocated memory.
 
const_iterator Dune::UncachedSelection< TS, TG, TL, N >::begin () const
 Get the index set we are a selection for. More...
 
const_iterator Dune::UncachedSelection< TS, TG, TL, N >::end () const
 Get an iterator over the selected indices. More...
 
void Dune::UncachedSelection< TS, TG, TL, N >::setIndexSet (const ParallelIndexSet &indexset)
 Set the index set of the selection. More...
 

Variables

int Dune::IndicesSyncer< T >::MessageInformation::publish
 The number of indices we publish for the other process.
 
int Dune::IndicesSyncer< T >::MessageInformation::pairs
 The number of pairs (attribute and process number) we publish to the neighbour process.
 

Detailed Description

Provides classes for syncing distributed indexed data structures.

In a parallel representation a container \(x\), e.g. a plain C-array, cannot be stored with all entries on each process because of limited memory and efficiency reasons. Therefore it is represented by individual pieces \(x_p\), \(p=0, \ldots, P-1\), where \(x_p\) is the piece stored on process \(p\) of the \(P\) processes participating in the calculation. Although the global representation of the container is not available on any process, a process \(p\) needs to know how the entries of it's local piece \(x_p\) correspond to the entries of the global container \(x\), which would be used in a sequential program. In this module we present classes describing the mapping of the local pieces to the global view and the communication interfaces.

Parallel Index Sets

Form an abstract point of view a random access container \(x: I \rightarrow K\) provides a mapping from an index set \(I \subset N_0\) onto a set of objects \(K\). Note that we do not require \(I\) to be consecutive. The piece \(x_p\) of the container \(x\) stored on process \(p\) is a mapping \(x_p:I_p \rightarrow K\), where \(I_p \subset I\). Due to efficiency the entries of \(x_p\) should be stored in consecutive memory.

This means that for the local computation the data must be addressable by a consecutive index starting from \(0\). When using adaptive discretisation methods there might be the need to reorder the indices after adding and/or deleting some of the discretisation points. Therefore this index does not have to be persistent. Further on we will call this index local index.

For the communication phases of our algorithms these locally stored entries must also be addressable by a global identifier to be able to store the received values tagged with the global identifiers at the correct local index in the consecutive local memory chunk. To ease the addition and removal of discretisation points this global identifier has to be persistent. Further on we will call this global identifier global index.

Classes to build the mapping are ParallelIndexSet and ParallelLocalIndex. As these just provide a mapping from the global index to the local index, the wrapper class GlobalLookupIndexSet facilitates the reverse lookup.

Remote Index Information

To setup communication between the processes every process needs to know what indices are also known to other processes and what attributes are attached to them on the remote side. This information is calculated and encapsulated in class RemoteIndices.

Communication

Based on the information about the distributed index sets, data independent interfaces between different sets of the index sets can be setup using the class Interface. For the actual communication data dependent communicators can be setup using BufferedCommunicator, DatatypeCommunicator VariableSizeCommunicator based on the interface information. In contrast to the former the latter is independent of the class Interface can work on a map from process number to a pair of index lists describing which local indices are send and received from that processs, respectively.

Typedef Documentation

◆ GlobalIndex

template<class TG , class TL >
typedef TG Dune::IndexPair< TG, TL >::GlobalIndex

the type of the global index.

This type has to provide at least a operator< for sorting.

◆ LocalIndex [1/2]

template<class TG , class TL >
typedef TL Dune::IndexPair< TG, TL >::LocalIndex

the type of the local index.

This class to provide the following functions:

LocalIndex operator=(int);
operator int() const;
LocalIndexState state() const;
void setState(LocalIndexState);
TL LocalIndex
the type of the local index.
Definition: indexset.hh:119
LocalIndexState
The states avaiable for the local indices.
Definition: localindex.hh:26

◆ LocalIndex [2/2]

template<typename TG , typename TL , int N = 100>
typedef TL Dune::ParallelIndexSet< TG, TL, N >::LocalIndex

The type of the local index, e.g. ParallelLocalIndex.

This class to provide the following functions:

LocalIndex operator=(int);
operator int() const;
void setState(LocalIndexState);
const ParallelIndexSetState & state()
Get the state the index set is in.
Definition: indexset.hh:318
TL LocalIndex
The type of the local index, e.g. ParallelLocalIndex.
Definition: indexset.hh:238

Enumeration Type Documentation

◆ anonymous enum

template<typename TG , typename TL , int N = 100>
anonymous enum
Enumerator
arraySize 

The size of the individual arrays in the underlying ArrayList.

The default value is 100.

See also
ArrayList::size

◆ LocalIndexState

The states avaiable for the local indices.

See also
LocalIndex::state()

◆ ParallelIndexSetState

The states the index set can be in.

See also
ParallelIndexSet::state_
Enumerator
GROUND 

The default mode. Indicates that the index set is ready to be used.

RESIZE 

Indicates that the index set is currently being resized.

Indicates that all previously deleted indices are now deleted.

CLEAN,

Indicates that the index set is currently being reordered.

REORDER

Function Documentation

◆ add() [1/2]

template<typename TG , typename TL , int N = 100>
void Dune::ParallelIndexSet< TG, TL, N >::add ( const GlobalIndex global)
inline

Add an new index to the set.

The local index is created by the default constructor.

Parameters
globalThe globally unique id of the index.
Exceptions
InvalidStateIf index set is not in ParallelIndexSetState::RESIZE mode.

Referenced by Dune::OwnerOverlapCopyCommunication< GlobalIdType, LocalIdType >::OwnerOverlapCopyCommunication().

◆ add() [2/2]

template<typename TG , typename TL , int N = 100>
void Dune::ParallelIndexSet< TG, TL, N >::add ( const GlobalIndex global,
const LocalIndex local 
)
inline

Add an new index to the set.

Parameters
globalThe globally unique id of the index.
localThe local index.
Exceptions
InvalidStateIf index set is not in ParallelIndexSetState::RESIZE mode.

◆ at() [1/2]

template<typename TG , typename TL , int N = 100>
IndexPair & Dune::ParallelIndexSet< TG, TL, N >::at ( const GlobalIndex global)
inline

Find the index pair with a specific global id.

This starts a binary search for the entry and therefore has complexity log(N).

Parameters
globalThe globally unique id of the pair.
Returns
The pair of indices for the id.
Exceptions
RangeErrorThrown if the global id is not known.

◆ at() [2/2]

template<typename TG , typename TL , int N = 100>
const IndexPair & Dune::ParallelIndexSet< TG, TL, N >::at ( const GlobalIndex global) const
inline

Find the index pair with a specific global id.

This starts a binary search for the entry and therefore has complexity log(N).

Parameters
globalThe globally unique id of the pair.
Returns
The pair of indices for the id.
Exceptions
RangeErrorThrown if the global id is not known.

◆ attribute()

template<class T >
const T Dune::ParallelLocalIndex< T >::attribute
inline

Get the attribute of the index.

Returns
The associated attribute.

◆ begin() [1/5]

◆ begin() [2/5]

template<typename TG , typename TL , int N = 100>
const_iterator Dune::ParallelIndexSet< TG, TL, N >::begin ( ) const
inline

Get an iterator over the indices positioned at the first index.

Returns
Iterator over the local indices.

◆ begin() [3/5]

template<class I >
const_iterator Dune::GlobalLookupIndexSet< I >::begin ( ) const
inline

Get an iterator over the indices positioned at the first index.

Returns
Iterator over the local indices.

◆ begin() [4/5]

template<typename TS , typename TG , typename TL , int N>
SelectionIterator< TS, TG, TL, N > Dune::UncachedSelection< TS, TG, TL, N >::begin

Get the index set we are a selection for.

Get an iterator over the selected indices.

Returns
An iterator positioned at the first selected index.

◆ begin() [5/5]

template<typename TS , typename TG , typename TL , int N>
uint32_t * Dune::Selection< TS, TG, TL, N >::begin

Get the index set we are a selection for.

Get an iterator over the selected indices.

Returns
An iterator positioned at the first selected index.

◆ beginResize()

template<typename TG , typename TL , int N = 100>
void Dune::ParallelIndexSet< TG, TL, N >::beginResize ( )

Indicate that the index set is to be resized.

Exceptions
InvalidStateIf index set was not in ParallelIndexSetState::GROUND mode.

Referenced by Dune::OwnerOverlapCopyCommunication< GlobalIdType, LocalIdType >::OwnerOverlapCopyCommunication().

◆ build()

template<typename R , typename T1 , typename T2 >
void Dune::Interface::build ( const R &  remoteIndices,
const T1 &  sourceFlags,
const T2 &  destFlags 
)
inline

Builds the interface.

The types T1 and T2 are classes representing a set of enumeration values of type Interface::Attribute. They have to provide a (static) method

bool contains(Attribute flag) const;

for checking whether the set contains a specific flag. This functionality is for example provided the classes EnumItem, EnumRange and Combine.

Parameters
remoteIndicesThe indices known to remote processes.
sourceFlagsThe set of flags marking indices we send from.
destFlagsThe set of flags marking indices we receive for.

References Dune::Interface::communicator_.

◆ buildInterface()

template<class R , class T1 , class T2 , class Op , bool send>
void Dune::InterfaceBuilder::buildInterface ( const R &  remoteIndices,
const T1 &  sourceFlags,
const T2 &  destFlags,
Op &  functor 
) const
protected

Builds the interface between remote processes.

The types T1 and T2 are classes representing a set of enumeration values of type InterfaceBuilder::Attribute. They have to provide a (static) method

bool contains(Attribute flag) const;

for checking whether the set contains a specific flag. This functionality is for example provided the classes EnumItem, EnumRange and Combine.

If the template parameter send is true the sending side of the interface will be built, otherwise the information for receiving will be built.

If the template parameter send is true we create interface for sending in a forward communication.

Parameters
remoteIndicesThe indices known to remote processes.
sourceFlagsThe set of flags marking source indices.
destFlagsThe setof flags markig destination indices.
functorA functor for callbacks. It should provide the following methods:
// Reserve memory for the interface to processor proc. The interface
// has to hold size entries
void reserve(int proc, int size);
// Add an entry to the interface
// We will send/receive size entries at index local to process proc
void add(int proc, int local);

References DUNE_THROW, and Dune::RemoteIndices< T, A >::end().

◆ end() [1/5]

◆ end() [2/5]

template<typename TG , typename TL , int N = 100>
const_iterator Dune::ParallelIndexSet< TG, TL, N >::end ( ) const
inline

Get an iterator over the indices positioned after the last index.

Returns
Iterator over the local indices.

◆ end() [3/5]

template<class I >
const_iterator Dune::GlobalLookupIndexSet< I >::end ( ) const
inline

Get an iterator over the indices positioned after the last index.

Returns
Iterator over the local indices.

◆ end() [4/5]

template<typename TS , typename TG , typename TL , int N>
SelectionIterator< TS, TG, TL, N > Dune::UncachedSelection< TS, TG, TL, N >::end

Get an iterator over the selected indices.

Returns
An iterator positioned at the first selected index.

◆ end() [5/5]

template<typename TS , typename TG , typename TL , int N>
uint32_t * Dune::Selection< TS, TG, TL, N >::end

Get an iterator over the selected indices.

Returns
An iterator positioned at the first selected index.

◆ endResize()

template<typename TG , typename TL , int N = 100>
void Dune::ParallelIndexSet< TG, TL, N >::endResize ( )

Indicate that the resizing finishes.

Warning
Invalidates all pointers stored to the elements of this index set. The local indices will be ordered according to the global indices: Let \((g_i,l_i)_{i=0}^N \) be the set of all indices then \(l_i < l_j\) if and only if \(g_i < g_j\) for arbitrary \(i \neq j\).
Exceptions
InvalidStateIf index set was not in ParallelIndexSetState::RESIZE mode.

Referenced by Dune::OwnerOverlapCopyCommunication< GlobalIdType, LocalIdType >::OwnerOverlapCopyCommunication().

◆ global()

template<class TG , class TL >
const GlobalIndex & Dune::IndexPair< TG, TL >::global ( ) const
inline

Get the global index.

Returns
The global index.

Referenced by Dune::RemoteIndexListModifier< T, A, mode >::insert().

◆ globalIndexPair()

template<typename T >
std::pair< typename IndicesSyncer< T >::GlobalIndex, typename IndicesSyncer< T >::Attribute > & Dune::IndicesSyncer< T >::Iterators::globalIndexPair
inline

Get the global index of the remote index at current position.

Returns
The current global index.

◆ GlobalLookupIndexSet() [1/2]

template<class I >
Dune::GlobalLookupIndexSet< I >::GlobalLookupIndexSet ( const ParallelIndexSet indexset)

Constructor.

Parameters
indexsetThe index set we want to be able to lookup the corresponding global index of a local index.

◆ GlobalLookupIndexSet() [2/2]

template<class I >
Dune::GlobalLookupIndexSet< I >::GlobalLookupIndexSet ( const ParallelIndexSet indexset,
std::size_t  size 
)

Constructor.

Parameters
indexsetThe index set we want to be able to lookup the corresponding global index of a local index.
sizeThe number of indices present, i.e. one more than the maximum local index.

◆ IndexPair() [1/2]

template<class TG , class TL >
Dune::IndexPair< TG, TL >::IndexPair ( const GlobalIndex global)

Constructs a new Pair.

The local index will be 0.

Parameters
globalThe global index.

◆ IndexPair() [2/2]

template<class TG , class TL >
Dune::IndexPair< TG, TL >::IndexPair ( const GlobalIndex global,
const LocalIndex local 
)

Constructs a new Pair.

Parameters
globalThe global index.
localThe local index.

◆ IndicesSyncer()

template<typename T >
Dune::IndicesSyncer< T >::IndicesSyncer ( ParallelIndexSet indexSet,
RemoteIndices remoteIndices 
)

Constructor.

The source as well as the target index set of the remote indices have to be the same as the provided index set.

Parameters
indexSetThe index set with the information of the locally present indices.
remoteIndicesThe remoteIndices.

References Dune::RemoteIndices< T, A >::communicator().

◆ insert()

template<typename T >
void Dune::IndicesSyncer< T >::Iterators::insert ( const RemoteIndex index,
const std::pair< GlobalIndex, Attribute > &  global 
)
inline

Insert a new remote index to the underlying remote index list.

Parameters
indexThe remote index.
globalThe global index corresponding to the remote index.

◆ interfaces() [1/2]

std::map< int, std::pair< InterfaceInformation, InterfaceInformation > > & Dune::Interface::interfaces ( )
inlineprotected

Get information about the interfaces.

Returns
Map of the interfaces. The key of the map is the process number and the value is the information pair (first the send and then the receive information).

◆ interfaces() [2/2]

const std::map< int, std::pair< InterfaceInformation, InterfaceInformation > > & Dune::Interface::interfaces ( ) const
inline

Get information about the interfaces.

Returns
Map of the interfaces. The key of the map is the process number and the value is the information pair (first the send and then the receive information).

◆ isAtEnd()

template<typename T >
bool Dune::IndicesSyncer< T >::Iterators::isAtEnd
inline

Are we at the end of the list?

Returns
True if the iterators are positioned at the end of the list and the tail of the list respectively.

◆ isNotAtEnd()

template<typename T >
bool Dune::IndicesSyncer< T >::Iterators::isNotAtEnd
inline

Are we not at the end of the list?

Returns
True if the iterators are not positioned at the end of the list and the tail of the list respectively.

◆ isOld()

template<typename T >
bool Dune::IndicesSyncer< T >::Iterators::isOld
inline

Was this entry already in the remote index list before the sync process?

Returns
True if the current index wasalready in the remote index list before the sync process.

◆ isPublic()

template<class T >
bool Dune::ParallelLocalIndex< T >::isPublic
inline

Check whether the index might also be known other processes.

Returns
True if the index might be known to other processors.

◆ Iterators()

template<typename T >
Dune::IndicesSyncer< T >::Iterators::Iterators ( RemoteIndexList remoteIndices,
GlobalIndexList globalIndices,
BoolList booleans 
)

Constructor.

Initializes all iterator to first entry and the one before the first entry, respectively.

Parameters
remoteIndicesThe list of the remote indices.
globalIndicesThe list of the coresponding global indices. This is needed because the the pointers to the local index will become invalid due to the merging of the index sets.
booleansWhether the remote index was there before the sync process started.

◆ local() [1/4]

template<class TG , class TL >
LocalIndex & Dune::IndexPair< TG, TL >::local ( )
inline

Get the local index.

Returns
The local index.

◆ local() [2/4]

template<class TG , class TL >
const LocalIndex & Dune::IndexPair< TG, TL >::local ( ) const
inline

Get the local index.

Returns
The local index.

◆ local() [3/4]

const std::size_t & Dune::LocalIndex::local ( ) const
inline

get the local index.

Returns
The local index.

◆ local() [4/4]

template<class T >
size_t Dune::ParallelLocalIndex< T >::local
inline

get the local index.

Returns
The local index.

◆ markAsDeleted()

template<typename TG , typename TL , int N = 100>
void Dune::ParallelIndexSet< TG, TL, N >::markAsDeleted ( const iterator position)
inline

Mark an index as deleted.

The index will be deleted during endResize().

Parameters
positionAn iterator at the position we want to delete.
Exceptions
InvalidStateIf index set is not in ParallelIndexSetState::RESIZE mode.

◆ operator()()

template<typename T >
std::size_t Dune::IndicesSyncer< T >::DefaultNumberer::operator() ( const GlobalIndex global)
inline

Provide the lcoal index, always std::numeric_limits<size_t>::max()

Parameters
globalThe global index (ignored).

References DUNE_UNUSED_PARAMETER.

◆ operator<<() [1/3]

template<class TG , class TL >
std::ostream & Dune::operator<< ( std::ostream &  os,
const IndexPair< TG, TL > &  pair 
)
inline

Print an index pair.

Parameters
osThe outputstream to print to.
pairThe index pair to print.

◆ operator<<() [2/3]

template<class TG , class TL , int N>
std::ostream & Dune::operator<< ( std::ostream &  os,
const ParallelIndexSet< TG, TL, N > &  indexSet 
)
inline

Print an index set.

Parameters
osThe outputstream to print to.
indexSetThe index set to print.

References Dune::ParallelIndexSet< TG, TL, N >::begin(), and Dune::ParallelIndexSet< TG, TL, N >::end().

◆ operator<<() [3/3]

template<class T >
std::ostream & Dune::operator<< ( std::ostream &  os,
const ParallelLocalIndex< T > &  index 
)

Print the local index to a stream.

Parameters
osThe output stream to print to.
indexThe index to print.

◆ operator=() [1/2]

template<class T >
ParallelLocalIndex< T > & Dune::ParallelLocalIndex< T >::operator= ( size_t  index)
inline

Assign a new local index.

Parameters
indexThe new local index.

◆ operator=() [2/2]

LocalIndex & Dune::LocalIndex::operator= ( std::size_t  index)
inline

Assign a new local index.

Parameters
indexThe new local index.

◆ operator[]() [1/3]

template<typename TG , typename TL , int N = 100>
IndexPair & Dune::ParallelIndexSet< TG, TL, N >::operator[] ( const GlobalIndex global)
inline

Find the index pair with a specific global id.

This starts a binary search for the entry and therefore has complexity log(N).

Parameters
globalThe globally unique id of the pair.
Returns
The pair of indices for the id.
Warning
If the global index is not in the set a wrong or even a null reference might be returned. To be save use the throwing alternative at.

◆ operator[]() [2/3]

template<typename TG , typename TL , int N = 100>
const IndexPair & Dune::ParallelIndexSet< TG, TL, N >::operator[] ( const GlobalIndex global) const
inline

Find the index pair with a specific global id.

This starts a binary search for the entry and therefore has complexity log(N).

Parameters
globalThe globally unique id of the pair.
Returns
The pair of indices for the id.
Warning
If the global index is not in the set a wrong or even a null reference might be returned. To be save use the throwing alternative at.

◆ operator[]() [3/3]

template<class I >
const IndexPair & Dune::GlobalLookupIndexSet< I >::operator[] ( const GlobalIndex global) const
inline

Find the index pair with a specific global id.

This starts a binary search for the entry and therefore has complexity log(N). This method is forwarded to the underlying index set.

Parameters
globalThe globally unique id of the pair.
Returns
The pair of indices for the id.
Exceptions
RangeErrorThrown if the global id is not known.

◆ ParallelLocalIndex() [1/3]

template<class T >
Dune::ParallelLocalIndex< T >::ParallelLocalIndex

Parameterless constructor.

Needed for use in container classes.

◆ ParallelLocalIndex() [2/3]

template<class T >
Dune::ParallelLocalIndex< T >::ParallelLocalIndex ( const Attribute attribute,
bool  isPublic 
)

Constructor.

The local index will be initialized to 0.

Parameters
attributeThe attribute of the index.
isPublicTrue if the index might also be known to other processes.

◆ ParallelLocalIndex() [3/3]

template<class T >
Dune::ParallelLocalIndex< T >::ParallelLocalIndex ( size_t  localIndex,
const Attribute attribute,
bool  isPublic = true 
)

Constructor.

Parameters
localIndexThe local index.
attributeThe attribute of the index.
isPublicTrue if the index might also be known to other processes.

◆ remoteIndex()

template<typename T >
IndicesSyncer< T >::RemoteIndex & Dune::IndicesSyncer< T >::Iterators::remoteIndex
inline

Get the remote index at current position.

Returns
The current remote index.

◆ renumberLocal()

template<typename TG , typename TL , int N = 100>
void Dune::ParallelIndexSet< TG, TL, N >::renumberLocal ( )
inline

Renumbers the local index numbers.

After this function returns the indices are consecutively numbered beginning from 0. Let $(g_i,l_i)$, $(g_j,l_j)$ be two arbitrary index pairs with $g_i<g_j$ then after renumbering $l_i<l_j$ will hold.

◆ repairLocalIndexPointers()

template<typename T , typename A , typename A1 >
void Dune::repairLocalIndexPointers ( std::map< int, SLList< std::pair< typename T::GlobalIndex, typename T::LocalIndex::Attribute >, A > > &  globalMap,
RemoteIndices< T, A1 > &  remoteIndices,
const T &  indexSet 
)
inline

Repair the pointers to the local indices in the remote indices.

Parameters
globalMapThe map of the process number to the list of global indices corresponding to the remote index list of the process.
remoteIndicesThe known remote indices.
indexSetThe set of local indices of the current process.

◆ reset()

template<typename T >
void Dune::IndicesSyncer< T >::Iterators::reset ( RemoteIndexList remoteIndices,
GlobalIndexList globalIndices,
BoolList booleans 
)
inline

Reset all the underlying iterators.

Position them to first list entry and the entry before the first entry respectively.

Parameters
remoteIndicesThe list of the remote indices.
globalIndicesThe list of the coresponding global indices. This is needed because the the pointers to the local index will become invalid due to the merging of the index sets.
booleansWhether the remote index was there before the sync process started.

References Dune::SLList< T, A >::beginModify().

◆ seqNo() [1/2]

template<typename TG , typename TL , int N = 100>
int Dune::ParallelIndexSet< TG, TL, N >::seqNo ( ) const
inline

Get the internal sequence number.

Is initially 0 is incremented for each resize.

Returns
The sequence number.

◆ seqNo() [2/2]

template<class I >
int Dune::GlobalLookupIndexSet< I >::seqNo ( ) const
inline

Get the internal sequence number.

Is initially 0 is incremented for each resize.

Returns
The sequence number.

◆ setAttribute()

template<class T >
void Dune::ParallelLocalIndex< T >::setAttribute ( const Attribute attribute)
inline

Set the attribute of the index.

Parameters
attributeThe associated attribute.

◆ setIndexSet() [1/2]

template<typename TS , typename TG , typename TL , int N>
void Dune::UncachedSelection< TS, TG, TL, N >::setIndexSet ( const ParallelIndexSet indexset)

Set the index set of the selection.

Parameters
indexsetThe index set to use.

◆ setIndexSet() [2/2]

template<typename TS , typename TG , typename TL , int N>
void Dune::Selection< TS, TG, TL, N >::setIndexSet ( const ParallelIndexSet indexset)
inline

Set the index set of the selection.

Parameters
indexsetThe index set to use.

References Dune::ParallelIndexSet< TG, TL, N >::begin(), and Dune::ParallelIndexSet< TG, TL, N >::end().

◆ setLocal()

template<class TG , class TL >
void Dune::IndexPair< TG, TL >::setLocal ( int  index)
inline

Set the local index.

Parameters
indexThe index to set it to.

◆ setState() [1/2]

template<class T >
void Dune::ParallelLocalIndex< T >::setState ( const LocalIndexState state)
inline

Set the state.

Parameters
stateThe state to set.

◆ setState() [2/2]

void Dune::LocalIndex::setState ( LocalIndexState  state)
inline

Set the state.

Parameters
stateThe state to set.

References Dune::LocalIndex::state().

◆ size() [1/2]

template<typename TG , typename TL , int N = 100>
size_t Dune::ParallelIndexSet< TG, TL, N >::size ( ) const
inline

Get the total number (public and nonpublic) indices.

Returns
The total number (public and nonpublic) indices.

Referenced by Dune::graphRepartition().

◆ size() [2/2]

template<class I >
size_t Dune::GlobalLookupIndexSet< I >::size ( ) const
inline

Get the total number (public and nonpublic) indices.

Returns
The total number (public and nonpublic) indices.

◆ state() [1/3]

template<typename TG , typename TL , int N = 100>
const ParallelIndexSetState & Dune::ParallelIndexSet< TG, TL, N >::state ( )
inline

Get the state the index set is in.

Returns
The state of the index set.

◆ state() [2/3]

LocalIndexState Dune::LocalIndex::state ( ) const
inline

Get the state.

Returns
The state.

Referenced by Dune::LocalIndex::setState().

◆ state() [3/3]

template<class T >
LocalIndexState Dune::ParallelLocalIndex< T >::state
inline

Get the state.

Returns
The state.

◆ storeGlobalIndicesOfRemoteIndices()

template<typename T , typename A , typename A1 >
void Dune::storeGlobalIndicesOfRemoteIndices ( std::map< int, SLList< std::pair< typename T::GlobalIndex, typename T::LocalIndex::Attribute >, A > > &  globalMap,
const RemoteIndices< T, A1 > &  remoteIndices 
)

Stores the corresponding global indices of the remote index information.

Whenever a ParallelIndexSet is resized all RemoteIndices that use it will be invalided as the pointers to the index set are invalid after calling ParallelIndexSet::Resize() One can rebuild them by storing the global indices in a map with this function and later repairing the pointers by calling repairLocalIndexPointers.

Warning
The RemoteIndices class has to be build with the same index set for both the sending and receiving side
Parameters
globalMapMap to store the corresponding global indices in.
remoteIndicesThe remote index information we need to store the corresponding global indices of.
indexSetThe index set that is for both the sending and receiving side of the remote index information.

◆ sync() [1/2]

template<typename T >
void Dune::IndicesSyncer< T >::sync
inline

Sync the index set.

Computes the missing indices in the local and the remote index list and adds them. No global communication is necessary! All indices added to the index will become the local index std::numeric_limits<size_t>::max()

References Dune::IndicesSyncer< T >::sync().

Referenced by Dune::IndicesSyncer< T >::sync().

◆ sync() [2/2]

template<typename T >
template<typename T1 >
void Dune::IndicesSyncer< T >::sync ( T1 &  numberer)

Synce the index set and assign local numbers to new indices.

Computes the missing indices in the local and the remote index list and adds them. No global communication is necessary!

Parameters
numbererFunctor providing the local indices for the added global indices. has to provide a function size_t operator()(const TG& global) that provides the local index to a global one. It will be called for ascending global indices.

References Dune::SLList< T, A >::begin(), Dune::RemoteIndices< T, A >::begin(), Dune::dverb, Dune::SLList< T, A >::end(), Dune::RemoteIndices< T, A >::end(), Dune::RemoteIndices< T, A >::neighbours(), and Dune::SLList< T, A >::push_back().

Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 26, 23:30, 2024)