imate
C++/CUDA Reference
cuCSRMatrix< DataType > Class Template Reference

#include <cu_csr_matrix.h>

Inheritance diagram for cuCSRMatrix< DataType >:
Collaboration diagram for cuCSRMatrix< DataType >:

Public Member Functions

 cuCSRMatrix ()
 Default constructor. More...
 
 cuCSRMatrix (const DataType *A_data_, const LongIndexType *A_indices_, const LongIndexType *A_index_pointer_, const LongIndexType num_rows_, const LongIndexType num_columns_, const int num_gpu_devices_)
 Constructor with arguments. More...
 
virtual ~cuCSRMatrix ()
 Virtual desructor. More...
 
virtual void dot (const DataType *device_vector, DataType *device_product)
 
virtual void dot_plus (const DataType *device_vector, const DataType alpha, DataType *device_product)
 
virtual void transpose_dot (const DataType *device_vector, DataType *device_product)
 
virtual void transpose_dot_plus (const DataType *device_vector, const DataType alpha, DataType *device_product)
 
- Public Member Functions inherited from cuMatrix< DataType >
 cuMatrix ()
 
 cuMatrix (int num_gpu_devices_)
 
virtual ~cuMatrix ()
 
- Public Member Functions inherited from cuLinearOperator< DataType >
 cuLinearOperator ()
 
 cuLinearOperator (int num_gpu_devices_)
 Constructor with setting num_rows and num_columns. More...
 
virtual ~cuLinearOperator ()
 
cublasHandle_t get_cublas_handle () const
 This function returns a reference to the cublasHandle_t object. The object will be created, if it is not created already. More...
 
- Public Member Functions inherited from cLinearOperator< DataType >
 cLinearOperator ()
 Default constructor. More...
 
 cLinearOperator (const LongIndexType num_rows_, const LongIndexType num_columns_)
 Constructor with setting num_rows and num_columns. More...
 
virtual ~cLinearOperator ()
 
LongIndexType get_num_rows () const
 
LongIndexType get_num_columns () const
 
void set_parameters (DataType *parameters_)
 Sets the scalar parameter this->parameters. Parameter is initialized to NULL. However, before calling dot or transpose_dot functions, the parameters must be set. More...
 
IndexType get_num_parameters () const
 
FlagType is_eigenvalue_relation_known () const
 Returns a flag that determines whether a relation between the parameters of the operator and its eigenvalue(s) is known. More...
 
- Public Member Functions inherited from cCSRMatrix< DataType >
 cCSRMatrix ()
 
 cCSRMatrix (const DataType *A_data_, const LongIndexType *A_indices_, const LongIndexType *A_index_pointer_, const LongIndexType num_rows_, const LongIndexType num_columns_)
 
virtual ~cCSRMatrix ()
 
virtual FlagType is_identity_matrix () const
 Checks whether the matrix is identity. More...
 
LongIndexType get_nnz () const
 Returns the number of non-zero elements of the sparse matrix. More...
 
- Public Member Functions inherited from cMatrix< DataType >
 cMatrix ()
 Default constructor. More...
 
virtual ~cMatrix ()
 
DataType get_eigenvalue (const DataType *known_parameters, const DataType known_eigenvalue, const DataType *inquiry_parameters) const
 This virtual function is implemented from its pure virtual function of the base class. In this class, this functio has no use and was only implemented so that this class be able to be instantiated (due to the pure virtual function). More...
 

Protected Member Functions

virtual void copy_host_to_device ()
 Copies the member data from the host memory to the device memory. More...
 
void allocate_buffer (const int device_id, cusparseOperation_t cusparse_operation, const DataType alpha, const DataType beta, cusparseDnVecDescr_t &cusparse_input_vector, cusparseDnVecDescr_t &cusparse_output_vector, cusparseSpMVAlg_t algorithm)
 Allocates an external buffer for matrix-vector multiplication using cusparseSpMV function. More...
 
