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

#include <cu_csc_matrix.h>

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

Public Member Functions

 cuCSCMatrix ()
 Default constructor. More...
 
 cuCSCMatrix (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 ~cuCSCMatrix ()
 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 cCSCMatrix< DataType >
 cCSCMatrix ()
 
 cCSCMatrix (const DataType *A_data_, const LongIndexType *A_indices_, const LongIndexType *A_index_pointer_, const LongIndexType num_rows_, const LongIndexType num_columns_)
 
virtual ~cCSCMatrix ()
 
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 cCSCMatrix< DataType >
const DataType * A_data
 
const LongIndexTypeA_indices
 
const LongIndexTypeA_index_pointer
 

Detailed Description

template<typename DataType>
class cuCSCMatrix< DataType >

Definition at line 30 of file cu_csc_matrix.h.

Constructor & Destructor Documentation

◆ cuCSCMatrix() [1/2]

template<typename DataType >
cuCSCMatrix< DataType >::cuCSCMatrix

Default constructor.

Definition at line 33 of file cu_csc_matrix.cu.

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

◆ cuCSCMatrix() [2/2]

template<typename DataType >
cuCSCMatrix< DataType >::cuCSCMatrix ( 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_csc_matrix.cu.

58  :
59 
60  // Base class constructor
61  cLinearOperator<DataType>(num_rows_, num_columns_),
62  cCSCMatrix<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 cuCSCMatrix< DataType >::copy_host_to_device(), cuCSCMatrix< DataType >::device_buffer, cuCSCMatrix< DataType >::device_buffer_num_bytes, cuLinearOperator< DataType >::initialize_cusparse_handle(), and cuLinearOperator< DataType >::num_gpu_devices.

Here is the call graph for this function:

◆ ~cuCSCMatrix()

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

Virtual desructor.

Definition at line 95 of file cu_csc_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 cuCSCMatrix< 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 250 of file cu_csc_matrix.cu.

258 {
259  // Find the buffer size needed for matrix-vector multiplication
260  size_t required_buffer_size;
262  this->cusparse_handle[device_id], cusparse_operation, alpha,
263  this->cusparse_matrix_A[device_id], cusparse_input_vector, beta,
264  cusparse_output_vector, algorithm, &required_buffer_size);
265 
266  if (this->device_buffer_num_bytes[device_id] != required_buffer_size)
267  {
268  // Update the buffer size
269  this->device_buffer_num_bytes[device_id] = required_buffer_size;
270 
271  // Delete buffer if it was allocated previously
273 
274  // Allocate (or reallocate) buffer on device.
276  this->device_buffer[device_id],
277  this->device_buffer_num_bytes[device_id]);
278  }
279 }
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 cuCSCMatrix< DataType >::copy_host_to_device
protectedvirtual

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

Note
Despite the input matrix is a CSC matrix, we treat it as a CSR matrix, since cusparse's interface is only for CSR matrices. For this, we swap the number of columns and rows from the input matrix to the cusparse matrix.

Implements cuMatrix< DataType >.

Definition at line 169 of file cu_csc_matrix.cu.

170 {
171  if (!this->copied_host_to_device)
172  {
173  // Set the number of threads
174  omp_set_num_threads(this->num_gpu_devices);
175 
176  // Array sizes
177  LongIndexType A_data_size = this->get_nnz();
178  LongIndexType A_indices_size = A_data_size;
179  LongIndexType A_index_pointer_size = this->num_rows + 1;
180  LongIndexType A_nnz = this->get_nnz();
181 
182  // Swapping the number of rows and columns to treat the input CSC
183  // matrix as a CSR matrix.
184  LongIndexType csc_num_rows = this->num_columns;
185  LongIndexType csc_num_columns = this->num_rows;
186 
187  // Create array of pointers for data on each gpu device
188  this->device_A_data = new DataType*[this->num_gpu_devices];
189  this->device_A_indices = new LongIndexType*[this->num_gpu_devices];
190  this->device_A_index_pointer = \
191  new LongIndexType*[this->num_gpu_devices];
192  this->cusparse_matrix_A = \
193  new cusparseSpMatDescr_t[this->num_gpu_devices];
194 
195  #pragma omp parallel
196  {
197  // Switch to a device with the same device id as the cpu thread id
198  unsigned int thread_id = omp_get_thread_num();
200 
201  // A_data
203  A_data_size);
205  this->A_data, A_data_size, this->device_A_data[thread_id]);
206 
207  // A_indices
209  this->device_A_indices[thread_id], A_indices_size);
211  this->A_indices, A_indices_size,
212  this->device_A_indices[thread_id]);
213 
214  // A_index_pointer
216  this->device_A_index_pointer[thread_id],
217  A_index_pointer_size);
219  this->A_index_pointer, A_index_pointer_size,
220  this->device_A_index_pointer[thread_id]);
221 
222  // Create cusparse matrix
224  this->cusparse_matrix_A[thread_id], csc_num_rows,
225  csc_num_columns, A_nnz, this->device_A_data[thread_id],
226  this->device_A_indices[thread_id],
227  this->device_A_index_pointer[thread_id]);
228  }
229 
230  // Flag to prevent reinitialization
231  this->copied_host_to_device = true;
232  }
233 }
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.
LongIndexType get_nnz() const
Returns the number of non-zero elements of the sparse matrix.
const DataType * A_data
Definition: c_csc_matrix.h:70
const LongIndexType * A_indices
Definition: c_csc_matrix.h:71
const LongIndexType * A_index_pointer
Definition: c_csc_matrix.h:72
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 cuCSCMatrix< DataType >::cuCSCMatrix().

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

