/
bandstructure.py
624 lines (552 loc) · 23.9 KB
/
bandstructure.py
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
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
"""
This module provides classes to define a phonon band structure.
"""
from __future__ import annotations
import numpy as np
from monty.json import MSONable
from pymatgen.core.lattice import Lattice
from pymatgen.core.structure import Structure
from pymatgen.electronic_structure.bandstructure import Kpoint
def get_reasonable_repetitions(n_atoms: int) -> tuple[int, int, int]:
"""
Choose the number of repetitions in a supercell
according to the number of atoms in the system
"""
if n_atoms < 4:
return (3, 3, 3)
if 4 <= n_atoms < 15:
return (2, 2, 2)
if 15 <= n_atoms < 50:
return (2, 2, 1)
return (1, 1, 1)
def eigenvectors_from_displacements(disp, masses):
"""
Calculate the eigenvectors from the atomic displacements
"""
sqrt_masses = np.sqrt(masses)
return np.einsum("nax,a->nax", disp, sqrt_masses)
def estimate_band_connection(prev_eigvecs, eigvecs, prev_band_order):
"""
A function to order the phonon eigenvectors taken from phonopy
"""
metric = np.abs(np.dot(prev_eigvecs.conjugate().T, eigvecs))
connection_order = []
for overlaps in metric:
max_val = 0
for i in reversed(range(len(metric))):
val = overlaps[i]
if i in connection_order:
continue
if val > max_val:
max_val = val
max_idx = i
connection_order.append(max_idx)
band_order = [connection_order[x] for x in prev_band_order]
return band_order
class PhononBandStructure(MSONable):
"""
This is the most generic phonon band structure data possible
it's defined by a list of qpoints + frequencies for each of them.
Additional information may be given for frequencies at Gamma, where
non-analytical contribution may be taken into account.
"""
def __init__(
self,
qpoints: list[Kpoint],
frequencies: np.ndarray,
lattice: Lattice,
nac_frequencies=None,
eigendisplacements=None,
nac_eigendisplacements=None,
labels_dict=None,
coords_are_cartesian=False,
structure: Structure | None = None,
):
"""
Args:
qpoints: list of qpoint as numpy arrays, in frac_coords of the
given lattice by default
frequencies: list of phonon frequencies in THz as a numpy array with shape
(3*len(structure), len(qpoints)). The First index of the array
refers to the band and the second to the index of the qpoint.
lattice: The reciprocal lattice as a pymatgen Lattice object.
Pymatgen uses the physics convention of reciprocal lattice vectors
WITH a 2*pi coefficient.
nac_frequencies: Frequencies with non-analytical contributions at Gamma in THz.
A list of tuples. The first element of each tuple should be a list
defining the direction (not necessarily a versor, will be normalized
internally). The second element containing the 3*len(structure)
phonon frequencies with non-analytical correction for that direction.
eigendisplacements: the phonon eigendisplacements associated to the
frequencies in Cartesian coordinates. A numpy array of complex
numbers with shape (3*len(structure), len(qpoints), len(structure), 3).
he First index of the array refers to the band, the second to the index
of the qpoint, the third to the atom in the structure and the fourth
to the Cartesian coordinates.
nac_eigendisplacements: the phonon eigendisplacements associated to the
non-analytical frequencies in nac_frequencies in Cartesian coordinates.
A list of tuples. The first element of each tuple should be a list
defining the direction. The second element containing a numpy array of
complex numbers with shape (3*len(structure), len(structure), 3).
labels_dict: (dict) of {} this links a qpoint (in frac coords or
Cartesian coordinates depending on the coords) to a label.
coords_are_cartesian: Whether the qpoint coordinates are Cartesian.
structure: The crystal structure (as a pymatgen Structure object)
associated with the band structure. This is needed if we
provide projections to the band structure
"""
self.lattice_rec = lattice
self.qpoints = []
self.labels_dict = {}
self.structure = structure
if eigendisplacements is None:
eigendisplacements = np.array([])
self.eigendisplacements = eigendisplacements
if labels_dict is None:
labels_dict = {}
for q_pt in qpoints:
# let see if this qpoint has been assigned a label
label = None
for c in labels_dict:
if np.linalg.norm(q_pt - np.array(labels_dict[c])) < 0.0001:
label = c
self.labels_dict[label] = Kpoint(
q_pt,
lattice,
label=label,
coords_are_cartesian=coords_are_cartesian,
)
self.qpoints.append(Kpoint(q_pt, lattice, label=label, coords_are_cartesian=coords_are_cartesian))
self.bands = frequencies
self.nb_bands = len(self.bands)
self.nb_qpoints = len(self.qpoints)
# normalize directions for nac_frequencies and nac_eigendisplacements
self.nac_frequencies = []
self.nac_eigendisplacements = []
if nac_frequencies is not None:
for t in nac_frequencies:
self.nac_frequencies.append(([i / np.linalg.norm(t[0]) for i in t[0]], t[1]))
if nac_eigendisplacements is not None:
for t in nac_eigendisplacements:
self.nac_eigendisplacements.append(([i / np.linalg.norm(t[0]) for i in t[0]], t[1]))
def min_freq(self) -> tuple[Kpoint, float]:
"""
Returns the point where the minimum frequency is reached and its value
"""
i = np.unravel_index(np.argmin(self.bands), self.bands.shape)
return self.qpoints[i[1]], self.bands[i]
def has_imaginary_freq(self, tol: float = 1e-5) -> bool:
"""
True if imaginary frequencies are present in the BS.
"""
return self.min_freq()[1] + tol < 0
@property
def has_nac(self) -> bool:
"""
True if nac_frequencies are present.
"""
return len(self.nac_frequencies) > 0
@property
def has_eigendisplacements(self) -> bool:
"""
True if eigendisplacements are present.
"""
return len(self.eigendisplacements) > 0
def get_nac_frequencies_along_dir(self, direction):
"""
Returns the nac_frequencies for the given direction (not necessarily a versor).
None if the direction is not present or nac_frequencies has not been calculated.
Args:
direction: the direction as a list of 3 elements
Returns:
the frequencies as a numpy array o(3*len(structure), len(qpoints)).
None if not found.
"""
versor = [i / np.linalg.norm(direction) for i in direction]
for d, f in self.nac_frequencies:
if np.allclose(versor, d):
return f
return None
def get_nac_eigendisplacements_along_dir(self, direction):
"""
Returns the nac_eigendisplacements for the given direction (not necessarily a versor).
None if the direction is not present or nac_eigendisplacements has not been calculated.
Args:
direction: the direction as a list of 3 elements
Returns:
the eigendisplacements as a numpy array of complex numbers with shape
(3*len(structure), len(structure), 3). None if not found.
"""
versor = [i / np.linalg.norm(direction) for i in direction]
for d, e in self.nac_eigendisplacements:
if np.allclose(versor, d):
return e
return None
def asr_breaking(self, tol_eigendisplacements=1e-5):
"""
Returns the breaking of the acoustic sum rule for the three acoustic modes,
if Gamma is present. None otherwise.
If eigendisplacements are available they are used to determine the acoustic
modes: selects the bands corresponding to the eigendisplacements that
represent to a translation within tol_eigendisplacements. If these are not
identified or eigendisplacements are missing the first 3 modes will be used
(indices [0:3]).
"""
for i in range(self.nb_qpoints):
if np.allclose(self.qpoints[i].frac_coords, (0, 0, 0)):
if self.has_eigendisplacements:
acoustic_modes_index = []
for j in range(self.nb_bands):
eig = self.eigendisplacements[j][i]
if np.max(np.abs(eig[1:] - eig[:1])) < tol_eigendisplacements:
acoustic_modes_index.append(j)
# if acoustic modes are not correctly identified return use
# the first three modes
if len(acoustic_modes_index) != 3:
acoustic_modes_index = [0, 1, 2]
return self.bands[acoustic_modes_index, i]
return self.bands[:3, i]
return None
def as_dict(self):
"""
:return: MSONable dict
"""
d = {
"@module": type(self).__module__,
"@class": type(self).__name__,
"lattice_rec": self.lattice_rec.as_dict(),
"qpoints": [],
}
# qpoints are not Kpoint objects dicts but are frac coords.Tthis makes
# the dict smaller and avoids the repetition of the lattice
for q in self.qpoints:
d["qpoints"].append(q.as_dict()["fcoords"])
d["bands"] = self.bands.tolist()
d["labels_dict"] = {}
for kpoint_letter, kpoint_object in self.labels_dict.items():
d["labels_dict"][kpoint_letter] = kpoint_object.as_dict()["fcoords"]
# split the eigendisplacements to real and imaginary part for serialization
d["eigendisplacements"] = {
"real": np.real(self.eigendisplacements).tolist(),
"imag": np.imag(self.eigendisplacements).tolist(),
}
d["nac_eigendisplacements"] = [
(direction, {"real": np.real(e).tolist(), "imag": np.imag(e).tolist()})
for direction, e in self.nac_eigendisplacements
]
d["nac_frequencies"] = [(direction, f.tolist()) for direction, f in self.nac_frequencies]
if self.structure:
d["structure"] = self.structure.as_dict()
return d
@classmethod
def from_dict(cls, d):
"""
:param d: Dict representation
:return: PhononBandStructure
"""
lattice_rec = Lattice(d["lattice_rec"]["matrix"])
eigendisplacements = np.array(d["eigendisplacements"]["real"]) + np.array(d["eigendisplacements"]["imag"]) * 1j
nac_eigendisplacements = [
(direction, np.array(e["real"]) + np.array(e["imag"]) * 1j) for direction, e in d["nac_eigendisplacements"]
]
nac_frequencies = [(direction, np.array(f)) for direction, f in d["nac_frequencies"]]
structure = Structure.from_dict(d["structure"]) if "structure" in d else None
return cls(
d["qpoints"],
np.array(d["bands"]),
lattice_rec,
nac_frequencies,
eigendisplacements,
nac_eigendisplacements,
d["labels_dict"],
structure=structure,
)
class PhononBandStructureSymmLine(PhononBandStructure):
r"""
This object stores phonon band structures along selected (symmetry) lines in the
Brillouin zone. We call the different symmetry lines (ex: \\Gamma to Z)
"branches".
"""
def __init__(
self,
qpoints,
frequencies,
lattice,
has_nac=False,
eigendisplacements=None,
labels_dict=None,
coords_are_cartesian=False,
structure=None,
):
"""
Args:
qpoints: list of qpoints as numpy arrays, in frac_coords of the
given lattice by default
frequencies: list of phonon frequencies in eV as a numpy array with shape
(3*len(structure), len(qpoints))
lattice: The reciprocal lattice as a pymatgen Lattice object.
Pymatgen uses the physics convention of reciprocal lattice vectors
WITH a 2*pi coefficient
has_nac: specify if the band structure has been produced taking into account
non-analytical corrections at Gamma. If True frequencies at Gamma from
different directions will be stored in naf. Default False.
eigendisplacements: the phonon eigendisplacements associated to the
frequencies in Cartesian coordinates. A numpy array of complex
numbers with shape (3*len(structure), len(qpoints), len(structure), 3).
he First index of the array refers to the band, the second to the index
of the qpoint, the third to the atom in the structure and the fourth
to the Cartesian coordinates.
labels_dict: (dict) of {} this links a qpoint (in frac coords or
Cartesian coordinates depending on the coords) to a label.
coords_are_cartesian: Whether the qpoint coordinates are cartesian.
structure: The crystal structure (as a pymatgen Structure object)
associated with the band structure. This is needed if we
provide projections to the band structure
"""
super().__init__(
qpoints=qpoints,
frequencies=frequencies,
lattice=lattice,
nac_frequencies=None,
eigendisplacements=eigendisplacements,
nac_eigendisplacements=None,
labels_dict=labels_dict,
coords_are_cartesian=coords_are_cartesian,
structure=structure,
)
self._reuse_init(eigendisplacements, frequencies, has_nac, qpoints)
def _reuse_init(self, eigendisplacements, frequencies, has_nac, qpoints):
self.distance = []
self.branches = []
one_group = []
branches_tmp = []
# get labels and distance for each qpoint
previous_qpoint = self.qpoints[0]
previous_distance = 0.0
previous_label = self.qpoints[0].label
for i in range(self.nb_qpoints):
label = self.qpoints[i].label
if label is not None and previous_label is not None:
self.distance.append(previous_distance)
else:
self.distance.append(
np.linalg.norm(self.qpoints[i].cart_coords - previous_qpoint.cart_coords) + previous_distance
)
previous_qpoint = self.qpoints[i]
previous_distance = self.distance[i]
if label and previous_label:
if len(one_group) != 0:
branches_tmp.append(one_group)
one_group = []
previous_label = label
one_group.append(i)
if len(one_group) != 0:
branches_tmp.append(one_group)
for b in branches_tmp:
self.branches.append(
{
"start_index": b[0],
"end_index": b[-1],
"name": str(self.qpoints[b[0]].label) + "-" + str(self.qpoints[b[-1]].label),
}
)
# extract the frequencies with non-analytical contribution at gamma
if has_nac:
naf = []
nac_eigendisplacements = []
for i in range(self.nb_qpoints):
# get directions with nac irrespectively of the label_dict. NB: with labels
# the gamma point is expected to appear twice consecutively.
if np.allclose(qpoints[i], (0, 0, 0)):
if i > 0 and not np.allclose(qpoints[i - 1], (0, 0, 0)):
q_dir = self.qpoints[i - 1]
direction = q_dir.frac_coords / np.linalg.norm(q_dir.frac_coords)
naf.append((direction, frequencies[:, i]))
if self.has_eigendisplacements:
nac_eigendisplacements.append((direction, eigendisplacements[:, i]))
if i < len(qpoints) - 1 and not np.allclose(qpoints[i + 1], (0, 0, 0)):
q_dir = self.qpoints[i + 1]
direction = q_dir.frac_coords / np.linalg.norm(q_dir.frac_coords)
naf.append((direction, frequencies[:, i]))
if self.has_eigendisplacements:
nac_eigendisplacements.append((direction, eigendisplacements[:, i]))
self.nac_frequencies = np.array(naf, dtype=object)
self.nac_eigendisplacements = np.array(nac_eigendisplacements, dtype=object)
def get_equivalent_qpoints(self, index):
"""
Returns the list of qpoint indices equivalent (meaning they are the
same frac coords) to the given one.
Args:
index: the qpoint index
Returns:
a list of equivalent indices
TODO: now it uses the label we might want to use coordinates instead
(in case there was a mislabel)
"""
# if the qpoint has no label it can"t have a repetition along the band
# structure line object
if self.qpoints[index].label is None:
return [index]
list_index_qpoints = []
for i in range(self.nb_qpoints):
if self.qpoints[i].label == self.qpoints[index].label:
list_index_qpoints.append(i)
return list_index_qpoints
def get_branch(self, index):
r"""
Returns in what branch(es) is the qpoint. There can be several
branches.
Args:
index: the qpoint index
Returns:
A list of dictionaries [{"name","start_index","end_index","index"}]
indicating all branches in which the qpoint is. It takes into
account the fact that one qpoint (e.g., \\Gamma) can be in several
branches
"""
to_return = []
for i in self.get_equivalent_qpoints(index):
for b in self.branches:
if b["start_index"] <= i <= b["end_index"]:
to_return.append(
{
"name": b["name"],
"start_index": b["start_index"],
"end_index": b["end_index"],
"index": i,
}
)
return to_return
def write_phononwebsite(self, filename):
"""
Write a json file for the phononwebsite:
http://henriquemiranda.github.io/phononwebsite
"""
import json
with open(filename, "w") as f:
json.dump(self.as_phononwebsite(), f)
def as_phononwebsite(self):
"""
Return a dictionary with the phononwebsite format:
http://henriquemiranda.github.io/phononwebsite
"""
d = {}
# define the lattice
d["lattice"] = self.structure.lattice._matrix.tolist()
# define atoms
atom_pos_car = []
atom_pos_red = []
atom_types = []
for site in self.structure.sites:
atom_pos_car.append(site.coords.tolist())
atom_pos_red.append(site.frac_coords.tolist())
atom_types.append(site.species_string)
# default for now
d["repetitions"] = get_reasonable_repetitions(len(atom_pos_car))
d["natoms"] = len(atom_pos_car)
d["atom_pos_car"] = atom_pos_car
d["atom_pos_red"] = atom_pos_red
d["atom_types"] = atom_types
d["atom_numbers"] = self.structure.atomic_numbers
d["formula"] = self.structure.formula
d["name"] = self.structure.formula
# get qpoints
qpoints = []
for q in self.qpoints:
qpoints.append(list(q.frac_coords))
d["qpoints"] = qpoints
# get labels
hsq_dict = {}
for nq, q in enumerate(self.qpoints):
if q.label is not None:
hsq_dict[nq] = q.label
# get distances
dist = 0
nqstart = 0
distances = [dist]
line_breaks = []
for nq in range(1, len(qpoints)):
q1 = np.array(qpoints[nq])
q2 = np.array(qpoints[nq - 1])
# detect jumps
if (nq in hsq_dict) and (nq - 1 in hsq_dict):
if hsq_dict[nq] != hsq_dict[nq - 1]:
hsq_dict[nq - 1] += "|" + hsq_dict[nq]
del hsq_dict[nq]
line_breaks.append((nqstart, nq))
nqstart = nq
else:
dist += np.linalg.norm(q1 - q2)
distances.append(dist)
line_breaks.append((nqstart, len(qpoints)))
d["distances"] = distances
d["line_breaks"] = line_breaks
d["highsym_qpts"] = list(hsq_dict.items())
# eigenvalues
thz2cm1 = 33.35641
bands = self.bands.copy() * thz2cm1
d["eigenvalues"] = bands.T.tolist()
# eigenvectors
eigenvectors = self.eigendisplacements.copy()
eigenvectors /= np.linalg.norm(eigenvectors[0, 0])
eigenvectors = eigenvectors.swapaxes(0, 1)
eigenvectors = np.array([eigenvectors.real, eigenvectors.imag])
eigenvectors = np.rollaxis(eigenvectors, 0, 5)
d["vectors"] = eigenvectors.tolist()
return d
def band_reorder(self):
"""
Re-order the eigenvalues according to the similarity of the eigenvectors
"""
eiv = self.eigendisplacements
eig = self.bands
nphonons, nqpoints = self.bands.shape
order = np.zeros([nqpoints, nphonons], dtype=int)
order[0] = np.array(range(nphonons))
# get the atomic masses
atomic_masses = [site.specie.atomic_mass for site in self.structure.sites]
# get order
for nq in range(1, nqpoints):
old_eiv = eigenvectors_from_displacements(eiv[:, nq - 1], atomic_masses)
new_eiv = eigenvectors_from_displacements(eiv[:, nq], atomic_masses)
order[nq] = estimate_band_connection(
old_eiv.reshape([nphonons, nphonons]).T,
new_eiv.reshape([nphonons, nphonons]).T,
order[nq - 1],
)
# reorder
for nq in range(1, nqpoints):
eivq = eiv[:, nq]
eigq = eig[:, nq]
eiv[:, nq] = eivq[order[nq]]
eig[:, nq] = eigq[order[nq]]
def as_dict(self):
"""
Returns: MSONable dict
"""
d = super().as_dict()
# remove nac_frequencies and nac_eigendisplacements as they are reconstructed
# in the __init__ when the dict is deserialized
nac_frequencies = d.pop("nac_frequencies")
d.pop("nac_eigendisplacements")
d["has_nac"] = len(nac_frequencies) > 0
return d
@classmethod
def from_dict(cls, d):
"""
Args:
d: Dict representation
Returns: PhononBandStructureSummLine
"""
lattice_rec = Lattice(d["lattice_rec"]["matrix"])
eigendisplacements = np.array(d["eigendisplacements"]["real"]) + np.array(d["eigendisplacements"]["imag"]) * 1j
structure = Structure.from_dict(d["structure"]) if "structure" in d else None
return cls(
d["qpoints"],
np.array(d["bands"]),
lattice_rec,
d["has_nac"],
eigendisplacements,
d["labels_dict"],
structure=structure,
)