xref: /petsc/src/mat/impls/aij/seq/seqcusparse/cusparsematimpl.h (revision 261a78b473f90dfbc8ba911c87f1bd46aa32a06a)
1 #if !defined(CUSPARSEMATIMPL)
2 #define CUSPARSEMATIMPL
3 
4 #include <petscpkg_version.h>
5 #include <petsc/private/cudavecimpl.h>
6 #include <petscaijdevice.h>
7 
8 #include <cusparse_v2.h>
9 
10 #include <algorithm>
11 #include <vector>
12 
13 #include <thrust/device_vector.h>
14 #include <thrust/device_ptr.h>
15 #include <thrust/device_malloc_allocator.h>
16 #include <thrust/transform.h>
17 #include <thrust/functional.h>
18 #include <thrust/sequence.h>
19 #include <thrust/system/system_error.h>
20 
21 #define PetscStackCallThrust(body) do {                                     \
22     try {                                                                   \
23       body;                                                                 \
24     } catch(thrust::system_error& e) {                                      \
25       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Thrust %s",e.what());\
26     }                                                                       \
27   } while (0)
28 
29 #if defined(PETSC_USE_COMPLEX)
30   #if defined(PETSC_USE_REAL_SINGLE)
31     const cuComplex PETSC_CUSPARSE_ONE        = {1.0f, 0.0f};
32     const cuComplex PETSC_CUSPARSE_ZERO       = {0.0f, 0.0f};
33   #elif defined(PETSC_USE_REAL_DOUBLE)
34     const cuDoubleComplex PETSC_CUSPARSE_ONE  = {1.0, 0.0};
35     const cuDoubleComplex PETSC_CUSPARSE_ZERO = {0.0, 0.0};
36   #endif
37 #else
38   const PetscScalar PETSC_CUSPARSE_ONE        = 1.0;
39   const PetscScalar PETSC_CUSPARSE_ZERO       = 0.0;
40 #endif
41 
42 #if PETSC_PKG_CUDA_VERSION_GE(9,0,0)
43   #define csrsvInfo_t              csrsv2Info_t
44   #define cusparseCreateCsrsvInfo  cusparseCreateCsrsv2Info
45   #define cusparseDestroyCsrsvInfo cusparseDestroyCsrsv2Info
46   #if defined(PETSC_USE_COMPLEX)
47     #if defined(PETSC_USE_REAL_SINGLE)
48       #define cusparseXcsrsv_buffsize(a,b,c,d,e,f,g,h,i,j)       cusparseCcsrsv2_bufferSize(a,b,c,d,e,(cuComplex*)(f),g,h,i,j)
49       #define cusparseXcsrsv_analysis(a,b,c,d,e,f,g,h,i,j,k)     cusparseCcsrsv2_analysis(a,b,c,d,e,(const cuComplex*)(f),g,h,i,j,k)
50       #define cusparseXcsrsv_solve(a,b,c,d,e,f,g,h,i,j,k,l,m,n)  cusparseCcsrsv2_solve(a,b,c,d,(const cuComplex*)(e),f,(const cuComplex*)(g),h,i,j,(const cuComplex*)(k),(cuComplex*)(l),m,n)
51     #elif defined(PETSC_USE_REAL_DOUBLE)
52       #define cusparseXcsrsv_buffsize(a,b,c,d,e,f,g,h,i,j)       cusparseZcsrsv2_bufferSize(a,b,c,d,e,(cuDoubleComplex*)(f),g,h,i,j)
53       #define cusparseXcsrsv_analysis(a,b,c,d,e,f,g,h,i,j,k)     cusparseZcsrsv2_analysis(a,b,c,d,e,(const cuDoubleComplex*)(f),g,h,i,j,k)
54       #define cusparseXcsrsv_solve(a,b,c,d,e,f,g,h,i,j,k,l,m,n)  cusparseZcsrsv2_solve(a,b,c,d,(const cuDoubleComplex*)(e),f,(const cuDoubleComplex*)(g),h,i,j,(const cuDoubleComplex*)(k),(cuDoubleComplex*)(l),m,n)
55     #endif
56   #else /* not complex */
57     #if defined(PETSC_USE_REAL_SINGLE)
58       #define cusparseXcsrsv_buffsize       cusparseScsrsv2_bufferSize
59       #define cusparseXcsrsv_analysis       cusparseScsrsv2_analysis
60       #define cusparseXcsrsv_solve          cusparseScsrsv2_solve
61     #elif defined(PETSC_USE_REAL_DOUBLE)
62       #define cusparseXcsrsv_buffsize       cusparseDcsrsv2_bufferSize
63       #define cusparseXcsrsv_analysis       cusparseDcsrsv2_analysis
64       #define cusparseXcsrsv_solve          cusparseDcsrsv2_solve
65     #endif
66   #endif
67 #else
68   #define csrsvInfo_t              cusparseSolveAnalysisInfo_t
69   #define cusparseCreateCsrsvInfo  cusparseCreateSolveAnalysisInfo
70   #define cusparseDestroyCsrsvInfo cusparseDestroySolveAnalysisInfo
71   #if defined(PETSC_USE_COMPLEX)
72     #if defined(PETSC_USE_REAL_SINGLE)
73       #define cusparseXcsrsv_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))
74       #define cusparseXcsrsv_analysis(a,b,c,d,e,f,g,h,i)  cusparseCcsrsv_analysis((a),(b),(c),(d),(e),(cuComplex*)(f),(g),(h),(i))
75     #elif defined(PETSC_USE_REAL_DOUBLE)
76       #define cusparseXcsrsv_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))
77       #define cusparseXcsrsv_analysis(a,b,c,d,e,f,g,h,i)  cusparseZcsrsv_analysis((a),(b),(c),(d),(e),(cuDoubleComplex*)(f),(g),(h),(i))
78     #endif
79   #else /* not complex */
80     #if defined(PETSC_USE_REAL_SINGLE)
81       #define cusparseXcsrsv_solve    cusparseScsrsv_solve
82       #define cusparseXcsrsv_analysis cusparseScsrsv_analysis
83     #elif defined(PETSC_USE_REAL_DOUBLE)
84       #define cusparseXcsrsv_solve    cusparseDcsrsv_solve
85       #define cusparseXcsrsv_analysis cusparseDcsrsv_analysis
86     #endif
87   #endif
88 #endif
89 
90 #if PETSC_PKG_CUDA_VERSION_GE(11,0,0)
91   #define cusparse_csr2csc cusparseCsr2cscEx2
92   #if defined(PETSC_USE_COMPLEX)
93     #if defined(PETSC_USE_REAL_SINGLE)
94       #define cusparse_scalartype CUDA_C_32F
95       #define cusparse_csr_spgeam(a,b,c,d,e,f,g,h,i,j,k,l,m,n,o,p,q,r,s,t)            cusparseCcsrgeam2(a,b,c,(cuComplex*)d,e,f,(cuComplex*)g,h,i,(cuComplex*)j,k,l,(cuComplex*)m,n,o,p,(cuComplex*)q,r,s,t)
96       #define cusparse_csr_spgeam_bufferSize(a,b,c,d,e,f,g,h,i,j,k,l,m,n,o,p,q,r,s,t) cusparseCcsrgeam2_bufferSizeExt(a,b,c,(cuComplex*)d,e,f,(cuComplex*)g,h,i,(cuComplex*)j,k,l,(cuComplex*)m,n,o,p,(cuComplex*)q,r,s,t)
97     #elif defined(PETSC_USE_REAL_DOUBLE)
98       #define cusparse_scalartype CUDA_C_64F
99       #define cusparse_csr_spgeam(a,b,c,d,e,f,g,h,i,j,k,l,m,n,o,p,q,r,s,t)            cusparseZcsrgeam2(a,b,c,(cuDoubleComplex*)d,e,f,(cuDoubleComplex*)g,h,i,(cuDoubleComplex*)j,k,l,(cuDoubleComplex*)m,n,o,p,(cuDoubleComplex*)q,r,s,t)
100       #define cusparse_csr_spgeam_bufferSize(a,b,c,d,e,f,g,h,i,j,k,l,m,n,o,p,q,r,s,t) cusparseZcsrgeam2_bufferSizeExt(a,b,c,(cuDoubleComplex*)d,e,f,(cuDoubleComplex*)g,h,i,(cuDoubleComplex*)j,k,l,(cuDoubleComplex*)m,n,o,p,(cuDoubleComplex*)q,r,s,t)
101     #endif
102   #else /* not complex */
103     #if defined(PETSC_USE_REAL_SINGLE)
104       #define cusparse_scalartype CUDA_R_32F
105       #define cusparse_csr_spgeam            cusparseScsrgeam2
106       #define cusparse_csr_spgeam_bufferSize cusparseScsrgeam2_bufferSizeExt
107     #elif defined(PETSC_USE_REAL_DOUBLE)
108       #define cusparse_scalartype CUDA_R_64F
109       #define cusparse_csr_spgeam            cusparseDcsrgeam2
110       #define cusparse_csr_spgeam_bufferSize cusparseDcsrgeam2_bufferSizeExt
111     #endif
112   #endif
113 #else
114   #if defined(PETSC_USE_COMPLEX)
115     #if defined(PETSC_USE_REAL_SINGLE)
116       #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))
117       #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))
118       #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))
119       #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))
120       #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))
121       #define cusparse_hyb2csr(a,b,c,d,e,f)                      cusparseChyb2csr((a),(b),(c),(cuComplex*)(d),(e),(f))
122       #define cusparse_csr_spgemm(a,b,c,d,e,f,g,h,i,j,k,l,m,n,o,p,q,r,s,t) cusparseCcsrgemm(a,b,c,d,e,f,g,h,(cuComplex*)i,j,k,l,m,(cuComplex*)n,o,p,q,(cuComplex*)r,s,t)
123       #define cusparse_csr_spgeam(a,b,c,d,e,f,g,h,i,j,k,l,m,n,o,p,q,r,s)   cusparseCcsrgeam(a,b,c,(cuComplex*)d,e,f,(cuComplex*)g,h,i,(cuComplex*)j,k,l,(cuComplex*)m,n,o,p,(cuComplex*)q,r,s)
124     #elif defined(PETSC_USE_REAL_DOUBLE)
125       #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))
126       #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))
127       #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))
128       #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))
129       #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))
130       #define cusparse_hyb2csr(a,b,c,d,e,f)                      cusparseZhyb2csr((a),(b),(c),(cuDoubleComplex*)(d),(e),(f))
131       #define cusparse_csr_spgemm(a,b,c,d,e,f,g,h,i,j,k,l,m,n,o,p,q,r,s,t) cusparseZcsrgemm(a,b,c,d,e,f,g,h,(cuDoubleComplex*)i,j,k,l,m,(cuDoubleComplex*)n,o,p,q,(cuDoubleComplex*)r,s,t)
132       #define cusparse_csr_spgeam(a,b,c,d,e,f,g,h,i,j,k,l,m,n,o,p,q,r,s)   cusparseZcsrgeam(a,b,c,(cuDoubleComplex*)d,e,f,(cuDoubleComplex*)g,h,i,(cuDoubleComplex*)j,k,l,(cuDoubleComplex*)m,n,o,p,(cuDoubleComplex*)q,r,s)
133     #endif
134   #else
135     #if defined(PETSC_USE_REAL_SINGLE)
136       #define cusparse_csr_spmv cusparseScsrmv
137       #define cusparse_csr_spmm cusparseScsrmm
138       #define cusparse_csr2csc  cusparseScsr2csc
139       #define cusparse_hyb_spmv cusparseShybmv
140       #define cusparse_csr2hyb  cusparseScsr2hyb
141       #define cusparse_hyb2csr  cusparseShyb2csr
142       #define cusparse_csr_spgemm cusparseScsrgemm
143       #define cusparse_csr_spgeam cusparseScsrgeam
144     #elif defined(PETSC_USE_REAL_DOUBLE)
145       #define cusparse_csr_spmv cusparseDcsrmv
146       #define cusparse_csr_spmm cusparseDcsrmm
147       #define cusparse_csr2csc  cusparseDcsr2csc
148       #define cusparse_hyb_spmv cusparseDhybmv
149       #define cusparse_csr2hyb  cusparseDcsr2hyb
150       #define cusparse_hyb2csr  cusparseDhyb2csr
151       #define cusparse_csr_spgemm cusparseDcsrgemm
152       #define cusparse_csr_spgeam cusparseDcsrgeam
153     #endif
154   #endif
155 #endif
156 
157 #define THRUSTINTARRAY32 thrust::device_vector<int>
158 #define THRUSTINTARRAY thrust::device_vector<PetscInt>
159 #define THRUSTARRAY thrust::device_vector<PetscScalar>
160 
161 /* A CSR matrix structure */
162 struct CsrMatrix {
163   PetscInt         num_rows;
164   PetscInt         num_cols;
165   PetscInt         num_entries;
166   THRUSTINTARRAY32 *row_offsets;
167   THRUSTINTARRAY32 *column_indices;
168   THRUSTARRAY      *values;
169 };
170 
171 /* This is struct holding the relevant data needed to a MatSolve */
172 struct Mat_SeqAIJCUSPARSETriFactorStruct {
173   /* Data needed for triangular solve */
174   cusparseMatDescr_t          descr;
175   cusparseOperation_t         solveOp;
176   CsrMatrix                   *csrMat;
177   csrsvInfo_t                 solveInfo;
178   cusparseSolvePolicy_t       solvePolicy;     /* whether level information is generated and used */
179   int                         solveBufferSize;
180   void                        *solveBuffer;
181   size_t                      csr2cscBufferSize; /* to transpose the triangular factor (only used for CUDA >= 11.0) */
182   void                        *csr2cscBuffer;
183   PetscScalar                 *AA_h; /* managed host buffer for moving values to the GPU */
184 };
185 
186 /* This is a larger struct holding all the triangular factors for a solve, transpose solve, and any indices used in a reordering */
187 struct Mat_SeqAIJCUSPARSETriFactors {
188   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorPtr; /* pointer for lower triangular (factored matrix) on GPU */
189   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorPtr; /* pointer for upper triangular (factored matrix) on GPU */
190   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorPtrTranspose; /* pointer for lower triangular (factored matrix) on GPU for the transpose (useful for BiCG) */
191   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorPtrTranspose; /* pointer for upper triangular (factored matrix) on GPU for the transpose (useful for BiCG)*/
192   THRUSTINTARRAY                    *rpermIndices;  /* indices used for any reordering */
193   THRUSTINTARRAY                    *cpermIndices;  /* indices used for any reordering */
194   THRUSTARRAY                       *workVector;
195   cusparseHandle_t                  handle;   /* a handle to the cusparse library */
196   PetscInt                          nnz;      /* number of nonzeros ... need this for accurate logging between ICC and ILU */
197   PetscScalar                       *a_band_d; /* GPU data for banded CSR LU factorization matrix diag(L)=1 */
198   int                               *i_band_d; /* this could be optimized away */
199   cudaDeviceProp                    dev_prop;
200   PetscBool                         init_dev_prop;
201 };
202 
203 struct Mat_CusparseSpMV {
204   PetscBool             initialized;    /* Don't rely on spmvBuffer != NULL to test if the struct is initialized, */
205   size_t                spmvBufferSize; /* since I'm not sure if smvBuffer can be NULL even after cusparseSpMV_bufferSize() */
206   void                  *spmvBuffer;
207  #if PETSC_PKG_CUDA_VERSION_GE(11,0,0)  /* these are present from CUDA 10.1, but PETSc code makes use of them from CUDA 11 on */
208   cusparseDnVecDescr_t  vecXDescr,vecYDescr; /* descriptor for the dense vectors in y=op(A)x */
209  #endif
210 };
211 
212 /* This is struct holding the relevant data needed to a MatMult */
213 struct Mat_SeqAIJCUSPARSEMultStruct {
214   void               *mat;  /* opaque pointer to a matrix. This could be either a cusparseHybMat_t or a CsrMatrix */
215   cusparseMatDescr_t descr; /* Data needed to describe the matrix for a multiply */
216   THRUSTINTARRAY     *cprowIndices;   /* compressed row indices used in the parallel SpMV */
217   PetscScalar        *alpha_one; /* pointer to a device "scalar" storing the alpha parameter in the SpMV */
218   PetscScalar        *beta_zero; /* pointer to a device "scalar" storing the beta parameter in the SpMV as zero*/
219   PetscScalar        *beta_one; /* pointer to a device "scalar" storing the beta parameter in the SpMV as one */
220  #if PETSC_PKG_CUDA_VERSION_GE(11,0,0)
221   cusparseSpMatDescr_t  matDescr;  /* descriptor for the matrix, used by SpMV and SpMM */
222   Mat_CusparseSpMV      cuSpMV[3]; /* different Mat_CusparseSpMV structs for non-transpose, transpose, conj-transpose */
223   Mat_SeqAIJCUSPARSEMultStruct() : matDescr(NULL) {
224     for (int i=0; i<3; i++) cuSpMV[i].initialized = PETSC_FALSE;
225   }
226  #endif
227 };
228 
229 /* This is a larger struct holding all the matrices for a SpMV, and SpMV Transpose */
230 struct Mat_SeqAIJCUSPARSE {
231   Mat_SeqAIJCUSPARSEMultStruct *mat;            /* pointer to the matrix on the GPU */
232   Mat_SeqAIJCUSPARSEMultStruct *matTranspose;   /* pointer to the matrix on the GPU (for the transpose ... useful for BiCG) */
233   THRUSTARRAY                  *workVector;     /* pointer to a workvector to which we can copy the relevant indices of a vector we want to multiply */
234   THRUSTINTARRAY32             *rowoffsets_gpu; /* rowoffsets on GPU in non-compressed-row format. It is used to convert CSR to CSC */
235   PetscInt                     nrows;           /* number of rows of the matrix seen by GPU */
236   MatCUSPARSEStorageFormat     format;          /* the storage format for the matrix on the device */
237   PetscBool                    use_cpu_solve;   /* Use AIJ_Seq (I)LU solve */
238   cudaStream_t                 stream;          /* a stream for the parallel SpMV ... this is not owned and should not be deleted */
239   cusparseHandle_t             handle;          /* a handle to the cusparse library ... this may not be owned (if we're working in parallel i.e. multiGPUs) */
240   PetscObjectState             nonzerostate;    /* track nonzero state to possibly recreate the GPU matrix */
241  #if PETSC_PKG_CUDA_VERSION_GE(11,0,0)
242   size_t                       csr2cscBufferSize; /* stuff used to compute the matTranspose above */
243   void                         *csr2cscBuffer;    /* This is used as a C struct and is calloc'ed by PetscNewLog() */
244   cusparseCsr2CscAlg_t         csr2cscAlg;        /* algorithms can be selected from command line options */
245   cusparseSpMVAlg_t            spmvAlg;
246   cusparseSpMMAlg_t            spmmAlg;
247  #endif
248   THRUSTINTARRAY               *csr2csc_i;
249   PetscSplitCSRDataStructure   deviceMat;       /* Matrix on device for, eg, assembly */
250 
251   /* Stuff for basic COO support */
252   THRUSTINTARRAY               *cooPerm;        /* permutation array that sorts the input coo entris by row and col */
253   THRUSTINTARRAY               *cooPerm_a;      /* ordered array that indicate i-th nonzero (after sorting) is the j-th unique nonzero */
254 
255   /* Stuff for extended COO support */
256   PetscBool                    use_extended_coo; /* Use extended COO format */
257   PetscCount                   *jmap_d; /* perm[disp+jmap[i]..disp+jmap[i+1]) gives indices of entries in v[] associated with i-th nonzero of the matrix */
258   PetscCount                   *perm_d;
259 
260   Mat_SeqAIJCUSPARSE() : use_extended_coo(PETSC_FALSE), perm_d(NULL), jmap_d(NULL) {}
261 };
262 
263 PETSC_INTERN PetscErrorCode MatSeqAIJCUSPARSECopyToGPU(Mat);
264 PETSC_INTERN PetscErrorCode MatSetPreallocationCOO_SeqAIJCUSPARSE_Basic(Mat,PetscCount,const PetscInt[],const PetscInt[]);
265 PETSC_INTERN PetscErrorCode MatSetValuesCOO_SeqAIJCUSPARSE_Basic(Mat,const PetscScalar[],InsertMode);
266 PETSC_INTERN PetscErrorCode MatSeqAIJCUSPARSEMergeMats(Mat,Mat,MatReuse,Mat*);
267 PETSC_INTERN PetscErrorCode MatSeqAIJCUSPARSETriFactors_Reset(Mat_SeqAIJCUSPARSETriFactors_p*);
268 
269 static inline bool isCudaMem(const void *data)
270 {
271   cudaError_t                  cerr;
272   struct cudaPointerAttributes attr;
273   enum cudaMemoryType          mtype;
274   cerr = cudaPointerGetAttributes(&attr,data); /* Do not check error since before CUDA 11.0, passing a host pointer returns cudaErrorInvalidValue */
275   cudaGetLastError(); /* Reset the last error */
276   #if (CUDART_VERSION < 10000)
277     mtype = attr.memoryType;
278   #else
279     mtype = attr.type;
280   #endif
281   if (cerr == cudaSuccess && mtype == cudaMemoryTypeDevice) return true;
282   else return false;
283 }
284 
285 #endif
286