Hide keyboard shortcuts

Hot-keys on this page

r m x p   toggle line displays

j k   next/prev highlighted chunk

0   (zero) top of page

1   (one) first highlighted chunk

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

625

626

627

628

629

630

631

632

633

634

635

636

637

638

639

640

641

642

643

644

645

646

647

648

649

650

651

652

653

654

655

656

657

658

659

660

661

662

663

664

665

666

667

668

669

670

671

672

673

674

675

""" 

This module provides functionality for storing and handling of force constants. 

""" 

 

import numpy as np 

import tarfile 

 

from abc import ABC, abstractmethod 

from itertools import product 

from string import ascii_lowercase, ascii_uppercase 

from scipy.linalg import eig 

from ase import units 

 

from .io import shengBTE as shengBTEIO 

from .io import phonopy as phonopyIO 

from .io.read_write_files import add_items_to_tarfile_pickle, read_items_pickle,\ 

add_ase_atoms_to_tarfile, read_ase_atoms_from_tarfile 

 

 

class ForceConstants(ABC): 

""" Base class for force constants """ 

 

def __init__(self, supercell): 

24 ↛ 25line 24 didn't jump to line 25, because the condition on line 24 was never true if not all(supercell.pbc): 

raise ValueError('configuration must have periodic boundary conditions') 

self._supercell = supercell.copy() 

 

@abstractmethod 

def __getitem__(self): 

pass 

 

@property 

def n_atoms(self): 

""" int : number of atoms """ 

return len(self.supercell) 

 

@property 

def supercell(self): 

""" ase.Atoms : supercell associated with force constants """ 

return self._supercell.copy() 

 

@property 

def clusters(self): 

""" list : sorted list of clusters """ 

return sorted(self._fc_dict.keys(), key=lambda c: (len(c), c)) 

 

def __len__(self): 

""" int : number of clusters (or force constants) """ 

return len(self._fc_dict) 

 

def __add__(self, fcs_other): 

""" ForceConstants : addition of two force constants """ 

53 ↛ 54line 53 didn't jump to line 54, because the condition on line 53 was never true if type(self) != type(fcs_other): 

raise ValueError('ForceConstants objects are of different types') 

 

# if they have overlapping orders, raise 

57 ↛ 58line 57 didn't jump to line 58, because the condition on line 57 was never true if any(order in self.orders for order in fcs_other.orders): 

raise ValueError('ForceConstants objects share at least one order') 

 

# check that they share the same supercell 

61 ↛ 62line 61 didn't jump to line 62, because the condition on line 61 was never true if self.supercell != fcs_other.supercell: 

raise ValueError('supercells do not match') 

 

# add force constants 

fc_merged = {**self._fc_dict, **fcs_other._fc_dict} 

return self.__class__(fc_merged, self.supercell) 

 

def get_fc_dict(self, order=None): 

""" Returns force constant dictionary for one specific order. 

 

The returned dict may be sparse or may be dense depending on the 

underlying force constants. 

 

Parameters 

---------- 

order : int 

force constants returned for this order 

 

Returns 

------- 

dict 

dictionary with keys corresponding to clusters and values to 

respective force constant 

""" 

if order is None: 

return self._fc_dict 

 

if order not in self.orders: 

raise ValueError('Order {} not in ForceConstants'.format(order)) 

 

fc_order = {} 

for c, fc in self._fc_dict.items(): 

if len(c) == order: 

fc_order[c] = fc 

return fc_order 

 

def get_fc_array(self, order, format='phonopy'): 

""" Returns force constants in array format for specified order. 

 

Parameters 

---------- 

order : int 

force constants for this order will be returned 

format : str 

specify which format (shape) the NumPy array should have, 

possible values are `phonopy` and `ase` 

 

Returns 

------- 

numpy.ndarray 

NumPy array with shape `(N,)*order + (3,)*order` where `N` is 

the number of atoms 

""" 

if order not in self.orders: 

