xref: /petsc/src/mat/impls/aij/seq/seqcusparse/cusparsematimpl.h (revision 82a78a4ef1c3dde9953e002d5a85008393775538)
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) {                                      \
257eaca502SStefano Zampini       SETERRQ1(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};
33afb2bd1cSJunchao Zhang   #elif defined(PETSC_USE_REAL_DOUBLE)
34afb2bd1cSJunchao Zhang     const cuDoubleComplex PETSC_CUSPARSE_ONE  = {1.0, 0.0};
35afb2bd1cSJunchao Zhang     const cuDoubleComplex PETSC_CUSPARSE_ZERO = {0.0, 0.0};
36afb2bd1cSJunchao Zhang   #endif
37afb2bd1cSJunchao Zhang #else
38afb2bd1cSJunchao Zhang   const PetscScalar PETSC_CUSPARSE_ONE        = 1.0;
39afb2bd1cSJunchao Zhang   const PetscScalar PETSC_CUSPARSE_ZERO       = 0.0;
40afb2bd1cSJunchao Zhang #endif
41afb2bd1cSJunchao Zhang 
421b0a6780SStefano Zampini #if PETSC_PKG_CUDA_VERSION_GE(9,0,0)
43afb2bd1cSJunchao Zhang   #define cusparse_create_analysis_info  cusparseCreateCsrsv2Info
44afb2bd1cSJunchao Zhang   #define cusparse_destroy_analysis_info cusparseDestroyCsrsv2Info
45afb2bd1cSJunchao Zhang   #if defined(PETSC_USE_COMPLEX)
46afb2bd1cSJunchao Zhang     #if defined(PETSC_USE_REAL_SINGLE)
47afb2bd1cSJunchao Zhang       #define cusparse_get_svbuffsize(a,b,c,d,e,f,g,h,i,j) cusparseCcsrsv2_bufferSize(a,b,c,d,e,(cuComplex*)(f),g,h,i,j)
48afb2bd1cSJunchao Zhang       #define cusparse_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)
49afb2bd1cSJunchao Zhang       #define cusparse_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)
50afb2bd1cSJunchao Zhang     #elif defined(PETSC_USE_REAL_DOUBLE)
51afb2bd1cSJunchao Zhang       #define cusparse_get_svbuffsize(a,b,c,d,e,f,g,h,i,j) cusparseZcsrsv2_bufferSize(a,b,c,d,e,(cuDoubleComplex*)(f),g,h,i,j)
52afb2bd1cSJunchao Zhang       #define cusparse_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)
53afb2bd1cSJunchao Zhang       #define cusparse_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)
54afb2bd1cSJunchao Zhang     #endif
55afb2bd1cSJunchao Zhang   #else /* not complex */
56afb2bd1cSJunchao Zhang     #if defined(PETSC_USE_REAL_SINGLE)
57afb2bd1cSJunchao Zhang       #define cusparse_get_svbuffsize cusparseScsrsv2_bufferSize
58afb2bd1cSJunchao Zhang       #define cusparse_analysis       cusparseScsrsv2_analysis
59afb2bd1cSJunchao Zhang       #define cusparse_solve          cusparseScsrsv2_solve
60afb2bd1cSJunchao Zhang     #elif defined(PETSC_USE_REAL_DOUBLE)
61afb2bd1cSJunchao Zhang       #define cusparse_get_svbuffsize cusparseDcsrsv2_bufferSize
62afb2bd1cSJunchao Zhang       #define cusparse_analysis       cusparseDcsrsv2_analysis
63afb2bd1cSJunchao Zhang       #define cusparse_solve          cusparseDcsrsv2_solve
64afb2bd1cSJunchao Zhang     #endif
65afb2bd1cSJunchao Zhang   #endif
66afb2bd1cSJunchao Zhang #else
67afb2bd1cSJunchao Zhang   #define cusparse_create_analysis_info  cusparseCreateSolveAnalysisInfo
68afb2bd1cSJunchao Zhang   #define cusparse_destroy_analysis_info cusparseDestroySolveAnalysisInfo
69afb2bd1cSJunchao Zhang   #if defined(PETSC_USE_COMPLEX)
70afb2bd1cSJunchao Zhang     #if defined(PETSC_USE_REAL_SINGLE)
71ec42abe4SAlejandro Lamas Daviña       #define cusparse_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))
72ec42abe4SAlejandro Lamas Daviña       #define cusparse_analysis(a,b,c,d,e,f,g,h,i)  cusparseCcsrsv_analysis((a),(b),(c),(d),(e),(cuComplex*)(f),(g),(h),(i))
731b0a6780SStefano Zampini     #elif defined(PETSC_USE_REAL_DOUBLE)
741b0a6780SStefano Zampini       #define cusparse_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))
751b0a6780SStefano Zampini       #define cusparse_analysis(a,b,c,d,e,f,g,h,i)  cusparseZcsrsv_analysis((a),(b),(c),(d),(e),(cuDoubleComplex*)(f),(g),(h),(i))
761b0a6780SStefano Zampini     #endif
771b0a6780SStefano Zampini   #else /* not complex */
781b0a6780SStefano Zampini     #if defined(PETSC_USE_REAL_SINGLE)
791b0a6780SStefano Zampini       #define cusparse_solve    cusparseScsrsv_solve
801b0a6780SStefano Zampini       #define cusparse_analysis cusparseScsrsv_analysis
811b0a6780SStefano Zampini     #elif defined(PETSC_USE_REAL_DOUBLE)
821b0a6780SStefano Zampini       #define cusparse_solve    cusparseDcsrsv_solve
831b0a6780SStefano Zampini       #define cusparse_analysis cusparseDcsrsv_analysis
841b0a6780SStefano Zampini     #endif
851b0a6780SStefano Zampini   #endif
861b0a6780SStefano Zampini #endif
871b0a6780SStefano Zampini 
881b0a6780SStefano Zampini #if PETSC_PKG_CUDA_VERSION_GE(11,0,0)
891b0a6780SStefano Zampini   #define cusparse_csr2csc cusparseCsr2cscEx2
901b0a6780SStefano Zampini   #if defined(PETSC_USE_COMPLEX)
911b0a6780SStefano Zampini     #if defined(PETSC_USE_REAL_SINGLE)
921b0a6780SStefano Zampini       #define cusparse_scalartype CUDA_C_32F
93039c6fbaSStefano 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)
94039c6fbaSStefano 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)
951b0a6780SStefano Zampini     #elif defined(PETSC_USE_REAL_DOUBLE)
961b0a6780SStefano Zampini       #define cusparse_scalartype CUDA_C_64F
97039c6fbaSStefano 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)
98039c6fbaSStefano 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)
991b0a6780SStefano Zampini     #endif
1001b0a6780SStefano Zampini   #else /* not complex */
1011b0a6780SStefano Zampini     #if defined(PETSC_USE_REAL_SINGLE)
1021b0a6780SStefano Zampini       #define cusparse_scalartype CUDA_R_32F
103039c6fbaSStefano Zampini       #define cusparse_csr_spgeam            cusparseScsrgeam2
104039c6fbaSStefano Zampini       #define cusparse_csr_spgeam_bufferSize cusparseScsrgeam2_bufferSizeExt
1051b0a6780SStefano Zampini     #elif defined(PETSC_USE_REAL_DOUBLE)
1061b0a6780SStefano Zampini       #define cusparse_scalartype CUDA_R_64F
107039c6fbaSStefano Zampini       #define cusparse_csr_spgeam            cusparseDcsrgeam2
108039c6fbaSStefano Zampini       #define cusparse_csr_spgeam_bufferSize cusparseDcsrgeam2_bufferSizeExt
1091b0a6780SStefano Zampini     #endif
1101b0a6780SStefano Zampini   #endif
1111b0a6780SStefano Zampini #else
1121b0a6780SStefano Zampini   #if defined(PETSC_USE_COMPLEX)
1131b0a6780SStefano Zampini     #if defined(PETSC_USE_REAL_SINGLE)
114ec42abe4SAlejandro 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))
115ccdfe979SStefano 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))
116ec42abe4SAlejandro 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))
117ec42abe4SAlejandro 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))
118ec42abe4SAlejandro 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))
119ec42abe4SAlejandro Lamas Daviña       #define cusparse_hyb2csr(a,b,c,d,e,f)                      cusparseChyb2csr((a),(b),(c),(cuComplex*)(d),(e),(f))
120fcdce8c4SStefano 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)
121039c6fbaSStefano 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)
122aa372e3fSPaul Mullowney     #elif defined(PETSC_USE_REAL_DOUBLE)
123ec42abe4SAlejandro 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))
124ccdfe979SStefano 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))
125ec42abe4SAlejandro 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))
126ec42abe4SAlejandro 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))
127ec42abe4SAlejandro 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))
128ec42abe4SAlejandro Lamas Daviña       #define cusparse_hyb2csr(a,b,c,d,e,f)                      cusparseZhyb2csr((a),(b),(c),(cuDoubleComplex*)(d),(e),(f))
129fcdce8c4SStefano 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)
130039c6fbaSStefano 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)
131aa372e3fSPaul Mullowney     #endif
132aa372e3fSPaul Mullowney   #else
133aa372e3fSPaul Mullowney     #if defined(PETSC_USE_REAL_SINGLE)
134aa372e3fSPaul Mullowney       #define cusparse_csr_spmv cusparseScsrmv
135ccdfe979SStefano Zampini       #define cusparse_csr_spmm cusparseScsrmm
136aa372e3fSPaul Mullowney       #define cusparse_csr2csc  cusparseScsr2csc
137aa372e3fSPaul Mullowney       #define cusparse_hyb_spmv cusparseShybmv
138aa372e3fSPaul Mullowney       #define cusparse_csr2hyb  cusparseScsr2hyb
139aa372e3fSPaul Mullowney       #define cusparse_hyb2csr  cusparseShyb2csr
140fcdce8c4SStefano Zampini       #define cusparse_csr_spgemm cusparseScsrgemm
141039c6fbaSStefano Zampini       #define cusparse_csr_spgeam cusparseScsrgeam
142aa372e3fSPaul Mullowney     #elif defined(PETSC_USE_REAL_DOUBLE)
143aa372e3fSPaul Mullowney       #define cusparse_csr_spmv cusparseDcsrmv
144ccdfe979SStefano Zampini       #define cusparse_csr_spmm cusparseDcsrmm
145aa372e3fSPaul Mullowney       #define cusparse_csr2csc  cusparseDcsr2csc
146aa372e3fSPaul Mullowney       #define cusparse_hyb_spmv cusparseDhybmv
147aa372e3fSPaul Mullowney       #define cusparse_csr2hyb  cusparseDcsr2hyb
148aa372e3fSPaul Mullowney       #define cusparse_hyb2csr  cusparseDhyb2csr
149fcdce8c4SStefano Zampini       #define cusparse_csr_spgemm cusparseDcsrgemm
150039c6fbaSStefano Zampini       #define cusparse_csr_spgeam cusparseDcsrgeam
151aa372e3fSPaul Mullowney     #endif
152aa372e3fSPaul Mullowney   #endif
153afb2bd1cSJunchao Zhang #endif
154aa372e3fSPaul Mullowney 
155aa372e3fSPaul Mullowney #define THRUSTINTARRAY32 thrust::device_vector<int>
156aa372e3fSPaul Mullowney #define THRUSTINTARRAY thrust::device_vector<PetscInt>
157aa372e3fSPaul Mullowney #define THRUSTARRAY thrust::device_vector<PetscScalar>
158aa372e3fSPaul Mullowney 
159aa372e3fSPaul Mullowney /* A CSR matrix structure */
160aa372e3fSPaul Mullowney struct CsrMatrix {
161aa372e3fSPaul Mullowney   PetscInt         num_rows;
162aa372e3fSPaul Mullowney   PetscInt         num_cols;
163aa372e3fSPaul Mullowney   PetscInt         num_entries;
164aa372e3fSPaul Mullowney   THRUSTINTARRAY32 *row_offsets;
165aa372e3fSPaul Mullowney   THRUSTINTARRAY32 *column_indices;
166aa372e3fSPaul Mullowney   THRUSTARRAY      *values;
1679ae82921SPaul Mullowney };
1689ae82921SPaul Mullowney 
169aa372e3fSPaul Mullowney /* This is struct holding the relevant data needed to a MatSolve */
170aa372e3fSPaul Mullowney struct Mat_SeqAIJCUSPARSETriFactorStruct {
171aa372e3fSPaul Mullowney   /* Data needed for triangular solve */
172aa372e3fSPaul Mullowney   cusparseMatDescr_t          descr;
173aa372e3fSPaul Mullowney   cusparseOperation_t         solveOp;
174aa372e3fSPaul Mullowney   CsrMatrix                   *csrMat;
1751b0a6780SStefano Zampini  #if PETSC_PKG_CUDA_VERSION_GE(9,0,0)
176afb2bd1cSJunchao Zhang   csrsv2Info_t                solveInfo;
177afb2bd1cSJunchao Zhang  #else
178afb2bd1cSJunchao Zhang   cusparseSolveAnalysisInfo_t solveInfo;
179afb2bd1cSJunchao Zhang  #endif
1801b0a6780SStefano Zampini   cusparseSolvePolicy_t       solvePolicy;     /* whether level information is generated and used */
1811b0a6780SStefano Zampini   int                         solveBufferSize;
1821b0a6780SStefano Zampini   void                        *solveBuffer;
1831b0a6780SStefano Zampini   size_t                      csr2cscBufferSize; /* to transpose the triangular factor (only used for CUDA >= 11.0) */
1841b0a6780SStefano Zampini   void                        *csr2cscBuffer;
1852cbc15d9SMark   PetscScalar                 *AA_h; /* managed host buffer for moving values to the GPU */
186aa372e3fSPaul Mullowney };
187aa372e3fSPaul Mullowney 
188afb2bd1cSJunchao Zhang /* This is a larger struct holding all the triangular factors for a solve, transpose solve, and any indices used in a reordering */
189aa372e3fSPaul Mullowney struct Mat_SeqAIJCUSPARSETriFactors {
190aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorPtr; /* pointer for lower triangular (factored matrix) on GPU */
191aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorPtr; /* pointer for upper triangular (factored matrix) on GPU */
192aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorPtrTranspose; /* pointer for lower triangular (factored matrix) on GPU for the transpose (useful for BiCG) */
193aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorPtrTranspose; /* pointer for upper triangular (factored matrix) on GPU for the transpose (useful for BiCG)*/
194aa372e3fSPaul Mullowney   THRUSTINTARRAY                    *rpermIndices;  /* indices used for any reordering */
195aa372e3fSPaul Mullowney   THRUSTINTARRAY                    *cpermIndices;  /* indices used for any reordering */
196aa372e3fSPaul Mullowney   THRUSTARRAY                       *workVector;
197aa372e3fSPaul Mullowney   cusparseHandle_t                  handle;   /* a handle to the cusparse library */
198aa372e3fSPaul Mullowney   PetscInt                          nnz;      /* number of nonzeros ... need this for accurate logging between ICC and ILU */
199bddcd29dSMark Adams   PetscScalar                       *a_band_d; /* GPU data for banded CSR LU factorization matrix diag(L)=1 */
200bddcd29dSMark Adams   int                               *i_band_d; /* this could be optimized away */
201e8d2b73aSMark Adams   cudaDeviceProp                    dev_prop;
202e8d2b73aSMark Adams   PetscBool                         init_dev_prop;
203aa372e3fSPaul Mullowney };
204aa372e3fSPaul Mullowney 
205afb2bd1cSJunchao Zhang struct Mat_CusparseSpMV {
206afb2bd1cSJunchao Zhang   PetscBool             initialized;    /* Don't rely on spmvBuffer != NULL to test if the struct is initialized, */
207afb2bd1cSJunchao Zhang   size_t                spmvBufferSize; /* since I'm not sure if smvBuffer can be NULL even after cusparseSpMV_bufferSize() */
208afb2bd1cSJunchao Zhang   void                  *spmvBuffer;
2099db3cbf9SStefano 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 */
210afb2bd1cSJunchao Zhang   cusparseDnVecDescr_t  vecXDescr,vecYDescr; /* descriptor for the dense vectors in y=op(A)x */
2119db3cbf9SStefano Zampini  #endif
212afb2bd1cSJunchao Zhang };
213afb2bd1cSJunchao Zhang 
214afb2bd1cSJunchao Zhang /* This is struct holding the relevant data needed to a MatMult */
215afb2bd1cSJunchao Zhang struct Mat_SeqAIJCUSPARSEMultStruct {
216afb2bd1cSJunchao Zhang   void               *mat;  /* opaque pointer to a matrix. This could be either a cusparseHybMat_t or a CsrMatrix */
217afb2bd1cSJunchao Zhang   cusparseMatDescr_t descr; /* Data needed to describe the matrix for a multiply */
218afb2bd1cSJunchao Zhang   THRUSTINTARRAY     *cprowIndices;   /* compressed row indices used in the parallel SpMV */
219afb2bd1cSJunchao Zhang   PetscScalar        *alpha_one; /* pointer to a device "scalar" storing the alpha parameter in the SpMV */
220afb2bd1cSJunchao Zhang   PetscScalar        *beta_zero; /* pointer to a device "scalar" storing the beta parameter in the SpMV as zero*/
221afb2bd1cSJunchao Zhang   PetscScalar        *beta_one; /* pointer to a device "scalar" storing the beta parameter in the SpMV as one */
222afb2bd1cSJunchao Zhang  #if PETSC_PKG_CUDA_VERSION_GE(11,0,0)
223afb2bd1cSJunchao Zhang   cusparseSpMatDescr_t  matDescr;  /* descriptor for the matrix, used by SpMV and SpMM */
224afb2bd1cSJunchao Zhang   Mat_CusparseSpMV      cuSpMV[3]; /* different Mat_CusparseSpMV structs for non-transpose, transpose, conj-transpose */
225afb2bd1cSJunchao Zhang   Mat_SeqAIJCUSPARSEMultStruct() : matDescr(NULL) {
226afb2bd1cSJunchao Zhang     for (int i=0; i<3; i++) cuSpMV[i].initialized = PETSC_FALSE;
227afb2bd1cSJunchao Zhang   }
228afb2bd1cSJunchao Zhang  #endif
229afb2bd1cSJunchao Zhang };
230afb2bd1cSJunchao Zhang 
231aa372e3fSPaul Mullowney /* This is a larger struct holding all the matrices for a SpMV, and SpMV Tranpose */
2329ae82921SPaul Mullowney struct Mat_SeqAIJCUSPARSE {
233aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSEMultStruct *mat;            /* pointer to the matrix on the GPU */
234aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSEMultStruct *matTranspose;   /* pointer to the matrix on the GPU (for the transpose ... useful for BiCG) */
235aa372e3fSPaul Mullowney   THRUSTARRAY                  *workVector;     /* pointer to a workvector to which we can copy the relevant indices of a vector we want to multiply */
236029808eeSJunchao Zhang   THRUSTINTARRAY32             *rowoffsets_gpu; /* rowoffsets on GPU in non-compressed-row format. It is used to convert CSR to CSC */
237213423ffSJunchao Zhang   PetscInt                     nrows;           /* number of rows of the matrix seen by GPU */
238e057df02SPaul Mullowney   MatCUSPARSEStorageFormat     format;          /* the storage format for the matrix on the device */
239365b711fSMark Adams   PetscBool                    use_cpu_solve;   /* Use AIJ_Seq (I)LU solve */
240aa372e3fSPaul Mullowney   cudaStream_t                 stream;          /* a stream for the parallel SpMV ... this is not owned and should not be deleted */
241aa372e3fSPaul Mullowney   cusparseHandle_t             handle;          /* a handle to the cusparse library ... this may not be owned (if we're working in parallel i.e. multiGPUs) */
24254da937aSStefano Zampini   PetscObjectState             nonzerostate;    /* track nonzero state to possibly recreate the GPU matrix */
243afb2bd1cSJunchao Zhang  #if PETSC_PKG_CUDA_VERSION_GE(11,0,0)
244afb2bd1cSJunchao Zhang   size_t                       csr2cscBufferSize; /* stuff used to compute the matTranspose above */
245afb2bd1cSJunchao Zhang   void                         *csr2cscBuffer;    /* This is used as a C struct and is calloc'ed by PetscNewLog() */
246afb2bd1cSJunchao Zhang   cusparseCsr2CscAlg_t         csr2cscAlg;        /* algorithms can be selected from command line options */
247afb2bd1cSJunchao Zhang   cusparseSpMVAlg_t            spmvAlg;
248afb2bd1cSJunchao Zhang   cusparseSpMMAlg_t            spmmAlg;
249afb2bd1cSJunchao Zhang  #endif
250a49f1ed0SStefano Zampini   THRUSTINTARRAY               *csr2csc_i;
251042217e8SBarry Smith   PetscSplitCSRDataStructure   deviceMat;       /* Matrix on device for, eg, assembly */
252ddea5d60SJunchao Zhang   THRUSTINTARRAY               *cooPerm;        /* permutation array that sorts the input coo entris by row and col */
253ddea5d60SJunchao Zhang   THRUSTINTARRAY               *cooPerm_a;      /* ordered array that indicate i-th nonzero (after sorting) is the j-th unique nonzero */
2549ae82921SPaul Mullowney };
2559ae82921SPaul Mullowney 
2565a576424SJed Brown PETSC_INTERN PetscErrorCode MatCUSPARSECopyToGPU(Mat);
257b06137fdSPaul Mullowney PETSC_INTERN PetscErrorCode MatCUSPARSESetStream(Mat, const cudaStream_t stream);
258b06137fdSPaul Mullowney PETSC_INTERN PetscErrorCode MatCUSPARSESetHandle(Mat, const cusparseHandle_t handle);
259b06137fdSPaul Mullowney PETSC_INTERN PetscErrorCode MatCUSPARSEClearHandle(Mat);
260*82a78a4eSJed Brown PETSC_INTERN PetscErrorCode MatSetPreallocationCOO_SeqAIJCUSPARSE(Mat,PetscCount,const PetscInt[],const PetscInt[]);
261e61fc153SStefano Zampini PETSC_INTERN PetscErrorCode MatSetValuesCOO_SeqAIJCUSPARSE(Mat,const PetscScalar[],InsertMode);
262ed502f03SStefano Zampini PETSC_INTERN PetscErrorCode MatSeqAIJCUSPARSEMergeMats(Mat,Mat,MatReuse,Mat*);
2635f101d05SStefano Zampini PETSC_INTERN PetscErrorCode MatSeqAIJCUSPARSETriFactors_Reset(Mat_SeqAIJCUSPARSETriFactors_p*);
264ed502f03SStefano Zampini 
26508391a17SStefano Zampini PETSC_STATIC_INLINE bool isCudaMem(const void *data)
26608391a17SStefano Zampini {
26708391a17SStefano Zampini   cudaError_t                  cerr;
26808391a17SStefano Zampini   struct cudaPointerAttributes attr;
26908391a17SStefano Zampini   enum cudaMemoryType          mtype;
27008391a17SStefano Zampini   cerr = cudaPointerGetAttributes(&attr,data); /* Do not check error since before CUDA 11.0, passing a host pointer returns cudaErrorInvalidValue */
27108391a17SStefano Zampini   cudaGetLastError(); /* Reset the last error */
27208391a17SStefano Zampini   #if (CUDART_VERSION < 10000)
27308391a17SStefano Zampini     mtype = attr.memoryType;
27408391a17SStefano Zampini   #else
27508391a17SStefano Zampini     mtype = attr.type;
27608391a17SStefano Zampini   #endif
27708391a17SStefano Zampini   if (cerr == cudaSuccess && mtype == cudaMemoryTypeDevice) return true;
27808391a17SStefano Zampini   else return false;
27908391a17SStefano Zampini }
28008391a17SStefano Zampini 
2819ae82921SPaul Mullowney #endif
282