Reference documentation for deal.II version 9.5.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Public Member Functions | Private Attributes | List of all members
Utilities::MPI::NoncontiguousPartitioner Class Reference

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

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

Public Member Functions

 NoncontiguousPartitioner ()=default
 
 NoncontiguousPartitioner (const IndexSet &indexset_locally_owned, const IndexSet &indexset_ghost, const MPI_Comm communicator)
 
 NoncontiguousPartitioner (const std::vector< types::global_dof_index > &indices_locally_owned, const std::vector< types::global_dof_index > &indices_ghost, const MPI_Comm communicator)
 
template<typename Number >
void export_to_ghosted_array (const ArrayView< const Number > &locally_owned_array, const ArrayView< Number > &ghost_array) const
 
template<typename Number >
void export_to_ghosted_array (const unsigned int communication_channel, const ArrayView< const Number > &locally_owned_array, const ArrayView< Number > &temporary_storage, const ArrayView< Number > &ghost_array, std::vector< MPI_Request > &requests) const
 
template<typename Number >
void export_to_ghosted_array_start (const unsigned int communication_channel, const ArrayView< const Number > &locally_owned_array, const ArrayView< Number > &temporary_storage, std::vector< MPI_Request > &requests) const
 
template<typename Number >
void export_to_ghosted_array_finish (const ArrayView< const Number > &temporary_storage, const ArrayView< Number > &ghost_array, std::vector< MPI_Request > &requests) const
 
std::pair< unsigned int, unsigned intn_targets () const
 
unsigned int temporary_storage_size () const
 
types::global_dof_index memory_consumption ()
 
MPI_Comm get_mpi_communicator () const override
 
void reinit (const IndexSet &locally_owned_indices, const IndexSet &ghost_indices, const MPI_Comm communicator) override
 
void reinit (const std::vector< types::global_dof_index > &locally_owned_indices, const std::vector< types::global_dof_index > &ghost_indices, const MPI_Comm communicator)
 

Private Attributes

MPI_Comm communicator
 
std::vector< unsigned intsend_ranks
 
std::vector< types::global_dof_indexsend_ptr
 
std::vector< types::global_dof_indexsend_indices
 
std::vector< unsigned intrecv_ranks
 
std::vector< types::global_dof_indexrecv_ptr
 
std::vector< types::global_dof_indexrecv_indices
 
std::vector< std::uint8_t > buffers
 
std::vector< MPI_Request > requests
 

Detailed Description

A flexible Partitioner class, which does not impose restrictions regarding the order of the underlying index sets. In other words, this class implements the interface of the Utilities::MPI::CommunicationPatternBase base class with no assumption that every process stores a contiguous part of the array of objects, but that indeed the locally owned indices can be an arbitrary subset of all indices of elements of the array to which they refer.

If you want to store only contiguous parts of these arrays on each process, take a look at Utilities::MPI::Partitioner.

Definition at line 47 of file mpi_noncontiguous_partitioner.h.

Constructor & Destructor Documentation

◆ NoncontiguousPartitioner() [1/3]

Utilities::MPI::NoncontiguousPartitioner::NoncontiguousPartitioner ( )
default

Default constructor. Requires calling one of the reinit() functions to create a valid object.

◆ NoncontiguousPartitioner() [2/3]

Utilities::MPI::NoncontiguousPartitioner::NoncontiguousPartitioner ( const IndexSet indexset_locally_owned,
const IndexSet indexset_ghost,
const MPI_Comm  communicator 
)

Constructor. Set up point-to-point communication pattern based on the IndexSets arguments indexset_locally_owned and indexset_ghost for the MPI communicator communicator.

Definition at line 32 of file mpi_noncontiguous_partitioner.cc.

◆ NoncontiguousPartitioner() [3/3]

Utilities::MPI::NoncontiguousPartitioner::NoncontiguousPartitioner ( const std::vector< types::global_dof_index > &  indices_locally_owned,
const std::vector< types::global_dof_index > &  indices_ghost,
const MPI_Comm  communicator 
)

