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

import itertools 

import numpy as np 

import spglib as spg 

from . import atoms as atoms_module 

from ..io.logging import logger 

 

logger = logger.getChild('relate_structures') 

 

 

def relate_structures(reference, target): 

"""Finds rotation and translation operations that align two structures with 

periodic boundary conditions. 

 

The rotation and translation in Cartesian coordinates will map the 

reference structure onto the target 

 

Parameters 

---------- 

reference : ASE Atoms 

The reference structure to be mapped 

target : ASE Atoms 

The target structure 

 

Returns 

------- 

R : NumPy (3, 3) array 

The rotation matrix in Cartesian coordinates 

T : NumPy (3) array 

The translation vector in Cartesian coordinates 

""" 

 

logger.debug('Reference atoms:') 

_debug_log_atoms(reference) 

 

reference_primitive_cell = get_primitive_cell(reference) 

 

logger.debug('Reference primitive cell') 

_debug_log_atoms(reference_primitive_cell) 

 

logger.debug('Target atoms:') 

_debug_log_atoms(target) 

 

target_primitive_cell = get_primitive_cell(target) 

 

logger.debug('Target primitive cell') 

_debug_log_atoms(target_primitive_cell) 

 

logger.debug('Sane check that primitive cells can match...') 

_assert_numbers_match(reference_primitive_cell.numbers, 

target_primitive_cell.numbers) 

_assert_volume_match(reference_primitive_cell.cell, 

target_primitive_cell.cell) 

 

logger.debug('Finding rotations...') 

rotations = _find_rotations(reference_primitive_cell.cell, 

target_primitive_cell.cell) 

 

logger.debug('Finding transformations...') 

59 ↛ 67line 59 didn't jump to line 67, because the loop on line 59 didn't complete for R in rotations: 

rotated_reference_primitive_cell = \ 

rotate_atoms(reference_primitive_cell, R) 

T = _find_translation(rotated_reference_primitive_cell, 

target_primitive_cell) 

if T is not None: 

break 

else: 

raise Exception(('Found no translation!\n' 

'Reference primitive cell basis:\n' 

'{}\n' 

'Target primitive cell basis:\n' 

'{}') 

.format(reference_primitive_cell.basis, 

target_primitive_cell.basis)) 

 

logger.debug(('Found rotation\n' 

'{}\n' 

'and translation\n' 

'{}') 

.format(R, T)) 

 

return R, T 

 

 

def is_rotation(R, cell_metric=None): 

"""Checks if rotation matrix is orthonormal 

 

A cell metric can be passed of the rotation matrix is in scaled coordinates 

 

Parameters 

---------- 

R : Numpy 3x3 matrix 

The rotation matrix 

cell_metric: Numpy 3x3 matrix 

Optinal cell metric if the rotation is in scaled coordinates 

""" 

96 ↛ 99line 96 didn't jump to line 99, because the condition on line 96 was never false if not cell_metric: 

cell_metric = np.eye(3) 

 

V = cell_metric 

V_inv = np.linalg.inv(V) 

lhs = np.linalg.multi_dot([V_inv, R.T, V, V.T, R, V_inv.T]) 

 

return np.allclose(lhs, np.eye(3), atol=1e-4) # TODO: tol 

 

 

def _find_rotations(reference_cell_metric, target_cell_metric): 

"""Generate all proper and improper rotations aligning two cell metrics""" 

 

rotations = [] 

V1 = reference_cell_metric 

for perm in itertools.permutations([0, 1, 2]): 

# Make sure the improper rotations are included 

for inv in itertools.product([1, -1], repeat=3): 

V2 = np.diag(inv) @ target_cell_metric[perm, :] 

R = np.linalg.solve(V1, V2).T 

# Make sure the rotation is orthonormal 

if is_rotation(R): 

for R_tmp in rotations: 

119 ↛ 120line 119 didn't jump to line 120, because the condition on line 119 was never true if np.allclose(R, R_tmp): # TODO: tol 

break 

else: 

rotations.append(R) 

 

assert rotations, ('Found no rotations! Reference cell metric:\n' 

'{}\n' 

'Target cell metric:\n' 

'{}').format(reference_cell_metric, target_cell_metric) 

 

logger.debug('Found {} rotations'.format(len(rotations))) 

 

return rotations 

 

 

def _assert_numbers_match(reference_numbers, target_numbers): 

""" Asserts that two configurations contain the same number and types of 

atoms. """ 

reference_numbers = sorted(reference_numbers) 

target_numbers = sorted(target_numbers) 

