Discontinuous Galerkin Library
#include "dg/algorithm.h"
dg::SurjectiveComm< Index, Vector > Struct Template Reference

Perform surjective gather and its transpose (scatter) operation across processes on distributed vectors using mpi. More...

Inheritance diagram for dg::SurjectiveComm< Index, Vector >:
[legend]

Public Member Functions

 SurjectiveComm ()
 no memory allocation; size 0 More...
 
 SurjectiveComm (unsigned local_size, const thrust::host_vector< int > &localIndexMap, const thrust::host_vector< int > &pidIndexMap, MPI_Comm comm)
 Construct from local indices and PIDs index map. More...
 
template<class ConversionPolicy >
 SurjectiveComm (const thrust::host_vector< int > &globalIndexMap, const ConversionPolicy &p)
 Construct from global indices index map. More...
 
template<class OtherIndex , class OtherVector >
 SurjectiveComm (const SurjectiveComm< OtherIndex, OtherVector > &src)
 reconstruct from another type; if src is empty, same as default constructor More...
 
const thrust::host_vector< int > & getLocalIndexMap () const
 read access to the local index index map given in constructor More...
 
const thrust::host_vector< int > & getPidIndexMap () const
 read access to the pid index map given in constructor More...
 
const Index & getSortedIndexMap () const
 
virtual SurjectiveCommclone () const override final
 Generic copy method. More...
 
bool isLocalBijective () const
 No reduction on this process? True: no reduction, False: need to reduce. More...
 
- Public Member Functions inherited from dg::aCommunicator< Vector >
Vector allocate_buffer () const
 Allocate a buffer object of size buffer_size() More...
 
void global_gather (const value_type *values, Vector &buffer) const
 \( w = G v\). Globally (across processes) gather data into a buffer More...
 
Vector global_gather (const value_type *values) const
 \( w = G v\). Globally (across processes) gather data into a buffer (memory allocating version) More...
 
void global_scatter_reduce (const Vector &toScatter, value_type *values) const
 \( v = G^\mathrm{T} w\). Globally (across processes) scatter data accross processes and reduce on multiple indices More...
 
unsigned buffer_size () const
 The local size of the buffer vector w = local map size. More...
 
unsigned local_size () const
 The local size of the source vector v = local size of the dg::MPI_Vector More...
 
bool isCommunicating () const
 True if the gather/scatter operation involves actual MPI communication. More...
 
MPI_Comm communicator () const
 The internal MPI communicator used. More...
 
virtual aCommunicatorclone () const=0
 Generic copy method. More...
 
virtual ~aCommunicator ()
 vritual destructor More...
 

Additional Inherited Members

- Public Types inherited from dg::aCommunicator< Vector >
using value_type = get_value_type< Vector >
 reveal value type More...
 
using container_type = Vector
 reveal local container type More...
 
- Protected Member Functions inherited from dg::aCommunicator< Vector >
 aCommunicator (unsigned local_size=0)
 only derived classes can construct More...
 
 aCommunicator (const aCommunicator &src)
 only derived classes can copy More...
 
aCommunicatoroperator= (const aCommunicator &src)
 only derived classes can assign More...
 
void set_local_size (unsigned new_size)
 Set the local size of the source vector v. More...
 

Detailed Description

template<class Index, class Vector>
struct dg::SurjectiveComm< Index, Vector >

Perform surjective gather and its transpose (scatter) operation across processes on distributed vectors using mpi.

This Communicator performs surjective global gather and scatter operations, which means that the index map is surjective: If the index map idx[i] is surjective, each element of the source vector v maps to at least one location in the buffer vector w. This means that the scatter matrix S can have more than one 1's in each line. (see aCommunicator for more details) Compared to BijectiveComm in the global_gather function there is an additional gather and in the global_scatter_reduce function a reduction needs to be performed.

Template Parameters
Indexan integer thrust Vector (needs to be int due to MPI interface)
Vectora thrust Vector

