/
structure.py
660 lines (549 loc) · 22.1 KB
/
structure.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
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
# coding: utf-8
# Copyright (c) Max-Planck-Institut für Eisenforschung GmbH - Computational Materials Design (CM) Department
# Distributed under the terms of "New BSD License", see the LICENSE file.
from __future__ import print_function
from collections import OrderedDict
import numpy as np
from pyiron_base import GenericParameters
import decimal as dec
try:
from ase.calculators.lammps import Prism
except ImportError:
try:
from ase.calculators.lammpsrun import Prism
except ImportError:
from ase.calculators.lammpsrun import prism as Prism
__author__ = "Joerg Neugebauer, Sudarsan Surendralal, Yury Lysogorskiy, Jan Janssen, Markus Tautschnig"
__copyright__ = (
"Copyright 2020, Max-Planck-Institut für Eisenforschung GmbH - "
"Computational Materials Design (CM) Department"
)
__version__ = "1.0"
__maintainer__ = "Sudarsan Surendralal"
__email__ = "surendralal@mpie.de"
__status__ = "production"
__date__ = "Sep 1, 2017"
class UnfoldingPrism(Prism):
"""
Create a lammps-style triclinic prism object from a cell
The main purpose of the prism-object is to create suitable
string representations of prism limits and atom positions
within the prism.
When creating the object, the digits parameter (default set to 10)
specify the precision to use.
lammps is picky about stuff being within semi-open intervals,
e.g. for atom positions (when using create_atom in the in-file),
x must be within [xlo, xhi).
Args:
cell:
pbc:
digits:
"""
def __init__(self, cell, pbc=(True, True, True), digits=10):
# Temporary fix. Since the arguments for the constructor have changed, try to see if it is compatible with
# the latest ase. If not, revert to the old __init__ parameters.
try:
super(UnfoldingPrism, self).__init__(
cell, pbc=pbc, tolerance=float("1e-{}".format(digits))
)
except TypeError:
super(UnfoldingPrism, self).__init__(cell, pbc=pbc, digits=digits)
a, b, c = cell
an, bn, cn = [np.linalg.norm(v) for v in cell]
alpha = np.arccos(np.dot(b, c) / (bn * cn))
beta = np.arccos(np.dot(a, c) / (an * cn))
gamma = np.arccos(np.dot(a, b) / (an * bn))
xhi = an
xyp = np.cos(gamma) * bn
yhi = np.sin(gamma) * bn
xzp = np.cos(beta) * cn
yzp = (bn * cn * np.cos(alpha) - xyp * xzp) / yhi
zhi = np.sqrt(cn ** 2 - xzp ** 2 - yzp ** 2)
# Set precision
self.car_prec = dec.Decimal("10.0") ** int(
np.floor(np.log10(max((xhi, yhi, zhi)))) - digits
)
self.dir_prec = dec.Decimal("10.0") ** (-digits)
self.acc = float(self.car_prec)
self.eps = np.finfo(xhi).eps
# For rotating positions from ase to lammps
apre = np.array(((xhi, 0, 0), (xyp, yhi, 0), (xzp, yzp, zhi)))
# np.linalg.inv(cell) ?= np.array([np.cross(b, c), np.cross(c, a), np.cross(a, b)]).T / np.linalg.det(cell)
self.R = np.dot(np.linalg.inv(cell), apre)
def fold(vec, pvec, i):
p = pvec[i]
x = vec[i] + 0.5 * p
n = (np.mod(x, p) - x) / p
return [float(self.f2qdec(vec_a)) for vec_a in (vec + n * pvec)], n
apre[1, :], n1 = fold(apre[1, :], apre[0, :], 0)
if np.abs(apre[1, 0] / apre[0, 0]) > 0.5:
apre[1, 0] -= np.sign(n1) * apre[0, 0]
n1 -= np.sign(n1)
apre[2, :], n2 = fold(apre[2, :], apre[1, :], 1)
if np.abs(apre[2, 1] / apre[1, 1]) > 0.5:
apre[2, 1] -= np.sign(n2) * apre[1, 1]
n2 -= np.sign(n2)
apre[2, :], n3 = fold(apre[2, :], apre[0, :], 0)
if np.abs(apre[2, 0] / apre[0, 0]) > 0.5:
apre[2, 0] -= np.sign(n3) * apre[0, 0]
n3 -= np.sign(n3)
self.ns = [n1, n2, n3]
d_a = apre[0, 0] / 2 - apre[1, 0]
if np.abs(d_a) < self.acc:
if d_a < 0:
print("debug: apply shift")
apre[1, 0] += 2 * d_a
apre[2, 0] += 2 * d_a
self.A = apre
if self.is_skewed() and (not (pbc[0] and pbc[1] and pbc[2])):
raise RuntimeError(
"Skewed lammps cells MUST have " "PBC == True in all directions!"
)
def unfold_cell(self, cell):
"""
Unfold LAMMPS cell to original
Let C be the pyiron cell and A be the Lammps cell, then define (in init) the rotation matrix between them as
R := C^inv.A
And recall that rotation matrices have the property
R^T == R^inv
Then left multiply the definition of R by C, and right multiply by R.T to get
C.R.R^T = C.C^inv.A.R^T
Then
C = A.R^T
After that, account for the folding process.
Args:
cell: LAMMPS cell,
Returns:
unfolded cell
"""
# Rotation
ucell = np.dot(cell, self.R.T)
# Folding
a = ucell[0]
bp = ucell[1]
cpp = ucell[2]
(n1, n2, n3) = self.ns
b = bp - n1 * a
c = cpp - n2 * bp - n3 * a
return np.array([a, b, c])
def pos_to_lammps(self, position):
"""
Rotate an ase-cell position to the lammps cell orientation
Args:
position:
Returns:
tuple of float.
"""
return tuple([x for x in np.dot(position, self.R)])
def f2qdec(self, f):
return dec.Decimal(repr(f)).quantize(self.car_prec, dec.ROUND_DOWN)
def f2s(self, f):
return str(dec.Decimal(repr(f)).quantize(self.car_prec, dec.ROUND_HALF_EVEN))
def get_lammps_prism_str(self):
"""Return a tuple of strings"""
p = self.get_lammps_prism()
return tuple([self.f2s(x) for x in p])
class LammpsStructure(GenericParameters):
"""
Args:
input_file_name:
"""
def __init__(self, input_file_name=None, bond_dict=None):
super(LammpsStructure, self).__init__(
input_file_name=input_file_name,
table_name="structure_inp",
comment_char="#",
val_only=True,
)
self._structure = None
self._potential = None
self._el_eam_lst = []
self.atom_type = None
self.cutoff_radius = None
self.digits = 10
self._bond_dict = bond_dict
@property
def potential(self):
return self._potential
@potential.setter
def potential(self, val):
self._potential = val
@property
def structure(self):
"""
Returns:
"""
return self._structure
@structure.setter
def structure(self, structure):
"""
Args:
structure:
Returns:
"""
self._structure = structure
if self.atom_type == "full":
input_str = self.structure_full()
elif self.atom_type == "bond":
input_str = self.structure_bond()
elif self.atom_type == "charge":
input_str = self.structure_charge()
else: # self.atom_type == 'atomic'
input_str = self.structure_atomic()
self.load_string(input_str)
@property
def el_eam_lst(self):
"""
Returns:
"""
return self._el_eam_lst
@el_eam_lst.setter
def el_eam_lst(self, el_eam_lst):
"""
Args:
el_eam_lst:
Returns:
"""
self._el_eam_lst = el_eam_lst
def load_default(self):
"""
Returns:
"""
input_str = ""
self.load_string(input_str)
def simulation_cell(self):
"""
Returns:
"""
self.prism = UnfoldingPrism(self._structure.cell, digits=15)
xhi, yhi, zhi, xy, xz, yz = self.prism.get_lammps_prism_str()
# Please, be carefull and not round xhi, yhi,..., otherwise you will get too skew cell from LAMMPS.
# These values are already checked in UnfoldingPrism to fullfill LAMMPS skewness criteria
simulation_cell = (
"0. {} xlo xhi\n".format(xhi)
+ "0. {} ylo yhi\n".format(yhi)
+ "0. {} zlo zhi\n".format(zhi)
)
if self.prism.is_skewed():
simulation_cell += "{0} {1} {2} xy xz yz\n".format(xy, xz, yz)
return simulation_cell
def structure_bond(self):
"""
Returns:
"""
# analyze structure to get molecule_ids, bonds, angles etc
coords = self.rotate_positions(self._structure)
elements = self._structure.get_chemical_symbols()
el_list = self._structure.get_species_symbols()
el_dict = OrderedDict()
for object_id, el in enumerate(el_list):
el_dict[el] = object_id
n_s = len(el_list)
bond_type = np.ones([n_s, n_s], dtype=np.int)
count = 0
for i in range(n_s):
for j in range(i, n_s):
count += 1
bond_type[i, j] = count
bond_type[j, i] = count
if self.structure.bonds is None:
if self.cutoff_radius is None:
bonds_lst = self.structure.get_bonds(max_shells=1)
else:
bonds_lst = self.structure.get_bonds(radius=self.cutoff_radius)
bonds = []
for ia, i_bonds in enumerate(bonds_lst):
el_i = el_dict[elements[ia]]
for el_j, b_lst in i_bonds.items():
b_type = bond_type[el_i][el_dict[el_j]]
for i_shell, ib_shell_lst in enumerate(b_lst):
for ib in np.unique(ib_shell_lst):
if ia < ib: # avoid double counting of bonds
bonds.append([ia + 1, ib + 1, b_type])
self.structure.bonds = np.array(bonds)
bonds = self.structure.bonds
atomtypes = (
" Start File for LAMMPS \n"
+ "{0:d} atoms".format(len(self._structure))
+ " \n"
+ "{0:d} bonds".format(len(bonds))
+ " \n"
+ "{0} atom types".format(self._structure.get_number_of_species())
+ " \n"
+ "{0} bond types".format(np.max(bond_type))
+ " \n"
)
cell_dimesions = self.simulation_cell()
masses = "Masses \n\n"
el_obj_list = self._structure.get_species_objects()
for object_id, el in enumerate(el_obj_list):
masses += "{0:3d} {1:f}".format(object_id + 1, el.AtomicMass) + "\n"
atoms = "Atoms \n\n"
# atom_style bond
# format: atom-ID, molecule-ID, atom_type, x, y, z
format_str = "{0:d} {1:d} {2:d} {3:f} {4:f} {5:f} "
if self._structure.dimension == 3:
for id_atom, (x, y, z) in enumerate(coords):
id_mol, id_species = 1, el_dict[elements[id_atom]]
atoms += (
format_str.format(id_atom + 1, id_mol, id_species + 1, x, y, z)
+ "\n"
)
elif self._structure.dimension == 2:
for id_atom, (x, y) in enumerate(coords):
id_mol, id_species = 1, el_dict[elements[id_atom]]
atoms += (
format_str.format(id_atom + 1, id_mol, id_species + 1, x, y, 0.0)
+ "\n"
)
else:
raise ValueError("dimension 1 not yet implemented")
bonds_str = "Bonds \n\n"
for i_bond, (i_a, i_b, b_type) in enumerate(bonds):
bonds_str += (
"{0:d} {1:d} {2:d} {3:d}".format(i_bond + 1, b_type, i_a, i_b) + "\n"
)
return (
atomtypes
+ "\n"
+ cell_dimesions
+ "\n"
+ masses
+ "\n"
+ atoms
+ "\n"
+ bonds_str
+ "\n"
)
def structure_full(self):
"""
Write routine to create atom structure static file for atom_type='full' that can be loaded by LAMMPS
Returns:
"""
coords = self.rotate_positions(self._structure)
# extract electric charges from potential file
q_dict = {}
species_translate_list = []
for species in self.structure.species:
species_name = species.Abbreviation
q_dict[species_name] = self.potential.get_charge(species_name)
species_translate_list.append(self.potential.get_element_id(species_name))
sorted_species_list = np.array(self._potential.get_element_lst())
molecule_lst, bonds_lst, angles_lst = [], [], []
bond_type_lst, angle_type_lst = [], []
# Using a cutoff distance to draw the bonds instead of the number of neighbors
cutoff_list = list()
for val in self._bond_dict.values():
cutoff_list.append(np.max(val["cutoff_list"]))
max_cutoff = np.max(cutoff_list)
# Calculate neighbors only once
neighbors = self._structure.get_neighbors_by_distance(cutoff_radius=max_cutoff)
id_mol = 0
indices = self._structure.indices
for id_el, id_species in enumerate(indices):
id_mol += 1
molecule_lst.append([id_el, id_mol, species_translate_list[id_species]])
# Draw bonds between atoms is defined in self._bond_dict
# Go through all elements for which bonds are defined
for element, val in self._bond_dict.items():
el_1_list = self._structure.select_index(element)
if el_1_list is not None:
if len(el_1_list) > 0:
for i, v in enumerate(val["element_list"]):
el_2_list = self._structure.select_index(v)
cutoff_dist = val["cutoff_list"][i]
for j, ind in enumerate(np.array(neighbors.indices)[el_1_list]):
# Only chose those indices within the cutoff distance and which belong
# to the species defined in the element_list
# i is the index of each bond type, and j is the element index
id_el = el_1_list[j]
bool_1 = np.array(neighbors.distances)[id_el] <= cutoff_dist
act_ind = ind[bool_1]
bool_2 = np.in1d(act_ind, el_2_list)
final_ind = act_ind[bool_2]
# Get the bond and angle type
bond_type = val["bond_type_list"][i]
angle_type = val["angle_type_list"][i]
# Draw only maximum allowed bonds
final_ind = final_ind[:val["max_bond_list"][i]]
for fi in final_ind:
bonds_lst.append([id_el + 1, fi + 1])
bond_type_lst.append(bond_type)
# Draw angles if at least 2 bonds are present and if an angle type is defined for this
# particular set of bonds
if len(final_ind) >= 2 and val["angle_type_list"][i] is not None:
angles_lst.append([final_ind[0] + 1, id_el + 1, final_ind[1] + 1])
angle_type_lst.append(angle_type)
m_lst = np.array(molecule_lst)
molecule_lst = m_lst[m_lst[:, 0].argsort()]
if len(bond_type_lst) == 0:
num_bond_types = 1
else:
num_bond_types = int(np.max(bond_type_lst))
if len(angle_type_lst) == 0:
num_angle_types = 1
else:
num_angle_types = int(np.max(angle_type_lst))
atomtypes = (
" Start File for LAMMPS \n"
+ "{0:d} atoms".format(len(self._structure))
+ " \n"
+ "{0:d} bonds".format(len(bonds_lst))
+ " \n"
+ "{0:d} angles".format(len(angles_lst))
+ " \n"
+ "{0} atom types".format(len(sorted_species_list))
+ " \n"
+ "{0} bond types".format(num_bond_types)
+ " \n"
+ "{0} angle types".format(num_angle_types)
+ " \n"
)
cell_dimensions = self.simulation_cell()
masses = "Masses" + "\n\n"
for ic, el_p in enumerate(sorted_species_list):
mass = self.structure._pse[el_p].AtomicMass
masses += "{0:3d} {1:f} # ({2}) \n".format(ic + 1, mass, el_p)
atoms = "Atoms \n\n"
# format: atom-ID, molecule-ID, atom_type, q, x, y, z
format_str = "{0:d} {1:d} {2:d} {3:f} {4:f} {5:f} {6:f}"
for atom in molecule_lst:
id_atom, id_mol, id_species = atom
x, y, z = coords[id_atom]
ind = np.argwhere(np.array(species_translate_list) == id_species).flatten()[0]
el_id = self._structure.species[ind].Abbreviation
atoms += (
format_str.format(
id_atom + 1, id_mol, id_species, q_dict[el_id], x, y, z
)
+ "\n"
)
if len(bonds_lst) > 0:
bonds_str = "Bonds \n\n"
for i_bond, id_vec in enumerate(bonds_lst):
bonds_str += (
"{0:d} {1:d} {2:d} {3:d}".format(
i_bond + 1, bond_type_lst[i_bond], id_vec[0], id_vec[1]
)
+ "\n"
)
else:
bonds_str = "\n"
if len(angles_lst) > 0:
angles_str = "Angles \n\n"
for i_angle, id_vec in enumerate(angles_lst):
angles_str += (
"{0:d} {1:d} {2:d} {3:d} {4:d}".format(
i_angle + 1, angle_type_lst[i_angle], id_vec[0], id_vec[1], id_vec[2]
)
+ "\n"
)
else:
angles_str = "\n"
return (
atomtypes
+ "\n"
+ cell_dimensions
+ "\n"
+ masses
+ "\n"
+ atoms
+ "\n"
+ bonds_str
+ "\n"
+ angles_str
+ "\n"
)
def structure_charge(self):
"""
Create atom structure including the atom charges.
By convention the LAMMPS atom type numbers are chose alphabetically for the chemical species.
Returns: LAMMPS readable structure.
"""
atomtypes = (
"Start File for LAMMPS \n"
+ "{0:d} atoms".format(self._structure.get_number_of_atoms())
+ " \n"
+ "{0} atom types".format(self._structure.get_number_of_species())
+ " \n"
)
cell_dimesions = self.simulation_cell()
masses = "Masses\n\n"
for ind, obj in enumerate(self._structure.get_species_objects()):
masses += "{0:3d} {1:f}".format(ind + 1, obj.AtomicMass) + "\n"
atoms = "Atoms\n\n"
coords = self.rotate_positions(self._structure)
el_charge_lst = self._structure.charge
el_lst = self._structure.get_chemical_symbols()
el_alphabet_dict = {}
for ind, el in enumerate(self._structure.get_species_symbols()):
el_alphabet_dict[el] = ind + 1
for id_atom, (el, coord) in enumerate(zip(el_lst, coords)):
id_el = el_alphabet_dict[el]
dim = self._structure.dimension
c = np.zeros(3)
c[:dim] = coord
atoms += (
"{0:d} {1:d} {2:f} {3:.15f} {4:.15f} {5:.15f}".format(
id_atom + 1, id_el, el_charge_lst[id_atom], c[0], c[1], c[2]
)
+ "\n"
)
return atomtypes + "\n" + cell_dimesions + "\n" + masses + "\n" + atoms + "\n"
def structure_atomic(self):
"""
Write routine to create atom structure static file that can be loaded by LAMMPS
Returns:
"""
atomtypes = (
"Start File for LAMMPS \n"
+ "{0:d} atoms".format(len(self._structure))
+ " \n"
+ "{0} atom types".format(len(self._el_eam_lst))
+ " \n"
) # '{0} atom types'.format(structure.get_number_of_species()) + ' \n'
cell_dimesions = self.simulation_cell()
masses = "Masses\n\n"
el_struct_lst = self._structure.get_species_symbols()
el_obj_lst = self._structure.get_species_objects()
el_dict = {}
for id_eam, el_eam in enumerate(self._el_eam_lst):
if el_eam in el_struct_lst:
id_el = list(el_struct_lst).index(el_eam)
el = el_obj_lst[id_el]
el_dict[el] = id_eam + 1
masses += "{0:3d} {1:f}".format(id_eam + 1, el.AtomicMass) + "\n"
else:
# element in EAM file but not used in structure, use dummy for atomic mass
masses += "{0:3d} {1:f}".format(id_eam + 1, 1.00) + "\n"
atoms = "Atoms\n\n"
coords = self.rotate_positions(self._structure)
el_lst = self._structure.get_chemical_elements()
for id_atom, (el, coord) in enumerate(zip(el_lst, coords)):
id_el = el_dict[el]
dim = self._structure.dimension
c = np.zeros(3)
c[:dim] = coord
atoms += (
"{0:d} {1:d} {2:.15f} {3:.15f} {4:.15f}".format(
id_atom + 1, id_el, c[0], c[1], c[2]
)
+ "\n"
)
return atomtypes + "\n" + cell_dimesions + "\n" + masses + "\n" + atoms + "\n"
def rotate_positions(self, structure):
"""
Rotate all atomic positions in given structure according to new Prism cell
Args:
structure: Atoms-like object. Should has .positions attribute
Returns:
(list): List of rotated coordinates
"""
prism = UnfoldingPrism(self._structure.cell)
coords = [prism.pos_to_lammps(position) for position in structure.positions]
return coords
def write_lammps_datafile(structure, file_name="lammps.data", cwd=None):
lammps_str = LammpsStructure()
lammps_str.el_eam_lst = structure.get_species_symbols()
lammps_str.structure = structure
lammps_str.write_file(file_name=file_name, cwd=cwd)