Common (dolfinx::common)

namespace common

Miscellaneous classes, functions and types.

This namespace provides utility type functions for managing subsystems, convenience classes and library-wide typedefs.

Functions

std::vector<int32_t> compute_owned_indices(const std::span<const std::int32_t> &indices, const IndexMap &map)

Given a vector of indices (local numbering, owned or ghost) and an index map, this function returns the indices owned by this process, including indices that might have been in the list of indices on another processes.

Parameters
  • indices[in] List of indices

  • map[in] The index map

Returns

Indices owned by the calling process

std::tuple<std::int64_t, std::vector<std::int32_t>, std::vector<std::vector<std::int64_t>>, std::vector<std::vector<int>>> stack_index_maps(const std::vector<std::pair<std::reference_wrapper<const common::IndexMap>, int>> &maps)

Compute layout data and ghost indices for a stacked (concatenated) index map, i.e. ‘splice’ multiple maps into one.

The input maps are concatenated, with indices in maps and owned by the caller remaining owned by the caller. Ghost data is stored at the end of the local range as normal, with the ghosts in blocks in the order of the index maps in maps.

Note

Index maps with a block size are unrolled in the data for the concatenated index map.

Note

Communication is required to compute the new ghost indices.

Parameters

maps[in] List of (index map, block size) pairs

Returns

The (0) global offset of a concatenated map for the calling rank, (1) local offset for the owned indices of each submap in the concatenated map, (2) new indices for the ghosts for each submap, and (3) owner rank of each ghost entry for each submap.

template<typename U, typename V>
std::pair<std::vector<typename U::value_type>, std::vector<typename V::value_type>> sort_unique(const U &indices, const V &values)

Sort two arrays based on the values in array indices. Any duplicate indices and the corresponding value are removed. In the case of duplicates, the entry with the smallest value is retained.

Parameters
  • indices[in] Array of indices

  • values[in] Array of values

Returns

Sorted (indices, values), with sorting based on indices

template<class T>
std::size_t hash_local(const T &x)

Compute a hash of a given object.

