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> 18*7eaca502SStefano Zampini #include <thrust/system/system_error.h> 19c41cb2e2SAlejandro Lamas Daviña 20185e5899SJunchao Zhang #if (CUSPARSE_VER_MAJOR > 10 || CUSPARSE_VER_MAJOR == 10 && CUSPARSE_VER_MINOR >= 2) /* According to cuda/10.1.168 on OLCF Summit */ 2157d48284SJunchao Zhang #define CHKERRCUSPARSE(stat) \ 2257d48284SJunchao Zhang do { \ 2357d48284SJunchao Zhang if (PetscUnlikely(stat)) { \ 2457d48284SJunchao Zhang const char *name = cusparseGetErrorName(stat); \ 2557d48284SJunchao Zhang const char *descr = cusparseGetErrorString(stat); \ 26589f383fSBarry Smith SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_GPU,"cuSPARSE error %d (%s) : %s",(int)stat,name,descr); \ 2757d48284SJunchao Zhang } \ 2857d48284SJunchao Zhang } while (0) 296d22fe3dSJunchao Zhang #else 30589f383fSBarry Smith #define CHKERRCUSPARSE(stat) do {if (PetscUnlikely(stat)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_GPU,"cusparse error %d",(int)stat);} while (0) 316d22fe3dSJunchao Zhang #endif 3257d48284SJunchao Zhang 33*7eaca502SStefano Zampini #define PetscStackCallThrust(body) do { \ 34*7eaca502SStefano Zampini try { \ 35*7eaca502SStefano Zampini body; \ 36*7eaca502SStefano Zampini } catch(thrust::system_error& e) { \ 37*7eaca502SStefano Zampini SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Thrust %s",e.what());\ 38*7eaca502SStefano Zampini } \ 39*7eaca502SStefano Zampini } while (0) 40*7eaca502SStefano Zampini 41aa372e3fSPaul Mullowney #if defined(PETSC_USE_COMPLEX) 42aa372e3fSPaul Mullowney #if defined(PETSC_USE_REAL_SINGLE) 43ccdfe979SStefano Zampini const cuComplex PETSC_CUSPARSE_ONE = {1.0f, 0.0f}; 44ccdfe979SStefano Zampini const cuComplex PETSC_CUSPARSE_ZERO = {0.0f, 0.0f}; 45afb2bd1cSJunchao Zhang #elif defined(PETSC_USE_REAL_DOUBLE) 46afb2bd1cSJunchao Zhang const cuDoubleComplex PETSC_CUSPARSE_ONE = {1.0, 0.0}; 47afb2bd1cSJunchao Zhang const cuDoubleComplex PETSC_CUSPARSE_ZERO = {0.0, 0.0}; 48afb2bd1cSJunchao Zhang #endif 49afb2bd1cSJunchao Zhang #else 50afb2bd1cSJunchao Zhang const PetscScalar PETSC_CUSPARSE_ONE = 1.0; 51afb2bd1cSJunchao Zhang const PetscScalar PETSC_CUSPARSE_ZERO = 0.0; 52afb2bd1cSJunchao Zhang #endif 53afb2bd1cSJunchao Zhang 541b0a6780SStefano Zampini #if PETSC_PKG_CUDA_VERSION_GE(9,0,0) 55afb2bd1cSJunchao Zhang #define cusparse_create_analysis_info cusparseCreateCsrsv2Info 56afb2bd1cSJunchao Zhang #define cusparse_destroy_analysis_info cusparseDestroyCsrsv2Info 57afb2bd1cSJunchao Zhang #if defined(PETSC_USE_COMPLEX) 58afb2bd1cSJunchao Zhang #if defined(PETSC_USE_REAL_SINGLE) 59afb2bd1cSJunchao 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) 60afb2bd1cSJunchao 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) 61afb2bd1cSJunchao 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) 62afb2bd1cSJunchao Zhang #elif defined(PETSC_USE_REAL_DOUBLE) 63afb2bd1cSJunchao 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) 64afb2bd1cSJunchao 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) 65afb2bd1cSJunchao 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) 66afb2bd1cSJunchao Zhang #endif 67afb2bd1cSJunchao Zhang #else /* not complex */ 68afb2bd1cSJunchao Zhang #if defined(PETSC_USE_REAL_SINGLE) 69afb2bd1cSJunchao Zhang #define cusparse_get_svbuffsize cusparseScsrsv2_bufferSize 70afb2bd1cSJunchao Zhang #define cusparse_analysis cusparseScsrsv2_analysis 71afb2bd1cSJunchao Zhang #define cusparse_solve cusparseScsrsv2_solve 72afb2bd1cSJunchao Zhang #elif defined(PETSC_USE_REAL_DOUBLE) 73afb2bd1cSJunchao Zhang #define cusparse_get_svbuffsize cusparseDcsrsv2_bufferSize 74afb2bd1cSJunchao Zhang #define cusparse_analysis cusparseDcsrsv2_analysis 75afb2bd1cSJunchao Zhang #define cusparse_solve cusparseDcsrsv2_solve 76afb2bd1cSJunchao Zhang #endif 77afb2bd1cSJunchao Zhang #endif 78afb2bd1cSJunchao Zhang #else 79afb2bd1cSJunchao Zhang #define cusparse_create_analysis_info cusparseCreateSolveAnalysisInfo 80afb2bd1cSJunchao Zhang #define cusparse_destroy_analysis_info cusparseDestroySolveAnalysisInfo 81afb2bd1cSJunchao Zhang #if defined(PETSC_USE_COMPLEX) 82afb2bd1cSJunchao Zhang #if defined(PETSC_USE_REAL_SINGLE) 83ec42abe4SAlejandro 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)) 84ec42abe4SAlejandro 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)) 851b0a6780SStefano Zampini #elif defined(PETSC_USE_REAL_DOUBLE) 861b0a6780SStefano 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)) 871b0a6780SStefano 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)) 881b0a6780SStefano Zampini #endif 891b0a6780SStefano Zampini #else /* not complex */ 901b0a6780SStefano Zampini #if defined(PETSC_USE_REAL_SINGLE) 911b0a6780SStefano Zampini #define cusparse_solve cusparseScsrsv_solve 921b0a6780SStefano Zampini #define cusparse_analysis cusparseScsrsv_analysis 931b0a6780SStefano Zampini #elif defined(PETSC_USE_REAL_DOUBLE) 941b0a6780SStefano Zampini #define cusparse_solve cusparseDcsrsv_solve 951b0a6780SStefano Zampini #define cusparse_analysis cusparseDcsrsv_analysis 961b0a6780SStefano Zampini #endif 971b0a6780SStefano Zampini #endif 981b0a6780SStefano Zampini #endif 991b0a6780SStefano Zampini 1001b0a6780SStefano Zampini #if PETSC_PKG_CUDA_VERSION_GE(11,0,0) 1011b0a6780SStefano Zampini #define cusparse_csr2csc cusparseCsr2cscEx2 1021b0a6780SStefano Zampini #if defined(PETSC_USE_COMPLEX) 1031b0a6780SStefano Zampini #if defined(PETSC_USE_REAL_SINGLE) 1041b0a6780SStefano Zampini #define cusparse_scalartype CUDA_C_32F 1051b0a6780SStefano Zampini #elif defined(PETSC_USE_REAL_DOUBLE) 1061b0a6780SStefano Zampini #define cusparse_scalartype CUDA_C_64F 1071b0a6780SStefano Zampini #endif 1081b0a6780SStefano Zampini #else /* not complex */ 1091b0a6780SStefano Zampini #if defined(PETSC_USE_REAL_SINGLE) 1101b0a6780SStefano Zampini #define cusparse_scalartype CUDA_R_32F 1111b0a6780SStefano Zampini #elif defined(PETSC_USE_REAL_DOUBLE) 1121b0a6780SStefano Zampini #define cusparse_scalartype CUDA_R_64F 1131b0a6780SStefano Zampini #endif 1141b0a6780SStefano Zampini #endif 1151b0a6780SStefano Zampini #else 1161b0a6780SStefano Zampini #if defined(PETSC_USE_COMPLEX) 1171b0a6780SStefano Zampini #if defined(PETSC_USE_REAL_SINGLE) 118ec42abe4SAlejandro 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)) 119ccdfe979SStefano 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)) 120ec42abe4SAlejandro 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)) 121ec42abe4SAlejandro 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)) 122ec42abe4SAlejandro 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)) 123ec42abe4SAlejandro Lamas Daviña #define cusparse_hyb2csr(a,b,c,d,e,f) cusparseChyb2csr((a),(b),(c),(cuComplex*)(d),(e),(f)) 124fcdce8c4SStefano 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) 125aa372e3fSPaul Mullowney #elif defined(PETSC_USE_REAL_DOUBLE) 126ec42abe4SAlejandro 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)) 127ccdfe979SStefano 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)) 128ec42abe4SAlejandro 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)) 129ec42abe4SAlejandro 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)) 130ec42abe4SAlejandro 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)) 131ec42abe4SAlejandro Lamas Daviña #define cusparse_hyb2csr(a,b,c,d,e,f) cusparseZhyb2csr((a),(b),(c),(cuDoubleComplex*)(d),(e),(f)) 132fcdce8c4SStefano 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) 133aa372e3fSPaul Mullowney #endif 134aa372e3fSPaul Mullowney #else 135aa372e3fSPaul Mullowney #if defined(PETSC_USE_REAL_SINGLE) 136aa372e3fSPaul Mullowney #define cusparse_csr_spmv cusparseScsrmv 137ccdfe979SStefano Zampini #define cusparse_csr_spmm cusparseScsrmm 138aa372e3fSPaul Mullowney #define cusparse_csr2csc cusparseScsr2csc 139aa372e3fSPaul Mullowney #define cusparse_hyb_spmv cusparseShybmv 140aa372e3fSPaul Mullowney #define cusparse_csr2hyb cusparseScsr2hyb 141aa372e3fSPaul Mullowney #define cusparse_hyb2csr cusparseShyb2csr 142fcdce8c4SStefano Zampini #define cusparse_csr_spgemm cusparseScsrgemm 143aa372e3fSPaul Mullowney #elif defined(PETSC_USE_REAL_DOUBLE) 144aa372e3fSPaul Mullowney #define cusparse_csr_spmv cusparseDcsrmv 145ccdfe979SStefano Zampini #define cusparse_csr_spmm cusparseDcsrmm 146aa372e3fSPaul Mullowney #define cusparse_csr2csc cusparseDcsr2csc 147aa372e3fSPaul Mullowney #define cusparse_hyb_spmv cusparseDhybmv 148aa372e3fSPaul Mullowney #define cusparse_csr2hyb cusparseDcsr2hyb 149aa372e3fSPaul Mullowney #define cusparse_hyb2csr cusparseDhyb2csr 150fcdce8c4SStefano Zampini #define cusparse_csr_spgemm cusparseDcsrgemm 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 */ 199aa372e3fSPaul Mullowney }; 200aa372e3fSPaul Mullowney 201afb2bd1cSJunchao Zhang struct Mat_CusparseSpMV { 202afb2bd1cSJunchao Zhang PetscBool initialized; /* Don't rely on spmvBuffer != NULL to test if the struct is initialized, */ 203afb2bd1cSJunchao Zhang size_t spmvBufferSize; /* since I'm not sure if smvBuffer can be NULL even after cusparseSpMV_bufferSize() */ 204afb2bd1cSJunchao Zhang void *spmvBuffer; 2059db3cbf9SStefano 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 */ 206afb2bd1cSJunchao Zhang cusparseDnVecDescr_t vecXDescr,vecYDescr; /* descriptor for the dense vectors in y=op(A)x */ 2079db3cbf9SStefano Zampini #endif 208afb2bd1cSJunchao Zhang }; 209afb2bd1cSJunchao Zhang 210afb2bd1cSJunchao Zhang /* This is struct holding the relevant data needed to a MatMult */ 211afb2bd1cSJunchao Zhang struct Mat_SeqAIJCUSPARSEMultStruct { 212afb2bd1cSJunchao Zhang void *mat; /* opaque pointer to a matrix. This could be either a cusparseHybMat_t or a CsrMatrix */ 213afb2bd1cSJunchao Zhang cusparseMatDescr_t descr; /* Data needed to describe the matrix for a multiply */ 214afb2bd1cSJunchao Zhang THRUSTINTARRAY *cprowIndices; /* compressed row indices used in the parallel SpMV */ 215afb2bd1cSJunchao Zhang PetscScalar *alpha_one; /* pointer to a device "scalar" storing the alpha parameter in the SpMV */ 216afb2bd1cSJunchao Zhang PetscScalar *beta_zero; /* pointer to a device "scalar" storing the beta parameter in the SpMV as zero*/ 217afb2bd1cSJunchao Zhang PetscScalar *beta_one; /* pointer to a device "scalar" storing the beta parameter in the SpMV as one */ 218afb2bd1cSJunchao Zhang #if PETSC_PKG_CUDA_VERSION_GE(11,0,0) 219afb2bd1cSJunchao Zhang cusparseSpMatDescr_t matDescr; /* descriptor for the matrix, used by SpMV and SpMM */ 220afb2bd1cSJunchao Zhang Mat_CusparseSpMV cuSpMV[3]; /* different Mat_CusparseSpMV structs for non-transpose, transpose, conj-transpose */ 221afb2bd1cSJunchao Zhang Mat_SeqAIJCUSPARSEMultStruct() : matDescr(NULL) { 222afb2bd1cSJunchao Zhang for (int i=0; i<3; i++) cuSpMV[i].initialized = PETSC_FALSE; 223afb2bd1cSJunchao Zhang } 224afb2bd1cSJunchao Zhang #endif 225afb2bd1cSJunchao Zhang }; 226afb2bd1cSJunchao Zhang 227aa372e3fSPaul Mullowney /* This is a larger struct holding all the matrices for a SpMV, and SpMV Tranpose */ 2289ae82921SPaul Mullowney struct Mat_SeqAIJCUSPARSE { 229aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSEMultStruct *mat; /* pointer to the matrix on the GPU */ 230aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSEMultStruct *matTranspose; /* pointer to the matrix on the GPU (for the transpose ... useful for BiCG) */ 231aa372e3fSPaul Mullowney THRUSTARRAY *workVector; /* pointer to a workvector to which we can copy the relevant indices of a vector we want to multiply */ 232029808eeSJunchao Zhang THRUSTINTARRAY32 *rowoffsets_gpu; /* rowoffsets on GPU in non-compressed-row format. It is used to convert CSR to CSC */ 233213423ffSJunchao Zhang PetscInt nrows; /* number of rows of the matrix seen by GPU */ 234e057df02SPaul Mullowney MatCUSPARSEStorageFormat format; /* the storage format for the matrix on the device */ 235aa372e3fSPaul Mullowney cudaStream_t stream; /* a stream for the parallel SpMV ... this is not owned and should not be deleted */ 236aa372e3fSPaul Mullowney cusparseHandle_t handle; /* a handle to the cusparse library ... this may not be owned (if we're working in parallel i.e. multiGPUs) */ 23754da937aSStefano Zampini PetscObjectState nonzerostate; /* track nonzero state to possibly recreate the GPU matrix */ 23854da937aSStefano Zampini PetscBool transgen; /* whether or not to generate explicit transpose for MatMultTranspose operations */ 239afb2bd1cSJunchao Zhang #if PETSC_PKG_CUDA_VERSION_GE(11,0,0) 240afb2bd1cSJunchao Zhang size_t csr2cscBufferSize; /* stuff used to compute the matTranspose above */ 241afb2bd1cSJunchao Zhang void *csr2cscBuffer; /* This is used as a C struct and is calloc'ed by PetscNewLog() */ 242afb2bd1cSJunchao Zhang cusparseCsr2CscAlg_t csr2cscAlg; /* algorithms can be selected from command line options */ 243afb2bd1cSJunchao Zhang cusparseSpMVAlg_t spmvAlg; 244afb2bd1cSJunchao Zhang cusparseSpMMAlg_t spmmAlg; 245afb2bd1cSJunchao Zhang #endif 2463fa6b06aSMark Adams PetscSplitCSRDataStructure *deviceMat; /* Matrix on device for, eg, assembly */ 2477e8381f9SStefano Zampini THRUSTINTARRAY *cooPerm; 2487e8381f9SStefano Zampini THRUSTINTARRAY *cooPerm_a; 2499ae82921SPaul Mullowney }; 2509ae82921SPaul Mullowney 2515a576424SJed Brown PETSC_INTERN PetscErrorCode MatCUSPARSECopyToGPU(Mat); 252b06137fdSPaul Mullowney PETSC_INTERN PetscErrorCode MatCUSPARSESetStream(Mat, const cudaStream_t stream); 253b06137fdSPaul Mullowney PETSC_INTERN PetscErrorCode MatCUSPARSESetHandle(Mat, const cusparseHandle_t handle); 254b06137fdSPaul Mullowney PETSC_INTERN PetscErrorCode MatCUSPARSEClearHandle(Mat); 255e61fc153SStefano Zampini PETSC_INTERN PetscErrorCode MatSetPreallocationCOO_SeqAIJCUSPARSE(Mat,PetscInt,const PetscInt[],const PetscInt[]); 256e61fc153SStefano Zampini PETSC_INTERN PetscErrorCode MatSetValuesCOO_SeqAIJCUSPARSE(Mat,const PetscScalar[],InsertMode); 257ed502f03SStefano Zampini PETSC_INTERN PetscErrorCode MatSeqAIJCUSPARSEGetArrayRead(Mat,const PetscScalar**); 258ed502f03SStefano Zampini PETSC_INTERN PetscErrorCode MatSeqAIJCUSPARSERestoreArrayRead(Mat,const PetscScalar**); 259ed502f03SStefano Zampini PETSC_INTERN PetscErrorCode MatSeqAIJCUSPARSEGetArrayWrite(Mat,PetscScalar**); 260ed502f03SStefano Zampini PETSC_INTERN PetscErrorCode MatSeqAIJCUSPARSERestoreArrayWrite(Mat,PetscScalar**); 261ed502f03SStefano Zampini PETSC_INTERN PetscErrorCode MatSeqAIJCUSPARSEMergeMats(Mat,Mat,MatReuse,Mat*); 262ed502f03SStefano Zampini 26308391a17SStefano Zampini PETSC_STATIC_INLINE bool isCudaMem(const void *data) 26408391a17SStefano Zampini { 26508391a17SStefano Zampini cudaError_t cerr; 26608391a17SStefano Zampini struct cudaPointerAttributes attr; 26708391a17SStefano Zampini enum cudaMemoryType mtype; 26808391a17SStefano Zampini cerr = cudaPointerGetAttributes(&attr,data); /* Do not check error since before CUDA 11.0, passing a host pointer returns cudaErrorInvalidValue */ 26908391a17SStefano Zampini cudaGetLastError(); /* Reset the last error */ 27008391a17SStefano Zampini #if (CUDART_VERSION < 10000) 27108391a17SStefano Zampini mtype = attr.memoryType; 27208391a17SStefano Zampini #else 27308391a17SStefano Zampini mtype = attr.type; 27408391a17SStefano Zampini #endif 27508391a17SStefano Zampini if (cerr == cudaSuccess && mtype == cudaMemoryTypeDevice) return true; 27608391a17SStefano Zampini else return false; 27708391a17SStefano Zampini } 27808391a17SStefano Zampini 2799ae82921SPaul Mullowney #endif 280