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

import logging 

from typing import List, Tuple 

 

import numpy as np 

from ase import Atoms 

from ase.build import cut 

 

from icet.io.logging import logger 

logger = logger.getChild('structure_mapping') 

 

 

def map_structure_to_reference(input_structure: Atoms, 

reference_structure: Atoms, 

tolerance_mapping: float, 

vacancy_type: str = None, 

inert_species: List[str] = None, 

tolerance_cell: float = 0.05, 

tolerance_positions: float = 0.01) \ 

-> Tuple[Atoms, float, float]: 

"""Maps a relaxed structure onto a reference structure. 

The function returns a tuple comprising 

 

* the ideal supercell most closely matching the input structure, 

* the largest deviation of any input coordinate from its ideal 

coordinate, and 

* the average deviation of the input coordinates from the ideal 

coordinates. 

 

Parameters 

---------- 

input_structure 

relaxed input structure 

reference_structure 

reference structure, which can but need not represent the primitive 

cell 

tolerance_mapping 

maximum allowed displacement for mapping an atom in the relaxed (but 

rescaled) structure to the reference supercell 

 

*Note*: A reasonable choice is up to 20-30% of the first 

nearest neighbor distance (`r1`). A value above 50% of `r1` 

will most likely lead to atoms being multiply assigned, 

whereby the mapping fails. 

vacancy_type 

If this parameter is set to a non-zero string unassigned sites in the 

reference structure will be assigned to this type. 

 

*Note 1*: By default (``None``) the method will fail if there 

are *any* unassigned sites in the reference structure. 

 

*Note 2*: ``vacancy_type`` must be a valid species as 

enforced by the :class:`ase.Atoms` class. 

inert_species 

List of chemical symbols (e.g., ``['Au', 'Pd']``) that are never 

substituted for a vacancy. Used to make an initial rescale of the cell 

and thus increases the probability for a successful mapping. Need not 

be specified if ``vacancy_type`` is ``None``. 

tolerance_cell 

tolerance factor applied when computing permutation matrix to generate 

supercell 

tolerance_positions 

tolerance factor applied when scanning for overlapping positions in 

Angstrom (forwarded to :func:`ase.build.cut`) 

 

Example 

------- 

The following code snippet illustrates the general usage. It first creates 

a primitive FCC cell, which is latter used as reference structure. To 

emulate a relaxed structure obtained from, e.g., a density functional 

theory calculation, the code then creates a 4x4x4 conventional FCC 

supercell, which is populated with two different atom types, has distorted 

cell vectors, and random displacements to the atoms. Finally, the present 

function is used to map the structure back the ideal lattice:: 

 

from ase.build import bulk 

reference = bulk('Au', a=4.09) 

structure = bulk('Au', cubic=True, a=4.09).repeat(4) 

structure.set_chemical_symbols(10 * ['Ag'] + (len(structure) - 10) * ['Au']) 

structure.set_cell(structure.cell * 1.02, scale_atoms=True) 

structure.rattle(0.1) 

mapped_structure = map_structure_to_reference(structure, reference, 1.0) 

 

""" 

assert np.all(input_structure.pbc == reference_structure.pbc), \ 

('The periodic boundary conditions differ' 

' between input and reference structure') 

 

88 ↛ 92line 88 didn't jump to line 92, because the condition on line 88 was never false if logger.isEnabledFor(logging.DEBUG): 

np.set_printoptions(suppress=True, precision=6) 

 

# Scale input cell and construct supercell of the reference structure 

scaled_cell = _get_scaled_cell(input_structure, reference_structure, 

vacancy_type=vacancy_type, 

inert_species=inert_species) 

P = _get_transformation_matrix(scaled_cell, 

reference_structure.cell, 

tolerance_cell=tolerance_cell) 

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

scaled_structure, ideal_supercell = \ 

_rescale_structures(input_structure, 

reference_structure, 

P, 

tolerance_positions=tolerance_positions) 

 

assert len(ideal_supercell) == len(scaled_structure) or \ 

vacancy_type is not None, \ 

('Number of atoms in ideal supercell does not match ' 

'input structure.\n' 

'ideal: {}\ninput: {}'.format(len(ideal_supercell), 

len(scaled_structure))) 

 

logger.debug('Number of atoms in reference structure:' 

' {}'.format(len(reference_structure))) 

logger.debug('Number of atoms in input structure:' 

' {}\n'.format(len(input_structure))) 

logger.debug('Reference cell metric:\n' 

'{}'.format(reference_structure.cell)) 

logger.debug('Input cell metric:\n' 

'{}\n'.format(input_structure.cell)) 

logger.debug('Transformation matrix connecting reference structure' 

' and idealized input structure:\n {}'.format(P)) 

logger.debug('Determinant of transformation matrix:' 

' {:.3f}\n'.format(np.linalg.det(P))) 

logger.debug('Cell metric of ideal supercell:\n' 

'{}'.format(ideal_supercell.cell)) 

logger.debug('Cell metric of rescaled input structure:\n' 

'{}\n'.format(scaled_structure.cell)) 

 

# per-atom-list for keeping track of mapped atoms 

