pymatgen.core.lattice module

class Lattice(matrix: Union[List[float], List[List[float]], numpy.ndarray])[source]

Bases: monty.json.MSONable

A lattice object. Essentially a matrix with conversion matrices. In general, it is assumed that length units are in Angstroms and angles are in degrees unless otherwise stated.

Create a lattice from any sequence of 9 numbers. Note that the sequence is assumed to be read one row at a time. Each row represents one lattice vector.

Parameters

matrix – Sequence of numbers in any form. Examples of acceptable input. i) An actual numpy array. ii) [[1, 0, 0], [0, 1, 0], [0, 0, 1]] iii) [1, 0, 0 , 0, 1, 0, 0, 0, 1] iv) (1, 0, 0, 0, 1, 0, 0, 0, 1) Each row should correspond to a lattice vector. E.g., [[10, 0, 0], [20, 10, 0], [0, 0, 30]] specifies a lattice with lattice vectors [10, 0, 0], [20, 10, 0] and [0, 0, 30].

a

a lattice parameter.

abc

Lengths of the lattice vectors, i.e. (a, b, c)

alpha

Angle alpha of lattice in degrees.

angles

Returns the angles (alpha, beta, gamma) of the lattice.

as_dict(verbosity: int = 0) → Dict[KT, VT][source]

“” Json-serialization dict representation of the Lattice.

Parameters

verbosity (int) – Verbosity level. Default of 0 only includes the matrix representation. Set to 1 for more details.

b

b lattice parameter.

beta

Angle beta of lattice in degrees.

c

c lattice parameter.

copy()[source]

Deep copy of self.

static cubic(a: float)[source]

Convenience constructor for a cubic lattice.

Parameters

a (float) – The a lattice parameter of the cubic cell.

Returns

Cubic lattice of dimensions a x a x a.

d_hkl(miller_index: Union[List[float], numpy.ndarray]) → float[source]

Returns the distance between the hkl plane and the origin

Parameters

miller_index ([h,k,l]) – Miller index of plane

Returns

d_hkl (float)

dot(coords_a: Union[List[float], numpy.ndarray], coords_b: Union[List[float], numpy.ndarray], frac_coords: bool = False) → numpy.ndarray[source]

Compute the scalar product of vector(s).

Parameters
  • coords_b (coords_a,) – Array-like objects with the coordinates.

  • frac_coords (bool) – Boolean stating whether the vector corresponds to fractional or cartesian coordinates.

Returns

one-dimensional numpy array.

find_all_mappings(other_lattice: pymatgen.core.lattice.Lattice, ltol: float = 1e-05, atol: float = 1, skip_rotation_matrix: bool = False) → Iterator[Tuple[pymatgen.core.lattice.Lattice, Optional[numpy.ndarray], numpy.ndarray]][source]

Finds all mappings between current lattice and another lattice.

Parameters
  • other_lattice (Lattice) – Another lattice that is equivalent to this one.

  • ltol (float) – Tolerance for matching lengths. Defaults to 1e-5.

  • atol (float) – Tolerance for matching angles. Defaults to 1.

  • skip_rotation_matrix (bool) – Whether to skip calculation of the rotation matrix

Yields

(aligned_lattice, rotation_matrix, scale_matrix) if a mapping is found. aligned_lattice is a rotated version of other_lattice that has the same lattice parameters, but which is aligned in the coordinate system of this lattice so that translational points match up in 3D. rotation_matrix is the rotation that has to be applied to other_lattice to obtain aligned_lattice, i.e., aligned_matrix = np.inner(other_lattice, rotation_matrix) and op = SymmOp.from_rotation_and_translation(rotation_matrix) aligned_matrix = op.operate_multi(latt.matrix) Finally, scale_matrix is the integer matrix that expresses aligned_matrix as a linear combination of this lattice, i.e., aligned_matrix = np.dot(scale_matrix, self.matrix)

None is returned if no matches are found.

find_mapping(other_lattice: pymatgen.core.lattice.Lattice, ltol: float = 1e-05, atol: float = 1, skip_rotation_matrix: bool = False) → Optional[Tuple[pymatgen.core.lattice.Lattice, Optional[numpy.ndarray], numpy.ndarray]][source]

Finds a mapping between current lattice and another lattice. There are an infinite number of choices of basis vectors for two entirely equivalent lattices. This method returns a mapping that maps other_lattice to this lattice.

