Reference documentation for deal.II version GIT relicensing-437-g81ec864850 2024-04-19 07:30:02+00:00
\(\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 | Public Attributes | Static Public Attributes | Private Member Functions | List of all members
Utilities::MPI::internal::ComputeIndexOwner::Dictionary Struct Reference

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

Public Member Functions

void reinit (const IndexSet &owned_indices, const MPI_Comm comm)
 
unsigned int dof_to_dict_rank (const types::global_dof_index i)
 
types::global_dof_index get_index_offset (const unsigned int rank)
 
unsigned int get_owning_rank_index (const unsigned int rank_in_owned_indices, const unsigned int guess=0)
 

Public Attributes

FlexibleIndexStorage actually_owning_ranks
 
std::vector< unsigned intactually_owning_rank_list
 
types::global_dof_index dofs_per_process
 
std::pair< types::global_dof_index, types::global_dof_indexlocal_range
 
types::global_dof_index locally_owned_size
 
types::global_dof_index size
 
unsigned int n_dict_procs_in_owned_indices
 
unsigned int stride_small_size
 

Static Public Attributes

static constexpr unsigned int range_minimum_grain_size = 64
 
static constexpr unsigned int sparsity_factor = 4
 

Private Member Functions

void partition (const IndexSet &owned_indices, const MPI_Comm comm)
 

Detailed Description

Dictionary class with basic partitioning in terms of a single interval of fixed size known to all MPI ranks for two-stage index lookup.

Definition at line 79 of file mpi_compute_index_owner_internal.h.

Member Function Documentation

◆ reinit()

void Utilities::MPI::internal::ComputeIndexOwner::Dictionary::reinit ( const IndexSet owned_indices,
const MPI_Comm  comm 
)

Set up the dictionary by computing the partitioning from the global size and sending the rank information on locally owned ranges to the owner of the dictionary part.

Definition at line 167 of file mpi_compute_index_owner_internal.cc.

◆ dof_to_dict_rank()

unsigned int Utilities::MPI::internal::ComputeIndexOwner::Dictionary::dof_to_dict_rank ( const types::global_dof_index  i)
inline

Translate a global dof index to the MPI rank in the dictionary using dofs_per_process. We multiply by stride_small_size to ensure a balance over the MPI ranks due to the grain size.

Definition at line 380 of file mpi_compute_index_owner_internal.h.

◆ get_index_offset()

types::global_dof_index Utilities::MPI::internal::ComputeIndexOwner::Dictionary::get_index_offset ( const unsigned int  rank)
inline

Given an MPI rank id of an arbitrary processor, return the index offset where the local range of that processor begins.

Definition at line 389 of file mpi_compute_index_owner_internal.h.

◆ get_owning_rank_index()

unsigned int Utilities::MPI::internal::ComputeIndexOwner::Dictionary::get_owning_rank_index ( const unsigned int  rank_in_owned_indices,
const unsigned int  guess = 0 
)
inline

Given the rank in the owned indices from actually_owning_ranks, this returns the index of the rank in the actually_owning_rank_list.

Definition at line 401 of file mpi_compute_index_owner_internal.h.

◆ partition()

void Utilities::MPI::internal::ComputeIndexOwner::Dictionary::partition ( const IndexSet owned_indices,
const MPI_Comm  comm 
)
private

Compute the partition from the global size of the index space and the number of ranks.

Definition at line 432 of file mpi_compute_index_owner_internal.cc.

Member Data Documentation

◆ range_minimum_grain_size

constexpr unsigned int Utilities::MPI::internal::ComputeIndexOwner::Dictionary::range_minimum_grain_size = 64
staticconstexpr

The minimum grain size for the intervals.

We choose to limit the smallest size an interval for the two-stage lookup can have with the following two conflicting goals in mind: On the one hand, we do not want intervals in the dictionary to become too short. For uneven distributions of unknowns (some ranks with several thousands of unknowns, others with none), the lookup DoFs -> dictionary then involves sending from one MPI rank to many other MPI ranks holding dictionary intervals, leading to an exceedingly high number of messages some ranks have to send. Also, fewer longer intervals are generally more efficient to look up. On the other hand, a range size too large leads to opposite effect of many messages that come into a particular dictionary owner in the lookup DoFs -> dictionary. With the current setting, we get at most 64 messages coming to a single MPI rank in case there is 1 dof per MPI rank, which is reasonably low. At the same time, uneven distributions up to factors of 4096 can be handled with at most 64 messages as well.

Definition at line 102 of file mpi_compute_index_owner_internal.h.

◆ sparsity_factor

constexpr unsigned int Utilities::MPI::internal::ComputeIndexOwner::Dictionary::sparsity_factor = 4
staticconstexpr

Factor that determines if an index set is sparse or not. An index set if sparse if less than 25% of the indices are owned by any process. If the index set is sparse, we switch the internal storage from a fast storage (vector) to a memory-efficient storage (map).

Definition at line 110 of file mpi_compute_index_owner_internal.h.

◆ actually_owning_ranks

FlexibleIndexStorage Utilities::MPI::internal::ComputeIndexOwner::Dictionary::actually_owning_ranks

A vector with as many entries as there are dofs in the dictionary of the current process, and each entry containing the rank of the owner of that dof in the IndexSet owned_indices. This is queried in the index lookup, so we keep an expanded list.

Definition at line 118 of file mpi_compute_index_owner_internal.h.

◆ actually_owning_rank_list

std::vector<unsigned int> Utilities::MPI::internal::ComputeIndexOwner::Dictionary::actually_owning_rank_list

A sorted vector containing the MPI ranks appearing in actually_owning_ranks.

Definition at line 124 of file mpi_compute_index_owner_internal.h.

◆ dofs_per_process

types::global_dof_index Utilities::MPI::internal::ComputeIndexOwner::Dictionary::dofs_per_process

The number of unknowns in the dictionary for on each MPI rank used for the index space splitting. For simplicity of index lookup without additional communication, this number is the same on all MPI ranks.

Definition at line 132 of file mpi_compute_index_owner_internal.h.

◆ local_range

std::pair<types::global_dof_index, types::global_dof_index> Utilities::MPI::internal::ComputeIndexOwner::Dictionary::local_range

The local range of the global index space that is represented in the dictionary, computed from dofs_per_process, the current MPI rank, and range_minimum_grain_size.

Definition at line 140 of file mpi_compute_index_owner_internal.h.

◆ locally_owned_size

types::global_dof_index Utilities::MPI::internal::ComputeIndexOwner::Dictionary::locally_owned_size

The actual size, computed as the minimum of dofs_per_process and the possible end of the index space. Equivalent to local_range.second - local_range.first.

Definition at line 147 of file mpi_compute_index_owner_internal.h.

◆ size

types::global_dof_index Utilities::MPI::internal::ComputeIndexOwner::Dictionary::size

The global size of the index space.

Definition at line 152 of file mpi_compute_index_owner_internal.h.

◆ n_dict_procs_in_owned_indices

unsigned int Utilities::MPI::internal::ComputeIndexOwner::Dictionary::n_dict_procs_in_owned_indices

The number of ranks the owned_indices IndexSet is distributed among.

Definition at line 158 of file mpi_compute_index_owner_internal.h.

◆ stride_small_size

unsigned int Utilities::MPI::internal::ComputeIndexOwner::Dictionary::stride_small_size

A stride to distribute the work more evenly over MPI ranks in case the grain size forces us to have fewer ranges than we have processors.

Definition at line 165 of file mpi_compute_index_owner_internal.h.


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