Skip to content

Python API

Data type Structure Description
cryst (lattice, species, positions) A crystal structure described by a primitive cell lattice, atomic species species, and fractional coordinates positions. Cell vectors are rows of lattice.
slm (hA, hB, q) A sublattice match (SLM) described by an integer-matrix triplet.

crystmatch

Enumerating and analyzing crystal-structure matches.

load_poscar(filename, to_primitive=True, tol=0.001, verbose=True)

Load the crystal structure from a POSCAR file.

Parameters:

Name Type Description Default
filename str

The name of the POSCAR file to be read.

required
to_primitive bool

If True, reduced the crystal structure to its primitive structure.

True
tol float

The tolerance for spglib symmetry detection.

0.001
verbose bool

If True, print verbose output during loading.

True

Returns:

Name Type Description
cryst cryst

The loaded crystal structure, represented as a tuple of (lattice, species, positions).

load_csmcar(filename, verbose=True)

Load the CSMCAR file, which contains crystmatch parameters.

Parameters:

Name Type Description Default
filename str

The name of the POSCAR file to be read.

required
verbose bool

If True, print verbose output during loading.

True

Returns:

Name Type Description
voigtA, voigtB : (6, 6) array

The loaded elastic tensor for the initial and final structure, in Voigt notation (ordered as XX, YY, ZZ, YZ, ZX, XY).

weight_func dict

The loaded weight function for the shuffle distance.

ori_rel (2, 2, 3) array

The two loaded parallelisms, representing the orientation relationship between the initial and final structure.

unique_filename(message, filename, ext=True)

Generate a unique filename by appending a counter to the base name if the file already exists.

Parameters:

Name Type Description Default
message str

A message to print when generating the unique filename.

required
filename str

The original filename.

required
ext bool

If True, assume that filename has an extension. This should be set to False if filename does not have an extension.

True

Returns:

Name Type Description
new_filename str

The unique filename, which is the original filename with a counter appended if necessary.

cryst_to_spglib(cryst, return_dict=False)

Convert a crystal structure to a spglib-compatible format.

Parameters:

Name Type Description Default
cryst cryst

The crystal structure, represented as a tuple of (lattice, species, positions).

required
return_dict bool

If True, return a dictionary mapping species to unique numbers.

False

Returns:

Name Type Description
spglib_cell 3 - tuple

The spglib-compatible format.

species_dict (dict, optional)

A dictionary mapping species to unique numbers. Only returned if return_dict is True.

spglib_to_cryst(spglib_cell, species_dict)

Convert a spglib-compatible cell to a crystal structure.

Parameters:

Name Type Description Default
spglib_cell 3 - tuple

The spglib-compatible cell, represented as (lattice, positions, species_numbers).

required
species_dict dict

A dictionary mapping species to unique numbers, usually obtained from cryst_to_spglib().

required

Returns:

Name Type Description
cryst cryst

The crystal structure, represented as a tuple of (lattice, species, positions).

primitive_cryst(cryst_sup, tol=0.001)

Reduce a crystal structure to its primitive structure using spglib.

Parameters:

Name Type Description Default
cryst_sup cryst

The crystal structure to be reduced, represented as a tuple of (lattice, species, positions).

required
tol float

The tolerance for symmetry detection.

0.001

Returns:

Name Type Description
cryst cryst

The primitive crystal structure, represented as a tuple of (lattice, species, positions).

check_stoichiometry(speciesA, speciesB)

Check if the stoichiometry of two species lists is the same.

Parameters:

Name Type Description Default
speciesA array - like

The species lists to be checked, which should be 1D arrays of atomic species.

required
speciesB array - like

The species lists to be checked, which should be 1D arrays of atomic species.

required

create_common_supercell(crystA, crystB, slm)

Return the initial, final, and half-distorted supercell, as well as the transformation matrices.

Parameters:

Name Type Description Default
crystA cryst

The initial and final structures.

required
crystB cryst

The initial and final structures.

required
slm slm

The SLM of the CSM.

required

Returns:

Name Type Description
crystA_sup cryst

The supercell structure of crystA.