raise ValueError('Order not in orders') 

116 ↛ 117line 116 didn't jump to line 117, because the condition on line 116 was never true if format not in ['ase', 'phonopy']: 

raise ValueError('Format must be either ase or phonopy') 

 

# generate fc array for phonopy format 

fc_array = np.zeros((self.n_atoms, ) * order + (3, ) * order) 

for cluster in product(range(self.n_atoms), repeat=order): 

fc_array[cluster] = self[cluster] 

 

if format == 'ase': 

125 ↛ 126line 125 didn't jump to line 126, because the condition on line 125 was never true if order != 2: 

raise ValueError('ASE format works only for order 2') 

return fc_array.transpose([0, 2, 1, 3]).reshape( 

self.n_atoms * 3, self.n_atoms * 3) 

else: 

return fc_array 

 

def compute_gamma_frequencies(self): 

""" Computes gamma frequencies using the second-order force constants. 

 

Returns 

------- 

gamma_frequencies : numpy.ndarray 

Gamma frequencies in THz 

""" 

fc2_array = self.get_fc_array(order=2) 

masses = self.supercell.get_masses() 

return _compute_gamma_frequencies(fc2_array, masses) 

 

def assert_acoustic_sum_rules(self, order=None, tol=1e-6): 

""" Asserts that acoustic sum rules are enforced for force constants. 

 

Parameters 

---------- 

order : int 

specifies which order to check, if None all are checked 

tol : float 

numeric tolerance for checking sum rules 

 

Raises 

------ 

AssertionError 

if acoustic sum rules are not enforced 

""" 

 

# set up orders 

if order is None: 

orders = self.orders 

else: 

164 ↛ 165line 164 didn't jump to line 165, because the condition on line 164 was never true if order not in self.orders: 

raise ValueError('Order not available in FCS') 

orders = [order] 

 

atomic_indices = range(self.n_atoms) 

for order in orders: 

assert_msg = 'Acoustic sum rule order {} not enforced for atom' 

assert_msg += ' {}' * (order - 1) + ' x' 

for ijk in product(atomic_indices, repeat=order-1): 

fc_sum_ijk = np.zeros((3, )*order) 

for l in atomic_indices: 

cluster = ijk + (l, ) 

fc_sum_ijk += self[cluster] 

assert np.linalg.norm(fc_sum_ijk) < tol, \ 

assert_msg.format(order, *ijk) 

 

def print_force_constant(self, cluster): 

""" 

Prints force constants for a cluster in a nice format. 

 

Parameters 

---------- 

cluster : tuple(int) 

sites belonging to the cluster 

""" 

print(self._repr_fc(cluster)) 

 

def __eq__(self, other): 

 

# check supercells are the same 

if not len(self.supercell) == len(other.supercell): 

return False 

196 ↛ 197line 196 didn't jump to line 197, because the condition on line 196 was never true if not np.allclose(self.supercell.positions, other.supercell.positions): 

return False 

198 ↛ 199line 198 didn't jump to line 199, because the condition on line 198 was never true if not np.allclose(self.supercell.cell, other.supercell.cell): 

return False 

200 ↛ 201line 200 didn't jump to line 201, because the condition on line 200 was never true if not all(self.supercell.numbers == other.supercell.numbers): 

return False 

 

# check orders and clusters are the same 

if self.orders != other.orders: 

return False 

206 ↛ 207line 206 didn't jump to line 207, because the condition on line 206 was never true if self.clusters != other.clusters: 

return False 

 

# check force constants 

for c in self.clusters: 

if not np.allclose(self[c], other[c]): 

return False 

return True 

 

def __str__(self): 

s = [] 

s.append(' ForceConstants '.center(54, '=')) 

s.append('Orders: {}'.format(self.orders)) 

s.append('Atoms: {}'.format(self.supercell)) 

s.append('') 

221 ↛ 228line 221 didn't jump to line 228, because the condition on line 221 was never false if len(self) > 10: 

