xref: /petsc/src/mat/impls/aij/seq/seqcusparse/cusparsematimpl.h (revision 9f7ba44d0dd30c921add853610c9bec7639b6901)
1042217e8SBarry Smith #if !defined(CUSPARSEMATIMPL)
2042217e8SBarry Smith   #define CUSPARSEMATIMPL
39ae82921SPaul Mullowney 
4afb2bd1cSJunchao Zhang   #include <petscpkg_version.h>
5303a667bSSatish Balay   #include <petsc/private/cudavecimpl.h>
6042217e8SBarry Smith   #include <petscaijdevice.h>
79ae82921SPaul Mullowney 
89ae82921SPaul Mullowney   #include <cusparse_v2.h>
99ae82921SPaul Mullowney 
109ae82921SPaul Mullowney   #include <algorithm>
119ae82921SPaul Mullowney   #include <vector>
129ae82921SPaul Mullowney 
13c41cb2e2SAlejandro Lamas Daviña   #include <thrust/device_vector.h>
14c41cb2e2SAlejandro Lamas Daviña   #include <thrust/device_ptr.h>
1550ab121bSSatish Balay   #include <thrust/device_malloc_allocator.h>
16c41cb2e2SAlejandro Lamas Daviña   #include <thrust/transform.h>
17c41cb2e2SAlejandro Lamas Daviña   #include <thrust/functional.h>
18554b8892SKarl Rupp   #include <thrust/sequence.h>
197eaca502SStefano Zampini   #include <thrust/system/system_error.h>
20c41cb2e2SAlejandro Lamas Daviña 
219371c9d4SSatish Balay   #define PetscCallThrust(body) \
229371c9d4SSatish Balay     do { \
237eaca502SStefano Zampini       try { \
247eaca502SStefano Zampini         body; \
25d71ae5a4SJacob Faibussowitsch       } catch (thrust::system_error & e) { \
26d71ae5a4SJacob Faibussowitsch         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in Thrust %s", e.what()); \
27d71ae5a4SJacob Faibussowitsch       } \
287eaca502SStefano Zampini     } while (0)
297eaca502SStefano Zampini 
30aa372e3fSPaul Mullowney   #if defined(PETSC_USE_COMPLEX)
31aa372e3fSPaul Mullowney     #if defined(PETSC_USE_REAL_SINGLE)
32ccdfe979SStefano Zampini const cuComplex PETSC_CUSPARSE_ONE  = {1.0f, 0.0f};
33ccdfe979SStefano Zampini const cuComplex PETSC_CUSPARSE_ZERO = {0.0f, 0.0f};
34da112707SJunchao Zhang       #define cusparseXcsrilu02_bufferSize(a, b, c, d, e, f, g, h, i)  cusparseCcsrilu02_bufferSize(a, b, c, d, (cuComplex *)e, f, g, h, i)
35da112707SJunchao Zhang       #define cusparseXcsrilu02_analysis(a, b, c, d, e, f, g, h, i, j) cusparseCcsrilu02_analysis(a, b, c, d, (cuComplex *)e, f, g, h, i, j)
36da112707SJunchao Zhang       #define cusparseXcsrilu02(a, b, c, d, e, f, g, h, i, j)          cusparseCcsrilu02(a, b, c, d, (cuComplex *)e, f, g, h, i, j)
37da112707SJunchao Zhang       #define cusparseXcsric02_bufferSize(a, b, c, d, e, f, g, h, i)   cusparseCcsric02_bufferSize(a, b, c, d, (cuComplex *)e, f, g, h, i)
38da112707SJunchao Zhang       #define cusparseXcsric02_analysis(a, b, c, d, e, f, g, h, i, j)  cusparseCcsric02_analysis(a, b, c, d, (cuComplex *)e, f, g, h, i, j)
39da112707SJunchao Zhang       #define cusparseXcsric02(a, b, c, d, e, f, g, h, i, j)           cusparseCcsric02(a, b, c, d, (cuComplex *)e, f, g, h, i, j)
40afb2bd1cSJunchao Zhang     #elif defined(PETSC_USE_REAL_DOUBLE)
41afb2bd1cSJunchao Zhang const cuDoubleComplex PETSC_CUSPARSE_ONE  = {1.0, 0.0};
42afb2bd1cSJunchao Zhang const cuDoubleComplex PETSC_CUSPARSE_ZERO = {0.0, 0.0};
43da112707SJunchao Zhang       #define cusparseXcsrilu02_bufferSize(a, b, c, d, e, f, g, h, i)  cusparseZcsrilu02_bufferSize(a, b, c, d, (cuDoubleComplex *)e, f, g, h, i)
44da112707SJunchao Zhang       #define cusparseXcsrilu02_analysis(a, b, c, d, e, f, g, h, i, j) cusparseZcsrilu02_analysis(a, b, c, d, (cuDoubleComplex *)e, f, g, h, i, j)
45da112707SJunchao Zhang       #define cusparseXcsrilu02(a, b, c, d, e, f, g, h, i, j)          cusparseZcsrilu02(a, b, c, d, (cuDoubleComplex *)e, f, g, h, i, j)
46da112707SJunchao Zhang       #define cusparseXcsric02_bufferSize(a, b, c, d, e, f, g, h, i)   cusparseZcsric02_bufferSize(a, b, c, d, (cuDoubleComplex *)e, f, g, h, i)
47da112707SJunchao Zhang       #define cusparseXcsric02_analysis(a, b, c, d, e, f, g, h, i, j)  cusparseZcsric02_analysis(a, b, c, d, (cuDoubleComplex *)e, f, g, h, i, j)
48da112707SJunchao Zhang       #define cusparseXcsric02(a, b, c, d, e, f, g, h, i, j)           cusparseZcsric02(a, b, c, d, (cuDoubleComplex *)e, f, g, h, i, j)
49afb2bd1cSJunchao Zhang     #endif
50afb2bd1cSJunchao Zhang   #else
51afb2bd1cSJunchao Zhang const PetscScalar PETSC_CUSPARSE_ONE  = 1.0;
52afb2bd1cSJunchao Zhang const PetscScalar PETSC_CUSPARSE_ZERO = 0.0;
53da112707SJunchao Zhang     #if defined(PETSC_USE_REAL_SINGLE)
54da112707SJunchao Zhang       #define cusparseXcsrilu02_bufferSize cusparseScsrilu02_bufferSize
55da112707SJunchao Zhang       #define cusparseXcsrilu02_analysis   cusparseScsrilu02_analysis
56da112707SJunchao Zhang       #define cusparseXcsrilu02            cusparseScsrilu02
57da112707SJunchao Zhang       #define cusparseXcsric02_bufferSize  cusparseScsric02_bufferSize
58da112707SJunchao Zhang       #define cusparseXcsric02_analysis    cusparseScsric02_analysis
59da112707SJunchao Zhang       #define cusparseXcsric02             cusparseScsric02
60da112707SJunchao Zhang     #elif defined(PETSC_USE_REAL_DOUBLE)
61da112707SJunchao Zhang       #define cusparseXcsrilu02_bufferSize cusparseDcsrilu02_bufferSize
62da112707SJunchao Zhang       #define cusparseXcsrilu02_analysis   cusparseDcsrilu02_analysis
63da112707SJunchao Zhang       #define cusparseXcsrilu02            cusparseDcsrilu02
64da112707SJunchao Zhang       #define cusparseXcsric02_bufferSize  cusparseDcsric02_bufferSize
65da112707SJunchao Zhang       #define cusparseXcsric02_analysis    cusparseDcsric02_analysis
66da112707SJunchao Zhang       #define cusparseXcsric02             cusparseDcsric02
67da112707SJunchao Zhang     #endif
68afb2bd1cSJunchao Zhang   #endif
69afb2bd1cSJunchao Zhang 
701b0a6780SStefano Zampini   #if PETSC_PKG_CUDA_VERSION_GE(9, 0, 0)
71261a78b4SJunchao Zhang     #define csrsvInfo_t              csrsv2Info_t
72261a78b4SJunchao Zhang     #define cusparseCreateCsrsvInfo  cusparseCreateCsrsv2Info
73261a78b4SJunchao Zhang     #define cusparseDestroyCsrsvInfo cusparseDestroyCsrsv2Info
74afb2bd1cSJunchao Zhang     #if defined(PETSC_USE_COMPLEX)
75afb2bd1cSJunchao Zhang       #if defined(PETSC_USE_REAL_SINGLE)
76261a78b4SJunchao Zhang         #define cusparseXcsrsv_buffsize(a, b, c, d, e, f, g, h, i, j)          cusparseCcsrsv2_bufferSize(a, b, c, d, e, (cuComplex *)(f), g, h, i, j)
77261a78b4SJunchao Zhang         #define cusparseXcsrsv_analysis(a, b, c, d, e, f, g, h, i, j, k)       cusparseCcsrsv2_analysis(a, b, c, d, e, (const cuComplex *)(f), g, h, i, j, k)
78261a78b4SJunchao Zhang         #define cusparseXcsrsv_solve(a, b, c, d, e, f, g, h, i, j, k, l, m, n) cusparseCcsrsv2_solve(a, b, c, d, (const cuComplex *)(e), f, (const cuComplex *)(g), h, i, j, (const cuComplex *)(k), (cuComplex *)(l), m, n)
79afb2bd1cSJunchao Zhang       #elif defined(PETSC_USE_REAL_DOUBLE)
80261a78b4SJunchao Zhang         #define cusparseXcsrsv_buffsize(a, b, c, d, e, f, g, h, i, j)          cusparseZcsrsv2_bufferSize(a, b, c, d, e, (cuDoubleComplex *)(f), g, h, i, j)
81261a78b4SJunchao Zhang         #define cusparseXcsrsv_analysis(a, b, c, d, e, f, g, h, i, j, k)       cusparseZcsrsv2_analysis(a, b, c, d, e, (const cuDoubleComplex *)(f), g, h, i, j, k)
82261a78b4SJunchao Zhang         #define cusparseXcsrsv_solve(a, b, c, d, e, f, g, h, i, j, k, l, m, n) cusparseZcsrsv2_solve(a, b, c, d, (const cuDoubleComplex *)(e), f, (const cuDoubleComplex *)(g), h, i, j, (const cuDoubleComplex *)(k), (cuDoubleComplex *)(l), m, n)
83afb2bd1cSJunchao Zhang       #endif
84afb2bd1cSJunchao Zhang     #else /* not complex */
85afb2bd1cSJunchao Zhang       #if defined(PETSC_USE_REAL_SINGLE)
86261a78b4SJunchao Zhang         #define cusparseXcsrsv_buffsize cusparseScsrsv2_bufferSize
87261a78b4SJunchao Zhang         #define cusparseXcsrsv_analysis cusparseScsrsv2_analysis
88261a78b4SJunchao Zhang         #define cusparseXcsrsv_solve    cusparseScsrsv2_solve
89afb2bd1cSJunchao Zhang       #elif defined(PETSC_USE_REAL_DOUBLE)
90261a78b4SJunchao Zhang         #define cusparseXcsrsv_buffsize cusparseDcsrsv2_bufferSize
91261a78b4SJunchao Zhang         #define cusparseXcsrsv_analysis cusparseDcsrsv2_analysis
92261a78b4SJunchao Zhang         #define cusparseXcsrsv_solve    cusparseDcsrsv2_solve
93afb2bd1cSJunchao Zhang       #endif
94afb2bd1cSJunchao Zhang     #endif
95*9f7ba44dSJacob Faibussowitsch   #else /* PETSC_PKG_CUDA_VERSION_GE(9, 0, 0) */
96261a78b4SJunchao Zhang     #define csrsvInfo_t              cusparseSolveAnalysisInfo_t
97261a78b4SJunchao Zhang     #define cusparseCreateCsrsvInfo  cusparseCreateSolveAnalysisInfo
98261a78b4SJunchao Zhang     #define cusparseDestroyCsrsvInfo cusparseDestroySolveAnalysisInfo
99afb2bd1cSJunchao Zhang     #if defined(PETSC_USE_COMPLEX)
100afb2bd1cSJunchao Zhang       #if defined(PETSC_USE_REAL_SINGLE)
101*9f7ba44dSJacob Faibussowitsch         #define cusparseXcsrsv_solve(a, b, c, d_IGNORED, e, f, g, h, i, j, k, l, m_IGNORED, n_IGNORED) cusparseCcsrsv_solve((a), (b), (c), (cuComplex *)(e), (f), (cuComplex *)(g), (h), (i), (j), (cuComplex *)(k), (cuComplex *)(l))
102*9f7ba44dSJacob Faibussowitsch         #define cusparseXcsrsv_analysis(a, b, c, d, e, f, g, h, i, j_IGNORED, k_IGNORED)               cusparseCcsrsv_analysis((a), (b), (c), (d), (e), (cuComplex *)(f), (g), (h), (i))
1031b0a6780SStefano Zampini       #elif defined(PETSC_USE_REAL_DOUBLE)
104*9f7ba44dSJacob Faibussowitsch         #define cusparseXcsrsv_solve(a, b, c, d_IGNORED, e, f, g, h, i, j, k, l, m_IGNORED, n_IGNORED) \
105*9f7ba44dSJacob Faibussowitsch           cusparseZcsrsv_solve((a), (b), (c), (cuDoubleComplex *)(e), (f), (cuDoubleComplex *)(g), (h), (i), (j), (cuDoubleComplex *)(k), (cuDoubleComplex *)(l))
106*9f7ba44dSJacob Faibussowitsch         #define cusparseXcsrsv_analysis(a, b, c, d, e, f, g, h, i, j_IGNORED, k_IGNORED) cusparseZcsrsv_analysis((a), (b), (c), (d), (e), (cuDoubleComplex *)(f), (g), (h), (i))
1071b0a6780SStefano Zampini       #endif
1081b0a6780SStefano Zampini     #else /* not complex */
1091b0a6780SStefano Zampini       #if defined(PETSC_USE_REAL_SINGLE)
110261a78b4SJunchao Zhang         #define cusparseXcsrsv_solve                                     cusparseScsrsv_solve
111*9f7ba44dSJacob Faibussowitsch         #define cusparseXcsrsv_analysis(a, b, c, d, e, f, g, h, i, j, k) cusparseScsrsv_analysis(a, b, c, d, e, f, g, h, i)
1121b0a6780SStefano Zampini       #elif defined(PETSC_USE_REAL_DOUBLE)
113261a78b4SJunchao Zhang         #define cusparseXcsrsv_solve                                     cusparseDcsrsv_solve
114*9f7ba44dSJacob Faibussowitsch         #define cusparseXcsrsv_analysis(a, b, c, d, e, f, g, h, i, j, k) cusparseDcsrsv_analysis(a, b, c, d, e, f, g, h, i)
1151b0a6780SStefano Zampini       #endif
1161b0a6780SStefano Zampini     #endif
117*9f7ba44dSJacob Faibussowitsch   #endif /* PETSC_PKG_CUDA_VERSION_GE(9, 0, 0) */
1181b0a6780SStefano Zampini 
1191b0a6780SStefano Zampini   #if PETSC_PKG_CUDA_VERSION_GE(11, 0, 0)
1201b0a6780SStefano Zampini     #define cusparse_csr2csc cusparseCsr2cscEx2
1211b0a6780SStefano Zampini     #if defined(PETSC_USE_COMPLEX)
1221b0a6780SStefano Zampini       #if defined(PETSC_USE_REAL_SINGLE)
1231b0a6780SStefano Zampini         #define cusparse_scalartype                                                             CUDA_C_32F
124039c6fbaSStefano Zampini         #define cusparse_csr_spgeam(a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q, r, s, t) cusparseCcsrgeam2(a, b, c, (cuComplex *)d, e, f, (cuComplex *)g, h, i, (cuComplex *)j, k, l, (cuComplex *)m, n, o, p, (cuComplex *)q, r, s, t)
1259371c9d4SSatish Balay         #define cusparse_csr_spgeam_bufferSize(a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q, r, s, t) \
1269371c9d4SSatish Balay           cusparseCcsrgeam2_bufferSizeExt(a, b, c, (cuComplex *)d, e, f, (cuComplex *)g, h, i, (cuComplex *)j, k, l, (cuComplex *)m, n, o, p, (cuComplex *)q, r, s, t)
1271b0a6780SStefano Zampini       #elif defined(PETSC_USE_REAL_DOUBLE)
1281b0a6780SStefano Zampini         #define cusparse_scalartype CUDA_C_64F
1299371c9d4SSatish Balay         #define cusparse_csr_spgeam(a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q, r, s, t) \
1309371c9d4SSatish Balay           cusparseZcsrgeam2(a, b, c, (cuDoubleComplex *)d, e, f, (cuDoubleComplex *)g, h, i, (cuDoubleComplex *)j, k, l, (cuDoubleComplex *)m, n, o, p, (cuDoubleComplex *)q, r, s, t)
1319371c9d4SSatish Balay         #define cusparse_csr_spgeam_bufferSize(a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q, r, s, t) \
1329371c9d4SSatish Balay           cusparseZcsrgeam2_bufferSizeExt(a, b, c, (cuDoubleComplex *)d, e, f, (cuDoubleComplex *)g, h, i, (cuDoubleComplex *)j, k, l, (cuDoubleComplex *)m, n, o, p, (cuDoubleComplex *)q, r, s, t)
1331b0a6780SStefano Zampini       #endif
1341b0a6780SStefano Zampini     #else /* not complex */
1351b0a6780SStefano Zampini       #if defined(PETSC_USE_REAL_SINGLE)
1361b0a6780SStefano Zampini         #define cusparse_scalartype            CUDA_R_32F
137039c6fbaSStefano Zampini         #define cusparse_csr_spgeam            cusparseScsrgeam2
138039c6fbaSStefano Zampini         #define cusparse_csr_spgeam_bufferSize cusparseScsrgeam2_bufferSizeExt
1391b0a6780SStefano Zampini       #elif defined(PETSC_USE_REAL_DOUBLE)
1401b0a6780SStefano Zampini         #define cusparse_scalartype            CUDA_R_64F
141039c6fbaSStefano Zampini         #define cusparse_csr_spgeam            cusparseDcsrgeam2
142039c6fbaSStefano Zampini         #define cusparse_csr_spgeam_bufferSize cusparseDcsrgeam2_bufferSizeExt
1431b0a6780SStefano Zampini       #endif
1441b0a6780SStefano Zampini     #endif
145*9f7ba44dSJacob Faibussowitsch   #else /* PETSC_PKG_CUDA_VERSION_GE(11, 0, 0) */
1461b0a6780SStefano Zampini     #if defined(PETSC_USE_COMPLEX)
1471b0a6780SStefano Zampini       #if defined(PETSC_USE_REAL_SINGLE)
148ec42abe4SAlejandro Lamas Daviña         #define cusparse_csr_spmv(a, b, c, d, e, f, g, h, i, j, k, l, m)                        cusparseCcsrmv((a), (b), (c), (d), (e), (cuComplex *)(f), (g), (cuComplex *)(h), (i), (j), (cuComplex *)(k), (cuComplex *)(l), (cuComplex *)(m))
149ccdfe979SStefano Zampini         #define cusparse_csr_spmm(a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p)               cusparseCcsrmm((a), (b), (c), (d), (e), (f), (cuComplex *)(g), (h), (cuComplex *)(i), (j), (k), (cuComplex *)(l), (m), (cuComplex *)(n), (cuComplex *)(o), (p))
150ec42abe4SAlejandro Lamas Daviña         #define cusparse_csr2csc(a, b, c, d, e, f, g, h, i, j, k, l)                            cusparseCcsr2csc((a), (b), (c), (d), (cuComplex *)(e), (f), (g), (cuComplex *)(h), (i), (j), (k), (l))
151ec42abe4SAlejandro Lamas Daviña         #define cusparse_hyb_spmv(a, b, c, d, e, f, g, h)                                       cusparseChybmv((a), (b), (cuComplex *)(c), (d), (e), (cuComplex *)(f), (cuComplex *)(g), (cuComplex *)(h))
152ec42abe4SAlejandro Lamas Daviña         #define cusparse_csr2hyb(a, b, c, d, e, f, g, h, i, j)                                  cusparseCcsr2hyb((a), (b), (c), (d), (cuComplex *)(e), (f), (g), (h), (i), (j))
153ec42abe4SAlejandro Lamas Daviña         #define cusparse_hyb2csr(a, b, c, d, e, f)                                              cusparseChyb2csr((a), (b), (c), (cuComplex *)(d), (e), (f))
154fcdce8c4SStefano Zampini         #define cusparse_csr_spgemm(a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q, r, s, t) cusparseCcsrgemm(a, b, c, d, e, f, g, h, (cuComplex *)i, j, k, l, m, (cuComplex *)n, o, p, q, (cuComplex *)r, s, t)
155039c6fbaSStefano Zampini         #define cusparse_csr_spgeam(a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q, r, s)    cusparseCcsrgeam(a, b, c, (cuComplex *)d, e, f, (cuComplex *)g, h, i, (cuComplex *)j, k, l, (cuComplex *)m, n, o, p, (cuComplex *)q, r, s)
156aa372e3fSPaul Mullowney       #elif defined(PETSC_USE_REAL_DOUBLE)
157ec42abe4SAlejandro Lamas Daviña         #define cusparse_csr_spmv(a, b, c, d, e, f, g, h, i, j, k, l, m) cusparseZcsrmv((a), (b), (c), (d), (e), (cuDoubleComplex *)(f), (g), (cuDoubleComplex *)(h), (i), (j), (cuDoubleComplex *)(k), (cuDoubleComplex *)(l), (cuDoubleComplex *)(m))
1589371c9d4SSatish Balay         #define cusparse_csr_spmm(a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p) \
1599371c9d4SSatish Balay           cusparseZcsrmm((a), (b), (c), (d), (e), (f), (cuDoubleComplex *)(g), (h), (cuDoubleComplex *)(i), (j), (k), (cuDoubleComplex *)(l), (m), (cuDoubleComplex *)(n), (cuDoubleComplex *)(o), (p))
160ec42abe4SAlejandro Lamas Daviña         #define cusparse_csr2csc(a, b, c, d, e, f, g, h, i, j, k, l)                            cusparseZcsr2csc((a), (b), (c), (d), (cuDoubleComplex *)(e), (f), (g), (cuDoubleComplex *)(h), (i), (j), (k), (l))
161ec42abe4SAlejandro Lamas Daviña         #define cusparse_hyb_spmv(a, b, c, d, e, f, g, h)                                       cusparseZhybmv((a), (b), (cuDoubleComplex *)(c), (d), (e), (cuDoubleComplex *)(f), (cuDoubleComplex *)(g), (cuDoubleComplex *)(h))
162ec42abe4SAlejandro Lamas Daviña         #define cusparse_csr2hyb(a, b, c, d, e, f, g, h, i, j)                                  cusparseZcsr2hyb((a), (b), (c), (d), (cuDoubleComplex *)(e), (f), (g), (h), (i), (j))
163ec42abe4SAlejandro Lamas Daviña         #define cusparse_hyb2csr(a, b, c, d, e, f)                                              cusparseZhyb2csr((a), (b), (c), (cuDoubleComplex *)(d), (e), (f))
164fcdce8c4SStefano Zampini         #define cusparse_csr_spgemm(a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q, r, s, t) cusparseZcsrgemm(a, b, c, d, e, f, g, h, (cuDoubleComplex *)i, j, k, l, m, (cuDoubleComplex *)n, o, p, q, (cuDoubleComplex *)r, s, t)
1659371c9d4SSatish Balay         #define cusparse_csr_spgeam(a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q, r, s) \
1669371c9d4SSatish Balay           cusparseZcsrgeam(a, b, c, (cuDoubleComplex *)d, e, f, (cuDoubleComplex *)g, h, i, (cuDoubleComplex *)j, k, l, (cuDoubleComplex *)m, n, o, p, (cuDoubleComplex *)q, r, s)
167aa372e3fSPaul Mullowney       #endif
168aa372e3fSPaul Mullowney     #else
169aa372e3fSPaul Mullowney       #if defined(PETSC_USE_REAL_SINGLE)
170aa372e3fSPaul Mullowney         #define cusparse_csr_spmv   cusparseScsrmv
171ccdfe979SStefano Zampini         #define cusparse_csr_spmm   cusparseScsrmm
172aa372e3fSPaul Mullowney         #define cusparse_csr2csc    cusparseScsr2csc
173aa372e3fSPaul Mullowney         #define cusparse_hyb_spmv   cusparseShybmv
174aa372e3fSPaul Mullowney         #define cusparse_csr2hyb    cusparseScsr2hyb
175aa372e3fSPaul Mullowney         #define cusparse_hyb2csr    cusparseShyb2csr
176fcdce8c4SStefano Zampini         #define cusparse_csr_spgemm cusparseScsrgemm
177039c6fbaSStefano Zampini         #define cusparse_csr_spgeam cusparseScsrgeam
178aa372e3fSPaul Mullowney       #elif defined(PETSC_USE_REAL_DOUBLE)
179aa372e3fSPaul Mullowney         #define cusparse_csr_spmv   cusparseDcsrmv
180ccdfe979SStefano Zampini         #define cusparse_csr_spmm   cusparseDcsrmm
181aa372e3fSPaul Mullowney         #define cusparse_csr2csc    cusparseDcsr2csc
182aa372e3fSPaul Mullowney         #define cusparse_hyb_spmv   cusparseDhybmv
183aa372e3fSPaul Mullowney         #define cusparse_csr2hyb    cusparseDcsr2hyb
184aa372e3fSPaul Mullowney         #define cusparse_hyb2csr    cusparseDhyb2csr
185fcdce8c4SStefano Zampini         #define cusparse_csr_spgemm cusparseDcsrgemm
186039c6fbaSStefano Zampini         #define cusparse_csr_spgeam cusparseDcsrgeam
187aa372e3fSPaul Mullowney       #endif
188aa372e3fSPaul Mullowney     #endif
189*9f7ba44dSJacob Faibussowitsch   #endif /* PETSC_PKG_CUDA_VERSION_GE(11, 0, 0) */
190aa372e3fSPaul Mullowney 
191aa372e3fSPaul Mullowney   #define THRUSTINTARRAY32 thrust::device_vector<int>
192aa372e3fSPaul Mullowney   #define THRUSTINTARRAY   thrust::device_vector<PetscInt>
193aa372e3fSPaul Mullowney   #define THRUSTARRAY      thrust::device_vector<PetscScalar>
194aa372e3fSPaul Mullowney 
195aa372e3fSPaul Mullowney /* A CSR matrix structure */
196aa372e3fSPaul Mullowney struct CsrMatrix {
197aa372e3fSPaul Mullowney   PetscInt          num_rows;
198aa372e3fSPaul Mullowney   PetscInt          num_cols;
199aa372e3fSPaul Mullowney   PetscInt          num_entries;
200aa372e3fSPaul Mullowney   THRUSTINTARRAY32 *row_offsets;
201aa372e3fSPaul Mullowney   THRUSTINTARRAY32 *column_indices;
202aa372e3fSPaul Mullowney   THRUSTARRAY      *values;
2039ae82921SPaul Mullowney };
2049ae82921SPaul Mullowney 
205aa372e3fSPaul Mullowney /* This is struct holding the relevant data needed to a MatSolve */
206aa372e3fSPaul Mullowney struct Mat_SeqAIJCUSPARSETriFactorStruct {
207aa372e3fSPaul Mullowney   /* Data needed for triangular solve */
208aa372e3fSPaul Mullowney   cusparseMatDescr_t    descr;
209aa372e3fSPaul Mullowney   cusparseOperation_t   solveOp;
210aa372e3fSPaul Mullowney   CsrMatrix            *csrMat;
211261a78b4SJunchao Zhang   csrsvInfo_t           solveInfo;
2121b0a6780SStefano Zampini   cusparseSolvePolicy_t solvePolicy; /* whether level information is generated and used */
2131b0a6780SStefano Zampini   int                   solveBufferSize;
2141b0a6780SStefano Zampini   void                 *solveBuffer;
2151b0a6780SStefano Zampini   size_t                csr2cscBufferSize; /* to transpose the triangular factor (only used for CUDA >= 11.0) */
2161b0a6780SStefano Zampini   void                 *csr2cscBuffer;
2172cbc15d9SMark   PetscScalar          *AA_h; /* managed host buffer for moving values to the GPU */
218aa372e3fSPaul Mullowney };
219aa372e3fSPaul Mullowney 
220afb2bd1cSJunchao Zhang /* This is a larger struct holding all the triangular factors for a solve, transpose solve, and any indices used in a reordering */
221aa372e3fSPaul Mullowney struct Mat_SeqAIJCUSPARSETriFactors {
222aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorPtr;          /* pointer for lower triangular (factored matrix) on GPU */
223aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorPtr;          /* pointer for upper triangular (factored matrix) on GPU */
224aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorPtrTranspose; /* pointer for lower triangular (factored matrix) on GPU for the transpose (useful for BiCG) */
225aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorPtrTranspose; /* pointer for upper triangular (factored matrix) on GPU for the transpose (useful for BiCG)*/
226aa372e3fSPaul Mullowney   THRUSTINTARRAY                    *rpermIndices;            /* indices used for any reordering */
227aa372e3fSPaul Mullowney   THRUSTINTARRAY                    *cpermIndices;            /* indices used for any reordering */
228aa372e3fSPaul Mullowney   THRUSTARRAY                       *workVector;
229aa372e3fSPaul Mullowney   cusparseHandle_t                   handle;   /* a handle to the cusparse library */
230aa372e3fSPaul Mullowney   PetscInt                           nnz;      /* number of nonzeros ... need this for accurate logging between ICC and ILU */
231bddcd29dSMark Adams   PetscScalar                       *a_band_d; /* GPU data for banded CSR LU factorization matrix diag(L)=1 */
232bddcd29dSMark Adams   int                               *i_band_d; /* this could be optimized away */
233e8d2b73aSMark Adams   cudaDeviceProp                     dev_prop;
234e8d2b73aSMark Adams   PetscBool                          init_dev_prop;
235da112707SJunchao Zhang 
236da112707SJunchao Zhang   /* csrilu0/csric0 appeared in cusparse-8.0, but we use it along with cusparseSpSV,
237da112707SJunchao Zhang      which first appeared in cusparse-11.5 with cuda-11.3.
238da112707SJunchao Zhang   */
239bc996fdcSJunchao Zhang   PetscBool factorizeOnDevice; /* Do factorization on device or not */
240da112707SJunchao Zhang   #if CUSPARSE_VERSION >= 11500
241da112707SJunchao Zhang   PetscScalar *csrVal;
242da112707SJunchao Zhang   int         *csrRowPtr, *csrColIdx; /* a,i,j of M. Using int since some cusparse APIs only support 32-bit indices */
243da112707SJunchao Zhang 
244da112707SJunchao Zhang   /* Mixed mat descriptor types? yes, different cusparse APIs use different types */
245da112707SJunchao Zhang   cusparseMatDescr_t   matDescr_M;
246da112707SJunchao Zhang   cusparseSpMatDescr_t spMatDescr_L, spMatDescr_U;
247da112707SJunchao Zhang   cusparseSpSVDescr_t  spsvDescr_L, spsvDescr_Lt, spsvDescr_U, spsvDescr_Ut;
248da112707SJunchao Zhang 
249da112707SJunchao Zhang   cusparseDnVecDescr_t dnVecDescr_X, dnVecDescr_Y;
250da112707SJunchao Zhang   PetscScalar         *X, *Y; /* data array of dnVec X and Y */
251da112707SJunchao Zhang 
252da112707SJunchao Zhang   /* Mixed size types? yes, CUDA-11.7.0 declared cusparseDcsrilu02_bufferSizeExt() that returns size_t but did not implement it! */
25312ba2bc6SJunchao Zhang   int    factBufferSize_M; /* M ~= LU or LLt */
254da112707SJunchao Zhang   size_t spsvBufferSize_L, spsvBufferSize_Lt, spsvBufferSize_U, spsvBufferSize_Ut;
25512ba2bc6SJunchao Zhang   /* cusparse needs various buffers for factorization and solve of L, U, Lt, or Ut.
25612ba2bc6SJunchao Zhang      So save memory, we share the factorization buffer with one of spsvBuffer_L/U.
25712ba2bc6SJunchao Zhang   */
25812ba2bc6SJunchao Zhang   void *factBuffer_M, *spsvBuffer_L, *spsvBuffer_U, *spsvBuffer_Lt, *spsvBuffer_Ut;
259da112707SJunchao Zhang 
260da112707SJunchao Zhang   csrilu02Info_t        ilu0Info_M;
261da112707SJunchao Zhang   csric02Info_t         ic0Info_M;
262da112707SJunchao Zhang   int                   structural_zero, numerical_zero;
263da112707SJunchao Zhang   cusparseSolvePolicy_t policy_M;
264da112707SJunchao Zhang 
26512ba2bc6SJunchao Zhang   /* In MatSolveTranspose() for ILU0, we use the two flags to do on-demand solve */
26612ba2bc6SJunchao Zhang   PetscBool createdTransposeSpSVDescr;    /* Have we created SpSV descriptors for Lt, Ut? */
26712ba2bc6SJunchao Zhang   PetscBool updatedTransposeSpSVAnalysis; /* Have we updated SpSV analysis with the latest L, U values? */
268da112707SJunchao Zhang 
269da112707SJunchao Zhang   PetscLogDouble numericFactFlops; /* Estimated FLOPs in ILU0/ICC0 numeric factorization */
270da112707SJunchao Zhang   #endif
271aa372e3fSPaul Mullowney };
272aa372e3fSPaul Mullowney 
273afb2bd1cSJunchao Zhang struct Mat_CusparseSpMV {
274afb2bd1cSJunchao Zhang   PetscBool initialized;    /* Don't rely on spmvBuffer != NULL to test if the struct is initialized, */
275afb2bd1cSJunchao Zhang   size_t    spmvBufferSize; /* since I'm not sure if smvBuffer can be NULL even after cusparseSpMV_bufferSize() */
276afb2bd1cSJunchao Zhang   void     *spmvBuffer;
2779db3cbf9SStefano Zampini   #if PETSC_PKG_CUDA_VERSION_GE(11, 0, 0)    /* these are present from CUDA 10.1, but PETSc code makes use of them from CUDA 11 on */
278afb2bd1cSJunchao Zhang   cusparseDnVecDescr_t vecXDescr, vecYDescr; /* descriptor for the dense vectors in y=op(A)x */
2799db3cbf9SStefano Zampini   #endif
280afb2bd1cSJunchao Zhang };
281afb2bd1cSJunchao Zhang 
282afb2bd1cSJunchao Zhang /* This is struct holding the relevant data needed to a MatMult */
283afb2bd1cSJunchao Zhang struct Mat_SeqAIJCUSPARSEMultStruct {
284afb2bd1cSJunchao Zhang   void              *mat;          /* opaque pointer to a matrix. This could be either a cusparseHybMat_t or a CsrMatrix */
285afb2bd1cSJunchao Zhang   cusparseMatDescr_t descr;        /* Data needed to describe the matrix for a multiply */
286afb2bd1cSJunchao Zhang   THRUSTINTARRAY    *cprowIndices; /* compressed row indices used in the parallel SpMV */
287afb2bd1cSJunchao Zhang   PetscScalar       *alpha_one;    /* pointer to a device "scalar" storing the alpha parameter in the SpMV */
288afb2bd1cSJunchao Zhang   PetscScalar       *beta_zero;    /* pointer to a device "scalar" storing the beta parameter in the SpMV as zero*/
289afb2bd1cSJunchao Zhang   PetscScalar       *beta_one;     /* pointer to a device "scalar" storing the beta parameter in the SpMV as one */
290afb2bd1cSJunchao Zhang   #if PETSC_PKG_CUDA_VERSION_GE(11, 0, 0)
291afb2bd1cSJunchao Zhang   cusparseSpMatDescr_t matDescr;  /* descriptor for the matrix, used by SpMV and SpMM */
292afb2bd1cSJunchao Zhang   Mat_CusparseSpMV     cuSpMV[3]; /* different Mat_CusparseSpMV structs for non-transpose, transpose, conj-transpose */
293d71ae5a4SJacob Faibussowitsch   Mat_SeqAIJCUSPARSEMultStruct() : matDescr(NULL)
294d71ae5a4SJacob Faibussowitsch   {
295afb2bd1cSJunchao Zhang     for (int i = 0; i < 3; i++) cuSpMV[i].initialized = PETSC_FALSE;
296afb2bd1cSJunchao Zhang   }
297afb2bd1cSJunchao Zhang   #endif
298afb2bd1cSJunchao Zhang };
299afb2bd1cSJunchao Zhang 
3007301b172SPierre Jolivet /* This is a larger struct holding all the matrices for a SpMV, and SpMV Transpose */
3019ae82921SPaul Mullowney struct Mat_SeqAIJCUSPARSE {
302aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSEMultStruct *mat;            /* pointer to the matrix on the GPU */
303aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSEMultStruct *matTranspose;   /* pointer to the matrix on the GPU (for the transpose ... useful for BiCG) */
304aa372e3fSPaul Mullowney   THRUSTARRAY                  *workVector;     /* pointer to a workvector to which we can copy the relevant indices of a vector we want to multiply */
305029808eeSJunchao Zhang   THRUSTINTARRAY32             *rowoffsets_gpu; /* rowoffsets on GPU in non-compressed-row format. It is used to convert CSR to CSC */
306213423ffSJunchao Zhang   PetscInt                      nrows;          /* number of rows of the matrix seen by GPU */
307e057df02SPaul Mullowney   MatCUSPARSEStorageFormat      format;         /* the storage format for the matrix on the device */
308365b711fSMark Adams   PetscBool                     use_cpu_solve;  /* Use AIJ_Seq (I)LU solve */
309aa372e3fSPaul Mullowney   cudaStream_t                  stream;         /* a stream for the parallel SpMV ... this is not owned and should not be deleted */
310aa372e3fSPaul Mullowney   cusparseHandle_t              handle;         /* a handle to the cusparse library ... this may not be owned (if we're working in parallel i.e. multiGPUs) */
31154da937aSStefano Zampini   PetscObjectState              nonzerostate;   /* track nonzero state to possibly recreate the GPU matrix */
312afb2bd1cSJunchao Zhang   #if PETSC_PKG_CUDA_VERSION_GE(11, 0, 0)
313afb2bd1cSJunchao Zhang   size_t               csr2cscBufferSize; /* stuff used to compute the matTranspose above */
314afb2bd1cSJunchao Zhang   void                *csr2cscBuffer;     /* This is used as a C struct and is calloc'ed by PetscNewLog() */
315afb2bd1cSJunchao Zhang   cusparseCsr2CscAlg_t csr2cscAlg;        /* algorithms can be selected from command line options */
316afb2bd1cSJunchao Zhang   cusparseSpMVAlg_t    spmvAlg;
317afb2bd1cSJunchao Zhang   cusparseSpMMAlg_t    spmmAlg;
318afb2bd1cSJunchao Zhang   #endif
319a49f1ed0SStefano Zampini   THRUSTINTARRAY            *csr2csc_i;
320042217e8SBarry Smith   PetscSplitCSRDataStructure deviceMat; /* Matrix on device for, eg, assembly */
321219fbbafSJunchao Zhang 
322219fbbafSJunchao Zhang   /* Stuff for basic COO support */
323ddea5d60SJunchao Zhang   THRUSTINTARRAY *cooPerm;   /* permutation array that sorts the input coo entris by row and col */
324ddea5d60SJunchao Zhang   THRUSTINTARRAY *cooPerm_a; /* ordered array that indicate i-th nonzero (after sorting) is the j-th unique nonzero */
325219fbbafSJunchao Zhang 
326219fbbafSJunchao Zhang   /* Stuff for extended COO support */
327219fbbafSJunchao Zhang   PetscBool   use_extended_coo; /* Use extended COO format */
328219fbbafSJunchao Zhang   PetscCount *jmap_d;           /* perm[disp+jmap[i]..disp+jmap[i+1]) gives indices of entries in v[] associated with i-th nonzero of the matrix */
329219fbbafSJunchao Zhang   PetscCount *perm_d;
330219fbbafSJunchao Zhang 
331219fbbafSJunchao Zhang   Mat_SeqAIJCUSPARSE() : use_extended_coo(PETSC_FALSE), perm_d(NULL), jmap_d(NULL) { }
3329ae82921SPaul Mullowney };
3339ae82921SPaul Mullowney 
334d3ecb6c7SJed Brown PETSC_INTERN PetscErrorCode MatSeqAIJCUSPARSECopyToGPU(Mat);
335e8729f6fSJunchao Zhang PETSC_INTERN PetscErrorCode MatSetPreallocationCOO_SeqAIJCUSPARSE_Basic(Mat, PetscCount, PetscInt[], PetscInt[]);
336219fbbafSJunchao Zhang PETSC_INTERN PetscErrorCode MatSetValuesCOO_SeqAIJCUSPARSE_Basic(Mat, const PetscScalar[], InsertMode);
337ed502f03SStefano Zampini PETSC_INTERN PetscErrorCode MatSeqAIJCUSPARSEMergeMats(Mat, Mat, MatReuse, Mat *);
3385f101d05SStefano Zampini PETSC_INTERN PetscErrorCode MatSeqAIJCUSPARSETriFactors_Reset(Mat_SeqAIJCUSPARSETriFactors_p *);
339ed502f03SStefano Zampini 
340d71ae5a4SJacob Faibussowitsch static inline bool isCudaMem(const void *data)
341d71ae5a4SJacob Faibussowitsch {
34208391a17SStefano Zampini   cudaError_t                  cerr;
34308391a17SStefano Zampini   struct cudaPointerAttributes attr;
34408391a17SStefano Zampini   enum cudaMemoryType          mtype;
34508391a17SStefano Zampini   cerr = cudaPointerGetAttributes(&attr, data); /* Do not check error since before CUDA 11.0, passing a host pointer returns cudaErrorInvalidValue */
34608391a17SStefano Zampini   cudaGetLastError();                           /* Reset the last error */
34708391a17SStefano Zampini   #if (CUDART_VERSION < 10000)
34808391a17SStefano Zampini   mtype = attr.memoryType;
34908391a17SStefano Zampini   #else
35008391a17SStefano Zampini   mtype = attr.type;
35108391a17SStefano Zampini   #endif
35208391a17SStefano Zampini   if (cerr == cudaSuccess && mtype == cudaMemoryTypeDevice) return true;
35308391a17SStefano Zampini   else return false;
35408391a17SStefano Zampini }
35508391a17SStefano Zampini 
3569ae82921SPaul Mullowney #endif
357