Tools for building things

ase.build.cut(atoms, a=(1, 0, 0), b=(0, 1, 0), c=None, clength=None, origo=(0, 0, 0), nlayers=None, extend=1.0, tolerance=0.01, maxatoms=None)[source]

Cuts out a cell defined by a, b, c and origo from a sufficiently repeated copy of atoms.

Typically, this function is used to create slabs of different sizes and orientations. The vectors a, b and c are in scaled coordinates and defines the returned cell and should normally be integer-valued in order to end up with a periodic structure. However, for systems with sub-translations, like fcc, integer multiples of 1/2 or 1/3 might also make sense for some directions (and will be treated correctly).

Parameters:

atoms: Atoms instance

This should correspond to a repeatable unit cell.

a: int | 3 floats

The a-vector in scaled coordinates of the cell to cut out. If integer, the a-vector will be the scaled vector from origo to the atom with index a.

b: int | 3 floats

The b-vector in scaled coordinates of the cell to cut out. If integer, the b-vector will be the scaled vector from origo to the atom with index b.

c: None | int | 3 floats

The c-vector in scaled coordinates of the cell to cut out. if integer, the c-vector will be the scaled vector from origo to the atom with index c. If None it will be along cross(a, b) converted to real space and normalised with the cube root of the volume. Note that this in general is not perpendicular to a and b for non-cubic systems. For cubic systems however, this is redused to c = cross(a, b).

clength: None | float

If not None, the length of the c-vector will be fixed to clength Angstroms. Should not be used together with nlayers.

origo: int | 3 floats

Position of origo of the new cell in scaled coordinates. If integer, the position of the atom with index origo is used.

nlayers: None | int

If nlayers is not None, the returned cell will have nlayers atomic layers in the c-direction.

extend: 1 or 3 floats

The extend argument scales the effective cell in which atoms will be included. It must either be three floats or a single float scaling all 3 directions. By setting to a value just above one, e.g. 1.05, it is possible to all the corner and edge atoms in the returned cell. This will of cause make the returned cell non-repeatable, but is very useful for visualisation.

tolerance: float

Determines what is defined as a plane. All atoms within tolerance Angstroms from a given plane will be considered to belong to that plane.

maxatoms: None | int

This option is used to auto-tune tolerance when nlayers is given for high zone axis systems. For high zone axis one needs to reduce tolerance in order to distinguise the atomic planes, resulting in the more atoms will be added and eventually MemoryError. A too small tolerance, on the other hand, might result in inproper splitting of atomic planes and that too few layers are returned. If maxatoms is not None, tolerance will automatically be gradually reduced until nlayers atomic layers is obtained, when the number of atoms exceeds maxatoms.

Example: Create an aluminium (111) slab with three layers.

>>> import ase
>>> from ase.spacegroup import crystal
>>> from ase.build.tools import cut

# First, a unit cell of Al >>> a = 4.05 >>> aluminium = crystal(‘Al’, [(0,0,0)], spacegroup=225, … cellpar=[a, a, a, 90, 90, 90])

# Then cut out the slab >>> al111 = cut(aluminium, (1,-1,0), (0,1,-1), nlayers=3)

Example: Visualisation of the skutterudite unit cell

>>> from ase.spacegroup import crystal
>>> from ase.build.tools import cut

# Again, create a skutterudite unit cell >>> a = 9.04 >>> skutterudite = crystal( … (‘Co’, ‘Sb’), … basis=[(0.25,0.25,0.25), (0.0, 0.335, 0.158)], … spacegroup=204, … cellpar=[a, a, a, 90, 90, 90])

# Then use origo to put ‘Co’ at the corners and extend to # include all corner and edge atoms. >>> s = cut(skutterudite, origo=(0.25, 0.25, 0.25), extend=1.01) >>> ase.view(s) # doctest:+SKIP

ase.build.stack(atoms1, atoms2, axis=2, cell=None, fix=0.5, maxstrain=0.5, distance=None, reorder=False, output_strained=False)[source]