The hash is computed using Boost container hash (https://www.boost.org/doc/libs/release/libs/container_hash/).

Parameters

x[in] The object to compute a hash of.

Returns

The hash values.

template<class T>
std::size_t hash_global(MPI_Comm comm, const T &x)

Compute a hash for a distributed (MPI) object.

A hash is computed on each process for the local part of the obejct. Then, a hash of the std::vector containing each local hash key in rank order is returned.

Note

Collective

Parameters
  • comm[in] The communicator on which to compute the hash.

  • x[in] The object to compute a hash of.

Returns

The hash values.

class IndexMap
#include <IndexMap.h>

This class represents the distribution index arrays across processes. An index array is a contiguous collection of N+1 indices [0, 1, …, N] that are distributed across 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 indices to a unique global index.

Public Functions

IndexMap(MPI_Comm comm, std::int32_t local_size)

Create an non-overlapping index map.

Note

Collective

Parameters
  • comm[in] The MPI communicator

  • local_size[in] Local size of the index map, i.e. the number of owned entries

IndexMap(MPI_Comm comm, std::int32_t local_size, const std::span<const std::int64_t> &ghosts, const std::span<const int> &owners)

Create an overlapping (ghosted) index map.

This constructor uses a ‘consensus’ algorithm to determine the ranks that ghost indices that are owned by the caller. This requires non-trivial MPI communication. If the ranks that ghost indices owned by the caller are known, it more efficient to use the constructor that takes these ranks as an argument.

Note

Collective

Parameters
  • comm[in] The MPI communicator

  • local_size[in] Local size of the index map, i.e. the number of owned entries

  • ghosts[in] The global indices of ghost entries

  • owners[in] Owner rank (on global communicator) of each entry in ghosts

IndexMap(MPI_Comm comm, std::int32_t local_size, const std::array<std::vector<int>, 2> &src_dest, const std::span<const std::int64_t> &ghosts, const std::span<const int> &owners)

Create an overlapping (ghosted) index map.

This constructor is optimised for the case where the ‘source’ (ranks that own indices ghosted by the caller) and ‘destination’ ranks (ranks that ghost indices owned by the caller) are already available. It allows the complex computation of the destination ranks from owners.

Note

Collective

Parameters
  • comm[in] The MPI communicator

  • local_size[in] Local size of the index map, i.e. the number

  • src_dest[in] Lists of [0] src and [1] dest ranks. The list in each must be sorted and not contain duplicates. src ranks are owners of the indices in ghosts. dest ranks are the rank that ghost indices owned by the caller.

  • ghosts[in] The global indices of ghost entries

  • owners[in] Owner rank (on global communicator) of each entry in ghosts

IndexMap(IndexMap &&map) = default

Move constructor.

~IndexMap() = default

Destructor.

IndexMap &operator=(IndexMap &&map) = default

Move assignment.

std::array<std::int64_t, 2> local_range() const noexcept

Range of indices (global) owned by this process.

std::int32_t num_ghosts() const noexcept

Number of ghost indices on this process.

std::int32_t size_local() const noexcept

Number of indices owned by on this process.

std::int64_t size_global() const noexcept

Number indices across communicator.

const std::vector<std::int64_t> &ghosts() const noexcept

Local-to-global map for ghosts (local indexing beyond end of local range)

MPI_Comm comm() const

Return the MPI communicator that the map is defined on.

Returns

Communicator

void local_to_global(const std::span<const std::int32_t> &local, const std::span<std::int64_t> &global) const

Compute global indices for array of local indices.

Parameters
  • local[in] Local indices

  • global[out] The global indices

void global_to_local(const std::span<const std::int64_t> &global, const std::span<std::int32_t> &local) const

Compute local indices for array of global indices.

Parameters
  • global[in] Global indices

  • local[out] The local of the corresponding global index in ‘global’. Returns -1 if the local index does not exist on this process.

std::vector<std::int64_t> global_indices() const

Build list of indices with global indexing.

Returns

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

inline const std::vector<int> &owners() const

The ranks that own each ghost index.

Returns

List of ghost owners. The owning rank of the ith ghost index is owners()[i].

std::pair<IndexMap, std::vector<std::int32_t>> create_submap(const std::span<const std::int32_t> &indices) const

Create new index map from a subset of indices in this index map.

The order of the owned indices is preserved, with new map effectively a ‘compressed’ map.

Parameters

indices[in] Local indices in the map that should appear in the new index map. All indices must be owned, i.e. indices must be less than this->size_local().

Pre

indices must be sorted and contain no duplicates.

Returns

The (i) new index map and (ii) a map from the ghost position in the new map to the ghost position in the original (this) map

graph::AdjacencyList<int> index_to_dest_ranks() const

Compute map from each local (owned) index to the set of ranks that have the index as a ghost.

Todo:

Aim to remove this function?

Returns

shared indices

std::vector<std::int32_t> shared_indices() const

Build a list of owned indices that are ghosted by another rank.

Returns

The local index of owned indices that are ghosts on other rank(s). The indicies are unique and sorted.

const std::vector<int> &src() const noexcept

Ordered set of MPI ranks that own caller’s ghost indices.

Typically used when creating neighbourhood communicators.

Returns

MPI ranks than own ghost indices. The ranks are unique and sorted.

const std::vector<int> &dest() const noexcept

Ordered set of MPI ranks that ghost indices owned by caller.

Typically used when creating neighbourhood communicators.

Returns

MPI ranks than own ghost indices. The ranks are unique and sorted.

bool overlapped() const noexcept

Check if index map has overlaps (ghosts on any rank).

The return value of this function is determined by which constructor was used to create the index map.

Returns

True if index map has overlaps on any ranks, otherwise false.

class Scatterer
#include <Scatterer.h>

A Scatterer supports the MPI scattering and gathering of data that is associated with a common::IndexMap.

Scatter and gather operations uses MPI neighbourhood collectives. The implementation is designed is for sparse communication patterns, as it typical of patterns based on and IndexMap.

Public Functions

Scatterer(const IndexMap &map, int bs)

Create a scatterer.

Parameters
  • map[in] The index map that describes the parallel layout of data

  • bs[in] The block size of data associated with each index in map that will be scattered/gathered

template<typename T>
inline void scatter_fwd_begin(const std::span<const T> &send_buffer, const std::span<T> &recv_buffer, MPI_Request &request) const

Start a non-blocking send of owned data to ranks that ghost the data.

The communication is completed by calling Scatterer::scatter_fwd_end. The send and receive buffer should not be changed until after Scatterer::scatter_fwd_end has been called.

Parameters
  • send_buffer[in] Local data associated with each owned local index to be sent to process where the data is ghosted. It must not be changed until after a call to Scatterer::scatter_fwd_end. The order of data in the buffer is given by Scatterer::local_indices.

  • recv_buffer – A buffer used for the received data. The position of ghost entries in the buffer is given by Scatterer::remote_indices. The buffer must not be accessed or changed until after a call to Scatterer::scatter_fwd_end.

  • request – The MPI request handle for tracking the status of the non-blocking communication

void scatter_fwd_end(MPI_Request &request) const

Complete a non-blocking send from the local owner to process ranks that have the index as a ghost.

This function completes the communication started by Scatterer::scatter_fwd_begin.

Parameters

request[in] The MPI request handle for tracking the status of the send

template<typename T, typename Functor>
inline void scatter_fwd_begin(const std::span<const T> &local_data, std::span<T> local_buffer, std::span<T> remote_buffer, Functor pack_fn, MPI_Request &request) const

Scatter data associated with owned indices to ghosting ranks.

Note

This function is intended for advanced usage, and in particular when using CUDA/device-aware MPI.

Parameters
  • local_data[in] All data associated with owned indices. Size is size_local() from the IndexMap used to create the scatterer, multiplied by the block size. The data for each index is blocked.

  • local_buffer[in] Working buffer. The required size is given by Scatterer::local_buffer_size.

  • remote_buffer[out] Working buffer. The required size is given by Scatterer::remote_buffer_size.

  • pack_fn[in] Function to pack data from local_data into the send buffer. It is passed as an argument to support CUDA/device-aware MPI.

  • request[in] The MPI request handle for tracking the status of the send

template<typename T, typename Functor>
inline void scatter_fwd_end(const std::span<const T> &remote_buffer, std::span<T> remote_data, Functor unpack_fn, MPI_Request &request) const

Complete a non-blocking send from the local owner to process ranks that have the index as a ghost, and unpack received buffer into remote data.

This function completes the communication started by Scatterer::scatter_fwd_begin.

Parameters
  • remote_buffer[in] Working buffer, same used in Scatterer::scatter_fwd_begin.

  • remote_data[out] Received data associated with the ghost indices. The order follows the order of the ghost indices in the IndexMap used to create the scatterer. The size equal to the number of ghosts in the index map multiplied by the block size. The data for each index is blocked.

  • unpack_fn[in] Function to unpack the received buffer into remote_data. It is passed as an argument to support CUDA/device-aware MPI.

  • request[in] The MPI request handle for tracking the status of the send

template<typename T>
inline void scatter_fwd(const std::span<const T> &local_data, std::span<T> remote_data) const

Scatter data associated with owned indices to ghosting ranks.

Parameters
  • local_data[in] All data associated with owned indices. Size is size_local() from the IndexMap used to create the scatterer, multiplied by the block size. The data for each index is blocked

  • remote_data[out] Received data associated with the ghost indices. The order follows the order of the ghost indices in the IndexMap used to create the scatterer. The size equal to the number of ghosts in the index map multiplied by the block size. The data for each index is blocked.

template<typename T>
inline void scatter_rev_begin(const std::span<const T> &send_buffer, const std::span<T> &recv_buffer, MPI_Request &request) const

Start a non-blocking send of ghost data to ranks that own the data.

The communication is completed by calling Scatterer::scatter_rev_end. The send and receive buffers should not be changed until after Scatterer::scatter_rev_end has been called.

The buffer must not be accessed or changed until after a call to Scatterer::scatter_fwd_end.

Parameters
  • send_buffer[in] Data associated with each ghost index. This data is sent to process that owns the index. It must not be changed until after a call to Scatterer::scatter_ref_end.

  • recv_buffer – Buffer used for the received data. The position of owned indices in the buffer is given by Scatterer::local_indices. Scatterer::local_displacements()[i] is the location of the first entry in recv_buffer received from neighbourhood rank i. The number of entries received from neighbourhood rank i is Scatterer::local_displacements()[i + 1] - Scatterer::local_displacements()[i]. recv_buffer[j] is the data associated with the index Scatterer::local_indices()[j] in the index map.

  • request – The MPI request handle for tracking the status of the non-blocking communication

void scatter_rev_end(MPI_Request &request) const

End the reverse scatter communication.

This function must be called after Scatterer::scatter_rev_begin. The buffers passed to Scatterer::scatter_rev_begin must not be modified until after the function has been called.

Parameters

request[in] The handle used when calling Scatterer::scatter_rev_begin

template<typename T, typename Functor>
inline void scatter_rev_begin(const std::span<const T> &remote_data, std::span<T> remote_buffer, std::span<T> local_buffer, Functor pack_fn, MPI_Request &request) const

Scatter data associated with ghost indices to owning ranks.

Note

This function is intended for advanced usage, and in particular when using CUDA/device-aware MPI.

Template Parameters
  • T – The data type to send

  • BinaryOp – The reduction to perform when reducing data received from ghosting ranks to the value associated with the index on the owner

  • Functor1 – The pack function

  • Functor2 – The unpack function

Parameters
  • remote_data[in] Received data associated with the ghost indices. The order follows the order of the ghost indices in the IndexMap used to create the scatterer. The size equal to the number of ghosts in the index map multiplied by the block size. The data for each index is blocked.

  • local_buffer[out] Working buffer. The requires size is given by Scatterer::local_buffer_size.

  • remote_buffer[out] Working buffer. The requires size is given by Scatterer::remote_buffer_size.

  • pack_fn[in] Function to pack data from local_data into the send buffer. It is passed as an argument to support CUDA/device-aware MPI.

  • request – The MPI request handle for tracking the status of the non-blocking communication

template<typename T, typename Functor, typename BinaryOp>
inline void scatter_rev_end(const std::span<const T> &local_buffer, std::span<T> local_data, Functor unpack_fn, BinaryOp op, MPI_Request &request)

End the reverse scatter communication, and unpack the received local buffer into local data.

This function must be called after Scatterer::scatter_rev_begin. The buffers passed to Scatterer::scatter_rev_begin must not be modified until after the function has been called.

Parameters
  • local_buffer[in] Working buffer. Same buffer should be used in Scatterer::scatter_rev_begin.

  • local_data[out] All data associated with owned indices. Size is size_local() from the IndexMap used to create the scatterer, multiplied by the block size. The data for each index is blocked.

  • unpack_fn[in] Function to unpack the receive buffer into local_data. It is passed as an argument to support CUDA/device-aware MPI.

  • op[in] The reduction operation when accumulating received values. To add the received values use std::plus<T>().

  • request[in] The handle used when calling Scatterer::scatter_rev_begin

template<typename T, typename BinaryOp>
inline void scatter_rev(std::span<T> local_data, const std::span<const T> &remote_data, BinaryOp op)

Scatter data associated with ghost indices to ranks that own the indices.

std::int32_t local_buffer_size() const noexcept

Size of buffer for local data (owned and shared) used in forward and reverse communication.

Returns

The required buffer size

std::int32_t remote_buffer_size() const noexcept

Buffer size for remote data (ghosts) used in forward and reverse communication.

Returns

The required buffer size

const std::vector<std::int32_t> &local_indices() const noexcept

Return a vector of local indices (owned) used to pack/unpack local data. These indices are grouped by neighbor process (process for which an index is a ghost).

const std::vector<std::int32_t> &remote_indices() const noexcept

Return a vector of remote indices (ghosts) used to pack/unpack ghost data. These indices are grouped by neighbor process (ghost owners).

int bs() const noexcept

The number values (block size) to send per index in the common::IndexMap use to create the scatterer.

Returns

The block size

class TimeLogger
#include <TimeLogger.h>

Timer logging.

Public Functions

TimeLogger() = default

Constructor.

~TimeLogger() = default

Destructor.

void register_timing(std::string task, double wall, double user, double system)

Register timing (for later summary)

Table timings(std::set<TimingType> type)

Return a summary of timings and tasks in a Table.

void list_timings(MPI_Comm comm, std::set<TimingType> type, Table::Reduction reduction)

List a summary of timings and tasks. Reduction type is printed.

Parameters
  • comm – MPI Communicator

  • type – Set of possible timings: wall, user or system

  • reduction – Reduction type (min, max or average)

std::tuple<int, double, double, double> timing(std::string task)

Return timing.

Parameters

task[in] The task name to retrieve the timing for

Returns

Values (count, total wall time, total user time, total system time) for given task.

class TimeLogManager
#include <TimeLogManager.h>

Logger initialisation.

Public Static Functions

static TimeLogger &logger()

Singleton instance of logger.

class Timer
#include <Timer.h>

A timer can be used for timing tasks. The basic usage is.

Timer timer(“Assembling over cells”);

The timer is started at construction and timing ends when the timer is destroyed (goes out of scope). It is also possible to start and stop a timer explicitly by

timer.start(); timer.stop();

Timings are stored globally and a summary may be printed by calling

list_timings();

Public Functions

Timer()

Create timer without logging.

Timer(const std::string &task)

Create timer with logging.

~Timer()

Destructor.

void start()

Zero and start timer.

void resume()

Resume timer. Not well-defined for logging timer.

double stop()

Stop timer, return wall time elapsed and store timing data into logger.

std::array<double, 3> elapsed() const

Return wall, user and system time in seconds.

namespace impl

Functions

template<int N, class InputIt, class OutputIt>
inline void copy_N(InputIt first, OutputIt result)

std::copy_n-type function with compile-time loop bound