for c in self.clusters[:3]: 

s.append(self._repr_fc(c)+'\n') 

s.append('...\n') 

for c in self.clusters[-3:]: 

s.append(self._repr_fc(c)+'\n') 

else: 

for c in self.clusters: 

s.append(self._repr_fc(c)+'\n') 

return '\n'.join(s) 

 

def __repr__(self): 

fc_dict_str = '{{{}: {}, ..., {}: {}}}'.format( 

self.clusters[0], self[self.clusters[0]].round(5), 

self.clusters[-1], self[self.clusters[-1]].round(5)) 

return ('ForceConstants(fc_dict={}, atoms={!r})' 

.format(fc_dict_str, self.supercell)) 

 

def _repr_fc(self, cluster, precision=5, suppress=True): 

""" 

Representation for single cluster and its force constant. 

 

Parameters 

---------- 

cluster : tuple 

tuple of ints indicating the sites belonging to the cluster 

""" 

s = [] 

s.append('Cluster: {}'.format(cluster)) 

for atom_index in cluster: 

s.append(str(self.supercell[atom_index])) 

s.append('Force constant:') 

s.append(np.array_str(self[cluster], precision=precision, 

suppress_small=suppress)) 

return '\n'.join(s) 

 

def _sanity_check_dict(self, fc_dict): 

""" Checks that all indices in clusters are between 0 and number of 

atoms. """ 

for cluster in fc_dict.keys(): 

if not all(0 <= i < self.n_atoms for i in cluster): 

raise ValueError('Cluster {} not in supercell'.format(cluster)) 

 

@classmethod 

def from_arrays(cls, supercell, fc2_array=None, fc3_array=None): 

""" Constructs FCs from numpy arrays. 

 

One or both of fc2_array and fc3_array must not be None 

 

Parameters 

---------- 

supercell : ase.Atoms 

supercell structure 

fc2_array : np.ndarray 

second-order force constant in phonopy format, i.e. must have shape (N, N, 3, 3) 

fc3_array : np.ndarray 

third-order force constant in phonopy format, i.e. must have shape (N, N, N, 3, 3, 3) 

""" 

279 ↛ 280line 279 didn't jump to line 280, because the condition on line 279 was never true if fc2_array is None and fc3_array is None: 

raise ValueError('Please provide force constant arrays.') 

 

n_atoms = len(supercell) 

if fc2_array is None: 

fc2_dict = dict() 

else: 

286 ↛ 287line 286 didn't jump to line 287, because the condition on line 286 was never true if fc2_array.shape != (n_atoms, n_atoms, 3, 3): 

raise ValueError('fc2 array has wrong shape') 

fc2_dict = array_to_dense_dict(fc2_array) 

 

if fc3_array is None: 

fc3_dict = dict() 

else: 

293 ↛ 294line 293 didn't jump to line 294, because the condition on line 293 was never true if fc3_array.shape != (n_atoms, n_atoms, n_atoms, 3, 3, 3): 

raise ValueError('fc2 array has wrong shape') 

fc3_dict = array_to_dense_dict(fc3_array) 

 

fc_dict = {**fc2_dict, **fc3_dict} 

return cls.from_dense_dict(fc_dict, supercell) 

 

@classmethod 

def from_sparse_dict(cls, fc_dict, supercell): 

""" Assumes label symmetries, meaning only one cluster for each 

permuation should be included 

 

Parameters 

---------- 

fc_dict : dict 

keys corresponding to clusters and values to the force constants 

supercell : ase.Atoms 

""" 

return SortedForceConstants(fc_dict, supercell=supercell) 

 

@classmethod 

def from_dense_dict(cls, fc_dict, supercell): 

""" All permutations of clusters that are not zero must be listed, 

if label symmetries are fullfilled will return a SortedForceConstants 

 

Parameters 

---------- 

fc_dict : dict 

keys corresponding to clusters and values to the force constants 

supercell : ase.Atoms 

""" 

