Skip to content

Commit d19ef07

Browse files
committed
Lanczos implementation
1 parent 76f5ff2 commit d19ef07

File tree

6 files changed

+287
-81
lines changed

6 files changed

+287
-81
lines changed

examples/09_Disorder_KPM.jl

Lines changed: 21 additions & 24 deletions
Original file line numberDiff line numberDiff line change
@@ -47,51 +47,46 @@ plot_intensities(res)
4747

4848
sys_inhom = to_inhomogeneous(repeat_periodically(sys, (10, 10, 1)))
4949

50-
# Use [`symmetry_equivalent_bonds`](@ref) to iterate over all nearest neighbor
51-
# bonds of the inhomogeneous system. Modify each AFM exchange with a noise term
52-
# that has variance of 1/3. The newly minimized energy configuration allows for
53-
# long wavelength modulations on top of the original 120° order.
50+
# Iterate over all [`symmetry_equivalent_bonds`](@ref) for nearest neighbor
51+
# sites in the inhomogeneous system. Set the AFM exchange to include a noise
52+
# term that has a standard deviation of 1/3. The newly minimized energy
53+
# configuration allows for long wavelength modulations on top of the original
54+
# 120° order.
5455

55-
for (site1, site2, offset) in symmetry_equivalent_bonds(sys_inhom, Bond(1,1,[1,0,0]))
56+
for (site1, site2, offset) in symmetry_equivalent_bonds(sys_inhom, Bond(1, 1, [1, 0, 0]))
5657
noise = randn()/3
5758
set_exchange_at!(sys_inhom, 1.0 + noise, site1, site2; offset)
5859
end
5960

6061
minimize_energy!(sys_inhom, maxiters=5_000)
6162
plot_spins(sys_inhom; color=[S[3] for S in sys_inhom.dipoles], ndims=2)
6263

63-
# Traditional spin wave theory calculations become impractical for large system
64-
# sizes. Significant acceleration is possible with the [kernel polynomial
65-
# method](https://arxiv.org/abs/2312.08349). Enable it by selecting
66-
# [`SpinWaveTheoryKPM`](@ref) in place of the traditional
67-
# [`SpinWaveTheory`](@ref). Using KPM, the cost of an [`intensities`](@ref)
68-
# calculation becomes linear in system size and scales inversely with the width
69-
# of the line broadening `kernel`. Error tolerance is controlled through the
70-
# dimensionless `tol` parameter. A relatively small value, `tol = 0.01`, helps
71-
# to resolve the large intensities near the ordering wavevector. The alternative
72-
# choice `tol = 0.1` would be twice faster, but would introduce significant
73-
# numerical artifacts.
64+
# Spin wave theory with exact diagonalization becomes impractical for large
65+
# system sizes. Significant acceleration is possible with an iterative Krylov
66+
# space solver. With [`SpinWaveTheoryKPM`](@ref), the cost of an
67+
# [`intensities`](@ref) calculation scales linearly in the system size and
68+
# inversely with the width of the line broadening `kernel`. Good choices for the
69+
# dimensionless error tolerance are `tol=0.05` (more speed) or `tol=0.01` (more
70+
# accuracy).
7471
#
7572
# Observe from the KPM calculation that disorder in the nearest-neighbor
7673
# exchange serves to broaden the discrete excitation bands into a continuum.
7774

78-
swt = SpinWaveTheoryKPM(sys_inhom; measure=ssf_perp(sys_inhom), tol=0.01)
75+
swt = SpinWaveTheoryKPM(sys_inhom; measure=ssf_perp(sys_inhom), tol=0.05)
7976
res = intensities(swt, path; energies, kernel)
8077
plot_intensities(res)
8178

8279
# Now apply a magnetic field of magnitude 7.5 (energy units) along the global
83-
# ``ẑ`` axis. This field fully polarizes the spins. Because gap opens, a larger
84-
# tolerance of `tol = 0.1` can be used to accelerate the KPM calculation without
85-
# sacrificing much accuracy. The resulting spin wave spectrum shows a sharp mode
86-
# at the Γ-point (zone center) that broadens into a continuum along the K and M
87-
# points (zone boundary).
80+
# ``ẑ`` axis. This field fully polarizes the spins and opens a gap. The spin
81+
# wave spectrum shows a sharp mode at the Γ-point (zone center) that broadens
82+
# into a continuum along the K and M points (zone boundary).
8883

8984
set_field!(sys_inhom, [0, 0, 7.5])
9085
randomize_spins!(sys_inhom)
9186
minimize_energy!(sys_inhom)
9287

9388
energies = range(0.0, 9.0, 150)
94-
swt = SpinWaveTheoryKPM(sys_inhom; measure=ssf_perp(sys_inhom), tol=0.1)
89+
swt = SpinWaveTheoryKPM(sys_inhom; measure=ssf_perp(sys_inhom), tol=0.05)
9590
res = intensities(swt, path; energies, kernel)
9691
plot_intensities(res)
9792

