This is a demo/instructions for solving equations in Laplacian and SDD matrices. The sections are:


In [1]:
using Laplacians

We first generate a SDD, Positive Definite, system, and solve it using a direct solver. This uses the amd ordering, and is very fast.


In [2]:
a = grid2(5)
la = lap(a)
la[1,1] = la[1,1] + 1
F = cholfact(la)


Out[2]:
Base.SparseMatrix.CHOLMOD.Factor{Float64}
type:          LLt
method: simplicial
maxnnz:        102
nnz:           102

We can now use F to solve systems in this matrix, la. It is a complex structure that encodes a cholesky factorization, but it is much more than that.


In [3]:
n = size(a)[1]
b = randn(n)
x = F \ b
norm(la*x-b)


Out[3]:
1.5873593020111658e-14

Let's poke around to see what F has inside it.


In [4]:
F


Out[4]:
Base.SparseMatrix.CHOLMOD.Factor{Float64}
type:          LLt
method: simplicial
maxnnz:        102
nnz:           102

In [5]:
super(typeof(F))


Out[5]:
Factorization{Float64}

In [6]:
isa(F,Factorization)


Out[6]:
true

In [7]:
fieldnames(F)


Out[7]:
1-element Array{Symbol,1}:
 :p

Testing the speed of that

Let's see how long this will take on biggish grids, and on random regular graphs (which will have a lot of fill)


In [8]:
n = 500;
a = grid2(n)
la = lap(a)
la[1,1] = la[1,1] + 1
@time F = cholfact(la);


  1.298336 seconds (205 allocations: 256.844 MB, 4.86% gc time)

For a 500-by-500 grid, it took 1.5 seconds. We will now see that to use the solver, it takes 0.4 seconds.


In [9]:
N = size(la)[1]
b = randn(N)
@time x = F \ b
norm(la*x-b)


  0.106618 seconds (19 allocations: 7.634 MB, 6.58% gc time)
Out[9]:
2.3477210680497606e-10

For a random regular graph, we hit 1.5 seconds at around 20k vertices.


In [10]:
N = 20000;
a = randRegular(N,3)
la = lap(a)
la[1,1] = la[1,1] + 1
@time F = cholfact(la);


  1.787637 seconds (59 allocations: 260.694 MB, 0.45% gc time)

The time required for the solve is then around 0.05 seconds.


In [11]:
b = randn(N)
@time x = F \ b
norm(la*x-b)


  0.034356 seconds (19 allocations: 648.680 KB)
Out[11]:
2.5930764480212705e-11

What about just using \ ?


In [12]:
n = 500;
a = grid2(n)
la = lap(a)
la[1,1] = la[1,1] + 1
N = size(la)[1]
b = randn(N)
@time x = la \ b
norm(la*x-b)


  1.429117 seconds (151.48 k allocations: 274.719 MB, 6.32% gc time)
Out[12]:
7.323784555975749e-10

Just using \ appears to be right, so it is probably using cholfact.

We solve Laplacian systems by solving a system in the induced submatrix. Here are the steps, which I will then put into a wrapper function. It works by solving in a submatrix, like this:


In [13]:
la = lap(grid2(500))
N = size(la)[1]
lasub = la[1:(N-1),1:(N-1)]
Fsub = cholfact(lasub);

In [14]:
b = randn(N)
b = b - mean(b)
bs = b[1:(N-1)]
xs = Fsub \ bs;
x = [xs;0]
x = x - mean(x)
norm(la*x-b)


Out[14]:
6.440573359592263e-11

The following wraps a solver for SDD systems into a solver for Laplacian systems. We really need to work on the types of solver, and actually for everything else inside.

Let's see this work.


In [15]:
la = lap(a);
f = lapWrapSolver(cholfact,la)
b = randn(size(a)[1]); b = b - mean(b);
norm(la*f(b) - b)


Out[15]:
3.3314581543536396e-9

We now make two more versions: one that just takes the solver, and one that takes b as well.


In [16]:
lapChol2 = lapWrapSolver(cholfact)


Out[16]:
(anonymous function)

In [17]:
f = lapChol(la)


Out[17]:
(anonymous function)

In [18]:
norm(la*f(b) - b)


Out[18]:
3.3314581543536396e-9

In [19]:
x = lapWrapSolver(cholfact,la,b)
norm(la*x - b)


Out[19]:
3.3314581543536396e-9

In [ ]:

I really like the fact that Julia lets me type the following. It still needs reasonable types.

Here are examples of how to solve systems using the Conjugate Gradient.


In [20]:
n = 50
a = randn(n,n); a = a * a';
b = randn(n)
x = cg(a,b,maxits=100)
norm(a*x - b)