Parameters
  • other_lattice (Lattice) – Another lattice that is equivalent to this one.

  • ltol (float) – Tolerance for matching lengths. Defaults to 1e-5.

  • atol (float) – Tolerance for matching angles. Defaults to 1.

Returns

(aligned_lattice, rotation_matrix, scale_matrix) if a mapping is found. aligned_lattice is a rotated version of other_lattice that has the same lattice parameters, but which is aligned in the coordinate system of this lattice so that translational points match up in 3D. rotation_matrix is the rotation that has to be applied to other_lattice to obtain aligned_lattice, i.e., aligned_matrix = np.inner(other_lattice, rotation_matrix) and op = SymmOp.from_rotation_and_translation(rotation_matrix) aligned_matrix = op.operate_multi(latt.matrix) Finally, scale_matrix is the integer matrix that expresses aligned_matrix as a linear combination of this lattice, i.e., aligned_matrix = np.dot(scale_matrix, self.matrix)

None is returned if no matches are found.

classmethod from_dict(d: Dict[KT, VT], fmt: str = None, **kwargs)[source]

Create a Lattice from a dictionary containing the a, b, c, alpha, beta, and gamma parameters if fmt is None.

If fmt == “abivars”, the function build a Lattice object from a dictionary with the Abinit variables acell and rprim in Bohr. If acell is not given, the Abinit default is used i.e. [1,1,1] Bohr

Example

Lattice.from_dict(fmt=”abivars”, acell=3*[10], rprim=np.eye(3))

static from_lengths_and_angles(abc: List[float], ang: List[float])[source]

Create a Lattice using unit cell lengths and angles (in degrees).

Parameters
  • abc (3x1 array) – Lattice parameters, e.g. (4, 4, 5).

  • ang (3x1 array) – Lattice angles in degrees, e.g., (90,90,120).

Returns

A Lattice with the specified lattice parameters.

classmethod from_parameters(a: float, b: float, c: float, alpha: float, beta: float, gamma: float, vesta: bool = False)[source]

Create a Lattice using unit cell lengths and angles (in degrees).

Parameters
  • a (float) – a lattice parameter.

  • b (float) – b lattice parameter.

  • c (float) – c lattice parameter.

  • alpha (float) – alpha angle in degrees.

  • beta (float) – beta angle in degrees.

  • gamma (float) – gamma angle in degrees.

  • vesta – True if you import Cartesian coordinates from VESTA.

Returns

Lattice with the specified lattice parameters.

gamma

Angle gamma of lattice in degrees.

get_all_distances(fcoords1: Union[List[float], numpy.ndarray, List[Union[List[float], numpy.ndarray]]], fcoords2: Union[List[float], numpy.ndarray, List[Union[List[float], numpy.ndarray]]]) → numpy.ndarray[source]

Returns the distances between two lists of coordinates taking into account periodic boundary conditions and the lattice. Note that this computes an MxN array of distances (i.e. the distance between each point in fcoords1 and every coordinate in fcoords2). This is different functionality from pbc_diff.

Parameters
  • fcoords1 – First set of fractional coordinates. e.g., [0.5, 0.6, 0.7] or [[1.1, 1.2, 4.3], [0.5, 0.6, 0.7]]. It can be a single coord or any array of coords.

  • fcoords2 – Second set of fractional coordinates.

Returns

2d array of cartesian distances. E.g the distance between fcoords1[i] and fcoords2[j] is distances[i,j]

get_brillouin_zone() → List[List[numpy.ndarray]][source]

Returns the Wigner-Seitz cell for the reciprocal lattice, aka the Brillouin Zone.

Returns

A list of list of coordinates. Each element in the list is a “facet” of the boundary of the Brillouin Zone. For instance, a list of four coordinates will represent a square facet.

get_cartesian_coords(fractional_coords: Union[List[float], numpy.ndarray]) → numpy.ndarray[source]

Returns the cartesian coordinates given fractional coordinates.

Parameters

fractional_coords (3x1 array) – Fractional coords.

Returns

Cartesian coordinates

get_distance_and_image(frac_coords1: Union[List[float], numpy.ndarray], frac_coords2: Union[List[float], numpy.ndarray], jimage: Union[List[int], numpy.ndarray, None] = None) → Tuple[float, numpy.ndarray][source]

Gets distance between two frac_coords assuming periodic boundary conditions. If the index jimage is not specified it selects the j image nearest to the i atom and returns the distance and jimage indices in terms of lattice vector translations. If the index jimage is specified it returns the distance between the frac_coords1 and the specified jimage of frac_coords2, and the given jimage is also returned.

