shakenbreak.analysis module#

Module containing functions to analyse rattled and bond-distorted defect structure relaxations.

shakenbreak.analysis.analyse_defect_site(structure: Structure, name: str | None = None, site_num: int | None = None, vac_site: list | None = None) tuple[source]#

Analyse coordination environment and bond distances to nearest neighbours of defect site.

Parameters:
  • structure (Structure) – pymatgen Structure object to analyse

  • name (str) – Defect name for printing. (Default: None)

  • site_num (int) – Defect site index in the structure, starting from 1 (VASP rather than python indexing). (Default: None)

  • vac_site (list) – For vacancies, the fractional coordinates of the vacant lattice site. (Default: None)

Returns:

Tuple of coordination analysis and bond length DataFrames, respectively.

Return type:

tuple

shakenbreak.analysis.analyse_structure(defect_species: str, structure: Structure, output_path: str = '.') tuple[source]#

Analyse the local distortion of the input defect structure. Requires access to the distortion_metadata.json file generated with ShakeNBreak to read info about defect site. If lacking this, can alternatively use analyse_defect_site().

Parameters:
  • defect_species (str) – Defect name including charge (e.g. ‘vac_1_Cd_0’)

  • structure (Structure) – Defect structure to analyse

  • output_path (str) – Path to directory containing distortion_metadata.json (Default: ‘.’, current directory)

Returns:

Tuple of coordination analysis and bond length DataFrames, respectively.

Return type:

tuple

shakenbreak.analysis.calculate_struct_comparison(defect_structures_dict: dict, metric: str = 'max_dist', ref_structure: str | float | Structure = 'Unperturbed', stol: float = 0.5, ltol: float = 0.3, angle_tol: float = 5, min_dist: float = 0.1, verbose: bool = False) dict[source]#

Calculate either the summed atomic displacement, with metric = “disp”, or the maximum distance between matched atoms, with metric = “max_dist”, (default) between each distorted structure in defect_struct_dict, and either ‘Unperturbed’ or a specified structure (ref_structure). For metric = “disp”, atomic displacements below min_dist (in Å, default 0.1 Å) are considered noise and omitted from the sum, to reduce supercell size dependence.

Parameters:
  • defect_structures_dict (dict) – Dictionary of bond distortions and corresponding (final) structures (as pymatgen Structure objects).

  • metric (str) – Structure comparison metric to use. Either summed atomic displacement normalised to the average free length per atom (‘disp’) or the maximum distance between matched atoms (‘max_dist’, default). (Default: “max_dist”)

  • ref_structure (str or float or Structure) – Structure to use as a reference for comparison (to compute atomic displacements). Either as a key from defect_structures_dict (e.g. ‘-0.4’ for “Bond_Distortion_-40.0%”) or a pymatgen Structure object (to compare with a specific external structure). (Default: “Unperturbed”)

  • stol (float) – Site tolerance used for structural comparison (via pymatgen’s StructureMatcher), as a fraction of the average free length per atom := ( V / Nsites ) ** (1/3). If output contains too many ‘NaN’ values, this likely needs to be increased. (Default: 0.5)

  • ltol (float) – Length tolerance used for structural comparison (via pymatgen’s StructureMatcher). (Default: 0.3)

  • angle_tol (float) – Angle tolerance used for structural comparison (via pymatgen’s StructureMatcher). (Default: 5)

  • min_dist (float) – Minimum atomic displacement threshold to include in atomic displacements sum (in Å, default 0.1 Å).

  • verbose (bool) – Whether to print information message about structures being compared.

Returns:

Dictionary matching bond distortions to structure comparison metric (disp or max_dist).

Return type:

dict

shakenbreak.analysis.compare_structures(defect_structures_dict: dict, defect_energies_dict: dict, ref_structure: str | float | Structure = 'Unperturbed', stol: float = 0.5, units: str = 'eV', min_dist: float = 0.1, display_df: bool = True, verbose: bool = True) None | DataFrame[source]#

Compare final bond-distorted structures with either ‘Unperturbed’ or a specified structure (ref_structure), and calculate the summed atomic displacement (in Å), and maximum distance between matched atomic sites (in Å).

Parameters:
  • defect_structures_dict (dict) – Dictionary mapping bond distortion to (relaxed) structure

  • defect_energies_dict (dict) – Dictionary matching distortion to final energy (eV), as produced by _organize_data().

  • ref_structure (str or float or Structure) – Structure to use as a reference for comparison (to compute atomic displacements). Either as a key from defect_structures_dict (e.g. ‘-0.4’ for ‘Bond_Distortion_-40.0%’) or a pymatgen Structure object (to compare with a specific external structure). (Default: “Unperturbed”)

  • stol (float) – Site tolerance used for structural comparison (via pymatgen’s StructureMatcher), as a fraction of the average free length per atom := ( V / Nsites ) ** (1/3). If structure comparison output contains too many ‘NaN’ values, this likely needs to be increased. (Default: 0.5)

  • units (str) – Energy units label for outputs (either ‘eV’ or ‘meV’). Should be the same as the units in defect_energies_dict, as this does not modify the supplied values. (Default: “eV”)

  • min_dist (float) – Minimum atomic displacement threshold to include in atomic displacements sum (in Å, default 0.1 Å).

  • display_df (bool) – Whether to display the structure comparison DataFrame interactively in Jupyter/Ipython (Default: True).

  • verbose (bool) – Whether to print information message about structures being compared.