- Protected Member Functions inherited from cuLinearOperator< DataType >
int query_gpu_devices () const
 Before any numerical computation, this method chechs if any gpu device is available on the machine, or notifies the user if nothing was found. More...
 
void initialize_cublas_handle ()
 Creates a cublasHandle_t object, if not created already. More...
 
void initialize_cusparse_handle ()
 Creates a cusparseHandle_t object, if not created already. More...
 

Protected Attributes

DataType ** device_A_data
 
LongIndexType ** device_A_indices
 
LongIndexType ** device_A_index_pointer
 
void ** device_buffer
 
size_t * device_buffer_num_bytes
 
cusparseSpMatDescr_t * cusparse_matrix_A
 
- Protected Attributes inherited from cuLinearOperator< DataType >
int num_gpu_devices
 
bool copied_host_to_device
 
cublasHandle_t * cublas_handle
 
cusparseHandle_t * cusparse_handle
 
- Protected Attributes inherited from cLinearOperator< DataType >
const LongIndexType num_rows
 
const LongIndexType num_columns
 
FlagType eigenvalue_relation_known
 
DataType * parameters
 
IndexType num_parameters
 
- Protected Attributes inherited from cCSRMatrix< DataType >
const DataType * A_data
 
const LongIndexTypeA_indices
 
const LongIndexTypeA_index_pointer
 

Detailed Description

template<typename DataType>
class cuCSRMatrix< DataType >

Definition at line 30 of file cu_csr_matrix.h.

Constructor & Destructor Documentation

◆ cuCSRMatrix() [1/2]

template<typename DataType >
cuCSRMatrix< DataType >::cuCSRMatrix

Default constructor.

Definition at line 33 of file cu_csr_matrix.cu.

33  :
34  device_A_data(NULL),
35  device_A_indices(NULL),
37  device_buffer(NULL),
39  cusparse_matrix_A(NULL)
40 {
41 }
size_t * device_buffer_num_bytes
Definition: cu_csr_matrix.h:86
cusparseSpMatDescr_t * cusparse_matrix_A
Definition: cu_csr_matrix.h:87
LongIndexType ** device_A_index_pointer
Definition: cu_csr_matrix.h:84
DataType ** device_A_data
Definition: cu_csr_matrix.h:82
LongIndexType ** device_A_indices
Definition: cu_csr_matrix.h:83
void ** device_buffer
Definition: cu_csr_matrix.h:85

◆ cuCSRMatrix() [2/2]

template<typename DataType >
cuCSRMatrix< DataType >::cuCSRMatrix ( const DataType *  A_data_,
const LongIndexType A_indices_,
const LongIndexType A_index_pointer_,
const LongIndexType  num_rows_,
const LongIndexType  num_columns_,
const int  num_gpu_devices_ 
)

Constructor with arguments.

Definition at line 52 of file cu_csr_matrix.cu.

58  :
59 
60  // Base class constructor
61  cLinearOperator<DataType>(num_rows_, num_columns_),
62  cCSRMatrix<DataType>(A_data_, A_indices_, A_index_pointer_, num_rows_,
63  num_columns_),
64  cuMatrix<DataType>(num_gpu_devices_),
65 
66  // Initializer list
67  device_A_data(NULL),
68  device_A_indices(NULL),
70  device_buffer(NULL),
71  cusparse_matrix_A(NULL)
72 {
74  this->copy_host_to_device();
75 
76  // Initialize device buffer
77  this->device_buffer = new void*[this->num_gpu_devices];
78  this->device_buffer_num_bytes = new size_t[this->num_gpu_devices];
79  for (int device_id=0; device_id < this->num_gpu_devices; ++device_id)
80  {
81  this->device_buffer[device_id] = NULL;
82  this->device_buffer_num_bytes[device_id] = 0;
83  }
84 }
Base class for linear operators. This class serves as interface for all derived classes.
virtual void copy_host_to_device()
Copies the member data from the host memory to the device memory.
void initialize_cusparse_handle()
Creates a cusparseHandle_t object, if not created already.
Base class for constant matrices.
Definition: cu_matrix.h:41