@@ -105,8 +100,10 @@ for site in eachsite(sys_inhom)
105100
noise = randn()/6
106101
sys_inhom.gs[site] = [1 0 0; 0 1 0; 0 0 1+noise]
107102
end
103+
randomize_spins!(sys_inhom)
104+
minimize_energy!(sys_inhom)
108105

109-
swt = SpinWaveTheoryKPM(sys_inhom; measure=ssf_perp(sys_inhom), tol=0.1)
106+
swt = SpinWaveTheoryKPM(sys_inhom; measure=ssf_perp(sys_inhom), tol=0.05)
110107
res = intensities(swt, path; energies, kernel)
111108
plot_intensities(res)
112109

src/KPM/Lanczos.jl

Lines changed: 105 additions & 28 deletions
Original file line numberDiff line numberDiff line change
@@ -1,13 +1,40 @@
1-
# Reference implementation, without consideration of speed.
1+
# Reference implementation of the generalized Lanczos method, without
2+
# consideration of speed. The input A can be any Hermitian matrix. The optional
3+
# argument S is an inner-product that must be both Hermitian and positive
4+
# definite. Intuitively, Lanczos finds a low-rank approximant A S ≈ Q T Q† S,
5+
# where T is tridiagonal and Q is orthogonal in the sense that Q† S Q = I. The
6+
# first column of Q is the initial vector v, which must be normalized. The
7+
# Lanczos method is useful in two ways. First, eigenvalues of T are
8+
# representative of extremal eigenvalues of A (in an appropriate S-measure
9+
# sense). Second, one obtains a very good approximation to the matrix-vector
10+
# product, f(A S) v ≈ Q f(T) e₁, valid for any function f [1].
11+
#
12+
# The generalization of Lanczos to an arbitrary measure S is implemented as
13+
# follows: Each matrix-vector product A v is replaced by A S v, and each vector
14+
# inner product w† v is replaced by w† S v. Similar generalizations of Lanczos
15+
# have been considered in [2] and [3].
16+
#
17+
# 1. N. Amsel, T. Chen, A. Greenbaum, C. Musco, C. Musco, Near-optimal
18+
# approximation of matrix functions by the Lanczos method, (2013)
19+
# https://arxiv.org/abs/2303.03358v2.
20+
# 2. H. A. van der Vorst, Math. Comp. 39, 559-561 (1982),
21+
# https://doi.org/10.1090/s0025-5718-1982-0669648-0
22+
# 3. M. Grüning, A. Marini, X. Gonze, Comput. Mater. Sci. 50, 2148-2156 (2011),
23+
# https://doi.org/10.1016/j.commatsci.2011.02.021.
224
function lanczos_ref(A, S, v; niters)
325
αs = Float64[]
426
βs = Float64[]
27+
vs = typeof(v)[]
28+
29+
c = sqrt(real(dot(v, S, v)))
30+
@assert c 1 "Initial vector not normalized"
31+
v /= c
532

6-
v /= sqrt(real(dot(v, S, v)))
733
w = A * S * v
834
α = real(dot(w, S, v))
935
w = w - α * v
1036
push!(αs, α)
37+
push!(vs, v)
1138

1239
for _ in 2:niters
1340
β = sqrt(real(dot(w, S, w)))
@@ -19,35 +46,54 @@ function lanczos_ref(A, S, v; niters)
1946
v = vp
2047
push!(αs, α)
2148
push!(βs, β)
49+
push!(vs, v)
2250
end
2351

24-
return SymTridiagonal(αs, βs)
52+
Q = reduce(hcat, vs)
53+
T = SymTridiagonal(αs, βs)
54+
return (; Q, T)
2555
end
2656

2757