crystB_sup cryst

The supercell structure of crystB.

c_sup_half (3, 3) array of floats

The half-distorted supercell.

mA (3, 3) array of ints

The matrix that transforms crystA to crystA_sup.

mB (3, 3) array of ints

The matrix that transforms crystB to crystB_sup.

frac_cell(mA, mB)

The primitive cell of the lattice generated by mA^{-1} and mB^{-1}.

Parameters:

Name Type Description Default
mA (3, 3) array of ints

The transformation matrix from the initial crystal structure to the supercell.

required
mB (3, 3) array of ints

The transformation matrix from the final crystal structure to the supercell.

required

Returns:

Name Type Description
fcell (3, 3) array of floats

The primitive cell of the lattice generated by mA^{-1} and mB^{-1}.

imt_multiplicity(crystA, crystB, slmlist)

Return multiplicities of elements in slmlist.

Parameters:

Name Type Description Default
crystA cryst

The initial crystal structure, usually obtained by load_poscar.

required
crystB cryst

The final crystal structure, usually obtained by load_poscar.

required
slmlist list of slm

A list of SLMs, each represented by a triplet of integer matrices like (hA, hB, q).

required

Returns:

Name Type Description
mu int or (...,) array of ints

Multiplicities of each SLM in slmlist.

deformation_gradient(crystA, crystB, slmlist)

Compute the deformation gradient matrices of given IMTs.

Parameters:

Name Type Description Default
crystA cryst

The initial and final structures.

required
crystB cryst

The initial and final structures.

required
slmlist list of slm

A list of uuLMs, each represented by a triplet of integer matrices like (hA, hB, q).

required

Returns:

Name Type Description
slist (..., 3, 3) array

A list of deformation gradient matrices.

rmss(slist)

Root-mean-square strains of given deformation gradient matrices.

Parameters:

Name Type Description Default
slist (..., 3, 3) array

A list of deformation gradient matrices.

required

Returns:

Name Type Description
rmss (...) array

Root-mean-square strains.

zip_pct(p, ks)

Zip the permutation p and the class-wise translations ks into a single PCT array.

Parameters:

Name Type Description Default
p (N,) array of ints

The permutation part of the PCT.

required
ks (3, N) array of ints

The class-wise translations of the PCT.

required

Returns:

Name Type Description
pct (N, 4) array of ints

The zipped PCT array.

unzip_pct(pct)

Unzip the PCT array into the permutation and class-wise translations.

Parameters:

Name Type Description Default
pct (N, 4) array of ints

The PCT array, where the first column is the permutation and the rest are class-wise translations.

required

Returns:

Name Type Description
p (N,) array of ints

The permutation part of the PCT.

ks (3, N) array of ints

The class-wise translations of the PCT.

get_pure_rotation(cryst, tol=0.001)

Find all pure rotations appeared in the space group of cryst.

Parameters:

Name Type Description Default
cryst cryst

The crystal structure, usually obtained by load_poscar.

required
tol float

The tolerance for spglib symmetry detection.

0.001

Returns:

Name Type Description
g (..., 3, 3) array of ints

A point group of the first kind, containing all pure rotations appeared in the space group of cryst, elements of which are integer matrices (using the columns of cryst[1] as basis vectors).

setdiff2d(arr1, arr2)

Find the set difference of two 2D integer arrays, arr1 and arr2. This function is provided by @norok2 on StackOverflow: https://stackoverflow.com/a/66674679.

Parameters:

Name Type Description Default
arr1 (N, M) array of ints

The first array, from which elements will be removed.

required
arr2 (P, M) array of ints

The second array, containing elements to be removed from arr1.

required

Returns:

Name Type Description
result (Q, M) array of ints

The set difference of arr1 and arr2.

int_vec_inside(c)

Integer vectors inside the cell c @ [0, 1)^3 whose elements are integers.

Parameters:

Name Type Description Default
c (3, 3) array of ints

A matrix whose columns are integer cell vectors.

required

Returns:

Name Type Description
v_int (3, ...) array of ints

Its columns are vectors satisfying v = c @ k, where k[0], k[1], k[2] are all in [0, 1).

