pyxtal package

main pyxtal module to create the pyxtal class

print out the logo and version information

class pyxtal.pyxtal(molecular=False)[source]

Bases: object

Class for handling atomic crystals based on symmetry constraints

Examples

To create a new structure instance

>>> from pyxtal import pyxtal
>>> struc = pyxtal()

one can either use pyxtal to generate a random symmetric structure

>>> struc.from_random(3, 227, ['C'], [8])

or load a structure from a file or ASE atoms or Pymatgen structure object

>>> struc.from_seed('diamond.cif')     # as a string
>>> struc.from_seed(diamond_ase)       # as a ase atoms object
>>> struc.from_seed(diamond_pymatgen)  # as a pymatgen structure object

as long as the struc is created, you can check their symmetry as follows

>>> struc.get_site_labels()
{'C': ['8a']}
>>> struc
------Crystal from random------
Dimension: 3
Composition: C8
Group: Fd-3m (227)
cubic lattice:   4.3529   4.3529   4.3529  90.0000  90.0000  90.0000
Wyckoff sites:
     C @ [0.1250 0.1250 0.1250], WP:  8a, Site symmetry: -4 3 m

The structure object can be easily manipulated via apply_perturbtation or subgroup function

>>> struc2 = struc.subgroup_once(H=141)
>>> struc2
------Crystal from Wyckoff Split------
Dimension: 3
Composition: C8
Group: I41/amd (141)
tetragonal lattice:   3.3535   3.3535   4.6461  90.0000  90.0000  90.0000
Wyckoff sites:
     C @ [0.0000 0.2500 0.3750], WP:  4b, Site symmetry: -4 m 2

Alternatively, one can easily compute its XRD via the pyxtal.XRD class

>>> xrd = struc.get_XRD()
>>> xrd
  2theta     d_hkl     hkl       Intensity  Multi
  32.706     2.738   [ 1  1  1]   100.00        8
  54.745     1.677   [ 2  2  0]    40.95       12
  65.249     1.430   [ 3  1  1]    20.65       24
  81.116     1.186   [ 4  0  0]     5.15        6
  90.236     1.088   [ 3  3  1]     8.24       24
 105.566     0.968   [ 4  2  2]    14.44       24
 115.271     0.913   [ 5  1  1]    10.03       24
 133.720     0.838   [ 4  4  0]     9.80       12
 148.177     0.802   [ 5  3  1]    28.27       48

One can also try to get the transition path between two pyxtals that are symmetry related via the get_transition function

>>> s1 = pyxtal()
>>> s2 = pyxtal()
>>> s1.from_seed("pyxtal/database/cifs/0-G62.cif") #structure with low symmetry
>>> s2.from_seed("pyxtal/database/cifs/2-G71.cif") #structure with high symmetry
>>> strucs, _, _, _, _ = s2.get_transition(s1) # get the transition from high to low
>>> strucs
[
------Crystal from Transition 0  0.000------
Dimension: 3
Composition: O24Mg4W4Pb8
Group: Pnma (62)
orthorhombic lattice:  11.6075   8.0526   5.8010  90.0000  90.0000  90.0000
Wyckoff sites:
    Mg @ [ 0.8750  0.2500  0.7500], WP [4c] Site [.m.]
    Pb @ [ 0.6406  0.0053  0.7856], WP [8d] Site [1]
     W @ [ 0.6119  0.2500  0.2483], WP [4c] Site [.m.]
     O @ [ 0.6292  0.0083  0.2235], WP [8d] Site [1]
     O @ [ 0.4966  0.2500  0.0093], WP [4c] Site [.m.]
     O @ [ 0.5055  0.2500  0.4897], WP [4c] Site [.m.]
     O @ [ 0.7308  0.2500  0.9717], WP [4c] Site [.m.]
     O @ [ 0.7467  0.2500  0.4570], WP [4c] Site [.m.],
------Crystal from Transition 1  0.323------
Dimension: 3
Composition: O24Mg4W4Pb8
Group: Pnma (62)
orthorhombic lattice:  11.6020   8.0526   5.8038  90.0000  90.0000  90.0000
Wyckoff sites:
    Mg @ [ 0.8750  0.2500  0.7500], WP [4c] Site [.m.]
    Pb @ [ 0.6250 -0.0053  0.7500], WP [8d] Site [1]
     W @ [ 0.6250  0.2500  0.2500], WP [4c] Site [.m.]
     O @ [ 0.6250  0.0083  0.2500], WP [8d] Site [1]
     O @ [ 0.5158  0.2500 -0.0068], WP [4c] Site [.m.]
     O @ [ 0.5158  0.2500  0.5068], WP [4c] Site [.m.]
     O @ [ 0.7342  0.2500  0.9932], WP [4c] Site [.m.]
     O @ [ 0.7342  0.2500  0.5068], WP [4c] Site [.m.]]