if check_label_symmetries(fc_dict): 

fc_sparse_dict = dense_dict_to_sparse_dict(fc_dict) 

return SortedForceConstants(fc_sparse_dict, supercell=supercell) 

else: 

return RawForceConstants(fc_dict, supercell=supercell) 

 

@classmethod 

def read_phonopy(cls, supercell, fname, format=None): 

""" Reads force constants from a phonopy calculation. 

 

Parameters 

---------- 

supercell : ase.Atoms 

supercell structure (`SPOSCAR`) 

fname : str 

name of second-order force constant file 

format : str 

format for second-order force constants 

""" 

fc2_array = phonopyIO.read_phonopy_fc2(fname, format=format) 

return cls.from_arrays(supercell, fc2_array) 

 

@classmethod 

def read_phono3py(cls, supercell, fname): 

""" Reads force constants from a phono3py calculation. 

 

Parameters 

---------- 

supercell : ase.Atoms 

supercell structure (`SPOSCAR`) 

fname : str 

name of third-order force constant file 

""" 

fc3_array = phonopyIO.read_phonopy_fc3(fname) 

return cls.from_arrays(supercell, fc3_array=fc3_array) 

 

@classmethod 

def read_shengBTE(cls, supercell, fname, prim): 

""" Reads third order force constants from a shengBTE calculation. 

 

shengBTE force constants will be mapped onto a supercell. 

 

Parameters 

---------- 

supercell : str 

supercell structure 

fname : str 

name of third-order force constant file 

prim : ase.Atoms 

primitive configuration (must be equivalent to structure used in 

the shengBTE calculation) 

""" 

fcs = shengBTEIO.read_shengBTE_fc3(fname, prim, supercell) 

return fcs 

 

@classmethod 

def read(cls, fname): 

""" Reads ForceConstants from file. 

 

Parameters 

---------- 

fname : str 

name of file from which to read 

""" 

tar_file = tarfile.open(mode='r', name=fname) 

items = read_items_pickle(tar_file, 'fc_dict') 

fc_dict = items['fc_dict'] 

fcs_type = items['fcs_type'] 

supercell = read_ase_atoms_from_tarfile(tar_file, 'supercell') 

tar_file.close() 

394 ↛ 396line 394 didn't jump to line 396, because the condition on line 394 was never false if fcs_type == 'SortedForceConstants': 

return SortedForceConstants(fc_dict, supercell) 

elif fcs_type == 'RawForceConstants': 

return RawForceConstants(fc_dict, supercell) 

else: 

raise ValueError('FCS type not recongnized') 

 

def write(self, fname): 

""" Writes entire ForceConstants object to file. 

 

Parameters 

---------- 

fname : str 

name of file to which to write 

""" 

tar_file = tarfile.open(name=fname, mode='w') 

items_pickle = dict(fc_dict=self._fc_dict, fcs_type=self.__class__.__name__) 

add_items_to_tarfile_pickle(tar_file, items_pickle, 'fc_dict') 

add_ase_atoms_to_tarfile(tar_file, self.supercell, 'supercell') 

tar_file.close() 

 

def write_to_phonopy(self, fname, format=None): 

""" 

Writes force constants in phonopy format. 

 

Parameters 

---------- 

fname : str 

name of file to which to write second-order force constant 

format : str 

format for second-order force constants 

""" 

phonopyIO.write_phonopy_fc2(fname, self, format=format) 

 

def write_to_phono3py(self, fname): 

""" 

Writes force constants in phono3py format. 

 

Parameters 

---------- 

fname : str 

name of file to which to write third-order force constant 

""" 

phonopyIO.write_phonopy_fc3(fname, self) 

 

def write_to_shengBTE(self, fname, prim, **kwargs): 