◆ dot()

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

Reimplemented from cCSCMatrix< DataType >.

Definition at line 287 of file cu_csc_matrix.cu.

290 {
291  assert(this->copied_host_to_device);
292 
293  // Create cusparse vector for the input vector
294  cusparseDnVecDescr_t cusparse_input_vector;
296  cusparse_input_vector, this->num_columns,
297  const_cast<DataType*>(device_vector));
298 
299  // Create cusparse vector for the output vector
300  cusparseDnVecDescr_t cusparse_output_vector;
302  cusparse_output_vector, this->num_rows, device_product);
303 
304  // Matrix vector settings
305  DataType alpha = 1.0;
306  DataType beta = 0.0;
307 
308  // Using transpose operation since we treat CSC matrix as CSR
309  cusparseOperation_t cusparse_operation = CUSPARSE_OPERATION_TRANSPOSE;
310  cusparseSpMVAlg_t algorithm = CUSPARSE_SPMV_ALG_DEFAULT;
311 
312  // Get device id
313  int device_id = CudaInterface<DataType>::get_device();
314 
315  // Allocate device buffer (or reallocation if needed)
316  this->allocate_buffer(device_id, cusparse_operation, alpha, beta,
317  cusparse_input_vector, cusparse_output_vector,
318  algorithm);
319 
320  // Matrix vector multiplication
322  this->cusparse_handle[device_id], cusparse_operation, alpha,
323  this->cusparse_matrix_A[device_id], cusparse_input_vector, beta,
324  cusparse_output_vector, algorithm, this->device_buffer[device_id]);
325 
326  // Destroy cusparse vectors
327  cusparse_interface::destroy_cusparse_vector(cusparse_input_vector);
328  cusparse_interface::destroy_cusparse_vector(cusparse_output_vector);
329 }
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 cuCSCMatrix< DataType >::dot_plus ( const DataType *  device_vector,
const DataType  alpha,
DataType *  device_product 
)
virtual

Reimplemented from cCSCMatrix< DataType >.

Definition at line 337 of file cu_csc_matrix.cu.

341 {
342  assert(this->copied_host_to_device);
343 
344  // Create cusparse vector for the input vector
345  cusparseDnVecDescr_t cusparse_input_vector;
347  cusparse_input_vector, this->num_columns,
348  const_cast<DataType*>(device_vector));
349 
350  // Create cusparse vector for the output vector
351  cusparseDnVecDescr_t cusparse_output_vector;
353  cusparse_output_vector, this->num_rows, device_product);
354 
355  // Matrix vector settings
356  DataType beta = 1.0;
357 
358  // Using transpose operation since we treat CSC matrix as CSR
359  cusparseOperation_t cusparse_operation = CUSPARSE_OPERATION_TRANSPOSE;
360  cusparseSpMVAlg_t algorithm = CUSPARSE_SPMV_ALG_DEFAULT;
361 
362  // Get device id
363  int device_id = CudaInterface<DataType>::get_device();
364 
365  // Allocate device buffer (or reallocation if needed)
366  this->allocate_buffer(device_id, cusparse_operation, alpha, beta,
367  cusparse_input_vector, cusparse_output_vector,
368  algorithm);
369 
370  // Matrix vector multiplication
372  this->cusparse_handle[device_id], cusparse_operation, alpha,
373  this->cusparse_matrix_A[device_id], cusparse_input_vector, beta,
374  cusparse_output_vector, algorithm, this->device_buffer[device_id]);
375 
376  // Destroy cusparse vectors
377  cusparse_interface::destroy_cusparse_vector(cusparse_input_vector);
378  cusparse_interface::destroy_cusparse_vector(cusparse_output_vector);
379 }

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 cuCSCMatrix< DataType >::transpose_dot ( const DataType *  device_vector,
DataType *  device_product 
)
virtual