Out[20]:
2.827848817086888e-6

In [21]:
bbig = convert(Array{BigFloat,1},b);
xbig = cg(a,bbig,maxits=100)
norm(a*xbig - bbig)


CG stopped after: 50 iterations with relative error 1.427647178875816150213595857463364575541693109394711014852702931527973366187885e-38.
Out[21]:
1.101247547814646235270580604305022613523442584566417362886403637881772618266604e-37

In [22]:
la = lap(grid2(200))
n = size(la)[1]
b = randn(n)
b = b - mean(b);

In [23]:
@time x = cg(la,b,maxits=1000)
norm(la*x-b)


  0.843329 seconds (51.00 k allocations: 201.458 MB, 4.29% gc time)
Out[23]:
0.00019683922955253898

In [24]:
a = mapweight(grid2(200),x->1/(rand(1)[1]));
la = lap(a)
n = size(la)[1]
b = randn(n)
b = b - mean(b);

In [25]:
@time x = cg(la,b,maxits=4000)
norm(la*x-b)


  4.633557 seconds (16.01 k allocations: 1.193 GB, 4.36% gc time)
Out[25]:
43.83572140320749

Now, let's try a diagonal preconditioner.


In [26]:
d = diag(la)
pre(x) = x ./ d
@time x = pcg(la,b,pre,tol=1e-1,maxits=10^5)
norm(la*x-b)


  1.465714 seconds (491.96 k allocations: 595.625 MB, 6.51% gc time)
Out[26]:
19.957461012625778

In [27]:
@time x = cg(la,b,tol=1e-1,maxits=10^5)
norm(la*x-b)


  6.533181 seconds (21.05 k allocations: 1.568 GB, 4.16% gc time)
Out[27]:
18.40084364386128

It is very different for a random regular graph of the same size


In [28]:
n = 1000000
la = lap(randRegular(n,3))
b = randn(n)
b = b - mean(b);

In [29]:
@time x = cg(la,b,maxits=100)
norm(la*x-b)


  2.986557 seconds (172 allocations: 328.068 MB, 4.27% gc time)
Out[29]:
0.0007762305877414395

The following is a test of our stretch computation code. I begin by creating a grid graph with random weights, using our stretch computation, and checking that it agrees with the trace computation.


In [30]:
a = grid2(3)
a = uniformWeight(a)
a = a + a';

In [31]:
mst = kruskal(a)


Out[31]:
9x9 sparse matrix with 16 Float64 entries:
	[2, 1]  =  0.298218
	[1, 2]  =  0.298218
	[3, 2]  =  1.16494
	[5, 2]  =  1.41178
	[2, 3]  =  1.16494
	[6, 3]  =  1.10545
	[5, 4]  =  1.76592
	[7, 4]  =  1.92661
	[2, 5]  =  1.41178
	[4, 5]  =  1.76592
	[8, 5]  =  1.59899
	[3, 6]  =  1.10545
	[9, 6]  =  1.31406
	[4, 7]  =  1.92661
	[5, 8]  =  1.59899
	[6, 9]  =  1.31406

The following computes a matrix with entries corresponding to the nonzeros of a. For each nonzero, it puts in the stretch. So, to find the total stretch, we should sum them all and then divide by 2.


In [32]:
st = compStretches(mst,a)


Out[32]:
9x9 sparse matrix with 24 Float64 entries:
	[2, 1]  =  1.0
	[4, 1]  =  0.115101
	[1, 2]  =  1.0
	[3, 2]  =  1.0
	[5, 2]  =  1.0
	[2, 3]  =  1.0
	[6, 3]  =  1.0
	[1, 4]  =  0.115101
	[5, 4]  =  1.0
	[7, 4]  =  1.0
	⋮
	[8, 5]  =  1.0
	[3, 6]  =  1.0
	[5, 6]  =  0.203503
	[9, 6]  =  1.0
	[4, 7]  =  1.0
	[8, 7]  =  1.3647
	[5, 8]  =  1.0
	[7, 8]  =  1.3647
	[9, 8]  =  2.72786
	[6, 9]  =  1.0
	[8, 9]  =  2.72786

In [33]:
sum(st)/2


Out[33]:
12.411163624965448

We now check that we got the right answer by using the algebraic formula.


In [34]:
trace( pinv( full(lap(mst))) * lap(a)  )


Out[34]:
12.41116362496545

Now, let's do a speed test on a randomly weighted grid of side 2000.


