Reference documentation for deal.II version 9.0.0
|
#include <deal.II/lac/cuda_sparse_matrix.h>
Public Types | |
typedef unsigned int | size_type |
typedef Number | value_type |
typedef Number | real_type |
Public Member Functions | |
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 () | |
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 |
Modifying entries | |
SparseMatrix & | operator*= (const Number factor) |
SparseMatrix & | operator/= (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 |
Public Member Functions inherited from Subscriptor | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
Subscriptor (Subscriptor &&) noexcept | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) noexcept |
void | subscribe (const char *identifier=nullptr) const |
void | unsubscribe (const char *identifier=nullptr) const |
unsigned int | n_subscriptions () const |
void | list_subscribers () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Access to underlying CUDA data | |
cusparseHandle_t | cusparse_handle |
int | nnz |
int | n_rows |
int | n_cols |
Number * | val_dev |
int * | column_index_dev |
int * | row_ptr_dev |
cusparseMatDescr_t | descr |
std::tuple< Number *, int *, int *, cusparseMatDescr_t > | get_cusparse_matrix () const |
Additional Inherited Members | |
Static Public Member Functions inherited from Subscriptor | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
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.
<float> and <double>
.Definition at line 45 of file cuda_sparse_matrix.h.
typedef unsigned int CUDAWrappers::SparseMatrix< Number >::size_type |
Declare type for container size.
Definition at line 51 of file cuda_sparse_matrix.h.
typedef Number CUDAWrappers::SparseMatrix< Number >::value_type |
Type of the matrix entries.
Definition at line 56 of file cuda_sparse_matrix.h.
typedef Number CUDAWrappers::SparseMatrix< Number >::real_type |
Declare a type that holds real-valued numbers with the same precision as the template argument to this class.
Definition at line 62 of file cuda_sparse_matrix.h.
CUDAWrappers::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.
CUDAWrappers::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.
CUDAWrappers::SparseMatrix< Number >::SparseMatrix | ( | CUDAWrappers::SparseMatrix< Number > && | ) |
Move constructor. Create a new SparseMatrix by stealing the internal data.
|
delete |
Copy constructor is deleted.
CUDAWrappers::SparseMatrix< Number >::~SparseMatrix | ( | ) |
Destructor. Free all memory.
void CUDAWrappers::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.
|
inline |
Return the dimension of the codomain (or range) space. Note that the matrix is of dimension \(m \times n\).
Definition at line 300 of file cuda_sparse_matrix.h.
|
inline |
Return the dimension of the domain space. Note that the matrix is of dimension \(m \times n\).
Definition at line 309 of file cuda_sparse_matrix.h.
|
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 318 of file cuda_sparse_matrix.h.
SparseMatrix& CUDAWrappers::SparseMatrix< Number >::operator*= | ( | const Number | factor | ) |
Multiply the entire matrix by a fixed factor.
SparseMatrix& CUDAWrappers::SparseMatrix< Number >::operator/= | ( | const Number | factor | ) |
Divide the entrie matrix by a fixed factor.
void CUDAWrappers::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.
void CUDAWrappers::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.
void CUDAWrappers::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.
void CUDAWrappers::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.
Number CUDAWrappers::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.
Number CUDAWrappers::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)\).
Number CUDAWrappers::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.
Number CUDAWrappers::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\).
Number CUDAWrappers::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\).
Number CUDAWrappers::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.
std::tuple<Number *, int *, int *, cusparseMatDescr_t> CUDAWrappers::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, and the cuSPARSE matrix description.
|
mutableprivate |
cuSPARSE handle used to call cuSPARSE functions. The cuSPARSE handle needs to be mutable to be called in a const function.
Definition at line 258 of file cuda_sparse_matrix.h.
|
private |
Number of non-zero elements in the sparse matrix.
Definition at line 263 of file cuda_sparse_matrix.h.
|
private |
Number of rows of the sparse matrix.
Definition at line 268 of file cuda_sparse_matrix.h.
|
private |
Number of columns of the sparse matrix.
Definition at line 273 of file cuda_sparse_matrix.h.
|
private |
Pointer to the values (on the device) of the sparse matrix.
Definition at line 278 of file cuda_sparse_matrix.h.
|
private |
Pointer to the column indices (on the device) of the sparse matrix.
Definition at line 283 of file cuda_sparse_matrix.h.
|
private |
Pointer to the row pointer (on the device) of the sparse matrix.
Definition at line 288 of file cuda_sparse_matrix.h.
|
private |
cuSPARSE description of the sparse matrix.
Definition at line 293 of file cuda_sparse_matrix.h.