DOLFIN-X
DOLFIN-X C++ interface
Public Types | Public Member Functions | List of all members
dolfinx::common::IndexMap Class Reference

This class represents the distribution index arrays across processes. An index array is a contiguous collection of N+1 block indices [0, 1, . . ., N] that are distributed across processes M processes. On a given process, the IndexMap stores a portion of the index set using local indices [0, 1, . . . , n], and a map from the local block indices to a unique global block index. More...

#include <IndexMap.h>

Public Types

enum  Mode { insert, add }
 Mode for reverse scatter operation.
 

Public Member Functions

 IndexMap (MPI_Comm mpi_comm, std::int32_t local_size, const std::vector< std::int64_t > &ghosts, int block_size)
 Create Index map with local_size owned blocks on this process, and blocks have size block_size. More...
 
 IndexMap (MPI_Comm mpi_comm, std::int32_t local_size, const Eigen::Ref< const Eigen::Array< std::int64_t, Eigen::Dynamic, 1 >> &ghosts, int block_size)
 Create Index map with local_size owned blocks on this process, and blocks have size block_size. More...
 
 IndexMap (const IndexMap &map)=delete
 Copy constructor.
 
 IndexMap (IndexMap &&map)=default
 Move constructor.
 
 ~IndexMap ()=default
 Destructor.
 
std::array< std::int64_t, 2 > local_range () const
 Range of indices (global) owned by this process.
 
int block_size () const
 Block size.
 
std::int32_t num_ghosts () const
 Number of ghost indices on this process.
 
std::int32_t size_local () const
 Number of indices owned by on this process.
 
std::int64_t size_global () const
 Number indices across communicator.
 
const Eigen::Array< std::int64_t, Eigen::Dynamic, 1 > & ghosts () const
 Local-to-global map for ghosts (local indexing beyond end of local range)
 
Eigen::Array< std::int64_t, Eigen::Dynamic, 1 > local_to_global (const Eigen::Ref< const Eigen::Array< std::int32_t, Eigen::Dynamic, 1 >> &indices, bool blocked=true) const
 Compute global indices for array of local indices. More...
 
std::vector< std::int64_t > local_to_global (const std::vector< std::int32_t > &indices, bool blocked=true) const
 
std::vector< std::int32_t > global_to_local (const std::vector< std::int64_t > &indices, bool blocked=true) const
 Compute local indices for array of global indices. More...
 
std::vector< std::int32_t > global_to_local (const Eigen::Ref< const Eigen::Array< std::int64_t, Eigen::Dynamic, 1 >> &indices, bool blocked=true) const
 Compute local indices for array of global indices. More...
 
std::vector< std::int64_t > global_indices (bool blocked=true) const
 Global indices. More...
 
std::int64_t local_to_global (std::int32_t local_index) const
 
const std::vector< std::int32_t > & forward_indices () const
 
Eigen::Array< std::int32_t, Eigen::Dynamic, 1 > ghost_owners () const
 Owner rank (on global communicator) of each ghost entry.
 
int owner (std::int64_t global_index) const
 Get process that owns index (global block index)
 
Eigen::Array< std::int64_t, Eigen::Dynamic, 1 > indices (bool unroll_block) const
 Return array of global indices for all indices on this process, including ghosts.
 
MPI_Comm mpi_comm () const
 Return MPI communicator. More...
 
const std::vector< std::int32_t > & neighbours () const
 Neighbors for neigborhood communicator.
 
std::map< std::int32_t, std::set< std::int32_t > > compute_shared_indices () const
 
void scatter_fwd (const std::vector< std::int64_t > &local_data, std::vector< std::int64_t > &remote_data, int n) const
 Send n values for each index that is owned to processes that have the index as a ghost. The size of the input array local_data must be the same as n * size_local(). More...
 
void scatter_fwd (const std::vector< std::int32_t > &local_data, std::vector< std::int32_t > &remote_data, int n) const
 Send n values for each index that is owned to processes that have the index as a ghost. The size of the input array local_data must be the same as n * size_local(). More...
 
std::vector< std::int64_t > scatter_fwd (const std::vector< std::int64_t > &local_data, int n) const
 Send n values for each index that is owned to processes that have the index as a ghost. The size of the input array local_data must be the same as n * size_local(). More...
 