In [35]:
a = grid2(2000)
@time a = uniformWeight(a)
@time a = a + a';
@time mst = kruskal(a);
@time st = compStretches(mst,a);


  3.110296 seconds (31.98 M allocations: 1.609 GB, 13.45% gc time)
  1.070956 seconds (35 allocations: 823.610 MB, 20.38% gc time)
  9.772786 seconds (95 allocations: 1.103 GB, 7.46% gc time)
  4.505141 seconds (115 allocations: 1.500 GB, 15.35% gc time)

Right now there are a multitude of solvers available in the repository. All of them have both calls for SDD systems and Laplacian systems.


In [36]:
SDDSolvers


Out[36]:
5-element Array{Function,1}:
 Laplacians.augTreeSolver    
 Laplacians.KMPSDDSolver     
 Laplacians.hybridSDDSolver  
 Laplacians.samplingSDDSolver
 Laplacians.AMGSolver        

In [37]:
LapSolvers


Out[37]:
5-element Array{Function,1}:
 Laplacians.augTreeLapSolver 
 Laplacians.KMPLapSolver     
 Laplacians.hybridLapSolver  
 Laplacians.samplingLapSolver
 Laplacians.AMGLapSolver     

In [47]:
a = mapweight(grid2(100),x->1/(rand(1)[1]));
n = a.n
la = lap(a); 
dval = zeros(n); dval[1] = dval[n] = 1e-3;
sdd = la + spdiagm(dval);
b = randn(n); b = b - mean(b);

The SDD solvers take in the sdd matrix alongside the tolerance, maxits and maxtime parameters.


In [49]:
for solver in SDDSolvers
    println("Solver ", solver)
    @time f = solver(sdd, maxits=1000, maxtime=10, tol=1e-4, verbose=false)
    @time x = f(b);
    println("Relative norm: ", norm(sdd * x - b) / norm(b), "\n")
end


Solver Laplacians.augTreeSolver
  0.067958 seconds (235.22 k allocations: 35.807 MB, 8.41% gc time)
  0.060869 seconds (2.72 k allocations: 70.339 MB, 17.35% gc time)
Relative norm: 9.250884271915595e-5

Solver Laplacians.KMPSDDSolver
  0.063388 seconds (254.46 k allocations: 44.168 MB, 14.65% gc time)
  0.102330 seconds (8.82 k allocations: 133.219 MB, 18.42% gc time)
Relative norm: 9.355785387633888e-5

Solver Laplacians.hybridSDDSolver
  0.064385 seconds (264.54 k allocations: 41.390 MB, 12.33% gc time)
  0.379060 seconds (2.07 M allocations: 266.961 MB, 9.90% gc time)
Relative norm: 9.832304243953024e-5

Solver Laplacians.samplingSDDSolver
  0.135522 seconds (706.74 k allocations: 101.318 MB, 13.69% gc time)
  0.179568 seconds (3.14 M allocations: 144.395 MB, 13.88% gc time)
Relative norm: 9.95716064856691e-5

Solver Laplacians.AMGSolver
  0.030372 seconds (166 allocations: 860.125 KB)
  0.116332 seconds (1.55 k allocations: 4.400 MB)
Relative norm: 9.519986283623698e-5

The SDD solvers take in the sdd matrix alongside the tolerance, maxits and maxtime parameters.


In [50]:
for solver in LapSolvers
    println("Solver ", solver)
    @time f = solver(a, maxits=1000, maxtime=10, tol=1e-4, verbose=false)
    @time x = f(b);
    println("Relative norm: ", norm(sdd * x - b) / norm(b), "\n")
end


Solver Laplacians.augTreeLapSolver
  0.072029 seconds (235.32 k allocations: 38.405 MB, 12.87% gc time)
  0.088307 seconds (6.51 k allocations: 98.563 MB, 18.56% gc time)
Relative norm: 0.00015469840172506333

Solver Laplacians.KMPLapSolver
  0.054602 seconds (254.18 k allocations: 39.445 MB, 15.07% gc time)
  0.107957 seconds (8.38 k allocations: 126.549 MB, 21.16% gc time)
Relative norm: 0.0001558714652838137

Solver Laplacians.hybridLapSolver
  0.053575 seconds (267.23 k allocations: 36.942 MB, 6.97% gc time)
  0.407590 seconds (3.02 M allocations: 286.329 MB, 10.58% gc time)
Relative norm: 0.00015420023270231535

Solver Laplacians.samplingLapSolver
  0.119625 seconds (706.35 k allocations: 90.641 MB, 21.37% gc time)
  0.218231 seconds (3.52 M allocations: 161.625 MB, 14.20% gc time)
Relative norm: 0.00015470313398451233

Solver Laplacians.AMGLapSolver
  0.043749 seconds (235 allocations: 2.744 MB)
  0.179676 seconds (1.57 k allocations: 4.402 MB)
Relative norm: 0.00014173955766018617


In [ ]: