ISTL Communication
[Iterative Solvers Template Library (ISTL)]

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

Collaboration diagram for ISTL Communication:

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...
class  Dune::InvalidIndexSetState
 Exception indicating that the index set is not in the expected state. 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::GlobalLookupIndexSet< I >
 Decorates an index set with the possibility to find a global index that is mapped to a specific local. More...
class  Dune::IndicesSyncer< T >
 Class for recomputing missing indices of a distributed index set. More...
struct  Dune::IndicesSyncer< T >::MessageInformation
 Information about the messages to send to a neighbouring process.
class  Dune::IndicesSyncer< T >::DefaultNumberer
 Default numberer for sync().
class  Dune::IndicesSyncer< T >::Iterators
 A tuple of iterators.
class  Dune::Interface
 Communication interface between remote and local indices. More...
class  Dune::LocalIndex
 An index present on the local process. More...
class  Dune::MPITraits< T >
 A traits class describing the mapping of types onto MPI_Datatypes. More...
struct  Dune::MPITraits< FieldVector< K, n > >
struct  Dune::MPITraits< bigunsignedint< k > >
class  Dune::MPITraits< std::pair< T1, T2 > >
struct  Dune::OwnerOverlapCopyAttributeSet
 Attribute set for overlapping schwarz. More...
class  Dune::IndexInfoFromGrid< G, L >
 Information about the index distribution. 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::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::RemoteIndex< T1, T2 >
 Information about an index residing on another processor. More...
class  Dune::RemoteIndices< T, A >
 The indices present on remote processes. 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 >
 An cached selection of indices. 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  owneroverlapcopy.hh
 

Classes providing communication interfaces for overlapping Schwarz methods.


file  plocalindex.hh
 

Provides classes for use as the local index in ParallelIndexSet for distributed computing.


file  remoteindices.hh
 

Classes discribing a distributed indexset.


file  selection.hh
 

Provides classes for selecting indices base on attribute flags.


Defines

#define ComposeMPITraits(p, m)

Typedefs

typedef TG Dune::IndexPair::GlobalIndex
 the type of the global index.
typedef TL Dune::IndexPair::LocalIndex
 the type of the local index.
typedef TG Dune::ParallelIndexSet::GlobalIndex
 the type of the global index. This type has to provide at least a operator< for sorting.
typedef TL Dune::ParallelIndexSet::LocalIndex
 The type of the local index, e.g. ParallelLocalIndex.
typedef Dune::IndexPair
< GlobalIndex, LocalIndex > 
Dune::ParallelIndexSet::IndexPair
 The type of the pair stored.
typedef ArrayList< IndexPair,
N >::const_iterator 
Dune::ParallelIndexSet::const_iterator
 The constant iterator over the pairs.
typedef
ParallelIndexSet::LocalIndex 
Dune::GlobalLookupIndexSet::LocalIndex
 The type of the local index.
typedef
ParallelIndexSet::GlobalIndex 
Dune::GlobalLookupIndexSet::GlobalIndex
 The type of the global index.
typedef
ParallelIndexSet::const_iterator 
Dune::GlobalLookupIndexSet::const_iterator
 The iterator over the index pairs.
typedef Dune::IndexPair
< typename I::GlobalIndex,
typename I::LocalIndex > 
Dune::GlobalLookupIndexSet::IndexPair
typedef ParallelIndexSet::IndexPair Dune::IndicesSyncer::IndexPair
 The type of the index pair.
typedef
ParallelIndexSet::GlobalIndex 
Dune::IndicesSyncer::GlobalIndex
 Type of the global index used in the index set.
typedef
ParallelIndexSet::LocalIndex::Attribute 
Dune::IndicesSyncer::Attribute
 Type of the attribute used in the index set.
typedef Dune::RemoteIndices
< ParallelIndexSet > 
Dune::IndicesSyncer::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::arraySize = (N>0)?N:1 }
enum  Dune::LocalIndexState { Dune::VALID, Dune::DELETED }
 

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.
template<class TG , class TL >
bool Dune::operator== (const IndexPair< TG, TL > &, const IndexPair< TG, TL > &)
template<class TG , class TL >
bool Dune::operator!= (const IndexPair< TG, TL > &, const IndexPair< TG, TL > &)
template<class TG , class TL >
bool Dune::operator< (const IndexPair< TG, TL > &, const IndexPair< TG, TL > &)
template<class TG , class TL >
bool Dune::operator> (const IndexPair< TG, TL > &, const IndexPair< TG, TL > &)
template<class TG , class TL >
bool Dune::operator<= (const IndexPair< TG, TL > &, const IndexPair< TG, TL > &)
template<class TG , class TL >
bool Dune::operator>= (const IndexPair< TG, TL > &, const IndexPair< TG, TL > &)
template<class TG , class TL >
bool Dune::operator== (const IndexPair< TG, TL > &, const TG &)
template<class TG , class TL >
bool Dune::operator!= (const IndexPair< TG, TL > &, const TG &)
template<class TG , class TL >
bool Dune::operator< (const IndexPair< TG, TL > &, const TG &)
template<class TG , class TL >
bool Dune::operator> (const IndexPair< TG, TL > &, const TG &)
template<class TG , class TL >
bool Dune::operator<= (const IndexPair< TG, TL > &, const TG &)
template<class TG , class TL >
bool Dune::operator>= (const IndexPair< TG, TL > &, const TG &)
 Dune::IndexPair::IndexPair (const GlobalIndex &global, const LocalIndex &local)
 Constructs a new Pair.
 Dune::IndexPair::IndexPair ()
 Construct a new Pair.
 Dune::IndexPair::IndexPair (const GlobalIndex &global)
 Constructs a new Pair.
const GlobalIndex & Dune::IndexPair::global () const
 Get the global index.
LocalIndex & Dune::IndexPair::local ()
 Get the local index.
