xref: /petsc/src/mat/impls/aij/seq/seqcusparse/cusparsematimpl.h (revision 3fa6b06a16b9e10b6e41dcaa1d4ab6dd90b29d60)
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); \
2557d48284SJunchao Zhang       SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_LIB,"cuSPARSE error %d (%s) : %s",(int)stat,name,descr); \
2657d48284SJunchao Zhang    } \
2757d48284SJunchao Zhang } while (0)
286d22fe3dSJunchao Zhang #else
29ce29e4ddSJunchao Zhang #define CHKERRCUSPARSE(stat) do {if (PetscUnlikely(stat)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"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 
45afb2bd1cSJunchao Zhang #if PETSC_PKG_CUDA_VERSION_GE(11,0,0)
46afb2bd1cSJunchao Zhang   #define cusparse_create_analysis_info   cusparseCreateCsrsv2Info
47afb2bd1cSJunchao Zhang   #define cusparse_destroy_analysis_info  cusparseDestroyCsrsv2Info
48afb2bd1cSJunchao Zhang   #define cusparse_csr2csc                cusparseCsr2cscEx2
49afb2bd1cSJunchao Zhang   #if defined(PETSC_USE_COMPLEX)
50afb2bd1cSJunchao Zhang     #if defined(PETSC_USE_REAL_SINGLE)
51afb2bd1cSJunchao Zhang       #define cusparse_scalartype                            CUDA_C_32F
52afb2bd1cSJunchao 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)
53afb2bd1cSJunchao 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)
54afb2bd1cSJunchao 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)
55afb2bd1cSJunchao Zhang     #elif defined(PETSC_USE_REAL_DOUBLE)
56afb2bd1cSJunchao Zhang       #define cusparse_scalartype                            CUDA_C_64F
57afb2bd1cSJunchao 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)
58afb2bd1cSJunchao 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)
59afb2bd1cSJunchao 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)
60afb2bd1cSJunchao Zhang     #endif
61afb2bd1cSJunchao Zhang   #else /* not complex */
62afb2bd1cSJunchao Zhang     #if defined(PETSC_USE_REAL_SINGLE)
63afb2bd1cSJunchao Zhang       #define cusparse_scalartype         CUDA_R_32F
64afb2bd1cSJunchao Zhang       #define cusparse_get_svbuffsize     cusparseScsrsv2_bufferSize
65afb2bd1cSJunchao Zhang       #define cusparse_analysis           cusparseScsrsv2_analysis
66afb2bd1cSJunchao Zhang       #define cusparse_solve              cusparseScsrsv2_solve
67afb2bd1cSJunchao Zhang     #elif defined(PETSC_USE_REAL_DOUBLE)
68afb2bd1cSJunchao Zhang       #define cusparse_scalartype         CUDA_R_64F
69afb2bd1cSJunchao Zhang       #define cusparse_get_svbuffsize     cusparseDcsrsv2_bufferSize
70afb2bd1cSJunchao Zhang       #define cusparse_analysis           cusparseDcsrsv2_analysis
71afb2bd1cSJunchao Zhang       #define cusparse_solve              cusparseDcsrsv2_solve
72afb2bd1cSJunchao Zhang     #endif
73afb2bd1cSJunchao Zhang   #endif
74afb2bd1cSJunchao Zhang #else
75afb2bd1cSJunchao Zhang   #define cusparse_create_analysis_info   cusparseCreateSolveAnalysisInfo
76afb2bd1cSJunchao Zhang   #define cusparse_destroy_analysis_info  cusparseDestroySolveAnalysisInfo
77afb2bd1cSJunchao Zhang   #if defined(PETSC_USE_COMPLEX)
78afb2bd1cSJunchao Zhang     #if defined(PETSC_USE_REAL_SINGLE)
79ec42abe4SAlejandro 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))
80ec42abe4SAlejandro 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))
81ec42abe4SAlejandro 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))
82ccdfe979SStefano 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))
83ec42abe4SAlejandro 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))
84ec42abe4SAlejandro 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))
85ec42abe4SAlejandro 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))
86ec42abe4SAlejandro Lamas Daviña       #define cusparse_hyb2csr(a,b,c,d,e,f)                      cusparseChyb2csr((a),(b),(c),(cuComplex*)(d),(e),(f))
87aa372e3fSPaul Mullowney     #elif defined(PETSC_USE_REAL_DOUBLE)
88ec42abe4SAlejandro Lamas Daviña       #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))
89ec42abe4SAlejandro Lamas Daviña       #define cusparse_analysis(a,b,c,d,e,f,g,h,i)               cusparseZcsrsv_analysis((a),(b),(c),(d),(e),(cuDoubleComplex*)(f),(g),(h),(i))
90ec42abe4SAlejandro 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))
91ccdfe979SStefano 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))
92ec42abe4SAlejandro 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))
93ec42abe4SAlejandro 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))
94ec42abe4SAlejandro 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))
95ec42abe4SAlejandro Lamas Daviña       #define cusparse_hyb2csr(a,b,c,d,e,f)                      cusparseZhyb2csr((a),(b),(c),(cuDoubleComplex*)(d),(e),(f))
96aa372e3fSPaul Mullowney     #endif
97aa372e3fSPaul Mullowney   #else
98aa372e3fSPaul Mullowney     #if defined(PETSC_USE_REAL_SINGLE)
99aa372e3fSPaul Mullowney       #define cusparse_solve    cusparseScsrsv_solve
100aa372e3fSPaul Mullowney       #define cusparse_analysis cusparseScsrsv_analysis
101aa372e3fSPaul Mullowney       #define cusparse_csr_spmv cusparseScsrmv
102ccdfe979SStefano Zampini       #define cusparse_csr_spmm cusparseScsrmm
103aa372e3fSPaul Mullowney       #define cusparse_csr2csc  cusparseScsr2csc
104aa372e3fSPaul Mullowney       #define cusparse_hyb_spmv cusparseShybmv
105aa372e3fSPaul Mullowney       #define cusparse_csr2hyb  cusparseScsr2hyb
106aa372e3fSPaul Mullowney       #define cusparse_hyb2csr  cusparseShyb2csr
107aa372e3fSPaul Mullowney     #elif defined(PETSC_USE_REAL_DOUBLE)
108aa372e3fSPaul Mullowney       #define cusparse_solve    cusparseDcsrsv_solve
109aa372e3fSPaul Mullowney       #define cusparse_analysis cusparseDcsrsv_analysis
110aa372e3fSPaul Mullowney       #define cusparse_csr_spmv cusparseDcsrmv
111ccdfe979SStefano Zampini       #define cusparse_csr_spmm cusparseDcsrmm
112aa372e3fSPaul Mullowney       #define cusparse_csr2csc  cusparseDcsr2csc
113aa372e3fSPaul Mullowney       #define cusparse_hyb_spmv cusparseDhybmv
114aa372e3fSPaul Mullowney       #define cusparse_csr2hyb  cusparseDcsr2hyb
115aa372e3fSPaul Mullowney       #define cusparse_hyb2csr  cusparseDhyb2csr
116aa372e3fSPaul Mullowney     #endif
117aa372e3fSPaul Mullowney   #endif
118afb2bd1cSJunchao Zhang #endif
119aa372e3fSPaul Mullowney 
120aa372e3fSPaul Mullowney #define THRUSTINTARRAY32 thrust::device_vector<int>
121aa372e3fSPaul Mullowney #define THRUSTINTARRAY thrust::device_vector<PetscInt>
122aa372e3fSPaul Mullowney #define THRUSTARRAY thrust::device_vector<PetscScalar>
123aa372e3fSPaul Mullowney 
124aa372e3fSPaul Mullowney /* A CSR matrix structure */
125aa372e3fSPaul Mullowney struct CsrMatrix {
126aa372e3fSPaul Mullowney   PetscInt         num_rows;
127aa372e3fSPaul Mullowney   PetscInt         num_cols;
128aa372e3fSPaul Mullowney   PetscInt         num_entries;
129aa372e3fSPaul Mullowney   THRUSTINTARRAY32 *row_offsets;
130aa372e3fSPaul Mullowney   THRUSTINTARRAY32 *column_indices;
131aa372e3fSPaul Mullowney   THRUSTARRAY      *values;
1329ae82921SPaul Mullowney };
1339ae82921SPaul Mullowney 
134aa372e3fSPaul Mullowney /* This is struct holding the relevant data needed to a MatSolve */
135aa372e3fSPaul Mullowney struct Mat_SeqAIJCUSPARSETriFactorStruct {
136aa372e3fSPaul Mullowney   /* Data needed for triangular solve */
137aa372e3fSPaul Mullowney   cusparseMatDescr_t          descr;
138aa372e3fSPaul Mullowney   cusparseOperation_t         solveOp;
139aa372e3fSPaul Mullowney   CsrMatrix                   *csrMat;
140afb2bd1cSJunchao Zhang  #if PETSC_PKG_CUDA_VERSION_GE(11,0,0)
141afb2bd1cSJunchao Zhang   csrsv2Info_t                solveInfo;
142afb2bd1cSJunchao Zhang   cusparseSolvePolicy_t       solvePolicy;     /* whether level information is generated and used */
143afb2bd1cSJunchao Zhang   int                         solveBufferSize;
144afb2bd1cSJunchao Zhang   void                        *solveBuffer;
145afb2bd1cSJunchao Zhang   size_t                      csr2cscBufferSize; /* to transpose the triangular factor */
146afb2bd1cSJunchao Zhang   void                        *csr2cscBuffer;
147afb2bd1cSJunchao Zhang   Mat_SeqAIJCUSPARSETriFactorStruct() :
148afb2bd1cSJunchao Zhang    csrMat(NULL),solveBuffer(NULL),csr2cscBuffer(NULL), solvePolicy(CUSPARSE_SOLVE_POLICY_NO_LEVEL){} /* TODO: allow options database for policy */
149afb2bd1cSJunchao Zhang  #else
150afb2bd1cSJunchao Zhang   cusparseSolveAnalysisInfo_t solveInfo;
151afb2bd1cSJunchao Zhang  #endif
152aa372e3fSPaul Mullowney };
153aa372e3fSPaul Mullowney 
154afb2bd1cSJunchao Zhang /* This is a larger struct holding all the triangular factors for a solve, transpose solve, and any indices used in a reordering */
155aa372e3fSPaul Mullowney struct Mat_SeqAIJCUSPARSETriFactors {
156aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorPtr; /* pointer for lower triangular (factored matrix) on GPU */
157aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorPtr; /* pointer for upper triangular (factored matrix) on GPU */
158aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorPtrTranspose; /* pointer for lower triangular (factored matrix) on GPU for the transpose (useful for BiCG) */
159aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorPtrTranspose; /* pointer for upper triangular (factored matrix) on GPU for the transpose (useful for BiCG)*/
160aa372e3fSPaul Mullowney   THRUSTINTARRAY                    *rpermIndices;  /* indices used for any reordering */
161aa372e3fSPaul Mullowney   THRUSTINTARRAY                    *cpermIndices;  /* indices used for any reordering */
162aa372e3fSPaul Mullowney   THRUSTARRAY                       *workVector;
163aa372e3fSPaul Mullowney   cusparseHandle_t                  handle;   /* a handle to the cusparse library */
164aa372e3fSPaul Mullowney   PetscInt                          nnz;      /* number of nonzeros ... need this for accurate logging between ICC and ILU */
165aa372e3fSPaul Mullowney };
166aa372e3fSPaul Mullowney 
167afb2bd1cSJunchao Zhang struct Mat_CusparseSpMV {
168afb2bd1cSJunchao Zhang   PetscBool             initialized;    /* Don't rely on spmvBuffer != NULL to test if the struct is initialized, */
169afb2bd1cSJunchao Zhang   size_t                spmvBufferSize; /* since I'm not sure if smvBuffer can be NULL even after cusparseSpMV_bufferSize() */
170afb2bd1cSJunchao Zhang   void                  *spmvBuffer;
171afb2bd1cSJunchao Zhang   cusparseDnVecDescr_t  vecXDescr,vecYDescr; /* descriptor for the dense vectors in y=op(A)x */
172afb2bd1cSJunchao Zhang };
173afb2bd1cSJunchao Zhang 
174afb2bd1cSJunchao Zhang /* This is struct holding the relevant data needed to a MatMult */
175afb2bd1cSJunchao Zhang struct Mat_SeqAIJCUSPARSEMultStruct {
176afb2bd1cSJunchao Zhang   void               *mat;  /* opaque pointer to a matrix. This could be either a cusparseHybMat_t or a CsrMatrix */
177afb2bd1cSJunchao Zhang   cusparseMatDescr_t descr; /* Data needed to describe the matrix for a multiply */
178afb2bd1cSJunchao Zhang   THRUSTINTARRAY     *cprowIndices;   /* compressed row indices used in the parallel SpMV */
179afb2bd1cSJunchao Zhang   PetscScalar        *alpha_one; /* pointer to a device "scalar" storing the alpha parameter in the SpMV */
180afb2bd1cSJunchao Zhang   PetscScalar        *beta_zero; /* pointer to a device "scalar" storing the beta parameter in the SpMV as zero*/
181afb2bd1cSJunchao Zhang   PetscScalar        *beta_one; /* pointer to a device "scalar" storing the beta parameter in the SpMV as one */
182afb2bd1cSJunchao Zhang  #if PETSC_PKG_CUDA_VERSION_GE(11,0,0)
183afb2bd1cSJunchao Zhang   cusparseSpMatDescr_t  matDescr;  /* descriptor for the matrix, used by SpMV and SpMM */
184afb2bd1cSJunchao Zhang   Mat_CusparseSpMV      cuSpMV[3]; /* different Mat_CusparseSpMV structs for non-transpose, transpose, conj-transpose */
185afb2bd1cSJunchao Zhang   Mat_SeqAIJCUSPARSEMultStruct() : matDescr(NULL) {
186afb2bd1cSJunchao Zhang     for (int i=0; i<3; i++) cuSpMV[i].initialized = PETSC_FALSE;
187afb2bd1cSJunchao Zhang   }
188afb2bd1cSJunchao Zhang  #endif
189afb2bd1cSJunchao Zhang };
190afb2bd1cSJunchao Zhang 
191aa372e3fSPaul Mullowney /* This is a larger struct holding all the matrices for a SpMV, and SpMV Tranpose */
1929ae82921SPaul Mullowney struct Mat_SeqAIJCUSPARSE {
193aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSEMultStruct *mat;            /* pointer to the matrix on the GPU */
194aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSEMultStruct *matTranspose;   /* pointer to the matrix on the GPU (for the transpose ... useful for BiCG) */
195aa372e3fSPaul Mullowney   THRUSTARRAY                  *workVector;     /* pointer to a workvector to which we can copy the relevant indices of a vector we want to multiply */
196029808eeSJunchao Zhang   THRUSTINTARRAY32             *rowoffsets_gpu; /* rowoffsets on GPU in non-compressed-row format. It is used to convert CSR to CSC */
197213423ffSJunchao Zhang   PetscInt                     nrows;           /* number of rows of the matrix seen by GPU */
198e057df02SPaul Mullowney   MatCUSPARSEStorageFormat     format;          /* the storage format for the matrix on the device */
199aa372e3fSPaul Mullowney   cudaStream_t                 stream;          /* a stream for the parallel SpMV ... this is not owned and should not be deleted */
200aa372e3fSPaul Mullowney   cusparseHandle_t             handle;          /* a handle to the cusparse library ... this may not be owned (if we're working in parallel i.e. multiGPUs) */
20154da937aSStefano Zampini   PetscObjectState             nonzerostate;    /* track nonzero state to possibly recreate the GPU matrix */
20254da937aSStefano Zampini   PetscBool                    transgen;        /* whether or not to generate explicit transpose for MatMultTranspose operations */
203afb2bd1cSJunchao Zhang  #if PETSC_PKG_CUDA_VERSION_GE(11,0,0)
204afb2bd1cSJunchao Zhang   size_t                       csr2cscBufferSize; /* stuff used to compute the matTranspose above */
205afb2bd1cSJunchao Zhang   void                         *csr2cscBuffer;    /* This is used as a C struct and is calloc'ed by PetscNewLog() */
206afb2bd1cSJunchao Zhang   cusparseCsr2CscAlg_t         csr2cscAlg;        /* algorithms can be selected from command line options */
207afb2bd1cSJunchao Zhang   cusparseSpMVAlg_t            spmvAlg;
208afb2bd1cSJunchao Zhang   cusparseSpMMAlg_t            spmmAlg;
209afb2bd1cSJunchao Zhang  #endif
210*3fa6b06aSMark Adams   PetscSplitCSRDataStructure   *deviceMat;       /* Matrix on device for, eg, assembly */
2119ae82921SPaul Mullowney };
2129ae82921SPaul Mullowney 
2135a576424SJed Brown PETSC_INTERN PetscErrorCode MatCUSPARSECopyToGPU(Mat);
214b06137fdSPaul Mullowney PETSC_INTERN PetscErrorCode MatCUSPARSESetStream(Mat, const cudaStream_t stream);
215b06137fdSPaul Mullowney PETSC_INTERN PetscErrorCode MatCUSPARSESetHandle(Mat, const cusparseHandle_t handle);
216b06137fdSPaul Mullowney PETSC_INTERN PetscErrorCode MatCUSPARSEClearHandle(Mat);
2179ae82921SPaul Mullowney #endif
218