COMPASS  5.4.4
End-to-end AO simulation tool using GPU acceleration
carma_sparse_obj.h
Go to the documentation of this file.
1 // -----------------------------------------------------------------------------
2 // This file is part of COMPASS <https://anr-compass.github.io/compass/>
3 //
4 // Copyright (C) 2011-2023 COMPASS Team <https://github.com/ANR-COMPASS>
5 // All rights reserved.
6 
7 // -----------------------------------------------------------------------------
8 
16 
17 
18 #ifndef CARMA_SPARSE_OBJ_H_
19 #define CARMA_SPARSE_OBJ_H_
20 #include <cusparse_v2.h>
21 #include "carma_obj.h"
22 
23 template <class T_data>
24 class CarmaSparseHostObj;
25 
26 #ifndef USE_MAGMA_SPARSE
27 typedef void *magma_d_sparse_matrix;
28 typedef void *magma_s_sparse_matrix;
29 #else
30 #include "magmasparse.h"
31 #endif
32 
33 template <class T_data>
35  public:
36  long dims_data[3];
37  int nz_elem;
38  int device;
40 
41  // ZERO-BASED INDEXING CSR-FORMAT
42  T_data *d_data; // nz_elem elements
43  int *d_rowind; // dim1+1 elements
44  int *d_colind; // nz_elem elements
45  cusparseMatDescr_t descr;
46 #if CUDA_VERSION >= 11000
47  cusparseSpMatDescr_t sp_descr; // New cuSPARSE API from CUDA 11
48  cusparseDnMatDescr_t dn_descr; // New cuSPARSE API from CUDA 11
49 #endif
50 
51  char major_dim;
52  std::string format;
53  int block_dim; // blockdim for BSR format
54 
55  // Magma stuff
56  union {
59  };
60 
61  private:
62  public:
68  CarmaSparseObj(CarmaContext *current_context, const long *dims, T_data *M,
69  bool load_from_host);
71  T_data *values, int *colind, int *rowind, int nz,
72  bool load_from_host);
73  virtual ~CarmaSparseObj();
74 
77 
78  void resize(int nnz_, int dim1_, int dim2_);
80  char get_major_dim() const { return major_dim; }
81  void set_majorDim(char c) { major_dim = c; }
82 
84  operator T_data *() { return d_data; }
85  T_data *operator[](int index) { return &d_data[index]; }
86  T_data *get_data() { return d_data; }
87  T_data *get_data(int index) { return &d_data[index]; }
88  const long *get_dims() { return dims_data; }
89  long get_dims(int i) { return dims_data[i]; }
90  int get_nonzero_elem() { return nz_elem; }
92 
93  int get_device() { return device; }
94 
95  void sparse_to_host(int *h_rowInd, int *h_colInd, T_data *h_data);
96 
97 #if CUDA_VERSION >= 11000
98  void transpose();
99  cudaDataType_t get_data_type() {
100  cudaDataType_t data_type;
101  if (std::is_same<T_data, float>::value)
102  data_type = CUDA_R_32F;
103  else if (std::is_same<T_data, double>::value)
104  data_type = CUDA_R_64F;
105  else
106  std::cerr << "Unsupported data type" << std::endl;
107  return data_type;
108  }
109 #endif
110 
111  void allocate(int nnz, int dim1, int dim2);
112 
113  private:
114  void _create(int nnz_, int dim1_, int dim2_);
115  void _clear();
116 
117 #if CUDA_VERSION < 12000
118  template <cusparseStatus_t CUSPARSEAPI (*ptr_nnz)(
119  cusparseHandle_t handle, cusparseDirection_t dirA, int m, int n,
120  const cusparseMatDescr_t descrA, const T_data *A, int lda,
121  int *nnzPerRowCol, int *nnzTotalDevHostPtr),
122  cusparseStatus_t CUSPARSEAPI (*ptr_dense2csr)(
123  cusparseHandle_t handle, int m, int n,
124  const cusparseMatDescr_t descrA, const T_data *A, int lda,
125  const int *nnzPerRow, T_data *csrValA, int *csrRowPtrA,
126  int *csrColIndA)>
127  void init_carma_sparse_obj(CarmaContext *current_context, const long *dims,
128  T_data *M, bool load_from_host);
129 # else
130  void init_carma_sparse_obj(CarmaContext *current_context, const long *dims,
131  T_data *M, bool load_from_host);
132 
133 #endif
134 };
135 
136 template <class T_data>
137 cusparseStatus_t carma_gemv(cusparseHandle_t handle, char op_A, T_data alpha,
138  CarmaSparseObj<T_data> *A, T_data *x, T_data beta,
139  T_data *y);
140 
141 template <class T_data>
142 cusparseStatus_t carma_gemm(cusparseHandle_t handle, char op_A, T_data alpha,
144  T_data beta, CarmaObj<T_data> *C);
145 
146 template <class T_data>
147 cusparseStatus_t carma_gemm(cusparseHandle_t handle, char op_A, char op_B,
151 
152 template <class T_data>
153 cusparseStatus_t carma_csr2dense(CarmaSparseObj<T_data> *src, T_data *dest);
154 
155 template <class T_data>
156 cusparseStatus_t carma_csr2bsr(CarmaSparseObj<T_data> *src, int block_dim,
157  CarmaSparseObj<T_data> *dest);
158 
159 template <class T_data>
160 cusparseStatus_t carma_bsr2csr(CarmaSparseObj<T_data> *src,
161  CarmaSparseObj<T_data> *dest);
162 
163 template <class T_data>
164 int carma_kgemv(CarmaSparseObj<T_data> *A, T_data alpha,
165  const T_data *__restrict x, T_data beta, T_data *y);
166 #endif /* CARMA_SPARSE_OBJ_H_ */
cusparseStatus_t carma_csr2bsr(CarmaSparseObj< T_data > *src, int block_dim, CarmaSparseObj< T_data > *dest)
cusparseStatus_t carma_csr2dense(CarmaSparseObj< T_data > *src, T_data *dest)
cusparseStatus_t carma_bsr2csr(CarmaSparseObj< T_data > *src, CarmaSparseObj< T_data > *dest)
void * magma_s_sparse_matrix
int carma_kgemv(CarmaSparseObj< T_data > *A, T_data alpha, const T_data *__restrict x, T_data beta, T_data *y)
void * magma_d_sparse_matrix
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_gemv(cusparseHandle_t handle, char op_A, T_data alpha, CarmaSparseObj< T_data > *A, T_data *x, T_data beta, T_data *y)
this class provides the context in which CarmaObj are created
Definition: carma_context.h:79
this class provides wrappers to the generic carma object
Definition: carma_obj.h:92
this class provides wrappers to the generic carma sparse host object
this class provides wrappers to the generic carma sparse object
cusparseMatDescr_t descr
CarmaSparseObj(CarmaObj< T_data > *M)
char get_major_dim() const
long get_dims(int i)
T_data * operator[](int index)
CarmaSparseObj(CarmaContext *current_context, CarmaSparseHostObj< T_data > *M)
void allocate(int nnz, int dim1, int dim2)
virtual ~CarmaSparseObj()
int device
device where the CarmaObj is allocate
void operator=(CarmaSparseObj< T_data > &M)
CarmaSparseObj(CarmaContext *current_context)
const long * get_dims()
long dims_data[3]
dimensions of the array
int nz_elem
number of elements in the array
magma_d_sparse_matrix d_sparse_mat
void set_majorDim(char c)
magma_s_sparse_matrix s_sparse_mat
void resize(int nnz_, int dim1_, int dim2_)
CarmaContext * current_context
CarmaSparseObj(CarmaSparseObj< T_data > *M)
void sparse_to_host(int *h_rowInd, int *h_colInd, T_data *h_data)
CarmaContext * get_context()
void operator=(CarmaSparseHostObj< T_data > &M)
T_data * get_data()
T_data * get_data(int index)
bool is_column_major()
CarmaSparseObj(CarmaContext *current_context, const long *dims, T_data *values, int *colind, int *rowind, int nz, bool load_from_host)
std::string format
CarmaSparseObj(CarmaContext *current_context, const long *dims, T_data *M, bool load_from_host)