all_hnfs(det)

Enumerate all (3, 3) Hermite normal forms (HNFs) with given determinant.

Parameters:

Name Type Description Default
det int

The determinant of HNFs.

required

Returns:

Name Type Description
l (..., 3, 3) array of ints

Contains all HNFs with determinant det.

hnf(m, return_q=False)

Compute the Hermite normal form of integer matrix m.

Parameters:

Name Type Description Default
m (M, N) array of ints

The integer matrix to decompose, with positive determinant and M <= N.

required
return_q bool

Whether to return the unimodular matrix q.

False

Returns:

Name Type Description
h (M, N) array of ints

The column-style Hermite normal form of m.

q (N, N) array of ints

The unimodular matrix satisfying m = h @ q. Only returned if return_q is True.

lll_reduce(c, delta=0.75)

Return the LLL-reduced cell cc and the unimodular matrix q such that cc = c @ q

Parameters:

Name Type Description Default
c (3, 3) array

The cell to be reduced, whose columns are cell vectors.

required
delta float

The parameter for the LLL algorithm.

0.75

Returns:

Name Type Description
cc (3, 3) array

The LLL-reduced cell.

q (3, 3) array of ints

The unimodular matrix satisfying cc = c @ q.

voigt_to_tensor(voigt_matrix, cryst=None, tol=0.001, verbose=True)

Convert a Voigt-notation tensor to a rank-4 tensor, and symmetrize it according to the symmetry of cryst (if provided).

Parameters:

Name Type Description Default
voigt_matrix (6, 6) array

The elastic tensor, in Voigt notation (ordered as XX, YY, ZZ, YZ, ZX, XY).

required
cryst Cryst

The crystal structure, whose symmetry is used to symmetrize the elastic tensor.

None
tol float

The tolerance for symmetry finding.

0.001
verbose bool

Whether to print information about the symmetrization.

True

Returns:

Name Type Description
tensor (3, 3, 3, 3) array

The rank-4 elastic tensor.

pct_distance(c, pA, pB, p, ks, weights=None, l=2, min_t0=True, return_t0=False)

Return the shuffle distance of a PCT.

Parameters:

Name Type Description Default
c (3, 3) array

The lattice vectors of the crystal structure.

required
pA (3, Z) array

The fractional coordinates of the atoms in the initial and final structures, respectively.

required
pB (3, Z) array

The fractional coordinates of the atoms in the initial and final structures, respectively.

required
p (Z, ) array of ints

The permutation of the atoms.

required
ks (3, Z) array of ints

The class-wise translations (fractional coordinates) of the atoms in pB.

required
weights (Z, ) array of floats

The weights of each atom. If None, all atoms have the same weight.

None
l float

The l-norm to be used for distance calculation, must not be less than 1.

2
min_t0 bool

Set to True to minimize the shuffle distance by translating the final structure.

True
return_t0 bool

Whether to return the best overall translation if min_t0 is True.

False

Returns:

Name Type Description
distance float

The shuffle distance.

t0 (3, 1) array

The best overall translation, reshaped as a 3x1 matrix. Only returned if return_t0 is True.

standardize_imt(slm, gA, gB)

The standard SLM of the congruence class of slm.

Parameters:

Name Type Description Default
slm slm

(hA, hB, q), representing a SLM.

required
gA (..., 3, 3) array of ints

The rotation group of the initial crystal structure, whose elements are integer matrices under fractional coordinates.

required
gB (..., 3, 3) array of ints

The rotation group of the final crystal structure, whose elements are integer matrices under fractional coordinates.

required

Returns:

Name Type Description
slm0 slm

The standardized SLM.

strain_energy_func(elasticA, elasticB)

Return the strain energy density function.

Parameters:

Name Type Description Default
elasticA (3, 3, 3, 3) arrays

The rank-4 elastic tensors of the initial and final crystal structures, respectively.

required
elasticB (3, 3, 3, 3) arrays

The rank-4 elastic tensors of the initial and final crystal structures, respectively.

required

Returns:

Name Type Description
strain Callable

A function that takes a deformation gradient s and returns the strain energy density.

