DOLFINx 0.8.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
Public Types | Public Member Functions | List of all members
MatrixCSR< Scalar, Container, ColContainer, RowPtrContainer > Class Template Reference

#include <MatrixCSR.h>

Public Types

using value_type = Scalar
 Scalar type.
 
using container_type = Container
 Matrix entries container type.
 
using column_container_type = ColContainer
 Column index container type.
 
using rowptr_container_type = RowPtrContainer
 Row pointer container type.
 

Public Member Functions

template<int BS0 = 1, int BS1 = 1>
auto mat_set_values ()
 Insertion functor for setting values in a matrix. It is typically used in finite element assembly functions.
 
template<int BS0 = 1, int BS1 = 1>
auto mat_add_values ()
 Insertion functor for adding values to a matrix. It is typically used in finite element assembly functions.
 
 MatrixCSR (const SparsityPattern &p, BlockMode mode=BlockMode::compact)
 Create a distributed matrix.
 
 MatrixCSR (MatrixCSR &&A)=default
 
void set (value_type x)
 Set all non-zero local entries to a value including entries in ghost rows.
 
template<int BS0, int BS1>
void set (std::span< const value_type > x, std::span< const std::int32_t > rows, std::span< const std::int32_t > cols)
 Set values in the matrix.
 
template<int BS0 = 1, int BS1 = 1>
void add (std::span< const value_type > x, std::span< const std::int32_t > rows, std::span< const std::int32_t > cols)
 Accumulate values in the matrix.
 
std::int32_t num_owned_rows () const
 Number of local rows excluding ghost rows.
 
std::int32_t num_all_rows () const
 Number of local rows including ghost rows.
 
std::vector< value_typeto_dense () const
 
void scatter_rev ()
 Transfer ghost row data to the owning ranks accumulating received values on the owned rows, and zeroing any existing data in ghost rows. This process is analogous to scatter_rev for Vector except that the values are always accumulated on the owning process.
 
void scatter_rev_begin ()
 Begin transfer of ghost row data to owning ranks, where it will be accumulated into existing owned rows.
 
void scatter_rev_end ()
 End transfer of ghost row data to owning ranks.
 
double squared_norm () const
 Compute the Frobenius norm squared across all processes.
 
std::shared_ptr< const common::IndexMapindex_map (int dim) const
 Index maps for the row and column space.
 
container_typevalues ()
 
const container_typevalues () const
 
const rowptr_container_typerow_ptr () const
 
const column_container_typecols () const
 
const rowptr_container_typeoff_diag_offset () const
 
std::array< int, 2 > block_size () const
 

Detailed Description

template<class Scalar, class Container = std::vector<Scalar>, class ColContainer = std::vector<std::int32_t>, class RowPtrContainer = std::vector<std::int64_t>>
class dolfinx::la::MatrixCSR< Scalar, Container, ColContainer, RowPtrContainer >

Distributed sparse matrix

The matrix storage format is compressed sparse row. The matrix is partitioned row-wise across MPI ranks.

Warning
Experimental storage of a matrix in CSR format which can be assembled into using the usual DOLFINx assembly routines. Matrix internal data can be accessed for interfacing with other code.
Template Parameters
ScalarScalar type of matrix entries
ContainerSequence container type to store matrix entries
ColContainerColumn index container type
RowPtrContainerRow pointer container type

Constructor & Destructor Documentation

◆ MatrixCSR() [1/2]

template<class U , class V , class W , class X >
MatrixCSR ( const SparsityPattern & p,
BlockMode mode = BlockMode::compact )

Create a distributed matrix.

The structure of the matrix depends entirely on the input SparsityPattern, which must be finalized. The matrix storage is distributed Compressed Sparse Row: the matrix is distributed by row across processes, and on each process, there is a list of column indices and matrix entries for each row stored. This exactly matches the layout of the SparsityPattern. There is some overlap of matrix rows between processes to allow for independent Finite Element assembly, after which, the ghost rows should be sent to the row owning processes by calling scatter_rev().

Note
The block size of the matrix is given by the block size of the input SparsityPattern.
Parameters
[in]pThe sparsity pattern which describes the parallel distribution and the non-zero structure.
[in]modeBlock mode. When the block size is greater than one, the storage can be "compact" where each matrix entry refers to a block of data (stored row major), or "expanded" where each matrix entry is individual. In the "expanded" case, the sparsity is expanded for every entry in the block, and the block size of the matrix is set to (1, 1).

