-
Notifications
You must be signed in to change notification settings - Fork 0
/
spinhalfsquare_small.jl
216 lines (190 loc) · 8.32 KB
/
spinhalfsquare_small.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
using SparseArrays
using LinearAlgebra
using ExactDiagonalization
using LatticeTools
using MinimalPerfectHash
@show Threads.nthreads()
## Set up lattice
n1, n2 = 3, 3
n_sites = n1 * n2;
unitcell = make_unitcell([1.0 0.0; 0.0 1.0]; SiteType=String)
addsite!(unitcell, "Spin", FractCoord([0, 0], [0.0, 0.0]))
lattice = make_lattice(unitcell, [n1 0; 0 n2])
tsym = TranslationSymmetry(lattice)
psym = project(PointSymmetryDatabase.get(13), [1 0 0; 0 1 0]) # inversion symmetry
tsymbed = embed(lattice, tsym)
psymbed = embed(lattice, psym)
## Setup Hilbert Space
(hs, σ) = ExactDiagonalization.Toolkit.spin_half_system(n_sites)
## Setup Operators
Sx = sum(σ(i,:x) for i in 1:n_sites)
Sy = sum(σ(i,:y) for i in 1:n_sites)
Sz = sum(σ(i,:z) for i in 1:n_sites)
spin_squared = simplify( Sx^2 + Sy^2 + Sz^2 )
j1 = NullOperator()
sub2ind(i1::Integer, i2::Integer) = (i1-1)%n1 + ((i2-1)%n2)*n1 + 1
ind2sub(i::Integer) = (mod(i, n1), mod(i ÷ n1, n2))
for i1 in 1:n1, i2 in 1:n2
global j1
j1 += sum(σ(sub2ind(i1, i2), μ) * σ(sub2ind(i1+1, i2), μ) for μ in [:x, :y, :z])
j1 += sum(σ(sub2ind(i1, i2), μ) * σ(sub2ind(i1, i2+1), μ) for μ in [:x, :y, :z])
end
j1 = simplify(j1)
## Full Solution
hsr = represent(hs)
m = Matrix(represent(hsr, j1))
alleigenvalues1 = eigvals(Hermitian(m))
@show length(alleigenvalues1)
## Using translation symmetry only
println("## Translation Symmetry Only")
if true
global alleigenvalues1
alleigenvalues2 = Float64[]
alleigenvalues3 = Float64[]
for qn in quantum_number_sectors(hs)
hss = HilbertSpaceSector(hs, qn)
hssr = represent_dict(hss)
# println(" quantum number: $qn")
# println(" number of irreps: $(num_irreps(tsym))")
for tsic in get_irrep_components(tsymbed)
#for tsym_irrep_index in 1:num_irreps(tsym)
#tsic = TranslationSymmetryIrrepComponent(tsym, tsym_irrep_index, 1)
rhssr = symmetry_reduce_serial(hssr, tsic)
rhssr2 = symmetry_reduce_parallel(hssr, tsic)
#j1_redrep = represent(rhssr, j1)
#m = Matrix(j1_redrep)
dimension(rhssr) == 0 && continue
# println("- QN: $qn\ttsym: $tsym_irrep_index/$(num_irreps(tsym))\tdimension: $(dimension(rhssr))")
# println(" momentum: $(tsym.hypercube.coordinates[tsym_irrep_index])")
# println(" hilbert dimension: $(dimension(rhssr))")
m = Matrix(represent(rhssr, j1))
m2 = Matrix(represent(rhssr2, j1))
append!(alleigenvalues2, eigvals(Hermitian(m)))
append!(alleigenvalues3, eigvals(Hermitian(m2)))
end
end
sort!(alleigenvalues2)
sort!(alleigenvalues3)
@show length(alleigenvalues2)
@show length(alleigenvalues3)
@show norm(alleigenvalues1 - alleigenvalues2)
@show norm(alleigenvalues1 - alleigenvalues3)
end
if true
println("## Point Symmetry Only")
alleigenvalues2 = Float64[]
alleigenvalues3 = Float64[]
for qn in quantum_number_sectors(hs)
hss = HilbertSpaceSector(hs, qn)
hssr = represent_dict(hss);
# println(" quantum number: $qn")
# println(" number of irreps: $(num_irreps(psym))")
# for psym_irrep_index in 1:num_irreps(psymbed)
for psic in get_irrep_components(psymbed)
# println(" psym irrep: $(psym.irreps[psym_irrep_index].name)")
# println(" psym irrep dimension: $(irrep_dimension(psym, psym_irrep_index))")
# for psym_irrep_compo in 1:irrep_dimension(psymbed, psym_irrep_index)
# psic = IrrepComponent(psymbed, psym_irrep_index, psym_irrep_compo)
rhssr = symmetry_reduce_serial(hssr, psic)
rhssr2 = symmetry_reduce_parallel(hssr, psic)
dimension(rhssr) == 0 && continue
# println("- QN: $qn\tpsym: $psym_irrep_index/$(num_irreps(psym))\t$psym_irrep_compo/$(irrep_dimension(psym, psym_irrep_index))\tdimension:$(dimension(rhssr))")
# println(" psym irrep component: $(psym_irrep_compo)")
# println(" hilbert dimension: $(dimension(rhssr))")
#
# for bvec in rhssr.basis_list
# println(string(bvec, base=2, pad=9))
# end
# println()
m = Matrix(represent(rhssr, j1))
m2 = Matrix(represent(rhssr2, j1))
append!(alleigenvalues2, eigvals(Hermitian(m)))
append!(alleigenvalues3, eigvals(Hermitian(m2)))
# end
end
end
sort!(alleigenvalues2)
sort!(alleigenvalues3)
@show length(alleigenvalues2)
@show length(alleigenvalues3)
@show norm(alleigenvalues1 - alleigenvalues2)
@show norm(alleigenvalues1 - alleigenvalues3)
end
## Use both
println("## Translation AND Point Symmetry")
let
alleigenvalues2 = Float64[]
alleigenvalues3 = Float64[]
for qn in quantum_number_sectors(hs)
hss = HilbertSpaceSector(hs, qn)
hssr = represent_dict(hss);
# println(" quantum number: $qn")
# println(" number of tsym_irreps: $(num_irreps(tsym))")
# # println(" momentum: $(tsym.hypercube.coordinates[tsym_irrep_index])")
# # println(" little_group: $(psym_little.hermann_mauguin)")
# # println(" number of irreps: $(num_irreps(psym_little))")
# # println(" psym_irrep_index: $psym_irrep_index")
# # println(" irrep dimension: $(irrep_dimension(psym_little, psym_irrep_index))")
for tsic in get_irrep_components(tsymbed)
psymbed_little = little_symmetry(tsic, psymbed)
for psic in get_irrep_components(psymbed_little)
ssic = SymmorphicIrrepComponent(tsic, psic)
rhssr = symmetry_reduce_serial(hssr, ssic)
rhssr2 = symmetry_reduce_parallel(hssr, ssic)
dimension(rhssr) == 0 && continue
# print("- QN: $qn")
# print("\ttsym: $tsym_irrep_index/$(num_irreps(tsym))")
# print("\tpsym: $(psymbed_little.symmetry.hermann_mauguin)")
# print("\t$psym_irrep_index/$(num_irreps(psymbed_little))")
# print("\tdimension:$(dimension(rhssr))")
# println()
# println(" psym_irrepo_compo: $(psym_irrep_compo)")
# println(" hilbert dimension: $(dimension(rhssr))")
m = Matrix(represent(rhssr, j1))
m2 = Matrix(represent(rhssr2, j1))
append!(alleigenvalues2, eigvals(Hermitian(m)))
append!(alleigenvalues3, eigvals(Hermitian(m2)))
# end
end
end
end
sort!(alleigenvalues2)
sort!(alleigenvalues3)
@show length(alleigenvalues2)
@show length(alleigenvalues3)
@show norm(alleigenvalues1 - alleigenvalues2)
@show norm(alleigenvalues1 - alleigenvalues3)
end
ssymbed = tsymbed ⋊ psymbed
println("## Symmorphic Space Symmetry (S = T ⋊ P)")
let
alleigenvalues2 = Float64[]
alleigenvalues3 = Float64[]
for qn in quantum_number_sectors(hs)
hss = HilbertSpaceSector(hs, qn)
hssr = represent_dict(hss);
# println(" quantum number: $qn")
# println(" number of tsym_irreps: $(num_irreps(tsym))")
for ssic in get_irrep_components(ssymbed)
rhssr = symmetry_reduce_serial(hssr, ssic)
rhssr2 = symmetry_reduce_parallel(hssr, ssic)
dimension(rhssr) == 0 && continue
# print("- QN: $qn")
# print("\ttsym: $(ssic.normal.irrep_index)/$(num_irreps(tsym))")
# print("\tpsym: $(symmetry_name(ssic.rest.symmetry))")
# print("\t$(ssic.rest.irrep_index)/$(num_irreps(ssic.rest.symmetry))")
# print("\tdimension:$(dimension(rhssr))")
# println()
m = Matrix(represent(rhssr, j1))
m2 = Matrix(represent(rhssr2, j1))
append!(alleigenvalues2, eigvals(Hermitian(m)))
append!(alleigenvalues3, eigvals(Hermitian(m2)))
end # for tsym_irrep_index
end
sort!(alleigenvalues2)
sort!(alleigenvalues3)
@show length(alleigenvalues2)
@show length(alleigenvalues3)
@show norm(alleigenvalues1 - alleigenvalues2)
@show norm(alleigenvalues1 - alleigenvalues3)
end