solvers

Linear Equation Solvers

For more, see the page on using solvers.

Laplacians.chol_lapFunction.
solver = chol_lap(A::AbstractArray)

Uses Cholesky Factorization to solve systems in Laplacians.

Laplacians.chol_sddmFunction.
solveSDDM = chol_sddm(sddm::AbstractMatrix; tol, maxits, maxtime, verbose, pcgIts=Int[])

This functions wraps cholfact so that it satsfies our interface. It ignores all the keyword arguments.

f = lapWrapSDDM(sddmSolver, A::AbstractArray; tol::Real=1e-6, maxits=Inf, maxtime=Inf, verbose=false, pcgIts=Int[], params...)
f = lapWrapSDDM(sddmSolver)

Uses a sddmSolver to solve systems of linear equations in Laplacian matrices.

Laplacians.cgFunction.
x = cg(mat, b; tol, maxits, maxtime, verbose, pcgIts)

solves a symmetric linear system mat x = b.

Arguments

  • tol is set to 1e-6 by default,
  • maxits defaults to Inf
  • maxtime defaults to Inf. It measures seconds.
  • verbose defaults to false
  • pcgIts is an array for returning the number of pcgIterations. Default is length 0, in which case nothing is returned.
x = cgLapSolver(A::AbstractMatrix; tol::Real=1e-6, maxits=Inf, maxtime=Inf, verbose=false, pcgIts=Int[])

Create a solver that uses cg to solve Laplacian systems in the laplacian of A. This just exists to satisfy our interface. It does nothing more than create the Laplacian and call cg on each connected component.

Laplacians.cgSolverFunction.
x = cgSolver(mat; tol, maxits, maxtime, verbose, pcgIts)

creates a solver for a PSD system mat. The parameters are as described in cg.

