Reference documentation for deal.II version GIT relicensing-487-ge9eb5ab491 2024-04-25 07:20: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 | Protected Attributes | List of all members
PETScWrappers::Partitioner Class Reference

#include <deal.II/lac/petsc_communication_pattern.h>

Inheritance diagram for PETScWrappers::Partitioner:
Inheritance graph
[legend]

Public Member Functions

 Partitioner ()
 
virtual ~Partitioner () override=default
 
virtual void reinit (const IndexSet &locally_owned_indices, const IndexSet &ghost_indices, const MPI_Comm communicator) override
 
void reinit (const IndexSet &locally_owned_indices, const IndexSet &ghost_indices, const IndexSet &larger_ghost_indices, const MPI_Comm communicator)
 
unsigned int n_ghost_indices () const
 
const IndexSetghost_indices () const
 
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_start (const ArrayView< const Number > &locally_owned_array, const ArrayView< Number > &ghost_array) const
 
template<typename Number >
void export_to_ghosted_array_finish (const ArrayView< const Number > &locally_owned_array, const ArrayView< Number > &ghost_array) const
 
template<typename Number >
void import_from_ghosted_array (const VectorOperation::values op, const ArrayView< const Number > &ghost_array, const ArrayView< Number > &locally_owned_array) const
 
template<typename Number >
void import_from_ghosted_array_start (const VectorOperation::values op, const ArrayView< const Number > &ghost_array, const ArrayView< Number > &locally_owned_array) const
 
template<typename Number >
void import_from_ghosted_array_finish (const VectorOperation::values op, const ArrayView< const Number > &ghost_array, const ArrayView< Number > &locally_owned_array) const
 
MPI_Comm get_mpi_communicator () const override
 

Protected Attributes

CommunicationPattern ghost
 
CommunicationPattern larger_ghost
 
IndexSet ghost_indices_data
 
types::global_dof_index n_ghost_indices_data
 
types::global_dof_index n_ghost_indices_larger
 

Detailed Description

Partitioner implementation based on the PetscSF object. This class implements the same communication patterns of Utilities::MPI::Partitioner, internally using PetscSF API calls. Differently from the Utilities::MPI::Partitioner implementation, here we don't need to specify the communication channel, the temporary storage, and the MPI requests within export and import functions. Moreover, the import API does not zero the input ghost array.

For additional information, see the paper [200]

Definition at line 219 of file petsc_communication_pattern.h.

Constructor & Destructor Documentation

◆ Partitioner()

PETScWrappers::Partitioner::Partitioner ( )

Default constructor.

◆ ~Partitioner()

virtual PETScWrappers::Partitioner::~Partitioner ( )
overridevirtualdefault

Destructor.

Member Function Documentation

◆ reinit() [1/2]

virtual void PETScWrappers::Partitioner::reinit ( const IndexSet locally_owned_indices,
const IndexSet ghost_indices,
const MPI_Comm  communicator 
)
overridevirtual

Reinitialize the partitioner. As for the Utilities::MPI::Partitioner, any entry of ghost_indices that is also present in locally_owned_indices is discarded.

Implements Utilities::MPI::CommunicationPatternBase.

◆ reinit() [2/2]

void PETScWrappers::Partitioner::reinit ( const IndexSet locally_owned_indices,
const IndexSet ghost_indices,
const IndexSet larger_ghost_indices,
const MPI_Comm  communicator 
)

Reinitialize the partitioner. As for the Utilities::MPI::Partitioner, any entry of ghost_indices that is also present in locally_owned_indices is discarded. This reinitialization will allow to perform communications either using a ghost data array of the size of ghost_indices or of larger_ghost_indices.

◆ n_ghost_indices()

unsigned int PETScWrappers::Partitioner::n_ghost_indices ( ) const
inline

Return the actual number of ghost indices.

Definition at line 259 of file petsc_communication_pattern.h.

◆ ghost_indices()

