xref: /petsc/src/mat/impls/aij/seq/seqcusparse/cusparsematimpl.h (revision 47f8145d78ec2f038a1e36d827f5e466125d1312)
1687625d7SJacob Faibussowitsch #ifndef PETSC_CUSPARSEMATIMPL_H
2687625d7SJacob Faibussowitsch #define PETSC_CUSPARSEMATIMPL_H
39ae82921SPaul Mullowney 
4afb2bd1cSJunchao Zhang #include <petscpkg_version.h>
56d54fb17SJacob Faibussowitsch #include <../src/vec/vec/impls/seq/cupm/vecseqcupm.hpp> /* for VecSeq_CUPM */
6*47f8145dSJacob Faibussowitsch #include <../src/sys/objects/device/impls/cupm/cupmthrustutility.hpp>
76d54fb17SJacob Faibussowitsch #include <petsc/private/petsclegacycupmblas.h>
89ae82921SPaul Mullowney 
99ae82921SPaul Mullowney #include <cusparse_v2.h>
109ae82921SPaul Mullowney 
119ae82921SPaul Mullowney #include <algorithm>
129ae82921SPaul Mullowney #include <vector>
139ae82921SPaul Mullowney 
14c41cb2e2SAlejandro Lamas Daviña #include <thrust/device_vector.h>
15c41cb2e2SAlejandro Lamas Daviña #include <thrust/device_ptr.h>
1650ab121bSSatish Balay #include <thrust/device_malloc_allocator.h>
17c41cb2e2SAlejandro Lamas Daviña #include <thrust/transform.h>
18c41cb2e2SAlejandro Lamas Daviña #include <thrust/functional.h>
19554b8892SKarl Rupp #include <thrust/sequence.h>
207eaca502SStefano Zampini #include <thrust/system/system_error.h>
21c41cb2e2SAlejandro Lamas Daviña 
22aa372e3fSPaul Mullowney #if defined(PETSC_USE_COMPLEX)
23aa372e3fSPaul Mullowney   #if defined(PETSC_USE_REAL_SINGLE)
24ccdfe979SStefano Zampini const cuComplex PETSC_CUSPARSE_ONE  = {1.0f, 0.0f};
25ccdfe979SStefano Zampini const cuComplex PETSC_CUSPARSE_ZERO = {0.0f, 0.0f};
26da112707SJunchao 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)
27da112707SJunchao 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)
28da112707SJunchao 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)
29da112707SJunchao 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)
30da112707SJunchao 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)
31da112707SJunchao 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)
32afb2bd1cSJunchao Zhang   #elif defined(PETSC_USE_REAL_DOUBLE)
33afb2bd1cSJunchao Zhang const cuDoubleComplex PETSC_CUSPARSE_ONE  = {1.0, 0.0};
34afb2bd1cSJunchao Zhang const cuDoubleComplex PETSC_CUSPARSE_ZERO = {0.0, 0.0};
35da112707SJunchao 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)
36da112707SJunchao 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)
37da112707SJunchao 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)
38da112707SJunchao 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)
39da112707SJunchao 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)
40da112707SJunchao 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)
41afb2bd1cSJunchao Zhang   #endif
42afb2bd1cSJunchao Zhang #else
43afb2bd1cSJunchao Zhang const PetscScalar PETSC_CUSPARSE_ONE  = 1.0;
44afb2bd1cSJunchao Zhang const PetscScalar PETSC_CUSPARSE_ZERO = 0.0;
45da112707SJunchao Zhang   #if defined(PETSC_USE_REAL_SINGLE)
46da112707SJunchao Zhang     #define cusparseXcsrilu02_bufferSize cusparseScsrilu02_bufferSize
47da112707SJunchao Zhang     #define cusparseXcsrilu02_analysis   cusparseScsrilu02_analysis
48da112707SJunchao Zhang     #define cusparseXcsrilu02            cusparseScsrilu02
49da112707SJunchao Zhang     #define cusparseXcsric02_bufferSize  cusparseScsric02_bufferSize
50da112707SJunchao Zhang     #define cusparseXcsric02_analysis    cusparseScsric02_analysis
51da112707SJunchao Zhang     #define cusparseXcsric02             cusparseScsric02
52da112707SJunchao Zhang   #elif defined(PETSC_USE_REAL_DOUBLE)
53da112707SJunchao Zhang     #define cusparseXcsrilu02_bufferSize cusparseDcsrilu02_bufferSize
54da112707SJunchao Zhang     #define cusparseXcsrilu02_analysis   cusparseDcsrilu02_analysis
55da112707SJunchao Zhang     #define cusparseXcsrilu02            cusparseDcsrilu02
56da112707SJunchao Zhang     #define cusparseXcsric02_bufferSize  cusparseDcsric02_bufferSize
57da112707SJunchao Zhang     #define cusparseXcsric02_analysis    cusparseDcsric02_analysis
58da112707SJunchao Zhang     #define cusparseXcsric02             cusparseDcsric02
59da112707SJunchao Zhang   #endif
60afb2bd1cSJunchao Zhang #endif
61afb2bd1cSJunchao Zhang 
621b0a6780SStefano Zampini #if PETSC_PKG_CUDA_VERSION_GE(9, 0, 0)
63261a78b4SJunchao Zhang   #define csrsvInfo_t              csrsv2Info_t
64261a78b4SJunchao Zhang   #define cusparseCreateCsrsvInfo  cusparseCreateCsrsv2Info
65261a78b4SJunchao Zhang   #define cusparseDestroyCsrsvInfo cusparseDestroyCsrsv2Info
66afb2bd1cSJunchao Zhang   #if defined(PETSC_USE_COMPLEX)
67afb2bd1cSJunchao Zhang     #if defined(PETSC_USE_REAL_SINGLE)
68261a78b4SJunchao 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)
69261a78b4SJunchao 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)
70261a78b4SJunchao 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)
71afb2bd1cSJunchao Zhang     #elif defined(PETSC_USE_REAL_DOUBLE)
72261a78b4SJunchao 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)
73261a78b4SJunchao 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)
74261a78b4SJunchao 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)
75afb2bd1cSJunchao Zhang     #endif
76afb2bd1cSJunchao Zhang   #else /* not complex */
77afb2bd1cSJunchao Zhang     #if defined(PETSC_USE_REAL_SINGLE)
78261a78b4SJunchao Zhang       #define cusparseXcsrsv_buffsize cusparseScsrsv2_bufferSize
79261a78b4SJunchao Zhang       #define cusparseXcsrsv_analysis cusparseScsrsv2_analysis
80261a78b4SJunchao Zhang       #define cusparseXcsrsv_solve    cusparseScsrsv2_solve
81afb2bd1cSJunchao Zhang     #elif defined(PETSC_USE_REAL_DOUBLE)
82261a78b4SJunchao Zhang       #define cusparseXcsrsv_buffsize cusparseDcsrsv2_bufferSize
83261a78b4SJunchao Zhang       #define cusparseXcsrsv_analysis cusparseDcsrsv2_analysis
84261a78b4SJunchao Zhang       #define cusparseXcsrsv_solve    cusparseDcsrsv2_solve
85afb2bd1cSJunchao Zhang     #endif
86afb2bd1cSJunchao Zhang   #endif
879f7ba44dSJacob Faibussowitsch #else /* PETSC_PKG_CUDA_VERSION_GE(9, 0, 0) */
88261a78b4SJunchao Zhang   #define csrsvInfo_t              cusparseSolveAnalysisInfo_t
89261a78b4SJunchao Zhang   #define cusparseCreateCsrsvInfo  cusparseCreateSolveAnalysisInfo
90261a78b4SJunchao Zhang   #define cusparseDestroyCsrsvInfo cusparseDestroySolveAnalysisInfo
91afb2bd1cSJunchao Zhang   #if defined(PETSC_USE_COMPLEX)
92afb2bd1cSJunchao Zhang     #if defined(PETSC_USE_REAL_SINGLE)
939f7ba44dSJacob 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))
949f7ba44dSJacob 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))
951b0a6780SStefano Zampini     #elif defined(PETSC_USE_REAL_DOUBLE)
969f7ba44dSJacob Faibussowitsch       #define cusparseXcsrsv_solve(a, b, c, d_IGNORED, e, f, g, h, i, j, k, l, m_IGNORED, n_IGNORED) \
979f7ba44dSJacob Faibussowitsch         cusparseZcsrsv_solve((a), (b), (c), (cuDoubleComplex *)(e), (f), (cuDoubleComplex *)(g), (h), (i), (j), (cuDoubleComplex *)(k), (cuDoubleComplex *)(l))
989f7ba44dSJacob 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))
991b0a6780SStefano Zampini     #endif
1001b0a6780SStefano Zampini   #else /* not complex */
1011b0a6780SStefano Zampini     #if defined(PETSC_USE_REAL_SINGLE)
102261a78b4SJunchao Zhang       #define cusparseXcsrsv_solve                                     cusparseScsrsv_solve
1039f7ba44dSJacob 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)
1041b0a6780SStefano Zampini     #elif defined(PETSC_USE_REAL_DOUBLE)
105261a78b4SJunchao Zhang       #define cusparseXcsrsv_solve                                     cusparseDcsrsv_solve
1069f7ba44dSJacob 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)
1071b0a6780SStefano Zampini     #endif
1081b0a6780SStefano Zampini   #endif
1099f7ba44dSJacob Faibussowitsch #endif /* PETSC_PKG_CUDA_VERSION_GE(9, 0, 0) */
1101b0a6780SStefano Zampini 
1111b0a6780SStefano Zampini #if PETSC_PKG_CUDA_VERSION_GE(11, 0, 0)
1121b0a6780SStefano Zampini   #define cusparse_csr2csc cusparseCsr2cscEx2
1131b0a6780SStefano Zampini   #if defined(PETSC_USE_COMPLEX)
1141b0a6780SStefano Zampini     #if defined(PETSC_USE_REAL_SINGLE)
1151b0a6780SStefano Zampini       #define cusparse_scalartype                                                             CUDA_C_32F
116039c6fbaSStefano 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)
1179371c9d4SSatish 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) \
1189371c9d4SSatish 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)
1191b0a6780SStefano Zampini     #elif defined(PETSC_USE_REAL_DOUBLE)
1201b0a6780SStefano Zampini       #define cusparse_scalartype CUDA_C_64F
1219371c9d4SSatish 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) \
1229371c9d4SSatish 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)
1239371c9d4SSatish 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) \
1249371c9d4SSatish 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)
1251b0a6780SStefano Zampini     #endif
1261b0a6780SStefano Zampini   #else /* not complex */
1271b0a6780SStefano Zampini     #if defined(PETSC_USE_REAL_SINGLE)
1281b0a6780SStefano Zampini       #define cusparse_scalartype            CUDA_R_32F
129039c6fbaSStefano Zampini       #define cusparse_csr_spgeam            cusparseScsrgeam2
130039c6fbaSStefano Zampini       #define cusparse_csr_spgeam_bufferSize cusparseScsrgeam2_bufferSizeExt
1311b0a6780SStefano Zampini     #elif defined(PETSC_USE_REAL_DOUBLE)
1321b0a6780SStefano Zampini       #define cusparse_scalartype            CUDA_R_64F
133039c6fbaSStefano Zampini       #define cusparse_csr_spgeam            cusparseDcsrgeam2
134039c6fbaSStefano Zampini       #define cusparse_csr_spgeam_bufferSize cusparseDcsrgeam2_bufferSizeExt
1351b0a6780SStefano Zampini     #endif
1361b0a6780SStefano Zampini   #endif
1379f7ba44dSJacob Faibussowitsch #else /* PETSC_PKG_CUDA_VERSION_GE(11, 0, 0) */
1381b0a6780SStefano Zampini   #if defined(PETSC_USE_COMPLEX)
1391b0a6780SStefano Zampini     #if defined(PETSC_USE_REAL_SINGLE)
140ec42abe4SAlejandro 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))
141ccdfe979SStefano 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))
142ec42abe4SAlejandro 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))
143ec42abe4SAlejandro 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))
144ec42abe4SAlejandro 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))
145ec42abe4SAlejandro Lamas Daviña       #define cusparse_hyb2csr(a, b, c, d, e, f)                                              cusparseChyb2csr((a), (b), (c), (cuComplex *)(d), (e), (f))
146fcdce8c4SStefano 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)
147039c6fbaSStefano 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)
148aa372e3fSPaul Mullowney     #elif defined(PETSC_USE_REAL_DOUBLE)
149ec42abe4SAlejandro 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))
1509371c9d4SSatish Balay       #define cusparse_csr_spmm(a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p) \
1519371c9d4SSatish Balay         cusparseZcsrmm((a), (b), (c), (d), (e), (f), (cuDoubleComplex *)(g), (h), (cuDoubleComplex *)(i), (j), (k), (cuDoubleComplex *)(l), (m), (cuDoubleComplex *)(n), (cuDoubleComplex *)(o), (p))
152ec42abe4SAlejandro 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))
153ec42abe4SAlejandro 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))
154ec42abe4SAlejandro 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))
155ec42abe4SAlejandro Lamas Daviña       #define cusparse_hyb2csr(a, b, c, d, e, f)                                              cusparseZhyb2csr((a), (b), (c), (cuDoubleComplex *)(d), (e), (f))
156fcdce8c4SStefano 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)
1579371c9d4SSatish Balay       #define cusparse_csr_spgeam(a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q, r, s) \
1589371c9d4SSatish 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)
159aa372e3fSPaul Mullowney     #endif
160aa372e3fSPaul Mullowney   #else
161aa372e3fSPaul Mullowney     #if defined(PETSC_USE_REAL_SINGLE)
162aa372e3fSPaul Mullowney       #define cusparse_csr_spmv   cusparseScsrmv
163ccdfe979SStefano Zampini       #define cusparse_csr_spmm   cusparseScsrmm
164aa372e3fSPaul Mullowney       #define cusparse_csr2csc    cusparseScsr2csc
165aa372e3fSPaul Mullowney       #define cusparse_hyb_spmv   cusparseShybmv
166aa372e3fSPaul Mullowney       #define cusparse_csr2hyb    cusparseScsr2hyb
167aa372e3fSPaul Mullowney       #define cusparse_hyb2csr    cusparseShyb2csr
168fcdce8c4SStefano Zampini       #define cusparse_csr_spgemm cusparseScsrgemm
169039c6fbaSStefano Zampini       #define cusparse_csr_spgeam cusparseScsrgeam
170aa372e3fSPaul Mullowney     #elif defined(PETSC_USE_REAL_DOUBLE)
171aa372e3fSPaul Mullowney       #define cusparse_csr_spmv   cusparseDcsrmv
172ccdfe979SStefano Zampini       #define cusparse_csr_spmm   cusparseDcsrmm
173aa372e3fSPaul Mullowney       #define cusparse_csr2csc    cusparseDcsr2csc
174aa372e3fSPaul Mullowney       #define cusparse_hyb_spmv   cusparseDhybmv
175aa372e3fSPaul Mullowney       #define cusparse_csr2hyb    cusparseDcsr2hyb
176aa372e3fSPaul Mullowney       #define cusparse_hyb2csr    cusparseDhyb2csr
177fcdce8c4SStefano Zampini       #define cusparse_csr_spgemm cusparseDcsrgemm
178039c6fbaSStefano Zampini       #define cusparse_csr_spgeam cusparseDcsrgeam
179aa372e3fSPaul Mullowney     #endif
180aa372e3fSPaul Mullowney   #endif
1819f7ba44dSJacob Faibussowitsch #endif /* PETSC_PKG_CUDA_VERSION_GE(11, 0, 0) */
182aa372e3fSPaul Mullowney 
183aa372e3fSPaul Mullowney #define THRUSTINTARRAY32 thrust::device_vector<int>
184aa372e3fSPaul Mullowney #define THRUSTINTARRAY   thrust::device_vector<PetscInt>
185aa372e3fSPaul Mullowney #define THRUSTARRAY      thrust::device_vector<PetscScalar>
186aa372e3fSPaul Mullowney 
187aa372e3fSPaul Mullowney /* A CSR matrix structure */
188aa372e3fSPaul Mullowney struct CsrMatrix {
189aa372e3fSPaul Mullowney   PetscInt          num_rows;
190aa372e3fSPaul Mullowney   PetscInt          num_cols;
191aa372e3fSPaul Mullowney   PetscInt          num_entries;
192aa372e3fSPaul Mullowney   THRUSTINTARRAY32 *row_offsets;
193aa372e3fSPaul Mullowney   THRUSTINTARRAY32 *column_indices;
194aa372e3fSPaul Mullowney   THRUSTARRAY      *values;
1959ae82921SPaul Mullowney };
1969ae82921SPaul Mullowney 
197aa372e3fSPaul Mullowney /* This is struct holding the relevant data needed to a MatSolve */
198aa372e3fSPaul Mullowney struct Mat_SeqAIJCUSPARSETriFactorStruct {
199aa372e3fSPaul Mullowney   /* Data needed for triangular solve */
200aa372e3fSPaul Mullowney   cusparseMatDescr_t  descr;
201aa372e3fSPaul Mullowney   cusparseOperation_t solveOp;
202aa372e3fSPaul Mullowney   CsrMatrix          *csrMat;
203b917901dSJunchao Zhang #if PETSC_PKG_CUDA_VERSION_LT(11, 4, 0)
204261a78b4SJunchao Zhang   csrsvInfo_t           solveInfo;
2051b0a6780SStefano Zampini   cusparseSolvePolicy_t solvePolicy; /* whether level information is generated and used */
206d460d7bfSJunchao Zhang #endif
2071b0a6780SStefano Zampini   int          solveBufferSize;
2081b0a6780SStefano Zampini   void        *solveBuffer;
2091b0a6780SStefano Zampini   size_t       csr2cscBufferSize; /* to transpose the triangular factor (only used for CUDA >= 11.0) */
2101b0a6780SStefano Zampini   void        *csr2cscBuffer;
2112cbc15d9SMark   PetscScalar *AA_h; /* managed host buffer for moving values to the GPU */
212aa372e3fSPaul Mullowney };
213aa372e3fSPaul Mullowney 
214afb2bd1cSJunchao Zhang /* This is a larger struct holding all the triangular factors for a solve, transpose solve, and any indices used in a reordering */
215aa372e3fSPaul Mullowney struct Mat_SeqAIJCUSPARSETriFactors {
216b917901dSJunchao Zhang #if PETSC_PKG_CUDA_VERSION_LT(11, 4, 0)
217aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorPtr;          /* pointer for lower triangular (factored matrix) on GPU */
218aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorPtr;          /* pointer for upper triangular (factored matrix) on GPU */
219aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorPtrTranspose; /* pointer for lower triangular (factored matrix) on GPU for the transpose (useful for BiCG) */
220aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorPtrTranspose; /* pointer for upper triangular (factored matrix) on GPU for the transpose (useful for BiCG)*/
221d460d7bfSJunchao Zhang #endif
222d460d7bfSJunchao Zhang 
223aa372e3fSPaul Mullowney   THRUSTINTARRAY  *rpermIndices; /* indices used for any reordering */
224aa372e3fSPaul Mullowney   THRUSTINTARRAY  *cpermIndices; /* indices used for any reordering */
225aa372e3fSPaul Mullowney   THRUSTARRAY     *workVector;
226aa372e3fSPaul Mullowney   cusparseHandle_t handle; /* a handle to the cusparse library */
227aa372e3fSPaul Mullowney   PetscInt         nnz;    /* number of nonzeros ... need this for accurate logging between ICC and ILU */
228e8d2b73aSMark Adams   cudaDeviceProp   dev_prop;
229e8d2b73aSMark Adams   PetscBool        init_dev_prop;
230da112707SJunchao Zhang 
231d460d7bfSJunchao Zhang   PetscBool factorizeOnDevice; /* Do factorization on device or not */
232b917901dSJunchao Zhang #if PETSC_PKG_CUDA_VERSION_GE(11, 4, 0)
233da112707SJunchao Zhang   /* csrilu0/csric0 appeared in cusparse-8.0, but we use it along with cusparseSpSV,
234da112707SJunchao Zhang      which first appeared in cusparse-11.5 with cuda-11.3.
235da112707SJunchao Zhang   */
236d460d7bfSJunchao Zhang   PetscScalar *csrVal, *diag;             // the diagonal D in UtDU of Cholesky
237d460d7bfSJunchao Zhang   int         *csrRowPtr32, *csrColIdx32; // i,j of M. cusparseScsrilu02/ic02() etc require 32-bit indices
238d460d7bfSJunchao Zhang 
239d460d7bfSJunchao Zhang   PetscInt    *csrRowPtr, *csrColIdx; // i, j of M on device for CUDA APIs that support 64-bit indices
240d460d7bfSJunchao Zhang   PetscScalar *csrVal_h, *diag_h;     // Since LU is done on host, we prepare a factored matrix in regular csr format on host and then copy it to device
241d460d7bfSJunchao Zhang   PetscInt    *csrRowPtr_h;           // csrColIdx_h is temporary, so it is not here
242da112707SJunchao Zhang 
243da112707SJunchao Zhang   /* Mixed mat descriptor types? yes, different cusparse APIs use different types */
244da112707SJunchao Zhang   cusparseMatDescr_t   matDescr_M;
245da112707SJunchao Zhang   cusparseSpMatDescr_t spMatDescr_L, spMatDescr_U;
246da112707SJunchao Zhang   cusparseSpSVDescr_t  spsvDescr_L, spsvDescr_Lt, spsvDescr_U, spsvDescr_Ut;
247da112707SJunchao Zhang 
248da112707SJunchao Zhang   cusparseDnVecDescr_t dnVecDescr_X, dnVecDescr_Y;
249da112707SJunchao Zhang   PetscScalar         *X, *Y; /* data array of dnVec X and Y */
250da112707SJunchao Zhang 
251da112707SJunchao Zhang   /* Mixed size types? yes, CUDA-11.7.0 declared cusparseDcsrilu02_bufferSizeExt() that returns size_t but did not implement it! */
25212ba2bc6SJunchao Zhang   int    factBufferSize_M; /* M ~= LU or LLt */
253da112707SJunchao Zhang   size_t spsvBufferSize_L, spsvBufferSize_Lt, spsvBufferSize_U, spsvBufferSize_Ut;
25412ba2bc6SJunchao Zhang   /* cusparse needs various buffers for factorization and solve of L, U, Lt, or Ut.
25512ba2bc6SJunchao Zhang      So save memory, we share the factorization buffer with one of spsvBuffer_L/U.
25612ba2bc6SJunchao Zhang   */
25712ba2bc6SJunchao Zhang   void *factBuffer_M, *spsvBuffer_L, *spsvBuffer_U, *spsvBuffer_Lt, *spsvBuffer_Ut;
258da112707SJunchao Zhang 
259da112707SJunchao Zhang   csrilu02Info_t        ilu0Info_M;
260da112707SJunchao Zhang   csric02Info_t         ic0Info_M;
261da112707SJunchao Zhang   int                   structural_zero, numerical_zero;
262da112707SJunchao Zhang   cusparseSolvePolicy_t policy_M;
263da112707SJunchao Zhang 
26412ba2bc6SJunchao Zhang   /* In MatSolveTranspose() for ILU0, we use the two flags to do on-demand solve */
26512ba2bc6SJunchao Zhang   PetscBool createdTransposeSpSVDescr;    /* Have we created SpSV descriptors for Lt, Ut? */
26612ba2bc6SJunchao Zhang   PetscBool updatedTransposeSpSVAnalysis; /* Have we updated SpSV analysis with the latest L, U values? */
267da112707SJunchao Zhang 
268da112707SJunchao Zhang   PetscLogDouble numericFactFlops; /* Estimated FLOPs in ILU0/ICC0 numeric factorization */
269da112707SJunchao Zhang #endif
270aa372e3fSPaul Mullowney };
271aa372e3fSPaul Mullowney 
272afb2bd1cSJunchao Zhang struct Mat_CusparseSpMV {
273afb2bd1cSJunchao Zhang   PetscBool initialized;    /* Don't rely on spmvBuffer != NULL to test if the struct is initialized, */
274afb2bd1cSJunchao Zhang   size_t    spmvBufferSize; /* since I'm not sure if smvBuffer can be NULL even after cusparseSpMV_bufferSize() */
275afb2bd1cSJunchao Zhang   void     *spmvBuffer;
2769db3cbf9SStefano 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 */
277afb2bd1cSJunchao Zhang   cusparseDnVecDescr_t vecXDescr, vecYDescr; /* descriptor for the dense vectors in y=op(A)x */
2789db3cbf9SStefano Zampini #endif
279afb2bd1cSJunchao Zhang };
280afb2bd1cSJunchao Zhang 
281afb2bd1cSJunchao Zhang /* This is struct holding the relevant data needed to a MatMult */
282afb2bd1cSJunchao Zhang struct Mat_SeqAIJCUSPARSEMultStruct {
283afb2bd1cSJunchao Zhang   void              *mat;          /* opaque pointer to a matrix. This could be either a cusparseHybMat_t or a CsrMatrix */
284afb2bd1cSJunchao Zhang   cusparseMatDescr_t descr;        /* Data needed to describe the matrix for a multiply */
285afb2bd1cSJunchao Zhang   THRUSTINTARRAY    *cprowIndices; /* compressed row indices used in the parallel SpMV */
286afb2bd1cSJunchao Zhang   PetscScalar       *alpha_one;    /* pointer to a device "scalar" storing the alpha parameter in the SpMV */
287afb2bd1cSJunchao Zhang   PetscScalar       *beta_zero;    /* pointer to a device "scalar" storing the beta parameter in the SpMV as zero*/
288afb2bd1cSJunchao Zhang   PetscScalar       *beta_one;     /* pointer to a device "scalar" storing the beta parameter in the SpMV as one */
289afb2bd1cSJunchao Zhang #if PETSC_PKG_CUDA_VERSION_GE(11, 0, 0)
290afb2bd1cSJunchao Zhang   cusparseSpMatDescr_t matDescr;  /* descriptor for the matrix, used by SpMV and SpMM */
291afb2bd1cSJunchao Zhang   Mat_CusparseSpMV     cuSpMV[3]; /* different Mat_CusparseSpMV structs for non-transpose, transpose, conj-transpose */
292d71ae5a4SJacob Faibussowitsch   Mat_SeqAIJCUSPARSEMultStruct() : matDescr(NULL)
293d71ae5a4SJacob Faibussowitsch   {
294afb2bd1cSJunchao Zhang     for (int i = 0; i < 3; i++) cuSpMV[i].initialized = PETSC_FALSE;
295afb2bd1cSJunchao Zhang   }
296afb2bd1cSJunchao Zhang #endif
297afb2bd1cSJunchao Zhang };
298afb2bd1cSJunchao Zhang 
2997301b172SPierre Jolivet /* This is a larger struct holding all the matrices for a SpMV, and SpMV Transpose */
3009ae82921SPaul Mullowney struct Mat_SeqAIJCUSPARSE {
301aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSEMultStruct *mat;            /* pointer to the matrix on the GPU */
302aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSEMultStruct *matTranspose;   /* pointer to the matrix on the GPU (for the transpose ... useful for BiCG) */
303aa372e3fSPaul Mullowney   THRUSTARRAY                  *workVector;     /* pointer to a workvector to which we can copy the relevant indices of a vector we want to multiply */
304029808eeSJunchao Zhang   THRUSTINTARRAY32             *rowoffsets_gpu; /* rowoffsets on GPU in non-compressed-row format. It is used to convert CSR to CSC */
305213423ffSJunchao Zhang   PetscInt                      nrows;          /* number of rows of the matrix seen by GPU */
306e057df02SPaul Mullowney   MatCUSPARSEStorageFormat      format;         /* the storage format for the matrix on the device */
307365b711fSMark Adams   PetscBool                     use_cpu_solve;  /* Use AIJ_Seq (I)LU solve */
308aa372e3fSPaul Mullowney   cudaStream_t                  stream;         /* a stream for the parallel SpMV ... this is not owned and should not be deleted */
309aa372e3fSPaul Mullowney   cusparseHandle_t              handle;         /* a handle to the cusparse library ... this may not be owned (if we're working in parallel i.e. multiGPUs) */
31054da937aSStefano Zampini   PetscObjectState              nonzerostate;   /* track nonzero state to possibly recreate the GPU matrix */
311afb2bd1cSJunchao Zhang #if PETSC_PKG_CUDA_VERSION_GE(11, 0, 0)
312afb2bd1cSJunchao Zhang   size_t               csr2cscBufferSize; /* stuff used to compute the matTranspose above */
313fea72eb4SStefano Zampini   void                *csr2cscBuffer;     /* This is used as a C struct and is calloc'ed by PetscNew() */
314afb2bd1cSJunchao Zhang   cusparseCsr2CscAlg_t csr2cscAlg;        /* algorithms can be selected from command line options */
315afb2bd1cSJunchao Zhang   cusparseSpMVAlg_t    spmvAlg;
316afb2bd1cSJunchao Zhang   cusparseSpMMAlg_t    spmmAlg;
317afb2bd1cSJunchao Zhang #endif
318a49f1ed0SStefano Zampini   THRUSTINTARRAY *csr2csc_i;
3192c4ab24aSJunchao Zhang   THRUSTINTARRAY *coords; /* permutation array used in MatSeqAIJCUSPARSEMergeMats */
3209ae82921SPaul Mullowney };
3219ae82921SPaul Mullowney 
32247d993e7Ssuyashtn typedef struct Mat_SeqAIJCUSPARSETriFactors *Mat_SeqAIJCUSPARSETriFactors_p;
32347d993e7Ssuyashtn 
324d3ecb6c7SJed Brown PETSC_INTERN PetscErrorCode MatSeqAIJCUSPARSECopyToGPU(Mat);
325ed502f03SStefano Zampini PETSC_INTERN PetscErrorCode MatSeqAIJCUSPARSEMergeMats(Mat, Mat, MatReuse, Mat *);
3265f101d05SStefano Zampini PETSC_INTERN PetscErrorCode MatSeqAIJCUSPARSETriFactors_Reset(Mat_SeqAIJCUSPARSETriFactors_p *);
327ed502f03SStefano Zampini 
3286d54fb17SJacob Faibussowitsch using VecSeq_CUDA = Petsc::vec::cupm::impl::VecSeq_CUPM<Petsc::device::cupm::DeviceType::CUDA>;
3296d54fb17SJacob Faibussowitsch 
330d71ae5a4SJacob Faibussowitsch static inline bool isCudaMem(const void *data)
331d71ae5a4SJacob Faibussowitsch {
3326d54fb17SJacob Faibussowitsch   using namespace Petsc::device::cupm;
3336d54fb17SJacob Faibussowitsch   auto mtype = PETSC_MEMTYPE_HOST;
3346d54fb17SJacob Faibussowitsch 
3356d54fb17SJacob Faibussowitsch   PetscFunctionBegin;
3366d54fb17SJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, impl::Interface<DeviceType::CUDA>::PetscCUPMGetMemType(data, &mtype));
3376d54fb17SJacob Faibussowitsch   PetscFunctionReturn(PetscMemTypeDevice(mtype));
33808391a17SStefano Zampini }
33908391a17SStefano Zampini 
340687625d7SJacob Faibussowitsch #endif // PETSC_CUSPARSEMATIMPL_H
341