References cuCSRMatrix< DataType >::copy_host_to_device(), cuCSRMatrix< DataType >::device_buffer, cuCSRMatrix< DataType >::device_buffer_num_bytes, cuLinearOperator< DataType >::initialize_cusparse_handle(), and cuLinearOperator< DataType >::num_gpu_devices.

Here is the call graph for this function:

◆ ~cuCSRMatrix()

template<typename DataType >
cuCSRMatrix< DataType >::~cuCSRMatrix
virtual

Virtual desructor.

Definition at line 95 of file cu_csr_matrix.cu.

96 {
97  // Member objects exist if the second constructor was called.
98  if (this->copied_host_to_device)
99  {
100  // Deallocate arrays of data on gpu
101  for (int device_id=0; device_id < this->num_gpu_devices; ++device_id)
102  {
103  // Switch to a device
105 
106  // Deallocate
109  this->device_A_indices[device_id]);
111  this->device_A_index_pointer[device_id]);
114  this->cusparse_matrix_A[device_id]);
115  }
116  }
117 
118  // Deallocate arrays of pointers on cpu
119  if (this->device_A_data != NULL)
120  {
121  delete[] this->device_A_data;
122  this->device_A_data = NULL;
123  }
124 
125  if (this->device_A_indices != NULL)
126  {
127  delete[] this->device_A_indices;
128  this->device_A_indices = NULL;
129  }
130 
131  if (this->device_A_index_pointer != NULL)
132  {
133  delete[] this->device_A_index_pointer;
134  this->device_A_index_pointer = NULL;
135  }
136 
137  if (this->device_buffer != NULL)
138  {
139  delete[] this->device_buffer;
140  this->device_buffer = NULL;
141  }
142 
143  if (this->device_buffer_num_bytes != NULL)
144  {
145  delete[] this->device_buffer_num_bytes;
146  this->device_buffer_num_bytes = NULL;
147  }
148 
149  if (this->cusparse_matrix_A != NULL)
150  {
151  delete[] this->cusparse_matrix_A;
152  this->cusparse_matrix_A = NULL;
153  }
154 }
static void del(void *device_array)
Deletes memory on gpu device if its pointer is not NULL, then sets the pointer to NULL.
static void set_device(int device_id)
Sets the current device in multi-gpu applications.
void destroy_cusparse_matrix(cusparseSpMatDescr_t &cusparse_matrix)
Destroys cusparse matrix.

References CudaInterface< ArrayType >::del(), cusparse_interface::destroy_cusparse_matrix(), and CudaInterface< ArrayType >::set_device().

Here is the call graph for this function:

Member Function Documentation

◆ allocate_buffer()

template<typename DataType >
void cuCSRMatrix< DataType >::allocate_buffer ( const int  device_id,
cusparseOperation_t  cusparse_operation,
const DataType  alpha,
const DataType  beta,
cusparseDnVecDescr_t &  cusparse_input_vector,
cusparseDnVecDescr_t &  cusparse_output_vector,
cusparseSpMVAlg_t  algorithm 
)
protected

Allocates an external buffer for matrix-vector multiplication using cusparseSpMV function.

If buffer size if not the same as required buffer size, allocate (or reallocate) memory. The allocation is always performed in the first call of this function since buffer size is initialized to zero in constructor. But for the next calls it might not be reallocated if the buffer size is the same.

Definition at line 241 of file cu_csr_matrix.cu.

