geometry__ASE_tosee
esta.geometry__ASE_tosee
¶
Utility tools for atoms/geometry manipulations. - convenient creation of slabs and interfaces of different orientations. - detection of duplicate atoms / atoms within cutoff radius
translate_pretty(fractional, pbc)
¶
Translates atoms such that fractional positions are minimized.
wrap_positions(positions, cell, pbc=True, center=(0.5, 0.5, 0.5), pretty_translation=False, eps=1e-07)
¶
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
:meth: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 ]])
get_layers(atoms, miller, tolerance=0.001)
¶
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)) # doctest: +ELLIPSIS (array([0, 1, 1, 0]...), array([ 0. , 2.025]))
naive_find_mic(v, cell)
¶
Finds the minimum-image representation of vector(s) v. Safe to use for (pbc.all() and (norm(v_mic) < 0.5 * min(cell.lengths()))). Can otherwise fail for non-orthorhombic cells. Described in: W. Smith, "The Minimum Image Convention in Non-Cubic MD Cells", 1989, http://citeseerx.ist.psu.edu/viewdoc/summary?doi=10.1.1.57.1696.
general_find_mic(v, cell, pbc=True)
¶
Finds the minimum-image representation of vector(s) v. Using the Minkowski reduction the algorithm is relatively slow but safe for any cell.
find_mic(v, cell, pbc=True)
¶
Finds the minimum-image representation of vector(s) v using either one of two find mic algorithms depending on the given cell, v and pbc.
conditional_find_mic(vectors, cell, pbc)
¶
Return list of vector arrays and corresponding list of vector lengths for a given list of vector arrays. The minimum image convention is applied if cell and pbc are set. Can be used like a simple version of get_distances.
get_angles(v0, v1, cell=None, pbc=None)
¶
Get angles formed by two lists of vectors.
Calculate angle in degrees between vectors v0 and v1
Set a cell and pbc to enable minimum image convention, otherwise angles are taken as-is.
get_angles_derivatives(v0, v1, cell=None, pbc=None)
¶
Get derivatives of angles formed by two lists of vectors (v0, v1) w.r.t. Cartesian coordinates in degrees.
Set a cell and pbc to enable minimum image convention, otherwise derivatives of angles are taken as-is.
There is a singularity in the derivatives for sin(angle) -> 0 for which a ZeroDivisionError is raised.
Derivative output format: [[dx_a0, dy_a0, dz_a0], [...], [..., dz_a2].
get_dihedrals(v0, v1, v2, cell=None, pbc=None)
¶
Get dihedral angles formed by three lists of vectors.
Calculate dihedral angle (in degrees) between the vectors a0->a1, a1->a2 and a2->a3, written as v0, v1 and v2.
Set a cell and pbc to enable minimum image convention, otherwise angles are taken as-is.
get_dihedrals_derivatives(v0, v1, v2, cell=None, pbc=None)
¶
Get derivatives of dihedrals formed by three lists of vectors (v0, v1, v2) w.r.t Cartesian coordinates in degrees.
Set a cell and pbc to enable minimum image convention, otherwise dihedrals are taken as-is.
Derivative output format: [[dx_a0, dy_a0, dz_a0], ..., [..., dz_a3]].
get_distances(p1, p2=None, cell=None, pbc=None)
¶
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.
get_distances_derivatives(v0, cell=None, pbc=None)
¶
Get derivatives of distances for all vectors in v0 w.r.t. Cartesian coordinates in Angstrom.
Set cell and pbc to use the minimum image convention.
There is a singularity for distances -> 0 for which a ZeroDivisionError is raised. Derivative output format: [[dx_a0, dy_a0, dz_a0], [dx_a1, dy_a1, dz_a1]].
get_duplicate_atoms(atoms, cutoff=0.1, delete=False)
¶
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.
permute_axes(atoms, permutation)
¶
Permute axes of unit cell and atom positions. Considers only cell and atomic positions. Other vector quantities such as momenta are not modified.