Parameters
  • fcoords1 (3x1 array) – Reference fcoords to get distance from.

  • fcoords2 (3x1 array) – fcoords to get distance from.

  • jimage (3x1 array) – Specific periodic image in terms of lattice translations, e.g., [1,0,0] implies to take periodic image that is one a-lattice vector away. If jimage is None, the image that is nearest to the site is found.

Returns

distance and periodic lattice translations of the other site for which the distance applies. This means that the distance between frac_coords1 and (jimage + frac_coords2) is equal to distance.

Return type

(distance, jimage)

get_frac_coords_from_lll(lll_frac_coords: Union[List[float], numpy.ndarray]) → numpy.ndarray[source]

Given fractional coordinates in the lll basis, returns corresponding fractional coordinates in the lattice basis.

get_fractional_coords(cart_coords: Union[List[float], numpy.ndarray]) → numpy.ndarray[source]

Returns the fractional coordinates given cartesian coordinates.

Parameters

cart_coords (3x1 array) – Cartesian coords.

Returns

Fractional coordinates.

get_lll_frac_coords(frac_coords: Union[List[float], numpy.ndarray]) → numpy.ndarray[source]

Given fractional coordinates in the lattice basis, returns corresponding fractional coordinates in the lll basis.

get_lll_reduced_lattice(delta: float = 0.75) → pymatgen.core.lattice.Lattice[source]
get_miller_index_from_coords(coords: Union[List[float], numpy.ndarray], coords_are_cartesian: bool = True, round_dp: int = 4, verbose: bool = True) → Tuple[int, int, int][source]

Get the Miller index of a plane from a list of site coordinates.

A minimum of 3 sets of coordinates are required. If more than 3 sets of coordinates are given, the best plane that minimises the distance to all points will be calculated.

Parameters
  • coords (iterable) – A list or numpy array of coordinates. Can be cartesian or fractional coordinates. If more than three sets of coordinates are provided, the best plane that minimises the distance to all sites will be calculated.

  • coords_are_cartesian (bool, optional) – Whether the coordinates are in cartesian space. If using fractional coordinates set to False.

  • round_dp (int, optional) – The number of decimal places to round the miller index to.

  • verbose (bool, optional) – Whether to print warnings.

Returns

The Miller index.

Return type

(tuple)

get_niggli_reduced_lattice(tol: float = 1e-05) → pymatgen.core.lattice.Lattice[source]

Get the Niggli reduced lattice using the numerically stable algo proposed by R. W. Grosse-Kunstleve, N. K. Sauter, & P. D. Adams, Acta Crystallographica Section A Foundations of Crystallography, 2003, 60(1), 1-6. doi:10.1107/S010876730302186X

Parameters

tol (float) – The numerical tolerance. The default of 1e-5 should result in stable behavior for most cases.

Returns

Niggli-reduced lattice.

get_points_in_sphere(frac_points: List[Union[List[float], numpy.ndarray]], center: Union[List[float], numpy.ndarray], r: float, zip_results=True) → Union[List[Tuple[numpy.ndarray, float, int, numpy.ndarray]], Tuple[List[numpy.ndarray], List[float], List[int], List[numpy.ndarray]]][source]

Find all points within a sphere from the point taking into account periodic boundary conditions. This includes sites in other periodic images.

Algorithm:

  1. place sphere of radius r in crystal and determine minimum supercell (parallelpiped) which would contain a sphere of radius r. for this we need the projection of a_1 on a unit vector perpendicular to a_2 & a_3 (i.e. the unit vector in the direction b_1) to determine how many a_1”s it will take to contain the sphere.

    Nxmax = r * length_of_b_1 / (2 Pi)

  2. keep points falling within r.

Parameters
  • frac_points – All points in the lattice in fractional coordinates.

  • center – Cartesian coordinates of center of sphere.

  • r – radius of sphere.

  • zip_results (bool) – Whether to zip the results together to group by point, or return the raw fcoord, dist, index arrays

Returns

[(fcoord, dist, index, supercell_image) …] since most of the time, subsequent

processing requires the distance, index number of the atom, or index of the image

else:

fcoords, dists, inds, image

Return type

if zip_results

get_vector_along_lattice_directions(cart_coords: Union[List[float], numpy.ndarray])[source]

Returns the coordinates along lattice directions given cartesian coordinates.

