cuSPARSE
Public
cuSPARSE.CuSparseMatrix — Type
Utility union type of CuSparseMatrixCSC, CuSparseMatrixCSR, CuSparseMatrixBSR, CuSparseMatrixCOO.
cuSPARSE.CuSparseMatrixBSR — Type
CuSparseMatrixBSRContainer to hold sparse matrices in block compressed sparse row (BSR) format on the GPU. BSR format is also used in Intel MKL, and is suited to matrices that are "block" sparse - rare blocks of non-sparse regions.
cuSPARSE.CuSparseMatrixCOO — Type
CuSparseMatrixCOOContainer to hold sparse matrices in coordinate (COO) format on the GPU. COO format is mainly useful to initially construct sparse matrices, afterwards switch to CuSparseMatrixCSR for more functionality.
cuSPARSE.CuSparseMatrixCSC — Type
CuSparseMatrixCSCContainer to hold sparse matrices in compressed sparse column (CSC) format on the GPU.
cuSPARSE.CuSparseMatrixCSR — Type
CuSparseMatrixCSR{Tv, Ti} <: AbstractCuSparseMatrix{Tv, Ti}Container to hold sparse matrices in compressed sparse row (CSR) format on the GPU.
cuSPARSE.axpby! — Method
axpby!(alpha::Number, X::CuSparseVector, beta::Number, Y::CuVector, index::SparseChar)Computes alpha * X + beta * Y for sparse X and dense Y.
cuSPARSE.axpby — Method
axpby(alpha::Number, x::CuSparseVector, beta::Number, y::CuSparseVector, index::SparseChar)Performs z = alpha * x + beta * y. x and y are sparse vectors.
cuSPARSE.color — Function
color(A::CuSparseMatrixCSC, index::SparseChar; percentage::Number=1.0)
color(A::CuSparseMatrixCSR, index::SparseChar; percentage::Number=1.0)This function performs the coloring of the adjacency graph associated with the matrix A. The coloring is an assignment of colors (integer numbers) to nodes, such that neighboring nodes have distinct colors. An approximate coloring algorithm is used in this routine, and is stopped when a certain percentage of nodes has been colored. The rest of the nodes are assigned distinct colors (an increasing sequence of integers numbers, starting from the last integer used previously). The reordering is such that nodes that have been assigned the same color are reordered to be next to each other.
The matrix A passed to this routine, must be stored as a general matrix and have a symmetric sparsity pattern. If the matrix is non-symmetric the user should pass A + Aᵀ as a parameter to this routine.
cuSPARSE.gather! — Method
gather!(X::CuSparseVector, Y::CuVector, index::SparseChar)Sets the nonzero elements of X equal to the nonzero elements of Y at the same indices.
cuSPARSE.geam — Method
geam(alpha::Number, A::CuSparseMatrix, beta::Number, B::CuSparseMatrix, index::SparseChar)Performs C = alpha * A + beta * B. A and B are sparse matrices defined in CSR or CSC storage formats.
cuSPARSE.gemv — Function
y = gemv(transa, alpha, A, x, index, [algo])Perform a product between a CuSparseMatrix and a CuSparseVector, returning a CuSparseVector. This function should only be used for highly sparse matrices and vectors, as the result is expected to have many non-zeros in practice. For this reason, high-level functions like mul! and * internally convert the sparse vector into a dense vector to use a more efficient CUSPARSE routine.
Supported formats for the sparse matrix are CuSparseMatrixCSC and CuSparseMatrixCSR.
cuSPARSE.gtsv2! — Function
gtsv2!(dl::CuVector, d::CuVector, du::CuVector, B::CuVecOrMat, index::SparseChar='O'; pivoting::Bool=true)Solve the linear system A * X = B where A is a tridiagonal matrix defined by three vectors corresponding to its lower (dl), main (d), and upper (du) diagonals. With pivoting, the solution is more accurate but also more expensive. Note that the solution X overwrites the right-hand side B.
cuSPARSE.ic02! — Function
ic02!(A::CuSparseMatrix, index::SparseChar='O')Incomplete Cholesky factorization with no pivoting. Preserves the sparse layout of matrix A.
cuSPARSE.ilu02! — Function
ilu02!(A::CuSparseMatrix, index::SparseChar='O')Incomplete LU factorization with no pivoting. Preserves the sparse layout of matrix A.
cuSPARSE.mm! — Function
mm!(transa::SparseChar, transb::SparseChar, alpha::Number, A::CuSparseMatrix, B::CuMatrix, beta::Number, C::CuMatrix, index::SparseChar)
mm!(transa::SparseChar, transb::SparseChar, alpha::Number, A::CuMatrix, B::Union{CuSparseMatrixCSC,CuSparseMatrixCSR,CuSparseMatrixCOO}, beta::Number, C::CuMatrix, index::SparseChar)Performs C = alpha * op(A) * op(B) + beta * C, where op can be nothing (transa = N), tranpose (transa = T) or conjugate transpose (transa = C).
cuSPARSE.mv! — Function
mv!(transa::SparseChar, alpha::Number, A::CuSparseMatrix, X::CuVector, beta::Number, Y::CuVector, index::SparseChar)Performs Y = alpha * op(A) * X + beta * Y, where op can be nothing (transa = N), tranpose (transa = T) or conjugate transpose (transa = C). X and Y are dense vectors.
cuSPARSE.rot! — Method
rot!(X::CuSparseVector, Y::CuVector, c::Number, s::Number, index::SparseChar)Performs the Givens rotation specified by c and s to sparse X and dense Y.
cuSPARSE.scatter! — Method
scatter!(Y::CuVector, X::CuSparseVector, index::SparseChar)Set Y[:] = X[:] for dense Y and sparse X.
cuSPARSE.sm2! — Method
sm2!(transa::SparseChar, transxy::SparseChar, uplo::SparseChar, diag::SparseChar, alpha::BlasFloat, A::CuSparseMatrixBSR, X::CuMatrix, index::SparseChar)Performs X = alpha * op(A) \ op(X), where op can be nothing (transa = N), tranpose (transa = T) or conjugate transpose (transa = C). X is a dense matrix, and uplo tells sm2! which triangle of the block sparse matrix A to reference. If the triangle has unit diagonal, set diag to 'U'.
cuSPARSE.sv2! — Method
sv2!(transa::SparseChar, uplo::SparseChar, diag::SparseChar, alpha::BlasFloat, A::CuSparseMatrixBSR, X::CuVector, index::SparseChar)Performs X = alpha * op(A) \ X, where op can be nothing (transa = N), tranpose (transa = T) or conjugate transpose (transa = C). X is a dense vector, and uplo tells sv2! which triangle of the block sparse matrix A to reference. If the triangle has unit diagonal, set diag to 'U'.
Private
SparseArrays.sparse — Method
sparse(x::DenseCuMatrix; fmt=:csc)
sparse(I::CuVector, J::CuVector, V::CuVector, [m, n]; fmt=:csc)Return a sparse cuda matrix, with type determined by fmt. Possible formats are :csc, :csr, :bsr, and :coo.
cuSPARSE.chkbmmdims — Method
check that the dimensions of arrays B and C make sense for a batched matrix-matrix multiplication
cuSPARSE.chkmmdims — Method
check that the dimensions of matrices B and C make sense for a multiplication
cuSPARSE.chkmvdims — Method
check that the dimensions of matrix X and vector Y make sense for a multiplication