std::vector< std::int32_t > scatter_fwd (const std::vector< std::int32_t > &local_data, int n) const
 Send n values for each index that is owned to processes that have the index as a ghost. More...
 
void scatter_rev (std::vector< std::int64_t > &local_data, const std::vector< std::int64_t > &remote_data, int n, IndexMap::Mode op) const
 Send n values for each ghost index to owning to the process. More...
 
void scatter_rev (std::vector< std::int32_t > &local_data, const std::vector< std::int32_t > &remote_data, int n, IndexMap::Mode op) const
 Send n values for each ghost index to owning to the process. More...
 

Detailed Description

This class represents the distribution index arrays across processes. An index array is a contiguous collection of N+1 block indices [0, 1, . . ., N] that are distributed across processes M processes. On a given process, the IndexMap stores a portion of the index set using local indices [0, 1, . . . , n], and a map from the local block indices to a unique global block index.

Constructor & Destructor Documentation

◆ IndexMap() [1/2]

IndexMap::IndexMap ( MPI_Comm  mpi_comm,
std::int32_t  local_size,
const std::vector< std::int64_t > &  ghosts,
int  block_size 
)

Create Index map with local_size owned blocks on this process, and blocks have size block_size.

Collective

Parameters
[in]mpi_commThe MPI communicator
[in]local_sizeLocal size of the IndexMap, i.e. the number of owned entries
[in]ghostsThe global indices of ghost entries
[in]block_sizeThe block size of the IndexMap

◆ IndexMap() [2/2]

IndexMap::IndexMap ( MPI_Comm  mpi_comm,
std::int32_t  local_size,
const Eigen::Ref< const Eigen::Array< std::int64_t, Eigen::Dynamic, 1 >> &  ghosts,
int  block_size 
)

Create Index map with local_size owned blocks on this process, and blocks have size block_size.

Collective

Parameters
[in]mpi_commThe MPI communicator
[in]local_sizeLocal size of the IndexMap, i.e. the number of owned entries
[in]ghostsThe global indices of ghost entries
[in]block_sizeThe block size of the IndexMap

Member Function Documentation

◆ compute_shared_indices()

std::map< int, std::set< int > > IndexMap::compute_shared_indices ( ) const
Todo:
Aim to remove this function

Compute map from each local index to the complete set of sharing processes for that index

Returns
shared indices

◆ forward_indices()

const std::vector<std::int32_t>& dolfinx::common::IndexMap::forward_indices ( ) const
inline
Todo:
Reconsider name Local (owned) indices shared with neighbour processes, i.e. are ghosts on other processes
Returns
List of indices that are ghosted on other processes

◆ global_indices()

std::vector< std::int64_t > IndexMap::global_indices ( bool  blocked = true) const

Global indices.

Returns
The global index for all local indices (0, 1, 2, ...) on this process, including ghosts

◆ global_to_local() [1/2]

std::vector< std::int32_t > IndexMap::global_to_local ( const Eigen::Ref< const Eigen::Array< std::int64_t, Eigen::Dynamic, 1 >> &  indices,
bool  blocked = true 
) const

Compute local indices for array of global indices.

Parameters
[in]indicesGlobal indices
[in]blockedIf true work with blocked indices. If false the input indices are not block-wise.
Returns
The local of the corresponding global index in indices. Return -1 if the local index does not exist on this process.

◆ global_to_local() [2/2]

std::vector< std::int32_t > IndexMap::global_to_local ( const std::vector< std::int64_t > &  indices,
bool  blocked = true 
) const

Compute local indices for array of global indices.

Parameters
[in]indicesGlobal indices
[in]blockedIf true work with blocked indices. If false the input indices are not block-wise.
Returns
The local of the corresponding global index in indices. Return -1 if the local index does not exist on this process.

◆ local_to_global() [1/3]

Eigen::Array< std::int64_t, Eigen::Dynamic, 1 > IndexMap::local_to_global ( const Eigen::Ref< const Eigen::Array< std::int32_t, Eigen::Dynamic, 1 >> &  indices,
bool  blocked = true 
) const

Compute global indices for array of local indices.