Return a new Atoms instance with atoms2 stacked on top of atoms1 along the given axis. Periodicity in all directions is ensured.

The size of the final cell is determined by cell, except that the length alongh axis will be the sum of atoms1.cell[axis] and atoms2.cell[axis]. If cell is None, it will be interpolated between atoms1 and atoms2, where fix determines their relative weight. Hence, if fix equals zero, the final cell will be determined purely from atoms1 and if fix equals one, it will be determined purely from atoms2.

An ase.geometry.IncompatibleCellError exception is raised if the cells of atoms1 and atoms2 are incompatible, e.g. if the far corner of the unit cell of either atoms1 or atoms2 is displaced more than maxstrain. Setting maxstrain to None disables this check.

If distance is not None, the size of the final cell, along the direction perpendicular to the interface, will be adjusted such that the distance between the closest atoms in atoms1 and atoms2 will be equal to distance. This option uses scipy.optimize.fmin() and hence require scipy to be installed.

If reorder is True, then the atoms will be reordered such that all atoms with the same symbol will follow sequencially after each other, eg: ‘Al2MnAl10Fe’ -> ‘Al12FeMn’.

If output_strained is True, then the strained versions of atoms1 and atoms2 are returned in addition to the stacked structure.

Example: Create an Ag(110)-Si(110) interface with three atomic layers on each side.

>>> import ase
>>> from ase.spacegroup import crystal
>>> from ase.build.tools import cut, stack
>>>
>>> a_ag = 4.09
>>> ag = crystal(['Ag'], basis=[(0,0,0)], spacegroup=225,
...              cellpar=[a_ag, a_ag, a_ag, 90., 90., 90.])
>>> ag110 = cut(ag, (0, 0, 3), (-1.5, 1.5, 0), nlayers=3)
>>>
>>> a_si = 5.43
>>> si = crystal(['Si'], basis=[(0,0,0)], spacegroup=227,
...              cellpar=[a_si, a_si, a_si, 90., 90., 90.])
>>> si110 = cut(si, (0, 0, 2), (-1, 1, 0), nlayers=3)
>>>
>>> interface = stack(ag110, si110, maxstrain=1)
>>> ase.view(interface)  
>>>
# Once more, this time adjusted such that the distance between
# the closest Ag and Si atoms will be 2.3 Angstrom (requires scipy).
>>> interface2 = stack(ag110, si110,
...                    maxstrain=1, distance=2.3)   
Optimization terminated successfully.
    ...
>>> ase.view(interface2)  
ase.build.sort(atoms, tags=None)[source]

Return a new Atoms object with sorted atomic order. The default is to order according to chemical symbols, but if tags is not None, it will be used instead. A stable sorting algorithm is used.

Example:

>>> from ase.build import bulk
>>> from ase.build.tools import sort
>>> # Two unit cells of NaCl:
>>> a = 5.64
>>> nacl = bulk('NaCl', 'rocksalt', a=a) * (2, 1, 1)
>>> nacl.get_chemical_symbols()
['Na', 'Cl', 'Na', 'Cl']
>>> nacl_sorted = sort(nacl)
>>> nacl_sorted.get_chemical_symbols()
['Cl', 'Cl', 'Na', 'Na']
>>> np.all(nacl_sorted.cell == nacl.cell)
True
ase.build.rotate(atoms, a1, a2, b1, b2, rotate_cell=True, center=(0, 0, 0))[source]

Rotate atoms, such that a1 will be rotated in the direction of a2 and b1 in the direction of b2. The point at center is fixed. Use center=’COM’ to fix the center of mass. If rotate_cell is true, the cell will be rotated together with the atoms.

Note that the 000-corner of the cell is by definition fixed at origo. Hence, setting center to something other than (0, 0, 0) will rotate the atoms out of the cell, even if rotate_cell is True.

ase.build.niggli_reduce(atoms)[source]

