The Atoms object is a collection of atoms. Here is how to define a CO molecule:
from ase import Atoms
d = 1.1
co = Atoms('CO', positions=[(0, 0, 0), (0, 0, d)])
Here, the first argument specifies the type of the atoms and we used the positions keywords to specify their positions. Other possible keywords are: numbers, tags, momenta, masses, magmoms and charges.
Here is how you could make an infinite gold wire with a bond length of 2.9 Å:
from ase import Atoms
d = 2.9
L = 10.0
wire = Atoms('Au',
positions=[(0, L / 2, L / 2)],
cell=(d, L, L),
pbc=(1, 0, 0))
Here, two more optional keyword arguments were used:
You can also use the following methods to work with the unit cell and the boundary conditions: set_pbc(), set_cell(), get_cell(), and get_pbc().
Like with a single Atom the properties of a collection of atoms can be accessed and changed with get- and set-methods. For example the positions of the atoms can be addressed as
>>> a = Atoms('N3', [(0, 0, 0), (1, 0, 0), (0, 0, 1)])
>>> a.get_positions()
array([[ 0., 0., 0.],
[ 1., 0., 0.],
[ 0., 0., 1.]])
>>> a.set_positions([(2, 0, 0), (0, 2, 2), (2, 2, 0)])
>>> a.get_positions()
array([[ 2., 0., 0.],
[ 0., 2., 2.],
[ 2., 2., 0.]])
Here is the full list of the get/set methods operating on all the atoms at once. The get methods return an array of quantities, one for each atom; the set methods take similar arrays. E.g. get_positions() return N * 3 numbers, get_atomic_numbers() return N integers.
These methods return copies of the internal arrays, it is thus safe to modify the returned arrays.
There are also a number of get/set methods that operate on quantities common to all the atoms or defined for the collection of atoms:
The Atoms object holds a unit cell which by default is the 3x3 unit matrix as can be seen from
>>> a.get_cell()
array([[ 1., 0., 0.],
[ 0., 1., 0.],
[ 0., 0., 1.]])
The cell can be defined or changed using the set_cell() method. Changing the unit cell does per default not move the atoms:
>>> a.set_cell(2 * identity(3))
>>> a.get_cell()
array([[ 2., 0., 0.],
[ 0., 2., 0.],
[ 0., 0., 2.]])
>>> a.get_positions()
array([[ 2., 0., 0.],
[ 1., 1., 0.],
[ 2., 2., 0.]])
However if we set scale_atoms=True the atomic positions are scaled with the unit cell:
>>> a.set_cell(identity(3), scale_atoms=True)
>>> a.get_positions()
array([[ 1. , 0. , 0. ],
[ 0.5, 0.5, 0. ],
[ 1. , 1. , 0. ]])
The set_pbc() method specifies whether periodic boundary conditions are to be used in the directions of the three vectors of the unit cell. A slab calculation with periodic boundary conditions in x and y directions and free boundary conditions in the z direction is obtained through
>>> a.set_pbc((True, True, False))
It is also possible to work directly with the attributes positions, numbers, pbc and cell. Here we change the position of the 2nd atom (which has count number 1 because Python starts counting at zero) and the type of the first atom:
>>> a.positions[1] = (1, 1, 0)
>>> a.get_positions()
array([[2., 0., 0.],
[1., 1., 0.],
[2., 2., 0.]])
>>> a.positions
array([[2., 0., 0.],
[1., 1., 0.],
[2., 2., 0.]])
>>> a.numbers
array([7, 7, 7])
>>> a.numbers[0] = 13
>>> a.get_chemical_symbols()
['Al', 'N', 'N']
Check for periodic boundary conditions:
>>> a.pbc # equivalent to a.get_pbc()
array([False, False, False], dtype=bool)
>>> a.pbc.any()
False
>>> a.pbc[2] = 1
>>> a.pbc
array([False, False, True], dtype=bool)
A calculator can be attached to the atoms with the purpose of calculating energies and forces on the atoms. ASE works with many different ase.calculators.
A calculator object calc is attached to the atoms like this:
>>> a.set_calculator(calc)
After the calculator has been appropriately setup the energy of the atoms can be obtained through
>>> a.get_potential_energy()
The term “potential energy” here means for example the total energy of a DFT calculation, which includes both kinetic, electrostatic, and exchange-correlation energy for the electrons. The reason it is called potential energy is that the atoms might also have a kinetic energy (from the moving nuclei) and that is obtained with
>>> a.get_kinetic_energy()
In case of a DFT calculator, it is up to the user to check exactly what the get_potential_energy() method returns. For example it may be the result of a calculation with a finite temperature smearing of the occupation numbers extrapolated to zero temperature. More about this can be found for the different ase.calculators.
The following methods can only be called if a calculator is present:
Not all of these methods are supported by all calculators.
method | example |
+ | wire2 = wire + co |
+=, extend() | wire += co wire.extend(co) |
append() | wire.append(Atom('H')) |
* | wire3 = wire * (3, 1, 1) |
*=, repeat() | wire *= (3, 1, 1) wire.repeat((3, 1, 1)) |
len | len(co) |
del | del wire3[0] del wire3[[1,3]] |
pop() | oxygen = wire2.pop() |
Note that the del method can be used with the more powerful numpy-style indexing, as in the second example above. This can be combined with python list comprehension in order to selectively delete atoms within an ASE Atoms object. For example, the below code creates an ethanol molecule and subsequently strips all the hydrogen atoms from it:
from ase.structure import molecule
atoms = molecule('CH3CH2OH')
del atoms[[atom.index for atom in atoms if atom.symbol=='H']]
Atoms object.
The Atoms object can represent an isolated molecule, or a periodically repeated structure. It has a unit cell and there may be periodic boundary conditions along any of the three unit cell axes.
Information about the atoms (atomic numbers and position) is stored in ndarrays. Optionally, there can be information about tags, momenta, masses, magnetic moments and charges.
In order to calculate energies, forces and stresses, a calculator object has to attached to the atoms object.
Parameters:
Dictionary of key-value pairs with additional information about the system. The following keys may be used by ase:
- spacegroup: Spacegroup instance
- unit_cell: ‘conventional’ | ‘primitive’ | int | 3 ints
- adsorbate_info:
Items in the info attribute survives copy and slicing and can be store to and retrieved from trajectory files given that the key is a string, the value is picklable and, if the value is a user-defined object, its base class is importable. One should not make any assumptions about the existence of keys.
Examples:
These three are equivalent:
>>> d = 1.104 # N2 bondlength
>>> a = Atoms('N2', [(0, 0, 0), (0, 0, d)])
>>> a = Atoms(numbers=[7, 7], positions=[(0, 0, 0), (0, 0, d)])
>>> a = Atoms([Atom('N', (0, 0, 0)), Atom('N', (0, 0, d)])
FCC gold:
>>> a = 4.05 # Gold lattice constant
>>> b = a / 2
>>> fcc = Atoms('Au',
... cell=[(0, b, b), (b, 0, b), (b, b, 0)],
... pbc=True)
Hydrogen wire:
>>> d = 0.9 # H-H distance
>>> L = 7.0
>>> h = Atoms('H', positions=[(0, L / 2, L / 2)],
... cell=(d, L, L),
... pbc=(1, 0, 0))
Calculator object.
Attribute for direct manipulation of the unit cell.
Center atoms in unit cell.
Centers the atoms in the unit cell, so there is the same amount of vacuum on all sides.
Constraints of the atoms.
Modify atoms interactively through ase-gui viewer.
Conflicts leading to undesirable behaviour might arise when matplotlib has been pre-imported with certain incompatible backends and while trying to use the plot feature inside the interactive ag. To circumvent, please set matplotlib.use(‘gtk’) before calling this method.
Return distances of all of the atoms with all of the atoms.
Use mic=True to use the Minimum Image Convention.
Get angle formed by three atoms.
calculate angle between the vectors list[1]->list[0] and list[1]->list[2], where list contains the atomic indexes in question.
Get an array.
Returns a copy unless the optional argument copy is false.
Get the center of mass.
If scaled=True the center of mass in scaled coordinates is returned.
Get the chemial formula as a string based on the chemical symbols.
Parameters:
There are three different modes available:
‘all’: The list of chemical symbols are contracted to at string, e.g. [‘C’, ‘H’, ‘H’, ‘H’, ‘O’, ‘H’] becomes ‘CHHHOH’.
‘reduce’: The same as ‘all’ where repeated elements are contracted to a single symbol and a number, e.g. ‘CHHHOCHHH’ is reduced to ‘CH3OCH3’.
‘hill’: The list of chemical symbols are contracted to a string following the Hill notation (alphabetical order with C and H first), e.g. ‘CHHHOCHHH’ is reduced to ‘C2H6O’ and ‘SOOHOHO’ to ‘H2O4S’. This is default.
Calculate dihedral angle.
Calculate dihedral angle between the vectors list[0]->list[1] and list[2]->list[3], where list contains the atomic indexes in question.
Calculate the electric dipole moment for the atoms object.
Only available for calculators which has a get_dipole_moment() method.
Return distance between two atoms.
Use mic=True to use the Minimum Image Convention. vector=True gives the distance vector (from a0 to a1).
Return distances of atom No.i with a list of atoms.
Use mic=True to use the Minimum Image Convention. vector=True gives the distance vector (from a to self[indices]).
Calculate atomic forces.
Ask the attached calculator to calculate the forces and apply constraints. Use apply_constraint=False to get the raw forces.
Get the moments of inertia along the principal axes.
The three principal moments of inertia are computed from the eigenvalues of the symmetric inertial tensor. Periodic boundary conditions are ignored. Units of the moments of inertia are amu*angstrom**2.
Returns the number of atoms.
Equivalent to len(atoms) in the standard ASE Atoms class.
Get array of positions. If wrap==True, wraps atoms back into unit cell.
Calculate the potential energies of all the atoms.
Only available with calculators supporting per-atom energies (e.g. classical potentials).
Calculate potential energy.
Ask the attached calculator to calculate the potential energy and apply constraints. Use apply_constraint=False to get the raw forces.
When supported by the calculator, either the energy extrapolated to zero Kelvin or the energy consistent with the forces (the free energy) can be returned.
Get the three reciprocal lattice vectors as a 3x3 ndarray.
Note that the commonly used factor of 2 pi for Fourier transforms is not included here.
Get positions relative to unit cell.
If wrap is True, atoms outside the unit cell will be wrapped into the cell in those directions with periodic boundary conditions so that the scaled coordinates are between zero and one.
Calculate stress tensor.
Returns an array of the six independent components of the symmetric stress tensor, in the traditional Voigt order (xx, yy, zz, yz, xz, xy) or as a 3x3 matrix. Default is Voigt order.
Calculate the stress-tensor of all the atoms.
Only available with calculators supporting per-atom energies and stresses (e.g. classical potentials). Even for such calculators there is a certain arbitrariness in defining per-atom stresses.
Get integer array of tags.
Check for existence of array.
name must be one of: ‘tags’, ‘momenta’, ‘masses’, ‘magmoms’, ‘charges’.
Add new array.
If shape is not None, the shape of a will be checked.
Attribute for direct manipulation of the atomic numbers.
Attribute for direct manipulation of the periodic boundary condition flags.
Attribute for direct manipulation of the positions.
Randomly displace atoms.
This method adds random displacements to the atomic positions, taking a possible constraint into account. The random numbers are drawn from a normal distribution of standard deviation stdev.
For a parallel calculation, it is important to use the same seed on all processors!
Create new repeated atoms object.
The rep argument should be a sequence of three positive integers like (2,3,1) or a single integer (r) equivalent to (r,r,r).
Rotate atoms based on a vector and an angle, or two vectors.
Parameters:
Examples:
Rotate 90 degrees around the z-axis, so that the x-axis is rotated into the y-axis:
>>> a = pi / 2
>>> atoms.rotate('z', a)
>>> atoms.rotate((0, 0, 1), a)
>>> atoms.rotate('-z', -a)
>>> atoms.rotate((0, 0, a))
>>> atoms.rotate('x', 'y')
Rotate dihedral angle.
Complementing the two routines above: rotate a group by a predefined dihedral angle, starting from its current configuration
Rotate atoms via Euler angles.
See e.g http://mathworld.wolfram.com/EulerAngles.html for explanation.
Parameters:
Set angle formed by three atoms.
Sets the angle between vectors list[1]->list[0] and list[1]->list[2].
Same usage as in set_dihedral.
Update array.
If shape is not None, the shape of a will be checked. If a is None, then the array is deleted.
Set unit cell vectors.
Parameters:
Examples:
Two equivalent ways to define an orthorhombic cell:
>>> a.set_cell([a, b, c])
>>> a.set_cell([(a, 0, 0), (0, b, 0), (0, 0, c)])
FCC unit cell:
>>> a.set_cell([(0, b, b), (b, 0, b), (b, b, 0)])
Apply one or more constrains.
The constraint argument must be one constraint object or a list of constraint objects.
Set the dihedral angle between vectors list[0]->list[1] and list[2]->list[3] by changing the atom indexed by list[3] if mask is not None, all the atoms described in mask (read: the entire subgroup) are moved. Alternatively to the mask, the indices of the atoms to be rotated can be supplied.
example: the following defines a very crude ethane-like molecule and twists one half of it by 30 degrees.
>>> atoms = Atoms('HHCCHH', [[-1, 1, 0], [-1, -1, 0], [0, 0, 0],
[1, 0, 0], [2, 1, 0], [2, -1, 0]])
>>> atoms.set_dihedral([1,2,3,4],7*pi/6,mask=[0,0,0,1,1,1])
Set the distance between two atoms.
Set the distance between atoms a0 and a1 to distance. By default, the center of the two atoms will be fixed. Use fix=0 to fix the first atom, fix=1 to fix the second atom and fix=0.5 (default) to fix the center of the bond.
Set the initial magnetic moments.
Use either one or three numbers for every atom (collinear or non-collinear spins).
Set atomic masses.
The array masses should contain a list of masses. In case the masses argument is not given or for those elements of the masses list that are None, standard values are set.
Set tags for all atoms. If only one tag is supplied, it is applied to all atoms.
Translate atomic positions.
The displacement argument can be a float an xyz vector or an nx3 array (where n is the number of atoms).
Wrap positions to unit cell.
Parameters:
See also the ase.utils.geometry.wrap_positions() function. Example:
>>> a = Atoms('H',
... [[-0.1, 1.01, -0.5]],
... cell=[[1, 0, 0], [0, 1, 0], [0, 0, 4]],
... pbc=[1, 1, 0])
>>> a.wrap()
>>> a.positions
array([[ 0.9 , 0.01, -0.5 ]])