const LocalIndex & Dune::IndexPair::local () const
 Get the local index.
void Dune::IndexPair::setLocal (int index)
 Set the local index.
 Dune::ParallelIndexSet::iterator::iterator (ParallelIndexSet< TG, TL, N > &indexSet, const Father &father)
 Dune::ParallelIndexSet::iterator::iterator (const iterator &other)
iterator & Dune::ParallelIndexSet::iterator::operator== (const iterator &other)
 Dune::ParallelIndexSet::ParallelIndexSet ()
 Constructor.
const ParallelIndexSetState & Dune::ParallelIndexSet::state ()
 Get the state the index set is in.
void Dune::ParallelIndexSet::beginResize () throw (InvalidIndexSetState)
 Indicate that the index set is to be resized.

Exceptions:
InvalidState If index set was not in ParallelIndexSetState::GROUND mode.

void Dune::ParallelIndexSet::add (const GlobalIndex &global) throw (InvalidIndexSetState)
 Add an new index to the set.
void Dune::ParallelIndexSet::add (const GlobalIndex &global, const LocalIndex &local) throw (InvalidIndexSetState)
 Add an new index to the set.
void Dune::ParallelIndexSet::markAsDeleted (const iterator &position) throw (InvalidIndexSetState)
 Mark an index as deleted.
void Dune::ParallelIndexSet::endResize () throw (InvalidIndexSetState)
 Indicate that the resizing finishes.
IndexPair & Dune::ParallelIndexSet::operator[] (const GlobalIndex &global)
 Find the index pair with a specific global id.
IndexPair & Dune::ParallelIndexSet::at (const GlobalIndex &global)
 Find the index pair with a specific global id.
const IndexPair & Dune::ParallelIndexSet::operator[] (const GlobalIndex &global) const
 Find the index pair with a specific global id.
const IndexPair & Dune::ParallelIndexSet::at (const GlobalIndex &global) const
 Find the index pair with a specific global id.
iterator Dune::ParallelIndexSet::begin ()
 Get an iterator over the indices positioned at the first index.
iterator Dune::ParallelIndexSet::end ()
 Get an iterator over the indices positioned after the last index.
const_iterator Dune::ParallelIndexSet::begin () const
 Get an iterator over the indices positioned at the first index.
const_iterator Dune::ParallelIndexSet::end () const
 Get an iterator over the indices positioned after the last index.
void Dune::ParallelIndexSet::renumberLocal ()
 Renumbers the local index numbers.
int Dune::ParallelIndexSet::seqNo () const
 Get the internal sequence number.
size_t Dune::ParallelIndexSet::size () const
 Get the total number (public and nonpublic) indices.
template<class TG , class TL , int N>
std::ostream & Dune::operator<< (std::ostream &os, const ParallelIndexSet< TG, TL, N > &indexSet)
 Print an index set.
 Dune::GlobalLookupIndexSet::GlobalLookupIndexSet (const ParallelIndexSet &indexset, std::size_t size)
 Constructor.
 Dune::GlobalLookupIndexSet::GlobalLookupIndexSet (const ParallelIndexSet &indexset)
 Constructor.
 Dune::GlobalLookupIndexSet::~GlobalLookupIndexSet ()
 Destructor.
const IndexPair & Dune::GlobalLookupIndexSet::operator[] (const GlobalIndex &global) const
 Find the index pair with a specific global id.
const IndexPair * Dune::GlobalLookupIndexSet::pair (const std::size_t &local) const
 Get the index pair corresponding to a local index.
const_iterator Dune::GlobalLookupIndexSet::begin () const
 Get an iterator over the indices positioned at the first index.
const_iterator Dune::GlobalLookupIndexSet::end () const
 Get an iterator over the indices positioned after the last index.
int Dune::GlobalLookupIndexSet::seqNo () const
 Get the internal sequence number.
size_t Dune::GlobalLookupIndexSet::size () const
 Get the total number (public and nonpublic) indices.
 Dune::IndicesSyncer::IndicesSyncer (ParallelIndexSet &indexSet, RemoteIndices &remoteIndices)
 Constructor.
void Dune::IndicesSyncer::sync ()
 Sync the index set.
template<typename T1 >
void Dune::IndicesSyncer::sync (T1 &numberer)
 Synce the index set and assign local numbers to new indices.
 Dune::IndicesSyncer::Iterators::Iterators (RemoteIndexList &remoteIndices, GlobalIndexList &globalIndices, BoolList &booleans)
 Constructor.
 Dune::IndicesSyncer::Iterators::Iterators ()
 Default constructor.
Iterators & Dune::IndicesSyncer::Iterators::operator++ ()
 Increment all iteraors.
void Dune::IndicesSyncer::Iterators::insert (const RemoteIndex &index, const GlobalIndex &global)
 Insert a new remote index to the underlying remote index list.
RemoteIndex & Dune::IndicesSyncer::Iterators::remoteIndex () const
 Get the remote index at current position.
GlobalIndex & Dune::IndicesSyncer::Iterators::globalIndex () const
 Get the global index of the remote index at current position.
bool Dune::IndicesSyncer::Iterators::isOld () const
 Was this entry already in the remote index list before the sync process?
void Dune::IndicesSyncer::Iterators::reset (RemoteIndexList &remoteIndices, GlobalIndexList &globalIndices, BoolList &booleans)
 Reset all the underlying iterators.
bool Dune::IndicesSyncer::Iterators::isNotAtEnd () const
 Are we not at the end of the list?
bool Dune::IndicesSyncer::Iterators::isAtEnd () const
 Are we at the end of the list?
template<typename T , typename A , typename A1 >
void Dune::storeGlobalIndicesOfRemoteIndices (std::map< int, SLList< typename T::GlobalIndex, A > > &globalMap, const RemoteIndices< T, A1 > &remoteIndices, const T &indexSet)
 Stores the corresponding global indices of the remote index information.