mapped = [-1] * len(ideal_supercell) 

# positions in scaled input structure 

original_positions = [(None, None, None)] * len(ideal_supercell) 

# displacement vectors 

displacements = [(None, None, None)] * len(ideal_supercell) 

# distances between ideal and input sites 

drs = [None] * len(ideal_supercell) 

for ideal_site in ideal_supercell: 

for atom in scaled_structure: 

# in order to compute the distance the current atom from 

# the input structure is temporarily added to the 

# ideal supercell. This allows one to simply use the ASE 

# Atoms method for computing the interatomic distance 

ideal_supercell.append(atom) 

dvec = ideal_supercell.get_distance(ideal_site.index, 

ideal_supercell[-1].index, 

mic=True, vector=True) 

dr = np.linalg.norm(dvec) 

del ideal_supercell[-1] 

if dr < tolerance_mapping: 

if mapped[ideal_site.index] >= 0: 

raise Exception('More than one atom from the relaxed' 

' (and rescaled) structure have been' 

' mapped onto the same ideal site.\n' 

' Try reducing `tolerance_mapping`.') 

mapped[ideal_site.index] = atom.index 

drs[ideal_site.index] = dr 

original_positions[ideal_site.index] = atom.position 

displacements[ideal_site.index] = dvec 

ideal_site.symbol = atom.symbol 

break 

else: 

assert vacancy_type is not None, \ 

('Failed to assign an atom from the relaxed (and' 

' rescaled) structure to the ideal lattice.' 

' Try increasing `tolerance_mapping`.\n' 

' {}'.format(ideal_site)) 

ideal_site.symbol = vacancy_type 

 

# store information concerning mapping in output structure 

ideal_supercell.new_array('displacements', displacements, 

float, (3,)) 

ideal_supercell.new_array('original_positions', original_positions, 

float, (3,)) 

for name, ar in reference_structure.arrays.items(): 

175 ↛ 177line 175 didn't jump to line 177, because the condition on line 175 was never false if name in ideal_supercell.arrays: 

continue 

ideal_supercell.new_array(name, ar) 

 

# maximum, average, and standard deviation of displacements 

dr_avg = np.average([d for d in drs if d is not None]) 

dr_sdv = np.std([d for d in drs if d is not None]) 

dr_max = np.max([d for d in drs if d is not None]) 

 

# check that not more than one atom was assigned to the same site 

for k in set(mapped): 

assert k < 0 or mapped.count(k) <= 1, \ 

('Site {} has been assigned more than once.'.format(k)) 

 

# check that the chemical composition of input and ideal supercell matches 

for symbol in set(input_structure.get_chemical_symbols()): 

n1 = input_structure.get_chemical_symbols().count(symbol) 

n2 = ideal_supercell.get_chemical_symbols().count(symbol) 

assert n1 == n2, ('Number of atoms of type {} differs between' 

' input structure ({}) and ideal' 

' supercell ({}).'.format(symbol, n1, n2)) 

 

197 ↛ 227line 197 didn't jump to line 227, because the condition on line 197 was never false if logger.isEnabledFor(logging.DEBUG): 

logger.debug('Maximum, average and standard deviation of atomic' 

' displacements: {} {} {}'.format(dr_max, dr_avg, dr_sdv)) 

 

logger.debug('{:52} {:}'.format('Input structure:', 

'Scaled structure:')) 

for k, (input_atom, scaled_atom) in enumerate(zip(input_structure, 

scaled_structure)): 

msg = '{:4} {:2}'.format(k, input_atom.symbol) 

msg += (3 * ' {:12.6f}').format(*input_atom.position) 

msg += ' -->' 

msg += (3 * ' {:12.6f}').format(*scaled_atom.position) 

logger.debug(msg) 

 

logger.debug('\n') 

 

logger.debug('{:52} {}'.format('Ideal supercell:', 

'Scaled structure:')) 

for ideal_atom, k, dr in zip(ideal_supercell, mapped, drs): 

msg = ' {:2}'.format(ideal_atom.symbol) 

msg += (3 * ' {:12.6f}').format(*ideal_atom.position) 

msg += ' -->' 

msg += ' {:4}'.format(k) 

 

if k >= 0: 

scaled_pos = scaled_structure[k].position 

msg += (3 * ' {:12.6f}').format(*scaled_pos) 

msg += ' --> {:.4f}'.format(dr) 

logger.debug(msg) 

 

return ideal_supercell, dr_max, dr_avg 

 

 

def _get_scaled_cell(input_structure: Atoms, 

reference_structure: Atoms, 

vacancy_type: str = None, 

inert_species: List[str] = None) -> np.ndarray: 

""" 

The input structure needs to be scaled in order to match the lattice 

structure of the reference structure. The reference structure can be a 

primitive cell, in which case the input structure would usually be a 

supercell thereof. Also, we need an ideal supercell that matches the input 

structure. 

 

Parameters 

---------- 

input_structure 

relaxed input structure 

reference_structure: ASE Atoms object 

reference structure, which can but need not represent the primitive 

cell 

vacancy_type 

if not None, the cell is scaled if and only if `inert_species` is not 

None 

inert_species 

list of chemical symbols(e.g., `['Au', 'Pd']`) that are never 

substituted for a vacancy. Needless if `vacancy_type` is `None` 

""" 