Laplacians.pcgFunction.
x = pcg(mat, b, pre; tol, maxits, maxtime, verbose, pcgIts, stag_test)`

solves a symmetric linear system using preconditioner pre.

Arguments

  • pre can be a function or a matrix. If a matrix, a function to solve it is created with cholFact.
  • tol is set to 1e-6 by default,
  • maxits defaults to Inf
  • maxtime defaults to Inf. It measures seconds.
  • verbose defaults to false
  • pcgIts is an array for returning the number of pcgIterations. Default is length 0, in which case nothing is returned.
  • stag_test=k stops the code if rho[it] > (1-1/k) rho[it-k]. Set to 0 to deactivate.
x = pcgLapSolver(A, B; tol::Real=1e-6, maxits=Inf, maxtime=Inf, verbose=false, pcgIts=Int[])

Create a solver that uses pcg to solve Laplacian systems in A Specialized for the case when the preconditioner the Laplacian matrix of B. It solves the preconditioner by Cholesky Factorization.

Laplacians.pcgSolverFunction.
x = pcgSolver(mat, pre; tol, maxits, maxtime, verbose, pcgIts)

creates a solver for a PSD system using preconditioner pre. The parameters are as described in pcg.

params = ApproxCholParams(order, output)

order can be one of

  • :deg (by degree, adaptive),
  • :wdeg (by original wted degree, nonadaptive),
  • :given
solver = approxchol_lap(a); x = solver(b);
solver = approxchol_lap(a; tol::Real=1e-6, maxits=1000, maxtime=Inf, verbose=false, pcgIts=Int[], params=ApproxCholParams())

A heuristic by Daniel Spielman inspired by the linear system solver in https://arxiv.org/abs/1605.02353 by Rasmus Kyng and Sushant Sachdeva. Whereas that paper eliminates vertices one at a time, this eliminates edges one at a time. It is probably possible to analyze it. The ApproxCholParams let you choose one of three orderings to perform the elimination.

  • ApproxCholParams(:given) - in the order given. This is the fastest for construction the preconditioner, but the slowest solve.
  • ApproxCholParams(:deg) - always eliminate the node of lowest degree. This is the slowest build, but the fastest solve.
  • ApproxCholParams(:wdeg) - go by a perturbed order of wted degree.

For more info, see http://danspielman.github.io/Laplacians.jl/latest/usingSolvers/index.html

solver = approxchol_sddm(sddm); x = solver(b);
solver = approxchol_sddm(sddm; tol=1e-6, maxits=1000, maxtime=Inf, verbose=false, pcgIts=Int[], params=ApproxCholParams())

Solves sddm systems by wrapping approxchol_lap. Not yet optimized directly for sddm.

For more info, see http://danspielman.github.io/Laplacians.jl/latest/usingSolvers/index.html

cn = condNumber(a, ldli; verbose=false)

Given an adjacency matrix a and an ldli computed by approxChol, this computes the condition number.

solver = augTreeLap(A; tol::Real=1e-6, maxits=Inf, maxtime=Inf, verbose=false, pcgIts=Int[], params=AugTreeParams())

An "augmented spanning tree" solver for Laplacians. It works by adding edges to a low stretch spanning tree. It calls augTreePrecon to form the preconditioner. params has entries

  • params.treeAlg default to akpw
  • params.opt if true, it interacts with cholmod to choose a good number of edges to add back. If false, it adds back 2*sqrt(n).
pre = augTreeLapPrecon{Tv,Ti}(A; params=AugTreeParams())

This is an augmented spanning tree preconditioner for Laplacians. It takes as optional input a tree growing algorithm. It adds back 2sqrt(n) edges via augmentTree: the sqrt(n) of highest stretch and another sqrt(n) sampled according to stretch. For most purposes, one should directly call augTreeLapSolver.

pre = augTreePrecon{Tv,Ti}(ddmat::SparseMatrixCSC{Tv,Ti}; params=AugTreeParams())

This is an augmented spanning tree preconditioner for diagonally dominant linear systems. It takes as optional input a tree growing algorithm. It adds back 2sqrt(n) edges via augmentTree: the sqrt(n) of highest stretch and another sqrt(n) sampled according to stretch. For most purposes, one should directly call augTreeSolver.

solver = augTreeSddm(sddm; tol::Real=1e-6, maxits=Inf, maxtime=Inf, verbose=false, pcgIts=Int[],  params=AugTreeParams())

An "augmented spanning tree" solver for positive definite diagonally dominant matrices. It works by adding edges to a low stretch spanning tree. It calls augTreePrecon to form the preconditioner. params has entries

  • params.treeAlg default to akpw
  • params.opt if true, it interacts with cholmod to choose a good number of edges to add back. If false, it adds back 2*sqrt(n).
B = augmentTree{Tv,Ti}(tree, A, k)

Takes as input a tree and an adjacency matrix of a graph. It then computes the stretch of every edge of the graph wrt the tree. It then adds back the k edges of highest stretch, and k edges sampled according to stretch.

This is the old alg. We now recommend using augmentTreeOpt.

Parameters for the KMP solver

lapSolver = KMPLapSolver(A; verbose, tol, maxits, maxtime, pcgIts, params::KMPParams)

Solves linear equations in the Laplacian of graph with adjacency matrix A.

Based on the paper "Approaching optimality for solving SDD systems" by Koutis, Miller, and Peng, <i>SIAM Journal on Computing</i>, 2014.

sddmSolver = KMPSDDMSolver(mat; verbose, tol, maxits, maxtime, pcgIts, params::KMPParams)

Solves linear equations in symmetric, diagonally dominant matrices with non-positive off-diagonals. Based on the paper "Approaching optimality for solving SDD systems" by Koutis, Miller, and Peng, <i>SIAM Journal on Computing</i>, 2014.

x = harmonic_interp(adj_mat, S, vals; tol=1e-6)

Interpolates a function on a graph, given by its adjacency matrix, by minizing the Laplacian quadratic form subject to the boundary conditions that x[S[i]] = vals[i] for i in S.

This is the algorithm sometimes known as Label Propagation, or Semi-Supervised Learning on Graphs. The idea comes from the paper "Semi-Supervised Learning Using Gaussian Fields and Harmonic Functions" by Zhu, Gharamani, and Lafferty from ICML 2003.

This version might fail for disconnected graphs. You can check if a graph is connected with isConnected(adj_mat).