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: object

Class for storing atomic Wyckoff positions with a single coordinate.

Parameters:
  • wp – a `Wyckoff_position <pyxtal.symmetry.Wyckoff_position.html> object

  • coordinate – a fractional 3-vector for the generating atom’s coordinate

  • specie – an Element, element name or symbol, or atomic number of the atom

  • search – to search for the optimum position for special wyckoff site

check_with_ws2(ws2, lattice, tm, same_group=True)[source]

Given two Wyckoff sites, checks the inter-atomic distances between them.

Parameters:
  • ws2 – a different Wyckoff_site object (will always return False if

  • provided) (two identical WS's are)

  • lattice – a 3x3 cell matrix

  • same_group – whether or not the two WS’s are in the same structure.

  • cost (Default value True reduces the calculation)

Returns:

True if all distances are greater than the allowed tolerances. False if any distance is smaller than the allowed tolerance

copy()[source]

Simply copy the structure

encode()[source]

transform dict to 1D vector [specie, wp.index, free x, y, z]

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]

return the displacement towards the reference positions

Parameters:
  • pos – reference position (1*3 vector)

  • lattice – 3*3 matrix

  • translation

  • axis

classmethod load_dict(dicts)[source]

load the sites from a dictionary

perturbate(lattice, magnitude=0.1)[source]

Random perturbation of the site

Parameters:
  • lattice – lattice vectors

  • magnitude – the magnitude of displacement (default: 0.1 A)

save_dict()[source]
search_position()[source]

Sometimes, the initial posiition 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 May occur for special WP in the I/A/B/C/F cases e.g., in space group 71 (Immm), the permutation 4j(1/2,0,0.2) -> (0.2,0,1/2) -> 4f(0.7,1/2,0) it requires a shift of (0.5,0.5,0.5)

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’

swap_axis(swap_id, shift=array([0., 0., 0.]))[source]

sometimes space groups like Pmm2 allows one to swap the a,b axes to get an alternative representation

to_mol_site(lattice, molecule, ori=[0, 0, 0], reflect=False, type_id=0)[source]

transform it to the mol_sites, i.e., to build a molecule on the current WP

update(pos=None, reset_wp=False)[source]

Used to generate coords from self.position

class pyxtal.wyckoff_site.mol_site(mol, position, orientation, wp, lattice=None, stype=0)[source]

Bases: object

Class for storing molecular Wyckoff positions and orientations within the molecular_crystal class. Each mol_site object represenents an entire Wyckoff position, not necessarily a single molecule. This is the molecular version of Wyckoff_site

Parameters:
  • mol – a pyxtal_molecule object

  • position – the 3-vector representing the generating molecule’s position

  • orientation – an Orientation object

  • wp – a Wyckoff_position object

  • lattice – a Lattice object

  • stype – integer number to specify the type of molecule

encode()[source]

transform dict to 1D vector [x, y, z, or1, or2, or3, rotor1, rotor2, .etc]

classmethod from_1D_dicts(dicts)[source]
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_distances(coord1, coord2, m2=None, center=True, ignore=False)[source]

Compute the distance matrix between the center molecule (m1 length) and neighbors (m2 length) within the PBC consideration (pbc)

Parameters:
  • coord1 – fractional coordinates of the center molecule

  • coord2 – fractional coordinates of the reference neighbors

  • m2 – the length of reference molecule

  • center – whether or not consider the self image for coord2

  • ignore

Returns:

[m1*m2*pbc, m1, m2] coord2 under PBC: [pbc, m2, 3]

Return type:

distance matrix

get_dists_WP(ignore=False, id=None)[source]

Compute the distances within the WP sites

Returns:

a distance matrix (M, N, N) list of molecular xyz (M, N, 3)

get_dists_auto(ignore=False)[source]

Compute the distances between the periodic images

Returns:

a distance matrix (M, N, N) list of molecular xyz (M, N, 3)

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()[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 the neigboring molecules

Parameters:
  • factor – volume factor

  • max_d – maximum intermolecular distance

  • ignore_E

  • detail – show detailed energies

Returns

min_ds: list of shortest distances neighs: list of neighboring molecular xyzs

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

is_valid_orientation()[source]
classmethod load_dict(dicts)[source]

load the sites from a dictionary

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

save_dict()[source]
short_dist()[source]

Check if the atoms are too close within the WP.

Returns:

True or False

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

show(id=None, **kwargs)[source]

display WP on the notebook

show_molecule_in_box(id=0)[source]

display molecule with box

Parameters:

id (int) – molecular id

to_1D_dicts()[source]

save the wp in 1D representation

to_atom_site(specie=1)[source]

transform it to the mol_sites, i.e., to build a molecule on the current WP

translate(disp=array([0., 0., 0.]), absolute=False)[source]

To translate the molecule

update(coords, 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.

update_lattice(lattice)[source]
update_molecule(mol)[source]
update_orientation(angles)[source]