cuSPARSE

Provides basic linear algebra operations for sparse matrices. See NVIDIA cuSPARSE for an in-depth description of the cuSPARSE library and its methods and data types. All functions are accessed through the pyculib.sparse.Sparse class:

class pyculib.sparse.Sparse(idxbase=0)

All cuSPARSE functions are available under the Sparse object.

Parameters:idxbase – The base for indexing, either 0 or 1. Optional, defaults to 0.

Similarly to the cuBLAS interface, no special naming convention is used for functions to operate on different datatypes - all datatypes are handled by each function, and dispatch of the corresponding library function is handled by Pyculib. However, it is often necessary to provide a matrix descriptor to functions, which provides some information about the format and properties of a matrix. A matrix descriptor can be obtained from the pyculib.sparse.Sparse.matdescr() method:

pyculib.sparse.Sparse.matdescr(indexbase, diagtype, fillmode, matrixtype)

Creates a matrix descriptor that describes a matrix with the given indexbase, diagtype, fillmode, and matrixtype. Note that not all of these options are relevant to every matrix storage format.

Parameters:
  • indexbase – Optional. 0 for 0-based indexing, or 1 for 1-based indexing. If not specified, the default given to the pyculib.sparse.Sparse constructor is used instead.
  • diagtype – Optional. Defaults to ‘N’. ‘N’ signifies that the matrix diagonal has non-unit elements. ‘U’ signifies that the matrix diagonal only contains unit elements.
  • fillmode – Optional. Defaults to ‘L’. ‘L’ indicates that the lower triangular part of the matrix is stored. ‘U’ indicates that the upper triangular part of the matrix is stored.
  • matrixtype – Optional. Defaults to ‘G’. ‘S’ indicates that the matrix is symmetric. ‘H’ indicates that the matrix is Hermitian. ‘T’ indicates that the matrix is triangular. ‘G’ is used for a general matrix, which is not symmetric, Hermitian, or triangular.
Returns:

A matrix descriptor.

Many of the methods of the pyculib.sparse.Sparse class accept the individual data structures that make up a sparse representation of a matrix (for example the values, the row pointers and the column indices for a CSR format matrix). However, some methods (such as pyculib.sparse.Sparse.csrgemm_ez()), accept an instance of the pyculib.sparse.CudaSparseMatrix class:

class pyculib.sparse.CudaSparseMatrix

Base class for a representation of a sparse matrix on a CUDA device. The constructor takes no arguments.

from_host_matrix(matrix, stream)

Initialise the matrix structure and values from an instance of a matrix on the host. The host matrix must be of the corresponding host type, which is documented for each subclass below.

copy_to_host(stream)

Create an instance of the corresponding host matrix type and copy the matrix structure and data into it from the device. See subclass documentation for an indication of the corresponding matrix type.

Subclasses of the sparse matrix type are:

class pyculib.sparse.CudaBSRMatrix

CUDA sparse matrix for which the corresponding type is a scipy.sparse.bsr_matrix.

class pyculib.sparse.CudaCSRMatrix

CUDA sparse matrix for which the corresponding type is a scipy.sparse.csr_matrix.

class pyculib.sparse.CudaCSCMatrix

CUDA sparse matrix for which the corresponding type is a scipy.sparse.csc_matrix.

There are also some convenience methods for constructing CUDA sparse matrices in a similar manner to Scipy sparse matrices:

sparse.bsr_matrix(**kws)

Takes the same arguments as scipy.sparse.bsr_matrix.

Returns a BSR CUDA matrix.

sparse.csr_matrix(**kws)

Takes the same arguments as scipy.sparse.csr_matrix.

Returns a CSR CUDA matrix.

sparse.csc_matrix(**kws)

Takes the same arguments as scipy.sparse.csc_matrix.

Returns a CSC CUDA matrix.

BLAS Level 1

pyculib.sparse.Sparse.axpyi(alpha, xVal, xInd, y)

