Geometry tools

class ase.geometry.Cell(array, pbc=None)[source]

Parallel epipedal unit cell of up to three dimensions.

This object resembles a 3x3 array whose [i, j]-th element is the jth Cartesian coordinate of the ith unit vector.

Cells of less than three dimensions are represented by placeholder unit vectors that are zero.

Create cell.

Parameters:

array: 3x3 arraylike object

The three cell vectors.

pbc: None or 3 booleans

For each cell vector, whether the system is periodic in the direction of that cell vector. If not given, the cell will be periodic along directions with nonzero cell vectors.

angles()[source]

Return an array with the three angles alpha, beta, and gamma.

classmethod ascell(cell)[source]

Return argument as a Cell object. See ase.cell.Cell.new().

A new Cell object is created if necessary.

bandpath(path=None, npoints=None, density=None, special_points=None, eps=0.0002)[source]

Build a BandPath for this cell.

If special points are None, determine the Bravais lattice of this cell and return a suitable Brillouin zone path with standard special points.

If special special points are given, interpolate the path directly from the available data.

Parameters:

path: string

String of special point names defining the path, e.g. ‘GXL’.

npoints: int

Number of points in total. Note that at least one point is added for each special point in the path.

density: float

density of kpoints along the path in Å⁻¹.

special_points: dict

Dictionary mapping special points to scaled kpoint coordinates. For example {'G': [0, 0, 0], 'X': [1, 0, 0]}.

eps: float

Tolerance for determining Bravais lattice.

Example

>>> cell = Cell.fromcellpar([4, 4, 4, 60, 60, 60])
>>> cell.bandpath('GXW', npoints=20)
BandPath(path='GXW', cell=[3x3], special_points={GKLUWX}, kpts=[20x3])
cartesian_positions(scaled_positions)[source]

Calculate Cartesian positions from scaled positions.

cellpar(radians=False)[source]

Get cell lengths and angles of this cell.

See also ase.geometry.cell.cell_to_cellpar().

complete()[source]

Convert missing cell vectors into orthogonal unit vectors.

copy()[source]

Return a copy of this cell.

classmethod fromcellpar(cellpar, ab_normal=(0, 0, 1), a_direction=None, pbc=None)[source]

Return new Cell from cell lengths and angles.

See also cellpar_to_cell().

get_bravais_lattice(eps=0.0002)[source]

Return BravaisLattice for this cell:

>>> cell = Cell.fromcellpar([4, 4, 4, 60, 60, 60])
>>> print(cell.get_bravais_lattice())
FCC(a=5.65685)

Note

The Bravais lattice object follows the AFlow conventions. cell.get_bravais_lattice().tocell() may differ from the original cell by a permutation or other operation which maps it to the AFlow convention. For example, the orthorhombic lattice enforces a < b < c.

To build a bandpath for a particular cell, use ase.cell.Cell.bandpath() instead of this method. This maps the kpoints back to the original input cell.

lengths()[source]

Return the length of each lattice vector as an array.

minkowski_reduce()[source]

Minkowski-reduce this cell, returning new cell and mapping.

See also ase.geometry.minkowski_reduction.minkowski_reduce().

classmethod new(cell=None, pbc=None)[source]

Create new cell from any parameters.

If cell is three numbers, assume three lengths with right angles.

If cell is six numbers, assume three lengths, then three angles.

If cell is 3x3, assume three cell vectors.

niggli_reduce(eps=1e-05)[source]

Niggli reduce this cell, returning a new cell and mapping.

See also ase.build.tools.niggli_reduce_cell().

oldbandpath(path=None, npoints=None, density=None, eps=0.0002)[source]

Legacy implementation, please ignore.

property orthorhombic

Return whether this cell is represented by a diagonal matrix.

property rank

“Return the dimension of the cell.

Equal to the number of nonzero lattice vectors.

reciprocal()[source]

Get reciprocal lattice as a 3x3 array.

Does not include factor of 2 pi.

scaled_positions(positions)[source]

Calculate scaled positions from Cartesian positions.

The scaled positions are the positions given in the basis of the cell vectors. For the purpose of defining the basis, cell vectors that are zero will be replaced by unit vectors as per complete().

property volume

Get the volume of this cell.

If there are less than 3 lattice vectors, return 0.

