-
Notifications
You must be signed in to change notification settings - Fork 89
/
symmetry.jl
237 lines (209 loc) · 9.19 KB
/
symmetry.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
## This file contains functions to handle the symetries
# A symmetry operation (symop) is a couple (Stilde, τtilde) of a
# unitary (in cartesian coordinates, but not in reduced coordinates)
# matrix Stilde and a translation τtilde such that, for each atom of
# type A at position a, Stilde a + τtilde is also an atom of type A.
# This induces a symmetry in the Brillouin zone that the Hamiltonian
# at S k is unitary equivalent to that at k, which we exploit to
# reduce computations. The relationship is
# S = Stilde'
# τ = -Stilde^-1 τtilde
# (valid both in reduced and cartesian coordinates)
# The full (reducible) Brillouin zone is implicitly represented by
# a set of (irreducible) kpoints (see explanation in docs). Each
# irreducible kpoint k comes with a list of symmetry operations
# (S,τ) (containing at least the trivial operation (I,0)), where S
# is a rotation matrix (/!\ not unitary in reduced coordinates)
# and τ a translation vector. The kpoint is then used to represent
# implicitly the information at all the kpoints Sk. The
# relationship between the Hamiltonians is
# H_{Sk} = U H_k U*, with
# (Uu)(x) = u(S^-1(x-τ))
# or in Fourier space
# (Uu)(G) = e^{-i G τ} u(S^-1 G)
# In particular, we can choose the eigenvectors at Sk as u_{Sk} = U u_k
# We represent then the BZ as a set of irreducible points `kpoints`, a
# list of symmetry operations `ksymops` allowing the reconstruction of
# the full (reducible) BZ, and a set of weights `kweights` (summing to
# 1). The value of an observable (eg energy) per unit cell is given as
# the value of that observable at each irreducible kpoint, weighted by
# kweight
# There is by decreasing cardinality
# - The group of symmetry operations of the lattice
# - The group of symmetry operations of the crystal (model.symops)
# - The group of symmetry operations of the crystal that preserves the BZ mesh (basis.symops)
# - The set of symmetry operations that we use to reduce the
# reducible Brillouin zone (RBZ) to the irreducible (IBZ) (basis.ksymops)
@doc raw"""
Return the ``k``-point symmetry operations associated to a lattice, model or basis.
Since the ``k``-point discretisations may break some of the symmetries, the latter
case will return a subset of the symmetries of the former two.
"""
function symmetry_operations(lattice, atoms; tol_symmetry=1e-5, kcoords=nothing)
symops = []
# Get symmetries from spglib
Stildes, τtildes = spglib_get_symmetry(lattice, atoms, tol_symmetry=tol_symmetry)
for isym = 1:length(Stildes)
S = Stildes[isym]' # in fractional reciprocal coordinates
τ = -Stildes[isym] \ τtildes[isym] # in fractional real-space coordinates
τ = τ .- floor.(τ)
@assert all(0 .≤ τ .< 1)
push!(symops, (S, τ))
end
symops = unique(symops)
symops_preserving_kgrid(symops)
end
function symops_preserving_kgrid(symops, kcoords=nothing)
kcoords === nothing && return symops
# filter only the operations that respect the symmetries of the discrete BZ grid
function preserves_grid(S)
all(normalize_kpoint_coordinate(S * k) in kcoords
for k in normalize_kpoint_coordinate.(kcoords))
end
filter(symop -> preserves_grid(symop[1]), symops)
end
"""
Implements a primitive search to find an irreducible subset of kpoints
amongst the provided kpoints.
"""
function find_irreducible_kpoints(kcoords, Stildes, τtildes)
# This function is required because spglib sometimes flags kpoints
# as reducible, where we cannot find a symmetry operation to
# generate them from the provided irreducible kpoints. This
# reimplements that part of spglib, with a possibly very slow
# algorithm.
# Flag which kpoints have already been mapped to another irred.
# kpoint or which have been decided to be irreducible.
kcoords_mapped = zeros(Bool, length(kcoords))
kirreds = empty(kcoords) # Container for irreducible kpoints
ksymops = Vector{Vector{SymOp}}() # Corresponding symops
while !all(kcoords_mapped)
# Select next not mapped kpoint as irreducible
ik = findfirst(isequal(false), kcoords_mapped)
push!(kirreds, kcoords[ik])
thisk_symops = [identity_symop()]
kcoords_mapped[ik] = true
for jk in findall(.!kcoords_mapped)
isym = findfirst(1:length(Stildes)) do isym
# If the difference between kred and Stilde' * k == Stilde^{-1} * k
# is only integer in fractional reciprocal-space coordinates, then
# kred and S' * k are equivalent k-Points
all(isinteger, kcoords[jk] - (Stildes[isym]' * kcoords[ik]))
end
if !isnothing(isym) # Found a reducible kpoint
kcoords_mapped[jk] = true
S = Stildes[isym]' # in fractional reciprocal coordinates
τ = -Stildes[isym] \ τtildes[isym] # in fractional real-space coordinates
τ = τ .- floor.(τ)
@assert all(0 .≤ τ .< 1)
push!(thisk_symops, (S, τ))
end
end # jk
push!(ksymops, thisk_symops)
end
kirreds, ksymops
end
"""
Apply a symmetry operation to eigenvectors `ψk` at a given `kpoint` to obtain an
equivalent point in [-0.5, 0.5)^3 and associated eigenvectors (expressed in the
basis of the new kpoint).
"""
function apply_ksymop(ksymop, basis, kpoint, ψk::AbstractVecOrMat)
S, τ = ksymop
S == I && iszero(τ) && return kpoint, ψk
# Apply S and reduce coordinates to interval [-0.5, 0.5)
# Doing this reduction is important because
# of the finite kinetic energy basis cutoff
@assert all(-0.5 .≤ kpoint.coordinate .< 0.5)
Sk_raw = S * kpoint.coordinate
Sk = normalize_kpoint_coordinate(Sk_raw)
kshift = convert.(Int, Sk - Sk_raw)
@assert all(-0.5 .≤ Sk .< 0.5)
# Check whether the resulting kpoint is in the basis:
ikfull = findfirst(1:length(basis.kpoints)) do idx
all(isinteger, basis.kpoints[idx].coordinate - Sk)
end
if isnothing(ikfull)
# Build a new kpoint datastructure:
Skpoint = build_kpoints(basis, [Sk])[1]
else
Skpoint = basis.kpoints[ikfull]
@assert Skpoint.coordinate ≈ Sk
end
# Since the eigenfunctions of the Hamiltonian at k and Sk satisfy
# u_{Sk}(x) = u_{k}(S^{-1} (x - τ))
# their Fourier transform is related as
# ̂u_{Sk}(G) = e^{-i G \cdot τ} ̂u_k(S^{-1} G)
invS = Mat3{Int}(inv(S))
Gs_full = [G + kshift for G in G_vectors(Skpoint)]
ψSk = zero(ψk)
for iband in 1:size(ψk, 2)
for (ig, G_full) in enumerate(Gs_full)
igired = index_G_vectors(basis, kpoint, invS * G_full)
@assert igired !== nothing
ψSk[ig, iband] = cis(-2π * dot(G_full, τ)) * ψk[igired, iband]
end
end
Skpoint, ψSk
end
"""
Apply a `k`-point symmetry operation (the tuple (S, τ)) to a partial density.
"""
function apply_ksymop(symop, ρin::RealFourierArray)
symop[1] == I && iszero(symop[2]) && return ρin
from_fourier(ρin.basis, symmetrize(ρin, [symop]))
end
# Accumulates the symmetrized versions of the density ρin into ρout (in Fourier space).
# No normalization is performed
function accumulate_over_symops!(ρaccu, ρin, basis, symops, Gs)
T = eltype(basis)
for (S, τ) in symops
invS = Mat3{Int}(inv(S))
# Common special case, where ρin does not need to be processed
if iszero(S - I) && iszero(τ)
ρaccu .+= ρin
continue
end
# Transform ρin -> to the partial density at S * k.
#
# Since the eigenfunctions of the Hamiltonian at k and Sk satisfy
# u_{Sk}(x) = u_{k}(S^{-1} (x - τ))
# with Fourier transform
# ̂u_{Sk}(G) = e^{-i G \cdot τ} ̂u_k(S^{-1} G)
# equivalently
# ̂ρ_{Sk}(G) = e^{-i G \cdot τ} ̂ρ_k(S^{-1} G)
for (ig, G) in enumerate(Gs)
igired = index_G_vectors(basis, invS * G)
if igired !== nothing
@inbounds ρaccu[ig] += cis(-2T(π) * dot(G, τ)) * ρin[igired]
end
end
end # (S, τ)
ρaccu
end
# Low-pass filters ρ (in Fourier) so that symmetry operations acting on it stay in the grid
function lowpass_for_symmetry!(ρ, basis; symops=basis.model.symops)
for (S, τ) in symops
for (ig, G) in enumerate(G_vectors(basis))
if index_G_vectors(basis, S * G) === nothing
ρ[ig] = 0
end
end
end
ρ
end
"""
Symmetrize a `RealFourierArray` by applying all the model symmetries (by default) and forming the average.
"""
function symmetrize(ρin::RealFourierArray; symops=ρin.basis.model.symops)
ρin_fourier = copy(ρin.fourier)
lowpass_for_symmetry!(ρin_fourier, ρin.basis; symops=symops)
ρout_fourier = accumulate_over_symops!(zero(ρin_fourier), ρin_fourier, ρin.basis, symops,
G_vectors(ρin.basis)) ./ length(symops)
from_fourier(ρin.basis, ρout_fourier)
end
function check_symmetric(ρin::RealFourierArray; tol=1e-10, symops=ρin.basis.model.symops)
for symop in symops
@assert norm(symmetrize(ρin, [symop]).fourier - ρin.fourier) < tol
end
end