Multiplies the sparse vector x by alpha and adds the result to the dense vector y.

Parameters:
  • alpha – scalar
  • xVal – vector of non-zero values of x
  • xInd – vector of indices of non-zero values of x
  • y – dense vector
Returns:

dense vector

pyculib.sparse.Sparse.doti(xVal, xInd, y)

Computes the dot product of the sparse vector x and dense vector y.

Parameters:
  • xVal – vector of non-zero values of x
  • xInd – vector of indices of non-zero values of x
  • y – dense vector
Returns:

scalar

pyculib.sparse.Sparse.dotci(xVal, xInd, y)

Computes the dot product of the complex conjugate of the sparse vector x and the dense vector y.

Parameters:
  • xVal – vector of non-zero values of x
  • xInd – vector of indices of non-zero values of x
  • y – dense vector
Returns:

scalar

pyculib.sparse.Sparse.gthr(y, xVal, xInd)

Gathers the elements of y at the indices xInd into the array xVal

Parameters:
  • xVal – vector of non-zero values of x
  • xInd – vector of indices of non-zero values of x
  • y – dense vector
Returns:

None

pyculib.sparse.Sparse.gthrz(y, xVal, xInd)

Gathers the elements of y at the indices xInd into the array xVal and zeroes out the gathered elements of y.

Parameters:
  • xVal – vector of non-zero values of x
  • xInd – vector of indices of non-zero values of x
  • y – dense vector
Returns:

None

pyculib.sparse.Sparse.roti(xVal, xInd, y, c, s)

Applies the Givens rotation matrix, G:

\[\begin{split}G = \left( \begin{array}{cc} C & S \\ -S & C \end{array}\right)\end{split}\]

to the sparse vector x and dense vector y.

Parameters:
  • xVal – vector of non-zero values of x
  • xInd – vector of indices of non-zero values of x
  • y – dense vector
  • c – cosine element of the rotation matrix
  • s – sine element of the rotation matrix
Returns:

None

pyculib.sparse.Sparse.sctr(xVal, xInd, y)

Scatters the elements of the sparse vector x into the dense vector y. Elements of y whose indices are not listed in xInd are unmodified.

Parameters:
  • xVal – vector of non-zero values of x
  • xInd – vector of indices of non-zero values of x
  • y – dense vector
Returns:

None

BLAS Level 2

All level 2 routines follow the following naming convention for the following arguments:

  • alpha, beta – (scalar) Can be real or complex numbers.

  • descr, descrA, descrB – (descriptor) Matrix descriptor. An appropriate descriptor may be obtained by calling pyculib.sparse.Sparse.matdescr(). descr only applies to the only matrix argument. descrA and descrB apply to matrix A and matrix B, respectively.

  • dir – (string) Can be ‘C’ to indicate column-major block storage or ‘R’ to indicate row-major block storage.

  • trans, transa, transb – (string)

    Select the operation op to apply to a matrix:

    • ‘N’: op(X) = X, the identity operation;
    • ‘T’: op(X) = X**T, the transpose;
    • ‘C’: op(X) = X**H, the conjugate transpose.

    trans only applies to the only matrix argument. transa and transb apply to matrix A and matrix B, respectively.

pyculib.sparse.Sparse.bsrmv_matrix(dir, trans, alpha, descr, bsrmat, x, beta, y)

Matrix-vector multiplication y = alpha * op(A) * x + beta * y with a BSR-format matrix.

Parameters:
  • dir – block storage direction
  • trans – operation to apply to the matrix
  • alpha – scalar
  • descr – matrix descriptor
  • bsrmat – the matrix A
  • x – dense vector
  • beta – scalar
  • y – dense vector
Returns:

None

pyculib.sparse.Sparse.bsrmv(dir, trans, mb, nb, nnzb, alpha, descr, bsrVal, bsrRowPtr, bsrColInd, blockDim, x, beta, y)