enumerate_imt(crystA, crystB, mu, max_strain, strain=rmss, tol=0.001, max_iter=1000, verbose=1)

Enumerating all IMTs of multiplicity mu with strain smaller than max_strain.

Parameters:

Name Type Description Default
crystA Cryst

The initial crystal structure, usually obtained by load_poscar.

required
crystB Cryst

The final crystal structure, usually obtained by load_poscar.

required
mu int

The multiplicity of SLMs to enumerate.

required
max_strain float

The maximum strain energy density, with the same units as strain.

required
strain Callable

How to quantify the strain, usually rmss or obtained via strain_energy_func().

rmss
tol float

The tolerance for determining the pure rotation group of the crystal structures.

0.001
max_iter int

The maximum number of consequtive iterations without finding any new SLMs.

1000
verbose int

The level of verbosity.

1

Returns:

Name Type Description
slmlist (N, 3, 3, 3) array of ints

The IMTs enumerated, where N is the number of noncongruent IMTs found.

optimize_pct_fixed(c, species, pA, pB, constraint=Constraint(set(), set()), weight_func=None, l=2)

Minimize the shuffle distance with variable PCT and fixed overall translation.

Parameters:

Name Type Description Default
c (3, 3) array

The base matrix of the shuffle lattice.

required
species (Z, ) array of str

The species of each atom.

required
pA (3, Z) array

The fractional coordinates of the atoms in the initial and final structures, respectively.

required
pB (3, Z) array

The fractional coordinates of the atoms in the initial and final structures, respectively.

required
constraint Constraint

The constraint on the permutation.

Constraint(set(), set())
weight_func dict

The weight function, with keys as atomic species. If None, all atoms have the same weight.

None
l float

The l-norm to be used for distance calculation, must not be less than 1.

2

Returns:

Name Type Description
d_hat float

The least shuffle distance.

p (Z, ) array of ints

The permutation part of the PCT with the least shuffle distance.

ks (3, Z) array of ints

The class-wise translation part of the PCT with the least shuffle distance.

optimize_pct(crystA, crystB, slm, constraint=Constraint(set(), set()), weight_func=None, l=2, t_grid=64)

Minimize the shuffle distance with variable PCT and overall translation.

Parameters:

Name Type Description Default
crystA cryst

The initial crystal structure, usually obtained by load_poscar.

required
crystB cryst

The final crystal structure, usually obtained by load_poscar.

required
slm slm

The SLM, represented by a triplet of integer matrices like (hA, hB, q).

required
constraint 2-tuple of sets

The permutation constraint used in the Murty's algorithm.

Constraint(set(), set())
weight_func dict
None

pct_fill(crystA, crystB, slm, max_d, p, ks0=None, weight_func=None, l=2, warning_threshold=5000)

Returns all class-wise translations with permutation p that has d <= max_d, including ks0.

murty_split(node, p)

cong_permutations(p, crystA, crystB, slm)

enumerate_pct(crystA, crystB, slm, max_d, weight_func=None, l=2, t_grid=16, non_cong=True, verbose=1, warning_threshold=5000)

optimize_ct_fixed(c, pA, pB, p, weights=None, l=2)

optimize_ct(crystA, crystB, slm, p, weight_func=None, l=2, t_grid=64)

enumerate_rep_csm(crystA, crystB, max_mu, max_strain, strain=rmss, tol=0.001, max_iter=1000, weight_func=None, l=2, t_grid=64, verbose=1)

Enumerating all representative CSMs with multiplicity and strain not larger than max_mu and max_strain.

Parameters:

Name Type Description Default
crystA Cryst

The initial crystal structure.

required
crystB Cryst

The final crystal structure.

required
max_mu int

The maximum multiplicity of the IMTs to be enumerated.

required
max_strain float

The maximum value of the function strain.

required
strain Callable

How to quantify the strain, usually rmss or obtained via strain_energy_func().

rmss
tol float

The tolerance for spglib symmetry detection.

0.001
max_iter int

The maximum iteration number for IMT generation.

1000
weight_func dict

