Sparse spectral clustering (SSC)

This example illustrates sparse spectral clustering using POGM applied to simulated data and (todo) hand-written digits.

Original version by Javier Salazar Cavazos.

This page comes from a single Julia file: ssc.jl.

You can access the source code for such Julia documentation using the 'Edit on GitHub' link in the top right. You can view the corresponding notebook in nbviewer here: ssc.ipynb, or open it in binder here: ssc.ipynb.

Setup

Add the Julia packages used in this demo. Change false to true in the following code block if you are using any of the following packages for the first time.

if false
    import Pkg
    Pkg.add([
        "Clustering"
        "InteractiveUtils"
        "LaTeXStrings"
        "LinearAlgebra"
        "MIRT"
        "MIRTjim"
        "MLDatasets"
        "Plots"
        "Random"
    ])
end

Tell Julia to use the following packages. Run Pkg.add() in the preceding code block first, if needed.

using Clustering: kmeans
using InteractiveUtils: versioninfo
using LaTeXStrings
using LinearAlgebra: Diagonal, eigen, I, opnorm
using MIRT: pogm_restart
using MIRTjim: jim, prompt
using Plots: default, gui, palette, plot, plot!, scatter, scatter!
using Random: randperm, seed!
default(); default(markersize=5, markerstrokecolor=:auto, label="")

The following line is helpful when running this file as a script; this way it will prompt user to hit a key after each figure is displayed.

isinteractive() ? jim(:prompt, true) : prompt(:draw);

Synthetic data

Generate synthetic data points in ℝ² that lie along $K = 2$ subspaces in the span of $(1,1)$ and $(1,-1)$.

seed!(3) # fix random generation for better debugging

xval = -9.5:1:9.5 # x locations before adding noise

N = length(xval)
σ = 0.5
x1 = xval + σ * randn(N) # noisy data that lies on union of subspaces
x2 = xval + σ * randn(N)
y1 = 1 * x1 .+ σ * randn(N) # y=x and y=-x are the 2 subspaces
y2 = -1 * x2 .+ σ * randn(N)

data = [ [x1';y1'] [x2';y2'] ] # gathering data to one matrix
clusters = [1*ones(Int,20); 2*ones(Int,20)]; # ground-truth clusters

if true # permute data points
    permuteOrder = randperm(40)
    data = data[:,permuteOrder]
    clusters = clusters[permuteOrder]
end
reord = invperm(permuteOrder)
40-element Vector{Int64}:
 24
 33
 27
 16
 36
 13
 20
 29
  6
  8
  ⋮
 14
 21
  5
 35
 34
 17
 30
 25
 40

plot subspaces and data points

p0 = plot(aspect_ratio = 1, size = (550, 500), xlabel=L"x_1", ylabel=L"x_2")
plot!(p0, xval, 1 .* xval)
plot!(p0, xval, -1 .* xval)
pd = deepcopy(p0)
scatter!(pd, x1, y1, color=1)
scatter!(pd, x2, y2, color=2)
plot!(pd, title = "Data and Subspaces")
Example block output

POGM for SSC

Solve the SSC problem with the self-representation cost function

\[\arg\min_C (1/2) ‖ Y (M ⊙ C) - Y ‖_{\mathrm{F}}² + λ ‖ C ‖_{1,1}\]

where $M$ is a mask matrix that is unity everywhere except 0 along the diagonal that forces each column of $Y$ to be represented as a (sparse) linear combination of other columns of $Y$. The regularizer encourages sparsity of $C$.

POGM is an optimal accelerated optimization method for convex composite cost functions.

λ = 0.001 # regularization parameter for sparsity term
Lf = opnorm(data,2)^2 # Lipschitz constant for ∇f: smooth term of obj function
npoint = size(data,2) # total # of points
x0 = zeros(npoint, npoint) # initialize solution
M = 1 .- I(npoint); # mask to force diag(C)=0

