Reference documentation for deal.II version GIT relicensing-233-g802318d791 2024-03-28 20:20:02+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Static Public Member Functions | Private Types | Private Member Functions | Private Attributes | Static Private Attributes | List of all members
CUDAWrappers::SparseMatrix< Number > Class Template Reference

#include <deal.II/lac/cuda_sparse_matrix.h>

Inheritance diagram for CUDAWrappers::SparseMatrix< Number >:
Inheritance graph
[legend]

Public Types

using size_type = int
 
using value_type = Number
 
using real_type = Number
 

Public Member Functions

template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 
Constructors and initialization
 SparseMatrix ()
 
 SparseMatrix (Utilities::CUDA::Handle &handle, const ::SparseMatrix< Number > &sparse_matrix_host)
 
 SparseMatrix (CUDAWrappers::SparseMatrix< Number > &&)
 
 SparseMatrix (const CUDAWrappers::SparseMatrix< Number > &)=delete
 
 ~SparseMatrix ()
 
SparseMatrixoperator= (CUDAWrappers::SparseMatrix< Number > &&)
 
SparseMatrixoperator= (const CUDAWrappers::SparseMatrix< Number > &)=delete
 
void reinit (Utilities::CUDA::Handle &handle, const ::SparseMatrix< Number > &sparse_matrix_host)
 
Information on the matrix
size_type m () const
 
size_type n () const
 
std::size_t n_nonzero_elements () const
 
template<typename StreamType >
void print (StreamType &out, const bool across=false, const bool diagonal_first=true) const
 
void print_formatted (std::ostream &out, const unsigned int precision=3, const bool scientific=true, const unsigned int width=0, const char *zero_string=" ", const double denominator=1., const char *separator=" ") const
 
Modifying entries
SparseMatrixoperator*= (const Number factor)
 
SparseMatrixoperator/= (const Number factor)
 
Multiplications
void vmult (LinearAlgebra::CUDAWrappers::Vector< Number > &dst, const LinearAlgebra::CUDAWrappers::Vector< Number > &src) const
 
void Tvmult (LinearAlgebra::CUDAWrappers::Vector< Number > &dst, const LinearAlgebra::CUDAWrappers::Vector< Number > &src) const
 
void vmult_add (LinearAlgebra::CUDAWrappers::Vector< Number > &dst, const LinearAlgebra::CUDAWrappers::Vector< Number > &src) const
 
void Tvmult_add (LinearAlgebra::CUDAWrappers::Vector< Number > &dst, const LinearAlgebra::CUDAWrappers::Vector< Number > &src) const
 
Number matrix_norm_square (const LinearAlgebra::CUDAWrappers::Vector< Number > &v) const
 
Number matrix_scalar_product (const LinearAlgebra::CUDAWrappers::Vector< Number > &u, const LinearAlgebra::CUDAWrappers::Vector< Number > &v) const
 
Number residual (LinearAlgebra::CUDAWrappers::Vector< Number > &dst, const LinearAlgebra::CUDAWrappers::Vector< Number > &x, const LinearAlgebra::CUDAWrappers::Vector< Number > &b) const
 
Matrix norms
Number l1_norm () const
 
Number linfty_norm () const
 
Number frobenius_norm () const
 
Access to underlying CUDA data
std::tuple< Number *, int *, int *, cusparseMatDescr_t, cusparseSpMatDescr_t > get_cusparse_matrix () const
 
Subscriptor functionality

Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer class.

void subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
void unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
unsigned int n_subscriptions () const
 
template<typename StreamType >
void list_subscribers (StreamType &stream) const
 
void list_subscribers () const
 

Static Public Member Functions

static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 

Private Types

using map_value_type = decltype(counter_map)::value_type
 
using map_iterator = decltype(counter_map)::iterator
 

Private Member Functions

void check_no_subscribers () const noexcept
 

Private Attributes

cusparseHandle_t cusparse_handle
 
