KASKADE 7 development version
Public Member Functions | List of all members
Kaskade::NumaCRSPattern< Index > Class Template Reference

A NUMA-aware compressed row storage sparsity pattern. More...

#include <threadedMatrix.hh>

Detailed Description

template<class Index = size_t>
class Kaskade::NumaCRSPattern< Index >

A NUMA-aware compressed row storage sparsity pattern.

Template Parameters
Indexintegral type for row/column indices

Definition at line 1847 of file threadedMatrix.hh.

Public Member Functions

 NumaCRSPattern ()
 Constructs an empty 0x0 pattern. More...
 
 NumaCRSPattern (NumaCRSPatternCreator< Index > const &creator)
 Constructor creating a sparsity pattern from the given creator. More...
 
template<class Expanded , class Condensed , class Matrix >
 NumaCRSPattern (Expanded const &eIndices, Condensed const &cIndices, Matrix const &mat)
 Constructor. More...
 
template<class Matrix >
 NumaCRSPattern (Matrix const &matrix, bool isSymmetric, bool isTransposed, bool symmetric)
 Constructor extracting the sparsity pattern of a given matrix (usually a Dune::BCRSMatrix). More...
 
Index N () const
 The number of rows. More...
 
Index M () const
 The number of columns. More...
 
size_t storage () const
 Returns the number of stored entries. More...
 
size_t nonzeroes () const
 Returns the number of structurally nonzero entries. More...
 
int nodes () const
 Returns the number of NUMA nodes/chunks used. More...
 
std::shared_ptr< ChunkPatternpattern (int i) const
 Returns the individual patterns. More...
 
bool isSymmetric () const
 Returns the symmetry status of the pattern. More...
 
int chunk (Index row) const
 Returns the number of the chunk containing the given row. More...
 
std::vector< Index > const & rowStart () const
 Returns the limiting row indices between the chunks. More...
 
bool exists (Index r, Index c) const
 queries whether an entry is present in the sparsity pattern More...
 

Constructor & Destructor Documentation

◆ NumaCRSPattern() [1/4]

template<class Index = size_t>
Kaskade::NumaCRSPattern< Index >::NumaCRSPattern ( )
inline

Constructs an empty 0x0 pattern.

Definition at line 1856 of file threadedMatrix.hh.

◆ NumaCRSPattern() [2/4]

template<class Index = size_t>
Kaskade::NumaCRSPattern< Index >::NumaCRSPattern ( NumaCRSPatternCreator< Index > const &  creator)
inline

Constructor creating a sparsity pattern from the given creator.

Definition at line 1863 of file threadedMatrix.hh.

◆ NumaCRSPattern() [3/4]

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

Constructor.

This works like Matlab A(idx,idx), where idx is given by eIndices.

Template Parameters
Expandedan array type with value type convertible to Index
Condensedan array type with value type convertible to Index
Parameters
eIndicessorted global array of expanded indices
cIndicessorted global array of condensed indices
mata matrix with BCRSMatrix interface

Definition at line 1888 of file threadedMatrix.hh.

◆ NumaCRSPattern() [4/4]

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

Constructor extracting the sparsity pattern of a given matrix (usually a Dune::BCRSMatrix).

Template Parameters
Matrixthe type of the supplied matrix to copy
Parameters
matrixthe matrix to be copied
isSymmetricwhether the supplied matrix is symmetric
isTransposedwhether the supplied matrix is transposed
symmetricwhether only the lower triangular part should be stored

This copies the sparsity pattern of the provided Dune ISTL matrix. The number of chunks created is the number of NUMA nodes as reported by the NumaThreadPool.

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

Definition at line 1923 of file threadedMatrix.hh.

Member Function Documentation

◆ chunk()

template<class Index = size_t>
int Kaskade::NumaCRSPattern< Index >::chunk ( Index  row) const
inline

Returns the number of the chunk containing the given row.

Parameters
rowthe index of the row in the range [0,N]

For row==N, the total number of chunks is returned (i.e. one behind the last chunk).

Definition at line 2003 of file threadedMatrix.hh.

Referenced by Kaskade::NumaCRSPattern< Index >::exists().

◆ exists()

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

queries whether an entry is present in the sparsity pattern

Definition at line 2024 of file threadedMatrix.hh.

◆ isSymmetric()

template<class Index = size_t>
bool Kaskade::NumaCRSPattern< Index >::isSymmetric ( ) const
inline

Returns the symmetry status of the pattern.

If true, the matrix is symmetric, and only the lower triangular entries are actually stored.

Definition at line 1995 of file threadedMatrix.hh.

Referenced by Kaskade::NumaCRSPattern< Index >::NumaCRSPattern(), and Kaskade::operator+().

◆ M()

template<class Index = size_t>
Index Kaskade::NumaCRSPattern< Index >::M ( ) const
inline

The number of columns.

Definition at line 1954 of file threadedMatrix.hh.

Referenced by Kaskade::operator+().

◆ N()

template<class Index = size_t>
Index Kaskade::NumaCRSPattern< Index >::N ( ) const
inline

The number of rows.

Definition at line 1949 of file threadedMatrix.hh.

Referenced by Kaskade::NumaCRSPattern< Index >::chunk(), and Kaskade::operator+().

◆ nodes()

template<class Index = size_t>
int Kaskade::NumaCRSPattern< Index >::nodes ( ) const
inline

Returns the number of NUMA nodes/chunks used.

Definition at line 1983 of file threadedMatrix.hh.

Referenced by Kaskade::NumaCRSPattern< Index >::NumaCRSPattern(), and Kaskade::operator+().

◆ nonzeroes()

template<class Index = size_t>
size_t Kaskade::NumaCRSPattern< Index >::nonzeroes ( ) const
inline

Returns the number of structurally nonzero entries.

Definition at line 1972 of file threadedMatrix.hh.

Referenced by Kaskade::operator+().

◆ pattern()

template<class Index = size_t>
std::shared_ptr< ChunkPattern > Kaskade::NumaCRSPattern< Index >::pattern ( int  i) const
inline

Returns the individual patterns.

Definition at line 1988 of file threadedMatrix.hh.

Referenced by Kaskade::operator+().

◆ rowStart()

template<class Index = size_t>
std::vector< Index > const & Kaskade::NumaCRSPattern< Index >::rowStart ( ) const
inline

Returns the limiting row indices between the chunks.

For \( 0\le i < n \), the chunk \( i \) contains the rows in the half-open range [ rowStart()[i], rowStart()[i+1] [.

Definition at line 2019 of file threadedMatrix.hh.

◆ storage()

template<class Index = size_t>
size_t Kaskade::NumaCRSPattern< Index >::storage ( ) const
inline

Returns the number of stored entries.

This differs from nonzeroes() for symmetric storage.

Definition at line 1961 of file threadedMatrix.hh.


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