Finally, the structure can be saved to different formats

>>> struc.to_file('my.cif')
>>> struc.to_file('my_poscar', fmt='poscar')

or to Pymatgen/ASE structure object

>>> pmg_struc = struc.to_pymatgen()
>>> ase_struc = struc.to_ase()

or to json file

>>> struc.to_json('1.json')
apply_perturbation(d_lat=0.05, d_coor=0.05, d_rot=1)[source]

perturb the structure without breaking the symmetry

Parameters:
  • d_coor – magnitude of perturbation on atomic coordinates (in A)

  • d_lat – magnitude of perturbation on lattice (in percentage)

build(group, species, numIons, lattice, sites, tol=0.01, use_hall=False)[source]

Build a atomic crystal based on the necessary input

Parameters:
  • group – 225

  • species – [‘Na’, ‘Cl’]

  • numIons – [4, 4]

  • lattice – lattice object

  • sites – [[{“4a”: [0.0, 0.0, 0.0]}], [{“4b”: [0.5, 0.5, 0.5]}]]

check_H_coordination(r=1.12)[source]

A function to check short if H is connected to more than one atom Mainly used for debug, powered by pymatgen

Parameters:

r – the given cutoff distances

Returns:

True or False

check_distance()[source]

Check intermolecular distance for molecular crystal

check_mapping(ref_struc)[source]

Compute the displacement w.r.t. the reference structure

Parameters:

ref_struc – reference pyxtal structure (assuming the same atom order)

Returns:

True or False

check_short_distances(r=0.7, exclude_H=True)[source]

A function to check short distance pairs Mainly used for debug, powered by pymatgen

Parameters:
  • r – the given cutoff distances

  • exclude_H – whether or not exclude the H atoms

Returns:

list of pairs within the cutoff

check_short_distances_by_dict(dicts)[source]

A function to check short distance pairs Mainly used for debug, powered by pymatgen

Parameters:

dicts – e.g., {“H-H”: 1.0, “O-O”: 2.0}

Returns:

number of atomic pairs within the cutoff

Return type:

N_pairs

copy()[source]

Simply copy the structure

find_matched_lattice(ref_struc, d_tol=2.0, f_tol=0.15)[source]

Compute the displacement w.r.t. the reference structure

Parameters:
  • ref_struc – reference pyxtal structure (assuming the same atomic ordering)

  • d_tol – tolerence of mismatch in the absolute scale

  • f_tol – tolerence of mismatch in the fractional scale

Returns:

ref_struc with matched lattice

from_CSD(csd_code)[source]

Download the crystal from CCDC if csd_code is given, return the single pyxtal object if csd_family is given, do group analysis and ignore high pressure form

Parameters:

csd_code – e.g., ACSALA01

from_json(filename)[source]

Load the model from a json file

from_random(dim=3, group=None, species=None, numIons=None, factor=1.1, thickness=None, area=None, lattice=None, sites=None, conventional=True, t_factor=1.0, max_count=10, torsions=None, force_pass=False, block=None, num_block=None, seed=None, tm=None, use_hall=False)[source]
from_seed(seed, molecules=None, tol=0.0001, a_tol=5.0, ignore_HH=True, add_H=False, backend='pymatgen', style='pyxtal', hn=None, standard=False)[source]

