API reference

Core

class IndexGenerator(*args, **kwargs)[source]

Bases: Protocol

matching_sites(structure, reference_structure, species=None)[source]

Returns a subset of sites from structure (as a list) where each site is the closest to one site in the reference structure.

Parameters:
  • structure (Structure) – The structure being analysed.

  • reference_structure (Structure) – A Structure object containing a set of reference sites.

  • species (list(str), optional) – Optional list of strings specifying a subset of species to return. Default is None, which specifies all species are returned.

Returns:

A list of length-2 tuples for each matching site. Each tuple

contains the corresponding pymatgen Site, and the site index.

Return type:

(list[tuple(Site,int)])

matching_site_indices(structure, reference_structure, species=None)[source]

Returns a subset of site indices from structure (as a list) where each site is the closest to one site in the reference structure.

Parameters:
  • structure (Structure) – The structure being analysed.

  • reference_structure (Structure) – A Structure object containing a set of reference sites.

  • species (list[str], optional) – A list of species labels. If this is set, only matching sites will be included in the returned set.

Returns:

(list[int])

Return type:

list[int]

create_matching_site_generator(reference_structure, species=None)[source]
Parameters:
  • reference_structure (Structure)

  • species (list[str], optional)

Return type:

(func)

generator_from_atom_argument(arg)[source]

Returns a generator function for selecting a subset of sites from a pymatgen Structure object.

Parameters:

arg (various) – Argument used to construct the generator function.

Returns:

Generator function that takes one argument (Structure) and

returns an appropriate list of site indices.

Return type:

(func)

polyhedra_from_distance_cutoff(central_atoms, vertex_atoms, cutoff, label=None)[source]

Construct polyhedra by including all vertex atoms within a cutoff distance.

Parameters:
  • central_atoms (list[Atom]) – Atoms to use as polyhedron centres.

  • vertex_atoms (list[Atom]) – Candidate vertex atoms.

  • cutoff (float) – Maximum distance for a vertex atom to be included.

  • label (str | None) – Optional label for the resulting polyhedra.

Returns:

A list of CoordinationPolyhedron objects.

Return type:

list[CoordinationPolyhedron]

polyhedra_from_nearest_neighbours(central_atoms, vertex_atoms, nn, label=None)[source]

Construct polyhedra from the n nearest vertex atoms to each centre.

Parameters:
  • central_atoms (list[Atom]) – Atoms to use as polyhedron centres.

  • vertex_atoms (list[Atom]) – Candidate vertex atoms.

  • nn (int) – Number of nearest neighbours to include.

  • label (str | None) – Optional label for the resulting polyhedra.

Returns:

A list of CoordinationPolyhedron objects.

Return type:

list[CoordinationPolyhedron]

polyhedra_from_closest_centre(central_atoms, vertex_atoms, label=None)[source]

Construct polyhedra by assigning each vertex atom to its closest centre.

Each vertex atom is assigned to the central atom it is nearest to, so every vertex belongs to exactly one polyhedron.

Parameters:
  • central_atoms (list[Atom]) – Atoms to use as polyhedron centres.

  • vertex_atoms (list[Atom]) – Candidate vertex atoms.

  • label (str | None) – Optional label for the resulting polyhedra.

Returns:

A list of CoordinationPolyhedron objects.

Return type:

list[CoordinationPolyhedron]

polyhedra_from_atom_indices(central_atoms, vertex_atoms, central_indices, vertex_indices, label=None)[source]

Construct a set of polyhedra from lists of atom indices for central and vertex atoms. c

Args:

central_atoms (list(Atom)): List of Atom objects describing the set of possible centre atoms. vertex_atoms (list(Atom)): List of Atom objects describing the set of possible vertex atoms. central_indices (list(int)): List of integer indices specifying each central atom. vertex_indices (list(list(int)): Nested list of integer indices for the vertex atoms in each polyhedron.

Returns:

list(CoordinationPolyhedron)

Raises:

ValueError: If the lengths of the central_indices and vertex_indices lists are unequal.

Parameters:
Return type:

list[CoordinationPolyhedron]

class CoordinationPolyhedron(central_atom, vertices, label=None)[source]

Bases: object

A coordination polyhedron defined by a central atom and its vertex atoms.

Provides geometric analysis including bond distances, angles, volume, edge connectivity, continuous symmetry measures, and neighbour-sharing relationships with other polyhedra.

Initialise a CoordinationPolyhedron object.

Parameters:
  • central_atom (Atom) – the central atom.

  • vertices (list(Atom)) – A list of atoms that define the coordination environment.

  • label (str, optional) – An optional string used to label this coordination polyhedron. if the label is not defined, the label of the central atom will be used.

Returns:

None

property label: str

The label for this polyhedron.

set_label(label)[source]

Set the label for this polyhedron.

Parameters:

label (str) – The new label string.

Return type:

None

update_vertex_neighbours()[source]

Update the neighbour lists on each vertex atom from the edge graph.

Return type:

None

intersection(other_polyhedron)[source]

Returns a tuple of atom indices for vertex atoms shared with another polyhedron.

Parameters:

other_polyhedron (CoordinationPolyhedron) – The other coordination polyhedron.

Returns:

Tuple of shared vertex indices.

Return type:

tuple(int)

property vertex_indices: list[int]

Global indices of the vertex atoms.

vertex_vectors(reference='centroid')[source]

Returns an array of polyhedron centre -> vertex vectors. The polyhedron centre can be set as either the centroid of the polyhedron (default) or the position of the central atom.

Parameters:

reference (str, optional) – The reference point used to compute the vertex vectors. Can either be ‘centroid’ (default) or ‘central_atom’.

Returns:

Array of vectors from the reference point to each vertex.

Return type:

np.ndarray

Raises:

ValueError – If an invalid reference point setting is passed.

property vertex_coords: ndarray

Cartesian coordinates of the vertex atoms as an Nx3 array.

property vertex_labels: list[str | None]

Species labels of the vertex atoms.

property coordination_number: int

Number of vertex atoms in this polyhedron.

property index: int

Global index of the central atom.

property edge_graph: dict[int, list[int]]

Edge connectivity graph for this polyhedron.

A dictionary mapping each vertex index to a list of indices of vertices connected to it by an edge. Computed lazily from the convex hull on first access.

edge_vertex_indices()[source]

Return all edges as sorted pairs of vertex indices.

Return type:

tuple[tuple[int, int], …]

property symmetry_measure: dict[str, float]

Continuous symmetry measures against all reference geometries for this coordination number.

Returns a dictionary mapping geometry names to CSM values. Computed lazily on first access.

Raises:

ValueError – If no reference geometries are defined for this coordination number.

property best_fit_geometry: dict[str, str | float]

The reference geometry with the lowest symmetry measure.

Returns a dictionary with keys 'geometry' (the name) and 'symmetry_measure' (the CSM value).

minimum_image_vertex_coordinates()[source]

Vertex coordinates under the minimum-image convention.

Returns Cartesian coordinates for each vertex, shifted to be the closest periodic image to the central atom.

Returns:

An Nx3 array of minimum-image vertex coordinates.

Return type:

np.ndarray

faces()[source]

Returns a nested set of tuples for each face of the polyhedron. Each per-face tuple returned contains the indices of the vertices that define that face, in sorted numerical order.

Parameters:

None

Returns:

tuple[tuple[int]]

Return type:

tuple[tuple[int, …], …]

Examples

>>> print(polyhedron)
Coordination Polyhedron 4c
255 [12.71362322 17.90999634 12.74490767] S
----------
31 [12.46919306 20.2317206  12.2641591 ] Li
55 [13.0016308  17.39863735 10.46318072] Li
71 [10.4034848  18.18407515 12.43873978] Li
103 [12.17924193 15.66932958 13.34077502] Li
159 [13.24242002 18.43469275 15.02193658] Li
175 [15.02830461 17.60091516 12.52079631] Li
>>> polyhedron.faces()
((31, 159, 175), (31, 55, 71), (31, 55, 175), (31, 71, 159)
 (55, 71, 103), (55, 103, 175), (71, 103, 159), (103, 159, 175))
convex_hull()[source]

Compute the convex hull of the minimum-image vertex coordinates.

Return type:

ConvexHull

construct_edge_graph()[source]

Build the edge connectivity graph from the convex hull.

For polyhedra with more than three vertices, edges are derived from the convex hull simplices, merging coplanar faces. For three or fewer vertices, all vertices are connected.

Returns:

A dictionary mapping each vertex index to the indices of its edge-connected neighbours.

Return type:

dict[int, list[int]]

property vertex_count: Counter[str | None]

Count of vertex atoms grouped by species label.

vertex_distances(reference='central_atom')[source]

Returns an array of distances from either the central atom or the centroid to the vertex atoms.

Parameters:

reference (str, optional) – The reference point for distance calculations. Can be either ‘central_atom’ (default) or ‘centroid’.

Returns:

An array of distances between each vertex and the reference point.

Return type:

np.ndarray[float, …]

Raises:

ValueError – If an invalid reference point is provided.

vertex_distances_and_labels(reference='central_atom')[source]

Returns a tuple of distances and species labels from the reference point to the vertex atoms.

Parameters:

reference (str, optional) – The reference point for distance calculations. Can be either ‘central_atom’ (default) or ‘centroid’.

Returns:

A tuple of length-2 tuples, containing for each vertex the

distance from the reference point and the species label of the vertex atom.

Return type:

tuple[tuple[float, str | None], …]

Raises:

ValueError – If an invalid reference point is provided.

equal_vertices(other)[source]

Test whether this CoordinationPolyhedron has vertices with the same labels as another CoordinationPolyhedron.

Parameters:

other (CoordinationPolyhedron) – The other CoordinationPolyhedron.

Returns:

True / False.

Return type:

bool

equal_edge_graph(other)[source]

Test whether this CoordinationPolyhedron has the same edge graph as another CoordinationPolyhedron.

Parameters:

other (CoordinationPolyhedron) – The other CoordinationPolyhedron.

Returns:

True or False.

Return type:

bool

equal_members(other)[source]

Test whether this CoordinationPolyhedron has the same member atoms as another CoordinationPolyhedron.

Parameters:

other (CoordinationPolyhedron) – The other CoordinationPolyhedron.

Returns:

True or False.

Return type:

bool

neighbours()[source]

Returns a tuple of neighbouring polyhedra. Two polyhedra are considered to be neighbours if they have one or more vertex atoms in common.

Parameters:

none

Returns:

Tuple of neighbouring polyhedra.

Return type:

tuple(CoordinationPolyhedron)

neighbours_by_index_and_shared_vertices()[source]

Returns a dictionary of vertex atoms shared with neighbouring polyhedra. Two polyhedra are considered to be neighbours if they have one of more vertex atoms in common.

Parameters:

None

Returns:

Dictionary where the keys are polyhedra

indices, and the values are the indices of common vertex atoms.

Return type:

dict(int, tuple(int, …))

corner_sharing_neighbour_list()[source]

Indices of neighbouring polyhedra that share exactly one vertex (corner-sharing).

Return type:

tuple[int, …]

edge_sharing_neighbour_list()[source]

Indices of neighbouring polyhedra that share exactly two vertices (edge-sharing).

Return type:

tuple[int, …]

face_sharing_neighbour_list()[source]

Indices of neighbouring polyhedra that share three or more vertices (face-sharing).

Return type:

tuple[int, …]

vertex_vector_projections(vectors, reference='centroid')[source]

Calculate the projection of each vertex vector on one or more input vectors.

Parameters:
  • vectors (np.ndarray) – A Nx3 numpy array, where each row is a vector used to calculate the projection. vectors can also be a single length 3 array.

  • reference (str, optional) – The reference point for vertex vectors. Can be either ‘centroid’ (default) or ‘central_atom’.

Returns:

A (N_vertex x N_vector) dimension numpy array.

Return type:

np.ndarray

coordination_distances()[source]

Returns a list of distances from the central atom to the vertex atoms.

This method is a wrapper around vertex_distances() that maintains backward compatibility and returns a list instead of a tuple.

Parameters:

None

Returns:

A list of distances between each vertex and the central atom.

Return type:

list[float]

angles()[source]

List of all vertex–centroid–vertex angles.

Parameters:

None

Returns:

List of angles.

Return type:

list

vertex_vector_orientations(units='degrees', return_distance=False, reference='centroid')[source]

Returns the angular orientations of each vertex vector.

The orientation is defined by two angles, theta and phi. Theta is the angle with respect to [0, 0, 1] and ranges from 0 to 180 degrees. Phi is the angle with respect to [1, 0, 0] and ranges from -180 to +180 degrees.

Parameters:
  • units (str, optional) – Optionally select the units for the calculated angles. Options are degrees or radians. Default is degrees.

  • return_distance (bool, optional) – Optionally also return the distance.

  • reference (str, optional) – The reference point for vertex vectors. Options are centroid or central_atom. Default is centroid.

Returns:

A list of (theta,phi) tuple pairs, or (theta,phi,distance) if

return_distance is True.

Return type:

list(tuple)

Raises:

ValueError – If an invalid reference point is provided.

property volume: float

Volume of this polyhedron.

Parameters:

None

Returns:

The volume.

Return type:

float

classmethod from_sites(central_site, vertex_sites, label=None)[source]

Create a CoordinationPolyhedron from a set of pymatgen PeriodicSite objects.

Parameters:
  • central_site (pymatgen.PeriodicSite) – A pymatgen PeriodicSite object describing an atom at the nominal centre of the polyhedron.

  • vertex_sites (list[pymatgen.PeriodicSite) – A list of pymatgen PeriodicSite objects describing the atoms at the vertices.

  • label (str, optional) – An optional string used to label this coordination polyhedron.

Returns:

The CoordinationPolyhedron object.

Return type:

CoordinationPolyhedron

vertices_by_indices(vertex_indices)[source]

Select a subset of vertices from this polyhedron with a list of vertex indices.

Parameters:

vertex_indices (list) – List of vertex indices (int).

Returns:

A list of Atom objects containing the matching vertices.

Return type:

list

vertex_internal_index_from_global_index(vertex_global_index)[source]

Returns the internal index for the vertex with a given global index.

Parameters:

vertex_global_index (int) – The global index for a given vertex.

Returns:

the internal index for the corresponding vertex.

Return type:

int

Raises:

ValueError – If this polyhedron does not have a vertex with a global index equal to vertex_global_index.

centroid()[source]

Calculate the centroid of the vertices, accounting for periodic boundary conditions.

Returns:

The 3D coordinates of the centroid.

Return type:

np.ndarray

centroid_to_central_atom_vector()[source]

Calculate the vector displacement from the centroid of the vertices to the central atom, accounting for periodic boundary conditions.

Returns:

A 3D vector representing the displacement from the vertex centroid to the central atom.

Return type:

np.ndarray

property off_centre_displacement: float

Returns the displacment from the centroid of the polyhedron vertices to the central atom, accounting for periodic boundary conditions.

Returns:

float

radial_distortion_parameter(reference='centroid', normalize=True, method='MSD')[source]

Calculate the radial distortion parameter for the coordination polyhedron.

The radial distortion parameter is defined as:

  • For MSD: Delta = mean(((d_i - d_mean) / (d_mean if normalize else 1))^2)

  • For MAD: Delta = mean(abs(d_i - d_mean) / (d_mean if normalize else 1))

where d_i are the N distances from the reference point to each vertex, and d_mean is their mean value.

Parameters:
  • reference (Literal['central_atom', 'centroid'], optional) – The reference point for calculating distances. Can be either ‘centroid’ (default) or ‘central_atom’.

  • normalize (bool, optional) – Whether to normalize the distortion parameter with respect to the mean vertex distance. Default is True.

  • method (Literal['MSD', 'MAD'], optional) – The method to use for calculating distortion. Can be either ‘MSD’ (Mean Squared Deviation, default) or ‘MAD’ (Mean Absolute Deviation).

Returns:

The radial distortion parameter.

Return type:

float

Raises:

ValueError – If an invalid reference point or method is provided.

vertex_angles(vertex_pairs, reference='central_atom')[source]

Calculate angles between specified pairs of vertices.

Parameters:
  • vertex_pairs (tuple[tuple[int, int], ...]) – Pairs of vertex global indices to calculate angles for.

  • reference (Literal['central_atom', 'centroid'], optional) – The reference point for vertex vectors. Defaults to ‘central_atom’.

Returns:

Array of calculated angles in degrees.

Return type:

np.ndarray

Raises:

ValueError – If an invalid vertex index is provided.

class Atom(index, site, label=None)[source]

Bases: object

Atom class

Initialise an Atom object.

Parameters:
  • index (int) – Numerical index identifying this atom.

  • site (pymatgen.Site) – The pymatgen Site (or PeriodicSite) object describing this atom.

  • label (str, optional) – An optional string labelling this atom. If no label is given the site species string will be used.

in_polyhedra

List of polyhedra that this atom is part of.

Type:

list

neighbours

(dict(int: list)): Dictionary of lists of neighbours for each polyhedron this atom is part of. Dictionary keys are the indexes of each respective polyhedron.

as_dict()[source]

json-serializable dict representation of Atom.

Return type:

dict[str, Any]

to(fmt=None, filename=None)[source]

Outputs the structure to a file or string.

Parameters:
  • fmt (str, optional) – Format to ouput to. Defaults to JSON unlees the filename is provided. If fmt is specified this overrides the filename extension in the filename.

  • filename (str, optional) – If provided, the output will be written to a file. If fmt is not specified, the format will be determined from the filename extension.

Returns:

(str) if filename is None. None otherwise.

Return type:

None | str

property frac_coords: ndarray

Fractional coordinates

property coords: ndarray

Cartesian coordinates

distance(other)[source]

Shortest distance using periodic boundary conditions to another Atom.

Parameters:

other (Atom) – The other atom.

Returns:

float

Return type:

float

Symmetry measures

class SymmetryMeasure(reference_points, name)[source]

Bases: object

Continuous symmetry measure calculator for a reference geometry.

Wraps a set of ideal reference points and computes the minimum continuous symmetry measure (CSM) over all symmetry-inequivalent vertex permutations. The permutations are computed lazily on first use via bsym.

Parameters:
  • reference_points (np.ndarray) – An Nx3 array of ideal vertex coordinates for the reference geometry.

  • name (str) – A human-readable name for this geometry (e.g. 'Octahedron').

property permutations: ndarray

Symmetry-inequivalent index permutations for this geometry.

Computed lazily on first access from the reference points using bsym.

minimum_symmetry_measure(distorted_points)[source]

Compute the minimum CSM of distorted points against this reference.

Iterates over all symmetry-inequivalent permutations and returns the smallest continuous symmetry measure.

Parameters:

distorted_points (ndarray) – An Nx3 array of vertex coordinates to compare against the reference geometry.

Returns:

The minimum symmetry measure value. A value of 0 indicates perfect agreement with the reference geometry.

Return type:

float

classmethod from_name(name)[source]

Construct a SymmetryMeasure from a named reference geometry.

Parameters:

name (str) – The name of the reference geometry (e.g. 'Octahedron', 'Tetrahedron'). See get_reference_geometry() for available names.

Returns:

A new SymmetryMeasure instance.

Return type:

SymmetryMeasure

Continuous symmetry measure (CSM) via SVD-based Procrustes analysis.

References

Pinsky & Avnir, Inorganic Chemistry 37, 5575 (1998).

class SymmetryMeasureResult(symmetry_measure, rotation_matrix, scaling_factor)[source]

Bases: NamedTuple

Result of a continuous symmetry measure calculation.

Create new instance of SymmetryMeasureResult(symmetry_measure, rotation_matrix, scaling_factor)

Parameters:
symmetry_measure: float

The CSM value as a percentage (0–100).

rotation_matrix: ndarray

The optimal 3x3 rotation matrix.

scaling_factor: float

The optimal scaling factor.

continuous_symmetry_measure(points_distorted, points_perfect)[source]

Compute the continuous symmetry measure (CSM) of a set of points.

Uses SVD-based orthogonal Procrustes analysis to find the optimal rotation and scaling that align the distorted points to the perfect reference, then computes the normalised residual as a percentage.

Parameters:
  • points_distorted (ndarray) – An Nx3 array of distorted vertex coordinates.

  • points_perfect (ndarray) – An Nx3 array of ideal vertex coordinates.

Returns:

A SymmetryMeasureResult with fields symmetry_measure (float, percentage 0–100), rotation_matrix (3x3 ndarray), and scaling_factor (float).

Return type:

SymmetryMeasureResult

Ideal vertex coordinates for reference coordination geometries.

Coordinates are centred on the origin and at unit distance, matching the conventions used by pymatgen’s chemenv module from which they were originally extracted.

get_reference_geometry(name)[source]

Return the ideal vertex coordinates for a named reference geometry.

Parameters:

name (str) – The geometry name, e.g. 'Octahedron', 'Tetrahedron'.

Returns:

An Nx3 array of ideal vertex coordinates centred on the origin.

Raises:

KeyError – If the name is not a recognised geometry.

Return type:

ndarray

Trajectories

class PolyhedronTrajectory(polyhedra, config_numbers=None)[source]

Bases: object

Collects a sequence of polyhedra from different configurations into a trajectory.

Initialise a PolyhedronTrajectory object.

Parameters:
  • polyhedra (list(CoordinationPolyhedron)) – A list of CoordinationPolyhedron objects.

  • ( (config_numbers) – obj:list(int), optional): An optional list of configuration numbers for each of the polyhedra.

  • config_numbers (list[int] | None)

Returns:

None

Analysis

check_octahedra(polyhedron)[source]

Check whether a polyhedron is considered an octahedron.

Parameters:

polyhedron (CoordinationPolyhedron) – The CoordinationPolyhedron to be tested.

Returns:

None

Raises:

ValueError – If the “best fit” geometry (lowest continuous symmetry measure) is not an octahedron, a ValueError is raised.

Return type:

None

adjacent_vertex_pairs(polyhedron, check=True)[source]

Find pairs of adjacent vertices in an octahedral polyhedron.

Parameters:
  • polyhedron (CoordinationPolyhedron) – The polyhedron to be analysed.

  • check (bool, optional) – Whether to check if the polyhedron is octahedral. Default is True.

Returns:

A list of tuples, each containing a pair of adjacent Atom objects.

Return type:

list[tuple[Atom, Atom]]

Raises:

ValueError – If the polyhedron is not recognized as an octahedron.

opposite_vertex_pairs(polyhedron, check=True)[source]

For an octahedral polyhedron, find the pairs of vertices opposite each other.

Parameters:
  • polyhedron (CoordinationPolyhedron) – The polyhedron to be analysed.

  • check (bool) – (Optional, bool): Optional flag to set whether to check that this polyhedron is an octahedron. Default is True.

Returns:

3 pairs of vertex atoms.

Return type:

(tuple)

opposite_vertex_distances(polyhedron, check=True)[source]

For an octahedral polyhedron, return the distances between pairs of opposite (trans) vertices.

Parameters:
  • polyhedron (CoordinationPolyhedron) – The polyhedron to be analysed.

  • check (bool) – (Optional, bool): Optional flag to set whether to check that this polyhedron is an octahedron. Default is True.

Returns:

a length-3 tuple of distances.

Return type:

(tuple)

trans_vertex_vectors(polyhedron, check=True)[source]

For an octahedral polyhedron, return the vectors between pairs of trans vertices.

Parameters:
  • polyhedron (CoordinationPolyhedron) – The polyhedron to be analysed.

  • check (Optional, bool) – Optional flag to set whether to check that this polyhedron is an octahedron. Default is True.

Returns:

A list of three numpy arrays, each representing a vector

between a pair of trans vertices.

Return type:

list[np.ndarray]

Raises:

ValueError – If the polyhedron is not recognized as an octahedron.

isomer_is_trans(polyhedron, check=True)[source]

For an octahedral polyhedron with 2+4 coordination, return True for trans coordination, and False for cis coordination.

Parameters:
  • polyhedron (CoordinationPolyhedron) – The polyhedron to be analysed.

  • check (bool) – (Optional, bool): Optional flag to set whether to check that this polyhedron is an octahedron. Default is True.

Returns:

True for trans, False for cis.

Return type:

(bool)

isomer_is_cis(polyhedron, check=True)[source]

For an octahedral polyhedron with 2+4 coordination, return True for cis coordination, and False for trans coordination.

Parameters:
  • polyhedron (CoordinationPolyhedron) – The polyhedron to be analysed.

  • check (bool) – (Optional, bool): Optional flag to set whether to check that this polyhedron is an octahedron. Default is True.

Returns:

True for cis, False for trans.

Return type:

(bool)

isomer_is_fac(polyhedron, check=True)[source]

For an octahedral polyhedron with 3+3 coordination, return True for fac coordination, and False for mer coordination.

Parameters:
  • polyhedron (CoordinationPolyhedron) – The polyhedron to be analysed.

  • check (bool) – (Optional, bool): Optional flag to set whether to check that this polyhedron is an octahedron. Default is True.

Returns:

True for fac, False for mer.

Return type:

(bool)

isomer_is_mer(polyhedron, check=True)[source]

For an octahedral polyhedron with 3+3 coordination, return True for mer coordination, and False for fac coordination.

Parameters:
  • polyhedron (CoordinationPolyhedron) – The polyhedron to be analysed.

  • check (bool) – (Optional, bool): Optional flag to set whether to check that this polyhedron is an octahedron. Default is True.

Returns:

True for fac, False for mer.

Return type:

(bool)

trans_vector_orthogonality(polyhedron, check=True)[source]

For each trans pair vector, calculate the smallest angle between the vector and the plane defined by the other two vectors.

Parameters:
  • polyhedron (CoordinationPolyhedron) – The polyhedron to be analysed.

  • check (bool, optional) – Optional flag to set whether to check that this polyhedron is an octahedron. Default is True.

Returns:

A tuple containing the three smallest angles (in degrees) between each vector and the plane defined by the other two vectors.

Return type:

tuple[float, float, float]

Raises:

ValueError – If the polyhedron is not recognized as an octahedron.

cos_theta(a, b)[source]

Calculate the cosine of the angle between two vectors.

Parameters:
  • a (np.ndarray) – Vector a.

  • b (np.ndarray) – Vector b.

Returns:

float

Return type:

float

projection_xyz(vec_in)[source]

Maximum projection score of a vector onto the Cartesian axes.

For each Cartesian axis, computes 2 * cos^2(theta) - 1 where theta is the angle between vec_in and that axis. The per-axis score ranges from 1.0 (parallel) to -1.0 (perpendicular). The function returns the maximum score across the three axes, so the overall range is 1.0 (aligned with an axis) to -1/3 (equally spaced between all three axes, e.g. along [1, 1, 1]).

Parameters:

vec_in (ndarray) – A length-3 vector.

Returns:

The maximum projection score across the three Cartesian axes.

Return type:

float

oct_rotational_order_parameter(points)[source]

Order parameter for rotational orientation of an octahedron with respect to <100> vectors.

Parameters:

points (ndarray) – A 6x3 array of centroid-centred vertex vectors.

Returns:

1.0 for a perfectly aligned octahedron. 0.33 for a 45 degree rotated octahedron around one axis.

Return type:

float

class OrientationDict[source]

Bases: TypedDict

class RotationAnalyser(reference_points)[source]

Bases: object

Class for analysing rotational orientation of polyhedra.

All reference orientations must be rotated copies of the same geometry (i.e. share the same point group). The symmetry analysis is performed once from the first reference orientation and reused for all others.

reference_points

An MxNx3 numpy array of points around the origin ([0.0, 0.0, 0.0]) that define sets of centre-vertex vectors for specific reference orientations, e.g. for a single orientation of a tetrahedron, reference points will be a size 1x4x3 array, e.g.:

[[[ 1.0, -1.0,  1.0],
  [-1.0, -1.0, -1.0],
  [ 1.0,  1.0, -1.0],
  [-1.0,  1.0,  1.0]]]
Type:

np.ndarray

Initialise a RotationAnalyser instance.

Parameters:

reference_points (np.ndarray) –

Either an Nx3 or MxNx3 numpy array of points around the origin ([0.0, 0.0, 0.0]) that define sets of centre-vertex vectors, and are used to classify the rotational orientation of polyhedra.

If an Nx3 array is passed in, this will be converted to a 1xNx3 array.

property reduced_permutations: ndarray

Symmetry-inequivalent permutation indices for the reference geometry.

Computed lazily on first access.

property proper_rotations: ndarray

Proper rotation matrices (Kx3x3) of the reference geometry’s point group.

Computed lazily on first access.

discrete_orientation(points)[source]

Find the discrete “closest orientation” for an input polyhedron of points.

For example, a tetrahedron has 12 proper rotation symmetry operations. A distorted tetrahedron with vertices approximately aligned along the unordered vectors:

[[ 1.0, -1.0,  1.0],
 [-1.0, -1.0, -1.0],
 [ 1.0,  1.0, -1.0],
 [-1.0,  1.0,  1.0]]

can be in one of 12 orientations. This method compares the input points to the reference orientations and returns the orientation that minimises the rotational distance.

The algorithm uses a two-phase approach to avoid the N! cost of evaluating all vertex permutations:

  1. Find the best vertex alignment using symmetry-reduced permutations of the reference points (N! / |G| CSM evaluations instead of N!).

  2. Compose the best-fit Procrustes rotation R_best with all proper rotation matrices S_i of the reference geometry’s point group to give the rotation matrix for each orientation: R_i = S_i @ R_best.

  3. Calculate rotational distances from the composed matrices.

  4. Return the orientation with the minimum rotational distance.

The composed rotation matrices in step 2 are exact for undistorted geometries. For distorted geometries they are approximate, since the true Procrustes rotation for each symmetry-equivalent permutation differs slightly from the composed matrix. In practice the results are numerically indistinguishable from the brute-force approach for moderate distortions.

If exact rotation matrices are needed for highly distorted geometries, an alternative is a two-phase orbit-expansion approach: use reduced permutations to identify the minimum-CSM orbit (phase 1), then expand that orbit to all |G| equivalent permutations and compute CSM for each to obtain exact rotation matrices (phase 2). This costs N! / |G| + |G| CSM evaluations per reference orientation rather than N! / |G|, but avoids the composition approximation.

Parameters:

points (ndarray) – Nx3 numpy array of points describing the coordinates of the input polyhedron, centred around zero.

Returns:

  • orientation_index (int): Index of this particular orientation.

  • reference_geometry_index (int): Index of the reference geometry the closest discrete orientation is equivalent to.

  • rotational_distance (float): Angle of rotation from the relevant reference orientation, in radians.

  • all_rotational_distances (np.ndarray): Array of angles of rotation from each reference orientation.

  • symmetry_measure (float): The continuous symmetry measure (CSM) for this polyhedron.

Return type:

Dictionary describing the orientation, with keys

polyhedron_orientation(polyhedron)[source]

Find the discrete orientation of a coordination polyhedron.

Convenience wrapper around discrete_orientation() that extracts vertex vectors from the polyhedron relative to its central atom.

Parameters:

polyhedron (CoordinationPolyhedron) – The coordination polyhedron to analyse.

Returns:

Dictionary describing the orientation. See discrete_orientation() for details.

Return type:

OrientationDict

Plotting

plot_orientation_distribution(orientations, title=None, figsize=(10, 8), cmap=<matplotlib.colors.ListedColormap object>, fontsize=None, plot=False)[source]

Create a contour map of the probability distribution for (phi, theta) orientations.

Parameters:
  • orientations (list[tuple[float, float]]) – List of (theta, phi) tuples.

  • title (str | None, optional) – Title for the plot. Defaults to None.

  • figsize (tuple[int, int], optional) – Figure size as (width, height). Defaults to (10, 8).

  • cmap (matplotlib.colors.Colormap, optional) – Colormap to use for the contour plot. Defaults to cm.lipari.

  • fontsize (int | None, optional) – Font size for labels and title. If None, uses matplotlib defaults.

  • plot (bool, optional) – Whether to display the plot immediately. Defaults to False.

Returns:

The generated figure object.

Return type:

matplotlib.figure.Figure

Raises:

ValueError – If orientations is empty or contains invalid data.

Utilities

class T

Utility functions

alias of TypeVar(‘T’)

flatten(this_list)[source]

Flattens a nested list.

Parameters:
Returns:

The flattened list.

Return type:

(list)

lattice_mc_string(polyhedron, neighbour_list)[source]

Returns a string representation of a polyhedron as a lattice_mc site-input formatted site.

Parameters:
Returns:

str

Return type:

str

Example

>>> nlist = {1: (3, 5, 8), ...}
>>> lattice_mc_string(polyhedron, neighbour_list=nlist
site: 1
centre: 2.94651000 1.70116834 1.20290767
neighbours: 3 5 8
label: oct
prune_neighbour_list(neighbours, indices)[source]

Prune a neighbour list to include only the specified indices.

Filters both the keys and the values of the neighbour dictionary, keeping only entries whose keys appear in indices and removing any neighbour references to indices not in the set.

Parameters:
  • neighbours (dict[int, tuple[int, ...]]) – A neighbour list mapping each polyhedron index to a tuple of its neighbour indices.

  • indices (list[int]) – The subset of indices to retain.

Returns:

The pruned neighbour list.

Return type:

dict[int, tuple[int, …]]