Constructor. Same as above but for vectors of indices indices_locally_owned and indices_ghost. This allows the indices to not be sorted and the values are read and written automatically at the right position of the vector during update_values(), update_values_start(), and update_values_finish(). It is allowed to include entries with the value numbers::invalid_dof_index which do not take part of the index exchange but are present in the data vectors as padding.

Definition at line 42 of file mpi_noncontiguous_partitioner.cc.

Member Function Documentation

◆ export_to_ghosted_array() [1/2]

template<typename Number >
void Utilities::MPI::NoncontiguousPartitioner::export_to_ghosted_array ( const ArrayView< const Number > &  locally_owned_array,
const ArrayView< Number > &  ghost_array 
) const

Fill the vector ghost_array according to the precomputed communication pattern with values from locally_owned_array.

Precondition
The vectors only have to provide a method begin(), which allows to access their raw data.
The size of both vectors must be at least as large as the number of entries in the index sets passed to the constructors or the reinit() functions.
Note
This function calls the methods update_values_start() and update_values_finish() in sequence. Users can call these two functions separately and hereby overlap communication and computation.

◆ export_to_ghosted_array() [2/2]

template<typename Number >
void Utilities::MPI::NoncontiguousPartitioner::export_to_ghosted_array ( const unsigned int  communication_channel,
const ArrayView< const Number > &  locally_owned_array,
const ArrayView< Number > &  temporary_storage,
const ArrayView< Number > &  ghost_array,
std::vector< MPI_Request > &  requests 
) const

Same as above but with an interface similar to Utilities::MPI::Partitioner::export_to_ghosted_array_start and Utilities::MPI::Partitioner::export_to_ghosted_array_finish. In this function, the user can provide the temporary data structures to be used.

Precondition
The size of the temporary_storage vector has to be at least temporary_storage_size. The reason for this is that this vector is used as buffer for both sending and receiving data.
Note
Any value less than 10 is a valid value of communication_channel.

◆ export_to_ghosted_array_start()

template<typename Number >
void Utilities::MPI::NoncontiguousPartitioner::export_to_ghosted_array_start ( const unsigned int  communication_channel,
const ArrayView< const Number > &  locally_owned_array,
const ArrayView< Number > &  temporary_storage,
std::vector< MPI_Request > &  requests 
) const

Start update: Data is packed, non-blocking send and receives are started.

Note
In contrast to the function Utilities::MPI::Partitioner::export_to_ghosted_array_start, the user does not pass a reference to the destination vector, since the data is received into a designated part of the buffer temporary_storage. This allows for padding and other post-processing of the received data.
Precondition
The required size of the vectors are the same as in the functions above.
Note
Any value less than 10 is a valid value of communication_channel.

◆ export_to_ghosted_array_finish()

template<typename Number >
void Utilities::MPI::NoncontiguousPartitioner::export_to_ghosted_array_finish ( const ArrayView< const Number > &  temporary_storage,
const ArrayView< Number > &  ghost_array,
std::vector< MPI_Request > &  requests 
) const

Finish update. The method waits until all data has been sent and received. Once data from any process is received it is processed and placed at the right position of the vector dst.

Note
In contrast to the function Utilities::MPI::Partitioner::export_to_ghosted_array_finish, the user also has to pass a reference to the buffer temporary_storage, since the data has been received into the buffer and not into the destination vector.
Precondition
The required size of the vectors are the same as in the functions above.

◆ n_targets()

std::pair< unsigned int, unsigned int > Utilities::MPI::NoncontiguousPartitioner::n_targets ( ) const

Returns the number of processes this process sends data to and the number of processes this process receives data from.

Definition at line 53 of file mpi_noncontiguous_partitioner.cc.

◆ temporary_storage_size()

unsigned int Utilities::MPI::NoncontiguousPartitioner::temporary_storage_size ( ) const

Return the size of the temporary storage needed by the export_to_ghosted_array() functions, if the temporary storage is handled by the user code.

Definition at line 61 of file mpi_noncontiguous_partitioner.cc.

◆ memory_consumption()

types::global_dof_index Utilities::MPI::NoncontiguousPartitioner::memory_consumption ( )

Return memory consumption in Byte.

Definition at line 69 of file mpi_noncontiguous_partitioner.cc.

◆ get_mpi_communicator()