Matrix-vector multiplication y = alpha * op(A) * x + beta * y with a BSR-format matrix. This function accepts the individual arrays that make up the structure of a BSR matrix - if a pyculib.sparse.CudaBSRMatrix instance is to hand, it is recommended to use the bsrmv_matrix() method instead.

Parameters:
  • dir – block storage direction
  • trans – operation to apply to the matrix
  • mb – Number of block rows of the matrix
  • nb – Number of block columns of the matrix
  • nnzb – Number of nonzero blocks of the matrix
  • alpha – scalar
  • descr – matrix descriptor
  • bsrVal – vector of nonzero values of the matrix
  • bsrRowPtr – vector of block row pointers of the matrix
  • bsrColInd – vector of block column indices of the matrix
  • blockDim – block dimension of the matrix
  • x – dense vector
  • beta – scalar
  • y – dense vector
Returns:

None

pyculib.sparse.Sparse.bsrxmv(dir, trans, sizeOfMask, mb, nb, nnzb, alpha, descr, bsrVal, bsrMaskPtr, bsrRowPtr, bsrEndPtr, bsrColInd, blockDim, x, beta, y)

Matrix-vector multiplication similar to bsrmv(), but including a mask operation: y(mask) = (alpha * op(A) * x + beta * y)(mask). The blocks of y to be updated are specified in bsrMaskPtr. Blocks whose indices are not specified in bsrMaskPtr are left unmodified.

Parameters:
  • dir – block storage direction
  • trans – operation to apply to the matrix
  • sizeOfMask – number of updated blocks of rows of y
  • mb – Number of block rows of the matrix
  • nb – Number of block columns of the matrix
  • nnzb – Number of nonzero blocks of the matrix
  • alpha – scalar
  • descr – matrix descriptor
  • bsrVal – vector of nonzero values of the matrix
  • bsrMaskPtr – vector of indices of the block elements to be updated
  • bsrRowPtr – vector of block row pointers of the matrix
  • bsrEndPtr – vector of pointers to the end of every block row plus one
  • bsrColInd – vector of block column indices of the matrix
  • blockDim – block dimension of the matrix
  • x – dense vector
  • beta – scalar
  • y – dense vector
Returns:

None

pyculib.sparse.Sparse.csrmv(trans, m, n, nnz, alpha, descr, csrVal, csrRowPtr, csrColInd, x, beta, y)

Matrix-vector multiplication y = alpha * op(A) * x + beta * y with a CSR-format matrix.

Parameters:
  • trans – operation to apply to the matrix
  • m – Number of rows of the matrix
  • n – Number of columns of the matrix
  • nnz – Number of nonzeroes of the matrix
  • alpha – scalar
  • descr – matrix descriptor
  • csrVal – vector of nonzero values of the matrix
  • csrRowPtr – vector of row pointers of the matrix
  • csrColInd – vector of column indices of the matrix
  • x – dense vector
  • beta – scalar
  • y – dense vector
Returns:

None

pyculib.sparse.Sparse.csrsv_analysis(trans, m, nnz, descr, csrVal, csrRowPtr, csrColInd)

Performs the analysis phase of the solution of the sparse triangular linear system op(A) * y = alpha * x. This needs to be executed only once for a given matrix and operation type.

Parameters:
  • trans – operation to apply to the matrix
  • m – number of rows of the matrix
  • nnz – number of nonzeroes of the matrix
  • descr – matrix descriptor
  • csrVal – vector of nonzero values of the matrix
  • csrRowPtr – vector of row pointers of the matrix
  • csrColInd – vector of column indices of the matrix
Returns:

the analysis result, which can be used as input to the solve phase

pyculib.sparse.Sparse.csrsv_solve(trans, m, alpha, descr, csrVal, csrRowPtr, csrColInd, info, x, y)

Performs the analysis phase of the solution of the sparse triangular linear system op(A) * y = alpha * x.

