1519f805aSKarl Rupp #if !defined(__CUSPARSEMATIMPL) 29ae82921SPaul Mullowney #define __CUSPARSEMATIMPL 39ae82921SPaul Mullowney 4c41cb2e2SAlejandro Lamas Daviña #include <../src/vec/vec/impls/seq/seqcuda/cudavecimpl.h> 59ae82921SPaul Mullowney 62692e278SPaul Mullowney #if CUDA_VERSION>=4020 79ae82921SPaul Mullowney #include <cusparse_v2.h> 82692e278SPaul Mullowney #else 92692e278SPaul Mullowney #include <cusparse.h> 102692e278SPaul Mullowney #endif 119ae82921SPaul Mullowney 129ae82921SPaul Mullowney #include <algorithm> 139ae82921SPaul Mullowney #include <vector> 149ae82921SPaul Mullowney 15c41cb2e2SAlejandro Lamas Daviña #include <thrust/device_vector.h> 16c41cb2e2SAlejandro Lamas Daviña #include <thrust/device_ptr.h> 17*50ab121bSSatish Balay #include <thrust/device_malloc_allocator.h> 18c41cb2e2SAlejandro Lamas Daviña #include <thrust/transform.h> 19c41cb2e2SAlejandro Lamas Daviña #include <thrust/functional.h> 20554b8892SKarl Rupp #include <thrust/sequence.h> 21c41cb2e2SAlejandro Lamas Daviña 22aa372e3fSPaul Mullowney #if defined(PETSC_USE_COMPLEX) 23aa372e3fSPaul Mullowney #if defined(PETSC_USE_REAL_SINGLE) 24ec42abe4SAlejandro 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)) 25ec42abe4SAlejandro 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)) 26ec42abe4SAlejandro 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)) 27ec42abe4SAlejandro 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)) 28ec42abe4SAlejandro 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)) 29ec42abe4SAlejandro 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)) 30ec42abe4SAlejandro Lamas Daviña #define cusparse_hyb2csr(a,b,c,d,e,f) cusparseChyb2csr((a),(b),(c),(cuComplex*)(d),(e),(f)) 317656d835SStefano Zampini cuFloatComplex PETSC_CUSPARSE_ONE = {1.0f, 0.0f}; 327656d835SStefano Zampini cuFloatComplex PETSC_CUSPARSE_ZERO = {0.0f, 0.0f}; 33aa372e3fSPaul Mullowney #elif defined(PETSC_USE_REAL_DOUBLE) 34ec42abe4SAlejandro 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)) 35ec42abe4SAlejandro 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)) 36ec42abe4SAlejandro 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)) 37ec42abe4SAlejandro 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)) 38ec42abe4SAlejandro 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)) 39ec42abe4SAlejandro 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)) 40ec42abe4SAlejandro Lamas Daviña #define cusparse_hyb2csr(a,b,c,d,e,f) cusparseZhyb2csr((a),(b),(c),(cuDoubleComplex*)(d),(e),(f)) 417656d835SStefano Zampini cuDoubleComplex PETSC_CUSPARSE_ONE = {1.0, 0.0}; 427656d835SStefano Zampini cuDoubleComplex PETSC_CUSPARSE_ZERO = {0.0, 0.0}; 43aa372e3fSPaul Mullowney #endif 44aa372e3fSPaul Mullowney #else 457656d835SStefano Zampini PetscScalar PETSC_CUSPARSE_ONE = 1.0; 467656d835SStefano Zampini PetscScalar PETSC_CUSPARSE_ZERO = 0.0; 47aa372e3fSPaul Mullowney #if defined(PETSC_USE_REAL_SINGLE) 48aa372e3fSPaul Mullowney #define cusparse_solve cusparseScsrsv_solve 49aa372e3fSPaul Mullowney #define cusparse_analysis cusparseScsrsv_analysis 50aa372e3fSPaul Mullowney #define cusparse_csr_spmv cusparseScsrmv 51aa372e3fSPaul Mullowney #define cusparse_csr2csc cusparseScsr2csc 52aa372e3fSPaul Mullowney #define cusparse_hyb_spmv cusparseShybmv 53aa372e3fSPaul Mullowney #define cusparse_csr2hyb cusparseScsr2hyb 54aa372e3fSPaul Mullowney #define cusparse_hyb2csr cusparseShyb2csr 55aa372e3fSPaul Mullowney #elif defined(PETSC_USE_REAL_DOUBLE) 56aa372e3fSPaul Mullowney #define cusparse_solve cusparseDcsrsv_solve 57aa372e3fSPaul Mullowney #define cusparse_analysis cusparseDcsrsv_analysis 58aa372e3fSPaul Mullowney #define cusparse_csr_spmv cusparseDcsrmv 59aa372e3fSPaul Mullowney #define cusparse_csr2csc cusparseDcsr2csc 60aa372e3fSPaul Mullowney #define cusparse_hyb_spmv cusparseDhybmv 61aa372e3fSPaul Mullowney #define cusparse_csr2hyb cusparseDcsr2hyb 62aa372e3fSPaul Mullowney #define cusparse_hyb2csr cusparseDhyb2csr 63aa372e3fSPaul Mullowney #endif 64aa372e3fSPaul Mullowney #endif 65aa372e3fSPaul Mullowney 66aa372e3fSPaul Mullowney #define THRUSTINTARRAY32 thrust::device_vector<int> 67aa372e3fSPaul Mullowney #define THRUSTINTARRAY thrust::device_vector<PetscInt> 68aa372e3fSPaul Mullowney #define THRUSTARRAY thrust::device_vector<PetscScalar> 69aa372e3fSPaul Mullowney 70aa372e3fSPaul Mullowney /* A CSR matrix structure */ 71aa372e3fSPaul Mullowney struct CsrMatrix { 72aa372e3fSPaul Mullowney PetscInt num_rows; 73aa372e3fSPaul Mullowney PetscInt num_cols; 74aa372e3fSPaul Mullowney PetscInt num_entries; 75aa372e3fSPaul Mullowney THRUSTINTARRAY32 *row_offsets; 76aa372e3fSPaul Mullowney THRUSTINTARRAY32 *column_indices; 77aa372e3fSPaul Mullowney THRUSTARRAY *values; 789ae82921SPaul Mullowney }; 799ae82921SPaul Mullowney 80aa372e3fSPaul Mullowney //#define CUSPMATRIXCSR32 cusp::csr_matrix<int,PetscScalar,cusp::device_memory> 81aa372e3fSPaul Mullowney 82aa372e3fSPaul Mullowney /* This is struct holding the relevant data needed to a MatSolve */ 83aa372e3fSPaul Mullowney struct Mat_SeqAIJCUSPARSETriFactorStruct { 84aa372e3fSPaul Mullowney /* Data needed for triangular solve */ 85aa372e3fSPaul Mullowney cusparseMatDescr_t descr; 86aa372e3fSPaul Mullowney cusparseSolveAnalysisInfo_t solveInfo; 87aa372e3fSPaul Mullowney cusparseOperation_t solveOp; 88aa372e3fSPaul Mullowney CsrMatrix *csrMat; 89aa372e3fSPaul Mullowney }; 90aa372e3fSPaul Mullowney 91aa372e3fSPaul Mullowney /* This is struct holding the relevant data needed to a MatMult */ 92aa372e3fSPaul Mullowney struct Mat_SeqAIJCUSPARSEMultStruct { 93aa372e3fSPaul Mullowney void *mat; /* opaque pointer to a matrix. This could be either a cusparseHybMat_t or a CsrMatrix */ 94aa372e3fSPaul Mullowney cusparseMatDescr_t descr; /* Data needed to describe the matrix for a multiply */ 95aa372e3fSPaul Mullowney THRUSTINTARRAY *cprowIndices; /* compressed row indices used in the parallel SpMV */ 96b06137fdSPaul Mullowney PetscScalar *alpha; /* pointer to a device "scalar" storing the alpha parameter in the SpMV */ 977656d835SStefano Zampini PetscScalar *beta_zero; /* pointer to a device "scalar" storing the beta parameter in the SpMV as zero*/ 987656d835SStefano Zampini PetscScalar *beta_one; /* pointer to a device "scalar" storing the beta parameter in the SpMV as one */ 99aa372e3fSPaul Mullowney }; 100aa372e3fSPaul Mullowney 101aa372e3fSPaul Mullowney /* This is a larger struct holding all the triangular factors for a solve, transpose solve, and 102aa372e3fSPaul Mullowney any indices used in a reordering */ 103aa372e3fSPaul Mullowney struct Mat_SeqAIJCUSPARSETriFactors { 104aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorPtr; /* pointer for lower triangular (factored matrix) on GPU */ 105aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorPtr; /* pointer for upper triangular (factored matrix) on GPU */ 106aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorPtrTranspose; /* pointer for lower triangular (factored matrix) on GPU for the transpose (useful for BiCG) */ 107aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorPtrTranspose; /* pointer for upper triangular (factored matrix) on GPU for the transpose (useful for BiCG)*/ 108aa372e3fSPaul Mullowney THRUSTINTARRAY *rpermIndices; /* indices used for any reordering */ 109aa372e3fSPaul Mullowney THRUSTINTARRAY *cpermIndices; /* indices used for any reordering */ 110aa372e3fSPaul Mullowney THRUSTARRAY *workVector; 111aa372e3fSPaul Mullowney cusparseHandle_t handle; /* a handle to the cusparse library */ 112aa372e3fSPaul Mullowney PetscInt nnz; /* number of nonzeros ... need this for accurate logging between ICC and ILU */ 113aa372e3fSPaul Mullowney }; 114aa372e3fSPaul Mullowney 115aa372e3fSPaul Mullowney /* This is a larger struct holding all the matrices for a SpMV, and SpMV Tranpose */ 1169ae82921SPaul Mullowney struct Mat_SeqAIJCUSPARSE { 117aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSEMultStruct *mat; /* pointer to the matrix on the GPU */ 118aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSEMultStruct *matTranspose; /* pointer to the matrix on the GPU (for the transpose ... useful for BiCG) */ 119aa372e3fSPaul Mullowney THRUSTARRAY *workVector; /*pointer to a workvector to which we can copy the relevant indices of a vector we want to multiply */ 1209ae82921SPaul Mullowney PetscInt nonzerorow; /* number of nonzero rows ... used in the flop calculations */ 121e057df02SPaul Mullowney MatCUSPARSEStorageFormat format; /* the storage format for the matrix on the device */ 122aa372e3fSPaul Mullowney cudaStream_t stream; /* a stream for the parallel SpMV ... this is not owned and should not be deleted */ 123aa372e3fSPaul Mullowney cusparseHandle_t handle; /* a handle to the cusparse library ... this may not be owned (if we're working in parallel i.e. multiGPUs) */ 12434d6c7a5SJose E. Roman PetscObjectState nonzerostate; 1259ae82921SPaul Mullowney }; 1269ae82921SPaul Mullowney 1275a576424SJed Brown PETSC_INTERN PetscErrorCode MatCUSPARSECopyToGPU(Mat); 128b06137fdSPaul Mullowney PETSC_INTERN PetscErrorCode MatCUSPARSESetStream(Mat, const cudaStream_t stream); 129b06137fdSPaul Mullowney PETSC_INTERN PetscErrorCode MatCUSPARSESetHandle(Mat, const cusparseHandle_t handle); 130b06137fdSPaul Mullowney PETSC_INTERN PetscErrorCode MatCUSPARSEClearHandle(Mat); 1319ae82921SPaul Mullowney #endif 132