◆ MatrixCSR() [2/2]

template<class Scalar , class Container = std::vector<Scalar>, class ColContainer = std::vector<std::int32_t>, class RowPtrContainer = std::vector<std::int64_t>>
MatrixCSR ( MatrixCSR< Scalar, Container, ColContainer, RowPtrContainer > && A)
default

Move constructor

Todo
Check handling of MPI_Request

Member Function Documentation

◆ add()

template<class Scalar , class Container = std::vector<Scalar>, class ColContainer = std::vector<std::int32_t>, class RowPtrContainer = std::vector<std::int64_t>>
template<int BS0 = 1, int BS1 = 1>
void add ( std::span< const value_type > x,
std::span< const std::int32_t > rows,
std::span< const std::int32_t > cols )
inline

Accumulate values in the matrix.

Note
Only entries included in the sparsity pattern used to initialize the matrix can be accumulated into.
All indices are local to the calling MPI rank and entries may go into ghost rows.
Use scatter_rev after all entries have been added to send ghost rows to owners. Adding more entries after scatter_rev is allowed, but another call to scatter_rev will then be required.
Template Parameters
BS0Row block size of data
BS1Column block size of data
Parameters
[in]xThe m by n dense block of values (row-major) to add to the matrix
[in]rowsThe row indices of x
[in]colsThe column indices of x

◆ block_size()

template<class Scalar , class Container = std::vector<Scalar>, class ColContainer = std::vector<std::int32_t>, class RowPtrContainer = std::vector<std::int64_t>>
std::array< int, 2 > block_size ( ) const
inline

Block size

Returns
block sizes for rows and columns

◆ cols()

template<class Scalar , class Container = std::vector<Scalar>, class ColContainer = std::vector<std::int32_t>, class RowPtrContainer = std::vector<std::int64_t>>
const column_container_type & cols ( ) const
inline

Get local column indices

Note
Includes columns in ghost rows

◆ index_map()

template<class Scalar , class Container = std::vector<Scalar>, class ColContainer = std::vector<std::int32_t>, class RowPtrContainer = std::vector<std::int64_t>>
std::shared_ptr< const common::IndexMap > index_map ( int dim) const
inline

Index maps for the row and column space.

The row IndexMap contains ghost entries for rows which may be inserted into and the column IndexMap contains all local and ghost columns that may exist in the owned rows.

Returns
Row (0) or column (1) index maps

◆ mat_add_values()

template<class Scalar , class Container = std::vector<Scalar>, class ColContainer = std::vector<std::int32_t>, class RowPtrContainer = std::vector<std::int64_t>>
template<int BS0 = 1, int BS1 = 1>
auto mat_add_values ( )
inline

Insertion functor for adding values to a matrix. It is typically used in finite element assembly functions.

Create a function to add values to a MatrixCSR. The function signature is int mat_add_fn(std::span<const std::int32_t rows, std::span<const std::int32_t cols, std::span<const value_type> data). The rows and columns use process local indexing, and the given rows and columns must pre-exist in the sparsity pattern of the matrix. Insertion into "ghost" rows (in the ghost region of the row IndexMap) is permitted, so long as there are correct entries in the sparsity pattern.

Note
Using rows or columns which are not in the sparsity will result in undefined behaviour (or an assert failure in Debug mode).
Matrix block size may be (1, 1) or (BS0, BS1)
Data block size may be (1, 1) or (BS0, BS1)
Template Parameters
BS0Row block size of data for insertion
BS1Column block size of data for insertion
Returns
Function for inserting values into A

◆ mat_set_values()

template<class Scalar , class Container = std::vector<Scalar>, class ColContainer = std::vector<std::int32_t>, class RowPtrContainer = std::vector<std::int64_t>>
template<int BS0 = 1, int BS1 = 1>
auto mat_set_values ( )
inline

Insertion functor for setting values in a matrix. It is typically used in finite element assembly functions.

Create a function to set values in a MatrixCSR. The function signature is int mat_set_fn(std::span<const std::int32_t rows, std::span<const std::int32_t cols, std::span<const value_type> data). The rows and columns use process local indexing, and the given rows and columns must pre-exist in the sparsity pattern of the matrix. Insertion into "ghost" rows (in the ghost region of the row IndexMap) is permitted, so long as there are correct entries in the sparsity pattern.

