pyxtal.wyckoff_site module
Module for handling Wyckoff sites for both atom and molecule
- class pyxtal.wyckoff_site.atom_site(wp=None, coordinate=None, specie=1, search=False)[source]
Bases:
objectClass for storing atomic Wyckoff positions with a single coordinate.
- Parameters:
wp – a WP <pyxtal.symmetry.Wyckoff_position.html> object
coordinate (float) – a fractional (x, y, z) coordinate
specie (str) – element name or symbol, or atomic number
search (bool) – whether or not search generator for special WP
- check_with_ws2(ws2, lattice, tol, same_group=True)[source]
Given two Wyckoff sites, checks the inter-atomic distances between them.
- Parameters:
ws2 (-) – a different WP object (will always return False if
provided) (two identical WS's are)
lattice (-) – a 3x3 cell matrix
tol (-) – tolerance value for distance checking
same_group (-) – whether or not two WS’s are in the same structure.
cost (Default value True reduces the calculation)
- Returns:
True or False
- equivalent_set(tran, indices)[source]
Transform the wp to another equivalent set. Needs to update both wp and positions
- Parameters:
tran – affine matrix
indices – the list of transformed wps
- get_disp(pos, lattice, translation)[source]
return the displacement towards the reference positions
- Parameters:
pos – reference position (1*3 vector)
lattice – 3*3 matrix
translation
- get_translations(pos, axis)[source]
Get the displacement towards the reference positions
- Parameters:
pos – reference position (1*3 vector)
lattice – 3*3 matrix
translation
axis
- perturbate(lattice, magnitude=0.1)[source]
Random perturbation of the site
- Parameters:
lattice – lattice vectors
magnitude – the magnitude of displacement (default: 0.1 A)
- search_position(tol=0.001)[source]
Sometimes, the initial position is not the proper generator Needs to find the proper generator
- shift_by_swap(swap_id)[source]
Check if a shift is needed during swap.
This may occur for special WP in the I/A/B/C/F cases. For example, in space group 71 (Immm), the permutation: 4j(1/2,0,0.2) -> (0.2,0,1/2) -> 4f(0.7,1/2,0) requires a shift of (0.5,0.5,0.5)
- Parameters:
swap_id – The index of the swap operation
- Returns:
The required shift vector
- Return type:
numpy array
- substitute_with_linear(eles, direction, lattice)[source]
chemical substitution with another linear building block, e.g. CN
- Parameters:
eles (str) – e.g., [‘C’, ‘N’]
neighbors – two nearest neiboring atom xyz
- substitute_with_single(ele)[source]
chemical substitution with another element
- Parameters:
ele (str) – e.g. ‘Zn’
- class pyxtal.wyckoff_site.mol_site(mol, position, orientation, wp, lattice=None, stype=0)[source]
Bases:
objectClass for storing molecular Wyckoff positions and orientations within the molecular_crystal class. Each mol_site object represents an entire Wyckoff position, not necessarily a single molecule. This is the molecular version of the Wyckoff_site class.
- mol
A pyxtal_molecule object representing the molecule at the site.
- Type:
- position
A 3-vector representing the generating molecule’s position in fractional coordinates.
- Type:
list or array
- orientation
An orientation object that describes the molecule’s orientation.
- Type:
- wp
A Wyckoff position object that holds symmetry information.
- Type:
- stype
An integer specifying the type of molecule. Default is 0.
- Type:
int
- symbols
List of atomic symbols for the atoms in the molecule.
- Type:
list
- numbers
List of atomic numbers for the atoms in the molecule.
- Type:
list
- PBC
Periodic boundary conditions inherited from the Wyckoff position.
- Type:
list
- radius
Radius of the molecule, typically used for collision detection.
- Type:
float
- tols_matrix
Tolerance matrix for the molecular structure.
- Type:
numpy array
- Parameters:
mol (pyxtal_molecule) – A pyxtal_molecule object that describes the molecule.
position (list or array) – The 3D fractional coordinates of mol_center in the unit cell.
orientation (Orientation) – The orientation object describing the molecule’s rotation.
wp (Wyckoff_position) – A Wyckoff_position object defining the symmetry of the site.
lattice (Lattice, optional) – The lattice of the crystal. Can be either a Lattice object or a matrix that will be converted into a Lattice object.
stype (int, optional) – Integer specifying the type of molecule. Default is 0.
- cut_lattice(ax, cut)[source]
Cut lattice length on the given direction
- Parameters:
ax (int) – 0, 1, 2
cut (float) – the cut
- get_coords_and_species(absolute=False, PBC=False, unitcell=False)[source]
Lazily generates and returns the atomic coordinate and species for the Wyckoff position. Plugs the molecule into the provided orientation (with angle=0), and calculates the new positions.
- Parameters:
absolute – return absolute or relative coordinates
PBC – whether or not to add coordinates in neighboring unit cells
unitcell – whether or not to move the molecule center to the unit cell
- Returns:
a np array of 3-vectors. species: a list of atomic symbols, e.g. [‘H’, ‘H’, ‘O’, ‘H’, ‘H’, ‘O’]
- Return type:
coords
- get_dist_angle_AD(d_total, ids_total, A_id, D_ids, H_ids, rA)[source]
Calculates distances/angles for H bonds between acceptor (A), donor (D), and hydrogen (H) atoms.
- Parameters:
d_total (numpy.ndarray) – Array containing atomic coordinates and distances
ids_total (numpy.ndarray) – Array containing atom indices and molecule IDs
A_id (int) – Index of the acceptor atom
D_ids (dict) – Dictionary mapping molecule IDs to donor atom indices
H_ids (dict) – Dictionary mapping molecule IDs to hydrogen atom indices
rA (numpy.ndarray) – Coordinates of the acceptor atom
- Returns:
A tuple containing three float values: - dAD (float): Distance between acceptor and donor atoms - angle (float): Angle between H-D and H-A vectors in radians - dHA (float): Distance between hydrogen and acceptor atoms If no donor atoms are found, returns (10.0, 0, 10.0)
- Return type:
tuple
- get_distances(coord1, coord2, m2=None, center=True, ignore=False)[source]
Compute the distance matrix between the central molecule (coord1) and neighboring molecules (coord2) under the periodic boundary condition.
- Parameters:
coord1 (numpy array) – Fractional coordinates of the central molecule. Shape: (m1, 3), where m1 is the number of atoms
coord2 (numpy array) – Fractional coordinates of the neighboring molecules. Shape: (N2*m2, 3), where N2 is the number of atoms and m2 is the number of atoms in each neighboring molecule.
m2 (int, optional) – N_atoms in each neighboring molecule. If not provided, it’s assumed to be equal m1.
center (bool, optional) – If True, count self-image of the reference molecule
ignore (bool, optional) – If True, ignores some periodic boundary conditions.
- Returns:
[m1*m2*pbc, m1, m2] coord2 under PBC: [pbc, m2, 3]
- Return type:
distance matrix
- get_dists_WP(matrix=None, ignore=False, idx=None, cutoff=None)[source]
Compute the distances within the WP site
- Returns:
a distance matrix (M, N, N) list of molecular xyz (M, N, 3)
- get_dists_WP2(wp2, matrix=None, ignore=False, cutoff=None)[source]
Compute the distances w.r.t to other WP site
- Returns:
a distance matrix (M, N, N) list of molecular xyz (M, N, 3)
- get_dists_auto(matrix=None, ignore=False, cutoff=None)[source]
Compute the distances between the periodic images N: number of atoms in the molecule M: number of periodic images
- Parameters:
ignore (bool, optional) – If True, ignores some periodic boundary conditions.
() (cutoff) – if not None, reduce the distance matrix
- Returns:
a distance matrix (M, N, N) list of molecular xyz (M, N, 3)
- get_energy(angle=None, wps=[], k1=1.0, r1=1.0, k2=2.0, k3=1.0, r2=6.0, req=2.8, r_cut=2.0, verbose=False)[source]
Compute the virtual energy based on Hbonding and repulsion
- Parameters:
angle (float) – rotation along a rotation axis
wps (list) – list of wps for consideration
- get_ijk_lists(value=None)[source]
Get the occupatation in the unit cell for the generating molecule This can be used to estimate the supercell size for finding neighbors
- Returns
PBC
- get_min_dist(angle=None)[source]
Compute the minimum interatomic distance within the WP.
- Returns:
minimum distance
- get_mol_object(id=0)[source]
make the pymatgen molecule object
- Parameters:
id – the index of molecules in the given site
- Returns:
a molecule object
- get_neighbors_auto(factor=1.1, max_d=4.0, ignore_E=True, detail=False, etol=-0.05)[source]
Find neighboring molecules within a given distance threshold.
The function identifies neighboring molecules within PBC and computes the shortest distances and (optionally) interaction energies. it returns detailed information about neighboring molecule pairs, distances, and energies.
- Parameters:
factor (float, optional) – Scaling factor for distance (default is 1.1).
max_d (float, optional) – Maximum distance for neighbors (default is 4.0 Å).
ignore_E (bool, optional) – Skips energy calculations (default is True).
detail (bool, optional) – Returns detailed eng, molecular pairs and distances instead of only shortest distances.
etol (float, optional) – Energy tolerance for filtering pairs in detailed mode (default is -5e-2).
- Returns:
- engs (list): List of interaction energies for valid molecular pairs.
pairs (list): List of tuples with neighbor molecules and positions. dists (list): List of distances between neighboring molecular pairs.
- If detail is False:
min_ds (list): List of shortest distances between central and neighbors. neighs (list): List of neighboring molecular coordinates with PBC. Ps (list): List of Wyckoff position multiplicities or translations. engs (list): List of energies, if energy calculation is not skipped.
- Return type:
If detail is True
- get_neighbors_wp2(wp2, factor=1.1, max_d=4.0, ignore_E=True, detail=False, etol=-0.05)[source]
Find the neigboring molecules from a 2nd wp site
- Returns
min_ds: list of shortest distances neighs: list of neighboring molecular xyzs
- optimize_orientation_by_dist(ori_attempts=10, verbose=False)[source]
Optimize the orientation based on the shortest distance
- optimize_orientation_by_energy(wps=[], max_ax=20, max_ori=5, early_quit=3.0, verbose=False)[source]
Iteratively optimize the orientation with the bisection method
- perturbate(lattice, trans=0.1, rot=5)[source]
Random perturbation of the molecular site
- Parameters:
lattice – lattice vectors
trans – magnitude of tranlation vectors (default: 0.1 A)
rot – magnitude of rotation degree (default: 5.0)
- rotate(ax_id=0, ax_vector=None, angle=180)[source]
To rotate the molecule :param ax_id: the principle axis id :param ax_vector: 3-vector to define the axis :type ax_vector: float :param angle: angle to rotate :type angle: float
- short_dist_with_wp2(wp2, tm=<pyxtal.tolerance.Tol_matrix object>)[source]
Check whether or not the molecules of two wp sites overlap. Uses ellipsoid overlapping approximation to check.
- Parameters:
wp2 – the 2nd wp sites
tm – a Tol_matrix object (or prototype string) for distance checking
- Returns:
True or False
- to_atom_site(specie=1)[source]
transform it to the mol_sites, i.e., to build a molecule on the current WP
- update(coords=None, lattice=None, absolute=False, update_mol=True)[source]
After the geometry relaxation, the returned atomic coordinates maybe rescaled to [0, 1] bound. In this case, we need to refind the molecular coordinates according to the original neighbor list. If the list does not change, we return the new coordinates otherwise, terminate the calculation.
- Parameters:
coords (float) – input xyz coordinates for the given molecule
absolute (bool) – whether the input is absolute coordinates
update_mol (bool) – whether update the molecule (used for substitution)