Constructor & Destructor Documentation

◆ SurjectiveComm() [1/4]

template<class Index , class Vector >
dg::SurjectiveComm< Index, Vector >::SurjectiveComm ( )
inline

no memory allocation; size 0

◆ SurjectiveComm() [2/4]

template<class Index , class Vector >
dg::SurjectiveComm< Index, Vector >::SurjectiveComm ( unsigned  local_size,
const thrust::host_vector< int > &  localIndexMap,
const thrust::host_vector< int > &  pidIndexMap,
MPI_Comm  comm 
)
inline

Construct from local indices and PIDs index map.

The indices in the index map are written with respect to the buffer vector. Each location in the source vector is uniquely specified by a local vector index and the process rank.

Parameters
local_sizelocal size of a dg::MPI_Vector (same for all processes)
localIndexMapEach element localIndexMap[i] represents a local vector index from (or to) where to take the value buffer[i]. There are local_buffer_size = localIndexMap.size() elements.
pidIndexMapEach element pidIndexMap[i] represents the pid/rank to which the corresponding index localIndexMap[i] is local. Same size as localIndexMap. The pid/rank needs to be element of the given communicator.
commThe MPI communicator participating in the scatter/gather operations
Note
we assume that the index map is surjective

◆ SurjectiveComm() [3/4]

template<class Index , class Vector >
template<class ConversionPolicy >
dg::SurjectiveComm< Index, Vector >::SurjectiveComm ( const thrust::host_vector< int > &  globalIndexMap,
const ConversionPolicy &  p 
)
inline

Construct from global indices index map.

Uses the global2localIdx() member of MPITopology to generate localIndexMap and pidIndexMap

Parameters
globalIndexMapEach element globalIndexMap[i] represents a global vector index from (or to) where to take the value buffer[i]. There are local_buffer_size = globalIndexMap.size() elements.
pthe conversion object
Template Parameters
ConversionPolicyhas to have the members:
  • bool global2localIdx(unsigned,unsigned&,unsigned&) const; where the first parameter is the global index and the other two are the output pair (localIdx, rank). return true if successful, false if global index is not part of the grid
  • MPI_Comm communicator() const; returns the communicator to use in the gather/scatter
  • local_size(); return the local vector size
See also
Topology base classes the MPI grids defined in Level 3 can all be used as a ConversionPolicy
Note
we assume that the index map is surjective

◆ SurjectiveComm() [4/4]

template<class Index , class Vector >
template<class OtherIndex , class OtherVector >
dg::SurjectiveComm< Index, Vector >::SurjectiveComm ( const SurjectiveComm< OtherIndex, OtherVector > &  src)
inline

reconstruct from another type; if src is empty, same as default constructor

Member Function Documentation

◆ clone()

template<class Index , class Vector >
virtual SurjectiveComm * dg::SurjectiveComm< Index, Vector >::clone ( ) const
inlinefinaloverridevirtual

Generic copy method.

Returns
pointer to allocated object

Implements dg::aCommunicator< Vector >.

◆ getLocalIndexMap()

template<class Index , class Vector >
const thrust::host_vector< int > & dg::SurjectiveComm< Index, Vector >::getLocalIndexMap ( ) const
inline

read access to the local index index map given in constructor

◆ getPidIndexMap()

template<class Index , class Vector >
const thrust::host_vector< int > & dg::SurjectiveComm< Index, Vector >::getPidIndexMap ( ) const
inline

read access to the pid index map given in constructor

◆ getSortedIndexMap()

template<class Index , class Vector >
const Index & dg::SurjectiveComm< Index, Vector >::getSortedIndexMap ( ) const
inline

◆ isLocalBijective()

template<class Index , class Vector >
bool dg::SurjectiveComm< Index, Vector >::isLocalBijective ( ) const
inline

No reduction on this process? True: no reduction, False: need to reduce.


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