Convert the supplied atoms object’s unit cell into its maximally-reduced Niggli unit cell. Even if the unit cell is already maximally reduced, it will be converted into its unique Niggli unit cell. This will also wrap all atoms into the new unit cell.

References:

Niggli, P. “Krystallographische und strukturtheoretische Grundbegriffe. Handbuch der Experimentalphysik”, 1928, Vol. 7, Part 1, 108-176.

Krivy, I. and Gruber, B., “A Unified Algorithm for Determining the Reduced (Niggli) Cell”, Acta Cryst. 1976, A32, 297-298.

Grosse-Kunstleve, R.W.; Sauter, N. K.; and Adams, P. D. “Numerically stable algorithms for the computation of reduced unit cells”, Acta Cryst. 2004, A60, 1-6.

ase.build.minimize_tilt(atoms, order=range(0, 3), fold_atoms=True)[source]

Minimize the tilt angles of the unit cell.

ase.build.minimize_rotation_and_translation(target, atoms)[source]

Minimize RMSD between atoms and target.

Rotate and translate atoms to best match target. Disregards rotation if PBC are found. Does not accound for changes in the cell. For more details, see:

Melander et al. J. Chem. Theory Comput., 2015, 11,1055
ase.build.find_optimal_cell_shape(cell, target_size, target_shape, lower_limit=-2, upper_limit=2, verbose=False)[source]

Returns the transformation matrix that produces a supercell corresponding to target_size unit cells with metric cell that most closely approximates the shape defined by target_shape.

Parameters:

cell: 2D array of floats

Metric given as a (3x3 matrix) of the input structure.

target_size: integer

Size of desired super cell in number of unit cells.

target_shape: str

Desired supercell shape. Can be ‘sc’ for simple cubic or ‘fcc’ for face-centered cubic.

lower_limit: int

Lower limit of search range.

upper_limit: int

Upper limit of search range.

verbose: bool

Set to True to obtain additional information regarding construction of transformation matrix.

ase.build.get_deviation_from_optimal_cell_shape(cell, target_shape='sc', norm=None)[source]

Calculates the deviation of the given cell metric from the ideal cell metric defining a certain shape. Specifically, the function evaluates the expression \(\Delta = || Q \mathbf{h} - \mathbf{h}_{target}||_2\), where \(\mathbf{h}\) is the input metric (cell) and \(Q\) is a normalization factor (norm) while the target metric \(\mathbf{h}_{target}\) (via target_shape) represent simple cubic (‘sc’) or face-centered cubic (‘fcc’) cell shapes.

Parameters:

cell: 2D array of floats

Metric given as a (3x3 matrix) of the input structure.

target_shape: str

Desired supercell shape. Can be ‘sc’ for simple cubic or ‘fcc’ for face-centered cubic.

norm: float

Specify the normalization factor. This is useful to avoid recomputing the normalization factor when computing the deviation for a series of P matrices.

ase.build.make_supercell(prim, P, *, wrap=True, order='cell-major', tol=1e-05)[source]

Generate a supercell by applying a general transformation (P) to the input configuration (prim).

The transformation is described by a 3x3 integer matrix \(\mathbf{P}\). Specifically, the new cell metric \(\mathbf{h}\) is given in terms of the metric of the input configuration \(\mathbf{h}_p\) by \(\mathbf{P h}_p = \mathbf{h}\).

Parameters:

prim: ASE Atoms object

Input configuration.

P: 3x3 integer matrix

Transformation matrix \(\mathbf{P}\).

wrap: bool

wrap in the end

order: str (default: “cell-major”)

how to order the atoms in the supercell

“cell-major”: [atom1_shift1, atom2_shift1, …, atom1_shift2, atom2_shift2, …] i.e. run first over all the atoms in cell1 and then move to cell2.

“atom-major”: [atom1_shift1, atom1_shift2, …, atom2_shift1, atom2_shift2, …] i.e. run first over atom1 in all the cells and then move to atom2. This may be the order preferred by most VASP users.

tol: float

tolerance for wrapping