template<typename T , typename A , typename A1 >
void Dune::repairLocalIndexPointers (std::map< int, SLList< typename T::GlobalIndex, A > > &globalMap, RemoteIndices< T, A1 > &remoteIndices, const T &indexSet)
 Repair the pointers to the local indices in the remote indices.
 Dune::ComposeMPITraits (char, MPI_CHAR)
 Dune::ComposeMPITraits (unsigned char, MPI_UNSIGNED_CHAR)
 Dune::ComposeMPITraits (short, MPI_SHORT)
 Dune::ComposeMPITraits (unsigned short, MPI_UNSIGNED_SHORT)
 Dune::ComposeMPITraits (int, MPI_INT)
 Dune::ComposeMPITraits (unsigned int, MPI_UNSIGNED)
 Dune::ComposeMPITraits (long, MPI_LONG)
 Dune::ComposeMPITraits (unsigned long, MPI_UNSIGNED_LONG)
 Dune::ComposeMPITraits (float, MPI_FLOAT)
 Dune::ComposeMPITraits (double, MPI_DOUBLE)
 Dune::ComposeMPITraits (long double, MPI_LONG_DOUBLE)
static MPI_Datatype Dune::MPITraits< FieldVector< K, n > >::getType ()
static MPI_Datatype Dune::MPITraits< bigunsignedint< k > >::getType ()
template<class T >
std::ostream & Dune::operator<< (std::ostream &os, const ParallelLocalIndex< T > &index)
 Print the local index to a stream.
template<typename T1 , typename T2 >
std::ostream & Dune::operator<< (std::ostream &os, const RemoteIndex< T1, T2 > &index)
template<class T , class A >
std::ostream & Dune::operator<< (std::ostream &os, const RemoteIndices< T, A > &indices)
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.
MPI_Comm Dune::Interface::communicator () const
 Get the MPI Communicator.
const InformationMap & Dune::Interface::interfaces () const
 Get information about the interfaces.
InformationMap & Dune::Interface::interfaces ()
 Get information about the interfaces.
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.
void Dune::Interface::strip ()
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.
 Dune::LocalIndex::operator std::size_t () const
 Convert to the local index represented by an int.
LocalIndex & Dune::LocalIndex::operator= (std::size_t index)
 Assign a new local index.
LocalIndexState Dune::LocalIndex::state () const
 Get the state.
void Dune::LocalIndex::setState (LocalIndexState state)
 Set the state.
static MPI_Datatype Dune::MPITraits< std::pair< T1, T2 > >::getType ()
 Dune::ParallelLocalIndex::ParallelLocalIndex (const Attribute &attribute, bool isPublic)
 Constructor.
 Dune::ParallelLocalIndex::ParallelLocalIndex (size_t localIndex, const Attribute &attribute, bool isPublic=true)
 Constructor.
 Dune::ParallelLocalIndex::ParallelLocalIndex ()
 Parameterless constructor.
const Attribute Dune::ParallelLocalIndex::attribute () const
 Get the attribute of the index.
void Dune::ParallelLocalIndex::setAttribute (const Attribute &attribute)
 Set the attribute of the index.
size_t Dune::ParallelLocalIndex::local () const
 get the local index.
 Dune::ParallelLocalIndex::operator size_t () const
 Convert to the local index represented by an int.
ParallelLocalIndex< Attribute > & Dune::ParallelLocalIndex::operator= (size_t index)
 Assign a new local index.
bool Dune::ParallelLocalIndex::isPublic () const
 Check whether the index might also be known other processes.
LocalIndexState Dune::ParallelLocalIndex::state () const
 Get the state.
void Dune::ParallelLocalIndex::setState (const LocalIndexState &state)
 Set the state.
static MPI_Datatype Dune::MPITraits< ParallelLocalIndex< T > >::getType ()
void Dune::Selection::setIndexSet (const ParallelIndexSet &indexset)
 Set the index set of the selection.
const_iterator Dune::Selection::begin () const
 Get the index set we are a selection for.
const_iterator Dune::Selection::end () const
 Get an iterator over the selected indices.
void Dune::Selection::free ()
 Free allocated memory.
 Dune::Selection::~Selection ()
const_iterator Dune::UncachedSelection::begin () const
 Get the index set we are a selection for.
const_iterator Dune::UncachedSelection::end () const
 Get an iterator over the selected indices.
void Dune::UncachedSelection::setIndexSet (const ParallelIndexSet &indexset)
 Set the index set of the selection.

Variables

int Dune::IndicesSyncer::MessageInformation::publish
 The number of indices we publish for the other process.
int Dune::IndicesSyncer::MessageInformation::pairs
 The number of pairs (attribute and process number) we publish to the neighbour process.
static MPI_Datatype Dune::MPITraits< FieldVector< K, n > >::vectortype = {MPI_DATATYPE_NULL}
static MPI_Datatype Dune::MPITraits< bigunsignedint< k > >::vectortype = MPI_DATATYPE_NULL
static MPI_Datatype Dune::MPITraits< FieldVector< K, n > >::datatype = MPI_DATATYPE_NULL

Friends