MPI_Comm Utilities::MPI::NoncontiguousPartitioner::get_mpi_communicator ( ) const
overridevirtual

Return the underlying MPI communicator.

Implements Utilities::MPI::CommunicationPatternBase.

Definition at line 84 of file mpi_noncontiguous_partitioner.cc.

◆ reinit() [1/2]

void Utilities::MPI::NoncontiguousPartitioner::reinit ( const IndexSet locally_owned_indices,
const IndexSet ghost_indices,
const MPI_Comm  communicator 
)
overridevirtual

Reinitialize the communication pattern.

Parameters
[in]locally_owned_indicesThe set of indices of elements in the array mentioned in the class documentation that are stored on the current process.
[in]ghost_indicesThe set of indices of elements in the array mentioned in the class documentation that the current process will need to be able to import.
[in]communicatorThe MPI communicator used to describe the entire set of processes that participate in the storage and access to elements of the array.

Implements Utilities::MPI::CommunicationPatternBase.

Definition at line 92 of file mpi_noncontiguous_partitioner.cc.

◆ reinit() [2/2]

void Utilities::MPI::NoncontiguousPartitioner::reinit ( const std::vector< types::global_dof_index > &  locally_owned_indices,
const std::vector< types::global_dof_index > &  ghost_indices,
const MPI_Comm  communicator 
)

Initialize the inner data structures using explicit sets of indices. See the documentation of the other reinit() function for what the function does.

Definition at line 169 of file mpi_noncontiguous_partitioner.cc.

Member Data Documentation

◆ communicator

MPI_Comm Utilities::MPI::NoncontiguousPartitioner::communicator
private

MPI communicator.

Definition at line 216 of file mpi_noncontiguous_partitioner.h.

◆ send_ranks

std::vector<unsigned int> Utilities::MPI::NoncontiguousPartitioner::send_ranks
private

The ranks this process sends data to.

Definition at line 221 of file mpi_noncontiguous_partitioner.h.

◆ send_ptr

std::vector<types::global_dof_index> Utilities::MPI::NoncontiguousPartitioner::send_ptr
private

Offset of each process within send_buffer.

Note
Together with send_indices this forms a CRS data structure.

Definition at line 228 of file mpi_noncontiguous_partitioner.h.

◆ send_indices

std::vector<types::global_dof_index> Utilities::MPI::NoncontiguousPartitioner::send_indices
private

Local index of each entry in send_buffer within the destination vector.

Note
Together with send_ptr this forms a CRS data structure.

Definition at line 236 of file mpi_noncontiguous_partitioner.h.

◆ recv_ranks

std::vector<unsigned int> Utilities::MPI::NoncontiguousPartitioner::recv_ranks
private

The ranks this process receives data from.

Definition at line 241 of file mpi_noncontiguous_partitioner.h.

◆ recv_ptr

std::vector<types::global_dof_index> Utilities::MPI::NoncontiguousPartitioner::recv_ptr
private

Offset of each process within recv_buffer.

Note
Together with recv_indices this forms a CRS data structure.

Definition at line 248 of file mpi_noncontiguous_partitioner.h.

◆ recv_indices

std::vector<types::global_dof_index> Utilities::MPI::NoncontiguousPartitioner::recv_indices
private

Local index of each entry in recv_buffer within the destination vector.

Note
Together with recv_ptr this forms a CRS data structure.

Definition at line 256 of file mpi_noncontiguous_partitioner.h.

◆ buffers

std::vector<std::uint8_t> Utilities::MPI::NoncontiguousPartitioner::buffers
mutableprivate

Buffer containing the values sorted by rank for sending and receiving.

Note
Only allocated if not provided externally by user.
At this place we do not know the type of the data to be sent. So we use an arbitrary type of size 1 byte. The type is cast to the requested type in the relevant functions.

Definition at line 267 of file mpi_noncontiguous_partitioner.h.

◆ requests

std::vector<MPI_Request> Utilities::MPI::NoncontiguousPartitioner::requests
mutableprivate

MPI requests for sending and receiving.

Note
Only allocated if not provided externally by user.

Definition at line 274 of file mpi_noncontiguous_partitioner.h.


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