Returns:

DataFrame containing structural comparison results (summed normalised atomic displacement and maximum distance between matched atomic sites), and relative energies.

Return type:

pd.DataFrame

shakenbreak.analysis.get_energies(defect_species: str, output_path: str = '.', units: str = 'eV', verbose: bool = True) dict[source]#

Parse final energies for each bond distortion and store them in a dictionary matching the bond distortion to the final energy in eV.

Parameters:
  • defect_species (str) – Defect name including charge (e.g. ‘vac_1_Cd_0’)

  • output_path (str) – Path to top-level directory containing defect_species subdirectories. (Default: current directory)

  • distortion_increment (float) – Bond distortion increment. Recommended values: 0.1-0.3 (Default: 0.1)

  • units (str) – Energy units for outputs (either ‘eV’ or ‘meV’). (Default: “eV”)

  • verbose (bool) – Whether to print information about energy lowering distortions, if found. (Default: True)

Returns:

Dictionary matching bond distortions to final energies in eV.

Return type:

dict

shakenbreak.analysis.get_gs_distortion(defect_energies_dict: dict) tuple[source]#

Calculate energy difference between Unperturbed structure and lowest energy distortion. Returns the energy (in eV) and bond distortion of the ground-state relative to Unperturbed. If Unperturbed not present, returns (None, ground-state bond distortion).

Parameters:

defect_energies_dict (dict) – Dictionary matching distortion to final energy, as produced by get_energies() or _sort_data.

Returns:

(Energy difference, ground state bond distortion)

Return type:

tuple

shakenbreak.analysis.get_homoionic_bonds(structure: Structure, elements: list, radius: float | None = 3.3, verbose: bool = True) dict[source]#

Returns a list of homoionic bonds for the given element. These bonds are often formed by the defect neighbouts to accomodate charge deficiency.

Parameters:
  • structure (Structure) – pymatgen Structure object to analyse

  • elements (list) – List of element symbols (wihout oxidation state) for which to find the homoionic bonds (e.g. [“Te”, “Se”]).

  • radius (float, optional) – Distance cutoff to look for homoionic bonds. Defaults to 3.3 A.

  • verbose (bool, optional) – Whether or not to print the list of homoionic bonds.

Returns:

dictionary with homoionic bonds, matching site to the homoionic neighbours and distances (A) (e.g. {‘O(1)’: {‘O(2)’: ‘2.0 A’, ‘O(3)’: ‘2.0 A’}})

Return type:

dict

shakenbreak.analysis.get_site_magnetizations(defect_species: str, distortions: list, output_path: str = '.', threshold: float = 0.1, defect_site: int | None = None, orbital_projections: bool | None = False, verbose: bool | None = True) dict | None[source]#

For given distortions, find sites with significant magnetization and return as dictionary. Only implemented for VASP calculations.

Parameters:
  • defect_species (str) – Name of defect including charge state (e.g. ‘vac_1_Cd_0’)

  • distortions (list) – List of distortions to analyse (e.g. [‘Unperturbed’, 0.1, -0.2])

  • output_path (str) – Path to top-level directory containing defect_species subdirectories. (Default: current directory)

  • threshold (float, optional) – Magnetization threshold to consider site. (Default: 0.1)

  • defect_site (int or list, optional) – Site index or fractional coordinates of the defect. If not specified, will try to read from distortion_metadata.json. If this file is not present, will not include the distance between the defect and the sites with significant magnetization in the dataframe.

  • orbital_projections (bool, optional) – Whether to print orbital projections. If not necessary, set to False (faster). (Default: False)

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

Returns:

Dictionary matching distortion to DataFrame containing magnetization info.

Return type:

dict

shakenbreak.analysis.get_structures(defect_species: str, output_path: str = '.', bond_distortions: list | None = None, code: str | None = 'vasp', structure_filename: str | None = 'CONTCAR') dict[source]#

Import all structures found with rattling & bond distortions, and store them in a dictionary matching the bond distortion to the final structure. By default, will read the structures from the distortion subdirectories present in each defect folder. If only certain distortions should be parsed, use the argument bond_distortions to specify them.

Parameters:
  • defect_species (str) – Defect name including charge (e.g. ‘vac_1_Cd_0’)

  • output_path (str) – Path to top-level directory containing defect_species subdirectories. (Default: current directory.

  • bond_distortions (list, optional) – List of bond distortions whose structures will be analysed (i.e. [-0.2, 0.0, 0.5]). By default, all distortions present in subdirectories are considered. (Default: None)

  • code (str, optional) – Code used for the geometry relaxations. (Default: vasp)

  • structure_filename (str, optional) – Name of the file containing the structure. (Default: CONTCAR)

Returns:

Dictionary of bond distortions and corresponding final structures.

Return type:

dict