pyxtal.operations module¶
Module for generating and analyzing transformation operations. Several functions for working with matrices are provided. The class OperationAnalyzer allows for comparison between pymatgen.core.operations.SymmOp objects, and can be used to identify conjugate operations. The orientation class can be used to identify degrees of freedom for molecules in Wyckoff positions with certain symmetry constraints.
-
class
pyxtal.operations.
OperationAnalyzer
(op)[source]¶ Bases:
pymatgen.core.operations.SymmOp
Class for comparing operations. Stores rotation axis and angle, as well as the type and order of operation (identity, inversion, rotation, or rotoinversion). By default, takes a SymmOp as argument. This information can be accessed by calling print(object). The class methods is_conjugate and are_conjugate can be used to determine if two operations are conjugate to each other. That is, whether they represent the same angle of rotation and are either both inverted or both uninverted.
- Note: rotoinversions with odd-order rotational parts will have an over-all
- even order. For example, the order of (-3) is 6.
Note: reflections are treated as rotoinversions of order 2.
Parameters: SymmOp – a pymatgen.core.structure.SymmOp object to analyze -
affine_matrix
= None¶ The 4x4 affine matrix of the op
-
are_conjugate
(op2)[source]¶ Returns whether or not two operations are conjugate (the same operation in a different reference frame). Rotations with the same order will not always return True. For example, a 5/12 and 1/12 rotation will not be considered conjugate.
Parameters: - op1 – a SymmOp or OperationAnalyzer object
- op2 – a SymmOp or OperationAnalyzer object to compare with op1
Returns: True if op2 is conjugate to op1, and False otherwise
-
det
= None¶ The determinant of self.m
-
is_conjugate
(op2)[source]¶ Returns whether or not another operation is conjugate (the same operation in a different reference frame). Rotations with the same order will not always return True. For example, a 5/12 and 1/12 rotation will not be considered conjugate.
Parameters: op2 – a SymmOp or OperationAnalyzer object to compare with Returns: True if op2 is conjugate to self.op, and False otherwise
-
m
= None¶ The 3x3 rotation (or rotoinversion) matrix, which ignores the translational part of self.op
-
op
= None¶ The original SymmOp object being analyzed
-
order
= None¶ The order of the operation. This is the number of times the operation must be applied consecutively to return to the identity operation. If no integer number if found, we set this to ‘irrational’.
-
rotation_order
= None¶ The order of the rotational (non-inversional) part of the operation. Must be used in conjunction with self.order to determine the properties of the operation.
-
tol
= None¶ The numerical tolerance associated with self.op
-
type
= None¶ The type of operation. Is one of ‘identity’, ‘inversion’, ‘rotation’, or ‘rotoinversion’.
-
pyxtal.operations.
aa2matrix
(axis, angle, radians=True, random=False)[source]¶ Given an axis and an angle, return a 3x3 rotation matrix. Based on: https://en.wikipedia.org/wiki/Rotation_matrix#Axis_and_angle
Parameters: - axis – a vector about which to perform a rotation
- angle – the angle of rotation
- radians – whether the supplied angle is in radians (True) or in degrees (False)
- random – whether or not to choose a random rotation matrix. If True, the axis and angle are ignored, and a random orientation is generated
Returns: a 3x3 numpy array representing a rotation matrix
-
pyxtal.operations.
angle
(v1, v2, radians=True)[source]¶ Calculate the angle (in radians) between two vectors.
Parameters: - v1 – a 1x3 vector
- v2 – a 1x3 vector
- radians – whether to return angle in radians (default) or degrees
Returns: the angle in radians between the two vectors
-
pyxtal.operations.
apply_ops
(coord, ops)[source]¶ Apply a list of SymmOps to a single 3-vector and return an array of the generated vectors. This is the inverse of SymmOp.operate_multi.
Parameters: - coord – a 3-vector (list or numpy array)
- ops – a list, tuple, or array of SymmOp objects
Returns: an np array of floating-point 3-vectors
-
pyxtal.operations.
apply_ops_diagonal
(coords, ops)[source]¶ Given a list of coordinates and SymmOps, apply the ith op to the ith coord and return the list of transformed coordinates
Parameters: coords – a list or array of 3-vectors Returns: a transformed numpy array of 3-vectors
-
pyxtal.operations.
are_equal
(op1, op2, PBC=[1, 1, 1], rtol=0.001, atol=0.001)[source]¶ Check whether two SymmOp objects are equal up to some numerical tolerance. Allows for optional consideration of periodic boundary conditions. This option is useful for handling the periodicity of crystals.
Parameters: - op1 – a SymmOp object
- op2 – another SymmOp object
- PBC – A periodic boundary condition list.
- rtol – the relative numerical tolerance for equivalence
- atol – the absolute numerical tolerance for equivalence
Returns: True if op1 and op2 are equivalent, False otherwise
-
pyxtal.operations.
check_distance
(coord1, coord2, species1, species2, lattice, PBC=[1, 1, 1], tm=<pyxtal.tolerance.Tol_matrix object>, d_factor=1.0)[source]¶ Check the distances between two atom set. Only distances between points from different sets are checked.
Parameters: - coord1 – a list of fractional coordinates e.g. [[.1,.6,.4] [.3,.8,.2]]
- coord2 – a list of new fractional coordinates e.g. [[.7,.8,.9], [.4,.5,.6]]
- species1 – a list of atomic species or numbers for coord1
- species2 – a list of atomic species or numbers for coord2
- lattice – matrix describing the unit cell vectors
- PBC – A periodic boundary condition list, where 1 means periodic, 0 means not periodic. [1,1,1] -> full 3d periodicity, [0,0,1] -> periodicity along the z axis
- tm – a Tol_matrix object, or a string representing Tol_matrix
- d_factor – the tolerance is multiplied by this amount. Larger values mean atoms must be farther apart
Returns: a bool for whether or not the atoms are sufficiently far enough apart
-
pyxtal.operations.
check_images
(coords, species, lattice, PBC=[1, 1, 1], tm=<pyxtal.tolerance.Tol_matrix object>, tol=None, d_factor=1.0)[source]¶ Given a set of (unfiltered) frac coordinates, checks if the periodic images are too close.
Parameters: - coords – a list of fractional coordinates
- species – the atomic species of each coordinate
- lattice – a 3x3 lattice matrix
- PBC – the periodic boundary conditions
- tm – a Tol_matrix object
- tol – a single override value for the distance tolerances
- d_factor – the tolerance is multiplied by this amount. Larger values mean atoms must be farther apart
Returns: False if distances are too close. True if distances are not too close
-
pyxtal.operations.
create_matrix
(PBC=[1, 1, 1], omit=False)[source]¶ Used for calculating distances in lattices with periodic boundary conditions. When multiplied with a set of points, generates additional points in cells adjacent to and diagonal to the original cell
Parameters: PBC – A periodic boundary condition list (1: periodic; 0: nonperiodic). Returns: A numpy array which can be multiplied by a set of coordinates
-
pyxtal.operations.
distance
(xyz, lattice, PBC=[1, 1, 1])[source]¶ Returns the Euclidean distance from the origin for a fractional displacement vector. Takes into account the lattice metric and periodic boundary conditions, including up to one non-periodic axis.
Parameters: - xyz – a fractional 3d displacement vector. Can be obtained by subtracting one fractional vector from another
- lattice – a 3x3 matrix describing a unit cell’s lattice vectors
- PBC – A periodic boundary condition list, where 1 means periodic, 0 means
- periodic. Ex (not) – [1,1,1] -> full 3d periodicity, [0,0,1] -> 1d
- along the z axis (periodicity) –
Returns: a scalar for the distance of the point from the origin
-
pyxtal.operations.
distance_matrix
(pts1, pts2, lattice, PBC=[1, 1, 1], single=False, metric='euclidean')[source]¶ Returns the distances between two sets of fractional coordinates. Takes into account the lattice metric and periodic boundary conditions.
Parameters: - pts1 – a list of fractional coordinates
- pts2 – another list of fractional coordinates
- lattice – a 3x3 matrix describing a unit cell’s lattice vectors
- PBC – A periodic boundary condition list, where 1 means periodic, 0 means not periodic. Ex: [1,1,1] -> full 3d periodicity, [0,0,1] -> 1d periodicity along the z axis
- single – return a scalor and matrix?
- metric – the metric to use with cdist, such as ‘euclidean’, ‘sqeuclidean’, ‘minkowski’, and others
Returns: a scalor or distance matrix
-
pyxtal.operations.
distance_matrix_no_PBC
(pts1, pts2, lattice, single=False, metric='euclidean')[source]¶ Returns the distances between two sets of fractional coordinates. Without periodic boundary conditions.
Parameters: - pts1 – a list of fractional coordinates (N1*3)
- pts2 – another list of fractional coordinates (N2*3)
- lattice – a 3x3 matrix describing a unit cell’s lattice vectors
- single – return the minimum distance or the matrix
- metric – the metric to use with cdist. e.g. euclidean, sqeuclidean, minkowski, and others
Returns: a scalor or distance matrix
-
pyxtal.operations.
filtered_coords
(coords, PBC=[1, 1, 1])[source]¶ Transform all coordinates to [0, 1] interval if PBC is allowed For example, [1.2, 1.6, -.4] becomes [0.2, 0.6, 0.6] when PBC=[1,1,1] [0.2, 1.6, 0.6] when PBC=[1,0,1]
Parameters: - coords – an array of real 3d vectors.
- PBC – A periodic boundary condition list (1: periodic; 0: nonperiodic).
Returns: an array of filtered coords with the same shape as coords
-
pyxtal.operations.
filtered_coords_euclidean
(coords, PBC=[1, 1, 1])[source]¶ Given an array of fractional 3-vectors, filters coordinates to between 0 and 1. Then, values which are greater than 0.5 are converted to 1 minus their value. This is used for converting displacement vectors with a Euclidean lattice.
Parameters: - coords – an array of real 3d vectors. The shape does not matter
- PBC – A periodic boundary condition list (1: periodic; 0: nonperiodic).
Returns: an array of filtered coords with the same shape as coords
-
pyxtal.operations.
get_best_match
(positions, ref, cell)[source]¶ find the best match with the reference from a set of positions
Parameters: - positions – N*3 array
- ref – 1*3 array
- cell – cell matrix 3*3 array
Returns: matched position id: matched id
Return type: position
-
pyxtal.operations.
get_inverse
(op)[source]¶ Given a SymmOp object, returns its inverse.
Parameters: op – a Symmop object Returns: the inverse
-
pyxtal.operations.
get_inverse_ops
(ops)[source]¶ Given a inverse list of SymmOp objects
Parameters: ops – a list of Symmop’s Returns: a list of equal shape to ops, with the inverse operations
-
pyxtal.operations.
is_orthogonal
(m, tol=0.001)[source]¶ Check whether or not a 3x3 matrix is orthogonal. An orthogonal matrix has the property that it is equal to its transpose.
Parameters: - m – a 3x3 matrix (list or numpy array)
- tol – the numerical tolerance for checking if two matrices are equal
Returns: True if the matrix is orthogonal, False if it is not
-
pyxtal.operations.
rotate_vector
(v1, v2, rtol=0.0001)[source]¶ Rotates a vector v1 to v2 about an axis perpendicular to both. Returns the 3x3 rotation matrix used to do so.
Parameters: - v1 – a 1x3 vector (list or array) of floats
- v2 – a 1x3 vector (list or array) of floats
Returns: a 3x3 matrix corresponding to a rotation which can be applied to v1 to obtain v2
-
pyxtal.operations.
verify_distances
(coordinates, species, lattice, factor=1.0, PBC=[1, 1, 1])[source]¶ Checks the inter-atomic distance between all pairs of atoms in a crystal.
Parameters: - coordinates – a 1x3 list of fractional coordinates
- species – a list of atomic symbols for each coordinate
- lattice – a 3x3 matrix representing the lattice vectors of the unit cell
- factor – a tolerance factor for checking distances. A larger value means atoms must be farther apart
- PBC – A periodic boundary condition list, where 1 means periodic, 0 means not periodic. Ex: [1,1,1] -> full 3d periodicity, [0,0,1] -> 1d periodicity along the z axis
Returns: True if no atoms are too close together, False if any pair is too close