Parameters
[in]indicesLocal indices
[in]blockedIf true work with blocked indices. If false the input indices are not block-wise.
Returns
The global index of the corresponding local index in indices.

◆ local_to_global() [2/3]

std::vector< std::int64_t > IndexMap::local_to_global ( const std::vector< std::int32_t > &  indices,
bool  blocked = true 
) const
Todo:
Consider removing this function in favour of the version that accepts an Eigen array.

Compute global indices for array of local indices

Parameters
[in]indicesLocal indices
[in]blockedIf true work with blocked indices. If false the input indices are not block-wise.
Returns
The global index of the corresponding local index in indices.

◆ local_to_global() [3/3]

std::int64_t dolfinx::common::IndexMap::local_to_global ( std::int32_t  local_index) const
inline
Todo:
Remove this function Get global index for local index i (index of the block)

◆ mpi_comm()

MPI_Comm IndexMap::mpi_comm ( ) const

Return MPI communicator.

Returns
The communicator on which the IndexMap is defined

◆ scatter_fwd() [1/4]

std::vector< std::int32_t > IndexMap::scatter_fwd ( const std::vector< std::int32_t > &  local_data,
int  n 
) const

Send n values for each index that is owned to processes that have the index as a ghost.

Parameters
[in]local_dataLocal data associated with each owned local index to be sent to process where the data is ghosted. Size must be n * size_local().
[in]nNumber of data items per index
Returns
Ghost data on this process received from the owning process. Size will be n * num_ghosts().

◆ scatter_fwd() [2/4]

void IndexMap::scatter_fwd ( const std::vector< std::int32_t > &  local_data,
std::vector< std::int32_t > &  remote_data,
int  n 
) const

Send n values for each index that is owned to processes that have the index as a ghost. The size of the input array local_data must be the same as n * size_local().

Parameters
[in]local_dataLocal data associated with each owned local index to be sent to process where the data is ghosted. Size must be n * size_local().
[in,out]remote_dataGhost data on this process received from the owning process. Size will be n * num_ghosts().
[in]nNumber of data items per index

◆ scatter_fwd() [3/4]

std::vector< std::int64_t > IndexMap::scatter_fwd ( const std::vector< std::int64_t > &  local_data,
int  n 
) const

Send n values for each index that is owned to processes that have the index as a ghost. The size of the input array local_data must be the same as n * size_local().

Parameters
[in]local_dataLocal data associated with each owned local index to be sent to process where the data is ghosted. Size must be n * size_local().
[in]nNumber of data items per index
Returns
Ghost data on this process received from the owning process. Size will be n * num_ghosts().

◆ scatter_fwd() [4/4]

void IndexMap::scatter_fwd ( const std::vector< std::int64_t > &  local_data,
std::vector< std::int64_t > &  remote_data,
int  n 
) const

Send n values for each index that is owned to processes that have the index as a ghost. The size of the input array local_data must be the same as n * size_local().

Parameters
[in]local_dataLocal data associated with each owned local index to be sent to process where the data is ghosted. Size must be n * size_local().
[in,out]remote_dataGhost data on this process received from the owning process. Size will be n * num_ghosts().
[in]nNumber of data items per index

◆ scatter_rev() [1/2]

void IndexMap::scatter_rev ( std::vector< std::int32_t > &  local_data,
const std::vector< std::int32_t > &  remote_data,
int  n,
IndexMap::Mode  op 
) const

Send n values for each ghost index to owning to the process.

Parameters
[in,out]local_dataLocal data associated with each owned local index to be sent to process where the data is ghosted. Size must be n * size_local().
[in]remote_dataGhost data on this process received from the owning process. Size will be n * num_ghosts().
[in]nNumber of data items per index
[in]opSum or set received values in local_data

◆ scatter_rev() [2/2]

void IndexMap::scatter_rev ( std::vector< std::int64_t > &  local_data,
const std::vector< std::int64_t > &  remote_data,
int  n,
IndexMap::Mode  op 
) const

Send n values for each ghost index to owning to the process.

Parameters
[in,out]local_dataLocal data associated with each owned local index to be sent to process where the data is ghosted. Size must be n * size_local().
[in]remote_dataGhost data on this process received from the owning process. Size will be n * num_ghosts().
[in]nNumber of data items per index
[in]opSum or set received values in local_data

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