Reference documentation for deal.II version GIT relicensing-437-g81ec864850 2024-04-19 07:30: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
Classes | Public Types | Public Member Functions | Private Attributes | List of all members
CUDAWrappers::PreconditionIC< Number > Class Template Reference

#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
 
PreconditionICoperator= (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 intP_row_ptr_dev
 
const intP_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
 

Detailed Description

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

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).

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

Definition at line 61 of file cuda_precondition.h.

Member Typedef Documentation

◆ size_type

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

Declare the type for container size.

Definition at line 67 of file cuda_precondition.h.

Constructor & Destructor Documentation

◆ PreconditionIC() [1/2]

template<typename Number >
CUDAWrappers::PreconditionIC< Number >::PreconditionIC ( const Utilities::CUDA::Handle handle)

Constructor.

Definition at line 1199 of file cuda_precondition.cc.

◆ PreconditionIC() [2/2]

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

The copy constructor is deleted.

◆ ~PreconditionIC()

template<typename Number >
CUDAWrappers::PreconditionIC< Number >::~PreconditionIC ( )

Destructor. Free all resources that were initialized in this class.

Definition at line 1248 of file cuda_precondition.cc.

Member Function Documentation

◆ operator=()

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

The copy assignment operator is deleted.

◆ initialize()

template<typename Number >
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.

Definition at line 1271 of file cuda_precondition.cc.

◆ vmult()

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

Apply the preconditioner.

Definition at line 1430 of file cuda_precondition.cc.

◆ Tvmult()

template<typename Number >
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().

Definition at line 1484 of file cuda_precondition.cc.

◆ m()

template<typename Number >
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\).

Note
This function should only be called if the preconditioner has been initialized.

◆ n()

template<typename Number >
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\).

Note
This function should only be called if the preconditioner has been initialized.

Member Data Documentation

◆ cusparse_handle

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

cuSPARSE handle used to call cuSPARSE functions.

Definition at line 165 of file cuda_precondition.h.

◆ descr_M

template<typename Number >
cusparseMatDescr_t CUDAWrappers::PreconditionIC< Number >::descr_M
private

cuSPARSE description of the sparse matrix \(M=LL^T\).

Definition at line 170 of file cuda_precondition.h.

◆ descr_L

template<typename Number >
cusparseMatDescr_t CUDAWrappers::PreconditionIC< Number >::descr_L
private

cuSPARSE description of the lower triangular matrix \(L\).

Definition at line 175 of file cuda_precondition.h.

◆ info_M

template<typename Number >
csric02Info_t CUDAWrappers::PreconditionIC< Number >::info_M
private

Solve and analysis structure for \(M=LL^T\).

Definition at line 180 of file cuda_precondition.h.

◆ info_L

template<typename Number >
csrsv2Info_t CUDAWrappers::PreconditionIC< Number >::info_L
private

Solve and analysis structure for the lower triangular matrix \(L\).

Definition at line 185 of file cuda_precondition.h.

◆ info_Lt

template<typename Number >
csrsv2Info_t CUDAWrappers::PreconditionIC< Number >::info_Lt
private

Solve and analysis structure for the upper triangular matrix \(L^T\).

Definition at line 190 of file cuda_precondition.h.

◆ matrix_pointer

template<typename Number >
SmartPointer<const SparseMatrix<Number> > CUDAWrappers::PreconditionIC< Number >::matrix_pointer
private

Pointer to the matrix this object was initialized with.

Definition at line 195 of file cuda_precondition.h.

◆ P_val_dev

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

Pointer to the values (on the device) of the computed preconditioning matrix.

Definition at line 201 of file cuda_precondition.h.

◆ P_row_ptr_dev

template<typename Number >
const int* CUDAWrappers::PreconditionIC< Number >::P_row_ptr_dev
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 207 of file cuda_precondition.h.

◆ P_column_index_dev

template<typename Number >
const int* CUDAWrappers::PreconditionIC< Number >::P_column_index_dev
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 213 of file cuda_precondition.h.

◆ tmp_dev

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

Pointer to the value (on the device) for a temporary (helper) vector used in vmult().

Definition at line 219 of file cuda_precondition.h.

◆ buffer_dev

template<typename Number >
std::unique_ptr<void, void (*)(void *)> CUDAWrappers::PreconditionIC< Number >::buffer_dev
private

Pointer to an internal buffer (on the device) that is used for computing the decomposition.

Definition at line 225 of file cuda_precondition.h.

◆ policy_L

template<typename Number >
cusparseSolvePolicy_t CUDAWrappers::PreconditionIC< Number >::policy_L
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 232 of file cuda_precondition.h.

◆ policy_Lt

template<typename Number >
cusparseSolvePolicy_t CUDAWrappers::PreconditionIC< Number >::policy_Lt
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 239 of file cuda_precondition.h.

◆ policy_M

template<typename Number >
cusparseSolvePolicy_t CUDAWrappers::PreconditionIC< Number >::policy_M
private

Determine if level information should be generated for \(M=LL^T\). This value can be modified through an AdditionalData object.

Definition at line 245 of file cuda_precondition.h.

◆ n_rows

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

The number of rows is the same as for the matrix this object has been initialized with.

Definition at line 251 of file cuda_precondition.h.

◆ n_nonzero_elements

template<typename Number >
int CUDAWrappers::PreconditionIC< Number >::n_nonzero_elements
private

The number of non-zero elements is the same as for the matrix this object has been initialized with.

Definition at line 257 of file cuda_precondition.h.


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