grad(x) = M .* (data' * (data*x - data)) # ∇f: gradient of smooth term
soft(x,t) = sign(x) * max(abs(x) - t, 0) # soft threshold at t
g_prox(z,c) = soft.(z, c * λ) # proximal operator for c*b*|x|_1

niter = 1000
A, _ = pogm_restart(x0, x->0, grad, Lf ; g_prox, niter) # POGM method
jim(A[reord,reord], "A")
Example block output

Spectral clustering

Cluster via a spectral method; see:

W = transpose(abs.(A)) + abs.(A) # Weight matrix, force Hermitian
jim(W[reord,reord], "W")
Example block output
D = vec(sum(W, dims=2)) # degree matrix of graph
D = D .^ (-1/2)
D = Diagonal(D) # normalized symmetric Laplacian formula
L = I - D * W * D
jim(L[reord,reord], "L")
Example block output

For $K=2$ subspaces we pick the bottom $K$ eigenvectors (smallest λ)

K = 2
E = eigen(L) # eigen value decomposition, really only need vectors
eigenVectors = E.vectors[:, 1:K];

seriescolor = palette([:orange, :skyblue], 2)
p4 = scatter(eigenVectors[:,1], eigenVectors[:,2],
 title="Spectral Embedding Plot",
 marker_z = clusters;
 seriescolor,
 colorbar = nothing,
)
Example block output

Since there are $K$ subspaces, we look for $K$ clusters in rows of eigenvectors using kmeans

results = kmeans(eigenVectors', K)
assign = results.assignments; # store assignments

Plot truth (on the left) and SSC results (on the right)

p1 = deepcopy(p0)
scatter!(p1, data[1,:], data[2,:];
 aspect_ratio = 1, size = (550, 450),
 xlims = (-11,11), ylims = (-11,11),
 marker_z = clusters,
 seriescolor,
 colorbar = nothing,
 title = "Truth",
);

p2 = deepcopy(p0)
scatter!(p2, data[1,:], data[2,:];
 aspect_ratio = 1, size = (550, 450),
 xlims = (-11,11), ylims = (-11,11),
 marker_z = assign,
 seriescolor,
 colorbar = nothing,
 title = "SSC (POGM)",
)
p12 = plot(p1, p2, layout = (1, 2), size=(1100, 450))
Example block output

Reproducibility

This page was generated with the following version of Julia:

using InteractiveUtils: versioninfo
io = IOBuffer(); versioninfo(io); split(String(take!(io)), '\n')
11-element Vector{SubString{String}}:
 "Julia Version 1.11.6"
 "Commit 9615af0f269 (2025-07-09 12:58 UTC)"
 "Build Info:"
 "  Official https://julialang.org/ release"
 "Platform Info:"
 "  OS: Linux (x86_64-linux-gnu)"
 "  CPU: 4 × AMD EPYC 7763 64-Core Processor"
 "  WORD_SIZE: 64"
 "  LLVM: libLLVM-16.0.6 (ORCJIT, znver3)"
 "Threads: 1 default, 0 interactive, 1 GC (on 4 virtual cores)"
 ""

And with the following package versions

import Pkg; Pkg.status()
Status `~/work/book-la-demo/book-la-demo/docs/Project.toml`
  [6e4b80f9] BenchmarkTools v1.6.0
  [aaaa29a8] Clustering v0.15.8
  [35d6a980] ColorSchemes v3.30.0
 [3da002f7] ColorTypes v0.11.5
 [c3611d14] ColorVectorSpace v0.10.0
  [717857b8] DSP v0.8.4
  [72c85766] Demos v0.1.0 `~/work/book-la-demo/book-la-demo`
  [e30172f5] Documenter v1.14.1
  [4f61f5a4] FFTViews v0.3.2
  [7a1cc6ca] FFTW v1.9.0
  [587475ba] Flux v0.16.4
  [a09fc81d] ImageCore v0.10.5
  [71a99df6] ImagePhantoms v0.8.1
  [b964fa9f] LaTeXStrings v1.4.0
  [7031d0ef] LazyGrids v1.1.0
  [599c1a8e] LinearMapsAA v0.12.0
  [98b081ad] Literate v2.20.1
  [7035ae7a] MIRT v0.18.2
  [170b2178] MIRTjim v0.25.0
  [eb30cadb] MLDatasets v0.7.18
  [efe261a4] NFFT v0.13.7
  [6ef6ca0d] NMF v1.0.3
  [15e1cf62] NPZ v0.4.3
  [0b1bfda6] OneHotArrays v0.2.10
  [429524aa] Optim v1.13.2
  [91a5bcdd] Plots v1.40.17
  [f27b6e38] Polynomials v4.1.0
  [2913bbd2] StatsBase v0.34.5
  [d6d074c3] VideoIO v1.1.1
  [b77e0a4c] InteractiveUtils v1.11.0
  [37e2e46d] LinearAlgebra v1.11.0
  [44cfe95a] Pkg v1.11.0
  [9a3f8284] Random v1.11.0
Info Packages marked with  and  have new versions available. Those with  may be upgradable, but those with  are restricted by compatibility constraints from upgrading. To see why use `status --outdated`

This page was generated using Literate.jl.