Geometry tools

ase.geometry.wrap_positions(positions, cell, pbc=True, center=(0.5, 0.5, 0.5), eps=1e-07)[source]

Wrap positions to unit cell.

Returns positions changed by a multiple of the unit cell vectors to fit inside the space spanned by these vectors. See also the ase.Atoms.wrap() method.

Parameters:

positions: float ndarray of shape (n, 3)
Positions of the atoms
cell: float ndarray of shape (3, 3)
Unit cell vectors.
pbc: one or 3 bool
For each axis in the unit cell decides whether the positions will be moved along this axis.
center: three float
The positons in fractional coordinates that the new positions will be nearest possible to.
eps: float
Small number to prevent slightly negative coordinates from being wrapped.

Example:

>>> from ase.geometry import wrap_positions
>>> wrap_positions([[-0.1, 1.01, -0.5]],
...                [[1, 0, 0], [0, 1, 0], [0, 0, 4]],
...                pbc=[1, 1, 0])
array([[ 0.9 ,  0.01, -0.5 ]])
ase.geometry.complete_cell(cell)[source]

Calculate complete cell with missing lattice vectors.

Returns a new 3x3 ndarray.

ase.geometry.is_orthorhombic(cell)[source]

Check that cell only has stuff in the diagonal.

ase.geometry.orthorhombic(cell)[source]

Return cell as three box dimensions or raise ValueError.

ase.geometry.get_layers(atoms, miller, tolerance=0.001)[source]

Returns two arrays describing which layer each atom belongs to and the distance between the layers and origo.

Parameters:

miller: 3 integers
The Miller indices of the planes. Actually, any direction in reciprocal space works, so if a and b are two float vectors spanning an atomic plane, you can get all layers parallel to this with miller=np.cross(a,b).
tolerance: float
The maximum distance in Angstrom along the plane normal for counting two atoms as belonging to the same plane.

Returns:

tags: array of integres
Array of layer indices for each atom.
levels: array of floats
Array of distances in Angstrom from each layer to origo.

Example:

>>> import numpy as np
>>> from ase.spacegroup import crystal
>>> atoms = crystal('Al', [(0,0,0)], spacegroup=225, cellpar=4.05)
>>> np.round(atoms.positions, decimals=5)
array([[ 0.   ,  0.   ,  0.   ],
       [ 0.   ,  2.025,  2.025],
       [ 2.025,  0.   ,  2.025],
       [ 2.025,  2.025,  0.   ]])
>>> get_layers(atoms, (0,0,1))  
(array([0, 1, 1, 0]...), array([ 0.   ,  2.025]))
ase.geometry.find_mic(D, cell, pbc=True)[source]

Finds the minimum-image representation of vector(s) D

ase.geometry.get_duplicate_atoms(atoms, cutoff=0.1, delete=False)[source]

Get list of duplicate atoms and delete them if requested.

Identify all atoms which lie within the cutoff radius of each other. Delete one set of them if delete == True.

ase.geometry.cell_to_cellpar(cell, radians=False)[source]

Returns the cell parameters [a, b, c, alpha, beta, gamma].

Angles are in degrees unless radian=True is used.

ase.geometry.cellpar_to_cell(cellpar, ab_normal=(0, 0, 1), a_direction=None)[source]

Return a 3x3 cell matrix from cellpar=[a,b,c,alpha,beta,gamma].

Angles must be in degrees.

The returned cell is orientated such that a and b are normal to \(ab_normal\) and a is parallel to the projection of \(a_direction\) in the a-b plane.

Default \(a_direction\) is (1,0,0), unless this is parallel to \(ab_normal\), in which case default \(a_direction\) is (0,0,1).

The returned cell has the vectors va, vb and vc along the rows. The cell will be oriented such that va and vb are normal to \(ab_normal\) and va will be along the projection of \(a_direction\) onto the a-b plane.

Example:

>>> cell = cellpar_to_cell([1, 2, 4, 10, 20, 30], (0, 1, 1), (1, 2, 3))
>>> np.round(cell, 3)
array([[ 0.816, -0.408,  0.408],
       [ 1.992, -0.13 ,  0.13 ],
       [ 3.859, -0.745,  0.745]])
ase.geometry.crystal_structure_from_cell(cell, eps=0.0002, niggli_reduce=True)[source]

Return the crystal structure as a string calculated from the cell.

Supply a cell (from atoms.get_cell()) and get a string representing the crystal structure returned. Works exactly the opposite way as ase.dft.kpoints.get_special_points().

Parameters:

cell : numpy.array or list
An array like atoms.get_cell()

Returns:

crystal structure : str
‘cubic’, ‘fcc’, ‘bcc’, ‘tetragonal’, ‘orthorhombic’, ‘hexagonal’ or ‘monoclinic’
ase.geometry.distance(s1, s2, permute=True)[source]

Get the distance between two structures s1 and s2.

The distance is defined by the Frobenius norm of the spatial distance between all coordinates (see numpy.linalg.norm for the definition).

permute: minimise the distance by ‘permuting’ same elements

ase.geometry.get_angles(v1, v2, cell=None, pbc=None)[source]

Get angles formed by two lists of vectors.

calculate angle in degrees between vectors v1 and v2

Set a cell and pbc to enable minimum image convention, otherwise angles are taken as-is.

ase.geometry.get_distances(p1, p2=None, cell=None, pbc=None)[source]

Return distance matrix of every position in p1 with every position in p2

if p2 is not set, it is assumed that distances between all positions in p1 are desired. p2 will be set to p1 in this case.

Use set cell and pbc to use the minimum image convention.