/
EDCore.jl
382 lines (330 loc) · 16.4 KB
/
EDCore.jl
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
module EDCore
using Arpack: eigs
using LinearAlgebra: Eigen, Factorization, norm
using QuantumLattices: plain, bonds, expand, id, idtype, reparameter
using QuantumLattices: AbelianNumber, AbstractLattice, Algorithm, Boundary, Frontend, Hilbert, Image, LinearTransformation, MatrixRepresentation, Metric, Neighbors, Operator, OperatorGenerator, OperatorPack, Operators, OperatorSum, OperatorUnit, Table, Term, VectorSpace, VectorSpaceEnumerative, VectorSpaceStyle
using SparseArrays: SparseMatrixCSC, spzeros
using TimerOutputs: TimerOutput, @timeit
import LinearAlgebra: eigen
import QuantumLattices: Parameters, ⊕, add!, contentnames, dtype, getcontent, kind, matrix, parameternames, statistics, update!
export edtimer, ED, EDEigen, EDKind, EDMatrix, EDMatrixRepresentation, Sector, SectorFilter, TargetSpace, productable, sumable
"""
const edtimer = TimerOutput()
The default shared timer for all exact diagonalization methods.
"""
const edtimer = TimerOutput()
"""
abstract type Sector <: OperatorUnit
A sector of the Hilbert space which forms the bases of an irreducible representation of the Hamiltonian of a quantum lattice system.
"""
abstract type Sector <: OperatorUnit end
@inline Base.hash(sector::Sector, h::UInt) = hash(id(sector), h)
@inline Base.:(==)(sector₁::Sector, sector₂::Sector) = isequal(id(sector₁), id(sector₂))
@inline Base.isequal(sector₁::Sector, sector₂::Sector) = isequal(id(sector₁), id(sector₂))
"""
sumable(sector₁::Sector, sector₂::Sector) -> Bool
Judge whether two sectors could be direct summed.
"""
function sumable end
"""
productable(sector₁::Sector, sector₂::Sector) -> Bool
Judge whether two sectors could be direct producted.
"""
function productable end
"""
matrix(ops::Operators, braket::NTuple{2, Sector}, table, dtype=valtype(eltype(ops))) -> SparseMatrixCSC{dtype, Int}
Get the CSC-formed sparse matrix representation of a set of operators.
Here, `table` specifies the order of the operator ids.
"""
function matrix(ops::Operators, braket::NTuple{2, Sector}, table, dtype=valtype(eltype(ops)))
result = spzeros(dtype, length(braket[1]), length(braket[2]))
for op in ops
result += matrix(op, braket, table, dtype)
end
return result
end
"""
TargetSpace{S<:Sector} <: VectorSpace{S}
The target Hilbert space in which the exact diagonalization method is performed, which could be the direct sum of several sectors.
"""
struct TargetSpace{S<:Sector} <: VectorSpace{S}
sectors::Vector{S}
end
@inline VectorSpaceStyle(::Type{<:TargetSpace}) = VectorSpaceEnumerative()
@inline contentnames(::Type{<:TargetSpace}) = (:contents,)
@inline getcontent(target::TargetSpace, ::Val{:contents}) = target.sectors
function add!(target::TargetSpace, sector::Sector)
@assert all(map(previous->sumable(previous, sector), target.sectors)) "add! error: could not be direct summed."
push!(target.sectors, sector)
return target
end
function add!(target::TargetSpace, another::TargetSpace)
for sector in another
add!(target, sector)
end
return target
end
"""
TargetSpace(sector::Sector, sectors::Sector...)
Construct a target space from sectors.
"""
@inline function TargetSpace(sector::Sector, sectors::Sector...)
result = TargetSpace([sector])
for sec in sectors
add!(result, sec)
end
return result
end
"""
TargetSpace(hilbert::Hilbert)
TargetSpace(hilbert::Hilbert, quantumnumbers::Union{AbelianNumber, Tuple{AbelianNumber, Vararg{AbelianNumber}}})
Construct a target space from the total Hilbert space and the associated quantum numbers.
"""
@inline TargetSpace(hilbert::Hilbert) = TargetSpace(Sector(hilbert))
@inline TargetSpace(hilbert::Hilbert, quantumnumbers::Union{AbelianNumber, Tuple{AbelianNumber, Vararg{AbelianNumber}}}) = TargetSpace(map(quantumnumber->Sector(hilbert, quantumnumber), wrapper(quantumnumbers))...)
"""
⊕(sector::Sector, sectors::Union{Sector, TargetSpace}...) -> TargetSpace
⊕(target::TargetSpace, sectors::Union{Sector, TargetSpace}...) -> TargetSpace
Get the direct sum of sectors and target spaces.
"""
@inline function ⊕(sector::Sector, sectors::Union{Sector, TargetSpace}...)
result = TargetSpace(sector)
map(op->add!(result, op), sectors)
return result
end
@inline function ⊕(target::TargetSpace, sectors::Union{Sector, TargetSpace}...)
result = TargetSpace(copy(target.sectors))
map(op->add!(result, op), sectors)
return result
end
"""
EDMatrix{S<:Sector, M<:SparseMatrixCSC} <: OperatorPack{M, Tuple{S, S}}
Matrix representation of quantum operators between a ket and a bra Hilbert space.
"""
struct EDMatrix{S<:Sector, M<:SparseMatrixCSC} <: OperatorPack{M, Tuple{S, S}}
bra::S
ket::S
matrix::M
end
@inline parameternames(::Type{<:EDMatrix}) = (:sector, :value)
@inline getcontent(m::EDMatrix, ::Val{:id}) = (m.bra, m.ket)
@inline getcontent(m::EDMatrix, ::Val{:value}) = m.matrix
@inline dtype(::Type{<:EDMatrix{<:Sector, M}}) where {M<:SparseMatrixCSC} = eltype(M)
@inline EDMatrix(m::SparseMatrixCSC, braket::NTuple{2, Sector}) = EDMatrix(braket[1], braket[2], m)
@inline Base.promote_rule(M::Type{<:EDMatrix}, N::Type{<:Number}) = reparameter(M, :value, reparameter(valtype(M), 1, promote_type(dtype(M), N)))
"""
EDMatrix(sector::Sector, m::SparseMatrixCSC)
EDMatrix(braket::NTuple{2, Sector}, m::SparseMatrixCSC)
Construct a matrix representation when
1) the ket and bra spaces share the same bases;
2-3) the ket and bra spaces may be different.
"""
@inline EDMatrix(sector::Sector, m::SparseMatrixCSC) = EDMatrix(sector, sector, m)
@inline EDMatrix(braket::NTuple{2, Sector}, m::SparseMatrixCSC) = EDMatrix(braket[1], braket[2], m)
"""
EDEigen{V<:Number, T<:Number, S<:Sector} <: Factorization{T}
Eigen decomposition in exact diagonalization method.
Compared to the usual eigen decomposition `Eigen`, `EDEigen` contains a `:sectors` attribute to store the sectors of Hilbert space in which the eigen values and eigen vectors are computed.
Furthermore, given that in different sectors the dimensions of the sub-Hilbert spaces can also be different, the `:vectors` attribute of `EDEigen` is a vector of vector instead of a matrix.
"""
struct EDEigen{V<:Number, T<:Number, S<:Sector} <: Factorization{T}
values::Vector{V}
vectors::Vector{Vector{T}}
sectors::Vector{S}
end
@inline Base.iterate(content::EDEigen) = (content.values, Val(:vectors))
@inline Base.iterate(content::EDEigen, ::Val{:vectors}) = (content.vectors, Val(:sectors))
@inline Base.iterate(content::EDEigen, ::Val{:sectors}) = (content.sectors, Val(:done))
@inline Base.iterate(content::EDEigen, ::Val{:done}) = nothing
"""
eigen(m::EDMatrix; nev=6, which=:SR, tol=0.0, maxiter=300, sigma=nothing, v₀=dtype(m)[]) -> Eigen
Solve the eigen problem by the restarted Lanczos method provided by the Arpack package.
"""
@inline function eigen(m::EDMatrix; nev=6, which=:SR, tol=0.0, maxiter=300, sigma=nothing, v₀=dtype(m)[])
@assert m.bra==m.ket "eigen error: eigen decomposition of an `EDMatrix` are only available for those with the same bra and ket spaces."
if size(m.matrix)[1] > 1
eigvals, eigvecs = eigs(m.matrix; nev=nev, which=which, tol=tol, maxiter=maxiter, sigma=sigma, ritzvec=true, v0=v₀)
@assert norm(imag(eigvals))<10^-14 "eigen error: non-vanishing imaginary parts of the eigen values."
eigvals = real(eigvals)
else
eigvals, eigvecs = eigen(collect(m.matrix))
end
return Eigen(eigvals, eigvecs)
end
"""
eigen(ms::OperatorSum{<:EDMatrix}; nev::Int=1, tol::Real=0.0, maxiter::Int=300, v₀::Union{AbstractVector, Dict{<:Sector, <:AbstractVector}, Dict{<:AbelianNumber, <:AbstractVector}}=dtype(eltype(ms))[], timer::TimerOutput=edtimer)
Solve the eigen problem by the restarted Lanczos method provided by the Arpack package.
"""
@inline function eigen(ms::OperatorSum{<:EDMatrix}; nev::Int=1, tol::Real=0.0, maxiter::Int=300, v₀::Union{AbstractVector, Dict{<:Sector, <:AbstractVector}, Dict{<:AbelianNumber, <:AbstractVector}}=dtype(eltype(ms))[], timer::TimerOutput=edtimer)
@timeit timer "eigen" begin
isa(v₀, AbstractVector) && (v₀ = Dict(m.ket=>v₀ for m in ms))
isa(v₀, Dict{<:AbelianNumber, <:AbstractVector}) && (v₀ = Dict(m.ket=>get(v₀, AbelianNumber(m.ket), dtype(eltype(ms))[]) for m in ms))
values, vectors, sectors = real(dtype(eltype(ms)))[], Vector{dtype(eltype(ms))}[], eltype(idtype(eltype(ms)))[]
for m in ms
@timeit timer string(AbelianNumber(m.ket)) begin
k = min(length(m.ket), nev)
eigensystem = eigen(m; nev=k, which=:SR, tol=tol, maxiter=maxiter, v₀=get(v₀, m.ket, dtype(eltype(ms))[]))
for i = 1:k
push!(values, eigensystem.values[i])
push!(vectors, eigensystem.vectors[:, i])
push!(sectors, m.ket)
end
end
end
perm = sortperm(values)[1:min(nev, length(values))]
result = EDEigen(values[perm], vectors[perm], sectors[perm])
end
return result
end
"""
EDMatrixRepresentation{D<:Number, S<:Sector, T} <: MatrixRepresentation
Exact matrix representation of a quantum lattice system on a target Hilbert space.
"""
struct EDMatrixRepresentation{D<:Number, S<:Sector, T} <: MatrixRepresentation
brakets::Vector{Tuple{S, S}}
table::T
EDMatrixRepresentation{D}(brakets::Vector{Tuple{S, S}}, table) where {D<:Number, S<:Sector} = new{D, S, typeof(table)}(brakets, table)
end
@inline function Base.valtype(::Type{<:EDMatrixRepresentation{D, S}}, M::Type{<:Operator}) where {D<:Number, S<:Sector}
@assert promote_type(D, valtype(M))==D "valtype error: convert $(valtype(M)) to $D is inexact."
E = EDMatrix{S, SparseMatrixCSC{D, Int}}
return OperatorSum{E, idtype(E)}
end
@inline Base.valtype(R::Type{<:EDMatrixRepresentation}, M::Type{<:Operators}) = valtype(R, eltype(M))
function (representation::EDMatrixRepresentation)(m::Operator; kwargs...)
result = zero(valtype(representation, m))
for braket in representation.brakets
add!(result, EDMatrix(braket, matrix(m, braket, representation.table, dtype(eltype(result)); kwargs...)))
end
return result
end
"""
EDMatrixRepresentation{D}(target::TargetSpace, table) where {D<:Number}
Construct a exact matrix representation.
"""
@inline EDMatrixRepresentation{D}(target::TargetSpace, table) where {D<:Number} = EDMatrixRepresentation{D}([(sector, sector) for sector in target], table)
"""
SectorFilter{S} <: LinearTransformation
Filter the target bra and ket Hilbert spaces.
"""
struct SectorFilter{S} <: LinearTransformation
brakets::Set{Tuple{S, S}}
end
@inline Base.valtype(::Type{<:SectorFilter}, M::Type{<:OperatorSum{<:EDMatrix}}) = M
@inline (sectorfileter::SectorFilter)(m::EDMatrix) = id(m)∈sectorfileter.brakets ? m : 0
@inline SectorFilter(sector::Sector, sectors::Sector...) = SectorFilter((sector, sector), map(op->(op, op), sectors)...)
@inline SectorFilter(braket::NTuple{2, Sector}, brakets::NTuple{2, Sector}...) = SectorFilter(push!(Set{typeof(braket)}(), braket, brakets...))
"""
EDKind{K}
The kind of the exact diagonalization method applied to a quantum lattice system.
"""
struct EDKind{K} end
@inline EDKind(K::Symbol) = EDKind{K}()
@inline EDKind(hilbert::Hilbert) = EDKind(typeof(hilbert))
@inline EDKind(::Type{H}) where {H<:Hilbert} = error("EDKind error: not defined.")
"""
ED{K<:EDKind, L<:AbstractLattice, G<:OperatorGenerator, M<:Image} <: Frontend
Exact diagonalization method of a quantum lattice system.
"""
struct ED{K<:EDKind, L<:AbstractLattice, G<:OperatorGenerator, M<:Image} <: Frontend
lattice::L
H::G
Hₘ::M
function ED{K}(lattice::AbstractLattice, H::OperatorGenerator, mr::EDMatrixRepresentation; timer::TimerOutput=edtimer) where {K<:EDKind}
@timeit timer "matrix" begin
@timeit timer "prepare" Hₘ = mr(H)
end
new{K, typeof(lattice), typeof(H), typeof(Hₘ)}(lattice, H, Hₘ)
end
end
@inline kind(ed::ED) = kind(typeof(ed))
@inline kind(::Type{<:ED{K}}) where K = K()
@inline Base.valtype(::Type{<:ED{<:EDKind, <:AbstractLattice, G}}) where {G<:OperatorGenerator} = valtype(eltype(G))
@inline statistics(ed::ED) = statistics(typeof(ed))
@inline statistics(::Type{<:ED{<:EDKind, <:AbstractLattice, G}}) where {G<:OperatorGenerator} = statistics(eltype(eltype(G)))
@inline function update!(ed::ED; kwargs...)
if length(kwargs)>0
update!(ed.H; kwargs...)
update!(ed.Hₘ, ed.H)
end
return ed
end
@inline Parameters(ed::ED) = Parameters(ed.H)
"""
matrix(ed::ED, sectors::Union{AbelianNumber, Sector}...; timer::TimerOutput=edtimer, kwargs...) -> OperatorSum{<:EDMatrix}
matrix(ed::Algorithm{<:ED}, sectors::Union{AbelianNumber, Sector}...; kwargs...) -> OperatorSum{<:EDMatrix}
Get the sparse matrix representation of a quantum lattice system in the target space.
"""
function matrix(ed::ED; timer::TimerOutput=edtimer, kwargs...)
@timeit timer "matrix" begin
@timeit timer "expand" (result = expand(ed.Hₘ))
end
return result
end
function matrix(ed::ED, sectors::Sector...; timer::TimerOutput=edtimer, kwargs...)
@timeit timer "matrix" begin
@timeit timer "expand" (result = expand(SectorFilter(sectors...)(ed.Hₘ)))
end
return result
end
function matrix(ed::ED, quantumnumbers::AbelianNumber...; timer::TimerOutput=edtimer, kwargs...)
sectors = [braket[1] for braket in ed.Hₘ.transformation.brakets if AbelianNumber(braket[1]) in quantumnumbers]
return matrix(ed, sectors...; timer=timer, kwargs...)
end
function matrix(ed::Algorithm{<:ED}, sectors::Union{AbelianNumber, Sector}...; kwargs...)
return matrix(ed.frontend, sectors...; timer=ed.timer, kwargs...)
end
"""
eigen(ed::ED, sectors::Union{AbelianNumber, Sector}...; timer::TimerOutput=edtimer, kwargs...) -> EDEigen
eigen(ed::Algorithm{<:ED}, sectors::Union{AbelianNumber, Sector}...; kwargs...) -> EDEigen
Solve the eigen problem by the restarted Lanczos method provided by the Arpack package.
"""
@inline eigen(ed::ED, sectors::Union{AbelianNumber, Sector}...; timer::TimerOutput=edtimer, kwargs...) = eigen(matrix(ed, sectors...; timer=timer); timer=timer, kwargs...)
@inline eigen(ed::Algorithm{<:ED}, sectors::Union{AbelianNumber, Sector}...; kwargs...) = eigen(matrix(ed, sectors...; timer=ed.timer); timer=ed.timer, kwargs...)
"""
ED(
lattice::AbstractLattice, hilbert::Hilbert, terms::Union{Term, Tuple{Term, Vararg{Term}}}, targetspace::TargetSpace=TargetSpace(hilbert), dtype::Type{<:Number}=commontype(terms), boundary::Boundary=plain;
neighbors::Union{Nothing, Int, Neighbors}=nothing, timer::TimerOutput=edtimer
)
Construct the exact diagonalization method for a quantum lattice system.
"""
function ED(
lattice::AbstractLattice, hilbert::Hilbert, terms::Union{Term, Tuple{Term, Vararg{Term}}}, targetspace::TargetSpace=TargetSpace(hilbert), dtype::Type{<:Number}=commontype(terms), boundary::Boundary=plain;
neighbors::Union{Nothing, Int, Neighbors}=nothing, timer::TimerOutput=edtimer
)
terms = wrapper(terms)
isnothing(neighbors) && (neighbors = maximum(term->term.bondkind, terms))
H = OperatorGenerator(terms, bonds(lattice, neighbors), hilbert, boundary; half=false)
mr = EDMatrixRepresentation{dtype}(targetspace, Table(hilbert, Metric(EDKind(hilbert), hilbert)))
return ED{typeof(EDKind(hilbert))}(lattice, H, mr; timer=timer)
end
@inline wrapper(x) = (x,)
@inline wrapper(x::Tuple) = x
@inline commontype(term::Term) = valtype(term)
@inline @generated commontype(terms::Tuple{Vararg{Term}}) = mapreduce(valtype, promote_type, fieldtypes(terms))
"""
ED(
lattice::AbstractLattice,
hilbert::Hilbert,
terms::Union{Term, Tuple{Term, Vararg{Term}}},
quantumnumbers::Union{AbelianNumber, Tuple{AbelianNumber, Vararg{AbelianNumber}}},
dtype::Type{<:Number}=commontype(terms),
boundary::Boundary=plain;
neighbors::Union{Nothing, Int, Neighbors}=nothing,
timer::TimerOutput=edtimer,
)
Construct the exact diagonalization method for a quantum lattice system.
"""
function ED(
lattice::AbstractLattice,
hilbert::Hilbert,
terms::Union{Term, Tuple{Term, Vararg{Term}}},
quantumnumbers::Union{AbelianNumber, Tuple{AbelianNumber, Vararg{AbelianNumber}}},
boundary::Boundary=plain,
dtype::Type{<:Number}=commontype(terms);
neighbors::Union{Nothing, Int, Neighbors}=nothing,
timer::TimerOutput=edtimer
)
return ED(lattice, hilbert, terms, TargetSpace(hilbert, wrapper(quantumnumbers)), dtype, boundary; neighbors=neighbors, timer=timer)
end
end # module