Note, this is different than a projection of the cartesian vector along the lattice parameters. It is simply the fractional coordinates multiplied by the lattice vector magnitudes.

For example, this method is helpful when analyzing the dipole moment (in units of electron Angstroms) of a ferroelectric crystal. See the Polarization class in pymatgen.analysis.ferroelectricity.polarization.

Parameters

cart_coords (3x1 array) – Cartesian coords.

Returns

Lattice coordinates.

get_wigner_seitz_cell() → List[List[numpy.ndarray]][source]

Returns the Wigner-Seitz cell for the given lattice.

Returns

A list of list of coordinates. Each element in the list is a “facet” of the boundary of the Wigner Seitz cell. For instance, a list of four coordinates will represent a square facet.

static hexagonal(a: float, c: float)[source]

Convenience constructor for a hexagonal lattice.

Parameters
  • a (float) – a lattice parameter of the hexagonal cell.

  • c (float) – c lattice parameter of the hexagonal cell.

Returns

Hexagonal lattice of dimensions a x a x c.

inv_matrix

Inverse of lattice matrix.

is_hexagonal(hex_angle_tol: float = 5, hex_length_tol: float = 0.01) → bool[source]
is_orthogonal
lengths
lengths_and_angles

Returns (lattice lengths, lattice angles).

lll_inverse
lll_mapping
lll_matrix
matrix

Copy of matrix representing the Lattice

metric_tensor

The metric tensor of the lattice.

static monoclinic(a: float, b: float, c: float, beta: float)[source]

Convenience constructor for a monoclinic lattice.

Parameters
  • a (float) – a lattice parameter of the monoclinc cell.

  • b (float) – b lattice parameter of the monoclinc cell.

  • c (float) – c lattice parameter of the monoclinc cell.

  • beta (float) – beta angle between lattice vectors b and c in degrees.

Returns

Monoclinic lattice of dimensions a x b x c with non right-angle beta between lattice vectors a and c.

norm(coords: Union[List[float], numpy.ndarray], frac_coords: bool = True) → float[source]

Compute the norm of vector(s).

Parameters
  • coords – Array-like object with the coordinates.

  • frac_coords – Boolean stating whether the vector corresponds to fractional or cartesian coordinates.

Returns

one-dimensional numpy array.

static orthorhombic(a: float, b: float, c: float)[source]

Convenience constructor for an orthorhombic lattice.

Parameters
  • a (float) – a lattice parameter of the orthorhombic cell.

  • b (float) – b lattice parameter of the orthorhombic cell.

  • c (float) – c lattice parameter of the orthorhombic cell.

Returns

Orthorhombic lattice of dimensions a x b x c.

reciprocal_lattice

Return the reciprocal lattice. Note that this is the standard reciprocal lattice used for solid state physics with a factor of 2 * pi. If you are looking for the crystallographic reciprocal lattice, use the reciprocal_lattice_crystallographic property. The property is lazily generated for efficiency.

reciprocal_lattice_crystallographic

Returns the crystallographic reciprocal lattice, i.e., no factor of 2 * pi.

static rhombohedral(a: float, alpha: float)[source]

Convenience constructor for a rhombohedral lattice.

Parameters
  • a (float) – a lattice parameter of the rhombohedral cell.

  • alpha (float) – Angle for the rhombohedral lattice in degrees.

Returns

Rhombohedral lattice of dimensions a x a x a.

scale(new_volume: float) → pymatgen.core.lattice.Lattice[source]

Return a new Lattice with volume new_volume by performing a scaling of the lattice vectors so that length proportions and angles are preserved.

Parameters

new_volume – New volume to scale to.

Returns

New lattice with desired volume.

static tetragonal(a: float, c: float)[source]

Convenience constructor for a tetragonal lattice.

Parameters
  • a (float) – a lattice parameter of the tetragonal cell.

  • c (float) – c lattice parameter of the tetragonal cell.

Returns

Tetragonal lattice of dimensions a x a x c.

volume

Volume of the unit cell.

get_integer_index(miller_index: bool, round_dp: int = 4, verbose: bool = True) → Tuple[int, int, int][source]

Attempt to convert a vector of floats to whole numbers.

Parameters
  • miller_index (list of float) – A list miller indexes.

  • round_dp (int, optional) – The number of decimal places to round the miller index to.

  • verbose (bool, optional) – Whether to print warnings.

Returns

The Miller index.

Return type

(tuple)