Parameters:
  • trans – operation to apply to the matrix
  • m – number of rows of the matrix
  • alpha – scalar
  • descr – matrix descriptor
  • csrVal – vector of nonzero values of the matrix
  • csrRowPtr – vector of row pointers of the matrix
  • csrColInd – vector of column indices of the matrix
  • info – the analysis result from csrsv_analysis()
  • x – dense vector
  • y – dense vector into which the solve result is stored
Returns:

None

BLAS Level 3

pyculib.sparse.Sparse.csrmm(transA, m, n, k, nnz, alpha, descrA, csrValA, csrRowPtrA, csrColIndA, B, ldb, beta, C, ldc)

Matrix-matrix multiplication C = alpha * op(A) * B + beta * C where A is a sparse matrix in CSR format and B and C are dense matrices.

Parameters:
  • transA – operation to apply to A
  • m – number of rows of A
  • n – number of columns of B and C
  • k – number of columns of A
  • nnz – number of nonzeroes in A
  • alpha – scalar
  • descrA – matrix descriptor
  • csrValA – vector of nonzero values of A
  • csrRowPtrA – vector of row pointers of A
  • csrColIndA – vector of column indices of A
  • B – dense matrix
  • ldb – leading dimension of B
  • beta – scalar
  • C – dense matrix
  • ldc – leading dimension of C
Returns:

None

pyculib.sparse.Sparse.csrmm2(transA, transB, m, n, k, nnz, alpha, descrA, csrValA, csrRowPtrA, csrColIndA, B, ldb, beta, C, ldc)

Matrix-matrix multiplication C = alpha * op(A) * op(B) + beta * C where A is a sparse matrix in CSR format and B and C are dense matrices.

Parameters:
  • transA – operation to apply to A
  • transB – operation to apply to B
  • m – number of rows of A
  • n – number of columns of B and C
  • k – number of columns of A
  • nnz – number of nonzeroes in A
  • alpha – scalar
  • descrA – matrix descriptor
  • csrValA – vector of nonzero values of A
  • csrRowPtrA – vector of row pointers of A
  • csrColIndA – vector of column indices of A
  • B – dense matrix
  • ldb – leading dimension of B
  • beta – scalar
  • C – dense matrix
  • ldc – leading dimension of C
Returns:

None

pyculib.sparse.Sparse.csrsm_analysis(transA, m, nnz, descrA, csrValA, csrRowPtrA, csrColIndA)

Performs the analysis phase of the solution of a sparse triangular linear system op(A) * Y = alpha * X with multiple right-hand sides where A is a sparse matrix in CSR format, and X and Y are dense matrices.

Parameters:
  • transA – operation to apply to A
  • m – number of rows of A
  • nnz – number of nonzeroes in A
  • descrA – matrix descriptor
  • csrValA – vector of nonzero values of A
  • csrRowPtrA – vector of row pointers of A
  • csrColIndA – vector of column indices of A
Returns:

the analysis result

pyculib.sparse.Sparse.csrsm_solve(transA, m, n, alpha, descrA, csrValA, csrRowPtrA, csrColIndA, info, X, ldx, Y, ldy)

Performs the analysis phase of the solution of a sparse triangular linear system op(A) * Y = alpha * X with multiple right-hand sides where A is a sparse matrix in CSR format, and X and Y are dense matrices.

Parameters:
  • transA – operation to apply to A
  • m – number of rows of A
  • n – number of columns of B and C
  • alpha – scalar
  • descrA – matrix descriptor
  • csrValA – vector of nonzero values of A
  • csrRowPtrA – vector of row pointers of A
  • csrColIndA – vector of column indices of A
  • info – the analysis result from csrsm_analysis()
  • X – dense matrix
  • ldx – leading dimension of X
  • Y – dense matrix
  • ldy – leading dimension of Y
Returns:

None

Extra Functions