class Dune::IndexPair::MPITraits< IndexPair< TG, TL > >
class Dune::ParallelIndexSet::iterator::ParallelIndexSet< GlobalIndex, LocalIndex, N >
bool Dune::IndexPair::operator== (const IndexPair< TG, TL > &, const IndexPair< TG, TL > &)
bool Dune::IndexPair::operator!= (const IndexPair< TG, TL > &, const IndexPair< TG, TL > &)
bool Dune::IndexPair::operator< (const IndexPair< TG, TL > &, const IndexPair< TG, TL > &)
bool Dune::IndexPair::operator> (const IndexPair< TG, TL > &, const IndexPair< TG, TL > &)
bool Dune::IndexPair::operator<= (const IndexPair< TG, TL > &, const IndexPair< TG, TL > &)
bool Dune::IndexPair::operator>= (const IndexPair< TG, TL > &, const IndexPair< TG, TL > &)
bool Dune::IndexPair::operator== (const IndexPair< TG, TL > &, const TG &)
bool Dune::IndexPair::operator!= (const IndexPair< TG, TL > &, const TG &)
bool Dune::IndexPair::operator< (const IndexPair< TG, TL > &, const TG &)
bool Dune::IndexPair::operator> (const IndexPair< TG, TL > &, const TG &)
bool Dune::IndexPair::operator<= (const IndexPair< TG, TL > &, const TG &)
bool Dune::IndexPair::operator>= (const IndexPair< TG, TL > &, const TG &)

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 it data dependant communicators can be setup using BufferedCommunicator or DatatypeCommunicator.


Define Documentation

#define ComposeMPITraits ( p,
 ) 
Value:
template<> \
  struct MPITraits<p>{ \
    static inline MPI_Datatype getType(){ \
      return m; \
    } \
  }

Typedef Documentation

template<typename T>
typedef ParallelIndexSet::LocalIndex::Attribute Dune::IndicesSyncer< T >::Attribute [inherited]

Type of the attribute used in the index set.

template<class I>
typedef ParallelIndexSet::const_iterator Dune::GlobalLookupIndexSet< I >::const_iterator [inherited]

The iterator over the index pairs.

template<typename TG, typename TL, int N = 100>
typedef ArrayList<IndexPair,N>::const_iterator Dune::ParallelIndexSet< TG, TL, N >::const_iterator [inherited]

The constant iterator over the pairs.

template<typename T>
typedef ParallelIndexSet::GlobalIndex Dune::IndicesSyncer< T >::GlobalIndex [inherited]

Type of the global index used in the index set.

template<class I>
typedef ParallelIndexSet::GlobalIndex Dune::GlobalLookupIndexSet< I >::GlobalIndex [inherited]

The type of the global index.

template<typename TG, typename TL, int N = 100>
typedef TG Dune::ParallelIndexSet< TG, TL, N >::GlobalIndex [inherited]

the type of the global index. This type has to provide at least a operator< for sorting.

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

the type of the global index.

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

template<typename T>
typedef ParallelIndexSet::IndexPair Dune::IndicesSyncer< T >::IndexPair [inherited]

The type of the index pair.

template<class I>
typedef Dune::IndexPair<typename I::GlobalIndex, typename I::LocalIndex> Dune::GlobalLookupIndexSet< I >::IndexPair [inherited]
template<typename TG, typename TL, int N = 100>
typedef Dune::IndexPair<GlobalIndex,LocalIndex> Dune::ParallelIndexSet< TG, TL, N >::IndexPair [inherited]

The type of the pair stored.

template<class I>
typedef ParallelIndexSet::LocalIndex Dune::GlobalLookupIndexSet< I >::LocalIndex [inherited]

The type of the local index.

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

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

This class to provide the following functions:

 LocalIndex operator=(int);
 operator int() const;
 LocalIndexState state() const;
 void setState(LocalIndexState);
template<class TG, class TL>
typedef TL Dune::IndexPair< TG, TL >::LocalIndex [inherited]

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);
template<typename T>
typedef Dune::RemoteIndices<ParallelIndexSet> Dune::IndicesSyncer< T >::RemoteIndices [inherited]

Type of the remote indices.


Enumeration Type Documentation

template<typename TG, typename TL, int N = 100>
anonymous enum [inherited]
Enumerator:
arraySize 

The size of the individual arrays in the underlying ArrayList.

The default value is 100.

See also:
ArrayList::size

The states avaiable for the local indices.

See also:
LocalIndex::state()
Enumerator:
VALID 
DELETED 

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

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

Add an new index to the set.

Parameters:
global The globally unique id of the index.
local The local index.
Exceptions:
InvalidState If index set is not in ParallelIndexSetState::RESIZE mode.
template<typename TG, typename TL, int N = 100>
void Dune::ParallelIndexSet< TG, TL, N >::add ( const GlobalIndex global  )  throw (InvalidIndexSetState) [inline, inherited]

Add an new index to the set.

The local index is created by the default constructor.

Parameters:
global The globally unique id of the index.
Exceptions:
InvalidState If index set is not in ParallelIndexSetState::RESIZE mode.

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

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

Find the index pair with a specific global id.

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

Parameters:
global The globally unique id of the pair.
Returns:
The pair of indices for the id.
Exceptions:
RangeError Thrown if the global id is not known.
template<typename TG, typename TL, int N = 100>
IndexPair& Dune::ParallelIndexSet< TG, TL, N >::at ( const GlobalIndex global  )  [inline, inherited]

Find the index pair with a specific global id.

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

Parameters:
global The globally unique id of the pair.
Returns:
The pair of indices for the id.
Exceptions:
RangeError Thrown if the global id is not known.
template<class T >
const T Dune::ParallelLocalIndex< T >::attribute (  )  const [inline, inherited]

Get the attribute of the index.

Returns:
The associated attribute.
template<typename TS , typename TG , typename TL , int N>
SelectionIterator< TS, TG, TL, N > Dune::UncachedSelection< TS, TG, TL, N >::begin (  )  const [inline, inherited]

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.

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

template<typename TS , typename TG , typename TL , int N>
uint32_t * Dune::Selection< TS, TG, TL, N >::begin (  )  const [inline, inherited]

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.
template<class I>
const_iterator Dune::GlobalLookupIndexSet< I >::begin (  )  const [inline, inherited]

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

Returns:
Iterator over the local indices.
template<typename TG, typename TL, int N = 100>
const_iterator Dune::ParallelIndexSet< TG, TL, N >::begin (  )  const [inline, inherited]

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

