xref: /petsc/src/mat/impls/aij/seq/seqhipsparse/hipsparsematimpl.h (revision a207d08edab748df642306d213e6e891ba48ee92)
1 /* Portions of this code are under:
2    Copyright (c) 2022 Advanced Micro Devices, Inc. All rights reserved.
3 */
4 #ifndef PETSC_HIPSPARSEMATIMPL_H
5 #define PETSC_HIPSPARSEMATIMPL_H
6 
7 #include <petscpkg_version.h>
8 #include <../src/vec/vec/impls/seq/cupm/vecseqcupm.hpp> /* for VecSeq_CUPM */
9 #include <petsc/private/petsclegacycupmblas.h>
10 #include <petscaijdevice.h>
11 
12 #if PETSC_PKG_HIP_VERSION_GE(5, 2, 0)
13   #include <hipsparse/hipsparse.h>
14 #else /* PETSC_PKG_HIP_VERSION_GE(5,2,0) */
15   #include <hipsparse.h>
16 #endif /* PETSC_PKG_HIP_VERSION_GE(5,2,0) */
17 #include "hip/hip_runtime.h"
18 
19 #include <algorithm>
20 #include <vector>
21 
22 #include <thrust/device_vector.h>
23 #include <thrust/device_ptr.h>
24 #include <thrust/device_malloc_allocator.h>
25 #include <thrust/transform.h>
26 #include <thrust/functional.h>
27 #include <thrust/sequence.h>
28 #include <thrust/system/system_error.h>
29 
30 #if defined(PETSC_USE_COMPLEX)
31   #if defined(PETSC_USE_REAL_SINGLE)
32 const hipComplex PETSC_HIPSPARSE_ONE  = {1.0f, 0.0f};
33 const hipComplex PETSC_HIPSPARSE_ZERO = {0.0f, 0.0f};
34     #define hipsparseXcsrilu02_bufferSize(a, b, c, d, e, f, g, h, i)  hipsparseCcsrilu02_bufferSize(a, b, c, d, (hipComplex *)e, f, g, h, i)
35     #define hipsparseXcsrilu02_analysis(a, b, c, d, e, f, g, h, i, j) hipsparseCcsrilu02_analysis(a, b, c, d, (hipComplex *)e, f, g, h, i, j)
36     #define hipsparseXcsrilu02(a, b, c, d, e, f, g, h, i, j)          hipsparseCcsrilu02(a, b, c, d, (hipComplex *)e, f, g, h, i, j)
37     #define hipsparseXcsric02_bufferSize(a, b, c, d, e, f, g, h, i)   hipsparseCcsric02_bufferSize(a, b, c, d, (hipComplex *)e, f, g, h, i)
38     #define hipsparseXcsric02_analysis(a, b, c, d, e, f, g, h, i, j)  hipsparseCcsric02_analysis(a, b, c, d, (hipComplex *)e, f, g, h, i, j)
39     #define hipsparseXcsric02(a, b, c, d, e, f, g, h, i, j)           hipsparseCcsric02(a, b, c, d, (hipComplex *)e, f, g, h, i, j)
40   #elif defined(PETSC_USE_REAL_DOUBLE)
41 const hipDoubleComplex PETSC_HIPSPARSE_ONE  = {1.0, 0.0};
42 const hipDoubleComplex PETSC_HIPSPARSE_ZERO = {0.0, 0.0};
43     #define hipsparseXcsrilu02_bufferSize(a, b, c, d, e, f, g, h, i)  hipsparseZcsrilu02_bufferSize(a, b, c, d, (hipDoubleComplex *)e, f, g, h, i)
44     #define hipsparseXcsrilu02_analysis(a, b, c, d, e, f, g, h, i, j) hipsparseZcsrilu02_analysis(a, b, c, d, (hipDoubleComplex *)e, f, g, h, i, j)
45     #define hipsparseXcsrilu02(a, b, c, d, e, f, g, h, i, j)          hipsparseZcsrilu02(a, b, c, d, (hipDoubleComplex *)e, f, g, h, i, j)
46     #define hipsparseXcsric02_bufferSize(a, b, c, d, e, f, g, h, i)   hipsparseZcsric02_bufferSize(a, b, c, d, (hipDoubleComplex *)e, f, g, h, i)
47     #define hipsparseXcsric02_analysis(a, b, c, d, e, f, g, h, i, j)  hipsparseZcsric02_analysis(a, b, c, d, (hipDoubleComplex *)e, f, g, h, i, j)
48     #define hipsparseXcsric02(a, b, c, d, e, f, g, h, i, j)           hipsparseZcsric02(a, b, c, d, (hipDoubleComplex *)e, f, g, h, i, j)
49   #endif /* Single or double */
50 #else    /* not complex */
51 const PetscScalar PETSC_HIPSPARSE_ONE  = 1.0;
52 const PetscScalar PETSC_HIPSPARSE_ZERO = 0.0;
53   #if defined(PETSC_USE_REAL_SINGLE)
54     #define hipsparseXcsrilu02_bufferSize hipsparseScsrilu02_bufferSize
55     #define hipsparseXcsrilu02_analysis   hipsparseScsrilu02_analysis
56     #define hipsparseXcsrilu02            hipsparseScsrilu02
57     #define hipsparseXcsric02_bufferSize  hipsparseScsric02_bufferSize
58     #define hipsparseXcsric02_analysis    hipsparseScsric02_analysis
59     #define hipsparseXcsric02             hipsparseScsric02
60   #elif defined(PETSC_USE_REAL_DOUBLE)
61     #define hipsparseXcsrilu02_bufferSize hipsparseDcsrilu02_bufferSize
62     #define hipsparseXcsrilu02_analysis   hipsparseDcsrilu02_analysis
63     #define hipsparseXcsrilu02            hipsparseDcsrilu02
64     #define hipsparseXcsric02_bufferSize  hipsparseDcsric02_bufferSize
65     #define hipsparseXcsric02_analysis    hipsparseDcsric02_analysis
66     #define hipsparseXcsric02             hipsparseDcsric02
67   #endif /* Single or double */
68 #endif   /* complex or not */
69 
70 #define csrsvInfo_t               csrsv2Info_t
71 #define hipsparseCreateCsrsvInfo  hipsparseCreateCsrsv2Info
72 #define hipsparseDestroyCsrsvInfo hipsparseDestroyCsrsv2Info
73 #if defined(PETSC_USE_COMPLEX)
74   #if defined(PETSC_USE_REAL_SINGLE)
75     #define hipsparseXcsrsv_buffsize(a, b, c, d, e, f, g, h, i, j)          hipsparseCcsrsv2_bufferSize(a, b, c, d, e, (hipComplex *)(f), g, h, i, j)
76     #define hipsparseXcsrsv_analysis(a, b, c, d, e, f, g, h, i, j, k)       hipsparseCcsrsv2_analysis(a, b, c, d, e, (const hipComplex *)(f), g, h, i, j, k)
77     #define hipsparseXcsrsv_solve(a, b, c, d, e, f, g, h, i, j, k, l, m, n) hipsparseCcsrsv2_solve(a, b, c, d, (const hipComplex *)(e), f, (const hipComplex *)(g), h, i, j, (const hipComplex *)(k), (hipComplex *)(l), m, n)
78   #elif defined(PETSC_USE_REAL_DOUBLE)
79     #define hipsparseXcsrsv_buffsize(a, b, c, d, e, f, g, h, i, j)          hipsparseZcsrsv2_bufferSize(a, b, c, d, e, (hipDoubleComplex *)(f), g, h, i, j)
80     #define hipsparseXcsrsv_analysis(a, b, c, d, e, f, g, h, i, j, k)       hipsparseZcsrsv2_analysis(a, b, c, d, e, (const hipDoubleComplex *)(f), g, h, i, j, k)
81     #define hipsparseXcsrsv_solve(a, b, c, d, e, f, g, h, i, j, k, l, m, n) hipsparseZcsrsv2_solve(a, b, c, d, (const hipDoubleComplex *)(e), f, (const hipDoubleComplex *)(g), h, i, j, (const hipDoubleComplex *)(k), (hipDoubleComplex *)(l), m, n)
82   #endif /* Single or double */
83 #else    /* not complex */
84   #if defined(PETSC_USE_REAL_SINGLE)
85     #define hipsparseXcsrsv_buffsize hipsparseScsrsv2_bufferSize
86     #define hipsparseXcsrsv_analysis hipsparseScsrsv2_analysis
87     #define hipsparseXcsrsv_solve    hipsparseScsrsv2_solve
88   #elif defined(PETSC_USE_REAL_DOUBLE)
89     #define hipsparseXcsrsv_buffsize hipsparseDcsrsv2_bufferSize
90     #define hipsparseXcsrsv_analysis hipsparseDcsrsv2_analysis
91     #define hipsparseXcsrsv_solve    hipsparseDcsrsv2_solve
92   #endif /* Single or double */
93 #endif   /* not complex */
94 
95 #if PETSC_PKG_HIP_VERSION_GE(4, 5, 0)
96   // #define cusparse_csr2csc cusparseCsr2cscEx2
97   #if defined(PETSC_USE_COMPLEX)
98     #if defined(PETSC_USE_REAL_SINGLE)
99       #define hipsparse_scalartype                                                             HIP_C_32F
100       #define hipsparse_csr_spgeam(a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q, r, s, t) hipsparseCcsrgeam2(a, b, c, (hipComplex *)d, e, f, (hipComplex *)g, h, i, (hipComplex *)j, k, l, (hipComplex *)m, n, o, p, (hipComplex *)q, r, s, t)
101       #define hipsparse_csr_spgeam_bufferSize(a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q, r, s, t) \
102         hipsparseCcsrgeam2_bufferSizeExt(a, b, c, (hipComplex *)d, e, f, (hipComplex *)g, h, i, (hipComplex *)j, k, l, (hipComplex *)m, n, o, p, (hipComplex *)q, r, s, t)
103     #elif defined(PETSC_USE_REAL_DOUBLE)
104       #define hipsparse_scalartype HIP_C_64F
105       #define hipsparse_csr_spgeam(a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q, r, s, t) \
106         hipsparseZcsrgeam2(a, b, c, (hipDoubleComplex *)d, e, f, (hipDoubleComplex *)g, h, i, (hipDoubleComplex *)j, k, l, (hipDoubleComplex *)m, n, o, p, (hipDoubleComplex *)q, r, s, t)
107       #define hipsparse_csr_spgeam_bufferSize(a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q, r, s, t) \
108         hipsparseZcsrgeam2_bufferSizeExt(a, b, c, (hipDoubleComplex *)d, e, f, (hipDoubleComplex *)g, h, i, (hipDoubleComplex *)j, k, l, (hipDoubleComplex *)m, n, o, p, (hipDoubleComplex *)q, r, s, t)
109     #endif /* Single or double */
110   #else    /* not complex */
111     #if defined(PETSC_USE_REAL_SINGLE)
112       #define hipsparse_scalartype            HIP_R_32F
113       #define hipsparse_csr_spgeam            hipsparseScsrgeam2
114       #define hipsparse_csr_spgeam_bufferSize hipsparseScsrgeam2_bufferSizeExt
115     #elif defined(PETSC_USE_REAL_DOUBLE)
116       #define hipsparse_scalartype            HIP_R_64F
117       #define hipsparse_csr_spgeam            hipsparseDcsrgeam2
118       #define hipsparse_csr_spgeam_bufferSize hipsparseDcsrgeam2_bufferSizeExt
119     #endif /* Single or double */
120   #endif   /* complex or not */
121 #endif     /* PETSC_PKG_HIP_VERSION_GE(4, 5, 0) */
122 
123 #if defined(PETSC_USE_COMPLEX)
124   #if defined(PETSC_USE_REAL_SINGLE)
125     #define hipsparse_scalartype                                                             HIP_C_32F
126     #define hipsparse_csr_spmv(a, b, c, d, e, f, g, h, i, j, k, l, m)                        hipsparseCcsrmv((a), (b), (c), (d), (e), (hipComplex *)(f), (g), (hipComplex *)(h), (i), (j), (hipComplex *)(k), (hipComplex *)(l), (hipComplex *)(m))
127     #define hipsparse_csr_spmm(a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p)               hipsparseCcsrmm((a), (b), (c), (d), (e), (f), (hipComplex *)(g), (h), (hipComplex *)(i), (j), (k), (hipComplex *)(l), (m), (hipComplex *)(n), (hipComplex *)(o), (p))
128     #define hipsparse_csr2csc(a, b, c, d, e, f, g, h, i, j, k, l)                            hipsparseCcsr2csc((a), (b), (c), (d), (hipComplex *)(e), (f), (g), (hipComplex *)(h), (i), (j), (k), (l))
129     #define hipsparse_hyb_spmv(a, b, c, d, e, f, g, h)                                       hipsparseChybmv((a), (b), (hipComplex *)(c), (d), (e), (hipComplex *)(f), (hipComplex *)(g), (hipComplex *)(h))
130     #define hipsparse_csr2hyb(a, b, c, d, e, f, g, h, i, j)                                  hipsparseCcsr2hyb((a), (b), (c), (d), (hipComplex *)(e), (f), (g), (h), (i), (j))
131     #define hipsparse_hyb2csr(a, b, c, d, e, f)                                              hipsparseChyb2csr((a), (b), (c), (hipComplex *)(d), (e), (f))
132     #define hipsparse_csr_spgemm(a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q, r, s, t) hipsparseCcsrgemm(a, b, c, d, e, f, g, h, (hipComplex *)i, j, k, l, m, (hipComplex *)n, o, p, q, (hipComplex *)r, s, t)
133   // #define hipsparse_csr_spgeam(a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q, r, s)    hipsparseCcsrgeam(a, b, c, (hipComplex *)d, e, f, (hipComplex *)g, h, i, (hipComplex *)j, k, l, (hipComplex *)m, n, o, p, (hipComplex *)q, r, s)
134   #elif defined(PETSC_USE_REAL_DOUBLE)
135     #define hipsparse_scalartype                                      HIP_C_64F
136     #define hipsparse_csr_spmv(a, b, c, d, e, f, g, h, i, j, k, l, m) hipsparseZcsrmv((a), (b), (c), (d), (e), (hipDoubleComplex *)(f), (g), (hipDoubleComplex *)(h), (i), (j), (hipDoubleComplex *)(k), (hipDoubleComplex *)(l), (hipDoubleComplex *)(m))
137     #define hipsparse_csr_spmm(a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p) \
138       hipsparseZcsrmm((a), (b), (c), (d), (e), (f), (hipDoubleComplex *)(g), (h), (hipDoubleComplex *)(i), (j), (k), (hipDoubleComplex *)(l), (m), (hipDoubleComplex *)(n), (hipDoubleComplex *)(o), (p))
139     #define hipsparse_csr2csc(a, b, c, d, e, f, g, h, i, j, k, l)                            hipsparseZcsr2csc((a), (b), (c), (d), (hipDoubleComplex *)(e), (f), (g), (hipDoubleComplex *)(h), (i), (j), (k), (l))
140     #define hipsparse_hyb_spmv(a, b, c, d, e, f, g, h)                                       hipsparseZhybmv((a), (b), (hipDoubleComplex *)(c), (d), (e), (hipDoubleComplex *)(f), (hipDoubleComplex *)(g), (hipDoubleComplex *)(h))
141     #define hipsparse_csr2hyb(a, b, c, d, e, f, g, h, i, j)                                  hipsparseZcsr2hyb((a), (b), (c), (d), (hipDoubleComplex *)(e), (f), (g), (h), (i), (j))
142     #define hipsparse_hyb2csr(a, b, c, d, e, f)                                              hipsparseZhyb2csr((a), (b), (c), (hipDoubleComplex *)(d), (e), (f))
143     #define hipsparse_csr_spgemm(a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q, r, s, t) hipsparseZcsrgemm(a, b, c, d, e, f, g, h, (hipDoubleComplex *)i, j, k, l, m, (hipDoubleComplex *)n, o, p, q, (hipDoubleComplex *)r, s, t)
144   // #define hipsparse_csr_spgeam(a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q, r, s)    hipsparseZcsrgeam(a, b, c, (hipDoubleComplex *)d, e, f, (hipDoubleComplex *)g, h, i, (hipDoubleComplex *)j, k, l, (hipDoubleComplex *)m, n, o, p, (hipDoubleComplex *)q, r, s)
145   #endif /* Single or double */
146 #else    /* not complex */
147   #if defined(PETSC_USE_REAL_SINGLE)
148     #define hipsparse_scalartype HIP_R_32F
149     #define hipsparse_csr_spmv   hipsparseScsrmv
150     #define hipsparse_csr_spmm   hipsparseScsrmm
151     #define hipsparse_csr2csc    hipsparseScsr2csc
152     #define hipsparse_hyb_spmv   hipsparseShybmv
153     #define hipsparse_csr2hyb    hipsparseScsr2hyb
154     #define hipsparse_hyb2csr    hipsparseShyb2csr
155     #define hipsparse_csr_spgemm hipsparseScsrgemm
156   // #define hipsparse_csr_spgeam hipsparseScsrgeam
157   #elif defined(PETSC_USE_REAL_DOUBLE)
158     #define hipsparse_scalartype HIP_R_64F
159     #define hipsparse_csr_spmv   hipsparseDcsrmv
160     #define hipsparse_csr_spmm   hipsparseDcsrmm
161     #define hipsparse_csr2csc    hipsparseDcsr2csc
162     #define hipsparse_hyb_spmv   hipsparseDhybmv
163     #define hipsparse_csr2hyb    hipsparseDcsr2hyb
164     #define hipsparse_hyb2csr    hipsparseDhyb2csr
165     #define hipsparse_csr_spgemm hipsparseDcsrgemm
166   // #define hipsparse_csr_spgeam hipsparseDcsrgeam
167   #endif /* Single or double */
168 #endif   /* complex or not */
169 
170 #define THRUSTINTARRAY32 thrust::device_vector<int>
171 #define THRUSTINTARRAY   thrust::device_vector<PetscInt>
172 #define THRUSTARRAY      thrust::device_vector<PetscScalar>
173 
174 /* A CSR matrix structure */
175 struct CsrMatrix {
176   PetscInt          num_rows;
177   PetscInt          num_cols;
178   PetscInt          num_entries;
179   THRUSTINTARRAY32 *row_offsets;
180   THRUSTINTARRAY32 *column_indices;
181   THRUSTARRAY      *values;
182 };
183 
184 /* This is struct holding the relevant data needed to a MatSolve */
185 struct Mat_SeqAIJHIPSPARSETriFactorStruct {
186   /* Data needed for triangular solve */
187   hipsparseMatDescr_t    descr;
188   hipsparseOperation_t   solveOp;
189   CsrMatrix             *csrMat;
190   csrsvInfo_t            solveInfo;
191   hipsparseSolvePolicy_t solvePolicy; /* whether level information is generated and used */
192   int                    solveBufferSize;
193   void                  *solveBuffer;
194   size_t                 csr2cscBufferSize; /* to transpose the triangular factor (only used for CUDA >= 11.0) */
195   void                  *csr2cscBuffer;
196   PetscScalar           *AA_h; /* managed host buffer for moving values to the GPU */
197 };
198 
199 /* This is a larger struct holding all the triangular factors for a solve, transpose solve, and any indices used in a reordering */
200 struct Mat_SeqAIJHIPSPARSETriFactors {
201   Mat_SeqAIJHIPSPARSETriFactorStruct *loTriFactorPtr;          /* pointer for lower triangular (factored matrix) on GPU */
202   Mat_SeqAIJHIPSPARSETriFactorStruct *upTriFactorPtr;          /* pointer for upper triangular (factored matrix) on GPU */
203   Mat_SeqAIJHIPSPARSETriFactorStruct *loTriFactorPtrTranspose; /* pointer for lower triangular (factored matrix) on GPU for the transpose (useful for BiCG) */
204   Mat_SeqAIJHIPSPARSETriFactorStruct *upTriFactorPtrTranspose; /* pointer for upper triangular (factored matrix) on GPU for the transpose (useful for BiCG)*/
205   THRUSTINTARRAY                     *rpermIndices;            /* indices used for any reordering */
206   THRUSTINTARRAY                     *cpermIndices;            /* indices used for any reordering */
207   THRUSTARRAY                        *workVector;
208   hipsparseHandle_t                   handle;   /* a handle to the hipsparse library */
209   PetscInt                            nnz;      /* number of nonzeros ... need this for accurate logging between ICC and ILU */
210   PetscScalar                        *a_band_d; /* GPU data for banded CSR LU factorization matrix diag(L)=1 */
211   int                                *i_band_d; /* this could be optimized away */
212   hipDeviceProp_t                     dev_prop;
213   PetscBool                           init_dev_prop;
214 
215   /* csrilu0/csric0 appeared in earlier versions of AMD ROCm^{TM}, but we use it along with hipsparseSpSV,
216      which first appeared in hipsparse with ROCm-4.5.0.
217   */
218   PetscBool factorizeOnDevice; /* Do factorization on device or not */
219 #if PETSC_PKG_HIP_VERSION_GE(4, 5, 0)
220   PetscScalar *csrVal;
221   int         *csrRowPtr, *csrColIdx; /* a,i,j of M. Using int since some hipsparse APIs only support 32-bit indices */
222 
223   /* Mixed mat descriptor types? yes, different hipsparse APIs use different types */
224   hipsparseMatDescr_t   matDescr_M;
225   hipsparseSpMatDescr_t spMatDescr_L, spMatDescr_U;
226   hipsparseSpSVDescr_t  spsvDescr_L, spsvDescr_Lt, spsvDescr_U, spsvDescr_Ut;
227 
228   hipsparseDnVecDescr_t dnVecDescr_X, dnVecDescr_Y;
229   PetscScalar          *X, *Y; /* data array of dnVec X and Y */
230 
231   /* Mixed size types? yes */
232   int    factBufferSize_M; /* M ~= LU or LLt */
233   size_t spsvBufferSize_L, spsvBufferSize_Lt, spsvBufferSize_U, spsvBufferSize_Ut;
234   /* hipsparse needs various buffers for factorization and solve of L, U, Lt, or Ut.
235      To save memory, we share the factorization buffer with one of spsvBuffer_L/U.
236   */
237   void *factBuffer_M, *spsvBuffer_L, *spsvBuffer_U, *spsvBuffer_Lt, *spsvBuffer_Ut;
238 
239   csrilu02Info_t         ilu0Info_M;
240   csric02Info_t          ic0Info_M;
241   int                    structural_zero, numerical_zero;
242   hipsparseSolvePolicy_t policy_M;
243 
244   /* In MatSolveTranspose() for ILU0, we use the two flags to do on-demand solve */
245   PetscBool createdTransposeSpSVDescr;    /* Have we created SpSV descriptors for Lt, Ut? */
246   PetscBool updatedTransposeSpSVAnalysis; /* Have we updated SpSV analysis with the latest L, U values? */
247 
248   PetscLogDouble numericFactFlops; /* Estimated FLOPs in ILU0/ICC0 numeric factorization */
249 #endif
250 };
251 
252 struct Mat_HipsparseSpMV {
253   PetscBool             initialized;    /* Don't rely on spmvBuffer != NULL to test if the struct is initialized, */
254   size_t                spmvBufferSize; /* since I'm not sure if smvBuffer can be NULL even after hipsparseSpMV_bufferSize() */
255   void                 *spmvBuffer;
256   hipsparseDnVecDescr_t vecXDescr, vecYDescr; /* descriptor for the dense vectors in y=op(A)x */
257 };
258 
259 /* This is struct holding the relevant data needed to a MatMult */
260 struct Mat_SeqAIJHIPSPARSEMultStruct {
261   void                 *mat;          /* opaque pointer to a matrix. This could be either a hipsparseHybMat_t or a CsrMatrix */
262   hipsparseMatDescr_t   descr;        /* Data needed to describe the matrix for a multiply */
263   THRUSTINTARRAY       *cprowIndices; /* compressed row indices used in the parallel SpMV */
264   PetscScalar          *alpha_one;    /* pointer to a device "scalar" storing the alpha parameter in the SpMV */
265   PetscScalar          *beta_zero;    /* pointer to a device "scalar" storing the beta parameter in the SpMV as zero*/
266   PetscScalar          *beta_one;     /* pointer to a device "scalar" storing the beta parameter in the SpMV as one */
267   hipsparseSpMatDescr_t matDescr;     /* descriptor for the matrix, used by SpMV and SpMM */
268   Mat_HipsparseSpMV     hipSpMV[3];   /* different Mat_CusparseSpMV structs for non-transpose, transpose, conj-transpose */
269   Mat_SeqAIJHIPSPARSEMultStruct() : matDescr(NULL)
270   {
271     for (int i = 0; i < 3; i++) hipSpMV[i].initialized = PETSC_FALSE;
272   }
273 };
274 
275 /* This is a larger struct holding all the matrices for a SpMV, and SpMV Transpose */
276 struct Mat_SeqAIJHIPSPARSE {
277   Mat_SeqAIJHIPSPARSEMultStruct *mat;               /* pointer to the matrix on the GPU */
278   Mat_SeqAIJHIPSPARSEMultStruct *matTranspose;      /* pointer to the matrix on the GPU (for the transpose ... useful for BiCG) */
279   THRUSTARRAY                   *workVector;        /* pointer to a workvector to which we can copy the relevant indices of a vector we want to multiply */
280   THRUSTINTARRAY32              *rowoffsets_gpu;    /* rowoffsets on GPU in non-compressed-row format. It is used to convert CSR to CSC */
281   PetscInt                       nrows;             /* number of rows of the matrix seen by GPU */
282   MatHIPSPARSEStorageFormat      format;            /* the storage format for the matrix on the device */
283   PetscBool                      use_cpu_solve;     /* Use AIJ_Seq (I)LU solve */
284   hipStream_t                    stream;            /* a stream for the parallel SpMV ... this is not owned and should not be deleted */
285   hipsparseHandle_t              handle;            /* a handle to the cusparse library ... this may not be owned (if we're working in parallel i.e. multiGPUs) */
286   PetscObjectState               nonzerostate;      /* track nonzero state to possibly recreate the GPU matrix */
287   size_t                         csr2cscBufferSize; /* stuff used to compute the matTranspose above */
288   void                          *csr2cscBuffer;     /* This is used as a C struct and is calloc'ed by PetscNewLog() */
289                                                     //  hipsparseCsr2CscAlg_t         csr2cscAlg; /* algorithms can be selected from command line options */
290   hipsparseSpMVAlg_t         spmvAlg;
291   hipsparseSpMMAlg_t         spmmAlg;
292   THRUSTINTARRAY            *csr2csc_i;
293   PetscSplitCSRDataStructure deviceMat; /* Matrix on device for, eg, assembly */
294   THRUSTINTARRAY            *cooPerm;   /* permutation array that sorts the input coo entris by row and col */
295   THRUSTINTARRAY            *cooPerm_a; /* ordered array that indicate i-th nonzero (after sorting) is the j-th unique nonzero */
296 
297   /* Stuff for extended COO support */
298   PetscBool   use_extended_coo; /* Use extended COO format */
299   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 */
300   PetscCount *perm_d;
301 
302   Mat_SeqAIJHIPSPARSE() : use_extended_coo(PETSC_FALSE), perm_d(NULL), jmap_d(NULL) { }
303 };
304 
305 typedef struct Mat_SeqAIJHIPSPARSETriFactors *Mat_SeqAIJHIPSPARSETriFactors_p;
306 
307 PETSC_INTERN PetscErrorCode MatSeqAIJHIPSPARSECopyToGPU(Mat);
308 PETSC_INTERN PetscErrorCode MatSetPreallocationCOO_SeqAIJHIPSPARSE_Basic(Mat, PetscCount, PetscInt[], PetscInt[]);
309 PETSC_INTERN PetscErrorCode MatSetValuesCOO_SeqAIJHIPSPARSE_Basic(Mat, const PetscScalar[], InsertMode);
310 PETSC_INTERN PetscErrorCode MatSeqAIJHIPSPARSEMergeMats(Mat, Mat, MatReuse, Mat *);
311 PETSC_INTERN PetscErrorCode MatSeqAIJHIPSPARSETriFactors_Reset(Mat_SeqAIJHIPSPARSETriFactors_p *);
312 
313 using VecSeq_HIP = Petsc::vec::cupm::impl::VecSeq_CUPM<Petsc::device::cupm::DeviceType::HIP>;
314 
315 static inline bool isHipMem(const void *data)
316 {
317   using namespace Petsc::device::cupm;
318   auto mtype = PETSC_MEMTYPE_HOST;
319 
320   PetscFunctionBegin;
321   PetscCallAbort(PETSC_COMM_SELF, impl::Interface<DeviceType::HIP>::PetscCUPMGetMemType(data, &mtype));
322   PetscFunctionReturn(PetscMemTypeDevice(mtype));
323 }
324 
325 #endif // PETSC_HIPSPARSEIMPL_H
326