Load the seed structure from Pymatgen/ASE/POSCAR/CIFs Internally they will be handled by Pymatgen

Parameters:
  • seed – cif/poscar file or a Pymatgen Structure object

  • molecules – a list of reference molecule (xyz file or Pyxtal molecule)

  • tol – scale factor for covalent bond distance

  • ignore_HH – whether or not ignore short H-H distance in molecules

  • add_H – whether or not add H atoms

  • backend – structure parser, default is pymatgen

  • style – pyxtal for spglib

  • standard – whether or not optimize lattice

get_1D_comp()[source]

Get composition for 1d rep of molecular xtal

get_1D_representation()[source]

Get the 1D representation class for molecular crystals

get_XRD(**kwargs)[source]

Compute the PXRD object.

** kwargs include
  • wavelength (1.54184)

  • thetas [0, 180]

  • preferred_orientation: False

  • march_parameter: None

get_alternatives(include_self=True, same_letters=False, ref_lat=None, d_tol=2.0, f_tol=0.15)[source]

Get alternative structure representations

Parameters:

include_self (bool) – return the original structure

Returns:

list of structures

get_density()[source]

Compute density

get_dimensionality(cutoff=None)[source]

A quick wrapper to compute dimensionality from pymatgen https://pymatgen.org/pymatgen.analysis.dimensionality.html The dimensionality of the structure can be 1/2/3

get_disps_optim(ref_struc, trans, d_tol)[source]
Parameters:
  • ref_struc – reference pyxtal structure (assuming the same atom order)

  • trans – translation vector

  • d_tol – tolerence of mismatch

Returns:

Atomic displacements in np.array translation:

get_disps_sets(ref_struc, d_tol, d_tol2=0.3, ld_tol=2.0, fd_tol=0.15, keep_lattice=False)[source]

Compute the displacement w.r.t. a reference structure (considering all wycsets)

Parameters:
  • ref_struc – reference pyxtal structure (assuming the same atomic ordering)

  • d_tol – maximally allowed atomic displacement

  • d_tol2 – displacement that allows early termination

  • kepp_lattice – whether or not change the WP sets

Returns:

Atomic displacements in np.array

get_disps_single(ref_struc, trans, d_tol=1.2)[source]

Compute the displacement w.r.t. the reference structure

Parameters:
  • ref_struc – reference pyxtal structure (assuming the same atom order)

  • trans – translation vector

  • d_tol – tolerence of mismatch

Returns:

Atomic displacements in np.array translation:

get_dof()[source]

Get the number of dof for the given structures:

get_forcefield(ff_style='openff', code='lammps', chargemethod='am1bcc', parameters=None)[source]

An interface to create forcefield for molecular simulation with - Charmm - LAMMPS Note that the current approach generates charge with its own conformer, need to provide the option to use the conformer from the given molecule

Parameters:
  • ff_style (-) – ‘gaff’ or ‘openff’

  • code (-) – ‘lammps’ or ‘charmm’

  • charge_method (-) – ‘am1bcc’, ‘am1-mulliken’, ‘mmff94’, ‘gasteiger’

  • parameters (-) – 1D-array of user-specified FF parameters

Returns:

An ase_atoms_objects with force field information

get_free_axis()[source]

Check if the a, b, c axis have free parameters

get_init_translations(ref_struc, tol=0.75)[source]

Compute the displacement w.r.t. the reference structure

Parameters:

ref_struc – reference pyxtal structure (assuming the same atom order)

Returns:

list of possible translations

get_intermolecular_energy(factor=2.0, max_d=10.0)[source]

For molecular crystals, get the intermolecular interactions from Gavezzotti, A., Filippini, G., J. Phys. Chem., 1994, 98 (18), 4831-4837

Returns:

Total energy

get_neighboring_dists(site_id=0, factor=1.5, max_d=5.0)[source]

For molecular crystals, get the neighboring molecules for a given WP

Parameters:
  • site_id – the index of reference site

  • factor – factor of vdw tolerance

  • max_d – maximum distances

Returns:

list of short contact pairs engs: list of energies from atom-atom potential

Return type:

pairs

get_neighboring_molecules(site_id=0, factor=1.5, max_d=5.0, ignore_E=True)[source]

For molecular crystals, get the neighboring molecules for a given WP

Parameters:
  • site_id – the index of reference site

  • factor – factor of vdw tolerance

  • max_d

  • ignore_E

Returns:

list of shortest distances neighs: list of neighboring molecular xyzs comps: list of molecular types Ps: list of [0, 1] to distinguish self and other molecules engs: list of energies from atom-atom potential

Return type:

min_ds

get_num_torsions()[source]

Get number of torsions for molecular xtal

get_site_labels()[source]

Get the site_labels as a dictionary

get_spherical_images(**kwargs)[source]

get the spherical image representation

Parameters:

model – either ‘molecule’ or ‘contacts’

Returns:

the sph class

get_std_representation(trans)[source]

Perform cell transformation so that the symmetry operations follow standard space group notation

get_structure_factor(hkl, coeffs=None)[source]

Compute the structure factor for a given hkl plane

Parameters:
  • hkl (-) – three indices hkl plane

  • coeffs (-) – Atomic scatering factors

get_transition(ref_struc, d_tol=1.0, d_tol2=0.3, N_images=2, max_path=30, both=False)[source]

Get the splitted wyckoff information along a given path:

Parameters:
  • ref_struc – structure with subgroup symmetry

  • d_tol – maximally allowed atomic displacement

  • d_tol2 – displacement that allows early termination

  • N_images – number of intermediate images

  • max_path – maximum number of paths

  • both – whether or not do interpolation along both sides

Returns:

  • displacements:

  • cell translation:

  • the list of space groups along the path

  • the list of splitters along the path

Return type:

  • strucs

get_transition_by_path(ref_struc, path, d_tol, d_tol2=0.5, N_images=2, both=False)[source]

Get the splitted wyckoff information along a given path:

Parameters:
  • ref_struc – structure with subgroup symmetry

  • path – a list of transition path

  • d_tol – maximally allowed atomic displacement

  • d_tol2 – displacement that allows early termination

  • N_images – number of intermediate images

  • both – interpolation on both sides

Returns:

  • displacements:

  • cell translation:

Return type:

  • strucs

get_zprime(integer=False)[source]

Get zprime for molecular xtal

has_special_site(species=None)[source]

Check if the crystal has a special site

is_duplicate(ref_strucs)[source]

check if the structure is exactly the same

load_dict(dict0)[source]

Load the structure from a dictionary

make_transitions(disps, lattice=None, translation=None, N_images=3, both=False)[source]

make the pyxtals by following the atomic displacements

Parameters:
  • disps – N*3 atomic displacements in fractional coordinates

  • lattice – 3*3 cell matrix (self.lattice.matrix if None)

  • translation – overall translation

  • N_images – number of images

  • both – whether or not interpolar on both sides

optimize_lattice(iterations=5, force=False, standard=False)[source]

Optimize the lattice if the cell has a bad inclination angles We first optimize the angle to some good-looking range. If standard, force the structure to have the standard setting This only applies to monoclinic systems

Parameters:
  • iterations – maximum number of iterations

  • force – whether or not do the early termination

  • standard – True or False

remove_species(species)[source]

To remove the atom_sites by specie name :param dicts: e.g., [“O”]

remove_water()[source]

Remove water from hydrates

resort()[source]

A short cut to resort the sites by self.molecules or self.species

resort_species(species)[source]

resort the atomic species

Parameters:

species – list of elements, e.g. [‘Si’, ‘O’]

save_dict()[source]

Save the model as a dictionary

set_cutoff(exclude_ii=False)[source]

get the cutoff dictionary

set_site_coordination(cutoff=None, verbose=False, exclude_ii=False)[source]

Compute the coordination number from each atomic site

show(**kwargs)[source]

display the crystal structure

show_mol_cluster(id, factor=1.5, max_d=4.0, plot=True, ignore_E=False, cmap='YlGn', **kwargs)[source]

