root / ase / atoms.py @ 18
Historique | Voir | Annoter | Télécharger (44,24 ko)
1 |
# Copyright 2008, 2009 CAMd
|
---|---|
2 |
# (see accompanying license files for details).
|
3 |
|
4 |
"""Definition of the Atoms class.
|
5 |
|
6 |
This module defines the central object in the ASE package: the Atoms
|
7 |
object.
|
8 |
"""
|
9 |
|
10 |
from math import cos, sin |
11 |
|
12 |
import numpy as np |
13 |
|
14 |
from ase.atom import Atom |
15 |
from ase.data import atomic_numbers, chemical_symbols, atomic_masses |
16 |
import ase.units as units |
17 |
|
18 |
class Atoms(object): |
19 |
"""Atoms object.
|
20 |
|
21 |
The Atoms object can represent an isolated molecule, or a
|
22 |
periodically repeated structure. It has a unit cell and
|
23 |
there may be periodic boundary conditions along any of the three
|
24 |
unit cell axes.
|
25 |
|
26 |
Information about the atoms (atomic numbers and position) is
|
27 |
stored in ndarrays. Optionally, there can be information about
|
28 |
tags, momenta, masses, magnetic moments and charges.
|
29 |
|
30 |
In order to calculate energies, forces and stresses, a calculator
|
31 |
object has to attached to the atoms object.
|
32 |
|
33 |
Parameters:
|
34 |
|
35 |
symbols: str (formula) or list of str
|
36 |
Can be a string formula, a list of symbols or a list of
|
37 |
Atom objects. Examples: 'H2O', 'COPt12', ['H', 'H', 'O'],
|
38 |
[Atom('Ne', (x, y, z)), ...].
|
39 |
positions: list of xyz-positions
|
40 |
Atomic positions. Anything that can be converted to an
|
41 |
ndarray of shape (n, 3) will do: [(x1,y1,z1), (x2,y2,z2),
|
42 |
...].
|
43 |
scaled_positions: list of scaled-positions
|
44 |
Like positions, but given in units of the unit cell.
|
45 |
Can not be set at the same time as positions.
|
46 |
numbers: list of int
|
47 |
Atomic numbers (use only one of symbols/numbers).
|
48 |
tags: list of int
|
49 |
Special purpose tags.
|
50 |
momenta: list of xyz-momenta
|
51 |
Momenta for all atoms.
|
52 |
masses: list of float
|
53 |
Atomic masses in atomic units.
|
54 |
magmoms: list of float
|
55 |
Magnetic moments.
|
56 |
charges: list of float
|
57 |
Atomic charges.
|
58 |
cell: 3x3 matrix
|
59 |
Unit cell vectors. Can also be given as just three
|
60 |
numbers for orthorhombic cells. Default value: [1, 1, 1].
|
61 |
pbc: one or three bool
|
62 |
Periodic boundary conditions flags. Examples: True,
|
63 |
False, 0, 1, (1, 1, 0), (True, False, False). Default
|
64 |
value: False.
|
65 |
constraint: constraint object(s)
|
66 |
Used for applying one or more constraints during structure
|
67 |
optimization.
|
68 |
calculator: calculator object
|
69 |
Used to attach a calculator for calulating energies and atomic
|
70 |
forces.
|
71 |
|
72 |
Examples:
|
73 |
|
74 |
These three are equivalent:
|
75 |
|
76 |
>>> d = 1.104 # N2 bondlength
|
77 |
>>> a = Atoms('N2', [(0, 0, 0), (0, 0, d)])
|
78 |
>>> a = Atoms(numbers=[7, 7], positions=[(0, 0, 0), (0, 0, d)])
|
79 |
>>> a = Atoms([Atom('N', (0, 0, 0)), Atom('N', (0, 0, d)])
|
80 |
|
81 |
FCC gold:
|
82 |
|
83 |
>>> a = 4.05 # Gold lattice constant
|
84 |
>>> b = a / 2
|
85 |
>>> fcc = Atoms('Au',
|
86 |
... cell=[(0, b, b), (b, 0, b), (b, b, 0)],
|
87 |
... pbc=True)
|
88 |
|
89 |
Hydrogen wire:
|
90 |
|
91 |
>>> d = 0.9 # H-H distance
|
92 |
>>> L = 7.0
|
93 |
>>> h = Atoms('H', positions=[(0, L / 2, L / 2)],
|
94 |
... cell=(d, L, L),
|
95 |
... pbc=(1, 0, 0))
|
96 |
"""
|
97 |
|
98 |
def __init__(self, symbols=None, |
99 |
positions=None, numbers=None, |
100 |
tags=None, momenta=None, masses=None, |
101 |
magmoms=None, charges=None, |
102 |
scaled_positions=None,
|
103 |
cell=None, pbc=None, |
104 |
constraint=None,
|
105 |
calculator=None):
|
106 |
|
107 |
atoms = None
|
108 |
|
109 |
if hasattr(symbols, 'GetUnitCell'): |
110 |
from ase.old import OldASEListOfAtomsWrapper |
111 |
atoms = OldASEListOfAtomsWrapper(symbols) |
112 |
symbols = None
|
113 |
elif hasattr(symbols, "get_positions"): |
114 |
atoms = symbols |
115 |
symbols = None
|
116 |
elif (isinstance(symbols, (list, tuple)) and |
117 |
len(symbols) > 0 and isinstance(symbols[0], Atom)): |
118 |
# Get data from a list or tuple of Atom objects:
|
119 |
data = zip(*[atom.get_data() for atom in symbols]) |
120 |
atoms = self.__class__(None, *data) |
121 |
symbols = None
|
122 |
|
123 |
if atoms is not None: |
124 |
# Get data from another Atoms object:
|
125 |
if scaled_positions is not None: |
126 |
raise NotImplementedError |
127 |
if symbols is None and numbers is None: |
128 |
numbers = atoms.get_atomic_numbers() |
129 |
if positions is None: |
130 |
positions = atoms.get_positions() |
131 |
if tags is None and atoms.has('tags'): |
132 |
tags = atoms.get_tags() |
133 |
if momenta is None and atoms.has('momenta'): |
134 |
momenta = atoms.get_momenta() |
135 |
if magmoms is None and atoms.has('magmoms'): |
136 |
magmoms = atoms.get_initial_magnetic_moments() |
137 |
if masses is None and atoms.has('masses'): |
138 |
masses = atoms.get_masses() |
139 |
if charges is None and atoms.has('charges'): |
140 |
charges = atoms.get_charges() |
141 |
if cell is None: |
142 |
cell = atoms.get_cell() |
143 |
if pbc is None: |
144 |
pbc = atoms.get_pbc() |
145 |
if constraint is None: |
146 |
constraint = [c.copy() for c in atoms.constraints] |
147 |
if calculator is None: |
148 |
calculator = atoms.get_calculator() |
149 |
|
150 |
self.arrays = {}
|
151 |
|
152 |
if symbols is None: |
153 |
if numbers is None: |
154 |
if positions is not None: |
155 |
natoms = len(positions)
|
156 |
elif scaled_positions is not None: |
157 |
natoms = len(scaled_positions)
|
158 |
else:
|
159 |
natoms = 0
|
160 |
numbers = np.zeros(natoms, int)
|
161 |
self.new_array('numbers', numbers, int) |
162 |
else:
|
163 |
if numbers is not None: |
164 |
raise ValueError( |
165 |
'Use only one of "symbols" and "numbers".')
|
166 |
else:
|
167 |
self.new_array('numbers', symbols2numbers(symbols), int) |
168 |
|
169 |
if cell is None: |
170 |
cell = np.eye(3)
|
171 |
self.set_cell(cell)
|
172 |
|
173 |
if positions is None: |
174 |
if scaled_positions is None: |
175 |
positions = np.zeros((len(self.arrays['numbers']), 3)) |
176 |
else:
|
177 |
positions = np.dot(scaled_positions, self._cell)
|
178 |
else:
|
179 |
if scaled_positions is not None: |
180 |
raise RuntimeError, 'Both scaled and cartesian positions set!' |
181 |
self.new_array('positions', positions, float, (3,)) |
182 |
|
183 |
self.set_constraint(constraint)
|
184 |
self.set_tags(default(tags, 0)) |
185 |
self.set_momenta(default(momenta, (0.0, 0.0, 0.0))) |
186 |
self.set_masses(default(masses, None)) |
187 |
self.set_initial_magnetic_moments(default(magmoms, 0.0)) |
188 |
self.set_charges(default(charges, 0.0)) |
189 |
if pbc is None: |
190 |
pbc = False
|
191 |
self.set_pbc(pbc)
|
192 |
|
193 |
self.adsorbate_info = {}
|
194 |
|
195 |
self.set_calculator(calculator)
|
196 |
|
197 |
def set_calculator(self, calc=None): |
198 |
"""Attach calculator object."""
|
199 |
if hasattr(calc, '_SetListOfAtoms'): |
200 |
from ase.old import OldASECalculatorWrapper |
201 |
calc = OldASECalculatorWrapper(calc, self)
|
202 |
if hasattr(calc, 'set_atoms'): |
203 |
calc.set_atoms(self)
|
204 |
self.calc = calc
|
205 |
|
206 |
def get_calculator(self): |
207 |
"""Get currently attached calculator object."""
|
208 |
return self.calc |
209 |
|
210 |
def set_constraint(self, constraint=None): |
211 |
"""Apply one or more constrains.
|
212 |
|
213 |
The *constraint* argument must be one constraint object or a
|
214 |
list of constraint objects."""
|
215 |
if constraint is None: |
216 |
self._constraints = []
|
217 |
else:
|
218 |
if isinstance(constraint, (list, tuple)): |
219 |
self._constraints = constraint
|
220 |
else:
|
221 |
self._constraints = [constraint]
|
222 |
|
223 |
def _get_constraints(self): |
224 |
return self._constraints |
225 |
|
226 |
def _del_constraints(self): |
227 |
self._constraints = []
|
228 |
|
229 |
constraints = property(_get_constraints, set_constraint, _del_constraints,
|
230 |
"Constraints of the atoms.")
|
231 |
|
232 |
def set_cell(self, cell, scale_atoms=False, fix=None): |
233 |
"""Set unit cell vectors.
|
234 |
|
235 |
Parameters:
|
236 |
|
237 |
cell :
|
238 |
Unit cell. A 3x3 matrix (the three unit cell vectors) or
|
239 |
just three numbers for an orthorhombic cell.
|
240 |
scale_atoms : bool
|
241 |
Fix atomic positions or move atoms with the unit cell?
|
242 |
Default behavior is to *not* move the atoms (scale_atoms=False).
|
243 |
|
244 |
Examples:
|
245 |
|
246 |
Two equivalent ways to define an orthorhombic cell:
|
247 |
|
248 |
>>> a.set_cell([a, b, c])
|
249 |
>>> a.set_cell([(a, 0, 0), (0, b, 0), (0, 0, c)])
|
250 |
|
251 |
FCC unit cell:
|
252 |
|
253 |
>>> a.set_cell([(0, b, b), (b, 0, b), (b, b, 0)])
|
254 |
"""
|
255 |
|
256 |
if fix is not None: |
257 |
raise TypeError('Please use scale_atoms=%s' % (not fix)) |
258 |
|
259 |
cell = np.array(cell, float)
|
260 |
if cell.shape == (3,): |
261 |
cell = np.diag(cell) |
262 |
elif cell.shape != (3, 3): |
263 |
raise ValueError('Cell must be length 3 sequence or ' |
264 |
'3x3 matrix!')
|
265 |
if scale_atoms:
|
266 |
M = np.linalg.solve(self._cell, cell)
|
267 |
self.arrays['positions'][:] = np.dot(self.arrays['positions'], M) |
268 |
self._cell = cell
|
269 |
|
270 |
def get_cell(self): |
271 |
"""Get the three unit cell vectors as a 3x3 ndarray."""
|
272 |
return self._cell.copy() |
273 |
|
274 |
def set_pbc(self, pbc): |
275 |
"""Set periodic boundary condition flags."""
|
276 |
if isinstance(pbc, int): |
277 |
pbc = (pbc,) * 3
|
278 |
self._pbc = np.array(pbc, bool) |
279 |
|
280 |
def get_pbc(self): |
281 |
"""Get periodic boundary condition flags."""
|
282 |
return self._pbc.copy() |
283 |
|
284 |
def new_array(self, name, a, dtype=None, shape=None): |
285 |
"""Add new array.
|
286 |
|
287 |
If *shape* is not *None*, the shape of *a* will be checked."""
|
288 |
if dtype is not None: |
289 |
a = np.array(a, dtype) |
290 |
else:
|
291 |
a = a.copy() |
292 |
|
293 |
if name in self.arrays: |
294 |
raise RuntimeError |
295 |
|
296 |
for b in self.arrays.values(): |
297 |
if len(a) != len(b): |
298 |
raise ValueError('Array has wrong length: %d != %d.' % |
299 |
(len(a), len(b))) |
300 |
break
|
301 |
|
302 |
if shape is not None and a.shape[1:] != shape: |
303 |
raise ValueError('Array has wrong shape %s != %s.' % |
304 |
(a.shape, (a.shape[0:1] + shape))) |
305 |
|
306 |
self.arrays[name] = a
|
307 |
|
308 |
def get_array(self, name, copy=True): |
309 |
"""Get an array.
|
310 |
|
311 |
Returns a copy unless the optional argument copy is false.
|
312 |
"""
|
313 |
|
314 |
if copy:
|
315 |
return self.arrays[name].copy() |
316 |
else:
|
317 |
return self.arrays[name] |
318 |
|
319 |
def set_array(self, name, a, dtype=None, shape=None): |
320 |
"""Update array.
|
321 |
|
322 |
If *shape* is not *None*, the shape of *a* will be checked.
|
323 |
If *a* is *None*, then the array is deleted."""
|
324 |
|
325 |
b = self.arrays.get(name)
|
326 |
if b is None: |
327 |
if a is not None: |
328 |
self.new_array(name, a, dtype, shape)
|
329 |
else:
|
330 |
if a is None: |
331 |
del self.arrays[name] |
332 |
else:
|
333 |
a = np.asarray(a) |
334 |
if a.shape != b.shape:
|
335 |
raise ValueError('Array has wrong shape %s != %s.' % |
336 |
(a.shape, b.shape)) |
337 |
b[:] = a |
338 |
|
339 |
def has(self, name): |
340 |
"""Check for existance of array.
|
341 |
|
342 |
name must be one of: 'tags', 'momenta', 'masses', 'magmoms',
|
343 |
'charges'."""
|
344 |
return name in self.arrays |
345 |
|
346 |
def set_atomic_numbers(self, numbers): |
347 |
"""Set atomic numbers."""
|
348 |
self.set_array('numbers', numbers, int, ()) |
349 |
|
350 |
def get_atomic_numbers(self): |
351 |
"""Get integer array of atomic numbers."""
|
352 |
return self.arrays['numbers'] |
353 |
|
354 |
def set_chemical_symbols(self, symbols): |
355 |
"""Set chemical symbols."""
|
356 |
self.set_array('numbers', symbols2numbers(symbols), int, ()) |
357 |
|
358 |
def get_chemical_symbols(self, reduce=False): |
359 |
"""Get list of chemical symbol strings.
|
360 |
|
361 |
If reduce is True, a single string is returned, where repeated
|
362 |
elements have been contracted to a single symbol and a number.
|
363 |
E.g. instead of ['C', 'O', 'O', 'H'], the string 'CO2H' is returned.
|
364 |
"""
|
365 |
if not reduce: |
366 |
return [chemical_symbols[Z] for Z in self.arrays['numbers']] |
367 |
else:
|
368 |
num = self.get_atomic_numbers()
|
369 |
N = len(num)
|
370 |
dis = np.concatenate(([0], np.arange(1, N)[num[1:] != num[:-1]])) |
371 |
repeat = np.append(dis[1:], N) - dis
|
372 |
symbols = ''.join([chemical_symbols[num[d]] + str(r) * (r != 1) |
373 |
for r, d in zip(repeat, dis)]) |
374 |
return symbols
|
375 |
|
376 |
def set_tags(self, tags): |
377 |
"""Set tags for all atoms."""
|
378 |
self.set_array('tags', tags, int, ()) |
379 |
|
380 |
def get_tags(self): |
381 |
"""Get integer array of tags."""
|
382 |
if 'tags' in self.arrays: |
383 |
return self.arrays['tags'].copy() |
384 |
else:
|
385 |
return np.zeros(len(self), int) |
386 |
|
387 |
def set_momenta(self, momenta): |
388 |
"""Set momenta."""
|
389 |
if len(self.constraints) > 0 and momenta is not None: |
390 |
momenta = np.array(momenta) # modify a copy
|
391 |
for constraint in self.constraints: |
392 |
constraint.adjust_forces(self.arrays['positions'], momenta) |
393 |
self.set_array('momenta', momenta, float, (3,)) |
394 |
|
395 |
def set_velocities(self, velocities): |
396 |
"""Set the momenta by specifying the velocities."""
|
397 |
self.set_momenta(self.get_masses()[:,np.newaxis] * velocities) |
398 |
|
399 |
def get_momenta(self): |
400 |
"""Get array of momenta."""
|
401 |
if 'momenta' in self.arrays: |
402 |
return self.arrays['momenta'].copy() |
403 |
else:
|
404 |
return np.zeros((len(self), 3)) |
405 |
|
406 |
def set_masses(self, masses='defaults'): |
407 |
"""Set atomic masses.
|
408 |
|
409 |
The array masses should contain the a list masses. In case
|
410 |
the masses argument is not given or for those elements of the
|
411 |
masses list that are None, standard values are set."""
|
412 |
|
413 |
if masses == 'defaults': |
414 |
masses = atomic_masses[self.arrays['numbers']] |
415 |
elif isinstance(masses, (list, tuple)): |
416 |
newmasses = [] |
417 |
for m, Z in zip(masses, self.arrays['numbers']): |
418 |
if m is None: |
419 |
newmasses.append(atomic_masses[Z]) |
420 |
else:
|
421 |
newmasses.append(m) |
422 |
masses = newmasses |
423 |
self.set_array('masses', masses, float, ()) |
424 |
|
425 |
def get_masses(self): |
426 |
"""Get array of masses."""
|
427 |
if 'masses' in self.arrays: |
428 |
return self.arrays['masses'].copy() |
429 |
else:
|
430 |
return atomic_masses[self.arrays['numbers']] |
431 |
|
432 |
def set_initial_magnetic_moments(self, magmoms): |
433 |
"""Set the initial magnetic moments."""
|
434 |
self.set_array('magmoms', magmoms, float, ()) |
435 |
|
436 |
def set_magnetic_moments(self, magmoms): |
437 |
print 'Please use set_initial_magnetic_moments() instead!' |
438 |
self.set_initial_magnetic_moments(magmoms)
|
439 |
|
440 |
def get_initial_magnetic_moments(self): |
441 |
"""Get array of initial magnetic moments."""
|
442 |
if 'magmoms' in self.arrays: |
443 |
return self.arrays['magmoms'].copy() |
444 |
else:
|
445 |
return np.zeros(len(self)) |
446 |
|
447 |
def get_magnetic_moments(self): |
448 |
"""Get calculated local magnetic moments."""
|
449 |
if self.calc is None: |
450 |
raise RuntimeError('Atoms object has no calculator.') |
451 |
if self.calc.get_spin_polarized(): |
452 |
return self.calc.get_magnetic_moments(self) |
453 |
else:
|
454 |
return np.zeros(len(self)) |
455 |
|
456 |
def get_magnetic_moment(self): |
457 |
"""Get calculated total magnetic moment."""
|
458 |
if self.calc is None: |
459 |
raise RuntimeError('Atoms object has no calculator.') |
460 |
if self.calc.get_spin_polarized(): |
461 |
return self.calc.get_magnetic_moment(self) |
462 |
else:
|
463 |
return 0.0 |
464 |
|
465 |
def set_charges(self, charges): |
466 |
"""Set charges."""
|
467 |
self.set_array('charges', charges, float, ()) |
468 |
|
469 |
def get_charges(self): |
470 |
"""Get array of charges."""
|
471 |
if 'charges' in self.arrays: |
472 |
return self.arrays['charges'].copy() |
473 |
else:
|
474 |
return np.zeros(len(self)) |
475 |
|
476 |
def set_positions(self, newpositions): |
477 |
"""Set positions."""
|
478 |
positions = self.arrays['positions'] |
479 |
if self.constraints: |
480 |
newpositions = np.asarray(newpositions, float)
|
481 |
for constraint in self.constraints: |
482 |
constraint.adjust_positions(positions, newpositions) |
483 |
|
484 |
self.set_array('positions', newpositions, shape=(3,)) |
485 |
|
486 |
def get_positions(self): |
487 |
"""Get array of positions."""
|
488 |
return self.arrays['positions'].copy() |
489 |
|
490 |
def get_potential_energy(self): |
491 |
"""Calculate potential energy."""
|
492 |
if self.calc is None: |
493 |
raise RuntimeError('Atoms object has no calculator.') |
494 |
return self.calc.get_potential_energy(self) |
495 |
|
496 |
def get_potential_energies(self): |
497 |
"""Calculate the potential energies of all the atoms.
|
498 |
|
499 |
Only available with calculators supporting per-atom energies
|
500 |
(e.g. classical potentials).
|
501 |
"""
|
502 |
if self.calc is None: |
503 |
raise RuntimeError('Atoms object has no calculator.') |
504 |
return self.calc.get_potential_energies(self) |
505 |
|
506 |
def get_kinetic_energy(self): |
507 |
"""Get the kinetic energy."""
|
508 |
momenta = self.arrays.get('momenta') |
509 |
if momenta is None: |
510 |
return 0.0 |
511 |
return 0.5 * np.vdot(momenta, self.get_velocities()) |
512 |
|
513 |
def get_velocities(self): |
514 |
"""Get array of velocities."""
|
515 |
momenta = self.arrays.get('momenta') |
516 |
if momenta is None: |
517 |
return None |
518 |
m = self.arrays.get('masses') |
519 |
if m is None: |
520 |
m = atomic_masses[self.arrays['numbers']] |
521 |
return momenta / m.reshape(-1, 1) |
522 |
|
523 |
def get_total_energy(self): |
524 |
"""Get the total energy - potential plus kinetic energy."""
|
525 |
return self.get_potential_energy() + self.get_kinetic_energy() |
526 |
|
527 |
def get_forces(self, apply_constraint=True): |
528 |
"""Calculate atomic forces.
|
529 |
|
530 |
Ask the attached calculator to calculate the forces and apply
|
531 |
constraints. Use *apply_constraint=False* to get the raw
|
532 |
forces."""
|
533 |
|
534 |
if self.calc is None: |
535 |
raise RuntimeError('Atoms object has no calculator.') |
536 |
forces = self.calc.get_forces(self) |
537 |
if apply_constraint:
|
538 |
for constraint in self.constraints: |
539 |
constraint.adjust_forces(self.arrays['positions'], forces) |
540 |
return forces
|
541 |
|
542 |
def get_stress(self): |
543 |
"""Calculate stress tensor.
|
544 |
|
545 |
Returns an array of the six independent components of the
|
546 |
symmetric stress tensor, in the traditional order
|
547 |
(s_xx, s_yy, s_zz, s_yz, s_xz, s_xy).
|
548 |
"""
|
549 |
if self.calc is None: |
550 |
raise RuntimeError('Atoms object has no calculator.') |
551 |
stress = self.calc.get_stress(self) |
552 |
shape = getattr(stress, "shape", None) |
553 |
if shape == (3,3): |
554 |
return np.array([stress[0,0], stress[1,1], stress[2,2], |
555 |
stress[1,2], stress[0,2], stress[0,1]]) |
556 |
else:
|
557 |
# Hopefully a 6-vector, but don't check in case some weird
|
558 |
# calculator does something else.
|
559 |
return stress
|
560 |
|
561 |
def get_stresses(self): |
562 |
"""Calculate the stress-tensor of all the atoms.
|
563 |
|
564 |
Only available with calculators supporting per-atom energies and
|
565 |
stresses (e.g. classical potentials). Even for such calculators
|
566 |
there is a certain arbitrariness in defining per-atom stresses.
|
567 |
"""
|
568 |
if self.calc is None: |
569 |
raise RuntimeError('Atoms object has no calculator.') |
570 |
return self.calc.get_stresses(self) |
571 |
|
572 |
def get_dipole_moment(self): |
573 |
"""Calculate the electric dipole moment for the atoms object.
|
574 |
|
575 |
Only available for calculators which has a get_dipole_moment() method."""
|
576 |
if self.calc is None: |
577 |
raise RuntimeError('Atoms object has no calculator.') |
578 |
try:
|
579 |
dipole = self.calc.get_dipole_moment(self) |
580 |
except AttributeError: |
581 |
raise AttributeError('Calculator object has no get_dipole_moment method.') |
582 |
return dipole
|
583 |
|
584 |
def copy(self): |
585 |
"""Return a copy."""
|
586 |
import copy |
587 |
atoms = Atoms(cell=self._cell, pbc=self._pbc) |
588 |
|
589 |
atoms.arrays = {} |
590 |
for name, a in self.arrays.items(): |
591 |
atoms.arrays[name] = a.copy() |
592 |
atoms.constraints = copy.deepcopy(self.constraints)
|
593 |
atoms.adsorbate_info = copy.deepcopy(self.adsorbate_info)
|
594 |
return atoms
|
595 |
|
596 |
def __len__(self): |
597 |
return len(self.arrays['positions']) |
598 |
|
599 |
def get_number_of_atoms(self): |
600 |
"""Returns the number of atoms.
|
601 |
|
602 |
Equivalent to len(atoms) in the standard ASE Atoms class.
|
603 |
"""
|
604 |
return len(self) |
605 |
|
606 |
def __repr__(self): |
607 |
num = self.get_atomic_numbers()
|
608 |
N = len(num)
|
609 |
if N == 0: |
610 |
symbols = ''
|
611 |
elif N <= 60: |
612 |
symbols = self.get_chemical_symbols(reduce=True) |
613 |
else:
|
614 |
symbols = ''.join([chemical_symbols[Z] for Z in num[:15]]) + '...' |
615 |
s = "%s(symbols='%s', " %(self.__class__.__name__, symbols) |
616 |
for name in self.arrays: |
617 |
if name == 'numbers': |
618 |
continue
|
619 |
s += '%s=..., ' % name
|
620 |
if (self._cell - np.diag(self._cell.diagonal())).any(): |
621 |
s += 'cell=%s, ' % self._cell.tolist() |
622 |
else:
|
623 |
s += 'cell=%s, ' % self._cell.diagonal().tolist() |
624 |
s += 'pbc=%s, ' % self._pbc.tolist() |
625 |
if len(self.constraints) == 1: |
626 |
s += 'constraint=%s, ' % repr(self.constraints[0]) |
627 |
if len(self.constraints) > 1: |
628 |
s += 'constraint=%s, ' % repr(self.constraints) |
629 |
if self.calc is not None: |
630 |
s += 'calculator=%s(...), ' % self.calc.__class__.__name__ |
631 |
return s[:-2] + ')' |
632 |
|
633 |
def __add__(self, other): |
634 |
atoms = self.copy()
|
635 |
atoms += other |
636 |
return atoms
|
637 |
|
638 |
def extend(self, other): |
639 |
"""Extend atoms object by appending atoms from *other*."""
|
640 |
if isinstance(other, Atom): |
641 |
other = self.__class__([other])
|
642 |
|
643 |
n1 = len(self) |
644 |
n2 = len(other)
|
645 |
|
646 |
for name, a1 in self.arrays.items(): |
647 |
a = np.zeros((n1 + n2,) + a1.shape[1:], a1.dtype)
|
648 |
a[:n1] = a1 |
649 |
a2 = other.arrays.get(name) |
650 |
if a2 is not None: |
651 |
a[n1:] = a2 |
652 |
self.arrays[name] = a
|
653 |
|
654 |
for name, a2 in other.arrays.items(): |
655 |
if name in self.arrays: |
656 |
continue
|
657 |
a = np.zeros((n1 + n2,) + a2.shape[1:], a2.dtype)
|
658 |
a[n1:] = a2 |
659 |
self.set_array(name, a)
|
660 |
|
661 |
return self |
662 |
|
663 |
__iadd__ = extend |
664 |
|
665 |
def append(self, atom): |
666 |
"""Append atom to end."""
|
667 |
self.extend(self.__class__([atom])) |
668 |
|
669 |
def __getitem__(self, i): |
670 |
"""Return a subset of the atoms.
|
671 |
|
672 |
i -- scalar integer, list of integers, or slice object
|
673 |
describing which atoms to return.
|
674 |
|
675 |
If i is a scalar, return an Atom object. If i is a list or a
|
676 |
slice, return an Atoms object with the same cell, pbc, and
|
677 |
other associated info as the original Atoms object. The
|
678 |
indices of the constraints will be shuffled so that they match
|
679 |
the indexing in the subset returned.
|
680 |
|
681 |
"""
|
682 |
if isinstance(i, int): |
683 |
natoms = len(self) |
684 |
if i < -natoms or i >= natoms: |
685 |
raise IndexError('Index out of range.') |
686 |
|
687 |
return Atom(atoms=self, index=i) |
688 |
|
689 |
import copy |
690 |
from ase.constraints import FixConstraint |
691 |
|
692 |
atoms = self.__class__(cell=self._cell, pbc=self._pbc) |
693 |
# TODO: Do we need to shuffle indices in adsorbate_info too?
|
694 |
atoms.adsorbate_info = self.adsorbate_info
|
695 |
|
696 |
atoms.arrays = {} |
697 |
for name, a in self.arrays.items(): |
698 |
atoms.arrays[name] = a[i].copy() |
699 |
|
700 |
# Constraints need to be deepcopied, since we need to shuffle
|
701 |
# the indices
|
702 |
atoms.constraints = copy.deepcopy(self.constraints)
|
703 |
condel = [] |
704 |
for con in atoms.constraints: |
705 |
if isinstance(con, FixConstraint): |
706 |
try:
|
707 |
con.index_shuffle(i) |
708 |
except IndexError: |
709 |
condel.append(con) |
710 |
for con in condel: |
711 |
atoms.constraints.remove(con) |
712 |
return atoms
|
713 |
|
714 |
def __delitem__(self, i): |
715 |
if len(self._constraints) > 0: |
716 |
raise RuntimeError('Remove constraint using set_constraint() ' + |
717 |
'before deleting atoms.')
|
718 |
mask = np.ones(len(self), bool) |
719 |
mask[i] = False
|
720 |
for name, a in self.arrays.items(): |
721 |
self.arrays[name] = a[mask]
|
722 |
|
723 |
def pop(self, i=-1): |
724 |
"""Remove and return atom at index *i* (default last)."""
|
725 |
atom = self[i]
|
726 |
atom.cut_reference_to_atoms() |
727 |
del self[i] |
728 |
return atom
|
729 |
|
730 |
def __imul__(self, m): |
731 |
"""In-place repeat of atoms."""
|
732 |
if isinstance(m, int): |
733 |
m = (m, m, m) |
734 |
|
735 |
M = np.product(m) |
736 |
n = len(self) |
737 |
|
738 |
for name, a in self.arrays.items(): |
739 |
self.arrays[name] = np.tile(a, (M,) + (1,) * (len(a.shape) - 1)) |
740 |
|
741 |
positions = self.arrays['positions'] |
742 |
i0 = 0
|
743 |
for m2 in range(m[2]): |
744 |
for m1 in range(m[1]): |
745 |
for m0 in range(m[0]): |
746 |
i1 = i0 + n |
747 |
positions[i0:i1] += np.dot((m0, m1, m2), self._cell)
|
748 |
i0 = i1 |
749 |
|
750 |
if self.constraints is not None: |
751 |
self.constraints = [c.repeat(m, n) for c in self.constraints] |
752 |
|
753 |
self._cell = np.array([m[c] * self._cell[c] for c in range(3)]) |
754 |
|
755 |
return self |
756 |
|
757 |
def repeat(self, rep): |
758 |
"""Create new repeated atoms object.
|
759 |
|
760 |
The *rep* argument should be a sequence of three positive
|
761 |
integers like *(2,3,1)* or a single integer (*r*) equivalent
|
762 |
to *(r,r,r)*."""
|
763 |
|
764 |
atoms = self.copy()
|
765 |
atoms *= rep |
766 |
return atoms
|
767 |
|
768 |
__mul__ = repeat |
769 |
|
770 |
def translate(self, displacement): |
771 |
"""Translate atomic positions.
|
772 |
|
773 |
The displacement argument can be a float an xyz vector or an
|
774 |
nx3 array (where n is the number of atoms)."""
|
775 |
|
776 |
self.arrays['positions'] += np.array(displacement) |
777 |
|
778 |
def center(self, vacuum=None, axis=None): |
779 |
"""Center atoms in unit cell.
|
780 |
|
781 |
Centers the atoms in the unit cell, so there is the same
|
782 |
amount of vacuum on all sides.
|
783 |
|
784 |
Parameters:
|
785 |
|
786 |
vacuum (default: None): If specified adjust the amount of
|
787 |
vacuum when centering. If vacuum=10.0 there will thus be 10
|
788 |
Angstrom of vacuum on each side.
|
789 |
|
790 |
axis (default: None): If specified, only act on the specified
|
791 |
axis. Default: Act on all axes.
|
792 |
"""
|
793 |
# Find the orientations of the faces of the unit cell
|
794 |
c = self.get_cell()
|
795 |
dirs = np.zeros_like(c) |
796 |
for i in range(3): |
797 |
dirs[i] = np.cross(c[i-1], c[i-2]) |
798 |
dirs[i] /= np.sqrt(np.dot(dirs[i], dirs[i])) #Normalize
|
799 |
if np.dot(dirs[i], c[i]) < 0.0: |
800 |
dirs[i] *= -1
|
801 |
|
802 |
# Now, decide how much each basis vector should be made longer
|
803 |
if axis is None: |
804 |
axes = (0,1,2) |
805 |
else:
|
806 |
axes = (axis,) |
807 |
p = self.arrays['positions'] |
808 |
longer = np.zeros(3)
|
809 |
shift = np.zeros(3)
|
810 |
for i in axes: |
811 |
p0 = np.dot(p, dirs[i]).min() |
812 |
p1 = np.dot(p, dirs[i]).max() |
813 |
height = np.dot(c[i], dirs[i]) |
814 |
if vacuum is not None: |
815 |
lng = (p1 - p0 + 2*vacuum) - height
|
816 |
else:
|
817 |
lng = 0.0 # Do not change unit cell size! |
818 |
top = lng + height - p1 |
819 |
shf = 0.5 * (top - p0)
|
820 |
cosphi = np.dot(c[i], dirs[i]) / np.sqrt(np.dot(c[i], c[i])) |
821 |
longer[i] = lng / cosphi |
822 |
shift[i] = shf / cosphi |
823 |
|
824 |
# Now, do it!
|
825 |
translation = np.zeros(3)
|
826 |
for i in axes: |
827 |
nowlen = np.sqrt(np.dot(c[i], c[i])) |
828 |
self._cell[i] *= 1 + longer[i] / nowlen |
829 |
translation += shift[i] * c[i] / nowlen |
830 |
self.arrays['positions'] += translation |
831 |
|
832 |
def get_center_of_mass(self, scaled = False): |
833 |
"""Get the center of mass.
|
834 |
|
835 |
If scaled=True the center of mass in scaled coordinates
|
836 |
is returned."""
|
837 |
m = self.arrays.get('masses') |
838 |
if m is None: |
839 |
m = atomic_masses[self.arrays['numbers']] |
840 |
com = np.dot(m, self.arrays['positions']) / m.sum() |
841 |
if scaled:
|
842 |
return np.dot(np.linalg.inv(self._cell), com) |
843 |
else:
|
844 |
return com
|
845 |
|
846 |
def get_moments_of_inertia(self): |
847 |
'''Get the moments of inertia
|
848 |
|
849 |
The three principal moments of inertia are computed from the
|
850 |
eigenvalues of the inertial tensor. periodic boundary
|
851 |
conditions are ignored. Units of the moments of inertia are
|
852 |
amu*angstrom**2.
|
853 |
'''
|
854 |
|
855 |
com = self.get_center_of_mass()
|
856 |
positions = self.get_positions()
|
857 |
positions -= com #translate center of mass to origin
|
858 |
masses = self.get_masses()
|
859 |
|
860 |
#initialize elements of the inertial tensor
|
861 |
I11 = I22 = I33 = I12 = I13 = I23 = 0.0
|
862 |
for i in range(len(self)): |
863 |
x,y,z = positions[i] |
864 |
m = masses[i] |
865 |
|
866 |
I11 += m*(y**2 + z**2) |
867 |
I22 += m*(x**2 + z**2) |
868 |
I33 += m*(x**2 + y**2) |
869 |
I12 += -m*x*y |
870 |
I13 += -m*x*z |
871 |
I23 += -m*y*z |
872 |
|
873 |
I = np.array([[I11, I12, I13], |
874 |
[I12, I22, I23], |
875 |
[I13, I23, I33]]) |
876 |
|
877 |
evals, evecs = np.linalg.eig(I) |
878 |
return evals
|
879 |
|
880 |
def rotate(self, v, a=None, center=(0, 0, 0), rotate_cell=False): |
881 |
"""Rotate atoms.
|
882 |
|
883 |
Rotate the angle *a* around the vector *v*. If *a* is not
|
884 |
given, the length of *v* is used as the angle. If *a* is a
|
885 |
vector, then *v* is rotated into *a*. The point at *center*
|
886 |
is fixed. Use *center='COM'* to fix the center of mass.
|
887 |
Vectors can also be strings: 'x', '-x', 'y', ... .
|
888 |
If both *v* and *a* are vectors: rotate *v* into *a*.
|
889 |
|
890 |
Examples:
|
891 |
|
892 |
Rotate 90 degrees around the z-axis, so that the x-axis is
|
893 |
rotated into the y-axis:
|
894 |
|
895 |
>>> a = pi / 2
|
896 |
>>> atoms.rotate('z', a)
|
897 |
>>> atoms.rotate((0, 0, 1), a)
|
898 |
>>> atoms.rotate('-z', -a)
|
899 |
>>> atoms.rotate((0, 0, a))
|
900 |
>>> atoms.rotate('x', 'y')
|
901 |
"""
|
902 |
|
903 |
norm = np.linalg.norm |
904 |
v = string2vector(v) |
905 |
if a is None: |
906 |
a = norm(v) |
907 |
if isinstance(a, (float, int)): |
908 |
v /= norm(v) |
909 |
c = cos(a) |
910 |
s = sin(a) |
911 |
else:
|
912 |
v2 = string2vector(a) |
913 |
v /= norm(v) |
914 |
v2 /= norm(v2) |
915 |
c = np.dot(v, v2) |
916 |
v = np.cross(v, v2) |
917 |
s = norm(v) |
918 |
# In case *v* and *a* are parallel, np.cross(v, v2) vanish
|
919 |
# and can't be used as a rotation axis. However, in this
|
920 |
# case any rotation axis perpendicular to v2 will do.
|
921 |
eps = 1e-7
|
922 |
if s < eps:
|
923 |
v = np.cross((0, 0, 1), v2) |
924 |
if norm(v) < eps:
|
925 |
v = np.cross((1, 0, 0), v2) |
926 |
assert norm(v) >= eps
|
927 |
if s > 0: v /= s |
928 |
|
929 |
if isinstance(center, str) and center.lower() == 'com': |
930 |
center = self.get_center_of_mass()
|
931 |
|
932 |
p = self.arrays['positions'] - center |
933 |
self.arrays['positions'][:] = (c * p - |
934 |
np.cross(p, s * v) + |
935 |
np.outer(np.dot(p, v), (1.0 - c) * v)+
|
936 |
center) |
937 |
if rotate_cell:
|
938 |
rotcell = self.get_cell()
|
939 |
rotcell[:] = (c * rotcell - |
940 |
np.cross(rotcell, s * v) + |
941 |
np.outer(np.dot(rotcell, v), (1.0 - c) * v))
|
942 |
self.set_cell(rotcell)
|
943 |
|
944 |
def rotate_euler(self, center=(0, 0, 0), phi=0.0, theta=0.0, psi=0.0): |
945 |
"""Rotate atoms via Euler angles.
|
946 |
|
947 |
See e.g http://mathworld.wolfram.com/EulerAngles.html for explanation.
|
948 |
|
949 |
Parameters:
|
950 |
|
951 |
center :
|
952 |
The point to rotate about. A sequence of length 3 with the
|
953 |
coordinates, or 'COM' to select the center of mass.
|
954 |
phi :
|
955 |
The 1st rotation angle around the z axis.
|
956 |
theta :
|
957 |
Rotation around the x axis.
|
958 |
psi :
|
959 |
2nd rotation around the z axis.
|
960 |
|
961 |
"""
|
962 |
if isinstance(center, str) and center.lower() == 'com': |
963 |
center = self.get_center_of_mass()
|
964 |
else:
|
965 |
center = np.array(center) |
966 |
# First move the molecule to the origin In contrast to MATLAB,
|
967 |
# numpy broadcasts the smaller array to the larger row-wise,
|
968 |
# so there is no need to play with the Kronecker product.
|
969 |
rcoords = self.positions - center
|
970 |
# First Euler rotation about z in matrix form
|
971 |
D = np.array(((cos(phi), sin(phi), 0.),
|
972 |
(-sin(phi), cos(phi), 0.),
|
973 |
(0., 0., 1.))) |
974 |
# Second Euler rotation about x:
|
975 |
C = np.array(((1., 0., 0.), |
976 |
(0., cos(theta), sin(theta)),
|
977 |
(0., -sin(theta), cos(theta))))
|
978 |
# Third Euler rotation, 2nd rotation about z:
|
979 |
B = np.array(((cos(psi), sin(psi), 0.),
|
980 |
(-sin(psi), cos(psi), 0.),
|
981 |
(0., 0., 1.))) |
982 |
# Total Euler rotation
|
983 |
A = np.dot(B, np.dot(C, D)) |
984 |
# Do the rotation
|
985 |
rcoords = np.dot(A, np.transpose(rcoords)) |
986 |
# Move back to the rotation point
|
987 |
self.positions = np.transpose(rcoords) + center
|
988 |
|
989 |
def get_dihedral(self,list): |
990 |
"""
|
991 |
calculate dihedral angle between the vectors list[0]->list[1] and list[2]->list[3],
|
992 |
where list contains the atomic indexes in question.
|
993 |
"""
|
994 |
# vector 0->1, 1->2, 2->3 and their normalized cross products:
|
995 |
a = self.positions[list[1]]-self.positions[list[0]] |
996 |
b = self.positions[list[2]]-self.positions[list[1]] |
997 |
c = self.positions[list[3]]-self.positions[list[2]] |
998 |
bxa = np.cross(b,a) |
999 |
bxa /= np.sqrt(np.vdot(bxa,bxa)) |
1000 |
cxb = np.cross(c,b) |
1001 |
cxb /= np.sqrt(np.vdot(cxb,cxb)) |
1002 |
angle = np.vdot(bxa,cxb) |
1003 |
# check for numerical trouble due to finite precision:
|
1004 |
if angle < -1: angle = -1 |
1005 |
if angle > 1: angle = 1 |
1006 |
angle = np.arccos(angle) |
1007 |
if (np.vdot(bxa,c)) > 0: angle = 2*np.pi-angle |
1008 |
return angle
|
1009 |
|
1010 |
def set_dihedral(self,list,angle,mask=None): |
1011 |
"""
|
1012 |
set the dihedral angle between vectors list[0]->list[1] and
|
1013 |
list[2]->list[3] by changing the atom indexed by list[3]
|
1014 |
if mask is not None, all the atoms described in mask
|
1015 |
(read: the entire subgroup) are moved
|
1016 |
|
1017 |
example: the following defines a very crude
|
1018 |
ethane-like molecule and twists one half of it by 30 degrees.
|
1019 |
|
1020 |
>>> atoms = Atoms('HHCCHH',[[-1,1,0],[-1,-1,0],[0,0,0],[1,0,0],[2,1,0],[2,-1,0]])
|
1021 |
>>> atoms.set_dihedral([1,2,3,4],7*pi/6,mask=[0,0,0,1,1,1])
|
1022 |
"""
|
1023 |
# if not provided, set mask to the last atom in the dihedral description
|
1024 |
if mask is None: |
1025 |
mask = np.zeros(len(self)) |
1026 |
mask[list[3]] = 1 |
1027 |
# compute necessary in dihedral change, from current value
|
1028 |
current =self.get_dihedral(list) |
1029 |
diff = angle - current |
1030 |
# do rotation of subgroup by copying it to temporary atoms object and then rotating that
|
1031 |
axis = self.positions[list[2]]-self.positions[list[1]] |
1032 |
center = self.positions[list[2]] |
1033 |
# recursive object definition might not be the most elegant thing, more generally useful might be a rotation function with a mask?
|
1034 |
group = self.__class__()
|
1035 |
for i in range(len(self)): |
1036 |
if mask[i]:
|
1037 |
group += self[i]
|
1038 |
group.translate(-center) |
1039 |
group.rotate(axis,diff) |
1040 |
group.translate(center) |
1041 |
# set positions in original atoms object
|
1042 |
j = 0
|
1043 |
for i in range(len(self)): |
1044 |
if mask[i]:
|
1045 |
self.positions[i] = group[j].get_position()
|
1046 |
j += 1
|
1047 |
|
1048 |
def rotate_dihedral(self,list,angle,mask=None): |
1049 |
""" complementing the two routines above: rotate a group by a predefined dihedral angle,
|
1050 |
starting from its current configuration
|
1051 |
"""
|
1052 |
start = self.get_dihedral(list) |
1053 |
self.set_dihedral(list,angle+start,mask) |
1054 |
|
1055 |
def rattle(self, stdev=0.001, seed=42): |
1056 |
"""Randomly displace atoms.
|
1057 |
|
1058 |
This method adds random displacements to the atomic positions,
|
1059 |
taking a possible constraint into account. The random numbers are
|
1060 |
drawn from a normal distribution of standard deviation stdev.
|
1061 |
|
1062 |
For a parallel calculation, it is important to use the same
|
1063 |
seed on all processors! """
|
1064 |
|
1065 |
rs = np.random.RandomState(seed) |
1066 |
positions = self.arrays['positions'] |
1067 |
self.set_positions(positions +
|
1068 |
rs.normal(scale=stdev, size=positions.shape)) |
1069 |
|
1070 |
def get_distance(self, a0, a1, mic=False): |
1071 |
"""Return distance between two atoms.
|
1072 |
|
1073 |
Use mic=True to use the Minimum Image Convention.
|
1074 |
"""
|
1075 |
|
1076 |
R = self.arrays['positions'] |
1077 |
D = R[a1] - R[a0] |
1078 |
if mic:
|
1079 |
raise NotImplemented # XXX |
1080 |
return np.linalg.norm(D)
|
1081 |
|
1082 |
def set_distance(self, a0, a1, distance, fix=0.5): |
1083 |
"""Set the distance between two atoms.
|
1084 |
|
1085 |
Set the distance between atoms *a0* and *a1* to *distance*.
|
1086 |
By default, the center of the two atoms will be fixed. Use
|
1087 |
*fix=0* to fix the first atom, *fix=1* to fix the second
|
1088 |
atom and *fix=0.5* (default) to fix the center of the bond."""
|
1089 |
|
1090 |
R = self.arrays['positions'] |
1091 |
D = R[a1] - R[a0] |
1092 |
x = 1.0 - distance / np.linalg.norm(D)
|
1093 |
R[a0] += (x * fix) * D |
1094 |
R[a1] -= (x * (1.0 - fix)) * D
|
1095 |
|
1096 |
def get_scaled_positions(self): |
1097 |
"""Get positions relative to unit cell.
|
1098 |
|
1099 |
Atoms outside the unit cell will be wrapped into the cell in
|
1100 |
those directions with periodic boundary conditions so that the
|
1101 |
scaled coordinates are beween zero and one."""
|
1102 |
|
1103 |
scaled = np.linalg.solve(self._cell.T, self.arrays['positions'].T).T |
1104 |
for i in range(3): |
1105 |
if self._pbc[i]: |
1106 |
# Yes, we need to do it twice.
|
1107 |
# See the scaled_positions.py test
|
1108 |
scaled[:, i] %= 1.0
|
1109 |
scaled[:, i] %= 1.0
|
1110 |
return scaled
|
1111 |
|
1112 |
def set_scaled_positions(self, scaled): |
1113 |
"""Set positions relative to unit cell."""
|
1114 |
self.arrays['positions'][:] = np.dot(scaled, self._cell) |
1115 |
|
1116 |
def get_temperature(self): |
1117 |
"""Get the temperature. in Kelvin"""
|
1118 |
ekin = self.get_kinetic_energy() / len(self) |
1119 |
return ekin /(1.5*units.kB) |
1120 |
|
1121 |
def get_isotropic_pressure(self, stress): |
1122 |
""" get the current calculated pressure, assume isotropic medium.
|
1123 |
in Bar
|
1124 |
"""
|
1125 |
if type(stress) == type(1.0) or type(stress) == type(1): |
1126 |
return -stress * 1e-5 / units.Pascal |
1127 |
elif stress.shape == (3,3): |
1128 |
return (-(stress[0, 0] + stress[1, 1] + stress[2, 2]) / 3.0) * \ |
1129 |
1e-5 / units.Pascal
|
1130 |
elif stress.shape == (6,): |
1131 |
return (-(stress[0] + stress[1] + stress[2]) / 3.0) * \ |
1132 |
1e-5 / units.Pascal
|
1133 |
else:
|
1134 |
raise ValueError, "The external stress has the wrong shape." |
1135 |
|
1136 |
|
1137 |
|
1138 |
def __eq__(self, other): |
1139 |
"""Check for identity of two atoms objects.
|
1140 |
|
1141 |
Identity means: same positions, atomic numbers, unit cell and
|
1142 |
periodic boundary conditions."""
|
1143 |
try:
|
1144 |
a = self.arrays
|
1145 |
b = other.arrays |
1146 |
return (len(self) == len(other) and |
1147 |
(a['positions'] == b['positions']).all() and |
1148 |
(a['numbers'] == b['numbers']).all() and |
1149 |
(self._cell == other.cell).all() and |
1150 |
(self._pbc == other.pbc).all())
|
1151 |
except AttributeError: |
1152 |
return NotImplemented |
1153 |
|
1154 |
def __ne__(self, other): |
1155 |
eq = self.__eq__(other)
|
1156 |
if eq is NotImplemented: |
1157 |
return eq
|
1158 |
else:
|
1159 |
return not eq |
1160 |
|
1161 |
__hash__ = None
|
1162 |
|
1163 |
def get_volume(self): |
1164 |
"""Get volume of unit cell."""
|
1165 |
return abs(np.linalg.det(self._cell)) |
1166 |
|
1167 |
def _get_positions(self): |
1168 |
"""Return reference to positions-array for inplace manipulations."""
|
1169 |
return self.arrays['positions'] |
1170 |
|
1171 |
def _set_positions(self, pos): |
1172 |
"""Set positions directly, bypassing constraints."""
|
1173 |
self.arrays['positions'][:] = pos |
1174 |
|
1175 |
positions = property(_get_positions, _set_positions,
|
1176 |
doc='Attribute for direct ' +
|
1177 |
'manipulation of the positions.')
|
1178 |
|
1179 |
def _get_atomic_numbers(self): |
1180 |
"""Return reference to atomic numbers for inplace manipulations."""
|
1181 |
return self.arrays['numbers'] |
1182 |
|
1183 |
numbers = property(_get_atomic_numbers, set_atomic_numbers,
|
1184 |
doc='Attribute for direct ' +
|
1185 |
'manipulation of the atomic numbers.')
|
1186 |
|
1187 |
def _get_cell(self): |
1188 |
"""Return reference to unit cell for inplace manipulations."""
|
1189 |
return self._cell |
1190 |
|
1191 |
cell = property(_get_cell, set_cell, doc='Attribute for direct ' + |
1192 |
'manipulation of the unit cell.')
|
1193 |
|
1194 |
def _get_pbc(self): |
1195 |
"""Return reference to pbc-flags for inplace manipulations."""
|
1196 |
return self._pbc |
1197 |
|
1198 |
pbc = property(_get_pbc, set_pbc,
|
1199 |
doc='Attribute for direct manipulation ' +
|
1200 |
'of the periodic boundary condition flags.')
|
1201 |
|
1202 |
def get_name(self): |
1203 |
"""Return a name extracted from the elements."""
|
1204 |
elements = {} |
1205 |
for a in self: |
1206 |
try:
|
1207 |
elements[a.symbol] += 1
|
1208 |
except:
|
1209 |
elements[a.symbol] = 1
|
1210 |
name = ''
|
1211 |
for element in elements: |
1212 |
name += element |
1213 |
if elements[element] > 1: |
1214 |
name += str(elements[element])
|
1215 |
return name
|
1216 |
|
1217 |
def write(self, filename, format=None): |
1218 |
"""Write yourself to a file."""
|
1219 |
from ase.io import write |
1220 |
write(filename, self, format)
|
1221 |
|
1222 |
def string2symbols(s): |
1223 |
"""Convert string to list of chemical symbols."""
|
1224 |
n = len(s)
|
1225 |
|
1226 |
if n == 0: |
1227 |
return []
|
1228 |
|
1229 |
c = s[0]
|
1230 |
|
1231 |
if c.isdigit():
|
1232 |
i = 1
|
1233 |
while i < n and s[i].isdigit(): |
1234 |
i += 1
|
1235 |
return int(s[:i]) * string2symbols(s[i:]) |
1236 |
|
1237 |
if c == '(': |
1238 |
p = 0
|
1239 |
for i, c in enumerate(s): |
1240 |
if c == '(': |
1241 |
p += 1
|
1242 |
elif c == ')': |
1243 |
p -= 1
|
1244 |
if p == 0: |
1245 |
break
|
1246 |
j = i + 1
|
1247 |
while j < n and s[j].isdigit(): |
1248 |
j += 1
|
1249 |
if j > i + 1: |
1250 |
m = int(s[i + 1:j]) |
1251 |
else:
|
1252 |
m = 1
|
1253 |
return m * string2symbols(s[1:i]) + string2symbols(s[j:]) |
1254 |
|
1255 |
if c.isupper():
|
1256 |
i = 1
|
1257 |
if 1 < n and s[1].islower(): |
1258 |
i += 1
|
1259 |
j = i |
1260 |
while j < n and s[j].isdigit(): |
1261 |
j += 1
|
1262 |
if j > i:
|
1263 |
m = int(s[i:j])
|
1264 |
else:
|
1265 |
m = 1
|
1266 |
return m * [s[:i]] + string2symbols(s[j:])
|
1267 |
else:
|
1268 |
raise ValueError |
1269 |
|
1270 |
def symbols2numbers(symbols): |
1271 |
if isinstance(symbols, str): |
1272 |
symbols = string2symbols(symbols) |
1273 |
numbers = [] |
1274 |
for s in symbols: |
1275 |
if isinstance(s, str): |
1276 |
numbers.append(atomic_numbers[s]) |
1277 |
else:
|
1278 |
numbers.append(s) |
1279 |
return numbers
|
1280 |
|
1281 |
def string2vector(v): |
1282 |
if isinstance(v, str): |
1283 |
if v[0] == '-': |
1284 |
return -string2vector(v[1:]) |
1285 |
w = np.zeros(3)
|
1286 |
w['xyz'.index(v)] = 1.0 |
1287 |
return w
|
1288 |
return np.array(v, float) |
1289 |
|
1290 |
def default(data, dflt): |
1291 |
"""Helper function for setting default values."""
|
1292 |
if data is None: |
1293 |
return None |
1294 |
elif isinstance(data, (list, tuple)): |
1295 |
newdata = [] |
1296 |
allnone = True
|
1297 |
for x in data: |
1298 |
if x is None: |
1299 |
newdata.append(dflt) |
1300 |
else:
|
1301 |
newdata.append(x) |
1302 |
allnone = False
|
1303 |
if allnone:
|
1304 |
return None |
1305 |
return newdata
|
1306 |
else:
|
1307 |
return data
|
1308 |
|