 |
COMPASS
5.0.0
End-to-end AO simulation tool using GPU acceleration
|
Go to the documentation of this file.
43 #ifndef CARMA_SPARSE_OBJ_H_
44 #define CARMA_SPARSE_OBJ_H_
45 #include <cusparse_v2.h>
48 template <
class T_data>
51 #ifndef USE_MAGMA_SPARSE
55 #include "magmasparse.h"
58 template <
class T_data>
92 T_data *values,
int *colind,
int *rowind,
int nz,
99 void resize(
int nnz_,
int dim1_,
int dim2_);
120 void _create(
int nnz_,
int dim1_,
int dim2_);
123 template <cusparseStatus_t CUSPARSEAPI (*ptr_nnz)(
124 cusparseHandle_t handle, cusparseDirection_t dirA,
int m,
int n,
125 const cusparseMatDescr_t descrA,
const T_data *A,
int lda,
126 int *nnzPerRowCol,
int *nnzTotalDevHostPtr),
127 cusparseStatus_t CUSPARSEAPI (*ptr_dense2csr)(
128 cusparseHandle_t handle,
int m,
int n,
129 const cusparseMatDescr_t descrA,
const T_data *A,
int lda,
130 const int *nnzPerRow, T_data *csrValA,
int *csrRowPtrA,
133 T_data *M,
bool load_from_host);
136 template <
class T_data>
137 cusparseStatus_t
carma_gemv(cusparseHandle_t handle,
char op_A, T_data alpha,
141 template <
class T_data>
142 cusparseStatus_t
carma_gemm(cusparseHandle_t handle,
char op_A, T_data alpha,
146 template <
class T_data>
147 cusparseStatus_t
carma_gemm(cusparseHandle_t handle,
char op_A,
char op_B,
152 template <
class T_data>
155 template <
class T_data>
159 template <
class T_data>
163 template <
class T_data>
165 const T_data *__restrict x, T_data beta, T_data *y);
CarmaContext * get_context()
magma_d_sparse_matrix d_sparse_mat
CarmaSparseObj(CarmaContext *current_context)
void operator=(CarmaSparseHostObj< T_data > &M)
void * magma_s_sparse_matrix
CarmaSparseObj(CarmaContext *current_context, const long *dims, T_data *values, int *colind, int *rowind, int nz, bool load_from_host)
void resize(int nnz_, int dim1_, int dim2_)
void * magma_d_sparse_matrix
cusparseStatus_t carma_csr2bsr(CarmaSparseObj< T_data > *src, int block_dim, CarmaSparseObj< T_data > *dest)
CarmaContext * current_context
long dims_data[3]
dimensions of the array
virtual ~CarmaSparseObj()
CarmaSparseObj(CarmaContext *current_context, CarmaSparseHostObj< T_data > *M)
void init_from_transpose(CarmaSparseObj< T_data > *M)
void set_majorDim(char c)
this class provides wrappers to the generic carma object
char get_major_dim() const
this class provides the context in which CarmaObj are created
this class provides wrappers to the generic carma sparse object
int nz_elem
number of elements in the array
CarmaSparseObj(CarmaObj< T_data > *M)
int carma_kgemv(CarmaSparseObj< T_data > *A, T_data alpha, const T_data *__restrict x, T_data beta, T_data *y)
cusparseStatus_t carma_gemm(cusparseHandle_t handle, char op_A, T_data alpha, CarmaSparseObj< T_data > *A, CarmaObj< T_data > *B, T_data beta, CarmaObj< T_data > *C)
cusparseStatus_t carma_csr2dense(CarmaSparseObj< T_data > *src, T_data *dest)
CarmaSparseObj(CarmaSparseObj< T_data > *M)
cusparseStatus_t carma_bsr2csr(CarmaSparseObj< T_data > *src, CarmaSparseObj< T_data > *dest)
void sparse_to_host(int *h_rowInd, int *h_colInd, T_data *h_data)
cusparseStatus_t carma_gemv(cusparseHandle_t handle, char op_A, T_data alpha, CarmaSparseObj< T_data > *A, T_data *x, T_data beta, T_data *y)
void operator=(CarmaSparseObj< T_data > &M)
T_data * operator[](int index)
this class provides wrappers to the generic carma sparse host object
int device
device where the CarmaObj is allocate
T_data * get_data(int index)
magma_s_sparse_matrix s_sparse_mat
CarmaSparseObj(CarmaContext *current_context, const long *dims, T_data *M, bool load_from_host)