pyculib.sparse.Sparse.XcsrgeamNnz(m, n, descrA, nnzA, csrRowPtrA, csrColIndA, descrB, nnzB, csrRowPtrB, csrColIndB, descrC, csrRowPtrC)

Set up the sparsity pattern for the matrix operation C = alpha * A + beta * B where A, B, and C are all sparse matrices in CSR format.

Parameters:
  • m – number of rows of all matrices
  • n – number of columns of all matrices
  • descrA – matrix descriptor for A
  • nnzA – number of nonzeroes in A
  • csrRowPtrA – vector of row pointers of A
  • csrColIndA – vector of column indices of A
  • descrB – matrix descriptor for B
  • nnzB – number of nonzeroes in B
  • csrRowPtrB – vector of row pointers of B
  • csrColIndB – vector of column indices of B
  • descrC – matrix descriptor for B
  • csrRowPtrC – vector of row pointers of C, written to by this method
Returns:

number of nonzeroes in C

pyculib.sparse.Sparse.csrgeam(m, n, alpha, descrA, nnzA, csrValA, csrRowPtrA, csrColIndA, beta, descrB, nnzB, csrValB, csrRowPtrB, csrColIndB, descrC, csrValC, csrRowPtrC, csrColIndC)

Performs the the matrix operation C = alpha * A + beta * B where A, B, and C are all sparse matrices in CSR format.

Parameters:
  • m – number of rows of all matrices
  • n – number of columns of all matrices
  • alpha – scalar
  • descrA – matrix descriptor for A
  • nnzA – number of nonzeroes in A
  • csrValA – vector of nonzero values of A
  • csrRowPtrA – vector of row pointers of A
  • csrColIndA – vector of column indices of A
  • beta – scalar
  • descrB – matrix descriptor for B
  • nnzB – number of nonzeroes in B
  • csrValB – vector of nonzero values of B
  • csrRowPtrB – vector of row pointers of B
  • csrColIndB – vector of column indices of B
  • descrC – matrix descriptor for B
  • csrValC – vector of nonzero values of C
  • csrRowPtrC – vector of row pointers of C
  • csrColIndC – vector of column indices of C
Returns:

None

pyculib.sparse.Sparse.XcsrgemmNnz(transA, transB, m, n, k, descrA, nnzA, csrRowPtrA, csrColIndA, descrB, nnzB, csrRowPtrB, csrColIndB, descrC, csrRowPtrC)

Set up the sparsity pattern for the matrix operation C = op(A) * op(B) where A, B, and C are all sparse matrices in CSR format.

Parameters:
  • transA – operation to apply to A
  • transB – operation to apply to B
  • m – number of rows of A and C
  • n – number of columns of B and C
  • k – number of columns/rows of A/B
  • descrA – matrix descriptor for A
  • nnzA – number of nonzeroes in A
  • csrRowPtrA – vector of row pointers of A
  • csrColIndA – vector of column indices of A
  • descrB – matrix descriptor for B
  • nnzB – number of nonzeroes in B
  • csrRowPtrB – vector of row pointers of B
  • csrColIndB – vector of column indices of B
  • descrC – matrix descriptor for C
  • csrRowPtrC – vector of row pointers of C, written by this function
Returns:

number of nonzeroes in C

pyculib.sparse.Sparse.csrgemm(transA, transB, m, n, k, descrA, nnzA, csrValA, csrRowPtrA, csrColIndA, descrB, nnzB, csrValB, csrRowPtrB, csrColIndB, descrC, csrValC, csrRowPtrC, csrColIndC)

Perform the matrix operation C = op(A) * op(B) where A, B, and C are all sparse matrices in CSR format.

