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