249 {
250  // Find the buffer size needed for matrix-vector multiplication
251  size_t required_buffer_size;
253  this->cusparse_handle[device_id], cusparse_operation, alpha,
254  this->cusparse_matrix_A[device_id], cusparse_input_vector, beta,
255  cusparse_output_vector, algorithm, &required_buffer_size);
256 
257  if (this->device_buffer_num_bytes[device_id] != required_buffer_size)
258  {
259  // Update the buffer size
260  this->device_buffer_num_bytes[device_id] = required_buffer_size;
261 
262  // Delete buffer if it was allocated previously
264 
265  // Allocate (or reallocate) buffer on device.
267  this->device_buffer[device_id],
268  this->device_buffer_num_bytes[device_id]);
269  }
270 }
static void alloc_bytes(void *&device_array, const size_t num_bytes)
Allocates memory on gpu device. This function uses an existing given pointer.
cusparseHandle_t * cusparse_handle
void cusparse_matrix_buffer_size(cusparseHandle_t cusparse_handle, cusparseOperation_t cusparse_operation, const DataType alpha, cusparseSpMatDescr_t cusparse_matrix, cusparseDnVecDescr_t cusparse_input_vector, const DataType beta, cusparseDnVecDescr_t cusparse_output_vector, cusparseSpMVAlg_t algorithm, size_t *buffer_size)

References CudaInterface< ArrayType >::alloc_bytes(), cusparse_interface::cusparse_matrix_buffer_size(), and CudaInterface< ArrayType >::del().

Here is the call graph for this function:

◆ copy_host_to_device()

template<typename DataType >
void cuCSRMatrix< DataType >::copy_host_to_device
protectedvirtual

Copies the member data from the host memory to the device memory.

Implements cuMatrix< DataType >.

Definition at line 165 of file cu_csr_matrix.cu.

166 {
167  if (!this->copied_host_to_device)
168  {
169  // Set the number of threads
170  omp_set_num_threads(this->num_gpu_devices);
171 
172  // Array sizes
173  LongIndexType A_data_size = this->get_nnz();
174  LongIndexType A_indices_size = A_data_size;
175  LongIndexType A_index_pointer_size = this->num_rows + 1;
176  LongIndexType A_nnz = this->get_nnz();
177 
178  // Create array of pointers for data on each gpu device
179  this->device_A_data = new DataType*[this->num_gpu_devices];
180  this->device_A_indices = new LongIndexType*[this->num_gpu_devices];
181  this->device_A_index_pointer = \
182  new LongIndexType*[this->num_gpu_devices];
183  this->cusparse_matrix_A = \
184  new cusparseSpMatDescr_t[this->num_gpu_devices];
185 
186  #pragma omp parallel
187  {
188  // Switch to a device with the same device id as the cpu thread id
189  unsigned int thread_id = omp_get_thread_num();
191 
192  // A_data
194  A_data_size);
196  this->A_data, A_data_size, this->device_A_data[thread_id]);
197 
198  // A_indices
200  this->device_A_indices[thread_id], A_indices_size);
202  this->A_indices, A_indices_size,
203  this->device_A_indices[thread_id]);
204 
205  // A_index_pointer
207  this->device_A_index_pointer[thread_id],
208  A_index_pointer_size);
210  this->A_index_pointer, A_index_pointer_size,
211  this->device_A_index_pointer[thread_id]);
212 
213  // Create cusparse matrix
215  this->cusparse_matrix_A[thread_id], this->num_rows,
216  this->num_columns, A_nnz, this->device_A_data[thread_id],
217  this->device_A_indices[thread_id],
218  this->device_A_index_pointer[thread_id]);
219  }
220 
221  // Flag to prevent reinitialization
222  this->copied_host_to_device = true;
223  }
224 }
static ArrayType * alloc(const LongIndexType array_size)
Allocates memory on gpu device. This function creates a pointer and returns it.
static void copy_to_device(const ArrayType *host_array, const LongIndexType array_size, ArrayType *device_array)
Copies memory on host to device memory.
const LongIndexType * A_index_pointer
Definition: c_csr_matrix.h:72
LongIndexType get_nnz() const
Returns the number of non-zero elements of the sparse matrix.
const DataType * A_data
Definition: c_csr_matrix.h:70
const LongIndexType * A_indices
Definition: c_csr_matrix.h:71
const LongIndexType num_rows
const LongIndexType num_columns
void create_cusparse_matrix(cusparseSpMatDescr_t &cusparse_matrix, const LongIndexType num_rows, const LongIndexType num_columns, const LongIndexType nnz, DataType *device_A_data, LongIndexType *device_A_indices, LongIndexType *device_A_index_pointer)
int LongIndexType
Definition: types.h:60