assert_text = ('Atom species do not match.\n' 

'Reference: {}\n' 

'Target: {}').format(reference_numbers, target_numbers) 

 

assert reference_numbers == target_numbers, assert_text 

 

 

def _assert_volume_match(reference_cell, target_cell): 

"""Assert cell sizes are equal""" 

reference_volume = np.linalg.det(reference_cell) 

target_volume = np.linalg.det(target_cell) 

assert_text = ('Cell sizes do not match. Check lattice constants!\n' 

'Reference volume: {}, primitive cell metric:\n{}\n' 

'Target colume: {}, primitive cell metric:\n{}' 

.format(reference_volume, reference_cell, 

target_volume, target_cell)) 

# TODO: tol 

assert np.isclose(reference_volume, target_volume), assert_text 

 

 

# TODO: tol 

def get_primitive_cell(atoms, to_primitive=True, no_idealize=True): 

"""Get primitive cell from spglib 

 

Parameters 

---------- 

atoms: ASE Atoms object 

to_primitive: bool 

passed to spglib 

no_idealize: bool 

passed to spglib 

""" 

171 ↛ 172line 171 didn't jump to line 172, because the condition on line 171 was never true if not all(atoms.pbc): 

raise ValueError('atoms must have pbc.') 

spg_primitive_cell = spg.standardize_cell(atoms, to_primitive=True, 

no_idealize=True) 

primitive_cell = atoms_module.Atoms(cell=spg_primitive_cell[0], 

scaled_positions=spg_primitive_cell[1], 

numbers=spg_primitive_cell[2], 

pbc=True) 

return primitive_cell 

 

 

def _debug_log_atoms(atoms): 

logger.debug('cell:\n{}'.format(atoms.cell)) 

logger.debug('spos:\n{}'.format(atoms.get_scaled_positions())) 

logger.debug('pos:\n{}'.format(atoms.positions)) 

logger.debug('numbers:\n{}'.format(atoms.numbers)) 

 

 

def rotate_atoms(atoms, rotation): 

"""Rotates the cell and positions of Atoms and returns a copy 

 

Parameters 

---------- 

atoms : ASE Atoms object 

rotation: Numpy 3x3 matrix 

""" 

 

cell = np.dot(rotation, atoms.cell.T).T 

positions = np.dot(rotation, atoms.positions.T).T 

return atoms_module.Atoms(cell=cell, positions=positions, 

numbers=atoms.numbers, pbc=atoms.pbc) 

 

 

def _find_translation(reference, target): 

"""Given two compatible atoms returns the translation 

 

I.e. the two atoms obj must describe the same structure when infinitely 

repeated but differ by a translation 

""" 

 

atoms = atoms_module.Atoms(cell=target.cell, 

positions=reference.positions, 

numbers=reference.numbers, 

pbc=True) 

atoms.wrap() 

 

atoms_atom_0 = atoms[0] 

for atom in target: 

if atoms_atom_0.symbol != atom.symbol: 

continue 

T = atom.position - atoms_atom_0.position 

atoms_copy = atoms.copy() 

atoms_copy.positions += T 

if atoms_equal(atoms_copy, target): 

return T 

return None 

 

 

def atoms_equal(atoms_1, atoms_2): 

"""Compares two configurations with periodic boundary conditions. 

 

The atoms should have the same 

* cell 

* elements 

* scaled positions (mod 1) 

* pbc 

 

Compared to Atoms.__eq__ the order of the atoms does not matter 

""" 

nAtoms = len(atoms_1) 

# TODO: tol 

242 ↛ 246line 242 didn't jump to line 246, because the condition on line 242 was never true if not (np.allclose(atoms_1.cell, atoms_2.cell, atol=1e-4) and 

nAtoms == len(atoms_2) and 

sorted(atoms_1.numbers) == sorted(atoms_2.numbers) and 

all(atoms_1.pbc == atoms_2.pbc)): 

return False 

new_cell = (atoms_1.cell + atoms_2.cell) / 2 

pos = [a.position for a in atoms_1] + [a.position for a in atoms_2] 

num = [a.number for a in atoms_1] + [a.number for a in atoms_2] 

s3 = atoms_module.Atoms(cell=new_cell, positions=pos, numbers=num, 

pbc=True) 

for i in range(nAtoms): 

for j in range(nAtoms, len(s3)): 

d = s3.get_distance(i, j, mic=True) 

if abs(d) < 1e-4: # TODO: tol 

256 ↛ 257line 256 didn't jump to line 257, because the condition on line 256 was never true if s3[i].number != s3[j].number: 

return False 

break 

else: 

return False 

return True