Reference documentation for deal.II version 8.5.1
Public Member Functions | Static Public Member Functions | Private Attributes | List of all members
Utilities::MPI::Partitioner Class Reference

#include <deal.II/base/partitioner.h>

Inheritance diagram for Utilities::MPI::Partitioner:
[legend]

Public Member Functions

 Partitioner ()
 
 Partitioner (const unsigned int size)
 
 Partitioner (const IndexSet &locally_owned_indices, const IndexSet &ghost_indices_in, const MPI_Comm communicator_in)
 
 Partitioner (const IndexSet &locally_owned_indices, const MPI_Comm communicator_in)
 
virtual void reinit (const IndexSet &vector_space_vector_index_set, const IndexSet &read_write_vector_index_set, const MPI_Comm &communicator)
 
void set_owned_indices (const IndexSet &locally_owned_indices)
 
void set_ghost_indices (const IndexSet &ghost_indices)
 
types::global_dof_index size () const
 
unsigned int local_size () const
 
const IndexSetlocally_owned_range () const
 
std::pair< types::global_dof_index, types::global_dof_indexlocal_range () const
 
bool in_local_range (const types::global_dof_index global_index) const
 
unsigned int global_to_local (const types::global_dof_index global_index) const
 
types::global_dof_index local_to_global (const unsigned int local_index) const
 
bool is_ghost_entry (const types::global_dof_index global_index) const
 
const IndexSetghost_indices () const
 
unsigned int n_ghost_indices () const
 
const std::vector< std::pair< unsigned int, unsigned int > > & ghost_targets () const
 
const std::vector< std::pair< unsigned int, unsigned int > > & import_indices () const
 
unsigned int n_import_indices () const
 
const std::vector< std::pair< unsigned int, unsigned int > > & import_targets () const
 
bool is_compatible (const Partitioner &part) const
 
bool is_globally_compatible (const Partitioner &part) const
 
unsigned int this_mpi_process () const
 
unsigned int n_mpi_processes () const
 
const MPI_Comm & get_communicator () const 1
 
virtual const MPI_Comm & get_mpi_communicator () const
 
bool ghost_indices_initialized () const
 
std::size_t memory_consumption () const
 
- Public Member Functions inherited from LinearAlgebra::CommunicationPatternBase
virtual ~CommunicationPatternBase ()
 

Static Public Member Functions

static ::ExceptionBaseExcIndexNotPresent (types::global_dof_index arg1, unsigned int arg2)
 

Private Attributes

types::global_dof_index global_size
 
IndexSet locally_owned_range_data
 
std::pair< types::global_dof_index, types::global_dof_indexlocal_range_data
 
IndexSet ghost_indices_data
 
unsigned int n_ghost_indices_data
 
std::vector< std::pair< unsigned int, unsigned int > > ghost_targets_data
 
std::vector< std::pair< unsigned int, unsigned int > > import_indices_data
 
unsigned int n_import_indices_data
 
std::vector< std::pair< unsigned int, unsigned int > > import_targets_data
 
unsigned int my_pid
 
unsigned int n_procs
 
MPI_Comm communicator
 
bool have_ghost_indices
 

Detailed Description

This class defines a model for the partitioning of a vector (or, in fact, any linear data structure) among processors using MPI.

The partitioner stores the global vector size and the locally owned range as a half-open interval [lower, upper). Furthermore, it includes a structure for the point-to-point communication patterns. It allows the inclusion of ghost indices (i.e. indices that a current processor needs to have access to, but are owned by another process) through an IndexSet. In addition, it also stores the other processors' ghost indices belonging to the current processor, which are the indices where other processors might require information from. In a sense, these import indices form the dual of the ghost indices. This information is gathered once when constructing the partitioner, which obviates subsequent global communication steps when exchanging data.

The partitioner includes a mechanism for converting global to local and local to global indices. Internally, this class stores vector elements using the convention as follows: The local range is associated with local indices [0,local_size), and ghost indices are stored consecutively in [local_size, local_size + n_ghost_indices). The ghost indices are sorted according to their global index.

Author
Katharina Kormann, Martin Kronbichler, 2010, 2011

Definition at line 62 of file partitioner.h.

Constructor & Destructor Documentation

◆ Partitioner() [1/4]

Utilities::MPI::Partitioner::Partitioner ( )

Empty Constructor.

Definition at line 24 of file partitioner.cc.

◆ Partitioner() [2/4]

Utilities::MPI::Partitioner::Partitioner ( const unsigned int  size)

Constructor with size argument. Creates an MPI_COMM_SELF structure where there is no real parallel layout.

Definition at line 39 of file partitioner.cc.

◆ Partitioner() [3/4]

Utilities::MPI::Partitioner::Partitioner ( const IndexSet locally_owned_indices,
const IndexSet ghost_indices_in,
const MPI_Comm  communicator_in 
)

Constructor with index set arguments. This constructor creates a distributed layout based on a given communicators, an IndexSet describing the locally owned range and another one for describing ghost indices that are owned by other processors, but we need to have read or write access to.

