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 Member Functions | Protected Attributes | List of all members
PETScWrappers::CommunicationPattern Class Reference

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

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

Public Member Functions

 CommunicationPattern ()
 
virtual ~CommunicationPattern () override
 
virtual 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 > &indices_locally_owned, const std::vector< types::global_dof_index > &indices_want, const MPI_Comm communicator)
 
void reinit (const types::global_dof_index local_size, const IndexSet &ghost_indices, 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_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
 
 operator const PetscSF & () const
 
void clear ()
 

Protected Member Functions

void do_reinit (const std::vector< PetscInt > &inidx, const std::vector< PetscInt > &inloc, const std::vector< PetscInt > &outidx, const std::vector< PetscInt > &outloc, const MPI_Comm communicator)
 

Protected Attributes

PetscSF sf
 

Detailed Description

CommunicationPattern implementation based on the PetscSF object. This class implements the same communication patterns of Utilities::MPI::NoncontiguousPartitioner, internally using PetscSF API calls.

For additional information, see the paper [200]

Definition at line 48 of file petsc_communication_pattern.h.

Constructor & Destructor Documentation

◆ CommunicationPattern()

PETScWrappers::CommunicationPattern::CommunicationPattern ( )

Default constructor.

Definition at line 35 of file petsc_communication_pattern.cc.

◆ ~CommunicationPattern()

PETScWrappers::CommunicationPattern::~CommunicationPattern ( )
overridevirtual

Destructor.

Definition at line 41 of file petsc_communication_pattern.cc.

Member Function Documentation

◆ reinit() [1/3]

void PETScWrappers::CommunicationPattern::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 87 of file petsc_communication_pattern.cc.

◆ reinit() [2/3]

void PETScWrappers::CommunicationPattern::reinit ( const std::vector< types::global_dof_index > &  indices_locally_owned,
const std::vector< types::global_dof_index > &  indices_want,
const MPI_Comm  communicator 
)

Reinitialize the communication pattern. The argument indices_locally_owned and indices_want indicates the owned and required dofs, respectively. This allows the indices to not be sorted and 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.

The export_to_ghost_array will populate an array with the values associated to the indices_want only.

This emulates the corresponding constructor in Utilities::MPI::NoncontiguousPartitioner.

Definition at line 107 of file petsc_communication_pattern.cc.

◆ reinit() [3/3]

void PETScWrappers::CommunicationPattern::reinit ( const types::global_dof_index  local_size,
const IndexSet ghost_indices,
const MPI_Comm  communicator 
)

Reinitialization that takes the number of locally-owned degrees of freedom local_size and an index set for the required ghost indices ghost_indices.

The local index range is translated to global indices in an ascending and one-to-one fashion, i.e., the indices of process \(p\) sit exactly between the indices of the processes \(p-1\) and \(p+1\), respectively.

The export_to_ghost_array will populate an array containing values from locally-owned AND ghost indices, as for the relevant set of dofs of a usual FEM simulation.

Definition at line 49 of file petsc_communication_pattern.cc.

◆ export_to_ghosted_array()

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

Definition at line 315 of file petsc_communication_pattern.cc.

◆ export_to_ghosted_array_start()

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

Definition at line 279 of file petsc_communication_pattern.cc.

◆ export_to_ghosted_array_finish()

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

Definition at line 297 of file petsc_communication_pattern.cc.

◆ import_from_ghosted_array()

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

Definition at line 358 of file petsc_communication_pattern.cc.

◆ import_from_ghosted_array_start()

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

Definition at line 327 of file petsc_communication_pattern.cc.

◆ import_from_ghosted_array_finish()

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

Definition at line 343 of file petsc_communication_pattern.cc.

◆ get_mpi_communicator()

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

Return the underlying MPI communicator.

Implements Utilities::MPI::CommunicationPatternBase.

Definition at line 270 of file petsc_communication_pattern.cc.

◆ operator const PetscSF &()

PETScWrappers::CommunicationPattern::operator const PetscSF & ( ) const
inline

Conversion operator to gain access to the underlying PETSc object.

Definition at line 179 of file petsc_communication_pattern.h.

◆ clear()

void PETScWrappers::CommunicationPattern::clear ( )

Reset the object.

Definition at line 262 of file petsc_communication_pattern.cc.

◆ do_reinit()

void PETScWrappers::CommunicationPattern::do_reinit ( const std::vector< PetscInt > &  inidx,
const std::vector< PetscInt > &  inloc,
const std::vector< PetscInt > &  outidx,
const std::vector< PetscInt > &  outloc,
const MPI_Comm  communicator 
)
protected

General setup

Definition at line 162 of file petsc_communication_pattern.cc.

Member Data Documentation

◆ sf

PetscSF PETScWrappers::CommunicationPattern::sf
protected

A generic PetscSF object that will perform the communication.

Definition at line 194 of file petsc_communication_pattern.h.


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