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 <../src/sys/objects/device/impls/cupm/cupmthrustutility.hpp> 10 #include <petsc/private/petsclegacycupmblas.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 hipDeviceProp_t dev_prop; 211 PetscBool init_dev_prop; 212 213 /* csrilu0/csric0 appeared in earlier versions of AMD ROCm^{TM}, but we use it along with hipsparseSpSV, 214 which first appeared in hipsparse with ROCm-4.5.0. 215 */ 216 PetscBool factorizeOnDevice; /* Do factorization on device or not */ 217 #if PETSC_PKG_HIP_VERSION_GE(4, 5, 0) 218 PetscScalar *csrVal; 219 int *csrRowPtr, *csrColIdx; /* a,i,j of M. Using int since some hipsparse APIs only support 32-bit indices */ 220 221 /* Mixed mat descriptor types? yes, different hipsparse APIs use different types */ 222 hipsparseMatDescr_t matDescr_M; 223 hipsparseSpMatDescr_t spMatDescr_L, spMatDescr_U; 224 hipsparseSpSVDescr_t spsvDescr_L, spsvDescr_Lt, spsvDescr_U, spsvDescr_Ut; 225 226 hipsparseDnVecDescr_t dnVecDescr_X, dnVecDescr_Y; 227 PetscScalar *X, *Y; /* data array of dnVec X and Y */ 228 229 /* Mixed size types? yes */ 230 int factBufferSize_M; /* M ~= LU or LLt */ 231 size_t spsvBufferSize_L, spsvBufferSize_Lt, spsvBufferSize_U, spsvBufferSize_Ut; 232 /* hipsparse needs various buffers for factorization and solve of L, U, Lt, or Ut. 233 To save memory, we share the factorization buffer with one of spsvBuffer_L/U. 234 */ 235 void *factBuffer_M, *spsvBuffer_L, *spsvBuffer_U, *spsvBuffer_Lt, *spsvBuffer_Ut; 236 237 csrilu02Info_t ilu0Info_M; 238 csric02Info_t ic0Info_M; 239 int structural_zero, numerical_zero; 240 hipsparseSolvePolicy_t policy_M; 241 242 /* In MatSolveTranspose() for ILU0, we use the two flags to do on-demand solve */ 243 PetscBool createdTransposeSpSVDescr; /* Have we created SpSV descriptors for Lt, Ut? */ 244 PetscBool updatedTransposeSpSVAnalysis; /* Have we updated SpSV analysis with the latest L, U values? */ 245 246 PetscLogDouble numericFactFlops; /* Estimated FLOPs in ILU0/ICC0 numeric factorization */ 247 #endif 248 }; 249 250 struct Mat_HipsparseSpMV { 251 PetscBool initialized; /* Don't rely on spmvBuffer != NULL to test if the struct is initialized, */ 252 size_t spmvBufferSize; /* since I'm not sure if smvBuffer can be NULL even after hipsparseSpMV_bufferSize() */ 253 void *spmvBuffer; 254 hipsparseDnVecDescr_t vecXDescr, vecYDescr; /* descriptor for the dense vectors in y=op(A)x */ 255 }; 256 257 /* This is struct holding the relevant data needed to a MatMult */ 258 struct Mat_SeqAIJHIPSPARSEMultStruct { 259 void *mat; /* opaque pointer to a matrix. This could be either a hipsparseHybMat_t or a CsrMatrix */ 260 hipsparseMatDescr_t descr; /* Data needed to describe the matrix for a multiply */ 261 THRUSTINTARRAY *cprowIndices; /* compressed row indices used in the parallel SpMV */ 262 PetscScalar *alpha_one; /* pointer to a device "scalar" storing the alpha parameter in the SpMV */ 263 PetscScalar *beta_zero; /* pointer to a device "scalar" storing the beta parameter in the SpMV as zero*/ 264 PetscScalar *beta_one; /* pointer to a device "scalar" storing the beta parameter in the SpMV as one */ 265 hipsparseSpMatDescr_t matDescr; /* descriptor for the matrix, used by SpMV and SpMM */ 266 Mat_HipsparseSpMV hipSpMV[3]; /* different Mat_CusparseSpMV structs for non-transpose, transpose, conj-transpose */ 267 Mat_SeqAIJHIPSPARSEMultStruct() : matDescr(NULL) 268 { 269 for (int i = 0; i < 3; i++) hipSpMV[i].initialized = PETSC_FALSE; 270 } 271 }; 272 273 /* This is a larger struct holding all the matrices for a SpMV, and SpMV Transpose */ 274 struct Mat_SeqAIJHIPSPARSE { 275 Mat_SeqAIJHIPSPARSEMultStruct *mat; /* pointer to the matrix on the GPU */ 276 Mat_SeqAIJHIPSPARSEMultStruct *matTranspose; /* pointer to the matrix on the GPU (for the transpose ... useful for BiCG) */ 277 THRUSTARRAY *workVector; /* pointer to a workvector to which we can copy the relevant indices of a vector we want to multiply */ 278 THRUSTINTARRAY32 *rowoffsets_gpu; /* rowoffsets on GPU in non-compressed-row format. It is used to convert CSR to CSC */ 279 PetscInt nrows; /* number of rows of the matrix seen by GPU */ 280 MatHIPSPARSEStorageFormat format; /* the storage format for the matrix on the device */ 281 PetscBool use_cpu_solve; /* Use AIJ_Seq (I)LU solve */ 282 hipStream_t stream; /* a stream for the parallel SpMV ... this is not owned and should not be deleted */ 283 hipsparseHandle_t handle; /* a handle to the cusparse library ... this may not be owned (if we're working in parallel i.e. multiGPUs) */ 284 PetscObjectState nonzerostate; /* track nonzero state to possibly recreate the GPU matrix */ 285 size_t csr2cscBufferSize; /* stuff used to compute the matTranspose above */ 286 void *csr2cscBuffer; /* This is used as a C struct and is calloc'ed by PetscNewLog() */ 287 // hipsparseCsr2CscAlg_t csr2cscAlg; /* algorithms can be selected from command line options */ 288 hipsparseSpMVAlg_t spmvAlg; 289 hipsparseSpMMAlg_t spmmAlg; 290 THRUSTINTARRAY *csr2csc_i; 291 THRUSTINTARRAY *coords; /* permutation array used in MatSeqAIJHIPSPARSEMergeMats */ 292 }; 293 294 typedef struct Mat_SeqAIJHIPSPARSETriFactors *Mat_SeqAIJHIPSPARSETriFactors_p; 295 296 PETSC_INTERN PetscErrorCode MatSeqAIJHIPSPARSECopyToGPU(Mat); 297 PETSC_INTERN PetscErrorCode MatSeqAIJHIPSPARSEMergeMats(Mat, Mat, MatReuse, Mat *); 298 PETSC_INTERN PetscErrorCode MatSeqAIJHIPSPARSETriFactors_Reset(Mat_SeqAIJHIPSPARSETriFactors_p *); 299 300 using VecSeq_HIP = Petsc::vec::cupm::impl::VecSeq_CUPM<Petsc::device::cupm::DeviceType::HIP>; 301 302 static inline bool isHipMem(const void *data) 303 { 304 using namespace Petsc::device::cupm; 305 auto mtype = PETSC_MEMTYPE_HOST; 306 307 PetscFunctionBegin; 308 PetscCallAbort(PETSC_COMM_SELF, impl::Interface<DeviceType::HIP>::PetscCUPMGetMemType(data, &mtype)); 309 PetscFunctionReturn(PetscMemTypeDevice(mtype)); 310 } 311 312 #endif // PETSC_HIPSPARSEIMPL_H 313