pymatgen.analysis.prototypes package

This module is intended to match crystal structures against known crystallographic “prototype” structures.

In this module, the AflowPrototypeMatcher uses the AFLOW LIBRARY OF CRYSTALLOGRAPHIC PROTOTYPES. If using this particular class, please cite their publication appropriately:

Mehl, M. J., Hicks, D., Toher, C., Levy, O., Hanson, R. M., Hart, G., & Curtarolo, S. (2017). The AFLOW library of crystallographic prototypes: part 1. Computational Materials Science, 136, S1-S828. https://doi.org/10.1016/j.commatsci.2017.01.017

The module also contains functions for getting the protostructure labels from a variety of symmetry detection libraries (spglib, moyopy, aflow-sym). The protostructure label is defined as the canonicalized aflow label with the alphabetically sorted chemical system appended - aflow_sym_label:chemsys.

The utilities for determining the protostructure label are upstreamed from the aviary package (https://github.com/CompRhys/aviary). If using these functions, please cite the following publications:

Goodall, R. E., Parackal, A. S., Faber, F. A., Armiento, R., & Lee, A. A. (2022). Rapid discovery of stable materials by coordinate-free coarse graining. Science advances, 8(30), eabn4117. https://doi.org/10.1126/sciadv.abn4117

Parackal, A. S., Goodall, R. E., Faber, F. A., & Armiento, R. (2024). Identifying crystal structures beyond known prototypes from x-ray powder diffraction spectra. Physical Review Materials, 8(10), 103801. https://doi.org/10.1103/PhysRevMaterials.8.103801

class AflowPrototypeMatcher(initial_ltol: float = 0.2, initial_stol: float = 0.3, initial_angle_tol: float = 5)[source]

Bases: object

This class will match structures to their crystal prototypes, and will attempt to group species together to match structures derived from prototypes (e.g. an A_xB_1-x_C from a binary prototype), and will give these the names the “-like” suffix.

This class uses data from the AFLOW LIBRARY OF CRYSTALLOGRAPHIC PROTOTYPES. If using this class, please cite their publication appropriately:

Mehl, M. J., Hicks, D., Toher, C., Levy, O., Hanson, R. M., Hart, G., & Curtarolo, S. (2017). The AFLOW library of crystallographic prototypes: part 1. Computational Materials Science, 136, S1-S828. https://doi.org/10.1016/j.commatsci.2017.01.017

Tolerances as defined in StructureMatcher. Tolerances will be gradually decreased until only a single match is found (if possible).

Parameters:
  • initial_ltol (float) – fractional length tolerance.

  • initial_stol (float) – site tolerance.

  • initial_angle_tol (float) – angle tolerance.

get_prototypes(structure: Structure) list[dict] | None[source]

Get prototype(s) structures for a given input structure. If you use this method in your work, please cite the appropriate AFLOW publication:

Mehl, M. J., Hicks, D., Toher, C., Levy, O., Hanson, R. M., Hart, G., & Curtarolo, S. (2017). The AFLOW library of crystallographic prototypes: part 1. Computational Materials Science, 136, S1-S828. https://doi.org/10.1016/j.commatsci.2017.01.017

Parameters:

structure (Structure) – structure to match

Returns:

A list of dicts with keys “snl” for the matched prototype and

”tags”, a dict of tags (“mineral”, “strukturbericht” and “aflow”) of that prototype. This should be a list containing just a single entry, but it is possible a material can match multiple prototypes.

Return type:

list[dict] | None

canonicalize_element_wyckoffs(element_wyckoffs: str, spg_num: int | str) str[source]

Given an element ordering, canonicalize the associated Wyckoff positions based on the alphabetical weight of equivalent choices of origin.

Parameters:
  • element_wyckoffs (str) – wyckoff substring section from aflow_label with the wyckoff letters for different elements separated by underscores.

  • spg_num (int | str) – International space group number.

Returns:

element_wyckoff string with canonical ordering of the wyckoff letters.

Return type:

str

count_crystal_dof(protostructure_label: str) int[source]

Count number of free parameters in coarse-grained protostructure_label representation: how many degrees of freedom would remain to optimize during a crystal structure relaxation.

Parameters:

protostructure_label (str) – label constructed as aflow_label:chemsys where aflow_label is an AFLOW-style prototype label chemsys is the alphabetically sorted chemical system.

Returns:

Number of free-parameters in given prototype

Return type:

int

count_crystal_sites(protostructure_label: str) int[source]

Count number of sites from protostructure_label.

Parameters:

protostructure_label (str) – label constructed as aflow_label:chemsys where aflow_label is an AFLOW-style prototype label chemsys is the alphabetically sorted chemical system.

Returns:

Number of free-parameters in given prototype

Return type:

int

count_distinct_wyckoff_letters(protostructure_label: str) int[source]

Count number of distinct Wyckoff letters in protostructure_label.

Parameters:

protostructure_label (str) – label constructed as aflow_label:chemsys where aflow_label is an AFLOW-style prototype label chemsys is the alphabetically sorted chemical system.

Returns:

number of distinct Wyckoff letters in protostructure_label

Return type:

int

count_values_for_wyckoff(element_wyckoffs: list[str], counts: list[str], spg_num: str, lookup_dict: dict[str, dict[str, int]])[source]

Count values from a lookup table and scale by wyckoff multiplicities.

count_wyckoff_positions(protostructure_label: str) int[source]

Count number of Wyckoff positions in protostructure_label.

Parameters:

protostructure_label (str) – label constructed as aflow_label:chemsys where aflow_label is an AFLOW-style prototype label chemsys is the alphabetically sorted chemical system.

Returns:

number of distinct Wyckoff positions in protostructure_label

Return type:

int

get_anonymous_formula_from_prototype_formula(prototype_formula: str) str[source]

Get an anonymous formula from a prototype formula.

get_formula_from_protostructure_label(protostructure_label: str) str[source]

Get a formula from a protostructure label.

get_protostructure_label(struct: Structure, method: Literal['aflow', 'spglib', 'moyopy'], raise_errors: bool = False, **kwargs) str | None[source]

Get protostructure label for a pymatgen Structure.

Parameters:
  • struct (Structure) – pymatgen Structure

  • method (Literal["aflow", "spglib", "moyopy"]) – Method to use for symmetry

  • detection

  • raise_errors (bool) – Whether to raise errors or annotate them. Defaults to False.

  • **kwargs – Additional arguments for the specific method

Returns:

protostructure_label which is constructed as aflow_label:chemsys or

explanation of failure if symmetry detection failed and raise_errors is False.

Return type:

str

get_protostructure_label_from_aflow(struct: Structure, raise_errors: bool = False, aflow_executable: str | None = None) str[source]

Get protostructure label for a pymatgen Structure. Make sure you’re running a recent version of the aflow CLI as there’s been several breaking changes. This code was tested under v3.2.12. The protostructure label is constructed as aflow_label:chemsys.

Install guide: https://aflow.org/install-aflow/#install_aflow

http://aflow.org/install-aflow/install-aflow.sh -o install-aflow.sh chmod 555 install-aflow.sh ./install-aflow.sh –slim

Parameters:
  • struct (Structure) – pymatgen Structure

  • aflow_executable (str) – path to aflow executable. Defaults to which(“aflow”).

  • raise_errors (bool) – Whether to raise errors or annotate them. Defaults to False.

Returns:

protostructure_label which is constructed as aflow_label:chemsys or

explanation of failure if symmetry detection failed and raise_errors is False.

Return type:

str

get_protostructure_label_from_moyopy(struct: Structure, raise_errors: bool = False, symprec: float = 0.1) str | None[source]

Get AFLOW prototype label using Moyopy for symmetry detection.

Parameters:
  • struct (Structure) – pymatgen Structure object.

  • raise_errors (bool) – Whether to raise errors or annotate them. Defaults to False.

  • symprec (float) – Initial symmetry precision for Moyopy. Defaults to 0.1.

Returns:

protostructure_label which is constructed as aflow_label:chemsys or

explanation of failure if symmetry detection failed and raise_errors is False.

Return type:

str

get_protostructure_label_from_spg_analyzer(spg_analyzer: SpacegroupAnalyzer, raise_errors: bool = False) str[source]

Get protostructure label for pymatgen SpacegroupAnalyzer.

Parameters:
  • spg_analyzer (SpacegroupAnalyzer) – pymatgen SpacegroupAnalyzer object.

  • raise_errors (bool) – Whether to raise errors or annotate them. Defaults to False.

Returns:

protostructure_label which is constructed as aflow_label:chemsys or

explanation of failure if symmetry detection failed and raise_errors is False.

Return type:

str

get_protostructure_label_from_spglib(struct: Structure, raise_errors: bool = False, init_symprec: float = 0.1, fallback_symprec: float | None = 1e-05) str[source]

Get AFLOW prototype label for pymatgen Structure.

Parameters:
  • struct (Structure) – pymatgen Structure object.

  • raise_errors (bool) – Whether to raise errors or annotate them. Defaults to False.

  • init_symprec (float) – Initial symmetry precision for spglib. Defaults to 0.1.

  • fallback_symprec (float) – Fallback symmetry precision for spglib if first symmetry detection failed. Defaults to 1e-5.

Returns:

protostructure_label which is constructed as aflow_label:chemsys or

explanation of failure if symmetry detection failed and raise_errors is False.

Return type:

str

get_protostructures_from_aflow_label_and_composition(aflow_label: str, composition: Composition) list[str][source]

Get a canonicalized string for the prototype.

Parameters:
  • aflow_label (str) – AFLOW-style prototype label

  • composition (Composition) – pymatgen Composition object

Returns:

List of possible protostructure labels that can be generated

from combinations of the input aflow_label and composition.

Return type:

list[str]

get_prototype_formula_from_composition(composition: Composition) str[source]

An anonymized formula. Unique species are arranged in alphabetical order and assigned ascending alphabets. This format is used in the aflow structure prototype labelling scheme.

Parameters:

composition (Composition) – Pymatgen Composition to process

Returns:

anonymized formula where the species are in alphabetical order

Return type:

str

get_prototype_from_protostructure(protostructure_label: str) str[source]

Get a canonicalized string for the prototype. This prototype should be the same for all isopointal protostructures.

Parameters:

protostructure_label (str) – label constructed as aflow_label:chemsys where aflow_label is an AFLOW-style prototype label chemsys is the alphabetically sorted chemical system.

Returns:

Canonicalized AFLOW-style prototype label

Return type:

str

get_random_structure_for_protostructure(protostructure_label: str, **kwargs) Structure[source]

Generate a random structure for a given prototype structure.

NOTE that due to the random nature of the generation, the output structure may be higher symmetry than the requested prototype structure.

Parameters:
  • protostructure_label (str) – label constructed as aflow_label:chemsys where aflow_label is an AFLOW-style prototype label chemsys is the alphabetically sorted chemical system.

  • **kwargs – Keyword arguments to pass to pyxtal().from_random()

sort_and_score_element_wyckoffs(element_wyckoffs: str) tuple[str, int][source]

Determines the order or Wyckoff positions when canonicalizing AFLOW labels.

Parameters:

element_wyckoffs (str) – wyckoff substring section from aflow_label with the wyckoff letters for different elements separated by underscores.

Returns:

containing - str: sorted Wyckoff position substring for AFLOW-style prototype label - int: integer score to rank order when canonicalizing

Return type:

tuple

split_alpha_numeric(s: str) dict[str, list[str]][source]

Split a string into separate lists of alpha and numeric groups.

Parameters:

s (str) – The input string to split.

Returns:

A dictionary with keys ‘alpha’ and ‘numeric’,

each containing a list of the respective groups.

Return type:

dict[str, list[str]]