A dictionary of atomic weights for each species. If None, all atoms have the same weight.

None
l float

The type of norm to be used for distance calculation.

2
t_grid int

The number of grid points for PCT optimization.

64
verbose int

The level of verbosity.

1

Returns:

Name Type Description
slmlist (N, 3, 3, 3) array of ints

The IMTs enumerated.

pct_arrs list of (..., 4) arrays of ints

The PCTs of representative CSMs, where ... is the number of CSMs with multiplicity mu.

mulist (N,) array of ints

The multiplicities of the representative CSMs.

strainlist (N,) array of floats

The strain values of the representative CSMs.

d0list (N, ) array of floats

The shuffle distances of the representative CSMs.

enumerate_all_csm(crystA, crystB, max_mu, max_strain, max_d, strain=rmss, tol=0.001, max_iter=1000, weight_func=None, l=2, t_grid=16, verbose=1)

Enumerating all CSMs with multiplicity, strain, and shuffle distance not larger than max_mu, max_strain, and max_d.

Parameters:

Name Type Description Default
crystA Cryst

The initial crystal structure.

required
crystB Cryst

The final crystal structure.

required
max_mu int

The maximum multiplicity.

required
max_strain float

The maximum value of the function strain.

required
max_d float

The maximum shuffle distance.

required
strain Callable

How to quantify the strain, usually rmss or obtained via strain_energy_func().

rmss
tol float

The tolerance for spglib symmetry detection.

0.001
max_iter int

The maximum iteration number for IMT generation.

1000
weight_func dict

A dictionary of atomic weights for each species. If None, all atoms have the same weight.

None
l float

The type of norm to be used for distance calculation.

2
t_grid int

The number of grid points for PCT optimization.

16
verbose int

The level of verbosity.

1

Returns:

Name Type Description
slmlist (N, 3, 3, 3) array of ints

The IMTs enumerated.

slm_ind (N,) array of ints

The indices of the IMTs of the CSMs.

pct_arrs list of (..., 4) arrays of ints

The PCTs of the CSMs, where ... is the number of CSMs with multiplicity mu.

mulist (N,) array of ints

The multiplicities of the CSMs.

strainlist (N,) array of floats

The strain values of the CSMs.

dlist (N, ) array of floats

The shuffle distances of the CSMs.

decompose_cryst(cryst_sup, tol=0.001)

Compute the primitive cryst_prim and integer matrix m such that c_sup = c_prim @ m and p_prim ∈ m @ p_sup (mod 1)

csm_to_cryst(crystA, crystB, slm, p, ks, orientation='norot', use_medium_cell=False, min_t0=False, weight_func=None, l=2)

Convert the IMT and PCT representation of a CSM to a pair of crysts.

Parameters:

Name Type Description Default
crystA cryst

The initial and final structures.

required
crystB cryst

The initial and final structures.

required
slm slm

The SLM of the CSM.

required
p (Z, ) array of ints

The permutaion of the shuffle.

required
ks (3, Z) array of ints

The lattice-vector translations of the shuffle.

required
orientation str

The orientation of the final structure, either 'norot' or 'uspfixed', which means that the deformation is rotation-free or fixing the uniformly scaled plane (USP). When 'uspfixed', two final structures are returned since there are two USPs.

'norot'
use_medium_cell bool

Whether to return the average cell of crystA_sup and crystB_sup.

False
min_t0 bool

Whether to use optimal translation of crystB_sup_norot to minimize the shuffle distance.

False
weight_func dict

The weight function used when minimizing the shuffle distance, with keys as atomic species. If None, all atoms have the same weight.

None
l float

The l-norm to be used for distance calculation, must not be less than 1.

2

Returns:

Name Type Description
crystA_sup cryst

The initial structure.

crystB_sup_norot cryst

The final structure, whose lattice vectors and atoms are matched to crystA_sup according to the CSM, with rotation-free orientation.

crystB_sup_uspfixed_1, crystB_sup_uspfixed_2 : cryst

The final structure, whose lattice vectors and atoms are matched to crystA_sup according to the CSM, with uniformly scaled plane fixed.