ase.geometry.wrap_positions(positions, cell, pbc=True, center=(0.5, 0.5, 0.5), pretty_translation=False, 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.

pretty_translation: bool

Translates atoms such that fractional coordinates are minimized.

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:

cellnumpy.array or list

An array like atoms.get_cell()

Returns:

crystal structurestr

‘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.

ase.geometry.minkowski_reduce(B)[source]

Calculate a Minkowski-reduced lattice basis. The reduced basis has the shortest possible vector lengths and has norm(a) <= norm(b) <= norm(c).

Implements the method described in:

Low-dimensional Lattice Basis Reduction Revisited Nguyen, Phong Q. and Stehlé, Damien, ACM Trans. Algorithms 5(4) 46:1–46:48, 2009 https://doi.org/10.1145/1597036.1597050

Parameters:

B: array

The lattice basis to reduce (in row-vector format).

Returns:

R: array

The reduced lattice basis.

H: array

The unimodular matrix transformation (R = H @ B).

Analysis tools

Provides the class Analysis for structural analysis of any Atoms object or list thereof (trajectories).

Example:

>>> import numpy as np
>>> from ase.build import molecule
>>> from ase.geometry.analysis import Analysis
>>> mol = molecule('C60')
>>> ana = Analysis(mol)
>>> CCBonds = ana.get_bonds('C', 'C', unique=True)
>>> CCCAngles = ana.get_angles('C', 'C', 'C', unique=True)
>>> print("There are {} C-C bonds in C60.".format(len(CCBonds[0])))
>>> print("There are {} C-C-C angles in C60.".format(len(CCCAngles[0])))
>>> CCBondValues = ana.get_values(CCBonds)
>>> CCCAngleValues = ana.get_values(CCCAngles)
>>> print("The average C-C bond length is {}.".format(np.average(CCBondValues)))
>>> print("The average C-C-C angle is {}.".format(np.average(CCCAngleValues)))

The Analysis class provides a getter and setter for the images. This allows you to use the same neighbourlist for different images, e.g. to analyze two MD simulations at different termperatures but constant bonding patterns. Using this approach saves the time to recalculate all bonds, angles and dihedrals and therefore speeds up your analysis.

Using the Analysis.clear_cache() function allows you to clear the calculated matrices/lists to reduce your memory usage.

The entire class can be used with few commands:

  • To retrieve tuples of bonds/angles/dihedrals (they are calculated the first time they are accessed) use instance.all_xxx where xxx is one of bonds/angles/dihedrals.

  • If you only want those one-way (meaning e.g. not bonds i-j and j-i but just i-j) use instance.unique_xxx.

  • To get selected bonds/angles/dihedrals use instance.get_xxx(A,B,...), see the API section for details on which arguments you can pass.

  • To get the actual value of a bond/angle/dihedral use instance.get_xxx_value(tuple).

  • To get a lot of bond/angle/dihedral values at once use Analysis.get_values().

  • There is also a wrapper to get radial distribution functions Analysis.get_rdf().

The main difference between properties (getters) and functions here is, that getters provide data that is cached. This means that getting information from Analysis.all_bonds more than once is instantanious, since the information is cached in Analysis._cache. If you call any Analysis.get_xxx() the information is calculated from the cached data, meaning each call will take the same amount of time.

API:

class ase.geometry.analysis.Analysis(images, nl=None, **kwargs)[source]

Analysis class

Parameters for initialization:

images: Atoms object or list of such

Images to analyze.

nl: None, NeighborList object or list of such

Neighborlist(s) for the given images. One or nImages, depending if bonding pattern changes or is constant. Using one Neigborlist greatly improves speed.

kwargs: options, dict

Arguments for constructing NeighborList object if nl is None.

The choice of bothways=True for the NeighborList object will not influence the amount of bonds/angles/dihedrals you get, all are reported in both directions. Use the unique-labeled properties to get lists without duplicates.

property adjacency_matrix

The adjacency/connectivity matrix.

If not already done, build a list of adjacency matrices for all nl.

No setter or deleter, only getter

property all_angles

All angles

A list with indices of atoms in angles for each neighborlist in self. Atom i forms an angle to the atoms inside the tuples in result[i]: i – result[i][x][0] – result[i][x][1] where x is in range(number of angles from i). See also unique_angles.

No setter or deleter, only getter

property all_bonds

All Bonds.

A list with indices of bonded atoms for each neighborlist in self. Atom i is connected to all atoms inside result[i]. Duplicates from PBCs are removed. See also unique_bonds.

No setter or deleter, only getter

property all_dihedrals

All dihedrals

Returns a list with indices of atoms in dihedrals for each neighborlist in this instance. Atom i forms a dihedral to the atoms inside the tuples in result[i]: i – result[i][x][0] – result[i][x][1] – result[i][x][2] where x is in range(number of dihedrals from i). See also unique_dihedrals.

No setter or deleter, only getter

clear_cache()[source]

Delete all cached information.

property distance_matrix

The distance matrix.

If not already done, build a list of distance matrices for all nl. See ase.neighborlist.get_distance_matrix().

No setter or deleter, only getter

get_angle_value(imIdx, idxs, mic=True, **kwargs)[source]

Get angle.

Parameters:

imIdx: int

Index of Image to get value from.

idxs: tuple or list of integers

Get angle between atoms idxs[0]-idxs[1]-idxs[2].

mic: bool

Passed on to ase.Atoms.get_angle() for retrieving the value, defaults to True. If the cell of the image is correctly set, there should be no reason to change this.

kwargs: options or dict

Passed on to ase.Atoms.get_angle().

Returns:

return: float

Value returned by image.get_angle.

get_angles(A, B, C, unique=True)[source]

Get angles from given elements A-B-C.

Parameters:

A, B, C: str

Get Angles between elements A, B and C. B will be the central atom.

unique: bool

Return the angles both ways or just one way (A-B-C and C-B-A or only A-B-C)

Returns:

return: list of lists of tuples

return[imageIdx][atomIdx][angleI], each tuple starts with atomIdx.

Use get_values() to convert the returned list to values.

get_bond_value(imIdx, idxs, mic=True, **kwargs)[source]

Get bond length.

Parameters:

imIdx: int

Index of Image to get value from.

idxs: tuple or list of integers

Get distance between atoms idxs[0]-idxs[1].

mic: bool

Passed on to ase.Atoms.get_distance() for retrieving the value, defaults to True. If the cell of the image is correctly set, there should be no reason to change this.

kwargs: options or dict

Passed on to ase.Atoms.get_distance().

Returns:

return: float

Value returned by image.get_distance.

get_bonds(A, B, unique=True)[source]

Get bonds from element A to element B.

Parameters:

A, B: str

Get Bonds between elements A and B

unique: bool

Return the bonds both ways or just one way (A-B and B-A or only A-B)

Returns:

return: list of lists of tuples

return[imageIdx][atomIdx][bondI], each tuple starts with atomIdx.

Use get_values() to convert the returned list to values.

get_dihedral_value(imIdx, idxs, mic=True, **kwargs)[source]

Get dihedral.

Parameters:

imIdx: int

Index of Image to get value from.

idxs: tuple or list of integers

Get angle between atoms idxs[0]-idxs[1]-idxs[2]-idxs[3].

mic: bool

Passed on to ase.Atoms.get_dihedral() for retrieving the value, defaults to True. If the cell of the image is correctly set, there should be no reason to change this.

kwargs: options or dict

Passed on to ase.Atoms.get_dihedral().

Returns:

return: float

Value returned by image.get_dihedral.

get_dihedrals(A, B, C, D, unique=True)[source]

Get dihedrals A-B-C-D.

Parameters:

A, B, C, D: str

Get Dihedralss between elements A, B, C and D. B-C will be the central axis.

unique: bool

Return the dihedrals both ways or just one way (A-B-C-D and D-C-B-A or only A-B-C-D)

Returns:

return: list of lists of tuples

return[imageIdx][atomIdx][dihedralI], each tuple starts with atomIdx.

Use get_values() to convert the returned list to values.

get_rdf(rmax, nbins, imageIdx=None, elements=None, return_dists=False)[source]

Get RDF.

Wrapper for ase.ga.utilities.get_rdf() with more selection possibilities.

Parameters:

rmax: float

Maximum distance of RDF.

nbins: int

Number of bins to devide RDF.

imageIdx: int/slice/None

Images to analyze, see _get_slice() for details.

elements: str/int/list/tuple

Make partial RDFs.

If elements is None, a full RDF is calculated. If elements is an integer or a list/tuple of integers, only those atoms will contribute to the RDF (like a mask). If elements is a string or a list/tuple of strings, only Atoms of those elements will contribute.

Returns:

return: list of lists / list of tuples of lists

If return_dists is True, the returned tuples contain (rdf, distances). Otherwise only rdfs for each image are returned.

get_values(inputList, imageIdx=None, mic=True, **kwargs)[source]

Get Bond/Angle/Dihedral values.

Parameters:

inputList: list of lists of tuples

Can be any list provided by get_bonds(), get_angles() or get_dihedrals().

imageIdx: integer or slice

The images from images to be analyzed. If None, all frames will be analyzed. See _get_slice() for details.

mic: bool

Passed on to Atoms for retrieving the values, defaults to True. If the cells of the images are correctly set, there should be no reason to change this.

kwargs: options or dict

Passed on to the Atoms classes functions for retrieving the values.

Returns:

return: list of lists of floats

return[imageIdx][valueIdx]. Has the same shape as the inputList, instead of each tuple there is a float with the value this tuple yields.

The type of value requested is determined from the length of the tuple inputList[0][0]. The methods from the Atoms class are used.

property images

Images.

Set during initialization but can also be set later.

property nImages

Number of Images in this instance.

Cannot be set, is determined automatically.

property nl

Neighbor Lists in this instance.

Set during initialization.

No setter or deleter, only getter

property unique_angles

Get Unique Angles.

all_angles i-j-k without k-j-i.

property unique_bonds

Get Unique Bonds.

all_bonds i-j without j-i. This is the upper triangle of the connectivity matrix (i,j), \(i < j\)

property unique_dihedrals

Get Unique Dihedrals.

all_dihedrals i-j-k-l without l-k-j-i.