-
Notifications
You must be signed in to change notification settings - Fork 10
/
magnetic.jl
354 lines (328 loc) · 14.1 KB
/
magnetic.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
export MagneticDataset,
MagneticSpacegroupType,
get_symmetry_with_collinear_spin,
get_symmetry_with_site_tensors,
get_magnetic_dataset,
get_magnetic_symmetry_from_database,
get_magnetic_spacegroup_type,
get_magnetic_spacegroup_type_from_symmetry
# Python version: https://github.com/spglib/spglib/blob/42527b0/python/spglib/spglib.py#L182-L319
"""
get_symmetry_with_collinear_spin(cell::SpglibCell, symprec=1e-5)
Find symmetry operations with collinear polarizations (spins) on atoms.
Except for the magmoms in the `cell`, the usage is basically the same as [`get_symmetry`](@ref).
But as an output, `equivalent_atoms` are obtained as the last returned value. The size of this
array is the same of number of atoms in the cell.
"""
function get_symmetry_with_collinear_spin(cell::SpglibCell, symprec=1e-5)
lattice, positions, atoms, magmoms = _unwrap_convert(cell)
num_atom = length(cell.magmoms)
# See https://github.com/spglib/spglib/blob/42527b0/python/spglib/spglib.py#L270
max_size = 96num_atom # 96 = 48 × 2 since we have spins
rotations = Array{Cint,3}(undef, 3, 3, max_size)
translations = Matrix{Cdouble}(undef, 3, max_size)
equivalent_atoms = Vector{Cint}(undef, num_atom)
num_sym = @ccall libsymspg.spg_get_symmetry_with_collinear_spin(
rotations::Ptr{Cint},
translations::Ptr{Cdouble},
equivalent_atoms::Ptr{Cint},
max_size::Cint,
lattice::Ptr{Cdouble},
positions::Ptr{Cdouble},
atoms::Ptr{Cint},
magmoms::Ptr{Cdouble},
num_atom::Cint,
symprec::Cdouble,
)::Cint
check_error()
rotations = map(
SMatrix{3,3,Int32,9} ∘ transpose, eachslice(rotations[:, :, begin:num_sym]; dims=3)
) # Remember to transpose, see https://github.com/singularitti/Spglib.jl/blob/8aed6e0/src/core.jl#L195-L198
translations = map(SVector{3,Float64}, eachcol(translations[:, begin:num_sym]))
return rotations, translations, equivalent_atoms .+ 1
end
"""
get_symmetry_with_site_tensors(cell::SpglibCell, symprec=1e-5; with_time_reversal=true)
Return magnetic symmetry operations represented by `rotation`, `translation`, and `spin_flips`.
Returned `spin_flips` represents sign of site tensors after applying time-reversal operations:
``1`` for non time reversal, and ``-1`` for time reversal.
"""
function get_symmetry_with_site_tensors(
cell::SpglibCell, symprec=1e-5; with_time_reversal=true
)
lattice, positions, atoms, magmoms = _unwrap_convert(cell)
num_atom = natoms(cell)
# See https://github.com/spglib/spglib/blob/42527b0/python/spglib/spglib.py#L270
max_size = 96num_atom # 96 = 48 × 2 since we have spins
rotations = Array{Cint,3}(undef, 3, 3, max_size)
translations = Matrix{Cdouble}(undef, 3, max_size)
equivalent_atoms = Vector{Cint}(undef, num_atom)
primitive_lattice = zeros(Cdouble, 3, 3)
spin_flips = zeros(Cint, max_size)
tensor_rank = ndims(magmoms) - 1 # See https://github.com/spglib/spglib/blob/v2.1.0/python/spglib/spglib.py#L275-L276 & https://github.com/spglib/spglib/blob/v2.1.0/python/spglib/spglib.py#L615
is_axial = iszero(tensor_rank) ? false : true # Collinear spin & non-collinear spin
num_sym = @ccall libsymspg.spg_get_symmetry_with_site_tensors(
rotations::Ptr{Cint},
translations::Ptr{Cdouble},
equivalent_atoms::Ptr{Cint},
primitive_lattice::Ptr{Cdouble},
spin_flips::Ptr{Cint},
max_size::Cint,
lattice::Ptr{Cdouble},
positions::Ptr{Cdouble},
atoms::Ptr{Cint},
magmoms::Ptr{Cdouble},
tensor_rank::Cint,
num_atom::Cint,
with_time_reversal::Cint,
is_axial::Cint,
symprec::Cdouble,
)::Cint
check_error()
rotations = map(
SMatrix{3,3,Int32,9} ∘ transpose, eachslice(rotations[:, :, begin:num_sym]; dims=3)
) # Remember to transpose, see https://github.com/singularitti/Spglib.jl/blob/8aed6e0/src/core.jl#L195-L198
translations = map(SVector{3,Float64}, eachcol(translations[:, begin:num_sym]))
return rotations, translations, spin_flips[begin:num_sym]
end
struct SpglibMagneticDataset <: AbstractDataset
uni_number::Cint
msg_type::Cint
hall_number::Cint
tensor_rank::Cint
n_operations::Cint
rotations::Ptr{NTuple{3,NTuple{3,Cint}}}
translations::Ptr{NTuple{3,Cdouble}}
time_reversals::Ptr{Cint}
n_atoms::Cint
equivalent_atoms::Ptr{Cint}
transformation_matrix::NTuple{3,NTuple{3,Cdouble}}
origin_shift::NTuple{3,Cdouble}
n_std_atoms::Cint
std_lattice::NTuple{3,NTuple{3,Cdouble}}
std_types::Ptr{Cint}
std_positions::Ptr{NTuple{3,Cdouble}}
std_tensors::Ptr{Cdouble}
std_rotation_matrix::NTuple{3,NTuple{3,Cdouble}}
primitive_lattice::NTuple{3,NTuple{3,Cdouble}}
end
"""
MagneticDataset(uni_number, msg_type, hall_number, tensor_rank, n_operations, rotations, translations, time_reversals, n_atoms, equivalent_atoms, transformation_matrix, origin_shift, n_std_atoms, std_lattice, std_types, std_positions, std_tensors, std_rotation_matrix, primitive_lattice)
Represent `MagneticDataset`, see its [official documentation](https://spglib.readthedocs.io/en/latest/magnetic_dataset.html).
See also [`get_magnetic_dataset`](@ref).
"""
@struct_hash_equal_isequal struct MagneticDataset <: AbstractDataset
uni_number::Int32
msg_type::Int32
hall_number::Int32
tensor_rank::Int32
n_operations::Int32
rotations::Vector{SMatrix{3,3,Int32,9}}
translations::Vector{SVector{3,Float64}}
time_reversals::BitVector
n_atoms::Int32
equivalent_atoms::Vector{Int32}
transformation_matrix::SMatrix{3,3,Float64,9}
origin_shift::SVector{3,Float64}
n_std_atoms::Int32
std_lattice::Lattice{Float64}
std_types::Vector{Int32}
std_positions::Vector{SVector{3,Float64}}
std_tensors::Vector{Union{Float64,SVector{3,Float64}}}
std_rotation_matrix::SMatrix{3,3,Float64,9}
primitive_lattice::Lattice{Float64}
end
function MagneticDataset(dataset::SpglibMagneticDataset)
rotations = transpose.(
_convert(SMatrix{3,3,Int32}, unsafe_load(dataset.rotations, i)) for
i in Base.OneTo(dataset.n_operations)
)
translations = SVector{3}.(
unsafe_load(dataset.translations, i) for i in Base.OneTo(dataset.n_operations)
)
time_reversals =
Bool.(unsafe_wrap(Vector{Int32}, dataset.time_reversals, dataset.n_operations))
equivalent_atoms = # Need to add 1 because of C-index starts from 0
unsafe_wrap(Vector{Int32}, dataset.equivalent_atoms, dataset.n_atoms) .+ 1
transformation_matrix = transpose(
_convert(SMatrix{3,3,Float64}, dataset.transformation_matrix)
)
std_lattice = Lattice(transpose(_convert(SMatrix{3,3,Float64}, dataset.std_lattice)))
std_types = unsafe_wrap(Vector{Int32}, dataset.std_types, dataset.n_std_atoms)
std_positions = SVector{3}.(
unsafe_load(dataset.std_positions, i) for i in Base.OneTo(dataset.n_std_atoms)
)
std_tensors = if iszero(dataset.tensor_rank) # Collinear spin
unsafe_wrap(Vector{Float64}, dataset.std_tensors, dataset.n_std_atoms)
else # Non-collinear spin
SVector{3}.(
eachcol(
unsafe_wrap(
Matrix{Float64},
dataset.std_tensors,
(3, Int(dataset.n_std_atoms)), # Issue to Julia community
),
),
)
end
std_rotation_matrix = transpose(
_convert(SMatrix{3,3,Float64}, dataset.std_rotation_matrix)
)
primitive_lattice = Lattice(
transpose(_convert(SMatrix{3,3,Float64}, dataset.primitive_lattice))
)
return MagneticDataset(
dataset.uni_number,
dataset.msg_type,
dataset.hall_number,
dataset.tensor_rank,
dataset.n_operations,
rotations,
translations,
time_reversals,
dataset.n_atoms,
equivalent_atoms,
transformation_matrix,
dataset.origin_shift,
dataset.n_std_atoms,
std_lattice,
std_types,
std_positions,
std_tensors,
std_rotation_matrix,
primitive_lattice,
)
end
"""
get_magnetic_dataset(cell::SpglibCell, symprec=1e-5)
Return magnetic symmetry operations and standardized structure of given structure with site tensors.
The description of returned dataset is given at [Magnetic dataset (experimental)](@ref).
"""
function get_magnetic_dataset(cell::SpglibCell, symprec=1e-5)
lattice, positions, atoms, magmoms = _unwrap_convert(cell)
tensor_rank = ndims(magmoms) - 1 # See https://github.com/spglib/spglib/blob/v2.1.0/python/spglib/spglib.py#L275-L276 & https://github.com/spglib/spglib/blob/v2.1.0/python/spglib/spglib.py#L615
is_axial = iszero(tensor_rank) ? false : true # Collinear spin & non-collinear spin
ptr = @ccall libsymspg.spg_get_magnetic_dataset(
lattice::Ptr{Cdouble},
positions::Ptr{Cdouble},
atoms::Ptr{Cint},
magmoms::Ptr{Cdouble},
tensor_rank::Cint,
natoms(cell)::Cint,
is_axial::Cint,
symprec::Cdouble,
)::Ptr{SpglibMagneticDataset}
if ptr == C_NULL
check_error()
else
raw_dataset = unsafe_load(ptr)
return MagneticDataset(raw_dataset)
end
end
# Python version: https://github.com/spglib/spglib/blob/v2.1.0/python/spglib/spglib.py#L1108-L1138
"""
get_magnetic_symmetry_from_database(uni_number, hall_number=0)
Accesses magnetic space-group operations in the built-in database using UNI number (from ``1`` to ``1651``).
Optionally alternative settings can be specified with hall_number. For type-I, type-II, and
type-III magnetic space groups, `hall_number` changes settings in family space group. For
type-IV, `hall_number` changes settings in maximal space group. When `hall_number=0`, the
smallest hall number corresponding to `uni_number` is used.
"""
function get_magnetic_symmetry_from_database(uni_number, hall_number=0)
@assert 1 <= uni_number <= 1651 # See https://github.com/spglib/spglib/blob/77a8e5d/src/spglib.h#L390
_rotations = Array{Cint,3}(undef, 3, 3, 384) # See https://github.com/spglib/spglib/blob/v2.1.0/src/spglib.h#L398
_translations = Matrix{Cdouble}(undef, 3, 384)
_time_reversals = Vector{Cint}(undef, 384) # It cannot be `Vector{Bool}`!
num_sym = @ccall libsymspg.spg_get_magnetic_symmetry_from_database(
_rotations::Ptr{Cint},
_translations::Ptr{Cdouble},
_time_reversals::Ptr{Cint},
uni_number::Cint,
hall_number::Cint,
)::Cint
check_error()
rotations = map(
SMatrix{3,3,Int32,9} ∘ transpose, eachslice(_rotations[:, :, begin:num_sym]; dims=3)
) # Remember to transpose, see https://github.com/singularitti/Spglib.jl/blob/8aed6e0/src/core.jl#L195-L198
translations = map(SVector{3,Float64}, eachcol(_translations[:, begin:num_sym]))
time_reversals = Bool.(_time_reversals[begin:num_sym])
return rotations, translations, time_reversals
end
struct SpglibMagneticSpacegroupType <: AbstractSpacegroupType
uni_number::Cint
litvin_number::Cint
bns_number::NTuple{8,Cchar}
og_number::NTuple{12,Cchar}
number::Cint
type::Cint
end
"""
MagneticSpacegroupType(uni_number, litvin_number, bns_number, og_number, number, type)
Represent `SpglibMagneticSpacegroupType`, see its [official documentation](https://spglib.readthedocs.io/en/latest/api.html#spg-get-magnetic-spacegroup-type).
# Arguments
- `uni_number::Int32`: Serial number of UNI (or BNS) symbols.
- `litvin_number::Int32`: Serial number in Litvin's [Magnetic Group Tables](https://www.iucr.org/publ/978-0-9553602-2-0).
- `bns_number::String`: BNS number, e.g. `"151.32"`.
- `og_number::String`: OG number, e.g. `"153.4.1270"`.
- `number::Int32`: ITA's serial number of space group for reference setting.
- `type::Int32`: Type of MSG, from ``1`` to ``4``.
See also [`get_magnetic_spacegroup_type`](@ref), [`get_magnetic_spacegroup_type_from_symmetry`](@ref).
"""
struct MagneticSpacegroupType <: AbstractSpacegroupType
uni_number::Int32
litvin_number::Int32
bns_number::String
og_number::String
number::Int32
type::Int32
end
function MagneticSpacegroupType(spgtype::SpglibMagneticSpacegroupType)
bns_number = tostring(spgtype.bns_number)
og_number = tostring(spgtype.og_number)
return MagneticSpacegroupType(
spgtype.uni_number,
spgtype.litvin_number,
bns_number,
og_number,
spgtype.number,
spgtype.type,
)
end
"""
get_magnetic_spacegroup_type(uni_number)
Accesses magnetic space-group type by serial number (from ``1`` to ``1651``) of UNI (or BNS) symbols.
"""
function get_magnetic_spacegroup_type(uni_number)
spgtype = @ccall libsymspg.spg_get_magnetic_spacegroup_type(
uni_number::Cint
)::SpglibMagneticSpacegroupType
check_error()
return MagneticSpacegroupType(spgtype)
end
"""
get_magnetic_spacegroup_type_from_symmetry(rotations, translations, time_reversals, lattice::Lattice, symprec=1e-5)
Determine magnetic space-group type from magnetic symmetry operations.
`time_reversals` takes `false` for ordinary operations and `true` for time-reversal operations.
"""
function get_magnetic_spacegroup_type_from_symmetry(
rotations, translations, time_reversals, lattice::Lattice, symprec=1e-5
)
if length(rotations) != length(translations)
throw(DimensionMismatch("the numbers of rotations and translations are different!"))
end
num_sym = length(translations)
rotations = Base.cconvert(Array{Cint,3}, cat(transpose.(rotations)...; dims=3))
translations = Base.cconvert(Matrix{Cdouble}, reduce(hcat, translations))
time_reversals = Base.cconvert(Vector{Cint}, time_reversals)
lattice = Base.cconvert(Matrix{Cdouble}, transpose(lattice)) # `transpose` must before `cconvert`!
spgtype = @ccall libsymspg.spg_get_magnetic_spacegroup_type_from_symmetry(
rotations::Ptr{Cint},
translations::Ptr{Cdouble},
time_reversals::Ptr{Cint},
num_sym::Cint,
lattice::Ptr{Cdouble},
symprec::Cdouble,
)::SpglibMagneticSpacegroupType
check_error()
return MagneticSpacegroupType(spgtype)
end