int nnz
 
int n_rows
 
int n_cols
 
std::unique_ptr< Number[], void(*)(Number *)> val_dev
 
std::unique_ptr< int[], void(*)(int *)> column_index_dev
 
std::unique_ptr< int[], void(*)(int *)> row_ptr_dev
 
cusparseMatDescr_t descr
 
cusparseSpMatDescr_t sp_descr
 
std::atomic< unsigned intcounter
 
std::map< std::string, unsigned intcounter_map
 
std::vector< std::atomic< bool > * > validity_pointers
 
const std::type_info * object_info
 

Static Private Attributes

static std::mutex mutex
 

Detailed Description

template<typename Number>
class CUDAWrappers::SparseMatrix< Number >

This class is a wrapper around cuSPARSE csr sparse matrix. Unlike deal.II's own SparseMatrix all elements within each row are stored in increasing column index order.

Note
Instantiations for this template are provided for <float> and <double>.

Definition at line 47 of file cuda_sparse_matrix.h.

Member Typedef Documentation

◆ size_type

template<typename Number >
using CUDAWrappers::SparseMatrix< Number >::size_type = int

Declare type for container size.

Definition at line 53 of file cuda_sparse_matrix.h.

◆ value_type

template<typename Number >
using CUDAWrappers::SparseMatrix< Number >::value_type = Number

Type of the matrix entries.

Definition at line 58 of file cuda_sparse_matrix.h.

◆ real_type

template<typename Number >
using CUDAWrappers::SparseMatrix< Number >::real_type = Number

Declare a type that holds real-valued numbers with the same precision as the template argument to this class.

Definition at line 64 of file cuda_sparse_matrix.h.

◆ map_value_type

using Subscriptor::map_value_type = decltype(counter_map)::value_type
privateinherited

The data type used in counter_map.

Definition at line 229 of file subscriptor.h.

◆ map_iterator

using Subscriptor::map_iterator = decltype(counter_map)::iterator
privateinherited

The iterator type used in counter_map.

Definition at line 234 of file subscriptor.h.

Constructor & Destructor Documentation

◆ SparseMatrix() [1/4]

template<typename Number >
SparseMatrix< Number >::SparseMatrix ( )

Constructor. Initialize the matrix to be empty, without any structure, i.e., the matrix is not usable at all. This constructor is therefore only useful for matrices which are members of a class.

You have to initialize the matrix before usage with reinit.

Definition at line 286 of file cuda_sparse_matrix.cc.

◆ SparseMatrix() [2/4]

template<typename Number >
SparseMatrix< Number >::SparseMatrix ( Utilities::CUDA::Handle handle,
const ::SparseMatrix< Number > &  sparse_matrix_host 
)

Constructor. Takes a Utilities::CUDA::Handle and a sparse matrix on the host. The sparse matrix on the host is copied on the device and the elements are reordered according to the format supported by cuSPARSE.

Definition at line 299 of file cuda_sparse_matrix.cc.

◆ SparseMatrix() [3/4]

template<typename Number >
SparseMatrix< Number >::SparseMatrix ( CUDAWrappers::SparseMatrix< Number > &&  other)

Move constructor. Create a new SparseMatrix by stealing the internal data.

Definition at line 314 of file cuda_sparse_matrix.cc.

◆ SparseMatrix() [4/4]

template<typename Number >
CUDAWrappers::SparseMatrix< Number >::SparseMatrix ( const CUDAWrappers::SparseMatrix< Number > &  )
delete

Copy constructor is deleted.

◆ ~SparseMatrix()

template<typename Number >
SparseMatrix< Number >::~SparseMatrix< Number > ( )

Destructor. Free all memory.

Definition at line 314 of file cuda_sparse_matrix.cc.

Member Function Documentation

◆ operator=() [1/2]

template<typename Number >
SparseMatrix< Number > & SparseMatrix< Number >::operator= ( CUDAWrappers::SparseMatrix< Number > &&  other)

