shakenbreak.distortions module#
Module containing functions for applying distortions to defect structures.
- shakenbreak.distortions.apply_dimer_distortion(structure: Structure, site_index: int | None = None, frac_coords: array | None = None, dimer_bond_length: float | None = None, verbose: bool | None = False) dict [source]#
Apply a dimer distortion to a defect structure.
The defect nearest neighbours are determined (using
CrystalNN
with default settings), from them the two closest in distance are selected, which are pushed towards each other so that their inter-atomic distance isdimer_bond_length
Å.If
dimer_bond_length
isNone
(default), then the dimer bond length is estimated usingget_dimer_bond_length
, which uses a database of known covalent bond lengths, or if this fails, using the sum of the covalent radii of the two species.- Parameters:
structure (Structure) – Defect structure.
site_index (Optional[int]) – Index of defect site (for non vacancy defects), using
python
/pymatgen
0-indexing. Defaults to None.frac_coords (Optional[np.array]) – Fractional coordinates of the defect site in the structure (for vacancies). Defaults to None.
dimer_bond_length (float) – The bond length to set the dimer to, in Å. If
None
(default), usesget_dimer_bond_length
to estimate the dimer bond length.verbose (Optional[bool]) – Print information about the dimer distortion. Defaults to False.
- Returns:
- Structure:
Distorted dimer structure
- Return type:
obj
- shakenbreak.distortions.distort(structure: Structure, num_nearest_neighbours: int, distortion_factor: float, site_index: int | None = None, frac_coords: array | None = None, distorted_element: str | list | None = None, distorted_atoms: list | None = None, verbose: bool | None = False) dict [source]#
Applies bond distortions to
num_nearest_neighbours
of the defect (specified bysite_index
(for substitutions or interstitials, counting from 1) orfrac_coords
(for vacancies)).The nearest neighbours to distort are chosen by taking all sites (or those matching
distorted_element
/distorted_atoms
, if provided), then sorting by distance to the defect site (rounded to 2 decimal places) and site index, and then taking the firstnum_nearest_neighbours
of these. If there are multiple non-degenerate combinations of (nearly) equidistant NNs to distort (e.g. cis vs trans when distorting 2 NNs in a 4 NN square coordination), then the combination with distorted NNs closest to each other is chosen.- Parameters:
structure (
Structure
) – Defect structure as apymatgen
Structure
object.num_nearest_neighbours (
int
) – Number of defect nearest neighbours to apply bond distortions todistortion_factor (
float
) – The distortion factor to apply to the bond distance between the defect and nearest neighbours. Typical choice is between 0.4 (-60%) and 1.6 (+60%).site_index (
int
, optional) – Index of defect site in structure (for substitutions or interstitials), usingpython
/pymatgen
0-indexing.frac_coords (
numpy.ndarray
, optional) – Fractional coordinates of the defect site in the structure (for vacancies).distorted_element (
str
, optional) – Neighbouring element(s) to distort, as a string or list of strings. If None, the closest neighbours to the defect will be chosen. (Default: None)distorted_atoms (
list
, optional) – List of atom indices to distort, using 0-indexing (i.e. python / pymatgen indexing). If None, the closest neighbours to the defect will be chosen. (Default: None)verbose (
bool
, optional) – Whether to print distortion information. (Default: False)
- Returns:
Dictionary with distorted defect structure and the distortion parameters.
- Return type:
dict
- shakenbreak.distortions.distort_and_rattle(structure: Structure, distortion_factor: float | str, num_nearest_neighbours: int = 0, site_index: int | None = None, frac_coords: array | None = None, local_rattle: bool = False, stdev: float | None = None, d_min: float | None = None, active_atoms: list | None = None, distorted_element: str | None = None, distorted_atoms: list | None = None, verbose: bool = False, dimer_bond_length: float | None = None, **mc_rattle_kwargs) dict [source]#
Applies bond distortions and rattling to
num_nearest_neighbours
of the defect (specified bysite_index
(for substitutions or interstitials, counting from 1) orfrac_coords
(for vacancies)).Note that by default, rattling is not applied to the defect site or distorted neighbours, but to all other atoms in the supercell, however this can be controlled with the
active_atoms
kwarg. Rattling is performed by applying random displacements to atomic positions, with the displacement distances randomly drawn from a Gaussian distribution of standard deviationstdev
, using a Monte Carlo algorithm which disfavours moves that bring atoms closer thand_min
. Rattling behaviour can be controlled withmc_rattle_kwargs
, see therattle()
docstring for more info. Note that dimer distortions can be generated by settingdistortion_factor
to “Dimer”; seeapply_dimer_distortion()
for more info.The nearest neighbours to distort are chosen by taking all sites (or those matching
distorted_element
/distorted_atoms
, if provided), then sorting by distance to the defect site (rounded to 2 decimal places) and site index, and then taking the firstnum_nearest_neighbours
of these. If there are multiple non-degenerate combinations of (nearly) equidistant NNs to distort (e.g. cis vs trans when distorting 2 NNs in a 4 NN square coordination), then the combination with distorted NNs closest to each other is chosen.- Parameters:
structure (
Structure
) – Defect structure as apymatgen
Structure
object.distortion_factor (
float
) – The distortion factor or distortion name (“Dimer”) to apply to the bond distance between the defect and nearest neighbours. Typical choice is between 0.4 (-60%) and 1.6 (+60%).num_nearest_neighbours (
int
) – Number of defect nearest neighbours to apply bond distortions to. Default is 0.site_index (
int
, optional) – Index of defect site in structure (for substitutions or interstitials), usingpython
/pymatgen
0-indexing.frac_coords (
numpy.ndarray
, optional) – Fractional coordinates of the defect site in the structure (for vacancies).local_rattle (
bool
) – Whether to apply random displacements that tail off as we move away from the defect site. If False, all supercell sites are rattled with the same amplitude. (Default: False)stdev (
float
) – Standard deviation (in Angstroms) of the Gaussian distribution from which random atomic displacement distances are drawn during rattling. Default is set to 10% of the bulk nearest neighbour distance.d_min (
float
) – Minimum interatomic distance (in Angstroms) in the rattled structure. Monte Carlo rattle moves that put atoms at distances less than this will be heavily penalised. Default is to set this to 80% of the nearest neighbour distance in the defect supercell (ignoring interstitials).active_atoms (
list
, optional) – List (or array) of which atomic indices should undergo Monte Carlo rattling. Default is to apply rattle to all atoms except the defect and the bond-distorted neighbours.distorted_element (
str
, optional) – Neighbouring element to distort. If None, the closest neighbours to the defect will be chosen. (Default: None)distorted_atoms (
list
, optional) – List of atom indices which can undergo bond distortions, using 0-indexing (i.e. python / pymatgen indexing). If None, the closest neighbours to the defect will be chosen. (Default: None)verbose (
bool
) – Whether to print distortion information. (Default: False)dimer_bond_length (
float
) – The bond length to set generated dimers in dimer distortions to, in Å. IfNone
(default), usesdistortions.get_dimer_bond_length
to estimate the dimer bond length.**mc_rattle_kwargs (
dict
) –Additional keyword arguments to pass to
hiphive
’smc_rattle
function. These include:- max_disp (
float
): Maximum atomic displacement (in Å) during Monte Carlo rattling. Rarely occurs and is used primarily as a safety net. (Default: 2.0)
- max_disp (
- max_attempts (
int
): Limit for how many attempted rattle moves are allowed a single atom.
- max_attempts (
- active_atoms (
list
): List of the atomic indices which should undergo Monte Carlo rattling. By default, all atoms are rattled. (Default: None)
- active_atoms (
- seed (
int
): Seed for setting up NumPy random state from which rattle random displacements are generated.
- seed (
- Returns:
Dictionary with distorted defect structure and the distortion parameters.
- Return type:
dict
- shakenbreak.distortions.get_dimer_bond_length(species_1: str | Element | Species | DummySpecies, species_2: str | Element | Species | DummySpecies)[source]#
Get the estimated dimer bond length between two species, using the
get_bond_length()
function frompymatgen
(which uses a database of known covalent bond lengths), or if this fails, using the sum of the covalent radii of the two species.- Parameters:
species_1 (SpeciesLike) – First species.
species_2 (SpeciesLike) – Second species.
- Returns:
The estimated dimer bond length between the two species.
- Return type:
float
- shakenbreak.distortions.local_mc_rattle(structure: Structure, site_index: int | None = None, frac_coords: array | None = None, stdev: float | None = None, d_min: float | None = None, verbose: bool | None = False, n_iter: int = 1, active_atoms: list | None = None, nbr_cutoff: float = 5, width: float = 0.1, max_attempts: int = 5000, max_disp: float = 2.0, seed: int = 42) Structure [source]#
Given a
pymatgen
Structure
object, apply random displacements to all atomic positions, with the displacement distances randomly drawn from a Gaussian distribution of standard deviationstdev
. The random displacements tail off as we move away from the defect site.- Parameters:
structure (
Structure
) – Structure as a pymatgen objectsite_index (
int
, optional) – Index of defect site in structure (for substitutions or interstitials), usingpython
/pymatgen
0-indexing.frac_coords (
numpy.ndarray
, optional) – Fractional coordinates of the defect site in the structure (for vacancies).stdev (
float
) – Standard deviation (in Angstroms) of the Gaussian distribution from which random atomic displacement distances are drawn during rattling. Default is set to 10% of the bulk nearest neighbour distance.d_min (
float
) – Minimum interatomic distance (in Angstroms) in the rattled structure. Monte Carlo rattle moves that put atoms at distances less than this will be heavily penalised. Default is to set this to 80% of the nearest neighbour distance in the defect supercell (ignoring interstitials).verbose (
bool
) – Whether to print information about the rattling process, if rattling initially fails with initiald_min
.n_iter (
int
) – Number of Monte Carlo cycles to perform. (Default: 1)active_atoms (
list
, optional) – List of which atomic indices should undergo Monte Carlo rattling. (Default: None)nbr_cutoff (
float
) – The radial cutoff distance (in Angstroms) used to construct the list of atomic neighbours for checking interatomic distances. (Default: 5)width (
float
) – Width of the Monte Carlo rattling error function, in Angstroms. (Default: 0.1)max_disp (
float
) – Maximum atomic displacement (in Angstroms) during Monte Carlo rattling. Rarely occurs and is used primarily as a safety net. (Default: 2.0)max_attempts (
int
) – Maximum Monte Carlo rattle move attempts allowed for a single atom; if this limit is reached anException
is raised. (Default: 5000)seed (
int
) – Seed for NumPy random state from which random rattle displacements are generated. (Default: 42)
- Returns:
Rattled
pymatgen
Structure
object- Return type:
Structure
- shakenbreak.distortions.rattle(structure: Structure, stdev: float | None = None, d_min: float | None = None, verbose: bool = False, n_iter: int = 1, active_atoms: list | None = None, nbr_cutoff: float = 5, width: float = 0.1, max_attempts: int = 5000, max_disp: float = 2.0, seed: int = 42) Structure [source]#
Given a
pymatgen
Structure
object, apply random displacements to all atomic positions, with the displacement distances randomly drawn from a Gaussian distribution of standard deviationstdev
, using a Monte Carlo algorithm which disfavours moves that bring atoms closer thand_min
.- Parameters:
structure (
Structure
) – Structure as a pymatgen objectstdev (
float
) – Standard deviation (in Angstroms) of the Gaussian distribution from which random atomic displacement distances are drawn during rattling. Default is set to 10% of the bulk nearest neighbour distance.d_min (
float
) – Minimum interatomic distance (in Angstroms) in the rattled structure. Monte Carlo rattle moves that put atoms at distances less than this will be heavily penalised. Default is to set this to 80% of the nearest neighbour distance in the defect supercell.verbose (
bool
) – Whether to print information about the rattling process, if rattling initially fails with initiald_min
.n_iter (
int
) – Number of Monte Carlo cycles to perform. (Default: 1)active_atoms (
list
, optional) – List of which atomic indices should undergo Monte Carlo rattling. If not set, rattles all atoms in the structure. (Default: None)nbr_cutoff (
float
) – The radial cutoff distance (in Angstroms) used to construct the list of atomic neighbours for checking interatomic distances. (Default: 5)width (
float
) – Width of the Monte Carlo rattling error function, in Angstroms. (Default: 0.1)max_disp (
float
) – Maximum atomic displacement (in Angstroms) during Monte Carlo rattling. Rarely occurs and is used primarily as a safety net. (Default: 2.0)max_attempts (
int
) – Maximum Monte Carlo rattle move attempts allowed for a single atom; if this limit is reached anException
is raised. (Default: 5000)seed (
int
) – Seed for NumPy random state from which random rattle displacements are generated. (Default: 42)
- Returns:
Rattled
pymatgen
Structure
object- Return type:
Structure