Returns:
Iterator over the local indices.
template<typename TG, typename TL, int N = 100>
iterator Dune::ParallelIndexSet< TG, TL, N >::begin (  )  [inline, inherited]
template<typename TG, typename TL, int N = 100>
void Dune::ParallelIndexSet< TG, TL, N >::beginResize (  )  throw (InvalidIndexSetState) [inherited]

Indicate that the index set is to be resized.

Exceptions:
InvalidState If index set was not in ParallelIndexSetState::GROUND mode.

Referenced by Dune::OwnerOverlapCopyCommunication< GlobalIdType, LocalIdType >::OwnerOverlapCopyCommunication(), and Dune::IndicesSyncer< T >::sync().

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

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 specfic flag. This functionality is for example provided the classes EnumItem, EnumRange and Combine.

Parameters:
remoteIndices The indices known to remote processes.
sourceFlags The set of flags marking indices we send from.
destFlags The set of flags marking indices we receive for.

References Dune::Interface::communicator_, and Dune::Interface::strip().

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

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 [inline, protected, inherited]

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 specfic 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:
remoteIndices The indices known to remote processes.
sourceFlags The set of flags marking source indices.
destFlags The setof flags markig destination indices.
functor A 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);
MPI_Comm Dune::Interface::communicator (  )  const [inline, inherited]

Get the MPI Communicator.

References Dune::Interface::communicator_.

Referenced by Dune::Interface::print().

Dune::ComposeMPITraits ( long  double,
MPI_LONG_DOUBLE   
)
Dune::ComposeMPITraits ( double  ,
MPI_DOUBLE   
)
Dune::ComposeMPITraits ( float  ,
MPI_FLOAT   
)
Dune::ComposeMPITraits ( unsigned  long,
MPI_UNSIGNED_LONG   
)
Dune::ComposeMPITraits ( long  ,
MPI_LONG   
)
Dune::ComposeMPITraits ( unsigned  int,
MPI_UNSIGNED   
)
Dune::ComposeMPITraits ( int  ,
MPI_INT   
)
Dune::ComposeMPITraits ( unsigned  short,
MPI_UNSIGNED_SHORT   
)
Dune::ComposeMPITraits ( short  ,
MPI_SHORT   
)
Dune::ComposeMPITraits ( unsigned  char,
MPI_UNSIGNED_CHAR   
)
Dune::ComposeMPITraits ( char  ,
MPI_CHAR   
)
template<typename TS , typename TG , typename TL , int N>
SelectionIterator< TS, TG, TL, N > Dune::UncachedSelection< TS, TG, TL, N >::end (  )  const [inline, inherited]

Get an iterator over the selected indices.

Returns:
An iterator positioned at the first selected index.

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

template<typename TS , typename TG , typename TL , int N>
uint32_t * Dune::Selection< TS, TG, TL, N >::end (  )  const [inline, inherited]

Get an iterator over the selected indices.

Returns:
An iterator positioned at the first selected index.

Referenced by Dune::Selection< TS, TG, TL, N >::setIndexSet().

template<class I>
const_iterator Dune::GlobalLookupIndexSet< I >::end (  )  const [inline, inherited]

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

Returns:
Iterator over the local indices.
template<typename TG, typename TL, int N = 100>
const_iterator Dune::ParallelIndexSet< TG, TL, N >::end (  )  const [inline, inherited]

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

Returns:
Iterator over the local indices.
template<typename TG, typename TL, int N = 100>
iterator Dune::ParallelIndexSet< TG, TL, N >::end (  )  [inline, inherited]
template<typename TG, typename TL, int N = 100>
void Dune::ParallelIndexSet< TG, TL, N >::endResize (  )  throw (InvalidIndexSetState) [inherited]

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:
InvalidState If index set was not in ParallelIndexSetState::RESIZE mode.

Referenced by Dune::OwnerOverlapCopyCommunication< GlobalIdType, LocalIdType >::OwnerOverlapCopyCommunication(), and Dune::IndicesSyncer< T >::sync().

template<typename TS , typename TG , typename TL , int N>
void Dune::Selection< TS, TG, TL, N >::free (  )  [inline, inherited]
void Dune::Interface::free (  )  [inline, inherited]
template<typename T >
MPI_Datatype Dune::MPITraits< ParallelLocalIndex< T > >::getType (  )  [inline, static, inherited]

References type.

template<typename T1 , typename T2 >
MPI_Datatype Dune::MPITraits< std::pair< T1, T2 > >::getType (  )  [inline, static, inherited]

References type.

template<int k>
static MPI_Datatype Dune::MPITraits< bigunsignedint< k > >::getType (  )  [inline, static, inherited]
template<class K , int n>
static MPI_Datatype Dune::MPITraits< FieldVector< K, n > >::getType (  )  [inline, static, inherited]
template<class TG, class TL>
const GlobalIndex& Dune::IndexPair< TG, TL >::global (  )  const [inline, inherited]
template<typename T >
IndicesSyncer< T >::GlobalIndex & Dune::IndicesSyncer< T >::Iterators::globalIndex (  )  const [inline, inherited]

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

Returns:
The current global index.
template<class I>
Dune::GlobalLookupIndexSet< I >::GlobalLookupIndexSet ( const ParallelIndexSet indexset  )  [inherited]

Constructor.

Parameters:
indexset The index set we want to be able to lookup the corresponding global index of a local index.
template<class I>
Dune::GlobalLookupIndexSet< I >::GlobalLookupIndexSet ( const ParallelIndexSet indexset,
std::size_t  size 
) [inherited]

Constructor.

Parameters:
indexset The index set we want to be able to lookup the corresponding global index of a local index.
size The number of indices present, i.e. one more than the maximum local index.
template<class TG, class TL>
Dune::IndexPair< TG, TL >::IndexPair ( const GlobalIndex global  )  [inherited]

