xref: /petsc/src/mat/impls/aij/seq/seqcusparse/cusparsematimpl.h (revision 54da937ac53361e7377cf6625207cd5b3a0dbe8c)
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)
33ccdfe979SStefano Zampini const cuComplex PETSC_CUSPARSE_ONE  = {1.0f, 0.0f};
34ccdfe979SStefano Zampini const cuComplex PETSC_CUSPARSE_ZERO = {0.0f, 0.0f};
35ec42abe4SAlejandro 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))
36ec42abe4SAlejandro 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))
37ec42abe4SAlejandro 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))
38ccdfe979SStefano 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))
39ec42abe4SAlejandro 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))
40ec42abe4SAlejandro 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))
41ec42abe4SAlejandro 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))
42ec42abe4SAlejandro Lamas Daviña #define cusparse_hyb2csr(a,b,c,d,e,f)                      cusparseChyb2csr((a),(b),(c),(cuComplex*)(d),(e),(f))
43aa372e3fSPaul Mullowney #elif defined(PETSC_USE_REAL_DOUBLE)
44ccdfe979SStefano Zampini const cuDoubleComplex PETSC_CUSPARSE_ONE  = {1.0, 0.0};
45ccdfe979SStefano Zampini const cuDoubleComplex PETSC_CUSPARSE_ZERO = {0.0, 0.0};
46ec42abe4SAlejandro 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))
47ec42abe4SAlejandro 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))
48ec42abe4SAlejandro 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))
49ccdfe979SStefano 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))
50ec42abe4SAlejandro 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))
51ec42abe4SAlejandro 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))
52ec42abe4SAlejandro 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))
53ec42abe4SAlejandro Lamas Daviña #define cusparse_hyb2csr(a,b,c,d,e,f)                      cusparseZhyb2csr((a),(b),(c),(cuDoubleComplex*)(d),(e),(f))
54aa372e3fSPaul Mullowney #endif
55aa372e3fSPaul Mullowney #else
5657d48284SJunchao Zhang const PetscScalar PETSC_CUSPARSE_ONE  = 1.0;
5757d48284SJunchao Zhang const PetscScalar PETSC_CUSPARSE_ZERO = 0.0;
58aa372e3fSPaul Mullowney #if defined(PETSC_USE_REAL_SINGLE)
59aa372e3fSPaul Mullowney #define cusparse_solve    cusparseScsrsv_solve
60aa372e3fSPaul Mullowney #define cusparse_analysis cusparseScsrsv_analysis
61aa372e3fSPaul Mullowney #define cusparse_csr_spmv cusparseScsrmv
62ccdfe979SStefano Zampini #define cusparse_csr_spmm cusparseScsrmm
63aa372e3fSPaul Mullowney #define cusparse_csr2csc  cusparseScsr2csc
64aa372e3fSPaul Mullowney #define cusparse_hyb_spmv cusparseShybmv
65aa372e3fSPaul Mullowney #define cusparse_csr2hyb  cusparseScsr2hyb
66aa372e3fSPaul Mullowney #define cusparse_hyb2csr  cusparseShyb2csr
67aa372e3fSPaul Mullowney #elif defined(PETSC_USE_REAL_DOUBLE)
68aa372e3fSPaul Mullowney #define cusparse_solve    cusparseDcsrsv_solve
69aa372e3fSPaul Mullowney #define cusparse_analysis cusparseDcsrsv_analysis
70aa372e3fSPaul Mullowney #define cusparse_csr_spmv cusparseDcsrmv
71ccdfe979SStefano Zampini #define cusparse_csr_spmm cusparseDcsrmm
72aa372e3fSPaul Mullowney #define cusparse_csr2csc  cusparseDcsr2csc
73aa372e3fSPaul Mullowney #define cusparse_hyb_spmv cusparseDhybmv
74aa372e3fSPaul Mullowney #define cusparse_csr2hyb  cusparseDcsr2hyb
75aa372e3fSPaul Mullowney #define cusparse_hyb2csr  cusparseDhyb2csr
76aa372e3fSPaul Mullowney #endif
77aa372e3fSPaul Mullowney #endif
78aa372e3fSPaul Mullowney 
79aa372e3fSPaul Mullowney #define THRUSTINTARRAY32 thrust::device_vector<int>
80aa372e3fSPaul Mullowney #define THRUSTINTARRAY thrust::device_vector<PetscInt>
81aa372e3fSPaul Mullowney #define THRUSTARRAY thrust::device_vector<PetscScalar>
82aa372e3fSPaul Mullowney 
83aa372e3fSPaul Mullowney /* A CSR matrix structure */
84aa372e3fSPaul Mullowney struct CsrMatrix {
85aa372e3fSPaul Mullowney   PetscInt         num_rows;
86aa372e3fSPaul Mullowney   PetscInt         num_cols;
87aa372e3fSPaul Mullowney   PetscInt         num_entries;
88aa372e3fSPaul Mullowney   THRUSTINTARRAY32 *row_offsets;
89aa372e3fSPaul Mullowney   THRUSTINTARRAY32 *column_indices;
90aa372e3fSPaul Mullowney   THRUSTARRAY      *values;
919ae82921SPaul Mullowney };
929ae82921SPaul Mullowney 
93aa372e3fSPaul Mullowney /* This is struct holding the relevant data needed to a MatSolve */
94aa372e3fSPaul Mullowney struct Mat_SeqAIJCUSPARSETriFactorStruct {
95aa372e3fSPaul Mullowney   /* Data needed for triangular solve */
96aa372e3fSPaul Mullowney   cusparseMatDescr_t          descr;
97aa372e3fSPaul Mullowney   cusparseSolveAnalysisInfo_t solveInfo;
98aa372e3fSPaul Mullowney   cusparseOperation_t         solveOp;
99aa372e3fSPaul Mullowney   CsrMatrix                   *csrMat;
100aa372e3fSPaul Mullowney };
101aa372e3fSPaul Mullowney 
102aa372e3fSPaul Mullowney /* This is struct holding the relevant data needed to a MatMult */
103aa372e3fSPaul Mullowney struct Mat_SeqAIJCUSPARSEMultStruct {
104aa372e3fSPaul Mullowney   void               *mat;  /* opaque pointer to a matrix. This could be either a cusparseHybMat_t or a CsrMatrix */
105aa372e3fSPaul Mullowney   cusparseMatDescr_t descr; /* Data needed to describe the matrix for a multiply */
106aa372e3fSPaul Mullowney   THRUSTINTARRAY     *cprowIndices;   /* compressed row indices used in the parallel SpMV */
107b06137fdSPaul Mullowney   PetscScalar        *alpha; /* pointer to a device "scalar" storing the alpha parameter in the SpMV */
1087656d835SStefano Zampini   PetscScalar        *beta_zero; /* pointer to a device "scalar" storing the beta parameter in the SpMV as zero*/
1097656d835SStefano Zampini   PetscScalar        *beta_one; /* pointer to a device "scalar" storing the beta parameter in the SpMV as one */
110aa372e3fSPaul Mullowney };
111aa372e3fSPaul Mullowney 
112aa372e3fSPaul Mullowney /* This is a larger struct holding all the triangular factors for a solve, transpose solve, and
113aa372e3fSPaul Mullowney  any indices used in a reordering */
114aa372e3fSPaul Mullowney struct Mat_SeqAIJCUSPARSETriFactors {
115aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorPtr; /* pointer for lower triangular (factored matrix) on GPU */
116aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorPtr; /* pointer for upper triangular (factored matrix) on GPU */
117aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorPtrTranspose; /* pointer for lower triangular (factored matrix) on GPU for the transpose (useful for BiCG) */
118aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorPtrTranspose; /* pointer for upper triangular (factored matrix) on GPU for the transpose (useful for BiCG)*/
119aa372e3fSPaul Mullowney   THRUSTINTARRAY                    *rpermIndices;  /* indices used for any reordering */
120aa372e3fSPaul Mullowney   THRUSTINTARRAY                    *cpermIndices;  /* indices used for any reordering */
121aa372e3fSPaul Mullowney   THRUSTARRAY                       *workVector;
122aa372e3fSPaul Mullowney   cusparseHandle_t                  handle;   /* a handle to the cusparse library */
123aa372e3fSPaul Mullowney   PetscInt                          nnz;      /* number of nonzeros ... need this for accurate logging between ICC and ILU */
124aa372e3fSPaul Mullowney };
125aa372e3fSPaul Mullowney 
126aa372e3fSPaul Mullowney /* This is a larger struct holding all the matrices for a SpMV, and SpMV Tranpose */
1279ae82921SPaul Mullowney struct Mat_SeqAIJCUSPARSE {
128aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSEMultStruct *mat;            /* pointer to the matrix on the GPU */
129aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSEMultStruct *matTranspose;   /* pointer to the matrix on the GPU (for the transpose ... useful for BiCG) */
130aa372e3fSPaul Mullowney   THRUSTARRAY                  *workVector;     /* pointer to a workvector to which we can copy the relevant indices of a vector we want to multiply */
131029808eeSJunchao Zhang   THRUSTINTARRAY32             *rowoffsets_gpu; /* rowoffsets on GPU in non-compressed-row format. It is used to convert CSR to CSC */
132213423ffSJunchao Zhang   PetscInt                     nrows;           /* number of rows of the matrix seen by GPU */
133e057df02SPaul Mullowney   MatCUSPARSEStorageFormat     format;          /* the storage format for the matrix on the device */
134aa372e3fSPaul Mullowney   cudaStream_t                 stream;          /* a stream for the parallel SpMV ... this is not owned and should not be deleted */
135aa372e3fSPaul Mullowney   cusparseHandle_t             handle;          /* a handle to the cusparse library ... this may not be owned (if we're working in parallel i.e. multiGPUs) */
136*54da937aSStefano Zampini   PetscObjectState             nonzerostate;    /* track nonzero state to possibly recreate the GPU matrix */
137*54da937aSStefano Zampini   PetscBool                    transgen;        /* whether or not to generate explicit transpose for MatMultTranspose operations */
1389ae82921SPaul Mullowney };
1399ae82921SPaul Mullowney 
1405a576424SJed Brown PETSC_INTERN PetscErrorCode MatCUSPARSECopyToGPU(Mat);
141b06137fdSPaul Mullowney PETSC_INTERN PetscErrorCode MatCUSPARSESetStream(Mat, const cudaStream_t stream);
142b06137fdSPaul Mullowney PETSC_INTERN PetscErrorCode MatCUSPARSESetHandle(Mat, const cusparseHandle_t handle);
143b06137fdSPaul Mullowney PETSC_INTERN PetscErrorCode MatCUSPARSEClearHandle(Mat);
1449ae82921SPaul Mullowney #endif
145