xref: /petsc/src/mat/impls/aij/seq/seqcusparse/cusparsematimpl.h (revision aa372e3fe8df1a25e86bd398820fb3017c3d0016)
1 #if !defined(__CUSPARSEMATIMPL)
2 #define __CUSPARSEMATIMPL
3 
4 #include <../src/vec/vec/impls/seq/seqcusp/cuspvecimpl.h>
5 
6 #include <cusparse_v2.h>
7 
8 #include <algorithm>
9 #include <vector>
10 #include <thrust/sort.h>
11 #include <thrust/fill.h>
12 #include <cusp/csr_matrix.h>
13 
14 /* Single instance of the cusparse handle for the class. */
15 MatCUSPARSEStorageFormat cusparseMatSolveStorageFormat=MAT_CUSPARSE_CSR;
16 
17 #if defined(PETSC_USE_COMPLEX)
18 #if defined(PETSC_USE_REAL_SINGLE)
19 #define cusparse_solve    cusparseCcsrsv_solve
20 #define cusparse_analysis cusparseCcsrsv_analysis
21 #define cusparse_csr_spmv cusparseCcsrmv
22 #define cusparse_csr2csc  cusparseCcsr2csc
23 #define cusparse_hyb_spmv cusparseChybmv
24 #define cusparse_csr2hyb  cusparseCcsr2hyb
25 #define cusparse_hyb2csr  cusparseChyb2csr
26 cuFloatComplex ALPHA = {1.0f, 0.0f};
27 cuFloatComplex BETA  = {0.0f, 0.0f};
28 #elif defined(PETSC_USE_REAL_DOUBLE)
29 #define cusparse_solve    cusparseZcsrsv_solve
30 #define cusparse_analysis cusparseZcsrsv_analysis
31 #define cusparse_csr_spmv cusparseZcsrmv
32 #define cusparse_csr2csc  cusparseZcsr2csc
33 #define cusparse_hyb_spmv cusparseZhybmv
34 #define cusparse_csr2hyb  cusparseZcsr2hyb
35 #define cusparse_hyb2csr  cusparseZhyb2csr
36 cuDoubleComplex ALPHA = {1.0, 0.0};
37 cuDoubleComplex BETA  = {0.0, 0.0};
38 #endif
39 #else
40 PetscScalar ALPHA = 1.0;
41 PetscScalar BETA  = 0.0;
42 #if defined(PETSC_USE_REAL_SINGLE)
43 #define cusparse_solve    cusparseScsrsv_solve
44 #define cusparse_analysis cusparseScsrsv_analysis
45 #define cusparse_csr_spmv cusparseScsrmv
46 #define cusparse_csr2csc  cusparseScsr2csc
47 #define cusparse_hyb_spmv cusparseShybmv
48 #define cusparse_csr2hyb  cusparseScsr2hyb
49 #define cusparse_hyb2csr  cusparseShyb2csr
50 #elif defined(PETSC_USE_REAL_DOUBLE)
51 #define cusparse_solve    cusparseDcsrsv_solve
52 #define cusparse_analysis cusparseDcsrsv_analysis
53 #define cusparse_csr_spmv cusparseDcsrmv
54 #define cusparse_csr2csc  cusparseDcsr2csc
55 #define cusparse_hyb_spmv cusparseDhybmv
56 #define cusparse_csr2hyb  cusparseDcsr2hyb
57 #define cusparse_hyb2csr  cusparseDhyb2csr
58 #endif
59 #endif
60 
61 #define THRUSTINTARRAY32 thrust::device_vector<int>
62 #define THRUSTINTARRAY thrust::device_vector<PetscInt>
63 #define THRUSTARRAY thrust::device_vector<PetscScalar>
64 
65 /* A CSR matrix structure */
66 struct CsrMatrix {
67   PetscInt         num_rows;
68   PetscInt         num_cols;
69   PetscInt         num_entries;
70   THRUSTINTARRAY32 *row_offsets;
71   THRUSTINTARRAY32 *column_indices;
72   THRUSTARRAY      *values;
73 };
74 
75 //#define CUSPMATRIXCSR32 cusp::csr_matrix<int,PetscScalar,cusp::device_memory>
76 
77 /* This is struct holding the relevant data needed to a MatSolve */
78 struct Mat_SeqAIJCUSPARSETriFactorStruct {
79   /* Data needed for triangular solve */
80   cusparseMatDescr_t          descr;
81   cusparseSolveAnalysisInfo_t solveInfo;
82   cusparseOperation_t         solveOp;
83   CsrMatrix                   *csrMat;
84 };
85 
86 /* This is struct holding the relevant data needed to a MatMult */
87 struct Mat_SeqAIJCUSPARSEMultStruct {
88   void               *mat;   /* opaque pointer to a matrix. This could be either a cusparseHybMat_t or a CsrMatrix */
89   cusparseMatDescr_t descr;   /* Data needed to describe the matrix for a multiply */
90   THRUSTINTARRAY     *cprowIndices;   /* compressed row indices used in the parallel SpMV */
91 };
92 
93 /* This is a larger struct holding all the triangular factors for a solve, transpose solve, and
94  any indices used in a reordering */
95 struct Mat_SeqAIJCUSPARSETriFactors {
96   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorPtr; /* pointer for lower triangular (factored matrix) on GPU */
97   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorPtr; /* pointer for upper triangular (factored matrix) on GPU */
98   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorPtrTranspose; /* pointer for lower triangular (factored matrix) on GPU for the transpose (useful for BiCG) */
99   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorPtrTranspose; /* pointer for upper triangular (factored matrix) on GPU for the transpose (useful for BiCG)*/
100   THRUSTINTARRAY                    *rpermIndices;  /* indices used for any reordering */
101   THRUSTINTARRAY                    *cpermIndices;  /* indices used for any reordering */
102   THRUSTARRAY                       *workVector;
103   MatCUSPARSEStorageFormat          format;   /* the storage format for the matrix on the device */
104   cusparseHandle_t                  handle;   /* a handle to the cusparse library */
105   PetscInt                          nnz;      /* number of nonzeros ... need this for accurate logging between ICC and ILU */
106 };
107 
108 /* This is a larger struct holding all the matrices for a SpMV, and SpMV Tranpose */
109 struct Mat_SeqAIJCUSPARSE {
110   Mat_SeqAIJCUSPARSEMultStruct *mat; /* pointer to the matrix on the GPU */
111   Mat_SeqAIJCUSPARSEMultStruct *matTranspose; /* pointer to the matrix on the GPU (for the transpose ... useful for BiCG) */
112   THRUSTARRAY                  *workVector; /*pointer to a workvector to which we can copy the relevant indices of a vector we want to multiply */
113   PetscInt                     nonzerorow; /* number of nonzero rows ... used in the flop calculations */
114   MatCUSPARSEStorageFormat     format;   /* the storage format for the matrix on the device */
115   cudaStream_t                 stream;   /* a stream for the parallel SpMV ... this is not owned and should not be deleted */
116   cusparseHandle_t             handle;   /* a handle to the cusparse library ... this may not be owned (if we're working in parallel i.e. multiGPUs) */
117 };
118 
119 PETSC_INTERN PetscErrorCode MatCUSPARSECopyToGPU(Mat);
120 #endif
121