Definition at line 58 of file partitioner.cc.

◆ Partitioner() [4/4]

Utilities::MPI::Partitioner::Partitioner ( const IndexSet locally_owned_indices,
const MPI_Comm  communicator_in 
)

Constructor with one index set argument. This constructor creates a distributed layout based on a given communicator, and an IndexSet describing the locally owned range. It allows to set the ghost indices at a later time. Apart from this, it is similar to the other constructor with two index sets.

Definition at line 76 of file partitioner.cc.

Member Function Documentation

◆ reinit()

void Utilities::MPI::Partitioner::reinit ( const IndexSet vector_space_vector_index_set,
const IndexSet read_write_vector_index_set,
const MPI_Comm &  communicator 
)
virtual

Reinitialize the communication pattern. The first argument vector_space_vector_index_set is the index set associated to a VectorSpaceVector object. The second argument read_write_vector_index_set is the index set associated to a ReadWriteVector object.

Implements LinearAlgebra::CommunicationPatternBase.

Definition at line 93 of file partitioner.cc.

◆ set_owned_indices()

void Utilities::MPI::Partitioner::set_owned_indices ( const IndexSet locally_owned_indices)

Set the locally owned indices. Used in the constructor.

Definition at line 106 of file partitioner.cc.

◆ set_ghost_indices()

void Utilities::MPI::Partitioner::set_ghost_indices ( const IndexSet ghost_indices)

Allows to set the ghost indices after the constructor has been called.

Definition at line 142 of file partitioner.cc.

◆ size()

types::global_dof_index Utilities::MPI::Partitioner::size ( ) const

Return the global size.

◆ local_size()

unsigned int Utilities::MPI::Partitioner::local_size ( ) const

Return the local size, i.e. local_range().second minus local_range().first.

◆ locally_owned_range()

const IndexSet& Utilities::MPI::Partitioner::locally_owned_range ( ) const

Return an IndexSet representation of the local range. This class only supports contiguous local ranges, so the IndexSet actually only consists of one single range of data, and is equivalent to the result of local_range().

◆ local_range()

std::pair<types::global_dof_index,types::global_dof_index> Utilities::MPI::Partitioner::local_range ( ) const

Return the local range. The returned pair consists of the index of the first element and the index of the element one past the last locally owned one.

◆ in_local_range()

bool Utilities::MPI::Partitioner::in_local_range ( const types::global_dof_index  global_index) const

Return true if the given global index is in the local range of this processor.

◆ global_to_local()

unsigned int Utilities::MPI::Partitioner::global_to_local ( const types::global_dof_index  global_index) const

Return the local index corresponding to the given global index. If the given global index is neither locally owned nor a ghost, an exception is thrown.

Note that the local index for locally owned indices is between 0 and local_size()-1, and the local index for ghosts is between local_size() and local_size()+n_ghost_indices()-1.

◆ local_to_global()

types::global_dof_index Utilities::MPI::Partitioner::local_to_global ( const unsigned int  local_index) const

Return the global index corresponding to the given local index.

Note that the local index for locally owned indices is between 0 and local_size()-1, and the local index for ghosts is between local_size() and local_size()+n_ghost_indices()-1.

◆ is_ghost_entry()

bool Utilities::MPI::Partitioner::is_ghost_entry ( const types::global_dof_index  global_index) const

Return whether the given global index is a ghost index on the present processor. Returns false for indices that are owned locally and for indices not present at all.

◆ ghost_indices()

const IndexSet& Utilities::MPI::Partitioner::ghost_indices ( ) const

Return an IndexSet representation of all ghost indices.

◆ n_ghost_indices()

unsigned int Utilities::MPI::Partitioner::n_ghost_indices ( ) const

Return the number of ghost indices. Same as ghost_indices().n_elements(), but cached for simpler access.

◆ ghost_targets()

const std::vector<std::pair<unsigned int, unsigned int> >& Utilities::MPI::Partitioner::ghost_targets ( ) const

Return a list of processors (first entry) and the number of degrees of freedom for the individual processor on the ghost elements present (second entry).

◆ import_indices()

const std::vector<std::pair<unsigned int, unsigned int> >& Utilities::MPI::Partitioner::import_indices ( ) const

The set of (local) indices that we are importing during compress(), i.e., others' ghosts that belong to the local range. Similar structure as in an IndexSet, but tailored to be iterated over, and some indices may be duplicates.

◆ n_import_indices()

unsigned int Utilities::MPI::Partitioner::n_import_indices ( ) const

Number of import indices, i.e., indices that are ghosts on other processors and we will receive data from.

◆ import_targets()

const std::vector<std::pair<unsigned int, unsigned int> >& Utilities::MPI::Partitioner::import_targets ( ) const

Return a list of processors (first entry) and the number of degrees of freedom for all the processors that data is obtained from (second entry), i.e., locally owned indices that are ghosts on other processors.

◆ is_compatible()

bool Utilities::MPI::Partitioner::is_compatible ( const Partitioner part) const