""" 

Writes third order force constants in shengBTE format. 

 

Parameters 

---------- 

fname : str 

name of file to which to write third-order force constant 

prim : ase.Atoms 

primitive configuration (must be equivalent to structure used in 

the shengBTE calculation) 

""" 

shengBTEIO.write_shengBTE_fc3(fname, self, prim, **kwargs) 

 

 

class SortedForceConstants(ForceConstants): 

""" Force constants with label symmetries. 

 

Parameters 

---------- 

fc_dict : dict 

keys corresponding to clusters and values to the force constants, 

should only contain sorted clusters 

supercell : ase.Atoms 

""" 

 

def __init__(self, fc_dict, supercell): 

super().__init__(supercell) 

self._sanity_check_dict(fc_dict) 

self._fc_dict = fc_dict 

self._orders = sorted(set(len(c) for c in self._fc_dict.keys())) 

 

def __getitem__(self, cluster): 

sorted_cluster = tuple(sorted(cluster)) 

 

# cluster not in fcs 

if sorted_cluster not in self._fc_dict.keys(): 

return np.zeros((3,)*len(cluster)) 

 

# return fc for the unsorted cluster 

inv_perm = np.argsort(np.argsort(cluster)) 

return self._fc_dict[sorted_cluster].transpose(inv_perm) 

 

@property 

def orders(self): 

""" list : orders for which force constants exist """ 

return self._orders.copy() 

 

def _sanity_check_dict(self, fc_dict): 

super()._sanity_check_dict(fc_dict) 

 

# also check clusters are sorted 

for cluster in fc_dict.keys(): 

if cluster != tuple(sorted(cluster)): 

raise ValueError('Found unsorted cluster {}'.format(cluster)) 

 

 

class RawForceConstants(ForceConstants): 

""" Force constants without label symmetries. 

 

Parameters 

---------- 

fc_dict : dict 

keys corresponding to clusters and values to the force constants, 

should contain all clusters with nonzero force constants 

supercell : ase.Atoms 

""" 

 

def __init__(self, fc_dict, supercell): 

super().__init__(supercell) 

self._sanity_check_dict(fc_dict) 

self._fc_dict = fc_dict 

self._orders = sorted(set(len(c) for c in self._fc_dict.keys())) 

 

def __getitem__(self, cluster): 

514 ↛ 515line 514 didn't jump to line 515, because the condition on line 514 was never true if cluster not in self._fc_dict.keys(): 

return np.zeros((3,)*len(cluster)) 

return self._fc_dict[cluster] 

 

@property 

def orders(self): 

""" list : orders for which force constants exist """ 

return self._orders.copy() 

 

 

# ====================== # 

# HELPER FUNCTIONS # 

# ====================== # 

 

 

def array_to_dense_dict(fc_array, fc_tol=1e-10): 

""" Constructs a dense dict from an fc array in phonopy format. 

 

Force constants with norm smaller than fc_tol will be considered zero and 

therefore not included in the fc_dict. 

 

Parameters 

---------- 

fc_array : np.ndarray 

force constant array in phonopy format 

fc_tol : float 

tolerance for considering force constants zero or not 

""" 

 

# sanity check 

n_atoms = fc_array.shape[0] 

order = int(len(fc_array.shape) / 2) 

546 ↛ 547line 546 didn't jump to line 547, because the condition on line 546 was never true if fc_array.shape != (n_atoms, ) * order + (3, ) * order: 

raise ValueError('fc array has bad shape') 

 

# construct dense dict 

fc_dict = dict() 

for cluster in product(range(n_atoms), repeat=order): 

fc = fc_array[cluster] 

if np.linalg.norm(fc) > fc_tol: 

fc_dict[cluster] = fc 

return fc_dict 

 

 

def check_label_symmetries(fc_dict): 

""" Checks label symmetries for dense fc dict. 

 

TODO 

---- 

tol, which one to use etc 

 

Parameters 

---------- 

fc_dict : dict 

keys corresponding to clusters and values to the force constants 

""" 

