pyxtal.lattice module

class pyxtal.lattice.Lattice(ltype, volume=None, matrix=None, PBC=[1, 1, 1], **kwargs)[source]

Bases: object

Class for storing and generating crystal lattices. Allows for specification of constraint values. Lattice types include triclinic, monoclinic, orthorhombic, tetragonal, trigonal, hexagonal, cubic, spherical, and ellipsoidal. The last two are used for generating point group structures, and do not actually represent a parallelepiped lattice.

Parameters:
  • ltype – a string representing the type of lattice (from the above list)

  • volume – the volume, in Angstroms cubed, of the lattice

  • matrix – matrix in 3*3 form

  • PBC – A periodic boundary condition list, where 1 is periodic, Ex: [1, 1, 1] -> 3d periodicity, [0, 0, 1] -> periodic at z axis

  • kwargs

    various values which may be defined. If none are defined, random ones will be generated. Values will be passed to generate_lattice. Options include: area: The cross-sectional area (in Ang^2). Only for 1D crystals thickness: The cell’s thickness (in Angstroms) for 2D crystals unique_axis: The unique axis for certain symmetry (and especially

    layer) groups. Because the symmetry operations are not also transformed, you should use the default values for random crystal generation

    random: If False, keeps the stored values for the lattice geometry

    even upon applying reset_matrix. To alter the matrix, use set_matrix() or set_para

    ’unique_axis’: the axis (‘a’, ‘b’, or ‘c’) which is not symmetrically

    equivalent to the other two

    ’min_l’: the smallest allowed cell vector. The smallest vector must

    be larger than this.

    ’mid_l’: the second smallest allowed cell vector. The second

    smallest vector must be larger than this.

    ’max_l’: the third smallest allowed cell vector. The largest cell

    vector must be larger than this.

    ’allow_volume_reset’: a bool stating whether or not the volume

    should be reset during each crystal generation attempt

add_vacuum(coor, frac=True, vacuum=15, PBC=[0, 0, 0])[source]

Adds space above and below a 2D or 1D crystal.

Parameters:
  • coor – the relative coordinates of the crystal

  • vacuum – the amount of space, in Angstroms, to add above and below

  • PBC – A periodic boundary condition list, Ex: [1,1,1] -> full 3d periodicity, [0,0,1] -> periodicity along the z axis

Returns:

The transformed lattice and coordinates after the vacuum is added

apply_transformation(trans)[source]

Optimize the lattice’s inclination angles

check_mismatch(trans, l_type, tol=1.0, a_tol=10)[source]

check if the lattice mismatch is big after a transformation This is mostly used in supergroup function QZ: to fix ===============

Parameters:
  • trans – 3*3 matrix

  • l_type – lattice_type like orthrhombic

  • tol – tolerance in a, b, c

  • a_tol – tolerance in alpha, beta, gamma

Returns:

True or False

copy()[source]

simply copy the structure

encode()[source]
find_transition_to_orthoslab(c=(0, 0, 1), a=(1, 0, 0), m=5)[source]

Create the slab model with an approximate orthogonal box shape

classmethod from_1d_representation(v, ltype)[source]
classmethod from_matrix(matrix, reset=True, shape='upper', ltype='triclinic', PBC=[1, 1, 1], **kwargs)[source]

Creates a Lattice object from a 3x3 cell matrix. Additional keywords are available. Unless specified by the keyword random=True, does not create a new matrix upon calling reset_matrix. This allows for a random crystals with a specific choice of unit cell.

