Sparse Linear Algebra

AppleAccelerate wraps Apple's Sparse Solvers library for sparse matrix operations and direct solvers.

AASparseMatrix

A wrapper around Apple's SparseMatrix format, constructed from Julia's SparseMatrixCSC.

using AppleAccelerate, SparseArrays
import AppleAccelerate: AASparseMatrix, muladd!

A_jl = sprandn(100, 100, 0.05)
A = AASparseMatrix(A_jl)

x = randn(100)
y = A * x  # Sparse matrix-vector multiply

The constructor automatically detects symmetric and triangular structure and sets the appropriate Apple Accelerate attributes.

Matrix operations

FunctionDescription
AASparseMatrixConstruct from Julia sparse matrix
A * xSparse matrix-vector or matrix-matrix multiply
alpha * A * xScaled sparse multiply
muladd!Multiply-add: y += A * x or y += alpha * A * x
transpose(A)Transpose (sets flag, no copy)

Query functions

FunctionDescription
size(A)Matrix dimensions
eltype(A)Element type
issymmetric(A)Check if symmetric
istriu(A)Check if upper triangular
istril(A)Check if lower triangular
A[i, j]Element access

AAFactorization

Wraps Apple's SparseOpaqueFactorization. Lazy factorization wrapper: the factorization is computed on the first call to solve or by explicitly calling factor!.

import AppleAccelerate: AAFactorization, solve, solve!, factor!

A = sprandn(100, 100, 0.1) + 20I
f = AAFactorization(A)

# Factorization computed lazily on first solve
b = randn(100)
x = solve(f, b)

# Or explicitly
factor!(f)
x = solve(f, b)

Factorization types

TypeUse case
SparseFactorizationQRDefault for non-symmetric matrices
SparseFactorizationCholeskyDefault for symmetric positive definite
SparseFactorizationLDLTSymmetric indefinite (default LDLT)
SparseFactorizationLDLTUnpivotedSymmetric indefinite, no pivoting
SparseFactorizationLDLTSBKSymmetric indefinite, Bunch-Kaufman
SparseFactorizationLDLTTPPSymmetric indefinite, threshold partial pivoting
SparseFactorizationCholeskyAtACholesky of A'A (for least squares)

Solve functions

FunctionDescription
AAFactorizationLazy factorization wrapper
solve(f, b)Solve Ax = b, returns new vector/matrix
solve!(f, xb)Solve in-place (xb is overwritten with solution)
f \ bEquivalent to solve(f, b)
ldiv!(f, xb)Equivalent to solve!(f, xb)
ldiv!(x, f, b)Solve Ax = b, store result in x
factor!(f)Explicitly compute the factorization
factor!(f, type)Compute factorization with specific type
factorize(A::AASparseMatrix)Create an AAFactorization from a sparse matrix
AppleAccelerate.AASparseMatrixType

Matrix wrapper, containing the Apple sparse matrix struct and the pointed-to data. Construct from a SparseMatrixCSC.

Multiplication (*) and multiply-add (muladd!) with both Vector and Matrix objects are working. transpose creates a new matrix structure with the opposite transpose flag, that references the same CSC data.

source
AppleAccelerate.AAFactorizationType

Factorization object.

Create via f = AAFactorization(A::SparseMatrixCSC{T, Int64}). Calls to solve, ldiv, and their in-place versions require explicitly passing in the factorization object as the first argument. On construction, the struct stores a placeholder yet-to-be-factored object: the factorization is computed upon the first call to solve, or by explicitly calling factor!. If the matrix is symmetric, it defaults to a Cholesky factorization; otherwise, it defaults to QR.

source
AppleAccelerate.muladd!Function

Computes y += A*x in place. Note that this modifies its LAST argument.

source

Computes y += alphaAx in place. Note that this modifies its LAST argument.

source
AppleAccelerate.factor!Function
factor!(f::AAFactorization, [type::SparseFactorization_t])

Explicitly compute the factorization. If type is not specified, Cholesky is used for symmetric matrices and QR for non-symmetric. Called automatically by solve if the factorization has not yet been computed.

source
AppleAccelerate.solveFunction
solve(f::AAFactorization, b::StridedVecOrMat)

Solve the linear system Ax = b using Apple's Sparse Solvers, returning the solution x. The factorization is computed lazily on the first call if not already factored. Equivalent to f \ b.

source
AppleAccelerate.solve!Function
solve!(f::AAFactorization, xb::StridedVecOrMat)

Solve the linear system Ax = b in-place, overwriting xb with the solution. On input xb contains the right-hand side b; on output it contains the solution x. Equivalent to ldiv!(f, xb).

source