for cluster, fc in fc_dict.items(): 

inv_perm = np.argsort(np.argsort(cluster)) 

sorted_cluster = tuple(sorted(cluster)) 

573 ↛ 574line 573 didn't jump to line 574, because the condition on line 573 was never true if sorted_cluster not in fc_dict: 

return False 

if not np.allclose(fc, fc_dict[sorted_cluster].transpose(inv_perm)): 

return False 

return True 

 

 

def dense_dict_to_sparse_dict(fc_dict): 

""" Converts dense dict to sparse dict. 

 

This does not check if label symmetry is True, but rather will just keep 

the sorted clusters and their force constants. 

 

Parameters 

---------- 

fc_dict : dict 

keys corresponding to clusters and values to the force constants 

""" 

fc_dict_sparse = dict() 

for cluster, fc in fc_dict.items(): 

if cluster == tuple(sorted(cluster)): 

fc_dict_sparse[cluster] = fc 

return fc_dict_sparse 

 

 

def symbolize_force_constant(fc, tol=1e-10): 

"""Carries out a symbolic symmetrization of a force constant tensor. 

 

Parameters 

---------- 

fc : numpy.ndarray 

force constant tensor 

tol : float 

tolerance used to decide whether two elements are identical 

 

Returns 

------- 

numpy.ndarray 

symbolic representation of force constant matrix 

""" 

fc_int = np.round(fc / tol).astype(int) 

fc_chars = np.empty(fc_int.shape, dtype=object).flatten() 

all_chars = ascii_lowercase + ascii_uppercase 

lookup_chars = {} 

for i, val in enumerate(fc_int.flatten()): 

if val == 0: 

fc_chars[i] = 0 

elif val in lookup_chars.keys(): 

fc_chars[i] = lookup_chars[val] 

elif -val in lookup_chars.keys(): 

fc_chars[i] = '-{:}'.format(lookup_chars[-val]) 

else: 

lookup_chars[val] = all_chars[len(lookup_chars.keys())] 

fc_chars[i] = lookup_chars[val] 

return fc_chars.reshape(fc_int.shape) 

 

 

def _compute_gamma_frequencies(fc2, masses): 

""" Computes Gamma frequencies using second-order force constants 

 

Assumes fc2 is in units of eV/A2. 

 

Parameters 

---------- 

fc2 : numpy.ndarray 

second-order force constants in phonopy format 

masses : list(float) 

mass of each atom 

 

Returns 

------- 

gamma_frequencies : numpy.ndarray 

Gamma frequencies in THz 

""" 

 

n_atoms = fc2.shape[0] 

649 ↛ 650line 649 didn't jump to line 650, because the condition on line 649 was never true if len(masses) != n_atoms: 

raise ValueError('Length of masses not compatible with fc2') 

mass_matrix = np.sqrt(np.outer(masses, masses)) 

 

# divide with mass matrix 

fc2_tmp = np.zeros((n_atoms, n_atoms, 3, 3)) 

for pair in product(range(n_atoms), repeat=2): 

fc2_tmp[pair] = fc2[pair] / mass_matrix[pair] 

 

# reshape into matrix and solve eigenvalues 

fc2_tmp = fc2_tmp.transpose([0, 2, 1, 3]).reshape(n_atoms * 3, n_atoms * 3) 

eigen_vals, _ = eig(fc2_tmp) 

eigen_vals *= 1e20 / units.J * units.kg # [1/s**2] 

eigen_vals.sort() 

 

# if negative eigenval, set frequency to negative 

gamma_frequencies = [] 

for val in eigen_vals.real: 

if val >= 0: 

gamma_frequencies.append(np.sqrt(val)) 

else: 

gamma_frequencies.append(-np.sqrt(np.abs(val))) 

 

# Sort and convert to THz 

gamma_frequencies = np.array(gamma_frequencies) / 1e12 / (2 * np.pi) 

gamma_frequencies.sort() 

return gamma_frequencies