28-
29-
# Use Lanczos iterations to find a truncated tridiagonal representation of A S,
30-
# up to similarity transformation. Here, A is any Hermitian matrix, while S is
31-
# both Hermitian and positive definite. Traditional Lanczos uses the identity
32-
# matrix, S = I. The extension to non-identity matrices S is as follows: Each
33-
# matrix-vector product A v becomes A S v, and each vector inner product w† v
34-
# becomes w† S v. The implementation below follows the notation of Wikipedia,
35-
# and is the most stable of the four variants considered by Paige [1]. Each
36-
# Lanczos iteration requires only one matrix-vector multiplication for A and S,
37-
# respectively.
38-
39-
# Similar generalizations of Lanczos have been considered in [2] and [3].
58+
# An optimized implementation of the generalized Lanczos method. See
59+
# `lanczos_ref` for explanation, and a reference implementation. This optimized
60+
# implementation follows "the most stable" of the four variants considered by
61+
# Paige [1], as listed on Wikipedia. Here, however, we allow for an arbitary
62+
# measure S. In practice, this means that each matrix-vector product A v is
63+
# replaced by A S v, and each inner product w† v is replaced by w† S v.
64+
#
65+
# In this implementation, each Lanczos iteration requires only one matrix-vector
66+
# multiplication involving A and S, respectively.
67+
#
68+
# Details:
69+
#
70+
# * The return value is a pair, (T, lhs† Q). The symbol `lhs` denotes "left-hand
71+
# side" column vectors, packed horizontally into a matrix. If the `lhs`
72+
# argument is ommitted, its default value will be a matrix of width 0.
73+
# * If the initial vector `v` is not normalized, then this normalization will be
74+
# performed automatically, and the scale `√v S v` will be absorbed into the
75+
# return value `lhs† Q`.
76+
# * The number of iterations will never exceed length(v). If this limit is hit
77+
# then, mathematically, the Krylov subspace should be complete, and the matrix
78+
# T should be exactly similar to the matrix A S. In practice, however,
79+
# numerical error at finite precision may be severe. Further testing is
80+
# required to determine whether the Lanczos method is appropriate in this
81+
# case, where the matrices have modest dimension. (Direct diagonalization may
82+
# be preferable.)
83+
# * After `min_iters` iterations, this procedure will estimate the spectral
84+
# bandwidth Δϵ between extreme eigenvalues of T. Then `Δϵ / resolution` will
85+
# be an upper bound on the total number of iterations (which includes the
86+
# initial `min_iters` iterations as well as subsequent ones).
4087
#
4188
# 1. C. C. Paige, IMA J. Appl. Math., 373-381 (1972),
42-
# https://doi.org/10.1093%2Fimamat%2F10.3.373.
43-
# 2. H. A. van der Vorst, Math. Comp. 39, 559-561 (1982),
44-
# https://doi.org/10.1090/s0025-5718-1982-0669648-0
45-
# 3. M. Grüning, A. Marini, X. Gonze, Comput. Mater. Sci. 50, 2148-2156 (2011),
46-
# https://doi.org/10.1016/j.commatsci.2011.02.021.
47-
function lanczos(mulA!, mulS!, v; niters)
89+
# https://doi.org/10.1093%2Fimamat%2F10.3.373.
90+
91+
function lanczos(mulA!, mulS!, v; min_iters, resolution=Inf, lhs=zeros(length(v), 0), verbose=false)
4892
αs = Float64[]
4993
βs = Float64[]
94+
lhs_adj_Q = Vector{ComplexF64}[]
5095

96+
v = copy(v)
5197
vp = zero(v)
5298
Sv = zero(v)
5399
Svp = zero(v)
@@ -56,6 +102,7 @@ function lanczos(mulA!, mulS!, v; niters)
56102

57103
mulS!(Sv, v)
58104
c = sqrt(real(dot(v, Sv)))
105+
@assert c 1 "Initial vector not normalized"
59106
@. v /= c
60107
@. Sv /= c
61108

@@ -64,10 +111,36 @@ function lanczos(mulA!, mulS!, v; niters)
64111
@. w = w - α * v
65112
mulS!(Sw, w)
66113
push!(αs, α)
67-
68-
for _ in 2:niters
69-
β = sqrt(real(dot(w, Sw)))
70-
iszero(β) && break
114+
push!(lhs_adj_Q, lhs' * v)
115+
116+
# The maximum number of iterations is length(v)-1, because that is the
117+
# dimension of the full vector space. Break out early if niters has been
118+
# reached, which will be set according to the spectral bounds Δϵ after
119+
# iteration i == min_iters.
120+
niters = typemax(Int)
121+
for i in 1:length(v)-1
122+
if i == min_iters
123+
T = SymTridiagonal(αs, βs)
124+
Δϵ = eigmax(T) - eigmin(T)
125+
niters = max(min_iters, fld(Δϵ, resolution))
126+
niters += mod(niters, 2) # Round up to nearest even integer
127+
if verbose
128+
println("Δϵ=$Δϵ, niters=$niters")
129+
end
130+
end
131+
132+
i >= niters && break
133+
134+
# Let β = w† S w. If β < 0 then S is not a valid positive definite
135+
# measure. If β = 0, this would formally indicate that v is a linear
136+
# combination of only a subset of the eigenvectors of A S. In this case,
137+
# it is valid to break early for purposes of approximating the
138+
# matrix-vector product f(A S) v.
139+
β² = real(dot(w, Sw))
140+
iszero(β²) && break
141+
β² < 0 && error("S is not a positive definite measure")
142+
143+
β = sqrt(β²)
71144
@. vp = w / β
72145
@. Svp = Sw / β
73146
mulA!(w, Svp)
@@ -77,9 +150,10 @@ function lanczos(mulA!, mulS!, v; niters)
77150
@. v = vp
78151
push!(αs, α)
79152
push!(βs, β)
153+
push!(lhs_adj_Q, lhs' * v)
80154
end
81155

82-
return SymTridiagonal(αs, βs)
156+
return SymTridiagonal(αs, βs), reduce(hcat, lhs_adj_Q)
83157
end
84158

85159

@@ -107,6 +181,9 @@ function eigbounds(swt, q_reshaped, niters)
107181
end
108182

109183
v = randn(ComplexF64, 2L)
110-
tridiag = lanczos(mulA!, mulS!, v; niters)
111-
return eigmin(tridiag), eigmax(tridiag)
184+
Sv = zero(v)
185+
mulS!(Sv, v)
186+
v /= sqrt(v' * Sv)
187+
T, _ = lanczos(mulA!, mulS!, v; min_iters=niters)
188+
return eigmin(T), eigmax(T)
112189
end

0 commit comments

Comments
 (0)