Note
Using rows or columns which are not in the sparsity will result in undefined behaviour (or an assert failure in Debug mode).
Matrix block size may be (1, 1) or (BS0, BS1)
Data block size may be (1, 1) or (BS0, BS1)
Template Parameters
BS0Row block size of data for insertion
BS1Column block size of data for insertion
Returns
Function for inserting values into A

◆ off_diag_offset()

template<class Scalar , class Container = std::vector<Scalar>, class ColContainer = std::vector<std::int32_t>, class RowPtrContainer = std::vector<std::int64_t>>
const rowptr_container_type & off_diag_offset ( ) const
inline

Get the start of off-diagonal (unowned columns) on each row, allowing the matrix to be split (virtually) into two parts. Operations (such as matrix-vector multiply) between the owned parts of the matrix and vector can then be performed separately from operations on the unowned parts.

Note
Includes ghost rows, which should be truncated manually if not required.

◆ row_ptr()

template<class Scalar , class Container = std::vector<Scalar>, class ColContainer = std::vector<std::int32_t>, class RowPtrContainer = std::vector<std::int64_t>>
const rowptr_container_type & row_ptr ( ) const
inline

Get local row pointers

Note
Includes pointers to ghost rows

◆ scatter_rev_begin()

template<typename U , typename V , typename W , typename X >
void scatter_rev_begin ( )

Begin transfer of ghost row data to owning ranks, where it will be accumulated into existing owned rows.

Note
Calls to this function must be followed by MatrixCSR::scatter_rev_end(). Between the two calls matrix values must not be changed.
This function does not change the matrix data. Data update only occurs with scatter_rev_end().

◆ scatter_rev_end()

template<typename U , typename V , typename W , typename X >
void scatter_rev_end ( )

End transfer of ghost row data to owning ranks.

Note
Must be preceded by MatrixCSR::scatter_rev_begin().
Matrix data received from other processes will be accumulated into locally owned rows, and ghost rows will be zeroed.

◆ set() [1/2]

template<class Scalar , class Container = std::vector<Scalar>, class ColContainer = std::vector<std::int32_t>, class RowPtrContainer = std::vector<std::int64_t>>
template<int BS0, int BS1>
void set ( std::span< const value_type > x,
std::span< const std::int32_t > rows,
std::span< const std::int32_t > cols )
inline

Set values in the matrix.

Note
Only entries included in the sparsity pattern used to initialize the matrix can be set.
All indices are local to the calling MPI rank and entries cannot be set in ghost rows.
This should be called after scatter_rev. Using before scatter_rev will set the values correctly, but incoming values may get added to them during a subsequent reverse scatter operation.
Template Parameters
BS0Data row block size
BS1Data column block size
Parameters
[in]xThe m by n dense block of values (row-major) to set in the matrix
[in]rowsThe row indices of x
[in]colsThe column indices of x

◆ set() [2/2]

template<class Scalar , class Container = std::vector<Scalar>, class ColContainer = std::vector<std::int32_t>, class RowPtrContainer = std::vector<std::int64_t>>
void set ( value_type x)
inline

Set all non-zero local entries to a value including entries in ghost rows.

Parameters
[in]xThe value to set non-zero matrix entries to

◆ squared_norm()

template<typename U , typename V , typename W , typename X >
double squared_norm ( ) const

Compute the Frobenius norm squared across all processes.

Note
MPI Collective

◆ to_dense()

template<typename U , typename V , typename W , typename X >
std::vector< typename MatrixCSR< U, V, W, X >::value_type > to_dense ( ) const

Copy to a dense matrix

Note
This function is typically used for debugging and not used in production
Ghost rows are also returned, and these can be truncated manually by using num_owned_rows() if required.
If the block size is greater than 1, the entries are expanded.
Returns
Dense copy of the part of the matrix on the calling rank. Storage is row-major.

◆ values() [1/2]

template<class Scalar , class Container = std::vector<Scalar>, class ColContainer = std::vector<std::int32_t>, class RowPtrContainer = std::vector<std::int64_t>>
container_type & values ( )
inline

Get local data values

Note
Includes ghost values

◆ values() [2/2]

template<class Scalar , class Container = std::vector<Scalar>, class ColContainer = std::vector<std::int32_t>, class RowPtrContainer = std::vector<std::int64_t>>
const container_type & values ( ) const
inline

Get local values (const version)

Note
Includes ghost values

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