Constructs a new Pair.

The local index will be 0.

Parameters:
global The global index.
template<class TG, class TL>
Dune::IndexPair< TG, TL >::IndexPair (  )  [inherited]

Construct a new Pair.

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

Constructs a new Pair.

Parameters:
global The global index.
local The local index.
template<typename T >
Dune::IndicesSyncer< T >::IndicesSyncer ( ParallelIndexSet indexSet,
RemoteIndices remoteIndices 
) [inline, inherited]

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:
indexSet The index set with the information of the locally present indices.
remoteIndices The remoteIndices.

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

template<typename T >
void Dune::IndicesSyncer< T >::Iterators::insert ( const RemoteIndex index,
const GlobalIndex global 
) [inline, inherited]

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

Parameters:
index The remote index.
global The global index corresponding to the remote index.
std::map< int, std::pair< InterfaceInformation, InterfaceInformation > > & Dune::Interface::interfaces (  )  [inline, protected, inherited]

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).
const std::map< int, std::pair< InterfaceInformation, InterfaceInformation > > & Dune::Interface::interfaces (  )  const [inline, inherited]

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).

Referenced by std::operator<<().

template<typename T >
bool Dune::IndicesSyncer< T >::Iterators::isAtEnd (  )  const [inline, inherited]

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.
template<typename T >
bool Dune::IndicesSyncer< T >::Iterators::isNotAtEnd (  )  const [inline, inherited]

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.
template<typename T >
bool Dune::IndicesSyncer< T >::Iterators::isOld (  )  const [inline, inherited]

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.
template<class T >
bool Dune::ParallelLocalIndex< T >::isPublic (  )  const [inline, inherited]

Check whether the index might also be known other processes.

Returns:
True if the index might be known to other processors.
template<typename TG, typename TL, int N = 100>
Dune::ParallelIndexSet< TG, TL, N >::iterator::iterator ( const iterator other  )  [inline, inherited]
template<typename TG, typename TL, int N = 100>
Dune::ParallelIndexSet< TG, TL, N >::iterator::iterator ( ParallelIndexSet< TG, TL, N > &  indexSet,
const Father &  father 
) [inline, inherited]
template<typename T >
Dune::IndicesSyncer< T >::Iterators::Iterators (  )  [inline, inherited]

Default constructor.

template<typename T >
Dune::IndicesSyncer< T >::Iterators::Iterators ( RemoteIndexList &  remoteIndices,
GlobalIndexList &  globalIndices,
BoolList &  booleans 
) [inline, inherited]

Constructor.

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

Parameters:
remoteIndices The list of the remote indices.
globalIndices The 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.
booleans Whether the remote index was there before the sync process started.
template<class T >
size_t Dune::ParallelLocalIndex< T >::local (  )  const [inline, inherited]

get the local index.

Returns:
The local index.
const std::size_t & Dune::LocalIndex::local (  )  const [inline, inherited]

get the local index.

Returns:
The local index.
template<class TG, class TL>
const LocalIndex& Dune::IndexPair< TG, TL >::local (  )  const [inline, inherited]

Get the local index.

Returns:
The local index.
template<class TG, class TL>
LocalIndex& Dune::IndexPair< TG, TL >::local (  )  [inline, inherited]
template<typename TG, typename TL, int N = 100>
void Dune::ParallelIndexSet< TG, TL, N >::markAsDeleted ( const iterator position  )  throw (InvalidIndexSetState) [inline, inherited]

Mark an index as deleted.

The index will be deleted during endResize().

Parameters:
position An iterator at the position we want to delete.
Exceptions:
InvalidState If index set is not in ParallelIndexSetState::RESIZE mode.
template<class T >
Dune::ParallelLocalIndex< T >::operator size_t (  )  const [inline, inherited]

Convert to the local index represented by an int.

Dune::LocalIndex::operator std::size_t (  )  const [inline, inherited]

Convert to the local index represented by an int.

template<class TG , class TL >
bool Dune::operator!= ( const IndexPair< TG, TL > &  a,
const TG &  b 
) [inline]
template<class TG , class TL >
bool Dune::operator!= ( const IndexPair< TG, TL > &  a,
const IndexPair< TG, TL > &  b 
) [inline]
template<typename T >
IndicesSyncer< T >::Iterators & Dune::IndicesSyncer< T >::Iterators::operator++ (  )  [inline, inherited]

Increment all iteraors.

template<class TG , class TL >
bool Dune::operator< ( const IndexPair< TG, TL > &  a,
const TG &  b 
) [inline]
template<class TG , class TL >
bool Dune::operator< ( const IndexPair< TG, TL > &  a,
const IndexPair< TG, TL > &  b 
) [inline]
template<class T , class A >
std::ostream & Dune::operator<< ( std::ostream &  os,
const RemoteIndices< T, A > &  indices 
) [inline]
template<typename T1 , typename T2 >
std::ostream& Dune::operator<< ( std::ostream &  os,
const RemoteIndex< T1, T2 > &  index 
) [inline]
template<class T >
std::ostream& Dune::operator<< ( std::ostream &  os,
const ParallelLocalIndex< T > &  index 
) [inline]

Print the local index to a stream.

Parameters:
os The output stream to print to.
index The index to print.

References index.

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:
os The outputstream to print to.
indexSet The index set to print.

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

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

Print an index pair.

Parameters:
os The outputstream to print to.
pair The index pair to print.
template<class TG , class TL >
bool Dune::operator<= ( const IndexPair< TG, TL > &  a,
const TG &  b 
) [inline]
template<class TG , class TL >
bool Dune::operator<= ( const IndexPair< TG, TL > &  a,
const IndexPair< TG, TL > &  b 
) [inline]
template<class T >
ParallelLocalIndex< T > & Dune::ParallelLocalIndex< T >::operator= ( size_t  index  )  [inline, inherited]