display the local packing environment for a selected molecule

sort_sites_by_mult()[source]
sort_sites_by_numIons(seq=None)[source]
subgroup(perms=None, H=None, eps=0.05, idx=None, group_type='t', max_cell=4, min_cell=0)[source]

Generate a structure with lower symmetry

Parameters:
  • perms – e.g., {“Si”: “C”}

  • H – space group number (int)

  • eps – pertubation term (float)

  • idx – list

  • group_typet, k or t+k

  • max_cell – maximum cell reconstruction (float)

Returns:

a list of pyxtal structures with lower symmetries

subgroup_by_path(gtypes, ids, eps=0, mut_lat=False)[source]

Generate a structure with lower symmetry (for atomic crystals only)

Parameters:
  • g_types – [‘t’, ‘k’]

  • idx – list of ids for the splitter

  • eps – degree of displacement

  • mut_lat – mutate the lattice of not

Returns:

a pyxtal structure with lower symmetries list of splitters

subgroup_once(eps=0.1, H=None, perms=None, group_type='t', max_cell=4, min_cell=0, mut_lat=True, ignore_special=False)[source]

Generate a structure with lower symmetry (for atomic crystals only)

Parameters:
  • perms – e.g., {“Si”: “C”}

  • H – space group number (int)

  • idx – list

  • group_typet or k

  • max_cell – maximum cell reconstruction (float)

Returns:

a pyxtal structure with lower symmetries

substitute(dicts)[source]

A quick function to apply substitution

Parameters:

dicts – e.g., {“F”: “Cl”}

substitute_linear(dicts)[source]

This is mainly designed for substitution between single atom and the linear block e.g., we want to make Zn(CN)2 from SiO2 :param dicts: e.g., {“Si”: [“Zn”], “O”: [“C”,”N”]}

supergroup(G=None, d_tol=1.0)[source]

Generate a structure with higher symmetry

Parameters:
  • G – super space group number (list of integers)

  • d_tol – maximum tolerance

Returns:

a list of pyxtal structures with minimum super group symmetries

supergroups(G=None, d_tol=1.0)[source]

Generate the structures with higher symmetry

Parameters:
  • G – super space group number (list of integers)

  • d_tol – maximum tolerance

Returns:

a list of pyxtal structures with minimum super group symmetries

to_ase(resort=True, center_only=False)[source]

Export to ase Atoms object.

to_atomic_xtal()[source]

Convert the molecular

to_file(filename=None, fmt=None, permission='w', sym_num=None, header='from_pyxtal')[source]

Creates a file with the given ame and type to store the structure. By default, creates cif files for crystals and xyz files for clusters. For other formats, Pymatgen is used

Parameters:
  • filename (string) – the file path

  • fmt (string) – the file type (cif, xyz, etc.)

  • permission (string) – w or a+

  • sym_num (int) – number of sym_ops, None means writing all symops

  • header (string) – header

Returns:

Nothing. Creates a file at the specified path

to_json(filename='pyxtal.json')[source]

Save the model as a dictionary

to_molecular_xtal(molecules, oris=None, reflects=None)[source]

Convert the atomic xtal to molecular xtal the input molecules must have the same length of the self.atom_sites

to_pymatgen(resort=True, shape='upper')[source]

Export to Pymatgen structure object.

to_pyxtal_center()[source]

Export to PyXtal object for molecular centers only.

to_standard_setting()[source]

A short cut to symmetrize the structure in the stardard setting

to_subgroup(path=None, t_only=True, iterate=False, species=None)[source]

Transform a crystal with speical sites to subgroup represenatation with general sites

Parameters:
  • Path – list of path to get the general sites

  • iterate (bool) – whether or not do it iteratively

transform(trans, lattice=None)[source]

Perform cell transformation and symmetry operation

Parameters:
  • trans – 3*3 matrix

  • lattice – pyxtal lattice object

  • update – whether or not update each wp

translate(trans, reset_wp=False)[source]

move the atomic sites along a translation Note that this may change the structure

Parameters:

trans – 1*3 vector

Subpackages

Submodules