Parameters:
  • transA – operation to apply to A
  • transB – operation to apply to B
  • m – number of rows of A and C
  • n – number of columns of B and C
  • k – number of columns/rows of A/B
  • descrA – matrix descriptor for A
  • nnzA – number of nonzeroes in A
  • csrValA – vector of nonzero values in A
  • csrRowPtrA – vector of row pointers of A
  • csrColIndA – vector of column indices of A
  • descrB – matrix descriptor for B
  • nnzB – number of nonzeroes in B
  • csrValB – vector of nonzero values in B
  • csrRowPtrB – vector of row pointers of B
  • csrColIndB – vector of column indices of B
  • descrC – matrix descriptor for C
  • csrValC – vector of nonzero values in C
  • csrRowPtrC – vector of row pointers of C
  • csrColIndC – vector of column indices of C
Returns:

None

pyculib.sparse.Sparse.csrgemm_ez(A, B, transA='N', transB='N', descrA=None, descrB=None, descrC=None)
Performs the matrix operation C = op(A) * op(B) where A, B and C are all sparse matrices in CSR format. This function accepts and returns pyculib.sparse.CudaCSRMatrix matrices, and makes calls to XcsrgemmNnz() and csrgemm().
Parameters:
  • Apyculib.sparse.CudaCSRMatrix
  • Bpyculib.sparse.CudaCSRMatrix
  • transA – optional, operation to apply to A
  • transB – optional, operation to apply to B
  • descrA – optional, matrix descriptor for A
  • descrB – optional, matrix descriptor for B
  • descrC – optional, matrix descriptor for C
Returns:

pyculib.sparse.CudaCSRMatrix

Preconditioners

pyculib.sparse.Sparse.csric0(trans, m, descr, csrValA, csrRowPtrA, csrColIndA, info)

Computes incomplete Cholesky factorization of a sparse matrix in CSR format with 0 fill-in and no pivoting: op(A) = R**T * R. This method must follow a call to csrsv_analysis(). The matrix A is overwritten with the upper or lower triangular factors R or R**T.

Parameters:
  • trans – operation to apply to the matrix
  • m – number of rows and columns of the matrix
  • descr – matrix descriptor
  • csrValA – vector of nonzero values in A
  • csrRowPtrA – vector of row pointers of A
  • csrColIndA – vector of column indices of A
  • info – analysis result
Returns:

None

pyculib.sparse.Sparse.csrilu0(trans, m, descr, csrValA, csrRowPtrA, csrColIndA, info)

Computes incomplete-LU factorization of a sparse matrix in CSR format with 0 fill-in and no pivoting: op(A) = L * U. This method must follow a call to csrsv_analysis(). The matrix A is overwritten with the lower and upper triangular factors L and U.

Parameters:
  • trans – operation to apply to the matrix
  • m – number of rows and columns of the matrix
  • descr – matrix descriptor
  • csrValA – vector of nonzero values in A
  • csrRowPtrA – vector of row pointers of A
  • csrColIndA – vector of column indices of A
  • info – analysis result
Returns:

None

pyculib.sparse.Sparse.gtsv(m, n, dl, d, du, B, ldb)

Computes the solution of a tridiagonal linear system with multiple right-hand sides: A * Y = alpha * X.

Parameters:
  • m – the size of the linear system
  • n – the number of right-hand sides in the system
  • dl – dense vector storing the lower-diagonal elements
  • d – dense vector storing the diagonal elements
  • du – dense vector storing the upper-diagonal elements
  • B – dense matrix holding the right-hand sides of the system
  • ldb – the leading dimension of B
Returns:

None

pyculib.sparse.Sparse.gtsv_nopivot(m, n, dl, d, du, B, ldb)

Similar to gtsv(), but computes the solution without performing any pivoting.

Parameters:
  • m – the size of the linear system
  • n – the number of right-hand sides in the system
  • dl – dense vector storing the lower-diagonal elements
  • d – dense vector storing the diagonal elements
  • du – dense vector storing the upper-diagonal elements
  • B – dense matrix holding the right-hand sides of the system
  • ldb – the leading dimension of B
Returns:

None

pyculib.sparse.Sparse.gtsvStridedBatch(m, dl, d, du, x, batchCount, batchStride)

