1 #if !defined(__CUSPARSEMATIMPL) 2 #define __CUSPARSEMATIMPL 3 4 #include <petscpkg_version.h> 5 #include <petsc/private/cudavecimpl.h> 6 7 #include <cusparse_v2.h> 8 9 #include <algorithm> 10 #include <vector> 11 12 #include <thrust/device_vector.h> 13 #include <thrust/device_ptr.h> 14 #include <thrust/device_malloc_allocator.h> 15 #include <thrust/transform.h> 16 #include <thrust/functional.h> 17 #include <thrust/sequence.h> 18 #include <thrust/system/system_error.h> 19 20 #if (CUSPARSE_VER_MAJOR > 10 || CUSPARSE_VER_MAJOR == 10 && CUSPARSE_VER_MINOR >= 2) /* According to cuda/10.1.168 on OLCF Summit */ 21 #define CHKERRCUSPARSE(stat) \ 22 do { \ 23 if (PetscUnlikely(stat)) { \ 24 const char *name = cusparseGetErrorName(stat); \ 25 const char *descr = cusparseGetErrorString(stat); \ 26 SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_GPU,"cuSPARSE error %d (%s) : %s",(int)stat,name,descr); \ 27 } \ 28 } while (0) 29 #else 30 #define CHKERRCUSPARSE(stat) do {if (PetscUnlikely(stat)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_GPU,"cusparse error %d",(int)stat);} while (0) 31 #endif 32 33 #define PetscStackCallThrust(body) do { \ 34 try { \ 35 body; \ 36 } catch(thrust::system_error& e) { \ 37 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Thrust %s",e.what());\ 38 } \ 39 } while (0) 40 41 #if defined(PETSC_USE_COMPLEX) 42 #if defined(PETSC_USE_REAL_SINGLE) 43 const cuComplex PETSC_CUSPARSE_ONE = {1.0f, 0.0f}; 44 const cuComplex PETSC_CUSPARSE_ZERO = {0.0f, 0.0f}; 45 #elif defined(PETSC_USE_REAL_DOUBLE) 46 const cuDoubleComplex PETSC_CUSPARSE_ONE = {1.0, 0.0}; 47 const cuDoubleComplex PETSC_CUSPARSE_ZERO = {0.0, 0.0}; 48 #endif 49 #else 50 const PetscScalar PETSC_CUSPARSE_ONE = 1.0; 51 const PetscScalar PETSC_CUSPARSE_ZERO = 0.0; 52 #endif 53 54 #if PETSC_PKG_CUDA_VERSION_GE(9,0,0) 55 #define cusparse_create_analysis_info cusparseCreateCsrsv2Info 56 #define cusparse_destroy_analysis_info cusparseDestroyCsrsv2Info 57 #if defined(PETSC_USE_COMPLEX) 58 #if defined(PETSC_USE_REAL_SINGLE) 59 #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) 60 #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) 61 #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) 62 #elif defined(PETSC_USE_REAL_DOUBLE) 63 #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) 64 #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) 65 #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) 66 #endif 67 #else /* not complex */ 68 #if defined(PETSC_USE_REAL_SINGLE) 69 #define cusparse_get_svbuffsize cusparseScsrsv2_bufferSize 70 #define cusparse_analysis cusparseScsrsv2_analysis 71 #define cusparse_solve cusparseScsrsv2_solve 72 #elif defined(PETSC_USE_REAL_DOUBLE) 73 #define cusparse_get_svbuffsize cusparseDcsrsv2_bufferSize 74 #define cusparse_analysis cusparseDcsrsv2_analysis 75 #define cusparse_solve cusparseDcsrsv2_solve 76 #endif 77 #endif 78 #else 79 #define cusparse_create_analysis_info cusparseCreateSolveAnalysisInfo 80 #define cusparse_destroy_analysis_info cusparseDestroySolveAnalysisInfo 81 #if defined(PETSC_USE_COMPLEX) 82 #if defined(PETSC_USE_REAL_SINGLE) 83 #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)) 84 #define cusparse_analysis(a,b,c,d,e,f,g,h,i) cusparseCcsrsv_analysis((a),(b),(c),(d),(e),(cuComplex*)(f),(g),(h),(i)) 85 #elif defined(PETSC_USE_REAL_DOUBLE) 86 #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)) 87 #define cusparse_analysis(a,b,c,d,e,f,g,h,i) cusparseZcsrsv_analysis((a),(b),(c),(d),(e),(cuDoubleComplex*)(f),(g),(h),(i)) 88 #endif 89 #else /* not complex */ 90 #if defined(PETSC_USE_REAL_SINGLE) 91 #define cusparse_solve cusparseScsrsv_solve 92 #define cusparse_analysis cusparseScsrsv_analysis 93 #elif defined(PETSC_USE_REAL_DOUBLE) 94 #define cusparse_solve cusparseDcsrsv_solve 95 #define cusparse_analysis cusparseDcsrsv_analysis 96 #endif 97 #endif 98 #endif 99 100 #if PETSC_PKG_CUDA_VERSION_GE(11,0,0) 101 #define cusparse_csr2csc cusparseCsr2cscEx2 102 #if defined(PETSC_USE_COMPLEX) 103 #if defined(PETSC_USE_REAL_SINGLE) 104 #define cusparse_scalartype CUDA_C_32F 105 #elif defined(PETSC_USE_REAL_DOUBLE) 106 #define cusparse_scalartype CUDA_C_64F 107 #endif 108 #else /* not complex */ 109 #if defined(PETSC_USE_REAL_SINGLE) 110 #define cusparse_scalartype CUDA_R_32F 111 #elif defined(PETSC_USE_REAL_DOUBLE) 112 #define cusparse_scalartype CUDA_R_64F 113 #endif 114 #endif 115 #else 116 #if defined(PETSC_USE_COMPLEX) 117 #if defined(PETSC_USE_REAL_SINGLE) 118 #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)) 119 #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)) 120 #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)) 121 #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)) 122 #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)) 123 #define cusparse_hyb2csr(a,b,c,d,e,f) cusparseChyb2csr((a),(b),(c),(cuComplex*)(d),(e),(f)) 124 #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) 125 #elif defined(PETSC_USE_REAL_DOUBLE) 126 #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)) 127 #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)) 128 #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)) 129 #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)) 130 #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)) 131 #define cusparse_hyb2csr(a,b,c,d,e,f) cusparseZhyb2csr((a),(b),(c),(cuDoubleComplex*)(d),(e),(f)) 132 #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) 133 #endif 134 #else 135 #if defined(PETSC_USE_REAL_SINGLE) 136 #define cusparse_csr_spmv cusparseScsrmv 137 #define cusparse_csr_spmm cusparseScsrmm 138 #define cusparse_csr2csc cusparseScsr2csc 139 #define cusparse_hyb_spmv cusparseShybmv 140 #define cusparse_csr2hyb cusparseScsr2hyb 141 #define cusparse_hyb2csr cusparseShyb2csr 142 #define cusparse_csr_spgemm cusparseScsrgemm 143 #elif defined(PETSC_USE_REAL_DOUBLE) 144 #define cusparse_csr_spmv cusparseDcsrmv 145 #define cusparse_csr_spmm cusparseDcsrmm 146 #define cusparse_csr2csc cusparseDcsr2csc 147 #define cusparse_hyb_spmv cusparseDhybmv 148 #define cusparse_csr2hyb cusparseDcsr2hyb 149 #define cusparse_hyb2csr cusparseDhyb2csr 150 #define cusparse_csr_spgemm cusparseDcsrgemm 151 #endif 152 #endif 153 #endif 154 155 #define THRUSTINTARRAY32 thrust::device_vector<int> 156 #define THRUSTINTARRAY thrust::device_vector<PetscInt> 157 #define THRUSTARRAY thrust::device_vector<PetscScalar> 158 159 /* A CSR matrix structure */ 160 struct CsrMatrix { 161 PetscInt num_rows; 162 PetscInt num_cols; 163 PetscInt num_entries; 164 THRUSTINTARRAY32 *row_offsets; 165 THRUSTINTARRAY32 *column_indices; 166 THRUSTARRAY *values; 167 }; 168 169 /* This is struct holding the relevant data needed to a MatSolve */ 170 struct Mat_SeqAIJCUSPARSETriFactorStruct { 171 /* Data needed for triangular solve */ 172 cusparseMatDescr_t descr; 173 cusparseOperation_t solveOp; 174 CsrMatrix *csrMat; 175 #if PETSC_PKG_CUDA_VERSION_GE(9,0,0) 176 csrsv2Info_t solveInfo; 177 #else 178 cusparseSolveAnalysisInfo_t solveInfo; 179 #endif 180 cusparseSolvePolicy_t solvePolicy; /* whether level information is generated and used */ 181 int solveBufferSize; 182 void *solveBuffer; 183 size_t csr2cscBufferSize; /* to transpose the triangular factor (only used for CUDA >= 11.0) */ 184 void *csr2cscBuffer; 185 PetscScalar *AA_h; /* managed host buffer for moving values to the GPU */ 186 }; 187 188 /* This is a larger struct holding all the triangular factors for a solve, transpose solve, and any indices used in a reordering */ 189 struct Mat_SeqAIJCUSPARSETriFactors { 190 Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorPtr; /* pointer for lower triangular (factored matrix) on GPU */ 191 Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorPtr; /* pointer for upper triangular (factored matrix) on GPU */ 192 Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorPtrTranspose; /* pointer for lower triangular (factored matrix) on GPU for the transpose (useful for BiCG) */ 193 Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorPtrTranspose; /* pointer for upper triangular (factored matrix) on GPU for the transpose (useful for BiCG)*/ 194 THRUSTINTARRAY *rpermIndices; /* indices used for any reordering */ 195 THRUSTINTARRAY *cpermIndices; /* indices used for any reordering */ 196 THRUSTARRAY *workVector; 197 cusparseHandle_t handle; /* a handle to the cusparse library */ 198 PetscInt nnz; /* number of nonzeros ... need this for accurate logging between ICC and ILU */ 199 }; 200 201 struct Mat_CusparseSpMV { 202 PetscBool initialized; /* Don't rely on spmvBuffer != NULL to test if the struct is initialized, */ 203 size_t spmvBufferSize; /* since I'm not sure if smvBuffer can be NULL even after cusparseSpMV_bufferSize() */ 204 void *spmvBuffer; 205 #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 */ 206 cusparseDnVecDescr_t vecXDescr,vecYDescr; /* descriptor for the dense vectors in y=op(A)x */ 207 #endif 208 }; 209 210 /* This is struct holding the relevant data needed to a MatMult */ 211 struct Mat_SeqAIJCUSPARSEMultStruct { 212 void *mat; /* opaque pointer to a matrix. This could be either a cusparseHybMat_t or a CsrMatrix */ 213 cusparseMatDescr_t descr; /* Data needed to describe the matrix for a multiply */ 214 THRUSTINTARRAY *cprowIndices; /* compressed row indices used in the parallel SpMV */ 215 PetscScalar *alpha_one; /* pointer to a device "scalar" storing the alpha parameter in the SpMV */ 216 PetscScalar *beta_zero; /* pointer to a device "scalar" storing the beta parameter in the SpMV as zero*/ 217 PetscScalar *beta_one; /* pointer to a device "scalar" storing the beta parameter in the SpMV as one */ 218 #if PETSC_PKG_CUDA_VERSION_GE(11,0,0) 219 cusparseSpMatDescr_t matDescr; /* descriptor for the matrix, used by SpMV and SpMM */ 220 Mat_CusparseSpMV cuSpMV[3]; /* different Mat_CusparseSpMV structs for non-transpose, transpose, conj-transpose */ 221 Mat_SeqAIJCUSPARSEMultStruct() : matDescr(NULL) { 222 for (int i=0; i<3; i++) cuSpMV[i].initialized = PETSC_FALSE; 223 } 224 #endif 225 }; 226 227 /* This is a larger struct holding all the matrices for a SpMV, and SpMV Tranpose */ 228 struct Mat_SeqAIJCUSPARSE { 229 Mat_SeqAIJCUSPARSEMultStruct *mat; /* pointer to the matrix on the GPU */ 230 Mat_SeqAIJCUSPARSEMultStruct *matTranspose; /* pointer to the matrix on the GPU (for the transpose ... useful for BiCG) */ 231 THRUSTARRAY *workVector; /* pointer to a workvector to which we can copy the relevant indices of a vector we want to multiply */ 232 THRUSTINTARRAY32 *rowoffsets_gpu; /* rowoffsets on GPU in non-compressed-row format. It is used to convert CSR to CSC */ 233 PetscInt nrows; /* number of rows of the matrix seen by GPU */ 234 MatCUSPARSEStorageFormat format; /* the storage format for the matrix on the device */ 235 cudaStream_t stream; /* a stream for the parallel SpMV ... this is not owned and should not be deleted */ 236 cusparseHandle_t handle; /* a handle to the cusparse library ... this may not be owned (if we're working in parallel i.e. multiGPUs) */ 237 PetscObjectState nonzerostate; /* track nonzero state to possibly recreate the GPU matrix */ 238 PetscBool transgen; /* whether or not to generate explicit transpose for MatMultTranspose operations */ 239 #if PETSC_PKG_CUDA_VERSION_GE(11,0,0) 240 size_t csr2cscBufferSize; /* stuff used to compute the matTranspose above */ 241 void *csr2cscBuffer; /* This is used as a C struct and is calloc'ed by PetscNewLog() */ 242 cusparseCsr2CscAlg_t csr2cscAlg; /* algorithms can be selected from command line options */ 243 cusparseSpMVAlg_t spmvAlg; 244 cusparseSpMMAlg_t spmmAlg; 245 #endif 246 PetscSplitCSRDataStructure *deviceMat; /* Matrix on device for, eg, assembly */ 247 THRUSTINTARRAY *cooPerm; 248 THRUSTINTARRAY *cooPerm_a; 249 }; 250 251 PETSC_INTERN PetscErrorCode MatCUSPARSECopyToGPU(Mat); 252 PETSC_INTERN PetscErrorCode MatCUSPARSESetStream(Mat, const cudaStream_t stream); 253 PETSC_INTERN PetscErrorCode MatCUSPARSESetHandle(Mat, const cusparseHandle_t handle); 254 PETSC_INTERN PetscErrorCode MatCUSPARSEClearHandle(Mat); 255 PETSC_INTERN PetscErrorCode MatSetPreallocationCOO_SeqAIJCUSPARSE(Mat,PetscInt,const PetscInt[],const PetscInt[]); 256 PETSC_INTERN PetscErrorCode MatSetValuesCOO_SeqAIJCUSPARSE(Mat,const PetscScalar[],InsertMode); 257 PETSC_INTERN PetscErrorCode MatSeqAIJCUSPARSEGetArrayRead(Mat,const PetscScalar**); 258 PETSC_INTERN PetscErrorCode MatSeqAIJCUSPARSERestoreArrayRead(Mat,const PetscScalar**); 259 PETSC_INTERN PetscErrorCode MatSeqAIJCUSPARSEGetArrayWrite(Mat,PetscScalar**); 260 PETSC_INTERN PetscErrorCode MatSeqAIJCUSPARSERestoreArrayWrite(Mat,PetscScalar**); 261 PETSC_INTERN PetscErrorCode MatSeqAIJCUSPARSEMergeMats(Mat,Mat,MatReuse,Mat*); 262 263 PETSC_STATIC_INLINE bool isCudaMem(const void *data) 264 { 265 cudaError_t cerr; 266 struct cudaPointerAttributes attr; 267 enum cudaMemoryType mtype; 268 cerr = cudaPointerGetAttributes(&attr,data); /* Do not check error since before CUDA 11.0, passing a host pointer returns cudaErrorInvalidValue */ 269 cudaGetLastError(); /* Reset the last error */ 270 #if (CUDART_VERSION < 10000) 271 mtype = attr.memoryType; 272 #else 273 mtype = attr.type; 274 #endif 275 if (cerr == cudaSuccess && mtype == cudaMemoryTypeDevice) return true; 276 else return false; 277 } 278 279 #endif 280