kgcnn.crystal package¶
Subpackages¶
Submodules¶
kgcnn.crystal.base module¶
-
class
kgcnn.crystal.base.
CrystalPreprocessor
(output_graph_as_dict: bool = False, lattice: str = 'graph_lattice', species: str = 'node_number', coords: str = 'node_coordinates', charge: str = 'charge')[source]¶ Bases:
Callable
[[pymatgen.core.structure.Structure
],networkx.classes.multidigraph.MultiDiGraph
]Base class for crystal preprocessors.
Concrete CrystalPreprocessors must be implemented as subclasses.
-
__call__
(structure: Union[pymatgen.core.structure.Structure, kgcnn.graph.base.GraphDict]) → Union[networkx.classes.multidigraph.MultiDiGraph, kgcnn.graph.base.GraphDict][source]¶ Function to process crystal structures. Executes
call
.- Parameters
structure (Structure) – Crystal for which the graph representation should be calculated.
- Raises
NotImplementedError – Should be implemented in a subclass.
- Returns
Graph representation of the crystal.
- Return type
MultiDiGraph
-
call
(structure: pymatgen.core.structure.Structure) → networkx.classes.multidigraph.MultiDiGraph[source]¶ Should be implemented in a subclass.
- Parameters
structure (Structure) – Crystal for which the graph representation should be calculated.
- Raises
NotImplementedError – Should be implemented in a subclass.
- Returns
Graph representation of the crystal.
- Return type
MultiDiGraph
-
edge_attributes
= []¶
-
get_config
() → dict[source]¶ Returns a dictionary uniquely identifying the CrystalPreprocessor and its configuration.
- Returns
A dictionary uniquely identifying the CrystalPreprocessor and its configuration.
- Return type
-
graph_attributes
= []¶
-
hash
() → str[source]¶ Generates a unique hash for the CrystalPreprocessor and its configuration.
- Returns
A unique hash for the CrystalPreprocessor and its configuration.
- Return type
-
node_attributes
= []¶
-
kgcnn.crystal.graph_builder module¶
-
kgcnn.crystal.graph_builder.
_estimate_nn_radius_from_density
(k: int, coordinates: numpy.ndarray, lattice: numpy.ndarray, empirical_tol_factor: float = 0.0)[source]¶ Rough estimate of the expected radius to find N nearest neighbours.
-
kgcnn.crystal.graph_builder.
_get_attr_from_graph
(graph: networkx.classes.multidigraph.MultiDiGraph, attr_name: str, make_copy: bool = False) → Union[Any, numpy.ndarray][source]¶ Utility function to obtain graph-level information of the underlying crystal.
-
kgcnn.crystal.graph_builder.
_get_cube
(dim: int) → numpy.ndarray[source]¶ Generate a cubic mesh.
- Parameters
dim (int) – Dimension for cubic mesh.
- Returns
Cubic mesh.
- Return type
np.ndarray
-
kgcnn.crystal.graph_builder.
_get_geometric_properties_of_unit_cell
(coordinates: numpy.ndarray, lattice: numpy.ndarray)[source]¶ Diameter of a 3D unit cell and other properties.
- Parameters
coordinates (np.ndarray) – Coordinates array.
lattice (np.ndarray) – Lattice matrix.
- Returns
(center_unit_cell, max_diameter_cell, volume_unit_cell, density_unit_cell)
- Return type
-
kgcnn.crystal.graph_builder.
_get_max_diameter
(lattice: numpy.ndarray) → Union[float, numpy.ndarray][source]¶ Determine the max diameter of a lattice.
- Parameters
lattice (np.ndarray) – Lattice matrix.
- Returns
Max diameter of the lattice.
- Return type
np.ndarray
-
kgcnn.crystal.graph_builder.
_get_mesh
(size: Union[int, list, tuple], dim: int) → numpy.ndarray[source]¶ Utility function to create a numpy mesh grid with indices at last dimension.
-
kgcnn.crystal.graph_builder.
_get_super_cell_grid_frac_coords
(lattice: numpy.ndarray, frac_coords: numpy.ndarray, size: Union[int, list, numpy.ndarray])[source]¶ Get frac coordinates for positions in a grid of unit cells that is a cubic super-cell.
..code - block:: python
import numpy as np from kgcnn.crystal.graph_builder import _get_super_cell_grid_frac_coords coordinates = _get_super_cell_grid_frac_coords(
np.array([[1.0, 0.0, 0.0], [0.0, 1.0, 0.5], [0.0, 0.5, 1.5]]), np.array([[0.0, 0.0, 0.0], [0.5, 0.5, 0.5]]), [3, 3, 3]
) print(coordinates.shape) # (7, 7, 7, 2, 3)
- Parameters
lattice (np.ndarray) – Lattice matrix.
frac_coords (np.ndarray) – Fractional coordinates of atoms in unit cell.
size (list) – Size of the super-cell in each dimension.
- Returns
List of fractional coordinates of atoms in the super-cell.
- Return type
np.ndarray
-
kgcnn.crystal.graph_builder.
_to_unit_cell
(frac_coords)[source]¶ Converts fractional coords to be within the \([0,1)\) interval.
- Parameters
frac_coords – Fractional coordinates to map into \([0,1)\) interval.
- Returns
Fractional coordinates within the [0,1) interval.
-
kgcnn.crystal.graph_builder.
add_edge_information
(graph: networkx.classes.multidigraph.MultiDiGraph, inplace=False, frac_offset=False, offset=True, distance=True) → networkx.classes.multidigraph.MultiDiGraph[source]¶ Adds edge information, such as offset ( frac_offset, offset ) and distances ( distance ) to edges.
- Parameters
graph (MultiDiGraph) – Graph for which to add edge information.
inplace (bool, optional) – Whether to add the edge information to the given graph or create a copy with added edges. Defaults to False.
frac_offset (bool, optional) – Whether to add fractional offsets (frac_offset attribute) to edges. Defaults to False.
offset (bool, optional) – Whether to add offsets (offset attribute) to edges. Defaults to True.
distance (bool, optional) – Whether to add distances (distance attribute) to edges. Defaults to True.
- Returns
The graph with added edge information.
- Return type
MultiDiGraph
-
kgcnn.crystal.graph_builder.
add_knn_bonds
(graph: networkx.classes.multidigraph.MultiDiGraph, k: int = 12, max_radius: float = 10.0, tolerance: Optional[float] = None, inplace: bool = False) → networkx.classes.multidigraph.MultiDiGraph[source]¶ Adds kNN-based edges to a unit cell graph.
- Parameters
graph (MultiDiGraph) – The unit cell graph to add kNN-based edges to.
k (int, optional) – How many neighbors to add for each node. Defaults to 12.
max_radius (float, optional) – This parameter has no effect on the outcome of the graph. It may only on the runtime. The algorithm starts the kNN search in the environment the radius of max_radius. If the kth neighbor is not within this radius the algorithm is called again with twice the initial radius. Defaults to 10.
tolerance (Optional[float], optional) – If tolerance is not None, edges with distances of the k-th nearest neighbor plus the tolerance value are included in the graph. Defaults to None.
inplace (bool, optional) – Whether to add the edges to the given graph or create a copy with added edges. Defaults to False.
- Returns
Graph with added edges.
- Return type
MultiDiGraph
-
kgcnn.crystal.graph_builder.
add_radius_bonds
(graph: networkx.classes.multidigraph.MultiDiGraph, radius: float = 5.0, inplace: bool = False) → networkx.classes.multidigraph.MultiDiGraph[source]¶ Adds radius-based edges to a unit cell graph.
- Parameters
- Returns
Graph with added edges.
- Return type
MultiDiGraph
-
kgcnn.crystal.graph_builder.
add_voronoi_bonds
(graph: networkx.classes.multidigraph.MultiDiGraph, min_ridge_area: Optional[float] = None, inplace: bool = False) → networkx.classes.multidigraph.MultiDiGraph[source]¶ Adds Voronoi-based edges to a unit cell graph.
- Parameters
graph (MultiDiGraph) – The unit cell graph to add radius-based edges to.
min_ridge_area (Optional[float], optional) – Threshold value for ridge area between two Voronoi cells. If a ridge area between two voronoi cells is smaller than this value the corresponding edge between the atoms of the cells is excluded from the graph. Defaults to None.
inplace (bool, optional) – Whether to add the edges to the given graph or create a copy with added edges. Defaults to False.
- Returns
Graph with added edges.
- Return type
MultiDiGraph
-
kgcnn.crystal.graph_builder.
get_ridge_area
(ridge_points)[source]¶ Computes the ridge area given ridge points.
Beware that this function, assumes that the ridge points are (roughly) within a flat subspace plane in the 3 dimensional space. It computes the area of the convex hull of the points in three dimensions and then divides it by two, since both sides of the flat convex hull are included.
- Parameters
ridge_points (np.ndarray) – Ridge points to calculate area for.
- Returns
Ridge area for the given points.
- Return type
-
kgcnn.crystal.graph_builder.
get_symmetrized_graph
(structure: Union[pymatgen.core.structure.Structure, pyxtal.pyxtal]) → networkx.classes.multidigraph.MultiDiGraph[source]¶ Builds a unit graph without any edges, but with symmetry information as node attributes.
Each node has a asymmetric_mapping attribute, which contains the id of the symmetry-equivalent atom in the asymmetric unit. Each node has a symmop attribute, which contains the affine matrix to generate the position (in fractional coordinates) of the atom, from its symmetry-equivalent atom position in the asymmetric unit. Each node has a multiplicity attribute, which contains the multiplicity of the atom (how many symmetry-equivalent atoms there are for this node). The resulting graph will have a spacegroup attribute, that specifies the spacegroup of the crystal.
- Parameters
structure (Union[Structure, pyxtal]) – Crystal structure to convert to a graph.
- Raises
ValueError – If the argument is not a pymatgen Structure or pyxtal object.
- Returns
Unit graph with symmetry information, but without any edges for the crystal.
- Return type
MultiDiGraph
-
kgcnn.crystal.graph_builder.
pairwise_diff
(coords1: numpy.ndarray, coords2: numpy.ndarray) → numpy.ndarray[source]¶ Get the pairwise offset difference between two vector sets.
- Parameters
coords1 (np.ndarray) – Coordinates of shape (…, n, 3)
coords2 (np.ndarray) – Coordinates of shape (…, m, 3)
- Returns
Difference values of shape (…, n, m, 3)
- Return type
np.ndarray
-
kgcnn.crystal.graph_builder.
prune_knn_bonds
(graph: networkx.classes.multidigraph.MultiDiGraph, k: int = 12, tolerance: Optional[float] = None, inplace: bool = False) → networkx.classes.multidigraph.MultiDiGraph[source]¶ Prunes edges of a graph to only the k with the smallest distance value.
- Parameters
graph (MultiDiGraph) – The unit cell graph with edges to prune.
k (int, optional) – How many neighbors each node should maximally have. Defaults to 12.
tolerance (Optional[float], optional) – If tolerance is not None, edges with distances of the k-th nearest neighbor plus the tolerance value are included in the graph. Defaults to None.
inplace (bool, optional) – Whether to add the edges to the given graph or create a copy with added edges. Defaults to False.
- Returns
The graph with pruned edges.
- Return type
MultiDiGraph
-
kgcnn.crystal.graph_builder.
prune_radius_bonds
(graph: networkx.classes.multidigraph.MultiDiGraph, radius: float = 4.0, inplace: bool = False) → networkx.classes.multidigraph.MultiDiGraph[source]¶ Prunes edges of a graph with larger distance than the specified radius.
- Parameters
graph (MultiDiGraph) – The unit cell graph with edges to prune.
radius (float, optional) – Distance threshold. Edges with larger distance than this value are removed from the graph. Defaults to 4.
inplace (bool, optional) – Whether to add the edges to the given graph or create a copy with added edges. Defaults to False.
- Returns
The graph with pruned edges.
- Return type
MultiDiGraph
-
kgcnn.crystal.graph_builder.
prune_voronoi_bonds
(graph: networkx.classes.multidigraph.MultiDiGraph, min_ridge_area: Optional[float] = None, inplace: bool = False) → networkx.classes.multidigraph.MultiDiGraph[source]¶ Prunes edges of a graph with a voronoi ridge are smaller then the specified min_ridge_area.
Only works for graphs with edges that contain voronoi_ridge_area as edge attributes. :param graph: The unit cell graph with edges to prune. :type graph: MultiDiGraph :param min_ridge_area: Threshold value for ridge area between two Voronoi cells.
If a ridge area between two voronoi cells is smaller than this value the corresponding edge between the atoms of the cells is excluded from the graph. Defaults to None.
- Parameters
inplace (bool, optional) – Whether to add the edges to the given graph or create a copy with added edges. Defaults to False.
- Returns
The graph with pruned edges.
- Return type
MultiDiGraph
-
kgcnn.crystal.graph_builder.
remove_duplicate_edges
(graph: networkx.classes.multidigraph.MultiDiGraph, inplace=False) → networkx.classes.multidigraph.MultiDiGraph[source]¶ Removes duplicate edges with same offset.
- Parameters
graph (MultiDiGraph) – The unit cell graph with edges to remove.
inplace (bool, optional) – Whether to add the edges to the given graph or create a copy with added edges. Defaults to False.
- Returns
The graph without duplicate edges.
- Return type
MultiDiGraph
-
kgcnn.crystal.graph_builder.
structure_to_empty_graph
(structure: Union[pymatgen.core.structure.Structure, pyxtal.pyxtal], symmetrize: bool = False) → networkx.classes.multidigraph.MultiDiGraph[source]¶ Builds an unit graph without any edges.
- Parameters
structure (Union[Structure, pyxtal]) – Crystal structure to convert to a graph. symmetrize (bool, optional): Whether to include symmetry information attributes (asymmetric_mapping, symmop, multiplicity attributes) in nodes and graph (spacegroup atribute). Defaults to False.
symmetrize (bool) – Whether to get symmetrized graph.
- Raises
ValueError – If the argument is not a pymatgen Structure or pyxtal object.
- Returns
Unit graph without any edges for the crystal.
- Return type
MultiDiGraph
-
kgcnn.crystal.graph_builder.
to_asymmetric_unit_graph
(graph: networkx.classes.multidigraph.MultiDiGraph) → networkx.classes.multidigraph.MultiDiGraph[source]¶ Generates super cell graph representation from unit cell graph representation.
- Parameters
graph (MultiDiGraph) – Unit cell graph to generate asymmetric unit graph for.
- Returns
Corresponding asymmetric unit graph for the given unit cell graph.
- Return type
MultiDiGraph
-
kgcnn.crystal.graph_builder.
to_non_periodic_unit_cell
(graph: networkx.classes.multidigraph.MultiDiGraph, add_reverse_edges: bool = True, inplace: bool = False) → networkx.classes.multidigraph.MultiDiGraph[source]¶ Generates non-periodic graph representation from unit cell graph representation.
- Parameters
graph (MultiDiGraph) – Unit cell graph to generate non-periodic graph for.
add_reverse_edges (bool, optional) – Whether to add incoming edges to atoms that lie outside the central unit cell. Defaults to True.
inplace (bool, optional) – Whether to add distances (distance attribute) to edges. Defaults to False.
- Returns
Corresponding non-periodic graph for the given unit cell graph.
- Return type
MultiDiGraph
-
kgcnn.crystal.graph_builder.
to_supercell_graph
(graph: networkx.classes.multidigraph.MultiDiGraph, size) → networkx.classes.multidigraph.MultiDiGraph[source]¶ Generates super-cell graph representation from unit cell graph representation.
- Parameters
graph (MultiDiGraph) – Unit cell graph to generate super cell graph for.
size (list) – How many cells the crystal will get expanded into each dimension.
- Returns
Corresponding super cell graph for the given unit cell graph.
- Return type
MultiDiGraph
kgcnn.crystal.preprocessor module¶
-
class
kgcnn.crystal.preprocessor.
AsymmetricUnitCell
(output_graph_as_dict: bool = False, lattice: str = 'graph_lattice', species: str = 'node_number', coords: str = 'node_coordinates', charge: str = 'charge')[source]¶ Bases:
Callable
[[pymatgen.core.structure.Structure
],networkx.classes.multidigraph.MultiDiGraph
]Preprocessor that builds asymmetric unit graphs without any edges for crystals.
-
call
(structure: pymatgen.core.structure.Structure) → networkx.classes.multidigraph.MultiDiGraph[source]¶ Builds the crystal graph (networkx.MultiDiGraph) for the pymatgen structure.
- Parameters
structure (Structure) – Structure to convert to a crystal graph.
- Returns
Crystal graph for the provided crystal.
- Return type
MultiDiGraph
-
edge_attributes
= []¶
-
graph_attributes
= ['lattice_matrix', 'spacegroup']¶
-
node_attributes
= ['atomic_number', 'frac_coords', 'coords', 'asymmetric_mapping', 'symmop', 'multiplicity']¶
-
-
class
kgcnn.crystal.preprocessor.
KNNAsymmetricUnitCell
(k: int = 12, tolerance: Optional[float] = 1e-09, **kwargs)[source]¶ Bases:
Callable
[[pymatgen.core.structure.Structure
],networkx.classes.multidigraph.MultiDiGraph
]Preprocessor that builds asymmetric unit graphs with kNN-based edges for crystals.
-
__init__
(k: int = 12, tolerance: Optional[float] = 1e-09, **kwargs)[source]¶ Initializes the crystal preprocessor.
-
call
(structure: pymatgen.core.structure.Structure) → networkx.classes.multidigraph.MultiDiGraph[source]¶ Builds the crystal graph (networkx.MultiDiGraph) for the pymatgen structure.
- Parameters
structure (Structure) – Structure to convert to a crystal graph.
- Returns
Crystal graph for the provided crystal.
- Return type
MultiDiGraph
-
edge_attributes
= ['cell_translation', 'distance', 'symmop', 'offset']¶
-
graph_attributes
= ['lattice_matrix', 'spacegroup']¶
-
node_attributes
= ['atomic_number', 'frac_coords', 'coords', 'multiplicity']¶
-
-
class
kgcnn.crystal.preprocessor.
KNNNonPeriodicUnitCell
(k=12, tolerance=1e-09, **kwargs)[source]¶ Bases:
Callable
[[pymatgen.core.structure.Structure
],networkx.classes.multidigraph.MultiDiGraph
]Preprocessor that builds non-periodic graphs with kNN-based edges for crystals.
-
call
(structure: pymatgen.core.structure.Structure) → networkx.classes.multidigraph.MultiDiGraph[source]¶ Builds the crystal graph (networkx.MultiDiGraph) for the pymatgen structure.
- Parameters
structure (Structure) – Structure to convert to a crystal graph.
- Returns
Crystal graph for the provided crystal.
- Return type
MultiDiGraph
-
edge_attributes
= ['cell_translation', 'distance', 'offset']¶
-
graph_attributes
= ['lattice_matrix']¶
-
node_attributes
= ['atomic_number', 'frac_coords', 'coords']¶
-
-
class
kgcnn.crystal.preprocessor.
KNNSuperCell
(k: int = 12, tolerance: Optional[float] = 1e-09, size: Optional[list] = None, **kwargs)[source]¶ Bases:
Callable
[[pymatgen.core.structure.Structure
],networkx.classes.multidigraph.MultiDiGraph
]Preprocessor that builds super-cell graphs with kNN-based edges for crystals.
-
__init__
(k: int = 12, tolerance: Optional[float] = 1e-09, size: Optional[list] = None, **kwargs)[source]¶ Initializes the crystal preprocessor.
- Parameters
k (int, optional) – How many nearest neighbours to consider for edge selection. Defaults to 12.
tolerance (Optional[float], optional) – If tolerance is not None, edges with distances of the k-th nearest neighbor plus the tolerance value are included in the graph. Defaults to 1e-9.
size (list, optional) – How many cells the crystal will get expanded into each dimension. Defaults to [3,3,3].
-
call
(structure: pymatgen.core.structure.Structure) → networkx.classes.multidigraph.MultiDiGraph[source]¶ Builds the crystal graph (networkx.MultiDiGraph) for the pymatgen structure.
- Parameters
structure (Structure) – Structure to convert to a crystal graph.
- Returns
Crystal graph for the provided crystal.
- Return type
MultiDiGraph
-
edge_attributes
= ['cell_translation', 'distance', 'offset']¶
-
graph_attributes
= ['lattice_matrix']¶
-
node_attributes
= ['atomic_number', 'frac_coords', 'coords']¶
-
-
class
kgcnn.crystal.preprocessor.
KNNUnitCell
(k: int = 12, tolerance: Optional[float] = 1e-09, **kwargs)[source]¶ Bases:
Callable
[[pymatgen.core.structure.Structure
],networkx.classes.multidigraph.MultiDiGraph
]Preprocessor that builds unit cell graphs with kNN-based edges for crystals.
-
__init__
(k: int = 12, tolerance: Optional[float] = 1e-09, **kwargs)[source]¶ Initializes the crystal preprocessor.
-
call
(structure: pymatgen.core.structure.Structure) → networkx.classes.multidigraph.MultiDiGraph[source]¶ Builds the crystal graph (networkx.MultiDiGraph) for the pymatgen structure.
- Parameters
structure (Structure) – Structure to convert to a crystal graph.
- Returns
Crystal graph for the provided crystal.
- Return type
MultiDiGraph
-
edge_attributes
= ['cell_translation', 'distance', 'offset']¶
-
graph_attributes
= ['lattice_matrix']¶
-
node_attributes
= ['atomic_number', 'frac_coords', 'coords']¶
-
-
class
kgcnn.crystal.preprocessor.
RadiusAsymmetricUnitCell
(radius: float = 3.0, **kwargs)[source]¶ Bases:
Callable
[[pymatgen.core.structure.Structure
],networkx.classes.multidigraph.MultiDiGraph
]Preprocessor that builds asymmetric unit graphs with radius-based edges for crystals.
-
__init__
(radius: float = 3.0, **kwargs)[source]¶ Initializes the crystal preprocessor.
- Parameters
radius (float, optional) – Cutoff radius for each atom in Angstrom units. Defaults to 3.0.
-
call
(structure: pymatgen.core.structure.Structure) → networkx.classes.multidigraph.MultiDiGraph[source]¶ Builds the crystal graph (networkx.MultiDiGraph) for the pymatgen structure.
- Parameters
structure (Structure) – Structure to convert to a crystal graph.
- Returns
Crystal graph for the provided crystal.
- Return type
MultiDiGraph
-
edge_attributes
= ['cell_translation', 'distance', 'symmop', 'offset']¶
-
graph_attributes
= ['lattice_matrix', 'spacegroup']¶
-
node_attributes
= ['atomic_number', 'frac_coords', 'coords', 'multiplicity']¶
-
-
class
kgcnn.crystal.preprocessor.
RadiusNonPeriodicUnitCell
(radius: float = 3.0, **kwargs)[source]¶ Bases:
Callable
[[pymatgen.core.structure.Structure
],networkx.classes.multidigraph.MultiDiGraph
]Preprocessor that builds non-periodic graphs with radius-based edges for crystals.
-
__init__
(radius: float = 3.0, **kwargs)[source]¶ Initializes the crystal preprocessor.
- Parameters
radius (float, optional) – Cutoff radius for each atom in Angstrom units. Defaults to 3.0.
-
call
(structure: pymatgen.core.structure.Structure) → networkx.classes.multidigraph.MultiDiGraph[source]¶ Builds the crystal graph (networkx.MultiDiGraph) for the pymatgen structure.
- Parameters
structure (Structure) – Structure to convert to a crystal graph.
- Returns
Crystal graph for the provided crystal.
- Return type
MultiDiGraph
-
edge_attributes
= ['cell_translation', 'distance', 'offset']¶
-
graph_attributes
= ['lattice_matrix']¶
-
node_attributes
= ['atomic_number', 'frac_coords', 'coords']¶
-
-
class
kgcnn.crystal.preprocessor.
RadiusSuperCell
(radius: float = 3.0, size: Optional[list] = None, **kwargs)[source]¶ Bases:
Callable
[[pymatgen.core.structure.Structure
],networkx.classes.multidigraph.MultiDiGraph
]Preprocessor that builds super-cell graphs with radius-based edges for crystals.
-
__init__
(radius: float = 3.0, size: Optional[list] = None, **kwargs)[source]¶ Initializes the crystal preprocessor.
-
call
(structure: pymatgen.core.structure.Structure) → networkx.classes.multidigraph.MultiDiGraph[source]¶ Builds the crystal graph (networkx.MultiDiGraph) for the pymatgen structure.
- Parameters
structure (Structure) – Structure to convert to a crystal graph.
- Returns
Crystal graph for the provided crystal.
- Return type
MultiDiGraph
-
edge_attributes
= ['cell_translation', 'distance', 'offset']¶
-
graph_attributes
= ['lattice_matrix']¶
-
node_attributes
= ['atomic_number', 'frac_coords', 'coords']¶
-
-
class
kgcnn.crystal.preprocessor.
RadiusUnitCell
(radius: float = 3.0, **kwargs)[source]¶ Bases:
Callable
[[pymatgen.core.structure.Structure
],networkx.classes.multidigraph.MultiDiGraph
]Preprocessor that builds unit cell graphs with radius-based edges for crystals.
-
__init__
(radius: float = 3.0, **kwargs)[source]¶ Initializes the crystal preprocessor.
- Parameters
radius (float, optional) – Cutoff radius for each atom in Angstrom units. Defaults to 3.0.
-
call
(structure: pymatgen.core.structure.Structure) → networkx.classes.multidigraph.MultiDiGraph[source]¶ Builds the crystal graph (networkx.MultiDiGraph) for the pymatgen structure.
- Parameters
structure (Structure) – Structure to convert to a crystal graph.
- Returns
Crystal graph for the provided crystal.
- Return type
MultiDiGraph
-
edge_attributes
= ['cell_translation', 'distance', 'offset']¶
-
graph_attributes
= ['lattice_matrix']¶
-
node_attributes
= ['atomic_number', 'frac_coords', 'coords']¶
-
-
class
kgcnn.crystal.preprocessor.
UnitCell
(output_graph_as_dict: bool = False, lattice: str = 'graph_lattice', species: str = 'node_number', coords: str = 'node_coordinates', charge: str = 'charge')[source]¶ Bases:
Callable
[[pymatgen.core.structure.Structure
],networkx.classes.multidigraph.MultiDiGraph
]Preprocessor that builds unit cell graphs without any edges for crystals.
-
call
(structure: pymatgen.core.structure.Structure) → networkx.classes.multidigraph.MultiDiGraph[source]¶ Builds the crystal graph (networkx.MultiDiGraph) for the pymatgen structure.
- Parameters
structure (Structure) – Structure to convert to a crystal graph.
- Returns
Crystal graph for the provided crystal.
- Return type
MultiDiGraph
-
edge_attributes
= []¶
-
graph_attributes
= ['lattice_matrix']¶
-
node_attributes
= ['atomic_number', 'frac_coords', 'coords']¶
-
-
class
kgcnn.crystal.preprocessor.
VoronoiAsymmetricUnitCell
(min_ridge_area: Optional[float] = 0.0, **kwargs)[source]¶ Bases:
Callable
[[pymatgen.core.structure.Structure
],networkx.classes.multidigraph.MultiDiGraph
]Preprocessor that builds asymmetric unit graphs with Voronoi-based edges for crystals.
-
__init__
(min_ridge_area: Optional[float] = 0.0, **kwargs)[source]¶ Initializes the crystal preprocessor.
- Parameters
min_ridge_area (Optional[float], optional) – Threshold value for ridge area between two Voronoi cells. If a ridge area between two voronoi cells is smaller than this value the corresponding edge between the atoms of the cells is excluded from the graph. Defaults to 0.0.
-
call
(structure: pymatgen.core.structure.Structure) → networkx.classes.multidigraph.MultiDiGraph[source]¶ Builds the crystal graph (networkx.MultiDiGraph) for the pymatgen structure.
- Parameters
structure (Structure) – Structure to convert to a crystal graph.
- Returns
Crystal graph for the provided crystal.
- Return type
MultiDiGraph
-
edge_attributes
= ['cell_translation', 'distance', 'symmop', 'offset', 'voronoi_ridge_area']¶
-
graph_attributes
= ['lattice_matrix', 'spacegroup']¶
-
node_attributes
= ['atomic_number', 'frac_coords', 'coords', 'multiplicity']¶
-
-
class
kgcnn.crystal.preprocessor.
VoronoiNonPeriodicUnitCell
(min_ridge_area: Optional[float] = 0.0, **kwargs)[source]¶ Bases:
Callable
[[pymatgen.core.structure.Structure
],networkx.classes.multidigraph.MultiDiGraph
]Preprocessor that builds non-periodic graphs with Voronoi-based edges for crystals.
-
__init__
(min_ridge_area: Optional[float] = 0.0, **kwargs)[source]¶ Initializes the crystal preprocessor.
- Parameters
min_ridge_area (Optional[float], optional) – Threshold value for ridge area between two Voronoi cells. If a ridge area between two voronoi cells is smaller than this value the corresponding edge between the atoms of the cells is excluded from the graph. Defaults to 0.0.
-
call
(structure: pymatgen.core.structure.Structure) → networkx.classes.multidigraph.MultiDiGraph[source]¶ Builds the crystal graph (networkx.MultiDiGraph) for the pymatgen structure.
- Parameters
structure (Structure) – Structure to convert to a crystal graph.
- Returns
Crystal graph for the provided crystal.
- Return type
MultiDiGraph
-
edge_attributes
= ['cell_translation', 'distance', 'offset', 'voronoi_ridge_area']¶
-
graph_attributes
= ['lattice_matrix']¶
-
node_attributes
= ['atomic_number', 'frac_coords', 'coords']¶
-
-
class
kgcnn.crystal.preprocessor.
VoronoiSuperCell
(min_ridge_area: Optional[float] = 0.0, size: Optional[list] = None, **kwargs)[source]¶ Bases:
Callable
[[pymatgen.core.structure.Structure
],networkx.classes.multidigraph.MultiDiGraph
]Preprocessor that builds super-cell graphs with Voronoi-based edges for crystals.
-
__init__
(min_ridge_area: Optional[float] = 0.0, size: Optional[list] = None, **kwargs)[source]¶ Initializes the crystal preprocessor.
- Parameters
min_ridge_area (Optional[float], optional) – Threshold value for ridge area between two Voronoi cells. If a ridge area between two voronoi cells is smaller than this value the corresponding edge between the atoms of the cells is excluded from the graph. Defaults to 0.0.
size (list, optional) – How many cells the crystal will get expanded into each dimension. Defaults to [3,3,3].
-
call
(structure: pymatgen.core.structure.Structure) → networkx.classes.multidigraph.MultiDiGraph[source]¶ Builds the crystal graph (networkx.MultiDiGraph) for the pymatgen structure.
- Parameters
structure (Structure) – Structure to convert to a crystal graph.
- Returns
Crystal graph for the provided crystal.
- Return type
MultiDiGraph
-
edge_attributes
= ['cell_translation', 'distance', 'offset', 'voronoi_ridge_area']¶
-
graph_attributes
= ['lattice_matrix']¶
-
node_attributes
= ['atomic_number', 'frac_coords', 'coords']¶
-
-
class
kgcnn.crystal.preprocessor.
VoronoiUnitCell
(min_ridge_area: Optional[float] = 0.0, **kwargs)[source]¶ Bases:
Callable
[[pymatgen.core.structure.Structure
],networkx.classes.multidigraph.MultiDiGraph
]Preprocessor that builds unit cell graphs with Voronoi-based edges for crystals.
-
__init__
(min_ridge_area: Optional[float] = 0.0, **kwargs)[source]¶ Initializes the crystal preprocessor.
- Parameters
min_ridge_area (Optional[float], optional) – Threshold value for ridge area between two Voronoi cells. If a ridge area between two voronoi cells is smaller than this value the corresponding edge between the atoms of the cells is excluded from the graph. Defaults to 0.0.
-
call
(structure: pymatgen.core.structure.Structure) → networkx.classes.multidigraph.MultiDiGraph[source]¶ Builds the crystal graph (networkx.MultiDiGraph) for the pymatgen structure.
- Parameters
structure (Structure) – Structure to convert to a crystal graph.
- Returns
Crystal graph for the provided crystal.
- Return type
MultiDiGraph
-
edge_attributes
= ['cell_translation', 'distance', 'offset', 'voronoi_ridge_area']¶
-
graph_attributes
= ['lattice_matrix']¶
-
node_attributes
= ['atomic_number', 'frac_coords', 'coords']¶
-