1 //////////////////////////////////////////////////////////////////////
2 // This file is distributed under the University of Illinois/NCSA Open Source
3 // License. See LICENSE file in top directory for details.
4 //
5 // Copyright (c) 2016 Jeongnim Kim and QMCPACK developers.
6 //
7 // File developed by:
8 // Lawrence Livermore National Laboratory
9 //
10 // File created by:
11 // Miguel A. Morales, moralessilva2@llnl.gov
12 // Lawrence Livermore National Laboratory
13 ////////////////////////////////////////////////////////////////////////////////
14
15 #ifndef CUSPARSE_FUNCTIONDEFS_DEPRECATED_H
16 #define CUSPARSE_FUNCTIONDEFS_DEPRECATED_H
17
18 #include <cassert>
19 #include <cuda_runtime.h>
20 #include "cusparse.h"
21 #include "AFQMC/Memory/CUDA/cuda_utilities.h"
22
23 #if CUSPARSE_VER_MAJOR > 10
24 #error
25 #endif
26
27 namespace cusparse
28 {
29 using qmc_cuda::cusparseOperation;
30
31 // Level-2
cusparse_csrmv(cusparseHandle_t handle,char Atrans,int m,int n,int nnz,const double alpha,const cusparseMatDescr_t & descrA,const double * csrValA,const int * csrRowPtrA,const int * csrColIndA,const double * x,const double beta,double * y)32 inline cusparseStatus_t cusparse_csrmv(cusparseHandle_t handle,
33 char Atrans,
34 int m,
35 int n,
36 int nnz,
37 const double alpha,
38 const cusparseMatDescr_t& descrA,
39 const double* csrValA,
40 const int* csrRowPtrA,
41 const int* csrColIndA,
42 const double* x,
43 const double beta,
44 double* y)
45
46 {
47 cusparseStatus_t sucess = cusparseDcsrmv(handle, cusparseOperation(Atrans), m, n, nnz, &alpha, descrA, csrValA,
48 csrRowPtrA, csrColIndA, x, &beta, y);
49 cudaDeviceSynchronize();
50 return sucess;
51 }
52
cusparse_csrmv(cusparseHandle_t handle,char Atrans,int m,int n,int nnz,const float alpha,const cusparseMatDescr_t & descrA,const float * csrValA,const int * csrRowPtrA,const int * csrColIndA,const float * x,const float beta,float * y)53 inline cusparseStatus_t cusparse_csrmv(cusparseHandle_t handle,
54 char Atrans,
55 int m,
56 int n,
57 int nnz,
58 const float alpha,
59 const cusparseMatDescr_t& descrA,
60 const float* csrValA,
61 const int* csrRowPtrA,
62 const int* csrColIndA,
63 const float* x,
64 const float beta,
65 float* y)
66
67 {
68 cusparseStatus_t sucess = cusparseScsrmv(handle, cusparseOperation(Atrans), m, n, nnz, &alpha, descrA, csrValA,
69 csrRowPtrA, csrColIndA, x, &beta, y);
70 cudaDeviceSynchronize();
71 return sucess;
72 }
73
cusparse_csrmv(cusparseHandle_t handle,char Atrans,int m,int n,int nnz,const std::complex<double> alpha,const cusparseMatDescr_t & descrA,const std::complex<double> * csrValA,const int * csrRowPtrA,const int * csrColIndA,const std::complex<double> * x,const std::complex<double> beta,std::complex<double> * y)74 inline cusparseStatus_t cusparse_csrmv(cusparseHandle_t handle,
75 char Atrans,
76 int m,
77 int n,
78 int nnz,
79 const std::complex<double> alpha,
80 const cusparseMatDescr_t& descrA,
81 const std::complex<double>* csrValA,
82 const int* csrRowPtrA,
83 const int* csrColIndA,
84 const std::complex<double>* x,
85 const std::complex<double> beta,
86 std::complex<double>* y)
87
88 {
89 cusparseStatus_t sucess =
90 cusparseZcsrmv(handle, cusparseOperation(Atrans), m, n, nnz, reinterpret_cast<cuDoubleComplex const*>(&alpha),
91 descrA, reinterpret_cast<cuDoubleComplex const*>(csrValA), csrRowPtrA, csrColIndA,
92 reinterpret_cast<cuDoubleComplex const*>(x), reinterpret_cast<cuDoubleComplex const*>(&beta),
93 reinterpret_cast<cuDoubleComplex*>(y));
94 cudaDeviceSynchronize();
95 return sucess;
96 }
97
98
cusparse_csrmv(cusparseHandle_t handle,char Atrans,int m,int n,int nnz,const std::complex<float> alpha,const cusparseMatDescr_t & descrA,const std::complex<float> * csrValA,const int * csrRowPtrA,const int * csrColIndA,const std::complex<float> * x,const std::complex<float> beta,std::complex<float> * y)99 inline cusparseStatus_t cusparse_csrmv(cusparseHandle_t handle,
100 char Atrans,
101 int m,
102 int n,
103 int nnz,
104 const std::complex<float> alpha,
105 const cusparseMatDescr_t& descrA,
106 const std::complex<float>* csrValA,
107 const int* csrRowPtrA,
108 const int* csrColIndA,
109 const std::complex<float>* x,
110 const std::complex<float> beta,
111 std::complex<float>* y)
112 {
113 cusparseStatus_t sucess =
114 cusparseCcsrmv(handle, cusparseOperation(Atrans), m, n, nnz, reinterpret_cast<cuComplex const*>(&alpha), descrA,
115 reinterpret_cast<cuComplex const*>(csrValA), csrRowPtrA, csrColIndA,
116 reinterpret_cast<cuComplex const*>(x), reinterpret_cast<cuComplex const*>(&beta),
117 reinterpret_cast<cuComplex*>(y));
118 cudaDeviceSynchronize();
119 return sucess;
120 }
121
122
cusparse_csrmm(cusparseHandle_t handle,char Atrans,int m,int n,int k,int nnz,const double alpha,const cusparseMatDescr_t & descrA,const double * csrValA,const int * csrRowPtrA,const int * csrColIndA,const double * B,const int ldb,const double beta,double * C,const int ldc)123 inline cusparseStatus_t cusparse_csrmm(cusparseHandle_t handle,
124 char Atrans,
125 int m,
126 int n,
127 int k,
128 int nnz,
129 const double alpha,
130 const cusparseMatDescr_t& descrA,
131 const double* csrValA,
132 const int* csrRowPtrA,
133 const int* csrColIndA,
134 const double* B,
135 const int ldb,
136 const double beta,
137 double* C,
138 const int ldc)
139
140 {
141 cusparseStatus_t sucess = cusparseDcsrmm(handle, cusparseOperation(Atrans), m, n, k, nnz, &alpha, descrA, csrValA,
142 csrRowPtrA, csrColIndA, B, ldb, &beta, C, ldc);
143 cudaDeviceSynchronize();
144 return sucess;
145 }
146
cusparse_csrmm(cusparseHandle_t handle,char Atrans,int m,int n,int k,int nnz,const float alpha,const cusparseMatDescr_t & descrA,const float * csrValA,const int * csrRowPtrA,const int * csrColIndA,const float * B,const int ldb,const float beta,float * C,const int ldc)147 inline cusparseStatus_t cusparse_csrmm(cusparseHandle_t handle,
148 char Atrans,
149 int m,
150 int n,
151 int k,
152 int nnz,
153 const float alpha,
154 const cusparseMatDescr_t& descrA,
155 const float* csrValA,
156 const int* csrRowPtrA,
157 const int* csrColIndA,
158 const float* B,
159 const int ldb,
160 const float beta,
161 float* C,
162 const int ldc)
163
164 {
165 cusparseStatus_t sucess = cusparseScsrmm(handle, cusparseOperation(Atrans), m, n, k, nnz, &alpha, descrA, csrValA,
166 csrRowPtrA, csrColIndA, B, ldb, &beta, C, ldc);
167 cudaDeviceSynchronize();
168 return sucess;
169 }
170
cusparse_csrmm(cusparseHandle_t handle,char Atrans,int m,int n,int k,int nnz,const std::complex<double> alpha,const cusparseMatDescr_t & descrA,const std::complex<double> * csrValA,const int * csrRowPtrA,const int * csrColIndA,const std::complex<double> * B,const int ldb,const std::complex<double> beta,std::complex<double> * C,const int ldc)171 inline cusparseStatus_t cusparse_csrmm(cusparseHandle_t handle,
172 char Atrans,
173 int m,
174 int n,
175 int k,
176 int nnz,
177 const std::complex<double> alpha,
178 const cusparseMatDescr_t& descrA,
179 const std::complex<double>* csrValA,
180 const int* csrRowPtrA,
181 const int* csrColIndA,
182 const std::complex<double>* B,
183 const int ldb,
184 const std::complex<double> beta,
185 std::complex<double>* C,
186 const int ldc)
187
188 {
189 cusparseStatus_t sucess =
190 cusparseZcsrmm(handle, cusparseOperation(Atrans), m, n, k, nnz, reinterpret_cast<cuDoubleComplex const*>(&alpha),
191 descrA, reinterpret_cast<cuDoubleComplex const*>(csrValA), csrRowPtrA, csrColIndA,
192 reinterpret_cast<cuDoubleComplex const*>(B), ldb, reinterpret_cast<cuDoubleComplex const*>(&beta),
193 reinterpret_cast<cuDoubleComplex*>(C), ldc);
194 cudaDeviceSynchronize();
195 return sucess;
196 }
197
cusparse_csrmm(cusparseHandle_t handle,char Atrans,int m,int n,int k,int nnz,const std::complex<float> alpha,const cusparseMatDescr_t & descrA,const std::complex<float> * csrValA,const int * csrRowPtrA,const int * csrColIndA,const std::complex<float> * B,const int ldb,const std::complex<float> beta,std::complex<float> * C,const int ldc)198 inline cusparseStatus_t cusparse_csrmm(cusparseHandle_t handle,
199 char Atrans,
200 int m,
201 int n,
202 int k,
203 int nnz,
204 const std::complex<float> alpha,
205 const cusparseMatDescr_t& descrA,
206 const std::complex<float>* csrValA,
207 const int* csrRowPtrA,
208 const int* csrColIndA,
209 const std::complex<float>* B,
210 const int ldb,
211 const std::complex<float> beta,
212 std::complex<float>* C,
213 const int ldc)
214
215 {
216 cusparseStatus_t sucess =
217 cusparseCcsrmm(handle, cusparseOperation(Atrans), m, n, k, nnz, reinterpret_cast<cuComplex const*>(&alpha),
218 descrA, reinterpret_cast<cuComplex const*>(csrValA), csrRowPtrA, csrColIndA,
219 reinterpret_cast<cuComplex const*>(B), ldb, reinterpret_cast<cuComplex const*>(&beta),
220 reinterpret_cast<cuComplex*>(C), ldc);
221 cudaDeviceSynchronize();
222 return sucess;
223 }
224
cusparse_csrmm2(cusparseHandle_t handle,char Atrans,char Btrans,int m,int n,int k,int nnz,const double alpha,const cusparseMatDescr_t & descrA,const double * csrValA,const int * csrRowPtrA,const int * csrColIndA,const double * B,const int ldb,const double beta,double * C,const int ldc)225 inline cusparseStatus_t cusparse_csrmm2(cusparseHandle_t handle,
226 char Atrans,
227 char Btrans,
228 int m,
229 int n,
230 int k,
231 int nnz,
232 const double alpha,
233 const cusparseMatDescr_t& descrA,
234 const double* csrValA,
235 const int* csrRowPtrA,
236 const int* csrColIndA,
237 const double* B,
238 const int ldb,
239 const double beta,
240 double* C,
241 const int ldc)
242
243 {
244 cusparseStatus_t sucess = cusparseDcsrmm2(handle, cusparseOperation(Atrans), cusparseOperation(Btrans), m, n, k, nnz,
245 &alpha, descrA, csrValA, csrRowPtrA, csrColIndA, B, ldb, &beta, C, ldc);
246 cudaDeviceSynchronize();
247 return sucess;
248 }
249
cusparse_csrmm2(cusparseHandle_t handle,char Atrans,char Btrans,int m,int n,int k,int nnz,const float alpha,const cusparseMatDescr_t & descrA,const float * csrValA,const int * csrRowPtrA,const int * csrColIndA,const float * B,const int ldb,const float beta,float * C,const int ldc)250 inline cusparseStatus_t cusparse_csrmm2(cusparseHandle_t handle,
251 char Atrans,
252 char Btrans,
253 int m,
254 int n,
255 int k,
256 int nnz,
257 const float alpha,
258 const cusparseMatDescr_t& descrA,
259 const float* csrValA,
260 const int* csrRowPtrA,
261 const int* csrColIndA,
262 const float* B,
263 const int ldb,
264 const float beta,
265 float* C,
266 const int ldc)
267
268 {
269 cusparseStatus_t sucess = cusparseScsrmm2(handle, cusparseOperation(Atrans), cusparseOperation(Btrans), m, n, k, nnz,
270 &alpha, descrA, csrValA, csrRowPtrA, csrColIndA, B, ldb, &beta, C, ldc);
271 cudaDeviceSynchronize();
272 return sucess;
273 }
274
cusparse_csrmm2(cusparseHandle_t handle,char Atrans,char Btrans,int m,int n,int k,int nnz,const std::complex<double> alpha,const cusparseMatDescr_t & descrA,const std::complex<double> * csrValA,const int * csrRowPtrA,const int * csrColIndA,const std::complex<double> * B,const int ldb,const std::complex<double> beta,std::complex<double> * C,const int ldc)275 inline cusparseStatus_t cusparse_csrmm2(cusparseHandle_t handle,
276 char Atrans,
277 char Btrans,
278 int m,
279 int n,
280 int k,
281 int nnz,
282 const std::complex<double> alpha,
283 const cusparseMatDescr_t& descrA,
284 const std::complex<double>* csrValA,
285 const int* csrRowPtrA,
286 const int* csrColIndA,
287 const std::complex<double>* B,
288 const int ldb,
289 const std::complex<double> beta,
290 std::complex<double>* C,
291 const int ldc)
292
293 {
294 cusparseStatus_t sucess =
295 cusparseZcsrmm2(handle, cusparseOperation(Atrans), cusparseOperation(Btrans), m, n, k, nnz,
296 reinterpret_cast<cuDoubleComplex const*>(&alpha), descrA,
297 reinterpret_cast<cuDoubleComplex const*>(csrValA), csrRowPtrA, csrColIndA,
298 reinterpret_cast<cuDoubleComplex const*>(B), ldb, reinterpret_cast<cuDoubleComplex const*>(&beta),
299 reinterpret_cast<cuDoubleComplex*>(C), ldc);
300 cudaDeviceSynchronize();
301 return sucess;
302 }
303
cusparse_csrmm2(cusparseHandle_t handle,char Atrans,char Btrans,int m,int n,int k,int nnz,const std::complex<float> alpha,const cusparseMatDescr_t & descrA,const std::complex<float> * csrValA,const int * csrRowPtrA,const int * csrColIndA,const std::complex<float> * B,const int ldb,const std::complex<float> beta,std::complex<float> * C,const int ldc)304 inline cusparseStatus_t cusparse_csrmm2(cusparseHandle_t handle,
305 char Atrans,
306 char Btrans,
307 int m,
308 int n,
309 int k,
310 int nnz,
311 const std::complex<float> alpha,
312 const cusparseMatDescr_t& descrA,
313 const std::complex<float>* csrValA,
314 const int* csrRowPtrA,
315 const int* csrColIndA,
316 const std::complex<float>* B,
317 const int ldb,
318 const std::complex<float> beta,
319 std::complex<float>* C,
320 const int ldc)
321
322 {
323 cusparseStatus_t sucess =
324 cusparseCcsrmm2(handle, cusparseOperation(Atrans), cusparseOperation(Btrans), m, n, k, nnz,
325 reinterpret_cast<cuComplex const*>(&alpha), descrA, reinterpret_cast<cuComplex const*>(csrValA),
326 csrRowPtrA, csrColIndA, reinterpret_cast<cuComplex const*>(B), ldb,
327 reinterpret_cast<cuComplex const*>(&beta), reinterpret_cast<cuComplex*>(C), ldc);
328 cudaDeviceSynchronize();
329 return sucess;
330 }
331
cusparse_gemmi(cusparseHandle_t handle,int m,int n,int k,int nnz,const double alpha,const double * A,const int lda,const double * cscValB,const int * cscColPtrB,const int * cscRowIndB,const double beta,double * C,const int ldc)332 inline cusparseStatus_t cusparse_gemmi(cusparseHandle_t handle,
333 int m,
334 int n,
335 int k,
336 int nnz,
337 const double alpha,
338 const double* A,
339 const int lda,
340 const double* cscValB,
341 const int* cscColPtrB,
342 const int* cscRowIndB,
343 const double beta,
344 double* C,
345 const int ldc)
346
347 {
348 cusparseStatus_t sucess =
349 cusparseDgemmi(handle, m, n, k, nnz, &alpha, A, lda, cscValB, cscColPtrB, cscRowIndB, &beta, C, ldc);
350 cudaDeviceSynchronize();
351 return sucess;
352 }
353
cusparse_gemmi(cusparseHandle_t handle,int m,int n,int k,int nnz,const float alpha,const float * A,const int lda,const float * cscValB,const int * cscColPtrB,const int * cscRowIndB,const float beta,float * C,const int ldc)354 inline cusparseStatus_t cusparse_gemmi(cusparseHandle_t handle,
355 int m,
356 int n,
357 int k,
358 int nnz,
359 const float alpha,
360 const float* A,
361 const int lda,
362 const float* cscValB,
363 const int* cscColPtrB,
364 const int* cscRowIndB,
365 const float beta,
366 float* C,
367 const int ldc)
368
369 {
370 cusparseStatus_t sucess =
371 cusparseSgemmi(handle, m, n, k, nnz, &alpha, A, lda, cscValB, cscColPtrB, cscRowIndB, &beta, C, ldc);
372 cudaDeviceSynchronize();
373 return sucess;
374 }
375
cusparse_gemmi(cusparseHandle_t handle,int m,int n,int k,int nnz,const std::complex<double> alpha,const std::complex<double> * A,const int lda,const std::complex<double> * cscValB,const int * cscColPtrB,const int * cscRowIndB,const std::complex<double> beta,std::complex<double> * C,const int ldc)376 inline cusparseStatus_t cusparse_gemmi(cusparseHandle_t handle,
377 int m,
378 int n,
379 int k,
380 int nnz,
381 const std::complex<double> alpha,
382 const std::complex<double>* A,
383 const int lda,
384 const std::complex<double>* cscValB,
385 const int* cscColPtrB,
386 const int* cscRowIndB,
387 const std::complex<double> beta,
388 std::complex<double>* C,
389 const int ldc)
390
391 {
392 cusparseStatus_t sucess =
393 cusparseZgemmi(handle, m, n, k, nnz, reinterpret_cast<cuDoubleComplex const*>(&alpha),
394 reinterpret_cast<cuDoubleComplex const*>(A), lda,
395 reinterpret_cast<cuDoubleComplex const*>(cscValB), cscColPtrB, cscRowIndB,
396 reinterpret_cast<cuDoubleComplex const*>(&beta), reinterpret_cast<cuDoubleComplex*>(C), ldc);
397 cudaDeviceSynchronize();
398 return sucess;
399 }
400
cusparse_gemmi(cusparseHandle_t handle,int m,int n,int k,int nnz,const std::complex<float> alpha,const std::complex<float> * A,const int lda,const std::complex<float> * cscValB,const int * cscColPtrB,const int * cscRowIndB,const std::complex<float> beta,std::complex<float> * C,const int ldc)401 inline cusparseStatus_t cusparse_gemmi(cusparseHandle_t handle,
402 int m,
403 int n,
404 int k,
405 int nnz,
406 const std::complex<float> alpha,
407 const std::complex<float>* A,
408 const int lda,
409 const std::complex<float>* cscValB,
410 const int* cscColPtrB,
411 const int* cscRowIndB,
412 const std::complex<float> beta,
413 std::complex<float>* C,
414 const int ldc)
415
416 {
417 cusparseStatus_t sucess =
418 cusparseCgemmi(handle, m, n, k, nnz, reinterpret_cast<cuComplex const*>(&alpha),
419 reinterpret_cast<cuComplex const*>(A), lda, reinterpret_cast<cuComplex const*>(cscValB),
420 cscColPtrB, cscRowIndB, reinterpret_cast<cuComplex const*>(&beta), reinterpret_cast<cuComplex*>(C),
421 ldc);
422 cudaDeviceSynchronize();
423 return sucess;
424 }
425
426 } // namespace cusparse
427
428 #endif
429