1687625d7SJacob Faibussowitsch #ifndef PETSC_CUSPARSEMATIMPL_H 2687625d7SJacob Faibussowitsch #define PETSC_CUSPARSEMATIMPL_H 39ae82921SPaul Mullowney 4afb2bd1cSJunchao Zhang #include <petscpkg_version.h> 56d54fb17SJacob Faibussowitsch #include <../src/vec/vec/impls/seq/cupm/vecseqcupm.hpp> /* for VecSeq_CUPM */ 66d54fb17SJacob Faibussowitsch #include <petsc/private/petsclegacycupmblas.h> 79ae82921SPaul Mullowney 89ae82921SPaul Mullowney #include <cusparse_v2.h> 99ae82921SPaul Mullowney 109ae82921SPaul Mullowney #include <algorithm> 119ae82921SPaul Mullowney #include <vector> 129ae82921SPaul Mullowney 13c41cb2e2SAlejandro Lamas Daviña #include <thrust/device_vector.h> 14c41cb2e2SAlejandro Lamas Daviña #include <thrust/device_ptr.h> 1550ab121bSSatish Balay #include <thrust/device_malloc_allocator.h> 16c41cb2e2SAlejandro Lamas Daviña #include <thrust/transform.h> 17c41cb2e2SAlejandro Lamas Daviña #include <thrust/functional.h> 18554b8892SKarl Rupp #include <thrust/sequence.h> 197eaca502SStefano Zampini #include <thrust/system/system_error.h> 20c41cb2e2SAlejandro Lamas Daviña 21aa372e3fSPaul Mullowney #if defined(PETSC_USE_COMPLEX) 22aa372e3fSPaul Mullowney #if defined(PETSC_USE_REAL_SINGLE) 23ccdfe979SStefano Zampini const cuComplex PETSC_CUSPARSE_ONE = {1.0f, 0.0f}; 24ccdfe979SStefano Zampini const cuComplex PETSC_CUSPARSE_ZERO = {0.0f, 0.0f}; 25da112707SJunchao Zhang #define cusparseXcsrilu02_bufferSize(a, b, c, d, e, f, g, h, i) cusparseCcsrilu02_bufferSize(a, b, c, d, (cuComplex *)e, f, g, h, i) 26da112707SJunchao Zhang #define cusparseXcsrilu02_analysis(a, b, c, d, e, f, g, h, i, j) cusparseCcsrilu02_analysis(a, b, c, d, (cuComplex *)e, f, g, h, i, j) 27da112707SJunchao Zhang #define cusparseXcsrilu02(a, b, c, d, e, f, g, h, i, j) cusparseCcsrilu02(a, b, c, d, (cuComplex *)e, f, g, h, i, j) 28da112707SJunchao Zhang #define cusparseXcsric02_bufferSize(a, b, c, d, e, f, g, h, i) cusparseCcsric02_bufferSize(a, b, c, d, (cuComplex *)e, f, g, h, i) 29da112707SJunchao Zhang #define cusparseXcsric02_analysis(a, b, c, d, e, f, g, h, i, j) cusparseCcsric02_analysis(a, b, c, d, (cuComplex *)e, f, g, h, i, j) 30da112707SJunchao Zhang #define cusparseXcsric02(a, b, c, d, e, f, g, h, i, j) cusparseCcsric02(a, b, c, d, (cuComplex *)e, f, g, h, i, j) 31afb2bd1cSJunchao Zhang #elif defined(PETSC_USE_REAL_DOUBLE) 32afb2bd1cSJunchao Zhang const cuDoubleComplex PETSC_CUSPARSE_ONE = {1.0, 0.0}; 33afb2bd1cSJunchao Zhang const cuDoubleComplex PETSC_CUSPARSE_ZERO = {0.0, 0.0}; 34da112707SJunchao Zhang #define cusparseXcsrilu02_bufferSize(a, b, c, d, e, f, g, h, i) cusparseZcsrilu02_bufferSize(a, b, c, d, (cuDoubleComplex *)e, f, g, h, i) 35da112707SJunchao Zhang #define cusparseXcsrilu02_analysis(a, b, c, d, e, f, g, h, i, j) cusparseZcsrilu02_analysis(a, b, c, d, (cuDoubleComplex *)e, f, g, h, i, j) 36da112707SJunchao Zhang #define cusparseXcsrilu02(a, b, c, d, e, f, g, h, i, j) cusparseZcsrilu02(a, b, c, d, (cuDoubleComplex *)e, f, g, h, i, j) 37da112707SJunchao Zhang #define cusparseXcsric02_bufferSize(a, b, c, d, e, f, g, h, i) cusparseZcsric02_bufferSize(a, b, c, d, (cuDoubleComplex *)e, f, g, h, i) 38da112707SJunchao Zhang #define cusparseXcsric02_analysis(a, b, c, d, e, f, g, h, i, j) cusparseZcsric02_analysis(a, b, c, d, (cuDoubleComplex *)e, f, g, h, i, j) 39da112707SJunchao Zhang #define cusparseXcsric02(a, b, c, d, e, f, g, h, i, j) cusparseZcsric02(a, b, c, d, (cuDoubleComplex *)e, f, g, h, i, j) 40afb2bd1cSJunchao Zhang #endif 41afb2bd1cSJunchao Zhang #else 42afb2bd1cSJunchao Zhang const PetscScalar PETSC_CUSPARSE_ONE = 1.0; 43afb2bd1cSJunchao Zhang const PetscScalar PETSC_CUSPARSE_ZERO = 0.0; 44da112707SJunchao Zhang #if defined(PETSC_USE_REAL_SINGLE) 45da112707SJunchao Zhang #define cusparseXcsrilu02_bufferSize cusparseScsrilu02_bufferSize 46da112707SJunchao Zhang #define cusparseXcsrilu02_analysis cusparseScsrilu02_analysis 47da112707SJunchao Zhang #define cusparseXcsrilu02 cusparseScsrilu02 48da112707SJunchao Zhang #define cusparseXcsric02_bufferSize cusparseScsric02_bufferSize 49da112707SJunchao Zhang #define cusparseXcsric02_analysis cusparseScsric02_analysis 50da112707SJunchao Zhang #define cusparseXcsric02 cusparseScsric02 51da112707SJunchao Zhang #elif defined(PETSC_USE_REAL_DOUBLE) 52da112707SJunchao Zhang #define cusparseXcsrilu02_bufferSize cusparseDcsrilu02_bufferSize 53da112707SJunchao Zhang #define cusparseXcsrilu02_analysis cusparseDcsrilu02_analysis 54da112707SJunchao Zhang #define cusparseXcsrilu02 cusparseDcsrilu02 55da112707SJunchao Zhang #define cusparseXcsric02_bufferSize cusparseDcsric02_bufferSize 56da112707SJunchao Zhang #define cusparseXcsric02_analysis cusparseDcsric02_analysis 57da112707SJunchao Zhang #define cusparseXcsric02 cusparseDcsric02 58da112707SJunchao Zhang #endif 59afb2bd1cSJunchao Zhang #endif 60afb2bd1cSJunchao Zhang 611b0a6780SStefano Zampini #if PETSC_PKG_CUDA_VERSION_GE(9, 0, 0) 62261a78b4SJunchao Zhang #define csrsvInfo_t csrsv2Info_t 63261a78b4SJunchao Zhang #define cusparseCreateCsrsvInfo cusparseCreateCsrsv2Info 64261a78b4SJunchao Zhang #define cusparseDestroyCsrsvInfo cusparseDestroyCsrsv2Info 65afb2bd1cSJunchao Zhang #if defined(PETSC_USE_COMPLEX) 66afb2bd1cSJunchao Zhang #if defined(PETSC_USE_REAL_SINGLE) 67261a78b4SJunchao Zhang #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) 68261a78b4SJunchao Zhang #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) 69261a78b4SJunchao Zhang #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) 70afb2bd1cSJunchao Zhang #elif defined(PETSC_USE_REAL_DOUBLE) 71261a78b4SJunchao Zhang #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) 72261a78b4SJunchao Zhang #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) 73261a78b4SJunchao Zhang #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) 74afb2bd1cSJunchao Zhang #endif 75afb2bd1cSJunchao Zhang #else /* not complex */ 76afb2bd1cSJunchao Zhang #if defined(PETSC_USE_REAL_SINGLE) 77261a78b4SJunchao Zhang #define cusparseXcsrsv_buffsize cusparseScsrsv2_bufferSize 78261a78b4SJunchao Zhang #define cusparseXcsrsv_analysis cusparseScsrsv2_analysis 79261a78b4SJunchao Zhang #define cusparseXcsrsv_solve cusparseScsrsv2_solve 80afb2bd1cSJunchao Zhang #elif defined(PETSC_USE_REAL_DOUBLE) 81261a78b4SJunchao Zhang #define cusparseXcsrsv_buffsize cusparseDcsrsv2_bufferSize 82261a78b4SJunchao Zhang #define cusparseXcsrsv_analysis cusparseDcsrsv2_analysis 83261a78b4SJunchao Zhang #define cusparseXcsrsv_solve cusparseDcsrsv2_solve 84afb2bd1cSJunchao Zhang #endif 85afb2bd1cSJunchao Zhang #endif 869f7ba44dSJacob Faibussowitsch #else /* PETSC_PKG_CUDA_VERSION_GE(9, 0, 0) */ 87261a78b4SJunchao Zhang #define csrsvInfo_t cusparseSolveAnalysisInfo_t 88261a78b4SJunchao Zhang #define cusparseCreateCsrsvInfo cusparseCreateSolveAnalysisInfo 89261a78b4SJunchao Zhang #define cusparseDestroyCsrsvInfo cusparseDestroySolveAnalysisInfo 90afb2bd1cSJunchao Zhang #if defined(PETSC_USE_COMPLEX) 91afb2bd1cSJunchao Zhang #if defined(PETSC_USE_REAL_SINGLE) 929f7ba44dSJacob Faibussowitsch #define cusparseXcsrsv_solve(a, b, c, d_IGNORED, e, f, g, h, i, j, k, l, m_IGNORED, n_IGNORED) cusparseCcsrsv_solve((a), (b), (c), (cuComplex *)(e), (f), (cuComplex *)(g), (h), (i), (j), (cuComplex *)(k), (cuComplex *)(l)) 939f7ba44dSJacob Faibussowitsch #define cusparseXcsrsv_analysis(a, b, c, d, e, f, g, h, i, j_IGNORED, k_IGNORED) cusparseCcsrsv_analysis((a), (b), (c), (d), (e), (cuComplex *)(f), (g), (h), (i)) 941b0a6780SStefano Zampini #elif defined(PETSC_USE_REAL_DOUBLE) 959f7ba44dSJacob Faibussowitsch #define cusparseXcsrsv_solve(a, b, c, d_IGNORED, e, f, g, h, i, j, k, l, m_IGNORED, n_IGNORED) \ 969f7ba44dSJacob Faibussowitsch cusparseZcsrsv_solve((a), (b), (c), (cuDoubleComplex *)(e), (f), (cuDoubleComplex *)(g), (h), (i), (j), (cuDoubleComplex *)(k), (cuDoubleComplex *)(l)) 979f7ba44dSJacob Faibussowitsch #define cusparseXcsrsv_analysis(a, b, c, d, e, f, g, h, i, j_IGNORED, k_IGNORED) cusparseZcsrsv_analysis((a), (b), (c), (d), (e), (cuDoubleComplex *)(f), (g), (h), (i)) 981b0a6780SStefano Zampini #endif 991b0a6780SStefano Zampini #else /* not complex */ 1001b0a6780SStefano Zampini #if defined(PETSC_USE_REAL_SINGLE) 101261a78b4SJunchao Zhang #define cusparseXcsrsv_solve cusparseScsrsv_solve 1029f7ba44dSJacob Faibussowitsch #define cusparseXcsrsv_analysis(a, b, c, d, e, f, g, h, i, j, k) cusparseScsrsv_analysis(a, b, c, d, e, f, g, h, i) 1031b0a6780SStefano Zampini #elif defined(PETSC_USE_REAL_DOUBLE) 104261a78b4SJunchao Zhang #define cusparseXcsrsv_solve cusparseDcsrsv_solve 1059f7ba44dSJacob Faibussowitsch #define cusparseXcsrsv_analysis(a, b, c, d, e, f, g, h, i, j, k) cusparseDcsrsv_analysis(a, b, c, d, e, f, g, h, i) 1061b0a6780SStefano Zampini #endif 1071b0a6780SStefano Zampini #endif 1089f7ba44dSJacob Faibussowitsch #endif /* PETSC_PKG_CUDA_VERSION_GE(9, 0, 0) */ 1091b0a6780SStefano Zampini 1101b0a6780SStefano Zampini #if PETSC_PKG_CUDA_VERSION_GE(11, 0, 0) 1111b0a6780SStefano Zampini #define cusparse_csr2csc cusparseCsr2cscEx2 1121b0a6780SStefano Zampini #if defined(PETSC_USE_COMPLEX) 1131b0a6780SStefano Zampini #if defined(PETSC_USE_REAL_SINGLE) 1141b0a6780SStefano Zampini #define cusparse_scalartype CUDA_C_32F 115039c6fbaSStefano Zampini #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) 1169371c9d4SSatish Balay #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) \ 1179371c9d4SSatish Balay 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) 1181b0a6780SStefano Zampini #elif defined(PETSC_USE_REAL_DOUBLE) 1191b0a6780SStefano Zampini #define cusparse_scalartype CUDA_C_64F 1209371c9d4SSatish Balay #define cusparse_csr_spgeam(a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q, r, s, t) \ 1219371c9d4SSatish Balay 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) 1229371c9d4SSatish Balay #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) \ 1239371c9d4SSatish Balay 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) 1241b0a6780SStefano Zampini #endif 1251b0a6780SStefano Zampini #else /* not complex */ 1261b0a6780SStefano Zampini #if defined(PETSC_USE_REAL_SINGLE) 1271b0a6780SStefano Zampini #define cusparse_scalartype CUDA_R_32F 128039c6fbaSStefano Zampini #define cusparse_csr_spgeam cusparseScsrgeam2 129039c6fbaSStefano Zampini #define cusparse_csr_spgeam_bufferSize cusparseScsrgeam2_bufferSizeExt 1301b0a6780SStefano Zampini #elif defined(PETSC_USE_REAL_DOUBLE) 1311b0a6780SStefano Zampini #define cusparse_scalartype CUDA_R_64F 132039c6fbaSStefano Zampini #define cusparse_csr_spgeam cusparseDcsrgeam2 133039c6fbaSStefano Zampini #define cusparse_csr_spgeam_bufferSize cusparseDcsrgeam2_bufferSizeExt 1341b0a6780SStefano Zampini #endif 1351b0a6780SStefano Zampini #endif 1369f7ba44dSJacob Faibussowitsch #else /* PETSC_PKG_CUDA_VERSION_GE(11, 0, 0) */ 1371b0a6780SStefano Zampini #if defined(PETSC_USE_COMPLEX) 1381b0a6780SStefano Zampini #if defined(PETSC_USE_REAL_SINGLE) 139ec42abe4SAlejandro Lamas Daviña #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)) 140ccdfe979SStefano Zampini #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)) 141ec42abe4SAlejandro Lamas Daviña #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)) 142ec42abe4SAlejandro Lamas Daviña #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)) 143ec42abe4SAlejandro Lamas Daviña #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)) 144ec42abe4SAlejandro Lamas Daviña #define cusparse_hyb2csr(a, b, c, d, e, f) cusparseChyb2csr((a), (b), (c), (cuComplex *)(d), (e), (f)) 145fcdce8c4SStefano Zampini #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) 146039c6fbaSStefano Zampini #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) 147aa372e3fSPaul Mullowney #elif defined(PETSC_USE_REAL_DOUBLE) 148ec42abe4SAlejandro Lamas Daviña #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)) 1499371c9d4SSatish Balay #define cusparse_csr_spmm(a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p) \ 1509371c9d4SSatish Balay cusparseZcsrmm((a), (b), (c), (d), (e), (f), (cuDoubleComplex *)(g), (h), (cuDoubleComplex *)(i), (j), (k), (cuDoubleComplex *)(l), (m), (cuDoubleComplex *)(n), (cuDoubleComplex *)(o), (p)) 151ec42abe4SAlejandro Lamas Daviña #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)) 152ec42abe4SAlejandro Lamas Daviña #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)) 153ec42abe4SAlejandro Lamas Daviña #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)) 154ec42abe4SAlejandro Lamas Daviña #define cusparse_hyb2csr(a, b, c, d, e, f) cusparseZhyb2csr((a), (b), (c), (cuDoubleComplex *)(d), (e), (f)) 155fcdce8c4SStefano Zampini #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) 1569371c9d4SSatish Balay #define cusparse_csr_spgeam(a, b, c, d, e, f, g, h, i, j, k, l, m, n, o, p, q, r, s) \ 1579371c9d4SSatish Balay 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) 158aa372e3fSPaul Mullowney #endif 159aa372e3fSPaul Mullowney #else 160aa372e3fSPaul Mullowney #if defined(PETSC_USE_REAL_SINGLE) 161aa372e3fSPaul Mullowney #define cusparse_csr_spmv cusparseScsrmv 162ccdfe979SStefano Zampini #define cusparse_csr_spmm cusparseScsrmm 163aa372e3fSPaul Mullowney #define cusparse_csr2csc cusparseScsr2csc 164aa372e3fSPaul Mullowney #define cusparse_hyb_spmv cusparseShybmv 165aa372e3fSPaul Mullowney #define cusparse_csr2hyb cusparseScsr2hyb 166aa372e3fSPaul Mullowney #define cusparse_hyb2csr cusparseShyb2csr 167fcdce8c4SStefano Zampini #define cusparse_csr_spgemm cusparseScsrgemm 168039c6fbaSStefano Zampini #define cusparse_csr_spgeam cusparseScsrgeam 169aa372e3fSPaul Mullowney #elif defined(PETSC_USE_REAL_DOUBLE) 170aa372e3fSPaul Mullowney #define cusparse_csr_spmv cusparseDcsrmv 171ccdfe979SStefano Zampini #define cusparse_csr_spmm cusparseDcsrmm 172aa372e3fSPaul Mullowney #define cusparse_csr2csc cusparseDcsr2csc 173aa372e3fSPaul Mullowney #define cusparse_hyb_spmv cusparseDhybmv 174aa372e3fSPaul Mullowney #define cusparse_csr2hyb cusparseDcsr2hyb 175aa372e3fSPaul Mullowney #define cusparse_hyb2csr cusparseDhyb2csr 176fcdce8c4SStefano Zampini #define cusparse_csr_spgemm cusparseDcsrgemm 177039c6fbaSStefano Zampini #define cusparse_csr_spgeam cusparseDcsrgeam 178aa372e3fSPaul Mullowney #endif 179aa372e3fSPaul Mullowney #endif 1809f7ba44dSJacob Faibussowitsch #endif /* PETSC_PKG_CUDA_VERSION_GE(11, 0, 0) */ 181aa372e3fSPaul Mullowney 182aa372e3fSPaul Mullowney #define THRUSTINTARRAY32 thrust::device_vector<int> 183aa372e3fSPaul Mullowney #define THRUSTINTARRAY thrust::device_vector<PetscInt> 184aa372e3fSPaul Mullowney #define THRUSTARRAY thrust::device_vector<PetscScalar> 185aa372e3fSPaul Mullowney 186aa372e3fSPaul Mullowney /* A CSR matrix structure */ 187aa372e3fSPaul Mullowney struct CsrMatrix { 188aa372e3fSPaul Mullowney PetscInt num_rows; 189aa372e3fSPaul Mullowney PetscInt num_cols; 190aa372e3fSPaul Mullowney PetscInt num_entries; 191aa372e3fSPaul Mullowney THRUSTINTARRAY32 *row_offsets; 192aa372e3fSPaul Mullowney THRUSTINTARRAY32 *column_indices; 193aa372e3fSPaul Mullowney THRUSTARRAY *values; 1949ae82921SPaul Mullowney }; 1959ae82921SPaul Mullowney 196aa372e3fSPaul Mullowney /* This is struct holding the relevant data needed to a MatSolve */ 197aa372e3fSPaul Mullowney struct Mat_SeqAIJCUSPARSETriFactorStruct { 198aa372e3fSPaul Mullowney /* Data needed for triangular solve */ 199aa372e3fSPaul Mullowney cusparseMatDescr_t descr; 200aa372e3fSPaul Mullowney cusparseOperation_t solveOp; 201aa372e3fSPaul Mullowney CsrMatrix *csrMat; 202b917901dSJunchao Zhang #if PETSC_PKG_CUDA_VERSION_LT(11, 4, 0) 203261a78b4SJunchao Zhang csrsvInfo_t solveInfo; 2041b0a6780SStefano Zampini cusparseSolvePolicy_t solvePolicy; /* whether level information is generated and used */ 205d460d7bfSJunchao Zhang #endif 2061b0a6780SStefano Zampini int solveBufferSize; 2071b0a6780SStefano Zampini void *solveBuffer; 2081b0a6780SStefano Zampini size_t csr2cscBufferSize; /* to transpose the triangular factor (only used for CUDA >= 11.0) */ 2091b0a6780SStefano Zampini void *csr2cscBuffer; 2102cbc15d9SMark PetscScalar *AA_h; /* managed host buffer for moving values to the GPU */ 211aa372e3fSPaul Mullowney }; 212aa372e3fSPaul Mullowney 213afb2bd1cSJunchao Zhang /* This is a larger struct holding all the triangular factors for a solve, transpose solve, and any indices used in a reordering */ 214aa372e3fSPaul Mullowney struct Mat_SeqAIJCUSPARSETriFactors { 215b917901dSJunchao Zhang #if PETSC_PKG_CUDA_VERSION_LT(11, 4, 0) 216aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorPtr; /* pointer for lower triangular (factored matrix) on GPU */ 217aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorPtr; /* pointer for upper triangular (factored matrix) on GPU */ 218aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorPtrTranspose; /* pointer for lower triangular (factored matrix) on GPU for the transpose (useful for BiCG) */ 219aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorPtrTranspose; /* pointer for upper triangular (factored matrix) on GPU for the transpose (useful for BiCG)*/ 220d460d7bfSJunchao Zhang #endif 221d460d7bfSJunchao Zhang 222aa372e3fSPaul Mullowney THRUSTINTARRAY *rpermIndices; /* indices used for any reordering */ 223aa372e3fSPaul Mullowney THRUSTINTARRAY *cpermIndices; /* indices used for any reordering */ 224aa372e3fSPaul Mullowney THRUSTARRAY *workVector; 225aa372e3fSPaul Mullowney cusparseHandle_t handle; /* a handle to the cusparse library */ 226aa372e3fSPaul Mullowney PetscInt nnz; /* number of nonzeros ... need this for accurate logging between ICC and ILU */ 227e8d2b73aSMark Adams cudaDeviceProp dev_prop; 228e8d2b73aSMark Adams PetscBool init_dev_prop; 229da112707SJunchao Zhang 230d460d7bfSJunchao Zhang PetscBool factorizeOnDevice; /* Do factorization on device or not */ 231b917901dSJunchao Zhang #if PETSC_PKG_CUDA_VERSION_GE(11, 4, 0) 232da112707SJunchao Zhang /* csrilu0/csric0 appeared in cusparse-8.0, but we use it along with cusparseSpSV, 233da112707SJunchao Zhang which first appeared in cusparse-11.5 with cuda-11.3. 234da112707SJunchao Zhang */ 235d460d7bfSJunchao Zhang PetscScalar *csrVal, *diag; // the diagonal D in UtDU of Cholesky 236d460d7bfSJunchao Zhang int *csrRowPtr32, *csrColIdx32; // i,j of M. cusparseScsrilu02/ic02() etc require 32-bit indices 237d460d7bfSJunchao Zhang 238d460d7bfSJunchao Zhang PetscInt *csrRowPtr, *csrColIdx; // i, j of M on device for CUDA APIs that support 64-bit indices 239d460d7bfSJunchao Zhang 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 240d460d7bfSJunchao Zhang PetscInt *csrRowPtr_h; // csrColIdx_h is temporary, so it is not here 241da112707SJunchao Zhang 242da112707SJunchao Zhang /* Mixed mat descriptor types? yes, different cusparse APIs use different types */ 243da112707SJunchao Zhang cusparseMatDescr_t matDescr_M; 244da112707SJunchao Zhang cusparseSpMatDescr_t spMatDescr_L, spMatDescr_U; 245da112707SJunchao Zhang cusparseSpSVDescr_t spsvDescr_L, spsvDescr_Lt, spsvDescr_U, spsvDescr_Ut; 246da112707SJunchao Zhang 247da112707SJunchao Zhang cusparseDnVecDescr_t dnVecDescr_X, dnVecDescr_Y; 248da112707SJunchao Zhang PetscScalar *X, *Y; /* data array of dnVec X and Y */ 249da112707SJunchao Zhang 250da112707SJunchao Zhang /* Mixed size types? yes, CUDA-11.7.0 declared cusparseDcsrilu02_bufferSizeExt() that returns size_t but did not implement it! */ 25112ba2bc6SJunchao Zhang int factBufferSize_M; /* M ~= LU or LLt */ 252da112707SJunchao Zhang size_t spsvBufferSize_L, spsvBufferSize_Lt, spsvBufferSize_U, spsvBufferSize_Ut; 25312ba2bc6SJunchao Zhang /* cusparse needs various buffers for factorization and solve of L, U, Lt, or Ut. 25412ba2bc6SJunchao Zhang So save memory, we share the factorization buffer with one of spsvBuffer_L/U. 25512ba2bc6SJunchao Zhang */ 25612ba2bc6SJunchao Zhang void *factBuffer_M, *spsvBuffer_L, *spsvBuffer_U, *spsvBuffer_Lt, *spsvBuffer_Ut; 257da112707SJunchao Zhang 258da112707SJunchao Zhang csrilu02Info_t ilu0Info_M; 259da112707SJunchao Zhang csric02Info_t ic0Info_M; 260da112707SJunchao Zhang int structural_zero, numerical_zero; 261da112707SJunchao Zhang cusparseSolvePolicy_t policy_M; 262da112707SJunchao Zhang 26312ba2bc6SJunchao Zhang /* In MatSolveTranspose() for ILU0, we use the two flags to do on-demand solve */ 26412ba2bc6SJunchao Zhang PetscBool createdTransposeSpSVDescr; /* Have we created SpSV descriptors for Lt, Ut? */ 26512ba2bc6SJunchao Zhang PetscBool updatedTransposeSpSVAnalysis; /* Have we updated SpSV analysis with the latest L, U values? */ 266da112707SJunchao Zhang 267da112707SJunchao Zhang PetscLogDouble numericFactFlops; /* Estimated FLOPs in ILU0/ICC0 numeric factorization */ 268da112707SJunchao Zhang #endif 269aa372e3fSPaul Mullowney }; 270aa372e3fSPaul Mullowney 271afb2bd1cSJunchao Zhang struct Mat_CusparseSpMV { 272afb2bd1cSJunchao Zhang PetscBool initialized; /* Don't rely on spmvBuffer != NULL to test if the struct is initialized, */ 273afb2bd1cSJunchao Zhang size_t spmvBufferSize; /* since I'm not sure if smvBuffer can be NULL even after cusparseSpMV_bufferSize() */ 274afb2bd1cSJunchao Zhang void *spmvBuffer; 2759db3cbf9SStefano Zampini #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 */ 276afb2bd1cSJunchao Zhang cusparseDnVecDescr_t vecXDescr, vecYDescr; /* descriptor for the dense vectors in y=op(A)x */ 2779db3cbf9SStefano Zampini #endif 278afb2bd1cSJunchao Zhang }; 279afb2bd1cSJunchao Zhang 280afb2bd1cSJunchao Zhang /* This is struct holding the relevant data needed to a MatMult */ 281afb2bd1cSJunchao Zhang struct Mat_SeqAIJCUSPARSEMultStruct { 282afb2bd1cSJunchao Zhang void *mat; /* opaque pointer to a matrix. This could be either a cusparseHybMat_t or a CsrMatrix */ 283afb2bd1cSJunchao Zhang cusparseMatDescr_t descr; /* Data needed to describe the matrix for a multiply */ 284afb2bd1cSJunchao Zhang THRUSTINTARRAY *cprowIndices; /* compressed row indices used in the parallel SpMV */ 285afb2bd1cSJunchao Zhang PetscScalar *alpha_one; /* pointer to a device "scalar" storing the alpha parameter in the SpMV */ 286afb2bd1cSJunchao Zhang PetscScalar *beta_zero; /* pointer to a device "scalar" storing the beta parameter in the SpMV as zero*/ 287afb2bd1cSJunchao Zhang PetscScalar *beta_one; /* pointer to a device "scalar" storing the beta parameter in the SpMV as one */ 288afb2bd1cSJunchao Zhang #if PETSC_PKG_CUDA_VERSION_GE(11, 0, 0) 289afb2bd1cSJunchao Zhang cusparseSpMatDescr_t matDescr; /* descriptor for the matrix, used by SpMV and SpMM */ 290afb2bd1cSJunchao Zhang Mat_CusparseSpMV cuSpMV[3]; /* different Mat_CusparseSpMV structs for non-transpose, transpose, conj-transpose */ 291d71ae5a4SJacob Faibussowitsch Mat_SeqAIJCUSPARSEMultStruct() : matDescr(NULL) 292d71ae5a4SJacob Faibussowitsch { 293afb2bd1cSJunchao Zhang for (int i = 0; i < 3; i++) cuSpMV[i].initialized = PETSC_FALSE; 294afb2bd1cSJunchao Zhang } 295afb2bd1cSJunchao Zhang #endif 296afb2bd1cSJunchao Zhang }; 297afb2bd1cSJunchao Zhang 2987301b172SPierre Jolivet /* This is a larger struct holding all the matrices for a SpMV, and SpMV Transpose */ 2999ae82921SPaul Mullowney struct Mat_SeqAIJCUSPARSE { 300aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSEMultStruct *mat; /* pointer to the matrix on the GPU */ 301aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSEMultStruct *matTranspose; /* pointer to the matrix on the GPU (for the transpose ... useful for BiCG) */ 302aa372e3fSPaul Mullowney THRUSTARRAY *workVector; /* pointer to a workvector to which we can copy the relevant indices of a vector we want to multiply */ 303029808eeSJunchao Zhang THRUSTINTARRAY32 *rowoffsets_gpu; /* rowoffsets on GPU in non-compressed-row format. It is used to convert CSR to CSC */ 304213423ffSJunchao Zhang PetscInt nrows; /* number of rows of the matrix seen by GPU */ 305e057df02SPaul Mullowney MatCUSPARSEStorageFormat format; /* the storage format for the matrix on the device */ 306365b711fSMark Adams PetscBool use_cpu_solve; /* Use AIJ_Seq (I)LU solve */ 307aa372e3fSPaul Mullowney cudaStream_t stream; /* a stream for the parallel SpMV ... this is not owned and should not be deleted */ 308aa372e3fSPaul Mullowney cusparseHandle_t handle; /* a handle to the cusparse library ... this may not be owned (if we're working in parallel i.e. multiGPUs) */ 30954da937aSStefano Zampini PetscObjectState nonzerostate; /* track nonzero state to possibly recreate the GPU matrix */ 310afb2bd1cSJunchao Zhang #if PETSC_PKG_CUDA_VERSION_GE(11, 0, 0) 311afb2bd1cSJunchao Zhang size_t csr2cscBufferSize; /* stuff used to compute the matTranspose above */ 312fea72eb4SStefano Zampini void *csr2cscBuffer; /* This is used as a C struct and is calloc'ed by PetscNew() */ 313afb2bd1cSJunchao Zhang cusparseCsr2CscAlg_t csr2cscAlg; /* algorithms can be selected from command line options */ 314afb2bd1cSJunchao Zhang cusparseSpMVAlg_t spmvAlg; 315afb2bd1cSJunchao Zhang cusparseSpMMAlg_t spmmAlg; 316afb2bd1cSJunchao Zhang #endif 317a49f1ed0SStefano Zampini THRUSTINTARRAY *csr2csc_i; 318*2c4ab24aSJunchao Zhang THRUSTINTARRAY *coords; /* permutation array used in MatSeqAIJCUSPARSEMergeMats */ 3199ae82921SPaul Mullowney }; 3209ae82921SPaul Mullowney 32147d993e7Ssuyashtn typedef struct Mat_SeqAIJCUSPARSETriFactors *Mat_SeqAIJCUSPARSETriFactors_p; 32247d993e7Ssuyashtn 323d3ecb6c7SJed Brown PETSC_INTERN PetscErrorCode MatSeqAIJCUSPARSECopyToGPU(Mat); 324ed502f03SStefano Zampini PETSC_INTERN PetscErrorCode MatSeqAIJCUSPARSEMergeMats(Mat, Mat, MatReuse, Mat *); 3255f101d05SStefano Zampini PETSC_INTERN PetscErrorCode MatSeqAIJCUSPARSETriFactors_Reset(Mat_SeqAIJCUSPARSETriFactors_p *); 326ed502f03SStefano Zampini 3276d54fb17SJacob Faibussowitsch using VecSeq_CUDA = Petsc::vec::cupm::impl::VecSeq_CUPM<Petsc::device::cupm::DeviceType::CUDA>; 3286d54fb17SJacob Faibussowitsch 329d71ae5a4SJacob Faibussowitsch static inline bool isCudaMem(const void *data) 330d71ae5a4SJacob Faibussowitsch { 3316d54fb17SJacob Faibussowitsch using namespace Petsc::device::cupm; 3326d54fb17SJacob Faibussowitsch auto mtype = PETSC_MEMTYPE_HOST; 3336d54fb17SJacob Faibussowitsch 3346d54fb17SJacob Faibussowitsch PetscFunctionBegin; 3356d54fb17SJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF, impl::Interface<DeviceType::CUDA>::PetscCUPMGetMemType(data, &mtype)); 3366d54fb17SJacob Faibussowitsch PetscFunctionReturn(PetscMemTypeDevice(mtype)); 33708391a17SStefano Zampini } 33808391a17SStefano Zampini 339687625d7SJacob Faibussowitsch #endif // PETSC_CUSPARSEMATIMPL_H 340