KASKADE 7 development version
Public Types | Related Functions | List of all members
Kaskade::NumaBCRSMatrix< Entry, Index > Class Template Reference

A NUMA-aware compressed row storage matrix adhering mostly to the Dune ISTL interface (to complete...) More...

#include <threadedMatrix.hh>

Detailed Description

template<class Entry, class Index = size_t>
class Kaskade::NumaBCRSMatrix< Entry, Index >

A NUMA-aware compressed row storage matrix adhering mostly to the Dune ISTL interface (to complete...)

This class distributes the matrix rows block-wise across available NUMA nodes to exploit the higher memory bandwidth in parallelized matrix-vector products and similar operations.

As the distribution of matrix entries to the NUMA nodes is optimized for equidistribution of storage size and computational effort in matrix vector products, transpose times vector multiplications are much slower than plain matrix-vector multiplications (roughly factor 6).

If the matrix is a priori specified as being symmetric, only the lower triangular elements are stored and can be accessed. Nevertheless, all linear algebra operations (i.e., matrix vector products) act on the complete symmetric matrix. Note that symmetric storage incurs a significant performance penalty for matrix-vector products (roughly factor 3), and should only be used if those operations are rare. In contrast, transpose times vector multiplications do not suffer.

For creating a NumaBCRSMatrix, a sparsity pattern (NumaCRSPattern) has to be provided. This in turn is created using a creator object (NumaCRSPatternCreator). This multi-stage creation allows efficient creation while avoiding stateful interfaces (incompletely created objects). The matrix creation can be done as follows:

// define the pattern creator for a 9x9 matrix
// add a diagonal entry
creator.addElement(2,2);
// create the matrix
// fill the entry
matrix[2][2] = 1.0;
A NUMA-aware compressed row storage matrix adhering mostly to the Dune ISTL interface (to complete....
A NUMA-aware creator for matrix sparsity patterns.

Concurrent calls to one NumaBCRSMatrix object are not safe.

Template Parameters
Entrythe type of matrix elements (in general Dune::FieldMatrix<double,n,m>)
Indexintegral type for row/column indices (usually int, long, or size_t), defaults to size_t

Definition at line 2114 of file threadedMatrix.hh.

Public Types

typedef ScalarType< Entry > Scalar
 
typedef Scalar field_type
 
typedef Entry block_type
 
using value_type = Entry
 
using size_type = Index
 
typedef ThreadedMatrixDetail::NumaBCRSMatrixIterator< Entry, Index > iterator
 iterator type stepping through the rows More...
 
typedef ThreadedMatrixDetail::NumaBCRSMatrixConstIterator< Entry, Index > const_iterator
 
typedef iterator row_type
 
typedef const_iterator const_row_type
 
typedef iterator ConstRowIterator
 
using ColIterator = typename iterator::Iterator
 column iterator stepping through the entries of a row More...
 
using ConstColIterator = typename iterator::ConstIterator
 column iterator stepping through the const entries of a row More...
 

Public Member Functions

Assignment and simple in-place modifications
Selfoperator= (Self const &mat)=default
 Copy assignment. More...
 
Selfoperator= (Self &&mat)=default
 Move assignment. More...
 
Selfoperator= (Entry const &a)
 Assigns the given value to each entry. More...
 
Selfoperator= (typename Entry::field_type const &a)
 Assigns the given scalar value to each entry. More...
 
template<class Arguments , class Operation >
Selfoperator= (ThreadedMatrixDetail::NumaBCRSMatrixExpression< Arguments, Operation > const &e)
 NOT YET IMPLEMENTED Assigns the given Numa matrix expression. More...
 
template<class EntryB , class IndexB >
Selfoperator+= (NumaBCRSMatrix< EntryB, IndexB > const &B)
 Adds a sparse matrix to this one. More...
 
template<class EntryB , class IndexB >
Selfoperator-= (NumaBCRSMatrix< EntryB, IndexB > const &B)
 Subtracts a sparse matrix from this one. More...
 
template<class Factor >
Selfoperator*= (Factor a)
 Multiplication with a "scalar". More...
 
Element access

Matrix entries can be accessed via iterators or by direct access. Iterators allow row-major scanning of existing matrix entries, with an iterator iterating over the rows and for each row an iterator iterating over the entries in that row.

Direct access is provided by the subscript operator [], which returns a row (reference), which in turn provides a subscript operator for selecting the column, e.g.,

A[i][j] = 42;

Note that the second (column) subscript has logarithmic complexity in the number of entries in the row.

iterator begin ()
 returns an iterator to the first row More...
 
const_iterator begin () const
 
iterator end ()
 returns an iterator to the first row More...
 
const_iterator end () const
 
row_type operator[] (Index r)
 Subscript operator allowing random access to rows. More...
 
const_row_type operator[] (Index r) const
 
Submatrix extraction
template<class RowIndices , class ColIndices >
Self operator() (RowIndices const &ri, ColIndices const &ci) const
 
General matrix information.
Index N () const
 The number of rows. More...
 
Index M () const
 The number of columns. More...
 
size_t nonzeroes () const
 Returns the number of structurally nonzero elements. More...
 
std::shared_ptr< NumaCRSPattern< Index > > getPattern () const
 Returns a pointer to the sparsity pattern. More...
 
Chunkchunk (int i)
 Obtains a reference to the given chunk. More...
 
bool exists (Index r, Index c) const
 returns true if (r,c) is structurally nonzero More...
 
field_type frobenius_norm () const
 Computes the Frobenius norm \( \|A\|_F = \sqrt{\mathrm{tr}(A^TA)} \). More...
 
field_type frobenius_norm2 () const
 Computes the square of the Frobenius norm \( \|A\|_F^2 = \mathrm{tr}(A^TA) \). More...
 
Matrix-vector operations
template<class X , class Y >
field_type mv (X const &x, Y &y) const
 Matrix-vector multiplication \( y \leftarrow Ax \) with computation of \( y^T x \) if A is square. More...
 
template<class X , class Y >
void mtv (X const &x, Y &y) const
 Matrix-vector multiplication \( y \leftarrow A^Tx \). More...
 
template<class X , class Y >
void mmv (X const &x, Y &y) const
 Matrix-vector multiplication \( y \leftarrow -Ax \). More...
 
template<class X , class Y >
field_type smv (field_type const &a, X const &x, Y &y) const
 Matrix-vector multiplication \( y \leftarrow aAx \). More...
 
template<class X , class Y >
void smtv (field_type const &a, X const &x, Y &y) const
 Matrix-vector multiplication \( y \leftarrow aA^Tx \). More...
 
template<class X , class Y >
field_type umv (X const &x, Y &y) const
 Matrix-vector multiplication \( y \leftarrow y + Ax \) and subsequent computation of \( y^T x \) if A is square. More...
 
template<class X , class Y >
void umtv (X const &x, Y &y) const
 Matrix-vector multiplication \( y \leftarrow y + A^Tx \). More...
 
template<class X , class Y >
field_type usmv (field_type const &a, X const &x, Y &y) const
 Matrix-vector multiplication \( y \leftarrow y + aAx \) and subsequent computation of \( y^T x \) if A is square. More...
 
template<class X , class Y >
void usmtv (field_type const &a, X const &x, Y &y) const
 Matrix-vector multiplication \( y \leftarrow y + aA^Tx \). More...
 

Related Functions

(Note that these are not member functions.)

template<class IndexP , class EntryP , class IndexA , class EntryA >
auto conjugation (NumaBCRSMatrix< EntryP, IndexP > const &P, NumaBCRSMatrix< EntryA, IndexA > const &A, bool onlyLowerTriangle=false, bool createDiagonal=false)
 Creates the conjugation product \( P^T A P\). More...
 
template<class Target , class Source , class RowIndices , class ColIndices >
Target submatrix (Source const &A, RowIndices const &ri, ColIndices const &ci)
 extracts a sparse or dense submatrix of \( A \) More...
 
template<class EntryA , class EntryB , class IndexA , class IndexB >
auto operator* (NumaBCRSMatrix< EntryA, IndexA > const &A, NumaBCRSMatrix< EntryB, IndexB > const &B)
 Computes the matrix-matrix product \( (A,B) \mapsto AB \). More...
 
template<class Entry , class Index >
std::ostream & operator<< (std::ostream &out, NumaBCRSMatrix< Entry, Index > const &A)
 Writes a NumaBCRSMatrix to an output stream. More...
 
template<class SparseMatrix >
DynamicMatrix< typename SparseMatrix::block_type > full (SparseMatrix const &A)
 Converts a sparse NumaBCRSMatrix to a dense matrix. More...
 
template<class SparseMatrix , class RowRange , class ColRange >
DynamicMatrix< typename SparseMatrix::block_type > full (SparseMatrix const &A, RowRange const &rows, ColRange const &cols)
 Converts a subrange of a sparse matrix to a dense matrix. More...
 
template<class Entry , class Index , class Index2 >
NumaBCRSMatrix< Entry, Index > operator+ (NumaBCRSMatrix< Entry, Index > const &A, NumaBCRSMatrix< Entry, Index2 > const &B)
 Matrix addition \( (A,B) \mapsto A+B \). The sparsity patterns of both matrices can be different. The size of the matrices have to be the same. Both must not have symmetric storage. More...
 
template<class Entry , class Index >
auto transpose (NumaBCRSMatrix< Entry, Index > const &A)
 Creates the transposed sparse matrix \( A^T \). More...
 
template<class Target , class Source , class RowIndices , class ColIndices >
Target eraseRowsNCols (Source const &A, RowIndices const &ri, ColIndices const &ci)
 "deletes" rows/columns by extracting a copy of the matrix keeping the non-deleted rows/columns. More...
 
template<class Target , class Source , class RowIndices >
Target eraseRows (Source const &A, RowIndices const &ri)
 "deletes" rows by extracting a copy of the matrix keeping only the non-deleted rows. More...
 
template<class Target , class Source , class ColIndices >
Target eraseCols (Source const &A, ColIndices const &ci)
 "deletes" columns by extracting a copy of the matrix keeping only the non-deleted columns. More...
 
template<class Source >
std::vector< size_t > nonZeroColumns (Source const &A)
 returns all indices of nonzero columns. More...
 
template<int row2, int col2, int row1, int col1, class Scalar , class Index >
NumaBCRSMatrix< Dune::FieldMatrix< Scalar, row2, col2 >, Index > reshapeBlocks (NumaBCRSMatrix< Dune::FieldMatrix< Scalar, row1, col1 >, Index > const &A)
 reshapes NumaBRCSMAtrix entry block structure More...
 
template<int blockrows, int blockcols, class Scalar , class Index >
NumaBCRSMatrix< Dune::FieldMatrix< Scalar, blockrows, blockcols >, Index > horzcat (NumaBCRSMatrix< Dune::FieldMatrix< Scalar, blockrows, blockcols >, Index > const &A, NumaBCRSMatrix< Dune::FieldMatrix< Scalar, blockrows, blockcols >, Index > const &B)
 concatenate two matrices horizontally More...
 
template<int blockrows, int blockcols, class Scalar , class Index >
NumaBCRSMatrix< Dune::FieldMatrix< Scalar, blockrows, blockcols >, Index > vertcat (NumaBCRSMatrix< Dune::FieldMatrix< Scalar, blockrows, blockcols >, Index > const &A, NumaBCRSMatrix< Dune::FieldMatrix< Scalar, blockrows, blockcols >, Index > const &B)
 concatenate two matrices vertically More...
 
template<int blockrows, int blockcols, class Index = size_t>
NumaBCRSMatrix< Dune::FieldMatrix< double, blockrows, blockcols >, Index > diagcat (NumaBCRSMatrix< Dune::FieldMatrix< double, blockrows, blockcols >, Index > const &A, NumaBCRSMatrix< Dune::FieldMatrix< double, blockrows, blockcols >, Index > const &B)
 concatenate two matrices diagonally, resulting matrix is zero on the off block-diagonal More...
 

Constructors

 NumaBCRSMatrix ()
 Constructs an empty 0x0 matrix. More...
 
 NumaBCRSMatrix (Self const &A)=default
 Copy constructor. More...
 
 NumaBCRSMatrix (Self &&A)=default
 Move constructor. More...
 
 NumaBCRSMatrix (NumaCRSPatternCreator< Index > const &creator, Entry const &init=Entry(0))
 Constructor creating a matrix from a given sparsity pattern creator. More...
 
 NumaBCRSMatrix (std::shared_ptr< NumaCRSPattern< Index > > const &pattern_, Entry const &init=Entry(0))
 Constructor creating a matrix from a given sparsity pattern. More...
 
template<class Expanded , class Condensed , class Matrix >
 NumaBCRSMatrix (Expanded const &eIndices, Condensed const &cIndices, Matrix const &mat)
 Deprecated. Use simpler version without cIndices instead. More...
 
template<class Expanded , class Matrix >
 NumaBCRSMatrix (Expanded const &eIndices, Matrix const &mat)
 Indexed submatrix constructor. More...
 
template<class Matrix >
 NumaBCRSMatrix (std::shared_ptr< NumaCRSPattern< Index > > const &pattern_, Matrix const &matrix, bool isSymmetric, bool isTransposed)
 Constructor copying a given matrix. More...
 
template<class OtherEntry >
 NumaBCRSMatrix (NumaBCRSMatrix< OtherEntry, Index > const &matrix)
 Constructor copying a given matrix. More...
 
template<class Matrix >
 NumaBCRSMatrix (Matrix const &matrix, bool isSymmetric, bool isTransposed=false, bool symmetric=false)
 Constructor. More...
 
 NumaBCRSMatrix (NumaBCRSMatrix< Entry, Index > const &matrix, bool isSymmetric, bool isTransposed=false, bool symmetric=false)
 Constructor. More...
 

Entering data

template<class LMIterator >
void scatter (LMIterator first, LMIterator last)
 Scatters given sub-matrices into the matrix by adding up their entries. More...
 
template<class RowIndices , class ColIndices , class BinaryOp = std::plus<Entry>>
void scatter (DynamicMatrix< Entry > const &B, RowIndices const &rows, ColIndices const &cols, BinaryOp const &op=BinaryOp())
 Scatters given submarix into the matrix. More...
 

Member Typedef Documentation

◆ block_type

template<class Entry , class Index = size_t>
typedef Entry Kaskade::NumaBCRSMatrix< Entry, Index >::block_type

Definition at line 2122 of file threadedMatrix.hh.

◆ ColIterator

template<class Entry , class Index = size_t>
using Kaskade::NumaBCRSMatrix< Entry, Index >::ColIterator = typename iterator::Iterator

column iterator stepping through the entries of a row

Definition at line 2137 of file threadedMatrix.hh.

◆ const_iterator

template<class Entry , class Index = size_t>
typedef ThreadedMatrixDetail::NumaBCRSMatrixConstIterator<Entry,Index> Kaskade::NumaBCRSMatrix< Entry, Index >::const_iterator

Definition at line 2130 of file threadedMatrix.hh.

◆ const_row_type

template<class Entry , class Index = size_t>
typedef const_iterator Kaskade::NumaBCRSMatrix< Entry, Index >::const_row_type

Definition at line 2132 of file threadedMatrix.hh.

◆ ConstColIterator

template<class Entry , class Index = size_t>
using Kaskade::NumaBCRSMatrix< Entry, Index >::ConstColIterator = typename iterator::ConstIterator

column iterator stepping through the const entries of a row

Definition at line 2140 of file threadedMatrix.hh.

◆ ConstRowIterator

template<class Entry , class Index = size_t>
typedef iterator Kaskade::NumaBCRSMatrix< Entry, Index >::ConstRowIterator

Definition at line 2134 of file threadedMatrix.hh.

◆ field_type

template<class Entry , class Index = size_t>
typedef Scalar Kaskade::NumaBCRSMatrix< Entry, Index >::field_type

Definition at line 2121 of file threadedMatrix.hh.

◆ iterator

template<class Entry , class Index = size_t>
typedef ThreadedMatrixDetail::NumaBCRSMatrixIterator<Entry,Index> Kaskade::NumaBCRSMatrix< Entry, Index >::iterator

iterator type stepping through the rows

Definition at line 2129 of file threadedMatrix.hh.

◆ row_type

template<class Entry , class Index = size_t>
typedef iterator Kaskade::NumaBCRSMatrix< Entry, Index >::row_type

Definition at line 2131 of file threadedMatrix.hh.

◆ Scalar

template<class Entry , class Index = size_t>
typedef ScalarType<Entry> Kaskade::NumaBCRSMatrix< Entry, Index >::Scalar

Definition at line 2120 of file threadedMatrix.hh.

◆ size_type

template<class Entry , class Index = size_t>
using Kaskade::NumaBCRSMatrix< Entry, Index >::size_type = Index

Definition at line 2124 of file threadedMatrix.hh.

◆ value_type

template<class Entry , class Index = size_t>
using Kaskade::NumaBCRSMatrix< Entry, Index >::value_type = Entry

Definition at line 2123 of file threadedMatrix.hh.

Constructor & Destructor Documentation

◆ NumaBCRSMatrix() [1/11]

template<class Entry , class Index = size_t>
Kaskade::NumaBCRSMatrix< Entry, Index >::NumaBCRSMatrix ( )
inline

Constructs an empty 0x0 matrix.

Definition at line 2151 of file threadedMatrix.hh.

◆ NumaBCRSMatrix() [2/11]

template<class Entry , class Index = size_t>
Kaskade::NumaBCRSMatrix< Entry, Index >::NumaBCRSMatrix ( Self const &  A)
default

Copy constructor.

◆ NumaBCRSMatrix() [3/11]

template<class Entry , class Index = size_t>
Kaskade::NumaBCRSMatrix< Entry, Index >::NumaBCRSMatrix ( Self &&  A)
default

Move constructor.

◆ NumaBCRSMatrix() [4/11]

template<class Entry , class Index = size_t>
Kaskade::NumaBCRSMatrix< Entry, Index >::NumaBCRSMatrix ( NumaCRSPatternCreator< Index > const &  creator,
Entry const &  init = Entry(0) 
)
inline

Constructor creating a matrix from a given sparsity pattern creator.

This is a convenience constructor that creates a new pattern under the hood.

Parameters
creatorthe sparsity pattern creatorof the matrix to be constructed
initinitialization value for the entries

Definition at line 2173 of file threadedMatrix.hh.

◆ NumaBCRSMatrix() [5/11]

template<class Entry , class Index = size_t>
Kaskade::NumaBCRSMatrix< Entry, Index >::NumaBCRSMatrix ( std::shared_ptr< NumaCRSPattern< Index > > const &  pattern_,
Entry const &  init = Entry(0) 
)
inline

Constructor creating a matrix from a given sparsity pattern.

Parameters
patternthe sparsity pattern of the matrix to be constructed
initinitialization value for the entries

Definition at line 2183 of file threadedMatrix.hh.

◆ NumaBCRSMatrix() [6/11]

template<class Entry , class Index = size_t>
template<class Expanded , class Condensed , class Matrix >
Kaskade::NumaBCRSMatrix< Entry, Index >::NumaBCRSMatrix ( Expanded const &  eIndices,
Condensed const &  cIndices,
Matrix const &  mat 
)
inline

Deprecated. Use simpler version without cIndices instead.

This constructor will be made private after 2022-12.

Definition at line 2203 of file threadedMatrix.hh.

◆ NumaBCRSMatrix() [7/11]

template<class Entry , class Index = size_t>
template<class Expanded , class Matrix >
Kaskade::NumaBCRSMatrix< Entry, Index >::NumaBCRSMatrix ( Expanded const &  eIndices,
Matrix const &  mat 
)
inline

Indexed submatrix constructor.

This works like Matlab A(idx,idx), where idx is given by eIndices. It creates a \( n \times n \) matrix \( a \) from the \( N \times \) N matrix \( A \), where \( a_{ij} = A_{e_ie_j} \). Of course, \( n \le N \) has to hold.

Template Parameters
Expandedan array type with value type convertible to Index
Parameters
eIndicessorted global array of size \( n \) of expanded indices
mata square matrix with BCRSMatrix interface

Definition at line 2242 of file threadedMatrix.hh.

◆ NumaBCRSMatrix() [8/11]

template<class Entry , class Index = size_t>
template<class Matrix >
Kaskade::NumaBCRSMatrix< Entry, Index >::NumaBCRSMatrix ( std::shared_ptr< NumaCRSPattern< Index > > const &  pattern_,
Matrix const &  matrix,
bool  isSymmetric,
bool  isTransposed 
)
inline

Constructor copying a given matrix.

Template Parameters
Matrixthe type of the supplied matrix to copy (usually a Dune::BCRSMatrix or a NumaBCRSMatrix with different scalar type).
Parameters
patternthe sparsity pattern to use
matrixthe matrix to be copied (shall have the given sparsity pattern)
isSymmetricwhether the supplied matrix is symmetric (then only its lower triangular part is accessed)
isTransposedwhether the supplied matrix is transposed

Definition at line 2258 of file threadedMatrix.hh.

◆ NumaBCRSMatrix() [9/11]

template<class Entry , class Index = size_t>
template<class OtherEntry >
Kaskade::NumaBCRSMatrix< Entry, Index >::NumaBCRSMatrix ( NumaBCRSMatrix< OtherEntry, Index > const &  matrix)
inline

Constructor copying a given matrix.

Template Parameters
OtherEntryThe entry type of the matrix to be copied. The entries shall have the same dimension (but possibly different scalar type).
Parameters
matrixthe matrix to be copied (shall conform to the Dune::BCRSMatrix interface)

Use this for switching between different entry types (e.g. float vs double). The sparsity pattern is shared between both matrices.

The number of chunks created is the number of NUMA nodes as reported by NumaThreadPool.

Definition at line 2281 of file threadedMatrix.hh.

◆ NumaBCRSMatrix() [10/11]

template<class Entry , class Index = size_t>
template<class Matrix >
Kaskade::NumaBCRSMatrix< Entry, Index >::NumaBCRSMatrix ( Matrix const &  matrix,
bool  isSymmetric,
bool  isTransposed = false,
bool  symmetric = false 
)
inline

Constructor.

Template Parameters
Matrixthe type of the supplied matrix to copy (shall conform to the Dune::BCRSMatrix interface)
Parameters
matrixthe matrix to be copied
isSymmetricwhether the supplied matrix is symmetric (then only its lower triangular part is accessed, requires square matrix and entries)
isTransposedwhether the supplied matrix is transposed
symmetricwhether only the lower triangular part should be stored (requires a square matrix and entries)

This is a convenience constructor that creates a new sparsity pattern from scratch. The number of chunks created is the number of NUMA nodes as reported by NumaThreadPool.

Note that !isSymmetric && symmetric is highly questionable and hence not allowed.

Definition at line 2302 of file threadedMatrix.hh.

◆ NumaBCRSMatrix() [11/11]

template<class Entry , class Index = size_t>
Kaskade::NumaBCRSMatrix< Entry, Index >::NumaBCRSMatrix ( NumaBCRSMatrix< Entry, Index > const &  matrix,
bool  isSymmetric,
bool  isTransposed = false,
bool  symmetric = false 
)
inline

Constructor.

Parameters
matrixthe matrix to be copied (shall conform to the Dune::BCRSMatrix interface)
isSymmetricwhether the supplied matrix is symmetric (then only its lower triangular part is accessed, requires square matrix and entries)
isTransposedwhether the supplied matrix is transposed
symmetricwhether only the lower triangular part should be stored (requires a square matrix and entries)

This is a convenience constructor that creates a new sparsity pattern from scratch. Use it for switching between symmetric/normal storage patterns and for creating transposed matrices.

The number of chunks created is the number of NUMA nodes as reported by NumaThreadPool.

Note that !isSymmetric && symmetric is highly questionable and hence not allowed.

Definition at line 2333 of file threadedMatrix.hh.

Member Function Documentation

◆ begin() [1/2]

template<class Entry , class Index = size_t>
iterator Kaskade::NumaBCRSMatrix< Entry, Index >::begin ( )
inline

◆ begin() [2/2]

template<class Entry , class Index = size_t>
const_iterator Kaskade::NumaBCRSMatrix< Entry, Index >::begin ( ) const
inline

Definition at line 2469 of file threadedMatrix.hh.

◆ chunk()

template<class Entry , class Index = size_t>
Chunk & Kaskade::NumaBCRSMatrix< Entry, Index >::chunk ( int  i)
inline

Obtains a reference to the given chunk.

Definition at line 2532 of file threadedMatrix.hh.

Referenced by Kaskade::ThreadedMatrixDetail::NumaBCRSRow< Entry, Index >::update().

◆ end() [1/2]

template<class Entry , class Index = size_t>
iterator Kaskade::NumaBCRSMatrix< Entry, Index >::end ( )
inline

◆ end() [2/2]

template<class Entry , class Index = size_t>
const_iterator Kaskade::NumaBCRSMatrix< Entry, Index >::end ( ) const
inline

Definition at line 2475 of file threadedMatrix.hh.

◆ exists()

template<class Entry , class Index = size_t>
bool Kaskade::NumaBCRSMatrix< Entry, Index >::exists ( Index  r,
Index  c 
) const
inline

returns true if (r,c) is structurally nonzero

Definition at line 2537 of file threadedMatrix.hh.

◆ frobenius_norm()

template<class Entry , class Index = size_t>
field_type Kaskade::NumaBCRSMatrix< Entry, Index >::frobenius_norm ( ) const
inline

Computes the Frobenius norm \( \|A\|_F = \sqrt{\mathrm{tr}(A^TA)} \).

Todo:
parallelize

Definition at line 2546 of file threadedMatrix.hh.

◆ frobenius_norm2()

template<class Entry , class Index = size_t>
field_type Kaskade::NumaBCRSMatrix< Entry, Index >::frobenius_norm2 ( ) const
inline

Computes the square of the Frobenius norm \( \|A\|_F^2 = \mathrm{tr}(A^TA) \).

Todo:
parallelize

Definition at line 2555 of file threadedMatrix.hh.

Referenced by Kaskade::NumaBCRSMatrix< Entry, Index >::frobenius_norm().

◆ getPattern()

template<class Entry , class Index = size_t>
std::shared_ptr< NumaCRSPattern< Index > > Kaskade::NumaBCRSMatrix< Entry, Index >::getPattern ( ) const
inline

◆ M()

template<class Entry , class Index = size_t>
Index Kaskade::NumaBCRSMatrix< Entry, Index >::M ( ) const
inline

◆ mmv()

template<class Entry , class Index = size_t>
template<class X , class Y >
void Kaskade::NumaBCRSMatrix< Entry, Index >::mmv ( X const &  x,
Y &  y 
) const
inline

Matrix-vector multiplication \( y \leftarrow -Ax \).

Definition at line 2594 of file threadedMatrix.hh.

◆ mtv()

template<class Entry , class Index = size_t>
template<class X , class Y >
void Kaskade::NumaBCRSMatrix< Entry, Index >::mtv ( X const &  x,
Y &  y 
) const
inline

Matrix-vector multiplication \( y \leftarrow A^Tx \).

Parameters
xa vector of size M()
ya vector of size N()

Note that transpose times vector operations are quite expensive (a factor of 6 slower than matrix-vector operations) and should be avoided. If you need transpose times vector products, consider storing a transposed matrix directly.

Definition at line 2588 of file threadedMatrix.hh.

◆ mv()

template<class Entry , class Index = size_t>
template<class X , class Y >
field_type Kaskade::NumaBCRSMatrix< Entry, Index >::mv ( X const &  x,
Y &  y 
) const
inline

Matrix-vector multiplication \( y \leftarrow Ax \) with computation of \( y^T x \) if A is square.

Parameters
xa vector of size M()
ya vector of size N()

Definition at line 2577 of file threadedMatrix.hh.

◆ N()

template<class Entry , class Index = size_t>
Index Kaskade::NumaBCRSMatrix< Entry, Index >::N ( ) const
inline

◆ nonzeroes()

template<class Entry , class Index = size_t>
size_t Kaskade::NumaBCRSMatrix< Entry, Index >::nonzeroes ( ) const
inline

Returns the number of structurally nonzero elements.

In symmetric storage, the number of actually stored entries is smaller, usually by a factor of almost two.

Definition at line 2522 of file threadedMatrix.hh.

◆ operator()()

template<class Entry , class Index = size_t>
template<class RowIndices , class ColIndices >
Self Kaskade::NumaBCRSMatrix< Entry, Index >::operator() ( RowIndices const &  ri,
ColIndices const &  ci 
) const
inline

Definition at line 2492 of file threadedMatrix.hh.

◆ operator*=()

template<class Entry , class Index = size_t>
template<class Factor >
Self & Kaskade::NumaBCRSMatrix< Entry, Index >::operator*= ( Factor  a)
inline

Multiplication with a "scalar".

The factor a can be a scalar or a quadratic matrix of compatible static size.

Definition at line 2435 of file threadedMatrix.hh.

◆ operator+=()

template<class Entry , class Index = size_t>
template<class EntryB , class IndexB >
Self & Kaskade::NumaBCRSMatrix< Entry, Index >::operator+= ( NumaBCRSMatrix< EntryB, IndexB > const &  B)
inline

Adds a sparse matrix to this one.

The sparsity pattern of B must be a subset of our sparsity pattern.

Definition at line 2411 of file threadedMatrix.hh.

◆ operator-=()

template<class Entry , class Index = size_t>
template<class EntryB , class IndexB >
Self & Kaskade::NumaBCRSMatrix< Entry, Index >::operator-= ( NumaBCRSMatrix< EntryB, IndexB > const &  B)
inline

Subtracts a sparse matrix from this one.

The sparsity pattern of B must be a subset of our sparsity pattern.

Definition at line 2423 of file threadedMatrix.hh.

◆ operator=() [1/5]

template<class Entry , class Index = size_t>
Self & Kaskade::NumaBCRSMatrix< Entry, Index >::operator= ( Entry const &  a)
inline

Assigns the given value to each entry.

Definition at line 2371 of file threadedMatrix.hh.

◆ operator=() [2/5]

template<class Entry , class Index = size_t>
Self & Kaskade::NumaBCRSMatrix< Entry, Index >::operator= ( Self &&  mat)
default

Move assignment.

The matrices need not be of same size or sparsity pattern.

◆ operator=() [3/5]

template<class Entry , class Index = size_t>
Self & Kaskade::NumaBCRSMatrix< Entry, Index >::operator= ( Self const &  mat)
default

Copy assignment.

The matrix becomes a copy of the given one, no matter what its previous structure or size was. The sparsity pattern is shared between both matrices.

◆ operator=() [4/5]

template<class Entry , class Index = size_t>
template<class Arguments , class Operation >
Self & Kaskade::NumaBCRSMatrix< Entry, Index >::operator= ( ThreadedMatrixDetail::NumaBCRSMatrixExpression< Arguments, Operation > const &  e)
inline

NOT YET IMPLEMENTED Assigns the given Numa matrix expression.

Definition at line 2394 of file threadedMatrix.hh.

◆ operator=() [5/5]

template<class Entry , class Index = size_t>
Self & Kaskade::NumaBCRSMatrix< Entry, Index >::operator= ( typename Entry::field_type const &  a)
inline

Assigns the given scalar value to each entry.

Definition at line 2385 of file threadedMatrix.hh.

◆ operator[]() [1/2]

template<class Entry , class Index = size_t>
row_type Kaskade::NumaBCRSMatrix< Entry, Index >::operator[] ( Index  r)
inline

Subscript operator allowing random access to rows.

This has constant complexity.

Definition at line 2482 of file threadedMatrix.hh.

◆ operator[]() [2/2]

template<class Entry , class Index = size_t>
const_row_type Kaskade::NumaBCRSMatrix< Entry, Index >::operator[] ( Index  r) const
inline

Definition at line 2483 of file threadedMatrix.hh.

◆ scatter()

template<class Entry , class Index = size_t>
template<class LMIterator >
void Kaskade::NumaBCRSMatrix< Entry, Index >::scatter ( LMIterator  first,
LMIterator  last 
)
inline

Scatters given sub-matrices into the matrix by adding up their entries.

Template Parameters
LMIteratora forward iterator with value type of LocalMatrix (requirements see below).

Local matrices \( a \) in the range [first,last[ are scattered into the global matrix as follows:

\[ A_{I_k,J_l} \leftarrow A_{I_k,J_l} + a_{i_k,j_l} \quad 0\le k < r, 0\le l < c \]

The possibly different ordering of entries in the local and global matrices is taken into account by the indirect indexing via \( (I_k,J_k) \) for the global indices and \( (i_k,j_l) \) for the local indices.

Objects a of LocalMatrix type have to provide the following access:

  • a(ik,jl) is the value of the local matrix entry \( a_{i_k,j_l} \)
  • ridx() returns a range of \( (I_k,i_k) \) pairs (of type std::pair), sorted ascendingly according to \( I_k \)
  • cidx() returns a range of \( (J_l,j_l) \) pairs (of type std::pair), sorted ascendingly according to \( J_l \)

In case LocalMatrix::lumped is true, then ridx()==cidx() shall hold and only the diagonal of \( a \) is to be scattered.

In contrast to the other methods, this method can (and should) be called concurrently from different threads.

Definition at line 2673 of file threadedMatrix.hh.

◆ smtv()

template<class Entry , class Index = size_t>
template<class X , class Y >
void Kaskade::NumaBCRSMatrix< Entry, Index >::smtv ( field_type const &  a,
X const &  x,
Y &  y 
) const
inline

Matrix-vector multiplication \( y \leftarrow aA^Tx \).

Note that transpose times vector operations are quite expensive (a factor of 6 slower than matrix-vector operations) and should be avoided. If you need transpose times vector products, consider storing a transposed matrix directly.

Definition at line 2609 of file threadedMatrix.hh.

◆ smv()

template<class Entry , class Index = size_t>
template<class X , class Y >
field_type Kaskade::NumaBCRSMatrix< Entry, Index >::smv ( field_type const &  a,
X const &  x,
Y &  y 
) const
inline

Matrix-vector multiplication \( y \leftarrow aAx \).

Definition at line 2600 of file threadedMatrix.hh.

◆ umtv()

template<class Entry , class Index = size_t>
template<class X , class Y >
void Kaskade::NumaBCRSMatrix< Entry, Index >::umtv ( X const &  x,
Y &  y 
) const
inline

Matrix-vector multiplication \( y \leftarrow y + A^Tx \).

Note that transpose times vector operations are quite expensive (a factor of 6 slower than matrix-vector operations) and should be avoided. If you need transpose times vector products, consider storing a transposed matrix directly.

Definition at line 2624 of file threadedMatrix.hh.

◆ umv()

template<class Entry , class Index = size_t>
template<class X , class Y >
field_type Kaskade::NumaBCRSMatrix< Entry, Index >::umv ( X const &  x,
Y &  y 
) const
inline

Matrix-vector multiplication \( y \leftarrow y + Ax \) and subsequent computation of \( y^T x \) if A is square.

Definition at line 2615 of file threadedMatrix.hh.

◆ usmtv()

template<class Entry , class Index = size_t>
template<class X , class Y >
void Kaskade::NumaBCRSMatrix< Entry, Index >::usmtv ( field_type const &  a,
X const &  x,
Y &  y 
) const
inline

Matrix-vector multiplication \( y \leftarrow y + aA^Tx \).

Note that transpose times vector operations are quite expensive (a factor of 6 slower than matrix-vector operations) and should be avoided. If you need transpose times vector products, consider storing a transposed matrix directly.

Definition at line 2639 of file threadedMatrix.hh.

◆ usmv()

template<class Entry , class Index = size_t>
template<class X , class Y >
field_type Kaskade::NumaBCRSMatrix< Entry, Index >::usmv ( field_type const &  a,
X const &  x,
Y &  y 
) const
inline

Matrix-vector multiplication \( y \leftarrow y + aAx \) and subsequent computation of \( y^T x \) if A is square.

Definition at line 2630 of file threadedMatrix.hh.

Friends And Related Function Documentation

◆ scatter()

template<class Entry , class Index = size_t>
template<class RowIndices , class ColIndices , class BinaryOp = std::plus<Entry>>
void scatter ( DynamicMatrix< Entry > const &  B,
RowIndices const &  rows,
ColIndices const &  cols,
BinaryOp const &  op = BinaryOp() 
)
related

Scatters given submarix into the matrix.

This performs \( A_{r_i c_j} = A_{r_i c_j} \circ B_{ij} \) for \( i,j \) indexing the provided matrix \( B \). Here, \( \circ \) denotes the given binary operation, which defaults to addition (the classical scatter operation).

Template Parameters
RowIndicesan STL sequence type with value_type convertible to Index
ColIndicesan STL sequence type with value_type convertible to Index
BinaryOpa binary callable type with op(Entry const, Entry const) convertible to Entry.
Parameters
Bthe given matrix data
rowsthe target row indices \( r_i \) for ourselves
colsthe target col indices \( c_i \) for ourselves
opan arbitrary binary operation on Entry, defaults to plus

Note that all entries \( A_{r_ic_i} \) must exist in this matrix.

As long as the row or column indices do not overlap, this method can safely be called concurrently.

Definition at line 2704 of file threadedMatrix.hh.


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