Parameters:
  • matrix – a 3x3 real matrix (numpy array or nested list) for the cell

  • ltype – the lattice type (“cubic, tetragonal, etc.”). Also can be - “spherical”, confines points to lie within a sphere, - “ellipsoidal”, points to lie within an ellipsoid (about z axis)

  • PBC – A periodic boundary condition list, where 1 is periodic Ex: [1,1,1] -> full 3d periodicity, [0,0,1] -> periodicity at z.

  • kwargs – various values which may be defined. Random ones if None Values will be passed to generate_lattice. Options include: area: The cross-sectional area (in Ang^2) for 1D crystals `thickness: The cell’s thickness (in Ang) for 2D crystals unique_axis: The unique axis for layer groups. random: If False, keeps the stored values for the lattice geometry even applying reset_matrix. To alter the matrix, use set_matrix() or set_para ‘unique_axis’: the axis (‘a’, ‘b’, or ‘c’) which is unique. ‘min_l’: the smallest allowed cell vector. ‘mid_l’: the second smallest allowed cell vector. ‘max_l’: the third smallest allowed cell vector.

Returns:

a Lattice object with the specified parameters

classmethod from_para(a, b, c, alpha, beta, gamma, ltype='triclinic', radians=False, PBC=[1, 1, 1], factor=1.0, **kwargs)[source]

Creates a Lattice object from 6 lattice parameters. Additional keyword arguments are available. Unless specified by the keyword random=True, does not create a new matrix upon calling reset_matrix. This allows for generation of random crystals with a specific choice of unit cell.

Parameters:
  • a – The length (in Angstroms) of the unit cell vectors

  • b – The length (in Angstroms) of the unit cell vectors

  • c – The length (in Angstroms) of the unit cell vectors

  • alpha – the angle (in degrees) between the b and c vectors

  • beta – the angle (in degrees) between the a and c vectors

  • gamma – the angle (in degrees) between the a and b vectors

  • ltype – the lattice type (“cubic, tetragonal, etc.”). Also available are “spherical”, which confines generated points to lie within a sphere, and “ellipsoidal”, which confines generated points to lie within an ellipse (oriented about the z axis)

  • radians – whether or not to use radians (instead of degrees) for the lattice angles

  • PBC – A periodic boundary condition list, where 1 means periodic, 0 means not periodic. Ex: [1,1,1] -> full 3d periodicity, [0,0,1] -> periodicity along the z axis

  • kwargs

    various values which may be defined. If none are defined, random ones will be generated. Values will be passed to generate_lattice. Options include: area: The cross-sectional area (in Angstroms squared). Only used

    to generate 1D crystals

    thickness: The unit cell’s non-periodic thickness (in Angstroms).

    Only used to generate 2D crystals

    unique_axis: The unique axis for certain symmetry (and especially

    layer) groups. Because the symmetry operations are not also transformed, you should use the default values for random crystal generation

    random: If False, keeps the stored values for the lattice geometry

    even upon applying reset_matrix. To alter the matrix, use set_matrix() or set_para

    ’unique_axis’: the axis (‘a’, ‘b’, or ‘c’) which is unique ‘min_l’: the smallest allowed cell vector. ‘mid_l’: the second smallest allowed cell vector. ‘max_l’: the third smallest allowed cell vector.

Returns:

a Lattice object with the specified parameters

generate_matrix()[source]

Generates a 3x3 matrix for a lattice based on the lattice type and volume

generate_para()[source]
generate_point()[source]
get_diff(l_ref)[source]

get the difference in length, angle, and check if switch is needed

classmethod get_dofs(ltype)[source]

get the number of degree of freedom

get_lengths()[source]
get_matrix(shape='upper')[source]

Returns a 3x3 numpy array representing the lattice vectors.

get_para(degree=False)[source]

Returns a tuple of lattice parameters.

get_permutation_matrices()[source]

Return the possible permutation matrices that donot violate the symmetry

get_transformation_matrices()[source]

Return possible transformation matrices that donot violate the symmetry

get_worst_angle()[source]

return the worst inclination angle difference w.r.t 90 degree

is_valid_matrix()[source]

check if the cell parameter is reasonable or not

mutate(degree=0.2, frozen=False)[source]

mutate the lattice object

optimize_multi(iterations=5)[source]

Optimize the lattice if the cell has a bad inclination angles

Parameters:
  • iterations – maximum number of iterations

  • force – whether or not do the early termination

Returns:

the optimized lattice

optimize_once(reset=False)[source]

Optimize the lattice’s inclination angles

reset_matrix(shape='upper')[source]
scale(factor=1.1)[source]
search_transformation(lat_ref, d_tol=1.0, f_tol=0.1)[source]

search the closest match to the reference lattice object

Parameters:
  • lat_ref – reference lattice object

  • d_tol – tolerance in angle

  • f_tol

  • a_tol

Returns:

a two steps of transformation matrix if the match is possible

search_transformations(lat_ref, d_tol=1.0, f_tol=0.1)[source]

search the closest match to the reference lattice object

Parameters:
  • lat_ref – reference lattice object

  • d_tol – tolerance in angle

  • f_tol

  • a_tol

Returns:

a two steps of transformation matrix if the match is possible

set_matrix(matrix=None)[source]
set_para(para=None, radians=False)[source]
set_volume(volume)[source]
standardize()[source]

Force the angle to be smaller than 90 degree

swap_angle(random=True, ids=None)[source]

If the angle is not 90. There will be two equivalent versions e.g., 80 and 100.

swap_axis(random=False, ids=None)[source]

For the lattice

transform(trans_mat=array([[1., 0., 0.], [0., 1., 0.], [0., 0., 1.]]), reset=False)[source]

Optimize the lattice’s inclination angles If reset is False, may return negative lattice

transform_multi(trans, reset=True)[source]

Optimize the lattice’s inclination angles

pyxtal.lattice.gaussian(min, max, sigma=3.0)[source]

Choose a random number from a Gaussian probability distribution centered between min and max. sigma is the number of standard deviations that min and max are away from the center. Thus, sigma is also the largest possible number of standard deviations corresponding to the returned value. sigma=2 corresponds to a 95.45% probability of choosing a number between min and max.

Parameters:
  • min – the minimum acceptable value

  • max – the maximum acceptable value

  • sigma – the number of standard deviations between the center and min/max

Returns:

a value chosen randomly between min and max

pyxtal.lattice.generate_cellpara(ltype, volume, minvec=1.2, minangle=0.5235987755982988, max_ratio=10.0, maxattempts=100, **kwargs)[source]

Generates the cell parameter (a, b, c, alpha, beta, gamma) according to the space group symmetry and number of atoms. If the spacegroup has centering, we will transform to conventional cell setting. If the generated lattice does not meet the minimum angle and vector requirements, we try to generate a new one, up to maxattempts times.

Parameters:
  • volume – volume of the conventional unit cell

  • minvec – minimum allowed lattice vector length (among a, b, and c)

  • minangle – minimum allowed lattice angle (among alpha, beta, and gamma)

  • max_ratio – largest allowed ratio of two lattice vector lengths

  • maxattempts – the maximum number of attempts for generating a lattice

  • kwargs – a dictionary of optional values. These include: ‘unique_axis’: the axis (‘a’, ‘b’, or ‘c’) which is unique. ‘min_l’: the smallest allowed cell vector. ‘mid_l’: the second smallest allowed cell vector. ‘max_l’: the third smallest allowed cell vector.

Returns:

a 6-length array representing the lattice of the unit cell. If generation fails, outputs a warning message and returns empty

pyxtal.lattice.generate_cellpara_0D(ltype, volume, area=None, minvec=1.2, max_ratio=10.0, maxattempts=100, **kwargs)[source]

Generates a cell parameter (a, b, c, alpha, beta, gamma) according to the point group symmetry and number of atoms. If the generated lattice does not meet the minimum angle and vector requirements, we try to generate a new one, up to maxattempts times.

Parameters:
  • num – number of the Rod group

  • volume – volume of the lattice

  • area – cross-sectional area of the unit cell in Angstroms squared. If set to None, a value is chosen automatically

  • minvec – minimum allowed lattice vector length (among a, b, and c)

  • max_ratio – largest allowed ratio of two lattice vector lengths

  • maxattempts – the maximum number of attempts for generating a lattice

  • kwargs – a dictionary of optional values. Only used for ellipsoidal lattices, which pass the value to generate_lattice. They include: ‘unique_axis’: the axis (‘a’, ‘b’, or ‘c’) which is unique. ‘min_l’: the smallest allowed cell vector. ‘mid_l’: the second smallest allowed cell vector. ‘max_l’: the third smallest allowed cell vector.

Returns:

a 3x3 matrix representing the lattice vectors of the unit cell. If generation fails, outputs a warning message and returns empty

pyxtal.lattice.generate_cellpara_1D(ltype, volume, area=None, minvec=1.2, minangle=0.5235987755982988, max_ratio=10.0, maxattempts=100, **kwargs)[source]

Generates a cell parameter (a, b, c, alpha, beta, gamma) according to the rod group symmetry and number of atoms. If the rod group has centering, we will transform to conventional cell setting. If the generated lattice does not meet the minimum angle and vector requirements, we try to generate a new one, up to maxattempts times.

Note: The monoclinic Rod groups have different unique axes. Groups 3-7

have unique axis a, while 8-12 have unique axis c. We use periodic axis c for all Rod groups.

Parameters:
  • num – number of the Rod group

  • volume – volume of the lattice

  • area – cross-sectional area of the unit cell in Angstroms squared. If set to None, a value is chosen automatically

  • minvec – minimum allowed lattice vector length (among a, b, and c)

  • minangle – minimum allowed lattice angle (among alpha, beta, and gamma)

  • max_ratio – largest allowed ratio of two lattice vector lengths

  • maxattempts – the maximum number of attempts for generating a lattice

  • kwargs – a dictionary of optional values. These include: ‘unique_axis’: the axis (‘a’, ‘b’, or ‘c’) which is unique. ‘min_l’: the smallest allowed cell vector. ‘mid_l’: the second smallest allowed cell vector. ‘max_l’: the third smallest allowed cell vector.

Returns:

a 6-length array representing the lattice of the unit cell. If generation fails, outputs a warning message and returns empty

pyxtal.lattice.generate_cellpara_2D(ltype, volume, thickness=None, minvec=1.2, minangle=0.5235987755982988, max_ratio=10.0, maxattempts=100, **kwargs)[source]

Generates the cell parameter (a, b, c, alpha, beta, gamma) according to the layer group symmetry and number of atoms. If the layer group has centering, we will transform to conventional cell setting. If the generated lattice does not meet the minimum angle and vector requirements, we try to generate a new one, up to maxattempts times.

Note: The monoclinic layer groups have different unique axes. Groups 3-7

have unique axis c, while 8-18 have unique axis a. We use non-periodic axis c for all layer groups.

Parameters:
  • num – International number of the space group

  • volume – volume of the lattice

  • thickness – 3rd-dimensional thickness of the unit cell. If set to None, a thickness is chosen automatically

  • minvec – minimum allowed lattice vector length (among a, b, and c)

  • minangle – minimum allowed lattice angle (among alpha, beta, and gamma)

  • max_ratio – largest allowed ratio of two lattice vector lengths

  • maxattempts – the maximum number of attempts for generating a lattice

  • kwargs – a dictionary of optional values. These include: ‘unique_axis’: the axis (‘a’, ‘b’, or ‘c’) which is unique. ‘min_l’: the smallest allowed cell vector. ‘mid_l’: the second smallest allowed cell vector. ‘max_l’: the third smallest allowed cell vector.

Returns:

a 6-length representing the lattice vectors of the unit cell. If generation fails, outputs a warning message and returns empty

pyxtal.lattice.matrix2para(matrix, radians=True)[source]

Given a 3x3 matrix representing a unit cell, outputs a list of lattice parameters.

Parameters:
  • matrix – a 3x3 array or list, where the first, second, and third rows represent the a, b, and c vectors respectively

  • radians – if True, outputs angles in radians. If False, outputs in degrees

Returns:

a 1x6 list of lattice parameters [a, b, c, alpha, beta, gamma]. a, b, and c are the length of the lattice vectos, and alpha, beta, and gamma are the angles between these vectors (in radians by default)

pyxtal.lattice.para2matrix(cell_para, radians=True, format='upper')[source]

Given a set of lattic parameters, generates a matrix representing the lattice vectors

Parameters:
  • cell_para – a 1x6 list of lattice parameters [a, b, c, alpha, beta, gamma]. a, b, and c are the length of the lattice vectos, and alpha, beta, and gamma are the angles between these vectors. Can be generated by matrix2para

  • radians – if True, lattice parameters should be in radians. If False, lattice angles should be in degrees

  • format – a string (‘lower’, ‘symmetric’, or ‘upper’) for the type of matrix to be output

Returns:

a 3x3 matrix representing the unit cell. By default (format=’lower’), the a vector is aligined along the x-axis, and the b vector is in the y-z plane

pyxtal.lattice.random_shear_matrix(width=1.0, unitary=False)[source]

Generate a random symmetric shear matrix with Gaussian elements. If unitary is True, normalize to determinant 1

Parameters:
  • width – the width of the normal distribution to use when choosing values. Passed to np.random.normal

  • unitary – whether or not to normalize the matrix to determinant 1

Returns:

a 3x3 numpy array of floats

pyxtal.lattice.random_vector(minvec=[0.0, 0.0, 0.0], maxvec=[1.0, 1.0, 1.0], width=0.35, unit=False)[source]

Generate a random vector for lattice constant generation. The ratios between x, y, and z of the returned vector correspond to the ratios between a, b, and c. Results in a Gaussian distribution of the natural log of the ratios.

Parameters:
  • minvec – the bottom-left-back minimum point which can be chosen

  • maxvec – the top-right-front maximum point which can be chosen

  • width – the width of the normal distribution to use when choosing values. Passed to np.random.normal

  • unit – whether or not to normalize the vector to determinant 1

Returns:

a 1x3 numpy array of floats