Move assignment operator.

Definition at line 361 of file cuda_sparse_matrix.cc.

◆ operator=() [2/2]

template<typename Number >
SparseMatrix & CUDAWrappers::SparseMatrix< Number >::operator= ( const CUDAWrappers::SparseMatrix< Number > &  )
delete

Copy assignment is deleted.

◆ reinit()

template<typename Number >
void SparseMatrix< Number >::reinit ( Utilities::CUDA::Handle handle,
const ::SparseMatrix< Number > &  sparse_matrix_host 
)

Reinitialize the sparse matrix. The sparse matrix on the host is copied to the device and the elementes are reordered according to the format supported by cuSPARSE.

Definition at line 386 of file cuda_sparse_matrix.cc.

◆ m()

template<typename Number >
SparseMatrix< Number >::size_type SparseMatrix< Number >::m ( ) const
inline

Return the dimension of the codomain (or range) space. Note that the matrix is of dimension \(m \times n\).

Definition at line 383 of file cuda_sparse_matrix.h.

◆ n()

template<typename Number >
SparseMatrix< Number >::size_type SparseMatrix< Number >::n ( ) const
inline

Return the dimension of the domain space. Note that the matrix is of dimension \(m \times n\).

Definition at line 392 of file cuda_sparse_matrix.h.

◆ n_nonzero_elements()

template<typename Number >
std::size_t SparseMatrix< Number >::n_nonzero_elements ( ) const
inline

Return the number of nonzero elements of this matrix. Actually, it returns the number of entries in the sparsity pattern; if any of the entries should happen to be zero, it is counted anyway.

Definition at line 401 of file cuda_sparse_matrix.h.

◆ print()

template<typename Number >
template<typename StreamType >
void SparseMatrix< Number >::print ( StreamType &  out,
const bool  across = false,
const bool  diagonal_first = true 
) const
inline

Print the matrix to the given stream, using the format (row,column) value, i.e. one nonzero entry of the matrix per line. If across is true, print all entries on a single line, using the format row,column:value.

If the argument diagonal_first is true, diagonal elements of quadratic matrices are printed first in their row. If it is false, the elements in a row are written in ascending column order.

Definition at line 411 of file cuda_sparse_matrix.h.

◆ print_formatted()

template<typename Number >
void SparseMatrix< Number >::print_formatted ( std::ostream &  out,
const unsigned int  precision = 3,
const bool  scientific = true,
const unsigned int  width = 0,
const char *  zero_string = " ",
const double  denominator = 1.,
const char *  separator = " " 
) const

Print the matrix in the usual format, i.e., as a matrix and not as a list of nonzero elements. For better readability, elements not in the matrix are displayed as empty space, while matrix elements which are explicitly set to zero are displayed as such.

The parameters allow for a flexible setting of the output format: precision and scientific are used to determine the number format, where scientific = false means fixed point notation. A zero entry for width makes the function compute a width, but it may be changed to a positive value, if output is crude.

Additionally, a character for an empty value may be specified in zero_string, and a character to separate row entries can be set in separator.

Finally, the whole matrix can be multiplied with a common denominator to produce more readable output, even integers.

Attention
This function may produce large amounts of output if applied to a large matrix!

Definition at line 467 of file cuda_sparse_matrix.h.

◆ operator*=()

template<typename Number >
SparseMatrix< Number > & SparseMatrix< Number >::operator*= ( const Number  factor)

Multiply the entire matrix by a fixed factor.

Definition at line 480 of file cuda_sparse_matrix.cc.

◆ operator/=()

template<typename Number >
SparseMatrix< Number > & SparseMatrix< Number >::operator/= ( const Number  factor)

Divide the entire matrix by a fixed factor.

Definition at line 495 of file cuda_sparse_matrix.cc.

◆ vmult()