Assign a new local index.

Parameters:
index The new local index.
LocalIndex & Dune::LocalIndex::operator= ( std::size_t  index  )  [inline, inherited]

Assign a new local index.

Parameters:
index The new local index.
template<typename TG, typename TL, int N = 100>
iterator& Dune::ParallelIndexSet< TG, TL, N >::iterator::operator== ( const iterator other  )  [inline, inherited]
template<class TG , class TL >
bool Dune::operator== ( const IndexPair< TG, TL > &  a,
const TG &  b 
) [inline]
template<class TG , class TL >
bool Dune::operator== ( const IndexPair< TG, TL > &  a,
const IndexPair< TG, TL > &  b 
) [inline]
template<class TG , class TL >
bool Dune::operator> ( const IndexPair< TG, TL > &  a,
const TG &  b 
) [inline]
template<class TG , class TL >
bool Dune::operator> ( const IndexPair< TG, TL > &  a,
const IndexPair< TG, TL > &  b 
) [inline]
template<class TG , class TL >
bool Dune::operator>= ( const IndexPair< TG, TL > &  a,
const TG &  b 
) [inline]
template<class TG , class TL >
bool Dune::operator>= ( const IndexPair< TG, TL > &  a,
const IndexPair< TG, TL > &  b 
) [inline]
template<class I>
const IndexPair& Dune::GlobalLookupIndexSet< I >::operator[] ( const GlobalIndex global  )  const [inline, inherited]

Find the index pair with a specific global id.

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

Parameters:
global The globally unique id of the pair.
Returns:
The pair of indices for the id.
Exceptions:
RangeError Thrown if the global id is not known.
template<typename TG, typename TL, int N = 100>
const IndexPair& Dune::ParallelIndexSet< TG, TL, N >::operator[] ( const GlobalIndex global  )  const [inline, inherited]

Find the index pair with a specific global id.

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

Parameters:
global The 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.
template<typename TG, typename TL, int N = 100>
IndexPair& Dune::ParallelIndexSet< TG, TL, N >::operator[] ( const GlobalIndex global  )  [inline, inherited]

Find the index pair with a specific global id.

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

Parameters:
global The 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.
template<class I>
const IndexPair* Dune::GlobalLookupIndexSet< I >::pair ( const std::size_t &  local  )  const [inline, inherited]
template<typename TG, typename TL, int N = 100>
Dune::ParallelIndexSet< TG, TL, N >::ParallelIndexSet (  )  [inherited]

Constructor.

template<class T >
Dune::ParallelLocalIndex< T >::ParallelLocalIndex (  )  [inline, inherited]

Parameterless constructor.

Needed for use in container classes.

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

Constructor.

Parameters:
localIndex The local index.
attribute The attribute of the index.
isPublic True if the index might also be known to other processes.
template<class T >
Dune::ParallelLocalIndex< T >::ParallelLocalIndex ( const Attribute attribute,
bool  isPublic 
) [inline, inherited]

Constructor.

The local index will be initialized to 0.

Parameters:
attribute The attribute of the index.
isPublic True if the index might also be known to other processes.
void Dune::Interface::print (  )  const [inline, inherited]

Print the interface to std::out for debugging.

References Dune::Interface::communicator(), and Dune::InterfaceInformation::size().

template<typename T >
IndicesSyncer< T >::RemoteIndex & Dune::IndicesSyncer< T >::Iterators::remoteIndex (  )  const [inline, inherited]

Get the remote index at current position.

Returns:
The current remote index.
template<typename TG, typename TL, int N = 100>
void Dune::ParallelIndexSet< TG, TL, N >::renumberLocal (  )  [inline, inherited]

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 arbituary index pairs with $g_i<g_j$ then after renumbering $l_i<l_j$ will hold.

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

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

Parameters:
globalMap The map of the process number to the list of global indices corresponding to the remote index list of the process.
remoteIndices The known remote indices.
indexSet The set of local indices of the current process.

References index, and Dune::ParallelIndexSet< TG, TL, N >::seqNo().

Referenced by Dune::fillIndexSetHoles().

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

Reset all the underlying iterators.

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

Parameters:
remoteIndices The list of the remote indices.
globalIndices The 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.
booleans Whether the remote index was there before the sync process started.
template<class I>
int Dune::GlobalLookupIndexSet< I >::seqNo (  )  const [inline, inherited]

Get the internal sequence number.

Is initially 0 is incremented for each resize.

Returns:
The sequence number.
template<typename TG, typename TL, int N = 100>
int Dune::ParallelIndexSet< TG, TL, N >::seqNo (  )  const [inline, inherited]

Get the internal sequence number.

Is initially 0 is incremented for each resize.

Returns:
The sequence number.

Referenced by Dune::OwnerOverlapCopyCommunication< GlobalIdType, LocalIdType >::buildGlobalLookup(), Dune::repairLocalIndexPointers(), and Dune::IndicesSyncer< T >::sync().

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

Set the attribute of the index.

Parameters:
attribute The associated attribute.
template<typename TS , typename TG , typename TL , int N>
void Dune::UncachedSelection< TS, TG, TL, N >::setIndexSet ( const ParallelIndexSet indexset  )  [inline, inherited]

Set the index set of the selection.

Parameters:
indexset The index set to use.
template<typename TS , typename TG , typename TL , int N>
void Dune::Selection< TS, TG, TL, N >::setIndexSet ( const ParallelIndexSet indexset  )  [inline, inherited]
template<class TG, class TL>
void Dune::IndexPair< TG, TL >::setLocal ( int  index  )  [inline, inherited]

Set the local index.

Parameters:
index The index to set it to.
template<class T >
void Dune::ParallelLocalIndex< T >::setState ( const LocalIndexState state  )  [inline, inherited]

