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