Computes the solution of i tridiagonal linear systems: A(i) * y(i) = alpha * x(i).

Parameters:
  • m – the size of the linear systems
  • dl – stacked dense vector storing the lower-diagonal elements of each system
  • d – stacked dense vector storing the diagonal elements of each system
  • du – stacked dense vector storing the upper-diagonal elements of each system
  • x – dense matrix holding the right-hand sides of the systems
  • batchCount – number of systems to solve
  • batchStride – number of elements separating the vectors of each system
Returns:

None

Format Conversion

pyculib.sparse.Sparse.bsr2csr(dirA, mb, nb, descrA, bsrValA, bsrRowPtrA, bsrColIndA, blockDim, descrC, csrValC, csrRowPtrC, csrColIndC)

Convert the sparse matrix A in BSR format to CSR format, stored in C.

Parameters:
  • dirA – row (‘R’) or column (‘C’) orientation of block storage
  • mb – number of block rows of A
  • nb – number of block columns of A
  • descrA – matrix descriptor for A
  • bsrValA – vector of nonzero values of A
  • bsrRowPtrA – vector of block row pointers of A
  • bsrColIndA – vector of block column indices of A
  • blockDim – block dimension of A
  • descrC – matrix descriptor for C
  • csrValA – vector of nonzero values in C
  • csrRowPtrA – vector of row pointers of C
  • csrColIndA – vector of column indices of C
Returns:

None

pyculib.sparse.Sparse.Xcoo2csr(cooRowInd, nnz, m, csrRowPtr)

Converts an array containing uncompressed row indices corresponding to the COO format into into an array of compressed row pointers corresponding to the CSR format.

Parameters:
  • cooRowInd – integer array of uncompressed row indices
  • nnz – number of nonzeroes
  • m – number of matrix rows
  • csrRowPtr – vector of row pointers to be written to
Returns:

None

pyculib.sparse.Sparse.csc2dense(m, n, descrA, cscValA, cscRowIndA, cscColPtrA, A, lda)

Convert the sparse matrix A in CSC format into a dense matrix.

Parameters:
  • m – number of rows of A
  • n – number of columns of A
  • descrA – matrix descriptor for A
  • cscValA – values in the CSC representation of A
  • cscRowIndA – row indices in the CSC representation of A
  • cscColPtrA – column pointers in the CSC representation of A
  • A – dense matrix representation of A to be written by this function.
  • lda – leading dimension of A
Returns:

None

pyculib.sparse.Sparse.Xcsr2bsrNnz(dirA, m, n, descrA, csrRowPtrA, csrColIndA, blockDim, descrC, bsrRowPtrC)

Performs the analysis necessary for converting a matrix in CSR format into BSR format.

Parameters:
  • dirA – row (‘R’) or column (‘C’) orientation of block storage
  • m – number of rows of matrix
  • n – number of columns of matrix
  • descrA – matrix descriptor for input matrix A
  • csrRowPtrA – row pointers of matrix
  • csrColIndA – column indices of matrix
  • blockDim – block dimension of output matrix C
  • descrC – matrix descriptor for output matrix C
Returns:

number of nonzeroes of matrix

pyculib.sparse.Sparse.csr2bsr(dirA, m, n, descrA, csrValA, csrRowPtrA, csrColIndA, blockDim, descrC, bsrValC, bsrRowPtrC, bsrColIndC)

Performs conversion of a matrix from CSR format into BSR format.

Parameters:
  • dirA – row (‘R’) or column (‘C’) orientation of block storage
  • m – number of rows of matrix
  • n – number of columns of matrix
  • descrA – matrix descriptor for input matrix A
  • csrValA – nonzero values of matrix
  • csrRowPtrA – row pointers of matrix
  • csrColIndA – column indices of matrix
  • blockDim – block dimension of output matrix C
  • descrC – matrix descriptor for output matrix C
  • bsrValC – nonzero values of output matrix C
  • bsrRowPtrC – block row pointers of output matrix C
  • bsrColIndC – block column indices of output matrix C