template<typename Number >
void SparseMatrix< Number >::vmult ( LinearAlgebra::CUDAWrappers::Vector< Number > &  dst,
const LinearAlgebra::CUDAWrappers::Vector< Number > &  src 
) const

Matrix-vector multiplication: let \(dst = M \cdot src\) with \(M\) being this matrix.

Definition at line 511 of file cuda_sparse_matrix.cc.

◆ Tvmult()

template<typename Number >
void SparseMatrix< Number >::Tvmult ( LinearAlgebra::CUDAWrappers::Vector< Number > &  dst,
const LinearAlgebra::CUDAWrappers::Vector< Number > &  src 
) const

Matrix-vector multiplication: let \(dst = M^T \cdot src\) with \(M\) being this matrix. This function does the same as vmult() but takes this transposed matrix.

Definition at line 529 of file cuda_sparse_matrix.cc.

◆ vmult_add()

template<typename Number >
void SparseMatrix< Number >::vmult_add ( LinearAlgebra::CUDAWrappers::Vector< Number > &  dst,
const LinearAlgebra::CUDAWrappers::Vector< Number > &  src 
) const

Adding matrix-vector multiplication. Add \(M \cdot src\) on \(dst\) with \(M\) being this matrix.

Definition at line 547 of file cuda_sparse_matrix.cc.

◆ Tvmult_add()

template<typename Number >
void SparseMatrix< Number >::Tvmult_add ( LinearAlgebra::CUDAWrappers::Vector< Number > &  dst,
const LinearAlgebra::CUDAWrappers::Vector< Number > &  src 
) const

Adding matrix-vector multiplication. Add \(M^T \cdot src\) to \(dst\) with \(M\) being this matrix. This function foes the same as vmult_add() but takes the transposed matrix.

Definition at line 565 of file cuda_sparse_matrix.cc.

◆ matrix_norm_square()

template<typename Number >
Number SparseMatrix< Number >::matrix_norm_square ( const LinearAlgebra::CUDAWrappers::Vector< Number > &  v) const

Return the square of the norm of the vector \(v\) with respect to the norm induced by this matrix, i.e., \(\left(v,Mv\right)\). This is useful, e.g., in the finite context, where the \(L_2\) norm of a function equals the matrix norm with respect to the mass matrix of the vector representing the nodal values of the finite element function.

Obviously, the matrix needs to be quadratic for this operation.

Definition at line 583 of file cuda_sparse_matrix.cc.

◆ matrix_scalar_product()

template<typename Number >
Number SparseMatrix< Number >::matrix_scalar_product ( const LinearAlgebra::CUDAWrappers::Vector< Number > &  u,
const LinearAlgebra::CUDAWrappers::Vector< Number > &  v 
) const

Compute the matrix scalar product \(\left(u,Mv\right)\).

Definition at line 596 of file cuda_sparse_matrix.cc.

◆ residual()

template<typename Number >
Number SparseMatrix< Number >::residual ( LinearAlgebra::CUDAWrappers::Vector< Number > &  dst,
const LinearAlgebra::CUDAWrappers::Vector< Number > &  x,
const LinearAlgebra::CUDAWrappers::Vector< Number > &  b 
) const

Compute the residual of an equation \(M \cdot x=b\), where the residual is defined to be \(r=b-M \cdot x\). Write the residual into \(dst\). The \(l_2\) norm of the residual vector is returned.

Source \(x\) and destination \(dst\) must not be the same vector.

Definition at line 610 of file cuda_sparse_matrix.cc.

◆ l1_norm()

template<typename Number >
Number SparseMatrix< Number >::l1_norm ( ) const

Return the \(l_1\)-norm of the matrix, that is \(|M|_1=\max_{\mathrm{all\ columns\ }j}\sum_{\mathrm{all\ rows\ }i} |M_{ij}|\), (max. sum of columns). This is the natural matrix norm that is compatible to the \(l_1\)-norm for vectors, i.e., \(|Mv|_1\leq |M|_1 |v|_1\).

