12#ifndef _CU_BASIC_ALGEBRA_CU_MATRIX_OPERATIONS_H_
13#define _CU_BASIC_ALGEBRA_CU_MATRIX_OPERATIONS_H_
20#include "../_definitions/types.h"
24 #pragma warning(push, 0)
25 #include <cublas_v2.h>
27#elif defined(__INTEL_LLVM_COMPILER) || defined(__INTEL_COMPILER)
28 #pragma warning(push, 0)
29 #include <cublas_v2.h>
31#elif defined(__GNUC__) || defined(__clang__)
32 #pragma GCC diagnostic push
33 #pragma GCC diagnostic ignored "-Wswitch-enum"
34 #include <cublas_v2.h>
35 #pragma GCC diagnostic pop
37 #include <cublas_v2.h>
42 #define RESTRICT __restrict
43#elif defined(__INTEL_COMPILER)
44 #define RESTRICT __restrict
45#elif defined(__CUDA__) || defined(__GNUC__) || defined(__clang__)
46 #define RESTRICT __restrict__
85template <
typename DataType>
92 cublasHandle_t cublas_handle,
102 cublasHandle_t cublas_handle,
105 const DataType alpha,
113 cublasHandle_t cublas_handle,
123 cublasHandle_t cublas_handle,
126 const DataType alpha,
134 cusparseHandle_t cusparse_handle,
144 cusparseHandle_t cusparse_handle,
149 const DataType alpha,
155 cusparseHandle_t cusparse_handle,
166 cusparseHandle_t cusparse_handle,
171 const DataType alpha,
178 cusparseHandle_t cusparse_handle,
189 cusparseHandle_t cusparse_handle,
194 const DataType alpha,
201 cusparseHandle_t cusparse_handle,
211 cusparseHandle_t cusparse_handle,
216 const DataType alpha,
222 cusparseHandle_t cublas_handle,
224 const DataType*
RESTRICT supdiagonals,
A static class for matrix-vector operations, which are similar to the level-2 operations of the BLAS ...
static void csr_transposed_matvec_plus(cusparseHandle_t cusparse_handle, const DataType *RESTRICT A_data, const LongIndexType *RESTRICT A_column_indices, const LongIndexType *RESTRICT A_index_pointer, const DataType *RESTRICT b, const DataType alpha, const LongIndexType num_rows, const LongIndexType num_columns, DataType *RESTRICT c)
Computes where is compressed sparse row (CSR) matrix and is a dense vector. The output is a dense...
static void create_band_matrix(cusparseHandle_t cublas_handle, const DataType *RESTRICT diagonals, const DataType *RESTRICT supdiagonals, const IndexType non_zero_size, const FlagType tridiagonal, DataType **matrix)
Creates bi-diagonal or symmetric tri-diagonal matrix from the diagonal array (diagonals) and off-diag...
static void csc_matvec_plus(cusparseHandle_t cusparse_handle, const DataType *RESTRICT A_data, const LongIndexType *RESTRICT A_row_indices, const LongIndexType *RESTRICT A_index_pointer, const DataType *RESTRICT b, const DataType alpha, const LongIndexType num_rows, const LongIndexType num_columns, DataType *RESTRICT c)
Computes where is compressed sparse column (CSC) matrix and is a dense vector. The output is a de...
static void dense_matvec(cublasHandle_t cublas_handle, const DataType *RESTRICT A, const DataType *RESTRICT b, const LongIndexType num_rows, const LongIndexType num_columns, const FlagType A_is_row_major, DataType *RESTRICT c)
Computes the matrix vector multiplication where is a dense matrix.
static void dense_transposed_matvec(cublasHandle_t cublas_handle, const DataType *RESTRICT A, const DataType *RESTRICT b, const LongIndexType num_rows, const LongIndexType num_columns, const FlagType A_is_row_major, DataType *RESTRICT c)
Computes matrix vector multiplication where is dense, and is the transpose of the matrix .
static void csc_transposed_matvec(cusparseHandle_t cusparse_handle, const DataType *RESTRICT A_data, const LongIndexType *RESTRICT A_row_indices, const LongIndexType *RESTRICT A_index_pointer, const DataType *RESTRICT b, const LongIndexType num_columns, DataType *RESTRICT c)
Computes where is compressed sparse column (CSC) matrix and is a dense vector. The output is a de...
static void csr_matvec(cusparseHandle_t cusparse_handle, const DataType *RESTRICT A_data, const LongIndexType *RESTRICT A_column_indices, const LongIndexType *RESTRICT A_index_pointer, const DataType *RESTRICT b, const LongIndexType num_rows, DataType *RESTRICT c)
Computes where is compressed sparse row (CSR) matrix and is a dense vector. The output is a dense...
static void csc_transposed_matvec_plus(cusparseHandle_t cusparse_handle, const DataType *RESTRICT A_data, const LongIndexType *RESTRICT A_row_indices, const LongIndexType *RESTRICT A_index_pointer, const DataType *RESTRICT b, const DataType alpha, const LongIndexType num_columns, DataType *RESTRICT c)
Computes where is compressed sparse column (CSC) matrix and is a dense vector. The output is a de...
static void csr_matvec_plus(cusparseHandle_t cusparse_handle, const DataType *RESTRICT A_data, const LongIndexType *RESTRICT A_column_indices, const LongIndexType *RESTRICT A_index_pointer, const DataType *RESTRICT b, const DataType alpha, const LongIndexType num_rows, DataType *RESTRICT c)
Computes where is compressed sparse row (CSR) matrix and is a dense vector. The output is a dense...
static void csr_transposed_matvec(cusparseHandle_t cusparse_handle, const DataType *RESTRICT A_data, const LongIndexType *RESTRICT A_column_indices, const LongIndexType *RESTRICT A_index_pointer, const DataType *RESTRICT b, const LongIndexType num_rows, const LongIndexType num_columns, DataType *RESTRICT c)
Computes where is compressed sparse row (CSR) matrix and is a dense vector. The output is a dense...
static void csc_matvec(cusparseHandle_t cusparse_handle, const DataType *RESTRICT A_data, const LongIndexType *RESTRICT A_row_indices, const LongIndexType *RESTRICT A_index_pointer, const DataType *RESTRICT b, const LongIndexType num_rows, const LongIndexType num_columns, DataType *RESTRICT c)
Computes where is compressed sparse column (CSC) matrix and is a dense vector. The output is a de...
static void dense_transposed_matvec_plus(cublasHandle_t cublas_handle, const DataType *RESTRICT A, const DataType *RESTRICT b, const DataType alpha, const LongIndexType num_rows, const LongIndexType num_columns, const FlagType A_is_row_major, DataType *RESTRICT c)
Computes where is dense, and is the transpose of the matrix .
static void dense_matvec_plus(cublasHandle_t cublas_handle, const DataType *RESTRICT A, const DataType *RESTRICT b, const DataType alpha, const LongIndexType num_rows, const LongIndexType num_columns, const FlagType A_is_row_major, DataType *RESTRICT c)
Computes the operation where is a dense matrix.