cryst_to_csm(crystA_sup, crystB_sup, tol=0.001)

Return the primitive crysts and IMT and PCT representation of a CSM determined by a pair of crysts.

csm_distance(crystA, crystB, slm, p, ks, weight_func=None, l=2, min_t0=True, return_t0=False)

Return the shuffle distance of a CSM.

primitive_shuffle(crystA, crystB, slm0, p0, ks0)

Identify the primitive (with minimal multiplicity) IMT and PCT representations of a CSM.

orient_matrix(vi, vf, wi, wf)

Rotation matrix r such that r @ vi || vf and r @ wi || wf.

Parameters:

Name Type Description Default
vi (3,) array_like

Vectors (cartesian coordinates) that satisfy r @ vi || vf.

required
vf (3,) array_like

Vectors (cartesian coordinates) that satisfy r @ vi || vf.

required
wi (3,) array_like

Vectors (cartesian coordinates) that satisfy r @ wi || wf.

required
wf (3,) array_like

Vectors (cartesian coordinates) that satisfy r @ wi || wf.

required

Returns:

Name Type Description
r (3, 3) array

A rotation matrix representing the given orientation relationship.

rot_usp(s)

The rotation that the uniformly scaled plane undergoes.

miller_to_vec(cryst, hkl, tol=0.001)

Convert Miller indices to cartesian coordinates.

Parameters:

Name Type Description Default
cryst cryst

The crystal structure, usually obtained by load_poscar.

required
hkl (3,) array_like

Miller indices.

required
tol float

Tolerance for determining the symmetry of the crystal.

0.001

Returns:

Name Type Description
vec (3,) array

Cartesian coordinates of the Miller index.

deviation_angle(crystA, crystB, slmlist, r, orientation)

Calculate how much each SLM in slmlist differ from a given orientation relationship.

Parameters:

Name Type Description Default
crystA cryst

The initial crystal structure, usually obtained by load_poscar.

required
crystB cryst

The final crystal structure, usually obtained by load_poscar.

required
slmlist (N, 3, 3, 3) array of ints

A list of SLMs, each represented by a triplet of integer matrices like (hA, hB, q).

required
r (3, 3) array

A rotation matrix representing the given orientation relationship.

required
orientation str

The orientation of the final structure, either 'norot' or 'uspfixed', which means that the deformation is rotation-free or fixing the uniformly scaled plane (USP). When 'uspfixed', two final structures are returned since there are two USPs.

required

Returns:

Name Type Description
anglelist (N,) array

Contains rotation angles that measure the difference of each SLM and the given orientation.

visualize_slmlist(filename, strainlist, dlist, colorlist=None, cmap=plt.get_cmap('viridis'), cbarlabel=None)

Scatter plot of the CSMs with colorbar.

visualize_csm(crystA, crystB, slm, p, ks, weight_func=None, l=2, tol=0.001, cluster_size=1.2, show_conventional=True, label=None)

Use with %matplotlib widget in Jupyter notebook (need to install ipympl) to interactively visualize the shuffling process.

save_poscar(filename, cryst, crystname=None)

Save the crystal structure to a POSCAR file.

Parameters:

Name Type Description Default
filename str

The name of the file to save, must not already exist in current directory. If filename = None, a string will be returned instead.

required
cryst cryst

The crystal structure to be saved, consisting of the lattice vectors, species, and positions.

required
crystname str

A system description to write to the comment line of the POSCAR file. If crystname = None, filename will be used.

None

save_xdatcar(filename, crystA_sup, crystB_sup, n_im=50, crystname=None)

Save the linear interpolation between crystA and crystB to a single XDATCAR file.

Parameters:

Name Type Description Default
filename str

The name of the file to save, must not already exist in current directory.

required
crystA_sup cryst

The initial and final crystal structures with specified atomic correspondence, usually obtained by minimize_rmsd.

required
crystB_sup cryst

The initial and final crystal structures with specified atomic correspondence, usually obtained by minimize_rmsd.

required
n_im int

Number of images to generate.

50
crystname str

A system description to write to the comment line of the POSCAR file. If crystname = None, filename will be used.

None