modcell = input_structure.get_cell() 

256 ↛ 258line 256 didn't jump to line 258, because the condition on line 256 was never true if vacancy_type is None: 

# Without scale factor we can just rescale with number of atoms 

atvol_in = input_structure.get_volume() / len(input_structure) 

atvol_ref = reference_structure.get_volume() / len(reference_structure) 

scale = atvol_in / atvol_ref 

 

262 ↛ 279line 262 didn't jump to line 279, because the condition on line 262 was never false if vacancy_type is not None: 

263 ↛ 264line 263 didn't jump to line 264, because the condition on line 263 was never true if inert_species is None: 

scale = 1.0 

else: 

# We cannot use the number of atoms since there may be vacancies 

# in the input_structure. Instead we count the species that we 

# know should always be present. 

n_in = 0 

n_ref = 0 

symbols_in = input_structure.get_chemical_symbols() 

symbols_ref = reference_structure.get_chemical_symbols() 

for species in inert_species: 

n_in += symbols_in.count(species) 

n_ref += symbols_ref.count(species) 

atvol_in = input_structure.get_volume() / n_in 

atvol_ref = reference_structure.get_volume() / n_ref 

scale = atvol_in / atvol_ref 

modcell *= (1.0 / scale) ** (1.0 / 3.0) 

return modcell 

 

 

def _get_transformation_matrix(input_cell: np.ndarray, 

reference_cell: np.ndarray, 

tolerance_cell: float = 0.05) -> np.ndarray: 

""" 

Obtains the (in general non-integer) transformation matrix connecting the 

input structure to the reference structure L=L_p.P - -> P=L_p ^ -1.L 

 

Parameters 

---------- 

input_cell 

cell metric of input structure(possibly scaled) 

reference_cell 

cell metric of reference structure 

tolerance_cell 

tolerance for how much the elements of P are allowed to deviate from 

the nearest integer before they are rounded 

 

Returns 

------- 

transformation matrix P of integers 

""" 

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

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

logger.debug('inv(reference_cell):\n {}' 

.format(np.linalg.inv(reference_cell))) 

P = np.dot(input_cell, np.linalg.inv(reference_cell)) 

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

 

# ensure that the transformation matrix does not deviate too 

# strongly from the nearest integer matrix 

313 ↛ 314line 313 didn't jump to line 314, because the condition on line 313 was never true if np.linalg.norm(P - np.around(P)) / 9 > tolerance_cell: 

print('reference:\n {}\n'.format(reference_cell)) 

print('input:\n {}\n'.format(input_cell)) 

print('P:\n {}\n'.format(P)) 

print('det P = {}\n'.format(np.linalg.det(P))) 

print('P_round:\n {}\n'.format(np.around(P))) 

print('Deviation: {}\n'.format(np.linalg.norm(P - np.around(P)) / 9)) 

raise Exception('Failed to map structure to reference structure' 

' (tolerance_cell exceeded). If there are vacancies,' 

' one can try specifying `inert_species`. Otherwise,' 

' one can try raising `tolerance_cell`.') 

 

# reduce the (real) transformation matrix to the nearest integer one 

P = np.around(P) 

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

return P 

 

 

def _rescale_structures(input_structure: Atoms, 

reference_structure: Atoms, 

P: np.ndarray, 

tolerance_positions: float = 0.01) \ 

-> Tuple[Atoms, Atoms]: 

""" 

Rescales `input_structure` with `P` so that it matches 

`reference_structure`, and creates a supercell of `reference_structure` 

using `P` 

 

Parameters 

---------- 

input_structure 

relaxed input structure 

reference_structure 

reference structure, which can but need not represent the primitive 

cell 

P 

transformation matrix of integers 

tolerance_positions 

tolerance factor applied when scanning for overlapping positions in 

Angstrom(forwarded to `ase.build.cut`) 

 

Returns 

------- 

a tuple with the scaled version of `input_structure` and the supercell of 

`reference_structure` matching cell metric of `scaled_structure` 

""" 

scaled_structure = input_structure.copy() 

scaled_structure.set_cell(np.dot(P, reference_structure.cell), 

scale_atoms=True) 

 

# generate supercell of (presumably primitive) reference structure 

ideal_supercell = cut(reference_structure, *P, 

tolerance=tolerance_positions) 

n_mapped = int(np.round(len(reference_structure) * np.linalg.det(P))) 

367 ↛ 368line 367 didn't jump to line 368, because the condition on line 367 was never true if len(ideal_supercell) != n_mapped: 

print('len(reference_structure)', len(reference_structure)) 

print('det(P)', np.linalg.det(P)) 

print('len(ideal_supercell)', len(ideal_supercell)) 

print('n_mapped', n_mapped) 

raise Exception('Supercell construction of reference structure failed' 

' (number of atoms do not match).\n' 

'Permutation matrix used:\n{}'.format(P) + 

'\nYou can try to change tolerance_positions.') 

 

return scaled_structure, ideal_supercell