xref: /petsc/src/mat/impls/aij/seq/seqcusparse/cusparsematimpl.h (revision 213423ffc686f176624d36698af1d5aff96e941b)
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 
69ae82921SPaul Mullowney #include <cusparse_v2.h>
79ae82921SPaul Mullowney 
89ae82921SPaul Mullowney #include <algorithm>
99ae82921SPaul Mullowney #include <vector>
109ae82921SPaul Mullowney 
11c41cb2e2SAlejandro Lamas Daviña #include <thrust/device_vector.h>
12c41cb2e2SAlejandro Lamas Daviña #include <thrust/device_ptr.h>
1350ab121bSSatish Balay #include <thrust/device_malloc_allocator.h>
14c41cb2e2SAlejandro Lamas Daviña #include <thrust/transform.h>
15c41cb2e2SAlejandro Lamas Daviña #include <thrust/functional.h>
16554b8892SKarl Rupp #include <thrust/sequence.h>
17c41cb2e2SAlejandro Lamas Daviña 
18185e5899SJunchao Zhang #if (CUSPARSE_VER_MAJOR > 10 || CUSPARSE_VER_MAJOR == 10 && CUSPARSE_VER_MINOR >= 2) /* According to cuda/10.1.168 on OLCF Summit */
1957d48284SJunchao Zhang #define CHKERRCUSPARSE(stat) \
2057d48284SJunchao Zhang do { \
2157d48284SJunchao Zhang    if (PetscUnlikely(stat)) { \
2257d48284SJunchao Zhang       const char *name  = cusparseGetErrorName(stat); \
2357d48284SJunchao Zhang       const char *descr = cusparseGetErrorString(stat); \
2457d48284SJunchao Zhang       SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_LIB,"cuSPARSE error %d (%s) : %s",(int)stat,name,descr); \
2557d48284SJunchao Zhang    } \
2657d48284SJunchao Zhang } while(0)
276d22fe3dSJunchao Zhang #else
28ce29e4ddSJunchao Zhang #define CHKERRCUSPARSE(stat) do {if (PetscUnlikely(stat)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"cusparse error %d",(int)stat);} while(0)
296d22fe3dSJunchao Zhang #endif
3057d48284SJunchao Zhang 
31aa372e3fSPaul Mullowney #if defined(PETSC_USE_COMPLEX)
32aa372e3fSPaul Mullowney #if defined(PETSC_USE_REAL_SINGLE)
33ec42abe4SAlejandro 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))
34ec42abe4SAlejandro 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))
35ec42abe4SAlejandro 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))
36ec42abe4SAlejandro 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))
37ec42abe4SAlejandro 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))
38ec42abe4SAlejandro 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))
39ec42abe4SAlejandro Lamas Daviña #define cusparse_hyb2csr(a,b,c,d,e,f)                cusparseChyb2csr((a),(b),(c),(cuComplex*)(d),(e),(f))
4057d48284SJunchao Zhang const cuFloatComplex PETSC_CUSPARSE_ONE  = {1.0f, 0.0f};
4157d48284SJunchao Zhang const cuFloatComplex PETSC_CUSPARSE_ZERO = {0.0f, 0.0f};
42aa372e3fSPaul Mullowney #elif defined(PETSC_USE_REAL_DOUBLE)
43ec42abe4SAlejandro 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))
44ec42abe4SAlejandro 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))
45ec42abe4SAlejandro 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))
46ec42abe4SAlejandro 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))
47ec42abe4SAlejandro 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))
48ec42abe4SAlejandro 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))
49ec42abe4SAlejandro Lamas Daviña #define cusparse_hyb2csr(a,b,c,d,e,f)                cusparseZhyb2csr((a),(b),(c),(cuDoubleComplex*)(d),(e),(f))
5057d48284SJunchao Zhang const cuDoubleComplex PETSC_CUSPARSE_ONE  = {1.0, 0.0};
5157d48284SJunchao Zhang const cuDoubleComplex PETSC_CUSPARSE_ZERO = {0.0, 0.0};
52aa372e3fSPaul Mullowney #endif
53aa372e3fSPaul Mullowney #else
5457d48284SJunchao Zhang const PetscScalar PETSC_CUSPARSE_ONE  = 1.0;
5557d48284SJunchao Zhang const PetscScalar PETSC_CUSPARSE_ZERO = 0.0;
56aa372e3fSPaul Mullowney #if defined(PETSC_USE_REAL_SINGLE)
57aa372e3fSPaul Mullowney #define cusparse_solve    cusparseScsrsv_solve
58aa372e3fSPaul Mullowney #define cusparse_analysis cusparseScsrsv_analysis
59aa372e3fSPaul Mullowney #define cusparse_csr_spmv cusparseScsrmv
60aa372e3fSPaul Mullowney #define cusparse_csr2csc  cusparseScsr2csc
61aa372e3fSPaul Mullowney #define cusparse_hyb_spmv cusparseShybmv
62aa372e3fSPaul Mullowney #define cusparse_csr2hyb  cusparseScsr2hyb
63aa372e3fSPaul Mullowney #define cusparse_hyb2csr  cusparseShyb2csr
64aa372e3fSPaul Mullowney #elif defined(PETSC_USE_REAL_DOUBLE)
65aa372e3fSPaul Mullowney #define cusparse_solve    cusparseDcsrsv_solve
66aa372e3fSPaul Mullowney #define cusparse_analysis cusparseDcsrsv_analysis
67aa372e3fSPaul Mullowney #define cusparse_csr_spmv cusparseDcsrmv
68aa372e3fSPaul Mullowney #define cusparse_csr2csc  cusparseDcsr2csc
69aa372e3fSPaul Mullowney #define cusparse_hyb_spmv cusparseDhybmv
70aa372e3fSPaul Mullowney #define cusparse_csr2hyb  cusparseDcsr2hyb
71aa372e3fSPaul Mullowney #define cusparse_hyb2csr  cusparseDhyb2csr
72aa372e3fSPaul Mullowney #endif
73aa372e3fSPaul Mullowney #endif
74aa372e3fSPaul Mullowney 
75aa372e3fSPaul Mullowney #define THRUSTINTARRAY32 thrust::device_vector<int>
76aa372e3fSPaul Mullowney #define THRUSTINTARRAY thrust::device_vector<PetscInt>
77aa372e3fSPaul Mullowney #define THRUSTARRAY thrust::device_vector<PetscScalar>
78aa372e3fSPaul Mullowney 
79aa372e3fSPaul Mullowney /* A CSR matrix structure */
80aa372e3fSPaul Mullowney struct CsrMatrix {
81aa372e3fSPaul Mullowney   PetscInt         num_rows;
82aa372e3fSPaul Mullowney   PetscInt         num_cols;
83aa372e3fSPaul Mullowney   PetscInt         num_entries;
84aa372e3fSPaul Mullowney   THRUSTINTARRAY32 *row_offsets;
85aa372e3fSPaul Mullowney   THRUSTINTARRAY32 *column_indices;
86aa372e3fSPaul Mullowney   THRUSTARRAY      *values;
879ae82921SPaul Mullowney };
889ae82921SPaul Mullowney 
89aa372e3fSPaul Mullowney //#define CUSPMATRIXCSR32 cusp::csr_matrix<int,PetscScalar,cusp::device_memory>
90aa372e3fSPaul Mullowney 
91aa372e3fSPaul Mullowney /* This is struct holding the relevant data needed to a MatSolve */
92aa372e3fSPaul Mullowney struct Mat_SeqAIJCUSPARSETriFactorStruct {
93aa372e3fSPaul Mullowney   /* Data needed for triangular solve */
94aa372e3fSPaul Mullowney   cusparseMatDescr_t          descr;
95aa372e3fSPaul Mullowney   cusparseSolveAnalysisInfo_t solveInfo;
96aa372e3fSPaul Mullowney   cusparseOperation_t         solveOp;
97aa372e3fSPaul Mullowney   CsrMatrix                   *csrMat;
98aa372e3fSPaul Mullowney };
99aa372e3fSPaul Mullowney 
100aa372e3fSPaul Mullowney /* This is struct holding the relevant data needed to a MatMult */
101aa372e3fSPaul Mullowney struct Mat_SeqAIJCUSPARSEMultStruct {
102aa372e3fSPaul Mullowney   void               *mat;  /* opaque pointer to a matrix. This could be either a cusparseHybMat_t or a CsrMatrix */
103aa372e3fSPaul Mullowney   cusparseMatDescr_t descr; /* Data needed to describe the matrix for a multiply */
104aa372e3fSPaul Mullowney   THRUSTINTARRAY     *cprowIndices;   /* compressed row indices used in the parallel SpMV */
105b06137fdSPaul Mullowney   PetscScalar        *alpha; /* pointer to a device "scalar" storing the alpha parameter in the SpMV */
1067656d835SStefano Zampini   PetscScalar        *beta_zero; /* pointer to a device "scalar" storing the beta parameter in the SpMV as zero*/
1077656d835SStefano Zampini   PetscScalar        *beta_one; /* pointer to a device "scalar" storing the beta parameter in the SpMV as one */
108aa372e3fSPaul Mullowney };
109aa372e3fSPaul Mullowney 
110aa372e3fSPaul Mullowney /* This is a larger struct holding all the triangular factors for a solve, transpose solve, and
111aa372e3fSPaul Mullowney  any indices used in a reordering */
112aa372e3fSPaul Mullowney struct Mat_SeqAIJCUSPARSETriFactors {
113aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorPtr; /* pointer for lower triangular (factored matrix) on GPU */
114aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorPtr; /* pointer for upper triangular (factored matrix) on GPU */
115aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorPtrTranspose; /* pointer for lower triangular (factored matrix) on GPU for the transpose (useful for BiCG) */
116aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorPtrTranspose; /* pointer for upper triangular (factored matrix) on GPU for the transpose (useful for BiCG)*/
117aa372e3fSPaul Mullowney   THRUSTINTARRAY                    *rpermIndices;  /* indices used for any reordering */
118aa372e3fSPaul Mullowney   THRUSTINTARRAY                    *cpermIndices;  /* indices used for any reordering */
119aa372e3fSPaul Mullowney   THRUSTARRAY                       *workVector;
120aa372e3fSPaul Mullowney   cusparseHandle_t                  handle;   /* a handle to the cusparse library */
121aa372e3fSPaul Mullowney   PetscInt                          nnz;      /* number of nonzeros ... need this for accurate logging between ICC and ILU */
122aa372e3fSPaul Mullowney };
123aa372e3fSPaul Mullowney 
124aa372e3fSPaul Mullowney /* This is a larger struct holding all the matrices for a SpMV, and SpMV Tranpose */
1259ae82921SPaul Mullowney struct Mat_SeqAIJCUSPARSE {
126aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSEMultStruct *mat; /* pointer to the matrix on the GPU */
127aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSEMultStruct *matTranspose; /* pointer to the matrix on the GPU (for the transpose ... useful for BiCG) */
128aa372e3fSPaul Mullowney   THRUSTARRAY                  *workVector; /*pointer to a workvector to which we can copy the relevant indices of a vector we want to multiply */
12981902715SJunchao Zhang   THRUSTINTARRAY               *rowoffsets_gpu; /* rowoffsets on GPU in non-compressed-row format. It is used to convert CSR to CSC */
130*213423ffSJunchao Zhang   PetscInt                     nrows;    /* number of rows of the matrix seen by GPU */
131e057df02SPaul Mullowney   MatCUSPARSEStorageFormat     format;   /* the storage format for the matrix on the device */
132aa372e3fSPaul Mullowney   cudaStream_t                 stream;   /* a stream for the parallel SpMV ... this is not owned and should not be deleted */
133aa372e3fSPaul Mullowney   cusparseHandle_t             handle;   /* a handle to the cusparse library ... this may not be owned (if we're working in parallel i.e. multiGPUs) */
13434d6c7a5SJose E. Roman   PetscObjectState             nonzerostate;
1359ae82921SPaul Mullowney };
1369ae82921SPaul Mullowney 
1375a576424SJed Brown PETSC_INTERN PetscErrorCode MatCUSPARSECopyToGPU(Mat);
138b06137fdSPaul Mullowney PETSC_INTERN PetscErrorCode MatCUSPARSESetStream(Mat, const cudaStream_t stream);
139b06137fdSPaul Mullowney PETSC_INTERN PetscErrorCode MatCUSPARSESetHandle(Mat, const cusparseHandle_t handle);
140b06137fdSPaul Mullowney PETSC_INTERN PetscErrorCode MatCUSPARSEClearHandle(Mat);
1419ae82921SPaul Mullowney #endif
142