Returns:

number of nonzeroes of matrix

pyculib.sparse.Sparse.Xcsr2coo(csrRowPtr, nnz, m, cooRowInd)

Converts an array of compressed row pointers corresponding to the CSR format into an array of uncompressed row indices corresponding to the COO format.

Parameters:
  • csrRowPtr – vector of row pointers
  • nnz – number of nonzeroes
  • m – number of rows of matrix
  • cooRowInd – vector of uncompressed row indices written by this function
Returns:

None

pyculib.sparse.Sparse.csr2csc(m, n, nnz, csrVal, csrRowPtr, csrColInd, cscVal, cscRowInd, cscColPtr, copyValues)

Converts a sparse matrix in CSR format into a sparse matrix in CSC format.

Parameters:
  • m – number of rows of matrix
  • n – number of columns of matrix
  • nnz – number of nonzeroes of the matrix
  • csrVal – values in the CSR representation
  • csrRowPtr – row indices in the CSR representation
  • csrColInd – column pointers in the CSR representation
  • cscVal – values in the CSC representation
  • cscRowInd – row indices in the CSC representation
  • cscColPtr – column pointers in the CSC representation
  • copyValues‘N’ or ‘S’ for symbolic or numeric copy of values
Returns:

None

pyculib.sparse.Sparse.csr2dense(m, n, descr, csrVal, csrRowPtr, csrColInd, A, lda)

Convert a sparse matrix in CSR format into dense format.

Parameters:
  • m – number of rows of matrix
  • n – number of columns of matrix
  • descr – matrix descriptor
  • csrVal – values in the CSR representation
  • csrRowPtr – row indices in the CSR representation
  • csrColInd – column pointers in the CSR representation
  • A – the dense representation, written to by this function
  • lda – leading dimension of the matrix
Returns:

None

pyculib.sparse.Sparse.dense2csc(m, n, descrA, A, lda, nnzPerCol, cscVal, cscRowInd, cscColPtr)

Convert a dense matrix into a sparse matrix in CSC format. The nnzPerCol parameter may be computed with a call to nnz().

Parameters:
  • m – number of rows of matrix
  • n – number of columns of matrix
  • descrA – matrix descriptor
  • A – the matrix in dense format
  • lda – leading dimension of the matrix
  • nnzPerCol – array containing the number of nonzero elements per column
  • cscVal – values in the CSC representation
  • cscRowInd – row indices in the CSC representation
  • cscColPtr – column pointers in the CSC representation
Returns:

None

pyculib.sparse.Sparse.dense2csr(m, n, descrA, A, lda, nnzPerRow, csrVal, csrRowPtr, csrColInd)

Convert a dense matrix into a sparse matrix in CSR format. The nnzPerRow parameter may be computed with a call to nnz().

Parameters:
  • m – number of rows of matrix
  • n – number of columns of matrix
  • descrA – matrix descriptor
  • A – the matrix in dense format
  • lda – leading dimension of the matrix
  • nnzPerRow – array containing the number of nonzero elements per row
  • csrVal – values in the CSR representation
  • csrRowPtr – row indices in the CSR representation
  • csrColInd – column pointers in the CSR representation
Returns:

None

pyculib.sparse.Sparse.nnz(dirA, m, n, descrA, A, lda, nnzPerRowCol)

Computes the number of nonzero elements per row or column of a dense matrix, and the total number of nonzero elements in the matrix.

Parameters:
  • dirA‘R’ for the number of nonzeroes per row, or ‘C’ for per column.
  • m – number of rows of matrix
  • n – number of columns of matrix
  • descrA – matrix descriptor
  • A – the matrix
  • lda – leading dimension of the matrix
  • nnzPerRowCol – array to contain the number of nonzeroes per row or column
Returns:

total number of nonzeroes in the matrix