References CudaInterface< ArrayType >::alloc(), CudaInterface< ArrayType >::copy_to_device(), cusparse_interface::create_cusparse_matrix(), and CudaInterface< ArrayType >::set_device().

Referenced by cuCSRMatrix< DataType >::cuCSRMatrix().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ dot()

template<typename DataType >
void cuCSRMatrix< DataType >::dot ( const DataType *  device_vector,
DataType *  device_product 
)
virtual

Reimplemented from cCSRMatrix< DataType >.

Definition at line 278 of file cu_csr_matrix.cu.

281 {
282  assert(this->copied_host_to_device);
283 
284  // Create cusparse vector for the input vector
285  cusparseDnVecDescr_t cusparse_input_vector;
287  cusparse_input_vector, this->num_columns,
288  const_cast<DataType*>(device_vector));
289 
290  // Create cusparse vector for the output vector
291  cusparseDnVecDescr_t cusparse_output_vector;
293  cusparse_output_vector, this->num_rows, device_product);
294 
295  // Matrix vector settings
296  DataType alpha = 1.0;
297  DataType beta = 0.0;
298  cusparseOperation_t cusparse_operation = CUSPARSE_OPERATION_NON_TRANSPOSE;
299  cusparseSpMVAlg_t algorithm = CUSPARSE_SPMV_ALG_DEFAULT;
300 
301  // Get device id
302  int device_id = CudaInterface<DataType>::get_device();
303 
304  // Allocate device buffer (or reallocation if needed)
305  this->allocate_buffer(device_id, cusparse_operation, alpha, beta,
306  cusparse_input_vector, cusparse_output_vector,
307  algorithm);
308 
309  // Matrix vector multiplication
311  this->cusparse_handle[device_id], cusparse_operation, alpha,
312  this->cusparse_matrix_A[device_id], cusparse_input_vector, beta,
313  cusparse_output_vector, algorithm, this->device_buffer[device_id]);
314 
315  // Destroy cusparse vectors
316  cusparse_interface::destroy_cusparse_vector(cusparse_input_vector);
317  cusparse_interface::destroy_cusparse_vector(cusparse_output_vector);
318 }
static int get_device()
Gets the current device in multi-gpu applications.
void allocate_buffer(const int device_id, cusparseOperation_t cusparse_operation, const DataType alpha, const DataType beta, cusparseDnVecDescr_t &cusparse_input_vector, cusparseDnVecDescr_t &cusparse_output_vector, cusparseSpMVAlg_t algorithm)
Allocates an external buffer for matrix-vector multiplication using cusparseSpMV function.
#define CUSPARSE_SPMV_ALG_DEFAULT
void destroy_cusparse_vector(cusparseDnVecDescr_t &cusparse_vector)
Destroys cusparse vector.
void cusparse_matvec(cusparseHandle_t cusparse_handle, cusparseOperation_t cusparse_operation, const DataType alpha, cusparseSpMatDescr_t cusparse_matrix, cusparseDnVecDescr_t cusparse_input_vector, const DataType beta, cusparseDnVecDescr_t cusparse_output_vector, cusparseSpMVAlg_t algorithm, void *external_buffer)
void create_cusparse_vector(cusparseDnVecDescr_t &cusparse_vector, const LongIndexType vector_size, DataType *device_vector)

References cusparse_interface::create_cusparse_vector(), cusparse_interface::cusparse_matvec(), CUSPARSE_SPMV_ALG_DEFAULT, cusparse_interface::destroy_cusparse_vector(), and CudaInterface< ArrayType >::get_device().

