Reference documentation for deal.II version 9.2.0
|
#include <deal.II/lac/cuda_precondition.h>
Classes | |
struct | AdditionalData |
Public Types | |
using | size_type = int |
Public Member Functions | |
PreconditionIC (const Utilities::CUDA::Handle &handle) | |
PreconditionIC (const PreconditionIC< Number > &)=delete | |
PreconditionIC & | operator= (const PreconditionIC< Number > &)=delete |
~PreconditionIC () | |
void | initialize (const SparseMatrix< Number > &matrix, const AdditionalData &additional_data=AdditionalData()) |
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 |
size_type | m () const |
size_type | n () const |
Private Attributes | |
cusparseHandle_t | cusparse_handle |
cusparseMatDescr_t | descr_M |
cusparseMatDescr_t | descr_L |
csric02Info_t | info_M |
csrsv2Info_t | info_L |
csrsv2Info_t | info_Lt |
SmartPointer< const SparseMatrix< Number > > | matrix_pointer |
std::unique_ptr< Number[], void(*)(Number *)> | P_val_dev |
const int * | P_row_ptr_dev |
const int * | P_column_index_dev |
std::unique_ptr< Number[], void(*)(Number *)> | tmp_dev |
std::unique_ptr< void, void(*)(void *)> | buffer_dev |
cusparseSolvePolicy_t | policy_L |
cusparseSolvePolicy_t | policy_Lt |
cusparseSolvePolicy_t | policy_M |
int | n_rows |
int | n_nonzero_elements |
This class implements an incomplete Cholesky factorization (IC) preconditioner for symmetric CUDAWrappers::SparseMatrix matrices.
The implementation closely follows the one documented in the cuSPARSE documentation (https://docs.nvidia.com/cuda/cusparse/index.html#cusparse-lt-t-gt-csric02).
<float> and <double>
.Definition at line 64 of file cuda_precondition.h.
using CUDAWrappers::PreconditionIC< Number >::size_type = int |
Declare the type for container size.
Definition at line 70 of file cuda_precondition.h.
CUDAWrappers::PreconditionIC< Number >::PreconditionIC | ( | const Utilities::CUDA::Handle & | handle | ) |
Constructor.
|
delete |
The copy constructor is deleted.
CUDAWrappers::PreconditionIC< Number >::~PreconditionIC | ( | ) |
Destructor. Free all resources that were initialized in this class.
|
delete |
The copy assignment operator is deleted.
void CUDAWrappers::PreconditionIC< Number >::initialize | ( | const SparseMatrix< Number > & | matrix, |
const AdditionalData & | additional_data = AdditionalData() |
||
) |
Initialize this object. In particular, the given matrix is copied to be modified in-place. For the underlying sparsity pattern pointers are stored. Specifically, this means that the current object can only be used reliably as long as matrix
is valid and has not been changed since calling this function.
The additional_data
determines if level information are used.
void CUDAWrappers::PreconditionIC< Number >::vmult | ( | LinearAlgebra::CUDAWrappers::Vector< Number > & | dst, |
const LinearAlgebra::CUDAWrappers::Vector< Number > & | src | ||
) | const |
Apply the preconditioner.
void CUDAWrappers::PreconditionIC< Number >::Tvmult | ( | LinearAlgebra::CUDAWrappers::Vector< Number > & | dst, |
const LinearAlgebra::CUDAWrappers::Vector< Number > & | src | ||
) | const |
Apply the preconditioner. Since the preconditioner is symmetric, this is the same as vmult().
size_type CUDAWrappers::PreconditionIC< Number >::m | ( | ) | const |
Return the dimension of the codomain (or range) space. Note that the matrix is square and has dimension \(m \times m\).
size_type CUDAWrappers::PreconditionIC< Number >::n | ( | ) | const |
Return the dimension of the codomain (or range) space. Note that the matrix is square and has dimension \(n \times n\).
|
private |
cuSPARSE handle used to call cuSPARSE functions.
Definition at line 168 of file cuda_precondition.h.
|
private |
cuSPARSE description of the sparse matrix \(M=LL^T\).
Definition at line 173 of file cuda_precondition.h.
|
private |
cuSPARSE description of the lower triangular matrix \(L\).
Definition at line 178 of file cuda_precondition.h.
|
private |
Solve and analysis structure for \(M=LL^T\).
Definition at line 183 of file cuda_precondition.h.
|
private |
Solve and analysis structure for the lower triangular matrix \(L\).
Definition at line 188 of file cuda_precondition.h.
|
private |
Solve and analysis structure for the upper triangular matrix \(L^T\).
Definition at line 193 of file cuda_precondition.h.
|
private |
Pointer to the matrix this object was initialized with.
Definition at line 198 of file cuda_precondition.h.
|
private |
Pointer to the values (on the device) of the computed preconditioning matrix.
Definition at line 204 of file cuda_precondition.h.
|
private |
Pointer to the row pointer (on the device) of the sparse matrix this object was initialized with. Guarded by matrix_pointer.
Definition at line 210 of file cuda_precondition.h.
|
private |
Pointer to the column indices (on the device) of the sparse matrix this object was initialized with. Guarded by matrix_pointer.
Definition at line 216 of file cuda_precondition.h.
|
private |
Pointer to the value (on the device) for a temporary (helper) vector used in vmult().
Definition at line 222 of file cuda_precondition.h.
|
private |
Pointer to an internal buffer (on the device) that is used for computing the decomposition.
Definition at line 228 of file cuda_precondition.h.
|
private |
Determine if level information should be generated for the lower triangular matrix \(L\). This value can be modified through an AdditionalData object.
Definition at line 235 of file cuda_precondition.h.
|
private |
Determine if level information should be generated for the upper triangular matrix \(L^T\). This value can be modified through an AdditionalData object.
Definition at line 242 of file cuda_precondition.h.
|
private |
Determine if level information should be generated for \(M=LL^T\). This value can be modified through an AdditionalData object.
Definition at line 248 of file cuda_precondition.h.
|
private |
The number of rows is the same as for the matrix this object has been initialized with.
Definition at line 254 of file cuda_precondition.h.
|
private |
The number of non-zero elements is the same as for the matrix this object has been initialized with.
Definition at line 260 of file cuda_precondition.h.