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 | Public Attributes | Private Member Functions | List of all members
Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload Class Reference

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

Inheritance diagram for Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload:
[legend]

Public Member Functions

 ConsensusAlgorithmsPayload (const IndexSet &owned_indices, const IndexSet &indices_to_look_up, const MPI_Comm comm, std::vector< unsigned int > &owning_ranks, const bool track_index_requests=false)
 
virtual void answer_request (const unsigned int other_rank, const std::vector< std::pair< types::global_dof_index, types::global_dof_index > > &buffer_recv, std::vector< unsigned int > &request_buffer) override
 
virtual std::vector< unsigned intcompute_targets () override
 
virtual void create_request (const unsigned int other_rank, std::vector< std::pair< types::global_dof_index, types::global_dof_index > > &send_buffer) override
 
virtual void read_answer (const unsigned int other_rank, const std::vector< unsigned int > &recv_buffer) override
 
std::map< unsigned int, IndexSetget_requesters ()
 

Public Attributes

const IndexSetowned_indices
 
const IndexSetindices_to_look_up
 
const MPI_Comm comm
 
const unsigned int my_rank
 
const unsigned int n_procs
 
const bool track_index_requests
 
std::vector< unsigned int > & owning_ranks
 
std::vector< std::vector< std::pair< unsigned int, std::vector< std::pair< unsigned int, unsigned int > > > > > requesters
 
Dictionary dict
 
std::map< unsigned int, std::pair< std::vector< types::global_dof_index >, std::vector< unsigned int > > > indices_to_look_up_by_dict_rank
 

Private Member Functions

void append_index_origin (const unsigned int index_within_dictionary, const unsigned int rank_of_request, const unsigned int rank_of_owner, unsigned int &owner_index_guess)
 

Detailed Description

Specialization of ConsensusAlgorithms::Process for the context of Utilities::MPI::compute_index_owner() and Utilities::MPI::Partitioner::set_ghost_indices() with additional payload.

Definition at line 217 of file mpi_compute_index_owner_internal.h.

Constructor & Destructor Documentation

◆ ConsensusAlgorithmsPayload()

Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload::ConsensusAlgorithmsPayload ( const IndexSet owned_indices,
const IndexSet indices_to_look_up,
const MPI_Comm  comm,
std::vector< unsigned int > &  owning_ranks,
const bool  track_index_requests = false 
)

Constructor.

Definition at line 462 of file mpi_compute_index_owner_internal.cc.

Member Function Documentation

◆ answer_request()

void Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload::answer_request ( const unsigned int  other_rank,
const std::vector< std::pair< types::global_dof_index, types::global_dof_index > > &  buffer_recv,
std::vector< unsigned int > &  request_buffer 
)
overridevirtual

Implementation of Utilities::MPI::ConsensusAlgorithms::Process::answer_request(), adding the owner of a particular index in request_buffer (and keeping track of who requested a particular index in case that information is also desired).

Reimplemented from Utilities::MPI::ConsensusAlgorithms::Process< std::vector< std::pair< types::global_dof_index, types::global_dof_index > >, std::vector< unsigned int > >.

Definition at line 483 of file mpi_compute_index_owner_internal.cc.

◆ compute_targets()

std::vector< unsigned int > Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload::compute_targets ( )
overridevirtual

◆ create_request()

void Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload::create_request ( const unsigned int  other_rank,
std::vector< std::pair< types::global_dof_index, types::global_dof_index > > &  send_buffer 
)
overridevirtual

◆ read_answer()

void Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload::read_answer ( const unsigned int  other_rank,
const std::vector< unsigned int > &  recv_buffer 
)
overridevirtual

◆ get_requesters()

std::map< unsigned int, IndexSet > Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload::get_requesters ( )

Resolve the origin of the requests by sending the information accumulated in terms of the dictionary owners during the run of the consensus algorithm back to the owner in the original IndexSet. This requires some point-to-point communication.

Returns
Map of processors and associated ranges of indices that are requested from the current rank

Definition at line 582 of file mpi_compute_index_owner_internal.cc.

◆ append_index_origin()

void Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload::append_index_origin ( const unsigned int  index_within_dictionary,
const unsigned int  rank_of_request,
const unsigned int  rank_of_owner,
unsigned int owner_index_guess 
)
private

Stores the index request in the requesters field. Given the rank of the owner, we start with a guess for the index at the owner's site. This is because we typically might look up on the same rank several times in a row, hence avoiding the binary search in Dictionary::get_owning_rank_index()). Once we know the index at the owner, we fill the vector entry with the rank of the request. Here, we utilize the fact that requests are processed rank-by-rank, so we can simply look at the end of the vector whether there is already some data stored or not. Finally, we build ranges, again using that the index list is sorted and we therefore only need to append at the end.

Definition at line 739 of file mpi_compute_index_owner_internal.cc.

Member Data Documentation

◆ owned_indices

const IndexSet& Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload::owned_indices

The index space which describes the locally owned space.

Definition at line 236 of file mpi_compute_index_owner_internal.h.

◆ indices_to_look_up

const IndexSet& Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload::indices_to_look_up

The indices which are "ghosts" on a given rank and should be looked up in terms of their owner rank from owned_indices.

Definition at line 242 of file mpi_compute_index_owner_internal.h.

◆ comm

const MPI_Comm Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload::comm

The underlying MPI communicator.

Definition at line 247 of file mpi_compute_index_owner_internal.h.

◆ my_rank

const unsigned int Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload::my_rank

The present MPI rank.

Definition at line 252 of file mpi_compute_index_owner_internal.h.

◆ n_procs

const unsigned int Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload::n_procs

The total number of ranks participating in the MPI communicator comm.

Definition at line 258 of file mpi_compute_index_owner_internal.h.

◆ track_index_requests

const bool Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload::track_index_requests

Controls whether the origin of ghost owner should also be stored. If true, it will be added into requesters and can be queried by get_requesters().

Definition at line 265 of file mpi_compute_index_owner_internal.h.

◆ owning_ranks

std::vector<unsigned int>& Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload::owning_ranks

The result of the index owner computation: To each index contained in indices_to_look_up, this vector contains the MPI rank of the owner in owned_indices.

Definition at line 272 of file mpi_compute_index_owner_internal.h.

◆ requesters

std::vector<std::vector< std::pair<unsigned int, std::vector<std::pair<unsigned int, unsigned int> > > > > Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload::requesters

Keeps track of the origin of the requests. The layout of the data structure is as follows: The outermost vector has as many entries as Dictionary::actually_owning_rank_list and represents the information we should send back to the owners from the present dictionary entry. The second vector then collects a list of MPI ranks that have requested data, using the rank in the first pair entry and a list of index ranges as the second entry.

Definition at line 286 of file mpi_compute_index_owner_internal.h.

◆ dict

Dictionary Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload::dict

The dictionary handling the requests.

Definition at line 291 of file mpi_compute_index_owner_internal.h.

◆ indices_to_look_up_by_dict_rank

std::map<unsigned int, std::pair<std::vector<types::global_dof_index>, std::vector<unsigned int> > > Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload::indices_to_look_up_by_dict_rank

Array to collect the indices to look up (first vector) and their local index among indices (second vector), sorted by the rank in the dictionary.

Definition at line 301 of file mpi_compute_index_owner_internal.h.


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