crystal_torture.pymatgen_interface

Functions for setting up a node, cluster and graph using pymatgen.

crystal_torture.pymatgen_interface.clusters_from_file(filename: str, rcut: float, elements: set[str]) set[Cluster][source]
crystal_torture.pymatgen_interface.clusters_from_structure(structure: Structure, rcut: float, elements: set[str]) set[Cluster][source]

Take a pymatgen structure and convert it to a graph object.

Parameters:
  • structure – Pymatgen structure object to set up graph from.

  • rcut – Cut-off radii for node-node connections in forming clusters.

  • elements – Set of element strings to include in setting up graph.

Returns:

Set of clusters.

crystal_torture.pymatgen_interface.create_halo(structure: Structure, neighbours: list[list[tuple]]) tuple[Structure, list[list[int]]][source]

Take a pymatgen structure object and set up a halo by making a 3x3x3 supercell.

Parameters:
  • structure – Pymatgen Structure object.

  • neighbours – List of neighbours for sites in structure (from get_all_neighbors_and_image).

Returns:

  • structure: 3x3x3 supercell pymatgen Structure object.

  • neighbours: New list of neighbours for sites in supercell structure.

Return type:

Tuple containing

crystal_torture.pymatgen_interface.filter_structure_by_species(structure: Structure, species_list: list[str]) Structure[source]

Filter structure to keep only specified species.

Creates a new Structure containing only the sites with species that are in the provided species list. The original structure is not modified.

Parameters:
  • structure – Pymatgen Structure object to filter.

  • species_list – List of element symbols to keep in the filtered structure.

Returns:

New Structure object containing only sites with species in species_list.

Raises:
  • ValueError – If species_list is empty.

  • ValueError – If species_list contains elements not present in the structure.

Example

>>> structure = Structure(lattice, ["Li", "Mg", "O"], coords)
>>> filtered_structure = filter_structure_by_species(structure, ["Li", "O"])
>>> filtered_structure.symbol_set  # {"Li", "O"}
crystal_torture.pymatgen_interface.get_all_neighbors_and_image(structure: Structure, r: float, include_index: bool = False) list[list[tuple]][source]

Get neighbours for each atom in the unit cell, out to a distance r.

Modified from pymatgen to return image (used for mapping to supercell), and to use the f2py wrapped OpenMP dist subroutine to get the distances (smaller memory footprint and faster than numpy).

Returns a list of list of neighbors for each site in structure. Use this method if you are planning on looping over all sites in the crystal. If you only want neighbors for a particular site, use the method get_neighbors as it may not have to build such a large supercell However if you are looping over all sites in the crystal, this method is more efficient since it only performs one pass over a large enough supercell to contain all possible atoms out to a distance r.

Parameters:
  • structure – Pymatgen Structure object.

  • r – Radius of sphere.

  • include_index – Whether to include the non-supercell site in the returned data.

Returns:

A list of a list of nearest neighbors for each site, i.e., [[(site, dist, index, image) …], ..]. Index only supplied if include_index = True. The index is the index of the site in the original (non-supercell) structure. This is needed for ewaldmatrix by keeping track of which sites contribute to the ewald sum.

crystal_torture.pymatgen_interface.graph_from_file(filename: str, rcut: float, elements: set[str]) Graph[source]
crystal_torture.pymatgen_interface.graph_from_structure(structure: Structure, rcut: float, elements: set[str]) Graph[source]

Create a graph from a pymatgen structure.

Creates a Graph object containing clusters of connected nodes and a filtered structure. The clusters are formed by connecting sites within the cutoff radius, and only sites with the specified elements are included.

Parameters:
  • structure – Pymatgen Structure object to create graph from.

  • rcut – Cutoff radius for node-node connections in forming clusters.

  • elements – Set of element strings to include in the graph.

Returns:

Graph object containing clusters and filtered structure.

Example

>>> structure = Structure(lattice, ["Li", "Mg", "O"], coords)
>>> graph = graph_from_structure(structure, 3.0, {"Li", "O"})
>>> graph.clusters  # Clusters containing only Li and O sites
crystal_torture.pymatgen_interface.map_index(uc_neighbours: list[list[int]], uc_index: list[int], x_d: int, y_d: int, z_d: int) list[list[int]][source]

Take a list of neighbour indices for sites in the original unit cell and map them onto all supercell sites.

Parameters:
  • uc_neighbours – List of lists containing neighbour indices for the nodes that are in the primitive cell.

  • uc_index – List of indices corresponding to the primitive cell nodes.

  • x_d – X dimension of supercell.

  • y_d – Y dimension of supercell.

  • z_d – Z dimension of supercell.

Returns:

List of neighbour indices for all nodes.

crystal_torture.pymatgen_interface.nodes_from_structure(structure: Structure, rcut: float, get_halo: bool = False) set[Node][source]

Take a pymatgen structure object and convert to Nodes for interrogation.

crystal_torture.pymatgen_interface.set_fort_nodes(nodes: set[Node]) None[source]

Set up a copy of the nodes and the neighbour indices in the tort.f90 Fortran module.

This allows access if using the Fortran tortuosity routines.

Parameters:

nodes – Set of Node objects to set up in Fortran module.

Sets:

tort.tort_mod.nodes: Allocates space to hold node indices for full graph. tort.tort_mod.uc_tort: Allocates space to hold unit cell node tortuosity for full graph.

Raises:

FortranNotAvailableError – If Fortran extensions are not available.