const IndexSet & PETScWrappers::Partitioner::ghost_indices ( ) const
inline

Return an IndexSet representation of the actual ghost indices.

Definition at line 268 of file petsc_communication_pattern.h.

◆ export_to_ghosted_array()

template<typename Number >
void PETScWrappers::Partitioner::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.

◆ export_to_ghosted_array_start()

template<typename Number >
void PETScWrappers::Partitioner::export_to_ghosted_array_start ( const ArrayView< const Number > &  locally_owned_array,
const ArrayView< Number > &  ghost_array 
) const

Start the communication round to fill the vector ghost_array according to the precomputed communication pattern with values from locally_owned_array. It can be overlapped with other communications. Differently from the Utilities::MPI::Partitioner implementation, here we don't need to specify the communication channel, the temporary storage, and the MPI requests.

◆ export_to_ghosted_array_finish()

template<typename Number >
void PETScWrappers::Partitioner::export_to_ghosted_array_finish ( const ArrayView< const Number > &  locally_owned_array,
const ArrayView< Number > &  ghost_array 
) const

Finish the communication round to fill the vector ghost_array according to the precomputed communication pattern with values from locally_owned_array. It can be overlapped with other communications. Differently from the Utilities::MPI::Partitioner implementation, here we don't need to specify the communication channel, the temporary storage, and the MPI requests.

◆ import_from_ghosted_array()

template<typename Number >
void PETScWrappers::Partitioner::import_from_ghosted_array ( const VectorOperation::values  op,
const ArrayView< const Number > &  ghost_array,
const ArrayView< Number > &  locally_owned_array 
) const

Modify the vector locally_owned_array according to the precomputed communication pattern and the operation op with values from ghost_array.

◆ import_from_ghosted_array_start()

template<typename Number >
void PETScWrappers::Partitioner::import_from_ghosted_array_start ( const VectorOperation::values  op,
const ArrayView< const Number > &  ghost_array,
const ArrayView< Number > &  locally_owned_array 
) const

Start the communication round to modify the vector locally_owned_array according to the precomputed communication pattern and the operation op with values from ghost_array. It can be overlapped with other communications. Differently from the Utilities::MPI::Partitioner implementation, here we don't need to specify the communication channel, the temporary storage, and the MPI requests.

◆ import_from_ghosted_array_finish()

template<typename Number >
void PETScWrappers::Partitioner::import_from_ghosted_array_finish ( const VectorOperation::values  op,
const ArrayView< const Number > &  ghost_array,
const ArrayView< Number > &  locally_owned_array 
) const

Finish the communication round to modify the vector locally_owned_array according to the precomputed communication pattern and the operation op with values from ghost_array. It can be overlapped with other communications. Differently from the Utilities::MPI::Partitioner implementation, here we don't need to specify the communication channel, the temporary storage, and the MPI requests.

◆ get_mpi_communicator()

MPI_Comm PETScWrappers::Partitioner::get_mpi_communicator ( ) const
overridevirtual

Return the underlying MPI communicator.

Implements Utilities::MPI::CommunicationPatternBase.

Member Data Documentation

◆ ghost

CommunicationPattern PETScWrappers::Partitioner::ghost
protected

Definition at line 360 of file petsc_communication_pattern.h.

◆ larger_ghost

CommunicationPattern PETScWrappers::Partitioner::larger_ghost
protected

Definition at line 360 of file petsc_communication_pattern.h.

◆ ghost_indices_data

IndexSet PETScWrappers::Partitioner::ghost_indices_data
protected

Definition at line 361 of file petsc_communication_pattern.h.

◆ n_ghost_indices_data

types::global_dof_index PETScWrappers::Partitioner::n_ghost_indices_data
protected

Definition at line 362 of file petsc_communication_pattern.h.

◆ n_ghost_indices_larger

types::global_dof_index PETScWrappers::Partitioner::n_ghost_indices_larger
protected

Definition at line 363 of file petsc_communication_pattern.h.


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