Here is the call graph for this function:

◆ dot_plus()

template<typename DataType >
void cuCSRMatrix< DataType >::dot_plus ( const DataType *  device_vector,
const DataType  alpha,
DataType *  device_product 
)
virtual

Reimplemented from cCSRMatrix< DataType >.

Definition at line 326 of file cu_csr_matrix.cu.

330 {
331  assert(this->copied_host_to_device);
332 
333  // Create cusparse vector for the input vector
334  cusparseDnVecDescr_t cusparse_input_vector;
336  cusparse_input_vector, this->num_columns,
337  const_cast<DataType*>(device_vector));
338 
339  // Create cusparse vector for the output vector
340  cusparseDnVecDescr_t cusparse_output_vector;
342  cusparse_output_vector, this->num_rows, device_product);
343 
344  // Matrix vector settings
345  DataType beta = 1.0;
346  cusparseOperation_t cusparse_operation = CUSPARSE_OPERATION_NON_TRANSPOSE;
347  cusparseSpMVAlg_t algorithm = CUSPARSE_SPMV_ALG_DEFAULT;
348 
349  // Get device id
350  int device_id = CudaInterface<DataType>::get_device();
351 
352  // Allocate device buffer (or reallocation if needed)
353  this->allocate_buffer(device_id, cusparse_operation, alpha, beta,
354  cusparse_input_vector, cusparse_output_vector,
355  algorithm);
356 
357  // Matrix vector multiplication
359  this->cusparse_handle[device_id], cusparse_operation, alpha,
360  this->cusparse_matrix_A[device_id], cusparse_input_vector, beta,
361  cusparse_output_vector, algorithm, this->device_buffer[device_id]);
362 
363  // Destroy cusparse vectors
364  cusparse_interface::destroy_cusparse_vector(cusparse_input_vector);
365  cusparse_interface::destroy_cusparse_vector(cusparse_output_vector);
366 }

References cusparse_interface::create_cusparse_vector(), cusparse_interface::cusparse_matvec(), CUSPARSE_SPMV_ALG_DEFAULT, cusparse_interface::destroy_cusparse_vector(), and CudaInterface< ArrayType >::get_device().

Here is the call graph for this function:

◆ transpose_dot()

template<typename DataType >
void cuCSRMatrix< DataType >::transpose_dot ( const DataType *  device_vector,
DataType *  device_product 
)
virtual

Reimplemented from cCSRMatrix< DataType >.

Definition at line 374 of file cu_csr_matrix.cu.

377 {
378  assert(this->copied_host_to_device);
379 
380  // Create cusparse vector for the input vector
381  cusparseDnVecDescr_t cusparse_input_vector;
383  cusparse_input_vector, this->num_columns,
384  const_cast<DataType*>(device_vector));
385 
386  // Create cusparse vector for the output vector
387  cusparseDnVecDescr_t cusparse_output_vector;
389  cusparse_output_vector, this->num_rows, device_product);
390 
391  // Matrix vector settings
392  DataType alpha = 1.0;
393  DataType beta = 0.0;
394  cusparseOperation_t cusparse_operation = CUSPARSE_OPERATION_TRANSPOSE;
395  cusparseSpMVAlg_t algorithm = CUSPARSE_SPMV_ALG_DEFAULT;
396 
397  // Get device id
398  int device_id = CudaInterface<DataType>::get_device();
399 
400  // Allocate device buffer (or reallocation if needed)
401  this->allocate_buffer(device_id, cusparse_operation, alpha, beta,
402  cusparse_input_vector, cusparse_output_vector,
403  algorithm);
404 
405  // Matrix vector multiplication
407  this->cusparse_handle[device_id], cusparse_operation, alpha,
408  this->cusparse_matrix_A[device_id], cusparse_input_vector, beta,
409  cusparse_output_vector, algorithm, this->device_buffer[device_id]);
410 
411  // Destroy cusparse vectors
412  cusparse_interface::destroy_cusparse_vector(cusparse_input_vector);
413  cusparse_interface::destroy_cusparse_vector(cusparse_output_vector);
414 }

