xref: /petsc/src/mat/impls/aij/seq/seqcusparse/cusparsematimpl.h (revision 12ba2bc62b18c623c4aa6aaa43d22510cfab07b4)
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 
217eaca502SStefano Zampini #define PetscStackCallThrust(body) do {                                     \
227eaca502SStefano Zampini     try {                                                                   \
237eaca502SStefano Zampini       body;                                                                 \
247eaca502SStefano Zampini     } catch(thrust::system_error& e) {                                      \
2598921bdaSJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Thrust %s",e.what());\
267eaca502SStefano Zampini     }                                                                       \
277eaca502SStefano Zampini   } while (0)
287eaca502SStefano Zampini 
29aa372e3fSPaul Mullowney #if defined(PETSC_USE_COMPLEX)
30aa372e3fSPaul Mullowney   #if defined(PETSC_USE_REAL_SINGLE)
31ccdfe979SStefano Zampini     const cuComplex PETSC_CUSPARSE_ONE        = {1.0f, 0.0f};
32ccdfe979SStefano Zampini     const cuComplex PETSC_CUSPARSE_ZERO       = {0.0f, 0.0f};
33da112707SJunchao 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)
34da112707SJunchao 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)
35da112707SJunchao 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)
36da112707SJunchao 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)
37da112707SJunchao 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)
38da112707SJunchao 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)
39afb2bd1cSJunchao Zhang   #elif defined(PETSC_USE_REAL_DOUBLE)
40afb2bd1cSJunchao Zhang     const cuDoubleComplex PETSC_CUSPARSE_ONE  = {1.0, 0.0};
41afb2bd1cSJunchao Zhang     const cuDoubleComplex PETSC_CUSPARSE_ZERO = {0.0, 0.0};
42da112707SJunchao 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)
43da112707SJunchao 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)
44da112707SJunchao 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)
45da112707SJunchao 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)
46da112707SJunchao 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)
47da112707SJunchao 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)
48afb2bd1cSJunchao Zhang   #endif
49afb2bd1cSJunchao Zhang #else
50afb2bd1cSJunchao Zhang   const PetscScalar PETSC_CUSPARSE_ONE        = 1.0;
51afb2bd1cSJunchao Zhang   const PetscScalar PETSC_CUSPARSE_ZERO       = 0.0;
52da112707SJunchao Zhang   #if defined(PETSC_USE_REAL_SINGLE)
53da112707SJunchao Zhang     #define cusparseXcsrilu02_bufferSize    cusparseScsrilu02_bufferSize
54da112707SJunchao Zhang     #define cusparseXcsrilu02_analysis      cusparseScsrilu02_analysis
55da112707SJunchao Zhang     #define cusparseXcsrilu02               cusparseScsrilu02
56da112707SJunchao Zhang     #define cusparseXcsric02_bufferSize     cusparseScsric02_bufferSize
57da112707SJunchao Zhang     #define cusparseXcsric02_analysis       cusparseScsric02_analysis
58da112707SJunchao Zhang     #define cusparseXcsric02                cusparseScsric02
59da112707SJunchao Zhang   #elif defined(PETSC_USE_REAL_DOUBLE)
60da112707SJunchao Zhang     #define cusparseXcsrilu02_bufferSize    cusparseDcsrilu02_bufferSize
61da112707SJunchao Zhang     #define cusparseXcsrilu02_analysis      cusparseDcsrilu02_analysis
62da112707SJunchao Zhang     #define cusparseXcsrilu02               cusparseDcsrilu02
63da112707SJunchao Zhang     #define cusparseXcsric02_bufferSize     cusparseDcsric02_bufferSize
64da112707SJunchao Zhang     #define cusparseXcsric02_analysis       cusparseDcsric02_analysis
65da112707SJunchao Zhang     #define cusparseXcsric02                cusparseDcsric02
66da112707SJunchao Zhang   #endif
67afb2bd1cSJunchao Zhang #endif
68afb2bd1cSJunchao Zhang 
691b0a6780SStefano Zampini #if PETSC_PKG_CUDA_VERSION_GE(9,0,0)
70261a78b4SJunchao Zhang   #define csrsvInfo_t              csrsv2Info_t
71261a78b4SJunchao Zhang   #define cusparseCreateCsrsvInfo  cusparseCreateCsrsv2Info
72261a78b4SJunchao Zhang   #define cusparseDestroyCsrsvInfo cusparseDestroyCsrsv2Info
73afb2bd1cSJunchao Zhang   #if defined(PETSC_USE_COMPLEX)
74afb2bd1cSJunchao Zhang     #if defined(PETSC_USE_REAL_SINGLE)
75261a78b4SJunchao 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)
76261a78b4SJunchao 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)
77261a78b4SJunchao 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)
78afb2bd1cSJunchao Zhang     #elif defined(PETSC_USE_REAL_DOUBLE)
79261a78b4SJunchao 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)
80261a78b4SJunchao 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)
81261a78b4SJunchao 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)
82afb2bd1cSJunchao Zhang     #endif
83afb2bd1cSJunchao Zhang   #else /* not complex */
84afb2bd1cSJunchao Zhang     #if defined(PETSC_USE_REAL_SINGLE)
85261a78b4SJunchao Zhang       #define cusparseXcsrsv_buffsize       cusparseScsrsv2_bufferSize
86261a78b4SJunchao Zhang       #define cusparseXcsrsv_analysis       cusparseScsrsv2_analysis
87261a78b4SJunchao Zhang       #define cusparseXcsrsv_solve          cusparseScsrsv2_solve
88afb2bd1cSJunchao Zhang     #elif defined(PETSC_USE_REAL_DOUBLE)
89261a78b4SJunchao Zhang       #define cusparseXcsrsv_buffsize       cusparseDcsrsv2_bufferSize
90261a78b4SJunchao Zhang       #define cusparseXcsrsv_analysis       cusparseDcsrsv2_analysis
91261a78b4SJunchao Zhang       #define cusparseXcsrsv_solve          cusparseDcsrsv2_solve
92afb2bd1cSJunchao Zhang     #endif
93afb2bd1cSJunchao Zhang   #endif
94afb2bd1cSJunchao Zhang #else
95261a78b4SJunchao Zhang   #define csrsvInfo_t              cusparseSolveAnalysisInfo_t
96261a78b4SJunchao Zhang   #define cusparseCreateCsrsvInfo  cusparseCreateSolveAnalysisInfo
97261a78b4SJunchao Zhang   #define cusparseDestroyCsrsvInfo cusparseDestroySolveAnalysisInfo
98afb2bd1cSJunchao Zhang   #if defined(PETSC_USE_COMPLEX)
99afb2bd1cSJunchao Zhang     #if defined(PETSC_USE_REAL_SINGLE)
100261a78b4SJunchao Zhang       #define cusparseXcsrsv_solve(a,b,c,d,e,f,g,h,i,j,k) cusparseCcsrsv_solve((a),(b),(c),(cuComplex*)(d),(e),(cuComplex*)(f),(g),(h),(i),(cuComplex*)(j),(cuComplex*)(k))
101261a78b4SJunchao Zhang       #define cusparseXcsrsv_analysis(a,b,c,d,e,f,g,h,i)  cusparseCcsrsv_analysis((a),(b),(c),(d),(e),(cuComplex*)(f),(g),(h),(i))
1021b0a6780SStefano Zampini     #elif defined(PETSC_USE_REAL_DOUBLE)
103261a78b4SJunchao Zhang       #define cusparseXcsrsv_solve(a,b,c,d,e,f,g,h,i,j,k) cusparseZcsrsv_solve((a),(b),(c),(cuDoubleComplex*)(d),(e),(cuDoubleComplex*)(f),(g),(h),(i),(cuDoubleComplex*)(j),(cuDoubleComplex*)(k))
104261a78b4SJunchao Zhang       #define cusparseXcsrsv_analysis(a,b,c,d,e,f,g,h,i)  cusparseZcsrsv_analysis((a),(b),(c),(d),(e),(cuDoubleComplex*)(f),(g),(h),(i))
1051b0a6780SStefano Zampini     #endif
1061b0a6780SStefano Zampini   #else /* not complex */
1071b0a6780SStefano Zampini     #if defined(PETSC_USE_REAL_SINGLE)
108261a78b4SJunchao Zhang       #define cusparseXcsrsv_solve    cusparseScsrsv_solve
109261a78b4SJunchao Zhang       #define cusparseXcsrsv_analysis cusparseScsrsv_analysis
1101b0a6780SStefano Zampini     #elif defined(PETSC_USE_REAL_DOUBLE)
111261a78b4SJunchao Zhang       #define cusparseXcsrsv_solve    cusparseDcsrsv_solve
112261a78b4SJunchao Zhang       #define cusparseXcsrsv_analysis cusparseDcsrsv_analysis
1131b0a6780SStefano Zampini     #endif
1141b0a6780SStefano Zampini   #endif
1151b0a6780SStefano Zampini #endif
1161b0a6780SStefano Zampini 
1171b0a6780SStefano Zampini #if PETSC_PKG_CUDA_VERSION_GE(11,0,0)
1181b0a6780SStefano Zampini   #define cusparse_csr2csc cusparseCsr2cscEx2
1191b0a6780SStefano Zampini   #if defined(PETSC_USE_COMPLEX)
1201b0a6780SStefano Zampini     #if defined(PETSC_USE_REAL_SINGLE)
1211b0a6780SStefano Zampini       #define cusparse_scalartype CUDA_C_32F
122039c6fbaSStefano 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)
123039c6fbaSStefano Zampini       #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) 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)
1241b0a6780SStefano Zampini     #elif defined(PETSC_USE_REAL_DOUBLE)
1251b0a6780SStefano Zampini       #define cusparse_scalartype CUDA_C_64F
126039c6fbaSStefano 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)            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)
127039c6fbaSStefano Zampini       #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) 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)
1281b0a6780SStefano Zampini     #endif
1291b0a6780SStefano Zampini   #else /* not complex */
1301b0a6780SStefano Zampini     #if defined(PETSC_USE_REAL_SINGLE)
1311b0a6780SStefano Zampini       #define cusparse_scalartype CUDA_R_32F
132039c6fbaSStefano Zampini       #define cusparse_csr_spgeam            cusparseScsrgeam2
133039c6fbaSStefano Zampini       #define cusparse_csr_spgeam_bufferSize cusparseScsrgeam2_bufferSizeExt
1341b0a6780SStefano Zampini     #elif defined(PETSC_USE_REAL_DOUBLE)
1351b0a6780SStefano Zampini       #define cusparse_scalartype CUDA_R_64F
136039c6fbaSStefano Zampini       #define cusparse_csr_spgeam            cusparseDcsrgeam2
137039c6fbaSStefano Zampini       #define cusparse_csr_spgeam_bufferSize cusparseDcsrgeam2_bufferSizeExt
1381b0a6780SStefano Zampini     #endif
1391b0a6780SStefano Zampini   #endif
1401b0a6780SStefano Zampini #else
1411b0a6780SStefano Zampini   #if defined(PETSC_USE_COMPLEX)
1421b0a6780SStefano Zampini     #if defined(PETSC_USE_REAL_SINGLE)
143ec42abe4SAlejandro 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))
144ccdfe979SStefano 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))
145ec42abe4SAlejandro 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))
146ec42abe4SAlejandro 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))
147ec42abe4SAlejandro 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))
148ec42abe4SAlejandro Lamas Daviña       #define cusparse_hyb2csr(a,b,c,d,e,f)                      cusparseChyb2csr((a),(b),(c),(cuComplex*)(d),(e),(f))
149fcdce8c4SStefano 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)
150039c6fbaSStefano 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)
151aa372e3fSPaul Mullowney     #elif defined(PETSC_USE_REAL_DOUBLE)
152ec42abe4SAlejandro 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))
153ccdfe979SStefano Zampini       #define cusparse_csr_spmm(a,b,c,d,e,f,g,h,i,j,k,l,m,n,o,p) cusparseZcsrmm((a),(b),(c),(d),(e),(f),(cuDoubleComplex*)(g),(h),(cuDoubleComplex*)(i),(j),(k),(cuDoubleComplex*)(l),(m),(cuDoubleComplex*)(n),(cuDoubleComplex*)(o),(p))
154ec42abe4SAlejandro 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))
155ec42abe4SAlejandro 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))
156ec42abe4SAlejandro 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))
157ec42abe4SAlejandro Lamas Daviña       #define cusparse_hyb2csr(a,b,c,d,e,f)                      cusparseZhyb2csr((a),(b),(c),(cuDoubleComplex*)(d),(e),(f))
158fcdce8c4SStefano 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)
159039c6fbaSStefano Zampini       #define cusparse_csr_spgeam(a,b,c,d,e,f,g,h,i,j,k,l,m,n,o,p,q,r,s)   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)
160aa372e3fSPaul Mullowney     #endif
161aa372e3fSPaul Mullowney   #else
162aa372e3fSPaul Mullowney     #if defined(PETSC_USE_REAL_SINGLE)
163aa372e3fSPaul Mullowney       #define cusparse_csr_spmv cusparseScsrmv
164ccdfe979SStefano Zampini       #define cusparse_csr_spmm cusparseScsrmm
165aa372e3fSPaul Mullowney       #define cusparse_csr2csc  cusparseScsr2csc
166aa372e3fSPaul Mullowney       #define cusparse_hyb_spmv cusparseShybmv
167aa372e3fSPaul Mullowney       #define cusparse_csr2hyb  cusparseScsr2hyb
168aa372e3fSPaul Mullowney       #define cusparse_hyb2csr  cusparseShyb2csr
169fcdce8c4SStefano Zampini       #define cusparse_csr_spgemm cusparseScsrgemm
170039c6fbaSStefano Zampini       #define cusparse_csr_spgeam cusparseScsrgeam
171aa372e3fSPaul Mullowney     #elif defined(PETSC_USE_REAL_DOUBLE)
172aa372e3fSPaul Mullowney       #define cusparse_csr_spmv cusparseDcsrmv
173ccdfe979SStefano Zampini       #define cusparse_csr_spmm cusparseDcsrmm
174aa372e3fSPaul Mullowney       #define cusparse_csr2csc  cusparseDcsr2csc
175aa372e3fSPaul Mullowney       #define cusparse_hyb_spmv cusparseDhybmv
176aa372e3fSPaul Mullowney       #define cusparse_csr2hyb  cusparseDcsr2hyb
177aa372e3fSPaul Mullowney       #define cusparse_hyb2csr  cusparseDhyb2csr
178fcdce8c4SStefano Zampini       #define cusparse_csr_spgemm cusparseDcsrgemm
179039c6fbaSStefano Zampini       #define cusparse_csr_spgeam cusparseDcsrgeam
180aa372e3fSPaul Mullowney     #endif
181aa372e3fSPaul Mullowney   #endif
182afb2bd1cSJunchao Zhang #endif
183aa372e3fSPaul Mullowney 
184aa372e3fSPaul Mullowney #define THRUSTINTARRAY32 thrust::device_vector<int>
185aa372e3fSPaul Mullowney #define THRUSTINTARRAY thrust::device_vector<PetscInt>
186aa372e3fSPaul Mullowney #define THRUSTARRAY thrust::device_vector<PetscScalar>
187aa372e3fSPaul Mullowney 
188aa372e3fSPaul Mullowney /* A CSR matrix structure */
189aa372e3fSPaul Mullowney struct CsrMatrix {
190aa372e3fSPaul Mullowney   PetscInt         num_rows;
191aa372e3fSPaul Mullowney   PetscInt         num_cols;
192aa372e3fSPaul Mullowney   PetscInt         num_entries;
193aa372e3fSPaul Mullowney   THRUSTINTARRAY32 *row_offsets;
194aa372e3fSPaul Mullowney   THRUSTINTARRAY32 *column_indices;
195aa372e3fSPaul Mullowney   THRUSTARRAY      *values;
1969ae82921SPaul Mullowney };
1979ae82921SPaul Mullowney 
198aa372e3fSPaul Mullowney /* This is struct holding the relevant data needed to a MatSolve */
199aa372e3fSPaul Mullowney struct Mat_SeqAIJCUSPARSETriFactorStruct {
200aa372e3fSPaul Mullowney   /* Data needed for triangular solve */
201aa372e3fSPaul Mullowney   cusparseMatDescr_t          descr;
202aa372e3fSPaul Mullowney   cusparseOperation_t         solveOp;
203aa372e3fSPaul Mullowney   CsrMatrix                   *csrMat;
204261a78b4SJunchao Zhang   csrsvInfo_t                 solveInfo;
2051b0a6780SStefano Zampini   cusparseSolvePolicy_t       solvePolicy;     /* whether level information is generated and used */
2061b0a6780SStefano Zampini   int                         solveBufferSize;
2071b0a6780SStefano Zampini   void                        *solveBuffer;
2081b0a6780SStefano Zampini   size_t                      csr2cscBufferSize; /* to transpose the triangular factor (only used for CUDA >= 11.0) */
2091b0a6780SStefano Zampini   void                        *csr2cscBuffer;
2102cbc15d9SMark   PetscScalar                 *AA_h; /* managed host buffer for moving values to the GPU */
211aa372e3fSPaul Mullowney };
212aa372e3fSPaul Mullowney 
213afb2bd1cSJunchao Zhang /* This is a larger struct holding all the triangular factors for a solve, transpose solve, and any indices used in a reordering */
214aa372e3fSPaul Mullowney struct Mat_SeqAIJCUSPARSETriFactors {
215aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorPtr; /* pointer for lower triangular (factored matrix) on GPU */
216aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorPtr; /* pointer for upper triangular (factored matrix) on GPU */
217aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorPtrTranspose; /* pointer for lower triangular (factored matrix) on GPU for the transpose (useful for BiCG) */
218aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorPtrTranspose; /* pointer for upper triangular (factored matrix) on GPU for the transpose (useful for BiCG)*/
219aa372e3fSPaul Mullowney   THRUSTINTARRAY                    *rpermIndices;  /* indices used for any reordering */
220aa372e3fSPaul Mullowney   THRUSTINTARRAY                    *cpermIndices;  /* indices used for any reordering */
221aa372e3fSPaul Mullowney   THRUSTARRAY                       *workVector;
222aa372e3fSPaul Mullowney   cusparseHandle_t                  handle;   /* a handle to the cusparse library */
223aa372e3fSPaul Mullowney   PetscInt                          nnz;      /* number of nonzeros ... need this for accurate logging between ICC and ILU */
224bddcd29dSMark Adams   PetscScalar                       *a_band_d; /* GPU data for banded CSR LU factorization matrix diag(L)=1 */
225bddcd29dSMark Adams   int                               *i_band_d; /* this could be optimized away */
226e8d2b73aSMark Adams   cudaDeviceProp                    dev_prop;
227e8d2b73aSMark Adams   PetscBool                         init_dev_prop;
228da112707SJunchao Zhang 
229da112707SJunchao Zhang   /* csrilu0/csric0 appeared in cusparse-8.0, but we use it along with cusparseSpSV,
230da112707SJunchao Zhang      which first appeared in cusparse-11.5 with cuda-11.3.
231da112707SJunchao Zhang   */
232da112707SJunchao Zhang  #if CUSPARSE_VERSION >= 11500
233da112707SJunchao Zhang   PetscScalar           *csrVal;
234da112707SJunchao Zhang   int                   *csrRowPtr,*csrColIdx; /* a,i,j of M. Using int since some cusparse APIs only support 32-bit indices */
235da112707SJunchao Zhang 
236da112707SJunchao Zhang   /* Mixed mat descriptor types? yes, different cusparse APIs use different types */
237da112707SJunchao Zhang   cusparseMatDescr_t    matDescr_M;
238da112707SJunchao Zhang   cusparseSpMatDescr_t  spMatDescr_L,spMatDescr_U;
239da112707SJunchao Zhang   cusparseSpSVDescr_t   spsvDescr_L,spsvDescr_Lt,spsvDescr_U,spsvDescr_Ut;
240da112707SJunchao Zhang 
241da112707SJunchao Zhang   cusparseDnVecDescr_t  dnVecDescr_X,dnVecDescr_Y;
242da112707SJunchao Zhang   PetscScalar           *X,*Y; /* data array of dnVec X and Y */
243da112707SJunchao Zhang 
244da112707SJunchao Zhang   /* Mixed size types? yes, CUDA-11.7.0 declared cusparseDcsrilu02_bufferSizeExt() that returns size_t but did not implement it! */
245*12ba2bc6SJunchao Zhang   int                   factBufferSize_M; /* M ~= LU or LLt */
246da112707SJunchao Zhang   size_t                spsvBufferSize_L,spsvBufferSize_Lt,spsvBufferSize_U,spsvBufferSize_Ut;
247*12ba2bc6SJunchao Zhang   /* cusparse needs various buffers for factorization and solve of L, U, Lt, or Ut.
248*12ba2bc6SJunchao Zhang      So save memory, we share the factorization buffer with one of spsvBuffer_L/U.
249*12ba2bc6SJunchao Zhang   */
250*12ba2bc6SJunchao Zhang   void                  *factBuffer_M,*spsvBuffer_L,*spsvBuffer_U,*spsvBuffer_Lt,*spsvBuffer_Ut;
251da112707SJunchao Zhang 
252da112707SJunchao Zhang   csrilu02Info_t        ilu0Info_M;
253da112707SJunchao Zhang   csric02Info_t         ic0Info_M;
254da112707SJunchao Zhang   int                   structural_zero,numerical_zero;
255da112707SJunchao Zhang   cusparseSolvePolicy_t policy_M;
256da112707SJunchao Zhang 
257*12ba2bc6SJunchao Zhang   /* In MatSolveTranspose() for ILU0, we use the two flags to do on-demand solve */
258*12ba2bc6SJunchao Zhang   PetscBool             createdTransposeSpSVDescr; /* Have we created SpSV descriptors for Lt, Ut? */
259*12ba2bc6SJunchao Zhang   PetscBool             updatedTransposeSpSVAnalysis; /* Have we updated SpSV analysis with the latest L, U values? */
260da112707SJunchao Zhang 
261da112707SJunchao Zhang   PetscLogDouble        numericFactFlops; /* Estimated FLOPs in ILU0/ICC0 numeric factorization */
262da112707SJunchao Zhang  #endif
263aa372e3fSPaul Mullowney };
264aa372e3fSPaul Mullowney 
265afb2bd1cSJunchao Zhang struct Mat_CusparseSpMV {
266afb2bd1cSJunchao Zhang   PetscBool             initialized;    /* Don't rely on spmvBuffer != NULL to test if the struct is initialized, */
267afb2bd1cSJunchao Zhang   size_t                spmvBufferSize; /* since I'm not sure if smvBuffer can be NULL even after cusparseSpMV_bufferSize() */
268afb2bd1cSJunchao Zhang   void                  *spmvBuffer;
2699db3cbf9SStefano 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 */
270afb2bd1cSJunchao Zhang   cusparseDnVecDescr_t  vecXDescr,vecYDescr; /* descriptor for the dense vectors in y=op(A)x */
2719db3cbf9SStefano Zampini  #endif
272afb2bd1cSJunchao Zhang };
273afb2bd1cSJunchao Zhang 
274afb2bd1cSJunchao Zhang /* This is struct holding the relevant data needed to a MatMult */
275afb2bd1cSJunchao Zhang struct Mat_SeqAIJCUSPARSEMultStruct {
276afb2bd1cSJunchao Zhang   void               *mat;  /* opaque pointer to a matrix. This could be either a cusparseHybMat_t or a CsrMatrix */
277afb2bd1cSJunchao Zhang   cusparseMatDescr_t descr; /* Data needed to describe the matrix for a multiply */
278afb2bd1cSJunchao Zhang   THRUSTINTARRAY     *cprowIndices;   /* compressed row indices used in the parallel SpMV */
279afb2bd1cSJunchao Zhang   PetscScalar        *alpha_one; /* pointer to a device "scalar" storing the alpha parameter in the SpMV */
280afb2bd1cSJunchao Zhang   PetscScalar        *beta_zero; /* pointer to a device "scalar" storing the beta parameter in the SpMV as zero*/
281afb2bd1cSJunchao Zhang   PetscScalar        *beta_one; /* pointer to a device "scalar" storing the beta parameter in the SpMV as one */
282afb2bd1cSJunchao Zhang  #if PETSC_PKG_CUDA_VERSION_GE(11,0,0)
283afb2bd1cSJunchao Zhang   cusparseSpMatDescr_t  matDescr;  /* descriptor for the matrix, used by SpMV and SpMM */
284afb2bd1cSJunchao Zhang   Mat_CusparseSpMV      cuSpMV[3]; /* different Mat_CusparseSpMV structs for non-transpose, transpose, conj-transpose */
285afb2bd1cSJunchao Zhang   Mat_SeqAIJCUSPARSEMultStruct() : matDescr(NULL) {
286afb2bd1cSJunchao Zhang     for (int i=0; i<3; i++) cuSpMV[i].initialized = PETSC_FALSE;
287afb2bd1cSJunchao Zhang   }
288afb2bd1cSJunchao Zhang  #endif
289afb2bd1cSJunchao Zhang };
290afb2bd1cSJunchao Zhang 
2917301b172SPierre Jolivet /* This is a larger struct holding all the matrices for a SpMV, and SpMV Transpose */
2929ae82921SPaul Mullowney struct Mat_SeqAIJCUSPARSE {
293aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSEMultStruct *mat;            /* pointer to the matrix on the GPU */
294aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSEMultStruct *matTranspose;   /* pointer to the matrix on the GPU (for the transpose ... useful for BiCG) */
295aa372e3fSPaul Mullowney   THRUSTARRAY                  *workVector;     /* pointer to a workvector to which we can copy the relevant indices of a vector we want to multiply */
296029808eeSJunchao Zhang   THRUSTINTARRAY32             *rowoffsets_gpu; /* rowoffsets on GPU in non-compressed-row format. It is used to convert CSR to CSC */
297213423ffSJunchao Zhang   PetscInt                     nrows;           /* number of rows of the matrix seen by GPU */
298e057df02SPaul Mullowney   MatCUSPARSEStorageFormat     format;          /* the storage format for the matrix on the device */
299365b711fSMark Adams   PetscBool                    use_cpu_solve;   /* Use AIJ_Seq (I)LU solve */
300aa372e3fSPaul Mullowney   cudaStream_t                 stream;          /* a stream for the parallel SpMV ... this is not owned and should not be deleted */
301aa372e3fSPaul Mullowney   cusparseHandle_t             handle;          /* a handle to the cusparse library ... this may not be owned (if we're working in parallel i.e. multiGPUs) */
30254da937aSStefano Zampini   PetscObjectState             nonzerostate;    /* track nonzero state to possibly recreate the GPU matrix */
303afb2bd1cSJunchao Zhang  #if PETSC_PKG_CUDA_VERSION_GE(11,0,0)
304afb2bd1cSJunchao Zhang   size_t                       csr2cscBufferSize; /* stuff used to compute the matTranspose above */
305afb2bd1cSJunchao Zhang   void                         *csr2cscBuffer;    /* This is used as a C struct and is calloc'ed by PetscNewLog() */
306afb2bd1cSJunchao Zhang   cusparseCsr2CscAlg_t         csr2cscAlg;        /* algorithms can be selected from command line options */
307afb2bd1cSJunchao Zhang   cusparseSpMVAlg_t            spmvAlg;
308afb2bd1cSJunchao Zhang   cusparseSpMMAlg_t            spmmAlg;
309afb2bd1cSJunchao Zhang  #endif
310a49f1ed0SStefano Zampini   THRUSTINTARRAY               *csr2csc_i;
311042217e8SBarry Smith   PetscSplitCSRDataStructure   deviceMat;       /* Matrix on device for, eg, assembly */
312219fbbafSJunchao Zhang 
313219fbbafSJunchao Zhang   /* Stuff for basic COO support */
314ddea5d60SJunchao Zhang   THRUSTINTARRAY               *cooPerm;        /* permutation array that sorts the input coo entris by row and col */
315ddea5d60SJunchao Zhang   THRUSTINTARRAY               *cooPerm_a;      /* ordered array that indicate i-th nonzero (after sorting) is the j-th unique nonzero */
316219fbbafSJunchao Zhang 
317219fbbafSJunchao Zhang   /* Stuff for extended COO support */
318219fbbafSJunchao Zhang   PetscBool                    use_extended_coo; /* Use extended COO format */
319219fbbafSJunchao 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 */
320219fbbafSJunchao Zhang   PetscCount                   *perm_d;
321219fbbafSJunchao Zhang 
322219fbbafSJunchao Zhang   Mat_SeqAIJCUSPARSE() : use_extended_coo(PETSC_FALSE), perm_d(NULL), jmap_d(NULL) {}
3239ae82921SPaul Mullowney };
3249ae82921SPaul Mullowney 
325d3ecb6c7SJed Brown PETSC_INTERN PetscErrorCode MatSeqAIJCUSPARSECopyToGPU(Mat);
326219fbbafSJunchao Zhang PETSC_INTERN PetscErrorCode MatSetPreallocationCOO_SeqAIJCUSPARSE_Basic(Mat,PetscCount,const PetscInt[],const PetscInt[]);
327219fbbafSJunchao Zhang PETSC_INTERN PetscErrorCode MatSetValuesCOO_SeqAIJCUSPARSE_Basic(Mat,const PetscScalar[],InsertMode);
328ed502f03SStefano Zampini PETSC_INTERN PetscErrorCode MatSeqAIJCUSPARSEMergeMats(Mat,Mat,MatReuse,Mat*);
3295f101d05SStefano Zampini PETSC_INTERN PetscErrorCode MatSeqAIJCUSPARSETriFactors_Reset(Mat_SeqAIJCUSPARSETriFactors_p*);
330ed502f03SStefano Zampini 
3319fbee547SJacob Faibussowitsch static inline bool isCudaMem(const void *data)
33208391a17SStefano Zampini {
33308391a17SStefano Zampini   cudaError_t                  cerr;
33408391a17SStefano Zampini   struct cudaPointerAttributes attr;
33508391a17SStefano Zampini   enum cudaMemoryType          mtype;
33608391a17SStefano Zampini   cerr = cudaPointerGetAttributes(&attr,data); /* Do not check error since before CUDA 11.0, passing a host pointer returns cudaErrorInvalidValue */
33708391a17SStefano Zampini   cudaGetLastError(); /* Reset the last error */
33808391a17SStefano Zampini   #if (CUDART_VERSION < 10000)
33908391a17SStefano Zampini     mtype = attr.memoryType;
34008391a17SStefano Zampini   #else
34108391a17SStefano Zampini     mtype = attr.type;
34208391a17SStefano Zampini   #endif
34308391a17SStefano Zampini   if (cerr == cudaSuccess && mtype == cudaMemoryTypeDevice) return true;
34408391a17SStefano Zampini   else return false;
34508391a17SStefano Zampini }
34608391a17SStefano Zampini 
3479ae82921SPaul Mullowney #endif
348