xref: /petsc/src/mat/impls/aij/seq/seqcusparse/cusparsematimpl.h (revision 1b0a6780d9463fdb93faab9d83a4082a63444bef)
1519f805aSKarl Rupp #if !defined(__CUSPARSEMATIMPL)
29ae82921SPaul Mullowney #define __CUSPARSEMATIMPL
39ae82921SPaul Mullowney 
4afb2bd1cSJunchao Zhang #include <petscpkg_version.h>
5303a667bSSatish Balay #include <petsc/private/cudavecimpl.h>
69ae82921SPaul Mullowney 
79ae82921SPaul Mullowney #include <cusparse_v2.h>
89ae82921SPaul Mullowney 
99ae82921SPaul Mullowney #include <algorithm>
109ae82921SPaul Mullowney #include <vector>
119ae82921SPaul Mullowney 
12c41cb2e2SAlejandro Lamas Daviña #include <thrust/device_vector.h>
13c41cb2e2SAlejandro Lamas Daviña #include <thrust/device_ptr.h>
1450ab121bSSatish Balay #include <thrust/device_malloc_allocator.h>
15c41cb2e2SAlejandro Lamas Daviña #include <thrust/transform.h>
16c41cb2e2SAlejandro Lamas Daviña #include <thrust/functional.h>
17554b8892SKarl Rupp #include <thrust/sequence.h>
18c41cb2e2SAlejandro Lamas Daviña 
19185e5899SJunchao Zhang #if (CUSPARSE_VER_MAJOR > 10 || CUSPARSE_VER_MAJOR == 10 && CUSPARSE_VER_MINOR >= 2) /* According to cuda/10.1.168 on OLCF Summit */
2057d48284SJunchao Zhang #define CHKERRCUSPARSE(stat) \
2157d48284SJunchao Zhang do { \
2257d48284SJunchao Zhang    if (PetscUnlikely(stat)) { \
2357d48284SJunchao Zhang       const char *name  = cusparseGetErrorName(stat); \
2457d48284SJunchao Zhang       const char *descr = cusparseGetErrorString(stat); \
25589f383fSBarry Smith       SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_GPU,"cuSPARSE error %d (%s) : %s",(int)stat,name,descr); \
2657d48284SJunchao Zhang    } \
2757d48284SJunchao Zhang } while (0)
286d22fe3dSJunchao Zhang #else
29589f383fSBarry Smith #define CHKERRCUSPARSE(stat) do {if (PetscUnlikely(stat)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_GPU,"cusparse error %d",(int)stat);} while (0)
306d22fe3dSJunchao Zhang #endif
3157d48284SJunchao Zhang 
32aa372e3fSPaul Mullowney #if defined(PETSC_USE_COMPLEX)
33aa372e3fSPaul Mullowney   #if defined(PETSC_USE_REAL_SINGLE)
34ccdfe979SStefano Zampini     const cuComplex PETSC_CUSPARSE_ONE        = {1.0f, 0.0f};
35ccdfe979SStefano Zampini     const cuComplex PETSC_CUSPARSE_ZERO       = {0.0f, 0.0f};
36afb2bd1cSJunchao Zhang   #elif defined(PETSC_USE_REAL_DOUBLE)
37afb2bd1cSJunchao Zhang     const cuDoubleComplex PETSC_CUSPARSE_ONE  = {1.0, 0.0};
38afb2bd1cSJunchao Zhang     const cuDoubleComplex PETSC_CUSPARSE_ZERO = {0.0, 0.0};
39afb2bd1cSJunchao Zhang   #endif
40afb2bd1cSJunchao Zhang #else
41afb2bd1cSJunchao Zhang   const PetscScalar PETSC_CUSPARSE_ONE        = 1.0;
42afb2bd1cSJunchao Zhang   const PetscScalar PETSC_CUSPARSE_ZERO       = 0.0;
43afb2bd1cSJunchao Zhang #endif
44afb2bd1cSJunchao Zhang 
45*1b0a6780SStefano Zampini #if PETSC_PKG_CUDA_VERSION_GE(9,0,0)
46afb2bd1cSJunchao Zhang   #define cusparse_create_analysis_info  cusparseCreateCsrsv2Info
47afb2bd1cSJunchao Zhang   #define cusparse_destroy_analysis_info cusparseDestroyCsrsv2Info
48afb2bd1cSJunchao Zhang   #if defined(PETSC_USE_COMPLEX)
49afb2bd1cSJunchao Zhang     #if defined(PETSC_USE_REAL_SINGLE)
50afb2bd1cSJunchao 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)
51afb2bd1cSJunchao 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)
52afb2bd1cSJunchao 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)
53afb2bd1cSJunchao Zhang     #elif defined(PETSC_USE_REAL_DOUBLE)
54afb2bd1cSJunchao 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)
55afb2bd1cSJunchao 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)
56afb2bd1cSJunchao 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)
57afb2bd1cSJunchao Zhang     #endif
58afb2bd1cSJunchao Zhang   #else /* not complex */
59afb2bd1cSJunchao Zhang     #if defined(PETSC_USE_REAL_SINGLE)
60afb2bd1cSJunchao Zhang       #define cusparse_get_svbuffsize cusparseScsrsv2_bufferSize
61afb2bd1cSJunchao Zhang       #define cusparse_analysis       cusparseScsrsv2_analysis
62afb2bd1cSJunchao Zhang       #define cusparse_solve          cusparseScsrsv2_solve
63afb2bd1cSJunchao Zhang     #elif defined(PETSC_USE_REAL_DOUBLE)
64afb2bd1cSJunchao Zhang       #define cusparse_get_svbuffsize cusparseDcsrsv2_bufferSize
65afb2bd1cSJunchao Zhang       #define cusparse_analysis       cusparseDcsrsv2_analysis
66afb2bd1cSJunchao Zhang       #define cusparse_solve          cusparseDcsrsv2_solve
67afb2bd1cSJunchao Zhang     #endif
68afb2bd1cSJunchao Zhang   #endif
69afb2bd1cSJunchao Zhang #else
70afb2bd1cSJunchao Zhang   #define cusparse_create_analysis_info  cusparseCreateSolveAnalysisInfo
71afb2bd1cSJunchao Zhang   #define cusparse_destroy_analysis_info cusparseDestroySolveAnalysisInfo
72afb2bd1cSJunchao Zhang   #if defined(PETSC_USE_COMPLEX)
73afb2bd1cSJunchao Zhang     #if defined(PETSC_USE_REAL_SINGLE)
74ec42abe4SAlejandro 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))
75ec42abe4SAlejandro 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))
76*1b0a6780SStefano Zampini     #elif defined(PETSC_USE_REAL_DOUBLE)
77*1b0a6780SStefano 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))
78*1b0a6780SStefano 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))
79*1b0a6780SStefano Zampini     #endif
80*1b0a6780SStefano Zampini   #else /* not complex */
81*1b0a6780SStefano Zampini     #if defined(PETSC_USE_REAL_SINGLE)
82*1b0a6780SStefano Zampini       #define cusparse_solve    cusparseScsrsv_solve
83*1b0a6780SStefano Zampini       #define cusparse_analysis cusparseScsrsv_analysis
84*1b0a6780SStefano Zampini     #elif defined(PETSC_USE_REAL_DOUBLE)
85*1b0a6780SStefano Zampini       #define cusparse_solve    cusparseDcsrsv_solve
86*1b0a6780SStefano Zampini       #define cusparse_analysis cusparseDcsrsv_analysis
87*1b0a6780SStefano Zampini     #endif
88*1b0a6780SStefano Zampini   #endif
89*1b0a6780SStefano Zampini #endif
90*1b0a6780SStefano Zampini 
91*1b0a6780SStefano Zampini #if PETSC_PKG_CUDA_VERSION_GE(11,0,0)
92*1b0a6780SStefano Zampini   #define cusparse_csr2csc cusparseCsr2cscEx2
93*1b0a6780SStefano Zampini   #if defined(PETSC_USE_COMPLEX)
94*1b0a6780SStefano Zampini     #if defined(PETSC_USE_REAL_SINGLE)
95*1b0a6780SStefano Zampini       #define cusparse_scalartype CUDA_C_32F
96*1b0a6780SStefano Zampini     #elif defined(PETSC_USE_REAL_DOUBLE)
97*1b0a6780SStefano Zampini       #define cusparse_scalartype CUDA_C_64F
98*1b0a6780SStefano Zampini     #endif
99*1b0a6780SStefano Zampini   #else /* not complex */
100*1b0a6780SStefano Zampini     #if defined(PETSC_USE_REAL_SINGLE)
101*1b0a6780SStefano Zampini       #define cusparse_scalartype CUDA_R_32F
102*1b0a6780SStefano Zampini     #elif defined(PETSC_USE_REAL_DOUBLE)
103*1b0a6780SStefano Zampini       #define cusparse_scalartype CUDA_R_64F
104*1b0a6780SStefano Zampini     #endif
105*1b0a6780SStefano Zampini   #endif
106*1b0a6780SStefano Zampini #else
107*1b0a6780SStefano Zampini   #if defined(PETSC_USE_COMPLEX)
108*1b0a6780SStefano Zampini     #if defined(PETSC_USE_REAL_SINGLE)
109ec42abe4SAlejandro 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))
110ccdfe979SStefano 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))
111ec42abe4SAlejandro 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))
112ec42abe4SAlejandro 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))
113ec42abe4SAlejandro 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))
114ec42abe4SAlejandro Lamas Daviña       #define cusparse_hyb2csr(a,b,c,d,e,f)                      cusparseChyb2csr((a),(b),(c),(cuComplex*)(d),(e),(f))
115aa372e3fSPaul Mullowney     #elif defined(PETSC_USE_REAL_DOUBLE)
116ec42abe4SAlejandro 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))
117ccdfe979SStefano 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))
118ec42abe4SAlejandro 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))
119ec42abe4SAlejandro 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))
120ec42abe4SAlejandro 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))
121ec42abe4SAlejandro Lamas Daviña       #define cusparse_hyb2csr(a,b,c,d,e,f)                      cusparseZhyb2csr((a),(b),(c),(cuDoubleComplex*)(d),(e),(f))
122aa372e3fSPaul Mullowney     #endif
123aa372e3fSPaul Mullowney   #else
124aa372e3fSPaul Mullowney     #if defined(PETSC_USE_REAL_SINGLE)
125aa372e3fSPaul Mullowney       #define cusparse_csr_spmv cusparseScsrmv
126ccdfe979SStefano Zampini       #define cusparse_csr_spmm cusparseScsrmm
127aa372e3fSPaul Mullowney       #define cusparse_csr2csc  cusparseScsr2csc
128aa372e3fSPaul Mullowney       #define cusparse_hyb_spmv cusparseShybmv
129aa372e3fSPaul Mullowney       #define cusparse_csr2hyb  cusparseScsr2hyb
130aa372e3fSPaul Mullowney       #define cusparse_hyb2csr  cusparseShyb2csr
131aa372e3fSPaul Mullowney     #elif defined(PETSC_USE_REAL_DOUBLE)
132aa372e3fSPaul Mullowney       #define cusparse_csr_spmv cusparseDcsrmv
133ccdfe979SStefano Zampini       #define cusparse_csr_spmm cusparseDcsrmm
134aa372e3fSPaul Mullowney       #define cusparse_csr2csc  cusparseDcsr2csc
135aa372e3fSPaul Mullowney       #define cusparse_hyb_spmv cusparseDhybmv
136aa372e3fSPaul Mullowney       #define cusparse_csr2hyb  cusparseDcsr2hyb
137aa372e3fSPaul Mullowney       #define cusparse_hyb2csr  cusparseDhyb2csr
138aa372e3fSPaul Mullowney     #endif
139aa372e3fSPaul Mullowney   #endif
140afb2bd1cSJunchao Zhang #endif
141aa372e3fSPaul Mullowney 
142aa372e3fSPaul Mullowney #define THRUSTINTARRAY32 thrust::device_vector<int>
143aa372e3fSPaul Mullowney #define THRUSTINTARRAY thrust::device_vector<PetscInt>
144aa372e3fSPaul Mullowney #define THRUSTARRAY thrust::device_vector<PetscScalar>
145aa372e3fSPaul Mullowney 
146aa372e3fSPaul Mullowney /* A CSR matrix structure */
147aa372e3fSPaul Mullowney struct CsrMatrix {
148aa372e3fSPaul Mullowney   PetscInt         num_rows;
149aa372e3fSPaul Mullowney   PetscInt         num_cols;
150aa372e3fSPaul Mullowney   PetscInt         num_entries;
151aa372e3fSPaul Mullowney   THRUSTINTARRAY32 *row_offsets;
152aa372e3fSPaul Mullowney   THRUSTINTARRAY32 *column_indices;
153aa372e3fSPaul Mullowney   THRUSTARRAY      *values;
1549ae82921SPaul Mullowney };
1559ae82921SPaul Mullowney 
156aa372e3fSPaul Mullowney /* This is struct holding the relevant data needed to a MatSolve */
157aa372e3fSPaul Mullowney struct Mat_SeqAIJCUSPARSETriFactorStruct {
158aa372e3fSPaul Mullowney   /* Data needed for triangular solve */
159aa372e3fSPaul Mullowney   cusparseMatDescr_t          descr;
160aa372e3fSPaul Mullowney   cusparseOperation_t         solveOp;
161aa372e3fSPaul Mullowney   CsrMatrix                   *csrMat;
162*1b0a6780SStefano Zampini  #if PETSC_PKG_CUDA_VERSION_GE(9,0,0)
163afb2bd1cSJunchao Zhang   csrsv2Info_t                solveInfo;
164afb2bd1cSJunchao Zhang  #else
165afb2bd1cSJunchao Zhang   cusparseSolveAnalysisInfo_t solveInfo;
166afb2bd1cSJunchao Zhang  #endif
167*1b0a6780SStefano Zampini   cusparseSolvePolicy_t       solvePolicy;     /* whether level information is generated and used */
168*1b0a6780SStefano Zampini   int                         solveBufferSize;
169*1b0a6780SStefano Zampini   void                        *solveBuffer;
170*1b0a6780SStefano Zampini   size_t                      csr2cscBufferSize; /* to transpose the triangular factor (only used for CUDA >= 11.0) */
171*1b0a6780SStefano Zampini   void                        *csr2cscBuffer;
172*1b0a6780SStefano Zampini   Mat_SeqAIJCUSPARSETriFactorStruct() :
173*1b0a6780SStefano Zampini    csrMat(NULL),solveBuffer(NULL),csr2cscBuffer(NULL),solvePolicy(CUSPARSE_SOLVE_POLICY_USE_LEVEL) {} /* TODO: allow options database for policy */
174aa372e3fSPaul Mullowney };
175aa372e3fSPaul Mullowney 
176afb2bd1cSJunchao Zhang /* This is a larger struct holding all the triangular factors for a solve, transpose solve, and any indices used in a reordering */
177aa372e3fSPaul Mullowney struct Mat_SeqAIJCUSPARSETriFactors {
178aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorPtr; /* pointer for lower triangular (factored matrix) on GPU */
179aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorPtr; /* pointer for upper triangular (factored matrix) on GPU */
180aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorPtrTranspose; /* pointer for lower triangular (factored matrix) on GPU for the transpose (useful for BiCG) */
181aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorPtrTranspose; /* pointer for upper triangular (factored matrix) on GPU for the transpose (useful for BiCG)*/
182aa372e3fSPaul Mullowney   THRUSTINTARRAY                    *rpermIndices;  /* indices used for any reordering */
183aa372e3fSPaul Mullowney   THRUSTINTARRAY                    *cpermIndices;  /* indices used for any reordering */
184aa372e3fSPaul Mullowney   THRUSTARRAY                       *workVector;
185aa372e3fSPaul Mullowney   cusparseHandle_t                  handle;   /* a handle to the cusparse library */
186aa372e3fSPaul Mullowney   PetscInt                          nnz;      /* number of nonzeros ... need this for accurate logging between ICC and ILU */
187aa372e3fSPaul Mullowney };
188aa372e3fSPaul Mullowney 
189afb2bd1cSJunchao Zhang struct Mat_CusparseSpMV {
190afb2bd1cSJunchao Zhang   PetscBool             initialized;    /* Don't rely on spmvBuffer != NULL to test if the struct is initialized, */
191afb2bd1cSJunchao Zhang   size_t                spmvBufferSize; /* since I'm not sure if smvBuffer can be NULL even after cusparseSpMV_bufferSize() */
192afb2bd1cSJunchao Zhang   void                  *spmvBuffer;
1939db3cbf9SStefano 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 */
194afb2bd1cSJunchao Zhang   cusparseDnVecDescr_t  vecXDescr,vecYDescr; /* descriptor for the dense vectors in y=op(A)x */
1959db3cbf9SStefano Zampini  #endif
196afb2bd1cSJunchao Zhang };
197afb2bd1cSJunchao Zhang 
198afb2bd1cSJunchao Zhang /* This is struct holding the relevant data needed to a MatMult */
199afb2bd1cSJunchao Zhang struct Mat_SeqAIJCUSPARSEMultStruct {
200afb2bd1cSJunchao Zhang   void               *mat;  /* opaque pointer to a matrix. This could be either a cusparseHybMat_t or a CsrMatrix */
201afb2bd1cSJunchao Zhang   cusparseMatDescr_t descr; /* Data needed to describe the matrix for a multiply */
202afb2bd1cSJunchao Zhang   THRUSTINTARRAY     *cprowIndices;   /* compressed row indices used in the parallel SpMV */
203afb2bd1cSJunchao Zhang   PetscScalar        *alpha_one; /* pointer to a device "scalar" storing the alpha parameter in the SpMV */
204afb2bd1cSJunchao Zhang   PetscScalar        *beta_zero; /* pointer to a device "scalar" storing the beta parameter in the SpMV as zero*/
205afb2bd1cSJunchao Zhang   PetscScalar        *beta_one; /* pointer to a device "scalar" storing the beta parameter in the SpMV as one */
206afb2bd1cSJunchao Zhang  #if PETSC_PKG_CUDA_VERSION_GE(11,0,0)
207afb2bd1cSJunchao Zhang   cusparseSpMatDescr_t  matDescr;  /* descriptor for the matrix, used by SpMV and SpMM */
208afb2bd1cSJunchao Zhang   Mat_CusparseSpMV      cuSpMV[3]; /* different Mat_CusparseSpMV structs for non-transpose, transpose, conj-transpose */
209afb2bd1cSJunchao Zhang   Mat_SeqAIJCUSPARSEMultStruct() : matDescr(NULL) {
210afb2bd1cSJunchao Zhang     for (int i=0; i<3; i++) cuSpMV[i].initialized = PETSC_FALSE;
211afb2bd1cSJunchao Zhang   }
212afb2bd1cSJunchao Zhang  #endif
213afb2bd1cSJunchao Zhang };
214afb2bd1cSJunchao Zhang 
215aa372e3fSPaul Mullowney /* This is a larger struct holding all the matrices for a SpMV, and SpMV Tranpose */
2169ae82921SPaul Mullowney struct Mat_SeqAIJCUSPARSE {
217aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSEMultStruct *mat;            /* pointer to the matrix on the GPU */
218aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSEMultStruct *matTranspose;   /* pointer to the matrix on the GPU (for the transpose ... useful for BiCG) */
219aa372e3fSPaul Mullowney   THRUSTARRAY                  *workVector;     /* pointer to a workvector to which we can copy the relevant indices of a vector we want to multiply */
220029808eeSJunchao Zhang   THRUSTINTARRAY32             *rowoffsets_gpu; /* rowoffsets on GPU in non-compressed-row format. It is used to convert CSR to CSC */
221213423ffSJunchao Zhang   PetscInt                     nrows;           /* number of rows of the matrix seen by GPU */
222e057df02SPaul Mullowney   MatCUSPARSEStorageFormat     format;          /* the storage format for the matrix on the device */
223aa372e3fSPaul Mullowney   cudaStream_t                 stream;          /* a stream for the parallel SpMV ... this is not owned and should not be deleted */
224aa372e3fSPaul Mullowney   cusparseHandle_t             handle;          /* a handle to the cusparse library ... this may not be owned (if we're working in parallel i.e. multiGPUs) */
22554da937aSStefano Zampini   PetscObjectState             nonzerostate;    /* track nonzero state to possibly recreate the GPU matrix */
22654da937aSStefano Zampini   PetscBool                    transgen;        /* whether or not to generate explicit transpose for MatMultTranspose operations */
227afb2bd1cSJunchao Zhang  #if PETSC_PKG_CUDA_VERSION_GE(11,0,0)
228afb2bd1cSJunchao Zhang   size_t                       csr2cscBufferSize; /* stuff used to compute the matTranspose above */
229afb2bd1cSJunchao Zhang   void                         *csr2cscBuffer;    /* This is used as a C struct and is calloc'ed by PetscNewLog() */
230afb2bd1cSJunchao Zhang   cusparseCsr2CscAlg_t         csr2cscAlg;        /* algorithms can be selected from command line options */
231afb2bd1cSJunchao Zhang   cusparseSpMVAlg_t            spmvAlg;
232afb2bd1cSJunchao Zhang   cusparseSpMMAlg_t            spmmAlg;
233afb2bd1cSJunchao Zhang  #endif
2343fa6b06aSMark Adams   PetscSplitCSRDataStructure   *deviceMat;       /* Matrix on device for, eg, assembly */
2357e8381f9SStefano Zampini   THRUSTINTARRAY               *cooPerm;
2367e8381f9SStefano Zampini   THRUSTINTARRAY               *cooPerm_a;
2377e8381f9SStefano Zampini   THRUSTARRAY                  *cooPerm_v;
2387e8381f9SStefano Zampini   THRUSTARRAY                  *cooPerm_w;
2399ae82921SPaul Mullowney };
2409ae82921SPaul Mullowney 
2415a576424SJed Brown PETSC_INTERN PetscErrorCode MatCUSPARSECopyToGPU(Mat);
242b06137fdSPaul Mullowney PETSC_INTERN PetscErrorCode MatCUSPARSESetStream(Mat, const cudaStream_t stream);
243b06137fdSPaul Mullowney PETSC_INTERN PetscErrorCode MatCUSPARSESetHandle(Mat, const cusparseHandle_t handle);
244b06137fdSPaul Mullowney PETSC_INTERN PetscErrorCode MatCUSPARSEClearHandle(Mat);
2459ae82921SPaul Mullowney #endif
246