4: #include <../src/vec/vec/impls/seq/seqcuda/cudavecimpl.h> 6: #include <cusparse_v2.h>
8: #include <algorithm>
9: #include <vector>
11: #include <thrust/device_vector.h>
12: #include <thrust/device_ptr.h>
13: #include <thrust/device_malloc_allocator.h>
14: #include <thrust/transform.h>
15: #include <thrust/functional.h>
16: #include <thrust/sequence.h>
18: #if defined(PETSC_USE_COMPLEX)
19: #if defined(PETSC_USE_REAL_SINGLE)
20: #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)) 21: #define cusparse_analysis(a,b,c,d,e,f,g,h,i) cusparseCcsrsv_analysis((a),(b),(c),(d),(e),(cuComplex*)(f),(g),(h),(i)) 22: #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)) 23: #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)) 24: #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)) 25: #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)) 26: #define cusparse_hyb2csr(a,b,c,d,e,f) cusparseChyb2csr((a),(b),(c),(cuComplex*)(d),(e),(f)) 27: cuFloatComplex PETSC_CUSPARSE_ONE = {1.0f, 0.0f};
28: cuFloatComplex PETSC_CUSPARSE_ZERO = {0.0f, 0.0f};
29: #elif defined(PETSC_USE_REAL_DOUBLE)
30: #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)) 31: #define cusparse_analysis(a,b,c,d,e,f,g,h,i) cusparseZcsrsv_analysis((a),(b),(c),(d),(e),(cuDoubleComplex*)(f),(g),(h),(i)) 32: #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)) 33: #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)) 34: #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)) 35: #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)) 36: #define cusparse_hyb2csr(a,b,c,d,e,f) cusparseZhyb2csr((a),(b),(c),(cuDoubleComplex*)(d),(e),(f)) 37: cuDoubleComplex PETSC_CUSPARSE_ONE = {1.0, 0.0};
38: cuDoubleComplex PETSC_CUSPARSE_ZERO = {0.0, 0.0};
39: #endif
40: #else
41: PetscScalar PETSC_CUSPARSE_ONE = 1.0;
42: PetscScalar PETSC_CUSPARSE_ZERO = 0.0;
43: #if defined(PETSC_USE_REAL_SINGLE)
44: #define cusparse_solve cusparseScsrsv_solve 45: #define cusparse_analysis cusparseScsrsv_analysis 46: #define cusparse_csr_spmv cusparseScsrmv 47: #define cusparse_csr2csc cusparseScsr2csc 48: #define cusparse_hyb_spmv cusparseShybmv 49: #define cusparse_csr2hyb cusparseScsr2hyb 50: #define cusparse_hyb2csr cusparseShyb2csr 51: #elif defined(PETSC_USE_REAL_DOUBLE)
52: #define cusparse_solve cusparseDcsrsv_solve 53: #define cusparse_analysis cusparseDcsrsv_analysis 54: #define cusparse_csr_spmv cusparseDcsrmv 55: #define cusparse_csr2csc cusparseDcsr2csc 56: #define cusparse_hyb_spmv cusparseDhybmv 57: #define cusparse_csr2hyb cusparseDcsr2hyb 58: #define cusparse_hyb2csr cusparseDhyb2csr 59: #endif
60: #endif
62: #define THRUSTINTARRAY32 thrust::device_vector<int> 63: #define THRUSTINTARRAY thrust::device_vector<PetscInt> 64: #define THRUSTARRAY thrust::device_vector<PetscScalar> 66: /* A CSR matrix structure */
67: struct CsrMatrix {
68: PetscInt num_rows;
69: PetscInt num_cols;
70: PetscInt num_entries;
71: THRUSTINTARRAY32 *row_offsets;
72: THRUSTINTARRAY32 *column_indices;
73: THRUSTARRAY *values;
74: };
76: //#define CUSPMATRIXCSR32 cusp::csr_matrix<int,PetscScalar,cusp::device_memory>
78: /* This is struct holding the relevant data needed to a MatSolve */
79: struct Mat_SeqAIJCUSPARSETriFactorStruct {
80: /* Data needed for triangular solve */
81: cusparseMatDescr_t descr;
82: cusparseSolveAnalysisInfo_t solveInfo;
83: cusparseOperation_t solveOp;
84: CsrMatrix *csrMat;
85: };
87: /* This is struct holding the relevant data needed to a MatMult */
88: struct Mat_SeqAIJCUSPARSEMultStruct {
89: void *mat; /* opaque pointer to a matrix. This could be either a cusparseHybMat_t or a CsrMatrix */
90: cusparseMatDescr_t descr; /* Data needed to describe the matrix for a multiply */
91: THRUSTINTARRAY *cprowIndices; /* compressed row indices used in the parallel SpMV */
92: PetscScalar *alpha; /* pointer to a device "scalar" storing the alpha parameter in the SpMV */
93: PetscScalar *beta_zero; /* pointer to a device "scalar" storing the beta parameter in the SpMV as zero*/
94: PetscScalar *beta_one; /* pointer to a device "scalar" storing the beta parameter in the SpMV as one */
95: };
97: /* This is a larger struct holding all the triangular factors for a solve, transpose solve, and
98: any indices used in a reordering */
99: struct Mat_SeqAIJCUSPARSETriFactors {
100: Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorPtr; /* pointer for lower triangular (factored matrix) on GPU */
101: Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorPtr; /* pointer for upper triangular (factored matrix) on GPU */
102: Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorPtrTranspose; /* pointer for lower triangular (factored matrix) on GPU for the transpose (useful for BiCG) */
103: Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorPtrTranspose; /* pointer for upper triangular (factored matrix) on GPU for the transpose (useful for BiCG)*/
104: THRUSTINTARRAY *rpermIndices; /* indices used for any reordering */
105: THRUSTINTARRAY *cpermIndices; /* indices used for any reordering */
106: THRUSTARRAY *workVector;
107: cusparseHandle_t handle; /* a handle to the cusparse library */
108: PetscInt nnz; /* number of nonzeros ... need this for accurate logging between ICC and ILU */
109: };
111: /* This is a larger struct holding all the matrices for a SpMV, and SpMV Tranpose */
112: struct Mat_SeqAIJCUSPARSE {
113: Mat_SeqAIJCUSPARSEMultStruct *mat; /* pointer to the matrix on the GPU */
114: Mat_SeqAIJCUSPARSEMultStruct *matTranspose; /* pointer to the matrix on the GPU (for the transpose ... useful for BiCG) */
115: THRUSTARRAY *workVector; /*pointer to a workvector to which we can copy the relevant indices of a vector we want to multiply */
116: THRUSTINTARRAY *rowoffsets_gpu; /* rowoffsets on GPU in non-compressed-row format. It is used to convert CSR to CSC */
117: PetscInt nonzerorow; /* number of nonzero rows ... used in the flop calculations */
118: MatCUSPARSEStorageFormat format; /* the storage format for the matrix on the device */
119: cudaStream_t stream; /* a stream for the parallel SpMV ... this is not owned and should not be deleted */
120: cusparseHandle_t handle; /* a handle to the cusparse library ... this may not be owned (if we're working in parallel i.e. multiGPUs) */
121: PetscObjectState nonzerostate;
122: };
124: PETSC_INTERN PetscErrorCode MatCUSPARSECopyToGPU(Mat);
125: PETSC_INTERN PetscErrorCode MatCUSPARSESetStream(Mat, const cudaStream_t stream);
126: PETSC_INTERN PetscErrorCode MatCUSPARSESetHandle(Mat, const cusparseHandle_t handle);
127: PETSC_INTERN PetscErrorCode MatCUSPARSEClearHandle(Mat);
128: #endif