Check whether the given partitioner is compatible with the partitioner used for this vector. Two partitioners are compatible if they have the same local size and the same ghost indices. They do not necessarily need to be the same data field. This is a local operation only, i.e., if only some processors decide that the partitioning is not compatible, only these processors will return false, whereas the other processors will return true.

Definition at line 366 of file partitioner.cc.

◆ is_globally_compatible()

bool Utilities::MPI::Partitioner::is_globally_compatible ( const Partitioner part) const

Check whether the given partitioner is compatible with the partitioner used for this vector. Two partitioners are compatible if they have the same local size and the same ghost indices. They do not necessarily need to be the same data field. As opposed to is_compatible(), this method checks for compatibility among all processors and the method only returns true if the partitioner is the same on all processors.

This method performs global communication, so make sure to use it only in a context where all processors call it the same number of times.

Definition at line 392 of file partitioner.cc.

◆ this_mpi_process()

unsigned int Utilities::MPI::Partitioner::this_mpi_process ( ) const

Return the MPI ID of the calling processor. Cached to have simple access.

◆ n_mpi_processes()

unsigned int Utilities::MPI::Partitioner::n_mpi_processes ( ) const

Return the total number of MPI processor participating in the given partitioner. Cached to have simple access.

◆ get_communicator()

const MPI_Comm& Utilities::MPI::Partitioner::get_communicator ( ) const

Return the MPI communicator underlying the partitioner object.

◆ get_mpi_communicator()

virtual const MPI_Comm& Utilities::MPI::Partitioner::get_mpi_communicator ( ) const
virtual

Return the MPI communicator underlying the partitioner object.

Implements LinearAlgebra::CommunicationPatternBase.

◆ ghost_indices_initialized()

bool Utilities::MPI::Partitioner::ghost_indices_initialized ( ) const

Return whether ghost indices have been explicitly added as a ghost_indices argument. Only true if a reinit call or constructor provided that argument.

◆ memory_consumption()

std::size_t Utilities::MPI::Partitioner::memory_consumption ( ) const

Compute the memory consumption of this structure.

Definition at line 401 of file partitioner.cc.

Member Data Documentation

◆ global_size

types::global_dof_index Utilities::MPI::Partitioner::global_size
private

The global size of the vector over all processors

Definition at line 297 of file partitioner.h.

◆ locally_owned_range_data

IndexSet Utilities::MPI::Partitioner::locally_owned_range_data
private

The range of the vector that is stored locally.

Definition at line 302 of file partitioner.h.

◆ local_range_data

std::pair<types::global_dof_index,types::global_dof_index> Utilities::MPI::Partitioner::local_range_data
private

The range of the vector that is stored locally. Extracted from locally_owned_range for performance reasons.

Definition at line 308 of file partitioner.h.

◆ ghost_indices_data

IndexSet Utilities::MPI::Partitioner::ghost_indices_data
private

The set of indices to which we need to have read access but that are not locally owned

Definition at line 314 of file partitioner.h.

◆ n_ghost_indices_data

unsigned int Utilities::MPI::Partitioner::n_ghost_indices_data
private

Caches the number of ghost indices. It would be expensive to use ghost_indices.n_elements() to compute this.

Definition at line 320 of file partitioner.h.

◆ ghost_targets_data

std::vector<std::pair<unsigned int, unsigned int> > Utilities::MPI::Partitioner::ghost_targets_data
private

Contains information which processors my ghost indices belong to and how many those indices are

Definition at line 326 of file partitioner.h.

◆ import_indices_data

std::vector<std::pair<unsigned int, unsigned int> > Utilities::MPI::Partitioner::import_indices_data
private

The set of (local) indices that we are importing during compress(), i.e., others' ghosts that belong to the local range. Similar structure as in an IndexSet, but tailored to be iterated over, and some indices may be duplicates.

Definition at line 334 of file partitioner.h.

◆ n_import_indices_data

unsigned int Utilities::MPI::Partitioner::n_import_indices_data
private

Caches the number of ghost indices. It would be expensive to compute it by iterating over the import indices and accumulate them.

Definition at line 340 of file partitioner.h.

◆ import_targets_data

std::vector<std::pair<unsigned int, unsigned int> > Utilities::MPI::Partitioner::import_targets_data
private

The set of processors and length of data field which send us their ghost data

Definition at line 346 of file partitioner.h.

◆ my_pid

unsigned int Utilities::MPI::Partitioner::my_pid
private

The ID of the current processor in the MPI network

Definition at line 351 of file partitioner.h.

◆ n_procs

unsigned int Utilities::MPI::Partitioner::n_procs
private

The total number of processors active in the problem

Definition at line 356 of file partitioner.h.

◆ communicator

MPI_Comm Utilities::MPI::Partitioner::communicator
private

The MPI communicator involved in the problem

Definition at line 361 of file partitioner.h.

◆ have_ghost_indices

bool Utilities::MPI::Partitioner::have_ghost_indices
private

Stores whether the ghost indices have been explicitly set.

Definition at line 366 of file partitioner.h.


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