Reimplemented from cCSCMatrix< DataType >.

Definition at line 387 of file cu_csc_matrix.cu.

390 {
391  assert(this->copied_host_to_device);
392 
393  // Create cusparse vector for the input vector
394  cusparseDnVecDescr_t cusparse_input_vector;
396  cusparse_input_vector, this->num_columns,
397  const_cast<DataType*>(device_vector));
398 
399  // Create cusparse vector for the output vector
400  cusparseDnVecDescr_t cusparse_output_vector;
402  cusparse_output_vector, this->num_rows, device_product);
403 
404  // Matrix vector settings
405  DataType alpha = 1.0;
406  DataType beta = 0.0;
407 
408  // Using non-transpose operation since we treat CSC matrix as CSR
409  cusparseOperation_t cusparse_operation = CUSPARSE_OPERATION_NON_TRANSPOSE;
410  cusparseSpMVAlg_t algorithm = CUSPARSE_SPMV_ALG_DEFAULT;
411 
412  // Get device id
413  int device_id = CudaInterface<DataType>::get_device();
414 
415  // Allocate device buffer (or reallocation if needed)
416  this->allocate_buffer(device_id, cusparse_operation, alpha, beta,
417  cusparse_input_vector, cusparse_output_vector,
418  algorithm);
419 
420  // Matrix vector multiplication
422  this->cusparse_handle[device_id], cusparse_operation, alpha,
423  this->cusparse_matrix_A[device_id], cusparse_input_vector, beta,
424  cusparse_output_vector, algorithm, this->device_buffer[device_id]);
425 
426  // Destroy cusparse vectors
427  cusparse_interface::destroy_cusparse_vector(cusparse_input_vector);
428  cusparse_interface::destroy_cusparse_vector(cusparse_output_vector);
429 }

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 cuCSCMatrix< DataType >::transpose_dot_plus ( const DataType *  device_vector,
const DataType  alpha,
DataType *  device_product 
)
virtual

Reimplemented from cCSCMatrix< DataType >.

Definition at line 437 of file cu_csc_matrix.cu.

441 {
442  assert(this->copied_host_to_device);
443 
444  // Create cusparse vector for the input vector
445  cusparseDnVecDescr_t cusparse_input_vector;
447  cusparse_input_vector, this->num_columns,
448  const_cast<DataType*>(device_vector));
449 
450  // Create cusparse vector for the output vector
451  cusparseDnVecDescr_t cusparse_output_vector;
453  cusparse_output_vector, this->num_rows, device_product);
454 
455  // Matrix vector settings
456  DataType beta = 1.0;
457 
458  // Using non-transpose operation since we treat CSC matrix as CSR
459  cusparseOperation_t cusparse_operation = CUSPARSE_OPERATION_NON_TRANSPOSE;
460  cusparseSpMVAlg_t algorithm = CUSPARSE_SPMV_ALG_DEFAULT;
461 
462  // Get device id
463  int device_id = CudaInterface<DataType>::get_device();
464 
465  // Allocate device buffer (or reallocation if needed)
466  this->allocate_buffer(device_id, cusparse_operation, alpha, beta,
467  cusparse_input_vector, cusparse_output_vector,
468  algorithm);
469 
470  // Matrix vector multiplication
472  this->cusparse_handle[device_id], cusparse_operation, alpha,
473  this->cusparse_matrix_A[device_id], cusparse_input_vector, beta,
474  cusparse_output_vector, algorithm, this->device_buffer[device_id]);
475 
476  // Destroy cusparse vectors
477  cusparse_interface::destroy_cusparse_vector(cusparse_input_vector);
478  cusparse_interface::destroy_cusparse_vector(cusparse_output_vector);
479 }

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* cuCSCMatrix< DataType >::cusparse_matrix_A
protected

Definition at line 87 of file cu_csc_matrix.h.

◆ device_A_data

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

Definition at line 82 of file cu_csc_matrix.h.

◆ device_A_index_pointer

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

Definition at line 84 of file cu_csc_matrix.h.

◆ device_A_indices

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

Definition at line 83 of file cu_csc_matrix.h.

◆ device_buffer

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

Definition at line 85 of file cu_csc_matrix.h.

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

◆ device_buffer_num_bytes

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

Definition at line 86 of file cu_csc_matrix.h.

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


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