Definition at line 625 of file cuda_sparse_matrix.cc.

◆ linfty_norm()

template<typename Number >
Number SparseMatrix< Number >::linfty_norm ( ) const

Return the \(l_\infty\)-norm of the matrix, that is \(|M|_\infty=\max_{\mathrm{all\ rows\ }i}\sum_{\mathrm{all\ columns\ }j} |M_{ij}|\), (max. sum of rows). This is the natural norm that is compatible to the \(l_\infty\)-norm of vectors, i.e., \(|Mv|_\infty \leq |M|_\infty |v|_\infty\).

Definition at line 644 of file cuda_sparse_matrix.cc.

◆ frobenius_norm()

template<typename Number >
Number SparseMatrix< Number >::frobenius_norm ( ) const

Return the frobenius norm of the matrix, i.e., the square root of the sum of squares of all entries in the matrix.

Definition at line 663 of file cuda_sparse_matrix.cc.

◆ get_cusparse_matrix()

template<typename Number >
std::tuple< Number *, int *, int *, cusparseMatDescr_t, cusparseSpMatDescr_t > SparseMatrix< Number >::get_cusparse_matrix ( ) const

Return a tuple containing the pointer to the values of matrix, the pointer to the columns indices, the pointer to the rows pointer, the cuSPARSE matrix description, and the cuSPARSE SP matrix description.

Definition at line 679 of file cuda_sparse_matrix.cc.

◆ subscribe()

void Subscriptor::subscribe ( std::atomic< bool > *const  validity,
const std::string &  identifier = "" 
) const
inherited

Subscribes a user of the object by storing the pointer validity. The subscriber may be identified by text supplied as identifier.

Definition at line 135 of file subscriptor.cc.

◆ unsubscribe()

void Subscriptor::unsubscribe ( std::atomic< bool > *const  validity,
const std::string &  identifier = "" 
) const
inherited

Unsubscribes a user from the object.

Note
The identifier and the validity pointer must be the same as the one supplied to subscribe().

Definition at line 155 of file subscriptor.cc.

◆ n_subscriptions()

unsigned int Subscriptor::n_subscriptions ( ) const
inlineinherited

Return the present number of subscriptions to this object. This allows to use this class for reference counted lifetime determination where the last one to unsubscribe also deletes the object.

Definition at line 300 of file subscriptor.h.

◆ list_subscribers() [1/2]

template<typename StreamType >
void Subscriptor::list_subscribers ( StreamType &  stream) const
inlineinherited

List the subscribers to the input stream.

Definition at line 317 of file subscriptor.h.

◆ list_subscribers() [2/2]

void Subscriptor::list_subscribers ( ) const
inherited

List the subscribers to deallog.

Definition at line 203 of file subscriptor.cc.

◆ serialize()

template<class Archive >
void Subscriptor::serialize ( Archive &  ar,
const unsigned int  version 
)
inlineinherited

Read or write the data of this object to or from a stream for the purpose of serialization using the BOOST serialization library.

This function does not actually serialize any of the member variables of this class. The reason is that what this class stores is only who subscribes to this object, but who does so at the time of storing the contents of this object does not necessarily have anything to do with who subscribes to the object when it is restored. Consequently, we do not want to overwrite the subscribers at the time of restoring, and then there is no reason to write the subscribers out in the first place.

Definition at line 309 of file subscriptor.h.

◆ check_no_subscribers()

void Subscriptor::check_no_subscribers ( ) const
privatenoexceptinherited

Check that there are no objects subscribing to this object. If this check passes then it is safe to destroy the current object. It this check fails then this function will either abort or print an error message to deallog (by using the AssertNothrow mechanism), but will not throw an exception.

Note
Since this function is just a consistency check it does nothing in release mode.
If this function is called when there is an uncaught exception then, rather than aborting, this function prints an error message to the standard error stream and returns.

Definition at line 52 of file subscriptor.cc.

Member Data Documentation