References cusparse_interface::create_cusparse_vector(), cusparse_interface::cusparse_matvec(), CUSPARSE_SPMV_ALG_DEFAULT, cusparse_interface::destroy_cusparse_vector(), and CudaInterface< ArrayType >::get_device().

Here is the call graph for this function:

◆ transpose_dot_plus()

template<typename DataType >
void cuCSRMatrix< DataType >::transpose_dot_plus ( const DataType *  device_vector,
const DataType  alpha,
DataType *  device_product 
)
virtual

Reimplemented from cCSRMatrix< DataType >.

Definition at line 422 of file cu_csr_matrix.cu.

426 {
427  assert(this->copied_host_to_device);
428 
429  // Create cusparse vector for the input vector
430  cusparseDnVecDescr_t cusparse_input_vector;
432  cusparse_input_vector, this->num_columns,
433  const_cast<DataType*>(device_vector));
434 
435  // Create cusparse vector for the output vector
436  cusparseDnVecDescr_t cusparse_output_vector;
438  cusparse_output_vector, this->num_rows, device_product);
439 
440  // Matrix vector settings
441  DataType beta = 1.0;
442  cusparseOperation_t cusparse_operation = CUSPARSE_OPERATION_TRANSPOSE;
443  cusparseSpMVAlg_t algorithm = CUSPARSE_SPMV_ALG_DEFAULT;
444 
445  // Get device id
446  int device_id = CudaInterface<DataType>::get_device();
447 
448  // Allocate device buffer (or reallocation if needed)
449  this->allocate_buffer(device_id, cusparse_operation, alpha, beta,
450  cusparse_input_vector, cusparse_output_vector,
451  algorithm);
452 
453  // Matrix vector multiplication
455  this->cusparse_handle[device_id], cusparse_operation, alpha,
456  this->cusparse_matrix_A[device_id], cusparse_input_vector, beta,
457  cusparse_output_vector, algorithm, this->device_buffer[device_id]);
458 
459  // Destroy cusparse vectors
460  cusparse_interface::destroy_cusparse_vector(cusparse_input_vector);
461  cusparse_interface::destroy_cusparse_vector(cusparse_output_vector);
462 }

References cusparse_interface::create_cusparse_vector(), cusparse_interface::cusparse_matvec(), CUSPARSE_SPMV_ALG_DEFAULT, cusparse_interface::destroy_cusparse_vector(), and CudaInterface< ArrayType >::get_device().

Here is the call graph for this function:

Member Data Documentation

◆ cusparse_matrix_A

template<typename DataType >
cusparseSpMatDescr_t* cuCSRMatrix< DataType >::cusparse_matrix_A
protected

Definition at line 87 of file cu_csr_matrix.h.

◆ device_A_data

template<typename DataType >
DataType** cuCSRMatrix< DataType >::device_A_data
protected

Definition at line 82 of file cu_csr_matrix.h.

◆ device_A_index_pointer

template<typename DataType >
LongIndexType** cuCSRMatrix< DataType >::device_A_index_pointer
protected

Definition at line 84 of file cu_csr_matrix.h.

◆ device_A_indices

template<typename DataType >
LongIndexType** cuCSRMatrix< DataType >::device_A_indices
protected

Definition at line 83 of file cu_csr_matrix.h.

◆ device_buffer

template<typename DataType >
void** cuCSRMatrix< DataType >::device_buffer
protected

Definition at line 85 of file cu_csr_matrix.h.

Referenced by cuCSRMatrix< DataType >::cuCSRMatrix().

◆ device_buffer_num_bytes

template<typename DataType >
size_t* cuCSRMatrix< DataType >::device_buffer_num_bytes
protected

Definition at line 86 of file cu_csr_matrix.h.

Referenced by cuCSRMatrix< DataType >::cuCSRMatrix().


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