Set the state.

Parameters:
state The state to set.
void Dune::LocalIndex::setState ( LocalIndexState  state  )  [inline, inherited]

Set the state.

Parameters:
state The state to set.
template<class I>
size_t Dune::GlobalLookupIndexSet< I >::size (  )  const [inline, inherited]

Get the total number (public and nonpublic) indices.

Returns:
The total number (public and nonpublic) indices.

Referenced by Dune::fillIndexSetHoles().

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

Get the total number (public and nonpublic) indices.

Returns:
The total number (public and nonpublic) indices.

Referenced by Dune::graphRepartition().

template<class T >
LocalIndexState Dune::ParallelLocalIndex< T >::state (  )  const [inline, inherited]

Get the state.

Returns:
The state.
LocalIndexState Dune::LocalIndex::state (  )  const [inline, inherited]

Get the state.

Returns:
The state.
template<typename TG, typename TL, int N = 100>
const ParallelIndexSetState& Dune::ParallelIndexSet< TG, TL, N >::state (  )  [inline, inherited]

Get the state the index set is in.

Returns:
The state of the index set.

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

template<typename T , typename A , typename A1 >
void Dune::storeGlobalIndicesOfRemoteIndices ( std::map< int, SLList< typename T::GlobalIndex, A > > &  globalMap,
const RemoteIndices< T, A1 > &  remoteIndices,
const T &  indexSet 
) [inline]

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:
globalMap Map to store the corresponding global indices in.
remoteIndices The remote index information we need to store the corresponding global indices of.
indexSet The index set that is for both the sending and receiving side of the remote index information.

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

Referenced by Dune::fillIndexSetHoles().

void Dune::Interface::strip (  )  [inline, inherited]

Referenced by Dune::Interface::build().

template<typename T >
template<typename T1 >
void Dune::IndicesSyncer< T >::sync ( T1 &  numberer  )  [inline, inherited]

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:
numberer Functor 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::RemoteIndices< T, A >::begin(), Dune::ParallelIndexSet< TG, TL, N >::beginResize(), Dune::RemoteIndices< T, A >::end(), Dune::ParallelIndexSet< TG, TL, N >::endResize(), index, Dune::RemoteIndices< T, A >::neighbours(), and Dune::ParallelIndexSet< TG, TL, N >::seqNo().

template<typename T >
void Dune::IndicesSyncer< T >::sync (  )  [inline, inherited]

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()

template<class I>
Dune::GlobalLookupIndexSet< I >::~GlobalLookupIndexSet (  )  [inherited]

Destructor.

Dune::Interface::~Interface (  )  [inline, virtual, inherited]

Destructor.

References Dune::Interface::free().

template<typename TS , typename TG , typename TL , int N>
Dune::Selection< TS, TG, TL, N >::~Selection (  )  [inline, inherited]

Variable Documentation

template<class K , int n>
MPI_Datatype Dune::MPITraits< FieldVector< K, n > >::datatype = MPI_DATATYPE_NULL [inline, static, inherited]
template<typename T>
int Dune::IndicesSyncer< T >::MessageInformation::pairs [inherited]

The number of pairs (attribute and process number) we publish to the neighbour process.

template<typename T>
int Dune::IndicesSyncer< T >::MessageInformation::publish [inherited]

The number of indices we publish for the other process.

template<int k>
MPI_Datatype Dune::MPITraits< bigunsignedint< k > >::vectortype = MPI_DATATYPE_NULL [inline, static, inherited]
template<class K , int n>
MPI_Datatype Dune::MPITraits< FieldVector< K, n > >::vectortype = {MPI_DATATYPE_NULL} [inline, static, inherited]

Friends

template<class TG, class TL>
friend class MPITraits< IndexPair< TG, TL > > [friend, inherited]
template<class TG, class TL>
bool operator!= ( const IndexPair< TG, TL > &  a,
const TG &  b 
) [friend, inherited]
template<class TG, class TL>
bool operator!= ( const IndexPair< TG, TL > &  a,
const IndexPair< TG, TL > &  b 
) [friend, inherited]
template<class TG, class TL>
bool operator< ( const IndexPair< TG, TL > &  a,
const TG &  b 
) [friend, inherited]
template<class TG, class TL>
bool operator< ( const IndexPair< TG, TL > &  a,
const IndexPair< TG, TL > &  b 
) [friend, inherited]
template<class TG, class TL>
bool operator<= ( const IndexPair< TG, TL > &  a,
const TG &  b 
) [friend, inherited]
template<class TG, class TL>
bool operator<= ( const IndexPair< TG, TL > &  a,
const IndexPair< TG, TL > &  b 
) [friend, inherited]
template<class TG, class TL>
bool operator== ( const IndexPair< TG, TL > &  a,
const TG &  b 
) [friend, inherited]
template<class TG, class TL>
bool operator== ( const IndexPair< TG, TL > &  a,
const IndexPair< TG, TL > &  b 
) [friend, inherited]
template<class TG, class TL>
bool operator> ( const IndexPair< TG, TL > &  a,
const TG &  b 
) [friend, inherited]
template<class TG, class TL>
bool operator> ( const IndexPair< TG, TL > &  a,
const IndexPair< TG, TL > &  b 
) [friend, inherited]
template<class TG, class TL>
bool operator>= ( const IndexPair< TG, TL > &  a,
const TG &  b 
) [friend, inherited]
template<class TG, class TL>
bool operator>= ( const IndexPair< TG, TL > &  a,
const IndexPair< TG, TL > &  b 
) [friend, inherited]
template<typename TG, typename TL, int N = 100>
friend class ParallelIndexSet< GlobalIndex, LocalIndex, N > [friend, inherited]
Generated on Sat Apr 24 11:13:48 2010 for dune-istl by  doxygen 1.6.3