◆ cusparse_handle

template<typename Number >
cusparseHandle_t CUDAWrappers::SparseMatrix< Number >::cusparse_handle
private

cuSPARSE handle used to call cuSPARSE functions.

Definition at line 336 of file cuda_sparse_matrix.h.

◆ nnz

template<typename Number >
int CUDAWrappers::SparseMatrix< Number >::nnz
private

Number of non-zero elements in the sparse matrix.

Definition at line 341 of file cuda_sparse_matrix.h.

◆ n_rows

template<typename Number >
int CUDAWrappers::SparseMatrix< Number >::n_rows
private

Number of rows of the sparse matrix.

Definition at line 346 of file cuda_sparse_matrix.h.

◆ n_cols

template<typename Number >
int CUDAWrappers::SparseMatrix< Number >::n_cols
private

Number of columns of the sparse matrix.

Definition at line 351 of file cuda_sparse_matrix.h.

◆ val_dev

template<typename Number >
std::unique_ptr<Number[], void (*)(Number *)> CUDAWrappers::SparseMatrix< Number >::val_dev
private

Pointer to the values (on the device) of the sparse matrix.

Definition at line 356 of file cuda_sparse_matrix.h.

◆ column_index_dev

template<typename Number >
std::unique_ptr<int[], void (*)(int *)> CUDAWrappers::SparseMatrix< Number >::column_index_dev
private

Pointer to the column indices (on the device) of the sparse matrix.

Definition at line 361 of file cuda_sparse_matrix.h.

◆ row_ptr_dev

template<typename Number >
std::unique_ptr<int[], void (*)(int *)> CUDAWrappers::SparseMatrix< Number >::row_ptr_dev
private

Pointer to the row pointer (on the device) of the sparse matrix.

Definition at line 366 of file cuda_sparse_matrix.h.

◆ descr

template<typename Number >
cusparseMatDescr_t CUDAWrappers::SparseMatrix< Number >::descr
private

cuSPARSE description of the matrix.

Definition at line 371 of file cuda_sparse_matrix.h.

◆ sp_descr

template<typename Number >
cusparseSpMatDescr_t CUDAWrappers::SparseMatrix< Number >::sp_descr
private

cuSPARSE description of the sparse matrix.

Definition at line 376 of file cuda_sparse_matrix.h.

◆ counter

std::atomic<unsigned int> Subscriptor::counter
mutableprivateinherited

Store the number of objects which subscribed to this object. Initially, this number is zero, and upon destruction it shall be zero again (i.e. all objects which subscribed should have unsubscribed again).

The creator (and owner) of an object is counted in the map below if HE manages to supply identification.

We use the mutable keyword in order to allow subscription to constant objects also.

This counter may be read from and written to concurrently in multithreaded code: hence we use the std::atomic class template.

Definition at line 218 of file subscriptor.h.

◆ counter_map

std::map<std::string, unsigned int> Subscriptor::counter_map
mutableprivateinherited

In this map, we count subscriptions for each different identification string supplied to subscribe().

Definition at line 224 of file subscriptor.h.

◆ validity_pointers

std::vector<std::atomic<bool> *> Subscriptor::validity_pointers
mutableprivateinherited

In this vector, we store pointers to the validity bool in the SmartPointer objects that subscribe to this class.

Definition at line 240 of file subscriptor.h.

◆ object_info

const std::type_info* Subscriptor::object_info
mutableprivateinherited

Pointer to the typeinfo object of this object, from which we can later deduce the class name. Since this information on the derived class is neither available in the destructor, nor in the constructor, we obtain it in between and store it here.

Definition at line 248 of file subscriptor.h.

◆ mutex

std::mutex Subscriptor::mutex
staticprivateinherited

A mutex used to ensure data consistency when accessing the mutable members of this class. This lock is used in the subscribe() and unsubscribe() functions, as well as in list_subscribers().

Definition at line 271 of file subscriptor.h.


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