doped.generation module

Code to generate Defect objects and supercell structures for ab-initio calculations.

class doped.generation.DefectsGenerator(structure: Structure | Atoms | str | Path, extrinsic: str | list | dict | None = None, interstitial_coords: list | None = None, generate_supercell: bool = True, charge_state_gen_kwargs: dict | None = None, supercell_gen_kwargs: dict[str, int | float | bool] | None = None, interstitial_gen_kwargs: dict | bool | None = None, target_frac_coords: list | None = None, processes: int | None = None, **kwargs)[source]

Bases: MSONable

Class for generating doped DefectEntry objects.

Generates doped DefectEntry objects for defects in the input host structure. By default, generates all intrinsic defects, but extrinsic defects (impurities) can also be created using the extrinsic argument.

Interstitial sites are generated using Voronoi tessellation by default (found to be the most reliable) using the get_interstitial_sites function, which also generates adsorbate sites if the structure is determined to be a slab/layer/rod. This can be controlled using the interstitial_gen_kwargs argument. Alternatively, a list of interstitial sites (or single interstitial site) can be manually specified using the interstitial_coords argument.

By default, supercells are generated for each defect using the doped get_ideal_supercell_matrix() function (see docstring), with default settings of min_image_distance = 10 (minimum distance between periodic images of 10 Å), min_atoms = 50 (minimum 50 atoms in the supercell) and ideal_threshold = 0.1 (allow up to 10% larger supercell if it is a diagonal expansion of the primitive or conventional cell). This uses a custom algorithm in doped to efficiently search over possible supercell transformations and identify that with the minimum number of atoms (hence computational cost) that satisfies the minimum image distance, number of atoms and ideal_threshold constraints. These settings can be controlled by specifying keyword arguments with supercell_gen_kwargs, which are passed to get_ideal_supercell_matrix() (e.g. for a minimum image distance of 15 Å with at least 100 atoms, supercell_gen_kwargs would be: {'min_image_distance': 15, 'min_atoms': 100}). If the input structure already satisfies these constraints (for the same number of atoms as the doped-generated supercell), then it will be used. Alternatively if generate_supercell = False, then no supercell is generated and the input structure is used as the defect & bulk supercell. (Note this may give a slightly different (but fully equivalent) set of coordinates).

The algorithm for determining defect entry names is to use the pymatgen defect name (e.g. v_Cd, Cd_Te etc.) for vacancies / antisites / substitutions, unless there are multiple inequivalent sites for the defect, in which case the point group of the defect site is appended (e.g. v_Cd_Td, Cd_Te_Td etc.), and if this is still not unique, then element identity and distance to the nearest neighbour of the defect site is appended (e.g. v_Cd_Td_Te2.83, Cd_Te_Td_Cd2.83 etc.). For interstitials, the same naming scheme is used, but the point group is always appended to the pymatgen name.

Possible charge states for the defects are estimated using the probability of the corresponding defect element oxidation state, the magnitude of the charge state, and the maximum magnitude of the host oxidation states (i.e. how ‘charged’ the host is). Large (absolute) charge states, low probability oxidation states and/or greater charge/oxidation state magnitudes than that of the host are disfavoured. This can be controlled using the probability_threshold (default = 0.0075) or padding (default = 1) keys in the charge_state_gen_kwargs parameter, which are passed to the guess_defect_charge_states() function. The input and computed values used to guess charge state probabilities are provided in the DefectEntry.charge_state_guessing_log attributes. See docs for examples of modifying the generated charge states, and the docstrings of charge_state_probability() & get_vacancy_charge_states() for more details on the charge state guessing algorithm. The doped algorithm was found to give optimal performance in terms of efficiency and completeness (see JOSS paper), but of course may not be perfect in all cases, so make sure to critically consider the estimated charge states for your system!

Note that Wyckoff letters can depend on the ordering of elements in the conventional standard structure, for which doped uses the spglib convention.

Parameters:
  • structure (Structure) – Structure of the host material, either as a pymatgen Structure, ASE Atoms or path to a structure file (e.g. CONTCAR). If this is not the primitive unit cell, it will be reduced to the primitive cell for defect generation, before supercell generation.

  • extrinsic (Union[str, list, dict]) – List or dict of elements (or string for single element) to be used for extrinsic defect generation (i.e. dopants/impurities). If a list is provided, all possible substitutional defects for each extrinsic element will be generated. If a dict is provided, the keys should be the host elements to be substituted, and the values the extrinsic element(s) to substitute in; as a string or list. In both cases, all possible extrinsic interstitials are generated.

  • interstitial_coords (list) – List of fractional coordinates (corresponding to the input structure), or a single set of fractional coordinates, to use as interstitial defect site(s). Default (when interstitial_coords not specified) is to automatically generate interstitial sites using Voronoi tessellation. The input interstitial_coords are converted to DefectsGenerator.prim_interstitial_coords, which are the corresponding fractional coordinates in DefectsGenerator.primitive_structure (which is used for defect generation), along with the multiplicity and equivalent coordinates, sorted according to the doped convention.

  • generate_supercell (bool) – Whether to generate a supercell for the output defect entries (using the custom algorithm in doped which efficiently searches over possible supercell transformations and identifies that with the minimum number of atoms (hence computational cost) that satisfies the minimum image distance, number of atoms and ideal_threshold constraints – which can be controlled with supercell_gen_kwargs). If False, then the input structure is used as the defect & bulk supercell (note this may give a slightly different (but fully equivalent) set of coordinates). Default is True.

  • charge_state_gen_kwargs (dict) – Keyword arguments to pass to the guess_defect_charge_states function (such as probability_threshold (default = 0.0075, used for substitutions and interstitials) and padding (default = 1, used for vacancies)) to control defect charge state generation.

  • supercell_gen_kwargs (dict) – Keyword arguments to pass to the get_ideal_supercell_matrix function (such as min_image_distance (default = 10), min_atoms (default = 50), ideal_threshold (default = 0.1), force_cubic – which enforces a (near-)cubic supercell output (default = False), or force_diagonal (default = False)).

  • interstitial_gen_kwargs (dict, bool) – Keyword arguments to be passed to get_interstitial_sites() (such as min_dist (0.9 Å), clustering_tol (0.8 Å), symmetry_preference (0.1 Å), stol (0.32), tight_stol (0.02), symprec (0.01), vacuum_radius (1.5 * bulk bond length) – see its docstring, parentheses indicate default values), or InterstitialGenerator if interstitial_coords is specified. If set to False, interstitial generation will be skipped entirely.

  • target_frac_coords (list) – Defects are placed at the closest equivalent site to these fractional coordinates in the generated supercells. Default is [0.5, 0.5, 0.5] (i.e. the supercell centre, to aid visualisation).

  • processes (int) – Number of processes to use for multiprocessing. If not set, defaults to one less than the number of CPUs available.

  • **kwargs – Additional keyword arguments for defect generation. Options: {defect}_elements where {defect} is vacancy, substitution, or interstitial, in which cases only those defects of the specified elements will be generated (where {defect}_elements is a list of element symbol strings). Setting {defect}_elements to an empty list will skip defect generation for that defect type entirely. {defect}_charge_states to specify the charge states to use for all defects of that type (as a list of integers). neutral_only to only generate neutral charge states.

Key attributes:
defect_entries (dict):

Dictionary of {defect_species: DefectEntry} for all defect entries (with charge state and supercell properties) generated.

defects (dict):

Dictionary of {defect_type: [Defect, ...]} for all defect objects generated.

primitive_structure (Structure):

Primitive cell structure of the host used to generate defects.

supercell_matrix (Matrix):

Matrix to generate defect/bulk supercells from the primitive cell structure.

bulk_supercell (Structure):

Supercell structure of the host (equal to self.primitive_structure * self.supercell_matrix).

conventional_structure (Structure):

Conventional cell structure of the host according to the Bilbao Crystallographic Server (BCS) definition, used to determine defect site Wyckoff labels and multiplicities.

prim_interstitial_coords (list):

List of interstitial coordinates in the primitive cell.

DefectsGenerator input parameters are also set as attributes.

add_charge_states(defect_entry_name: str, charge_states: list | int)[source]

Add additional DefectEntrys with the specified charge states to self.defect_entries.

Parameters:
  • defect_entry_name (str) – Name of defect entry to add charge states to. Doesn’t need to include the charge state.

  • charge_states (list) – List of charge states to add to defect entry (e.g. [-2, -3]), or a single charge state (e.g. -2).

as_dict()[source]

JSON-serializable dict representation of DefectsGenerator.

defect_generator_info()[source]

Prints information about the defects that have been generated.

classmethod from_dict(d)[source]

Reconstructs DefectsGenerator object from a dict representation created using DefectsGenerator.as_dict().

Parameters:

d (dict) – dict representation of DefectsGenerator.

Returns:

DefectsGenerator object

classmethod from_json(filename: str | Path)[source]

Load a DefectsGenerator object from a json(.gz) file.

Note that .json.gz files can be loaded directly.

Parameters:

filename (PathLike) – Filename of json file to load DefectsGenerator object from.

Returns:

DefectsGenerator object

remove_charge_states(defect_entry_name: str, charge_states: list | int)[source]

Remove DefectEntrys with the specified charge states from self.defect_entries.

Parameters:
  • defect_entry_name (str) – Name of defect entry to remove charge states from. Doesn’t need to include the charge state.

  • charge_states (list) – List of charge states to add to defect entry (e.g. [-2, -3]), or a single charge state (e.g. -2).

to_json(filename: str | Path | None = None)[source]

Save the DefectsGenerator object as a json file, which can be reloaded with the DefectsGenerator.from_json() class method.

Note that file extensions with “.gz” will be automatically compressed (recommended to save space)!

Parameters:

filename (PathLike) – Filename to save json file as. If None, the filename will be set as {Chemical Formula}_defects_generator.json.gz where {Chemical Formula} is the chemical formula of the host material.

doped.generation.charge_state_probability(charge_state: int, defect_el_oxi_state: int, defect_el_oxi_probability: float, max_host_oxi_magnitude: int, return_log: bool = False) float | dict[source]

Function to estimate the probability of a given defect charge state, using the probability of the corresponding defect element oxidation state, the magnitude of the charge state, and the maximum magnitude of the host oxidation states (i.e. how ‘charged’ the host is).

Disfavours large (absolute) charge states, low probability oxidation states and greater charge/oxidation state magnitudes than that of the host. This charge state probability function is primarily intended for substitutions and interstitials, while the get_vacancy_charge_states() function is used for vacancies.

Specifically, the overall probability is given by the product of these probability factors:

  • The probability of the corresponding oxidation state of the defect element (e.g. Na_i^+1 has Na in the +1 oxidation state), as given by its prevalence in the ICSD.

  • The magnitude of the charge state; with a probability function: 1/|charge_state|^(2/3)

  • The magnitude of the charge state relative to the max host oxidation state (i.e. how ‘charged’ the host is); with a probability function: 1/(2*|charge_state - max_host_oxi_magnitude|)^(2/3) if charge_state > max_host_oxi_magnitude, otherwise 1.

  • The magnitude of the defect element oxidation state relative to the max host oxidation state; with a probability function: 1/(2*|defect_el_oxi_state - max_host_oxi_magnitude|)^(2/3) if defect_el_oxi_state > max_host_oxi_magnitude, otherwise 1.

Note that neutral charge states are always included.

This probability function was found to give optimal performance in terms of efficiency and completeness when tested against other approaches (see the doped JOSS paper: https://doi.org/10.21105/joss.06433) but of course may not be perfect in all cases, so make sure to critically consider the estimated charge states for your system!

Parameters:
  • charge_state (int) – Charge state of defect.

  • defect_el_oxi_state (int) – Oxidation state of defect element.

  • defect_el_oxi_probability (float) – Probability of oxidation state of defect element.

  • max_host_oxi_magnitude (int) – Maximum host oxidation state magnitude.

  • return_log (bool) – If true, returns a dictionary of input & computed values used to determine charge state probability. Default is False.

Returns:

Probability of the defect charge state (between 0 and 1) if return_log is False, otherwise a dictionary of input & computed values used to determine charge state probability.

doped.generation.closest_site_info(defect_entry_or_defect: DefectEntry | Defect, n: int = 1, element_list: list[str] | None = None)[source]

Return the element and distance (rounded to 2 decimal places) of the nth closest site to the defect site in the input DefectEntry or Defect object.

If there are multiple elements with the same distance (to within ~0.01 Å), then the preference for which element to return is controlled by element_list. If element_list is not provided, then it is set to match the order of appearance in the structure composition.

If n > 1, then it returns the nth closest site, where the nth site must be at least 0.02 Å further away than the (n-1)th site.

Parameters:
  • defect_entry_or_defect (Union[DefectEntry, Defect]) – DefectEntry or Defect object, to get neighbour info.

  • n (int) – Return the element symbol and distance for the nth closest site. Default is 1, corresponding to the (1st) closest neighbour.

  • element_list (list) – Sorted list of elements in the host structure to govern the preference of elemental symbols to return, when the distance to multiple neighbours with different elements is the same. Default is to use _get_element_list(), which follows the order of appearance of elements in the Structure composition.

Returns:

Element symbol and distance (rounded to 2 decimal places) of the nth closest site to the defect site.

Return type:

str

doped.generation.get_defect_entry_from_defect(defect: ~doped.core.Defect, defect_supercell: ~pymatgen.core.structure.Structure, charge_state: int, dummy_species: ~pymatgen.core.periodic_table.DummySpecies = DummySpecies X0+, sc_defect_frac_coords: ~numpy.ndarray | None = None)[source]

Generate a doped DefectEntry object from a doped Defect object.

This is used to describe a Defect with a specified simulation cell. To set the sc_defect_frac_coords attribute for DefectEntry (fractional coordinates of the defect in the defect_supercell), either dummy_species must be present in the defect_supercell (which is taken as the defect site), or sc_defect_frac_coords must be set.

Parameters:
  • defect (Defect) – doped/pymatgen Defect object.

  • defect_supercell (Structure) – Defect supercell structure.

  • charge_state (int) – Charge state of the defect.

  • dummy_species (DummySpecies) – Dummy species present in the defect_supercell structure, used to determine sc_defect_frac_coords. If not found and sc_defect_frac_coords is not set, DefectEntry.sc_defect_frac_coords is set to None. Default is DummySpecies("X").

  • sc_defect_frac_coords (np.ndarray) – Fractional coordinates of the defect in the defect supercell. If not set and dummy_species is not found in the defect_supercell, DefectEntry.sc_defect_frac_coords is set to None. Default is None.

Returns:

doped DefectEntry object.

Return type:

DefectEntry

doped.generation.get_defect_name_from_defect(defect: Defect, element_list: list[str] | None = None, symprec: float = 0.01)[source]

Get the doped/SnB defect name from a Defect object.

Parameters:
  • defect (Defect) – Defect object.

  • element_list (list) – Sorted list of elements in the host structure, so that closest_site_info() returns deterministic results (in case two different elements located at the same distance from defect site). Default is None.

  • symprec (float) – Symmetry tolerance for spglib. Default is 0.01.

Returns:

Defect name.

Return type:

str

doped.generation.get_defect_name_from_entry(defect_entry: DefectEntry, element_list: list | None = None, symprec: float | None = None, relaxed: bool = True)[source]

Get the doped/SnB defect name from a DefectEntry object.

Note: If relaxed = True (default), then this tries to use the defect_entry.defect_supercell to determine the site symmetry. This will thus give the relaxed defect point symmetry if this is a DefectEntry created from parsed defect calculations. However, it should be noted that this is not guaranteed to work in all cases; namely for non-diagonal supercell expansions, or sometimes for non-scalar supercell expansion matrices (e.g. a 2x1x2 expansion)(particularly with high-symmetry materials) which can mess up the periodicity of the cell. doped tries to automatically check if this is the case, and will warn you if so.

This can also be checked by using this function on your doped generated defects:

from doped.generation import get_defect_name_from_entry
for defect_name, defect_entry in defect_gen.items():
    print(defect_name,
          get_defect_name_from_entry(defect_entry, relaxed=False),
          get_defect_name_from_entry(defect_entry), "\n")

And if the point symmetries match in each case, then using this function on your parsed relaxed DefectEntry objects should correctly determine the final relaxed defect symmetry (and closest site info) – otherwise periodicity-breaking prevents this.

Parameters:
  • defect_entry (DefectEntry) – DefectEntry object.

  • element_list (list) – Sorted list of elements in the host structure, so that closest_site_info() returns deterministic results (in case two different elements located at the same distance from defect site). Default is None.

  • symprec (float) – Symmetry tolerance for spglib. Default is 0.01 for unrelaxed structures, 0.2 for relaxed (to account for residual structural noise). You may want to adjust for your system (e.g. if there are very slight octahedral distortions etc).

  • relaxed (bool) – If False, determines the site symmetry using the defect site in the unrelaxed bulk supercell, otherwise tries to determine the point symmetry of the relaxed defect in the defect supercell). Default is True.

Returns:

Defect name.

Return type:

str

doped.generation.get_ideal_supercell_matrix(structure: Structure, min_image_distance: float = 10.0, min_atoms: int = 50, force_cubic: bool = False, force_diagonal: bool = False, ideal_threshold: float = 0.1, pbar: tqdm | None = None) ndarray | None[source]

Determine the ideal supercell matrix for a given structure, based on the minimum image distance, minimum number of atoms and ideal_threshold for further expanding if a diagonal expansion of the primitive/conventional cell is possible.

The ideal supercell is the smallest possible supercell which has a minimum image distance (i.e. minimum distance between periodic images of atoms/sites in a lattice) greater than min_image_distance (default = 10 Å – which is a typical threshold value used in DFT defect supercell calculations) and a number of atoms greater than min_atoms (default = 50). Once these criteria have been reached, doped will then continue searching up to supercell sizes (numbers of atoms) 1 + ideal_threshold times larger (rounded up) to see if they return a diagonal expansion of the primitive/conventional cell (which can make later visualisation and analysis much easier) – if so, this larger supercell will be returned.

This search for the ideal supercell transformation matrix is performed using the find_ideal_supercell function from doped.utils.supercells (see its docstring for more details), which efficiently scans over possible supercell matrices and identifies that with the minimum image distance and most cubic-like supercell shape. The advantage of this over that in pymatgen-analysis-defects is that it avoids the find_optimal_cell_shape function from ASE (which currently does not work for rotated matrices, is inefficient, and optimises based on cubic- like shape rather than minimum image distance), giving greatly reduced supercell sizes for a given minimum image distance.

If force_cubic or force_diagonal are True, then the CubicSupercellTransformation from pymatgen is used to identify any simple near-cubic supercell transformations which satisfy the minimum image distance and atom number criteria.

Parameters:
  • structure (Structure) – Primitive unit cell structure to generate supercell for.

  • min_image_distance (float) – Minimum image distance in Å of the supercell (i.e. minimum distance between periodic images of atoms/sites in the lattice). (Default = 10.0)

  • min_atoms (int) – Minimum number of atoms allowed in the supercell. (Default = 50)

  • force_cubic (bool) – Enforce usage of CubicSupercellTransformation from pymatgen for supercell generation. (Default = False)

  • force_diagonal (bool) – If True, return a transformation with a diagonal transformation matrix. (Default = False)

  • ideal_threshold (float) – Threshold for increasing supercell size (beyond that which satisfies min_image_distance and min_atoms`) to achieve an ideal supercell matrix (i.e. a diagonal expansion of the primitive or conventional cell). Supercells up to 1 + perfect_cell_threshold times larger (rounded up) are trialled, and will instead be returned if they yield an ideal transformation matrix. (Default = 0.1; i.e. 10% larger than the minimum size)

  • pbar (tqdm) – tqdm progress bar object to update (for internal doped usage). Default is None.

Returns:

Ideal supercell matrix (np.ndarray) or None if no suitable supercell could be found.

doped.generation.get_interstitial_sites(host_structure: Structure, interstitial_gen_kwargs: dict[str, Any] | None = None) list[source]

Generate candidate interstitial sites using Voronoi analysis.

This function uses a similar approach to that in VoronoiInterstitialGenerator from pymatgen-analysis-defects, but with modifications to make interstitial generation much faster, fix bugs with interstitial grouping (which could lead to undesired dropping of candidate sites) and achieve better control over site placement in order to favour sites which are higher-symmetry and furthest from the host lattice atoms (typically the most favourable interstitial sites).

If the structure is detected to have a significant vacuum volume (determined using vacuum_radius), thus corresponding to a slab/layer/line structure, then the pymatgen adsorbate finder will be used to generate further candidate interstitial sites.

The logic for picking interstitial sites is as follows:

  • Generate all candidate sites using (efficient) Voronoi analysis

  • Remove any sites which are within min_dist of any host atoms

  • Cluster the remaining sites using a tolerance of clustering_tol and symmetry-preference of symmetry_preference (see _doped_cluster_frac_coords)

  • Determine the multiplicities and symmetry-equivalent coordinates of the clustered sites using doped symmetry functions.

  • Group the clustered sites by symmetry using (looser) site matching as controlled by stol.

  • From each group, pick the site with the highest symmetry and furthest distance from the host atoms, if its min_dist is no more than symmetry_preference (0.1 Å by default) smaller than the site with the largest min_dist (to the host atoms).

(Parameters mentioned here can be supplied via interstitial_gen_kwargs as noted in the args section below.)

The motivation for favouring high symmetry interstitial sites and then distance to host atoms is because higher symmetry interstitial sites are typically the more intuitive sites for placement, cleaner, easier for analysis etc., and interstitials are often lowest-energy when furthest from host atoms (i.e. in the largest interstitial voids – particularly for fully-ionised charge states), and so this approach tries to strike a balance between these two goals.

One caveat to this preference for high symmetry interstitial sites, is that they can also be slightly more prone to being stuck in local minima on the PES, and so as always it is highly recommended to use ShakeNBreak or another structure-searching technique to account for symmetry-breaking when performing defect relaxations!

You can see what Cartesian distance the chosen stol corresponds to using the get_stol_equiv_dist function.

Parameters:
  • host_structure (Structure) – Host structure.

  • interstitial_gen_kwargs (dict) –

    Keyword arguments for interstitial generation. Supported kwargs:

    • min_dist (float):

      Minimum distance from host atoms for interstitial sites. Defaults to 0.9 Å.

    • clustering_tol (float):

      Tolerance for clustering interstitial sites. Defaults to 0.8 Å.

    • symmetry_preference (float):

      Symmetry preference for interstitial site selection. Defaults to 0.1 Å.

    • stol (float):

      Structure matcher tolerance for looser site matching. Defaults to 0.32.

    • tight_stol (float):

      Structure matcher tolerance for tighter site matching. Defaults to 0.02.

    • symprec (float):

      Symmetry precision for (symmetry-)equivalent site determination. Defaults to 0.01.

    • vacuum_radius (float):

      Tolerance radius for determining if a significant vacuum volume is present in the structure (if any Voronoi node has a minimum distance to a host atom greater than this value), in which case adsorbate sites will also be included for interstitial generation. Defaults to 1.5 times the bulk bond length.

Returns:

List of interstitial sites as fractional coordinates

Return type:

list

doped.generation.get_neighbour_distances_and_symbols(site: PeriodicSite, structure: Structure, n: int = 1, element_list: list[str] | None = None, dist_tol_prefactor: float = 3.0, min_dist: float = 0.02)[source]

Get a list of sorted tuples of (distance, element) for the nth closest sites to the input site in structure, where each consecutive neighbour in the list must be at least 0.02 Å further away than the previous neighbour.

If there are multiple elements with the same distance (to within ~0.01 Å), then the preference for which element to return is controlled by element_list. If element_list is not provided, then it is set to match the order of appearance in the structure composition.

For efficiency, this function dynamically increases the radius when searching for nearest neighbours as necessary, using the input dist_tol_prefactor and n to estimate a reasonable starting range.

Parameters:
  • site (PeriodicSite) – Site to get neighbour info.

  • structure (Structure) – Structure containing the site and neighbours.

  • n (int) – Return the element symbol and distance tuples for the nth closest neighbours. Default is 1, corresponding to only the nearest neighbour.

  • element_list (list) – Sorted list of elements in the host structure to govern the preference of elemental symbols to return, when the distance to multiple neighbours with different elements is the same. Default is to use the order of appearance of elements in the Structure composition.

  • dist_tol_prefactor (float) – Initial distance tolerance prefactor to use when searching for neighbours, where the initial search radius is set to dist_tol_prefactor * sqrt(n) Å. This value is dynamically updated as needed, and so should only be adjusted if providing odd systems with very small/large bond lengths for which efficiency improvements are possible. Default is 3.0.

  • min_dist (float) – Minimum distance in Å of neighbours from site. Intended as a distance tolerance to exclude the site itself from the neighbour list, for which the default value (0.02 Å) should be perfectly fine in most cases. Set to 0 to include the site itself in the output list (in which case it counts toward n).

Returns:

Sorted list of tuples of (distance, element) for the nth closest neighbours to the input site.

Return type:

list

doped.generation.get_oxi_probabilities(element_symbol: str) dict[source]

Get a dictionary of oxidation states and their probabilities for an element.

Tries to get the probabilities from the pymatgen tabulated ICSD oxidation state probabilities, and if not available, uses the common oxidation states of the element.

Parameters:

element_symbol (str) – Element symbol.

Returns:

Dictionary of oxidation states (ints) and their probabilities (floats).

Return type:

dict

doped.generation.get_stol_equiv_dist(stol: float, structure: Structure) float[source]

Get the equivalent Cartesian distance of a given stol value for a given Structure.

stol is a site tolerance parameter used in pymatgen StructureMatcher functions, defined as the fraction of the average free length per atom := ( V / Nsites ) ** (1/3).

Parameters:
  • stol (float) – Site tolerance parameter.

  • structure (Structure) – Structure to get equivalent distance for.

Returns:

Equivalent Cartesian distance for the given stol value.

Return type:

float

doped.generation.get_vacancy_charge_states(vacancy: Vacancy, padding: int = 1) list[int][source]

Get the estimated charge states for a vacancy defect, which is from +/-padding to the fully-ionised vacancy charge state (a.k.a. the vacancy oxidation state).

e.g. for vacancies in Sb2O5 (doi.org/10.1021/acs.chemmater.3c03257), the fully-ionised charge states for V_Sb and V_O are -5 and +2 respectively (i.e. the negative of the elemental oxidation states in Sb2O5), so the estimated charge states would be from +1 to -5 for V_Sb and from +2 to -1 for V_O for the default padding of 1.

This probability function was found to give optimal performance in terms of efficiency and completeness when tested against other approaches (see the doped JOSS paper: https://doi.org/10.21105/joss.06433) but of course may not be perfect in all cases, so make sure to critically consider the estimated charge states for your system!

Parameters:
  • vacancy (Defect) – A doped Vacancy object.

  • padding (int) – Padding for vacancy charge states, such that the vacancy charge states are set to range(vacancy oxi state, padding) if vacancy oxidation state is negative, or to range(-padding, vacancy oxi state) if positive. Default is 1.

Returns:

A list of estimated charge states for the defect.

Return type:

list[int]

doped.generation.guess_defect_charge_states(defect: Defect, probability_threshold: float = 0.0075, padding: int = 1, return_log: bool = False) list[int] | tuple[list[int], list[dict]][source]

Guess the possible stable charge states of a defect.

This function estimates the probabilities of the charge states of a defect, using the probability of the corresponding defect element oxidation states, the magnitudes of the charge states, and the maximum magnitude of the host oxidation states (i.e. how ‘charged’ the host is), and returns a list of charge states that have an estimated probability greater than the probability_threshold.

Disfavours large (absolute) charge states, low probability oxidation states and greater charge/oxidation state magnitudes than that of the host. Note that neutral charge states are always included.

For specific details on the probability functions employed, see the charge_state_probability (for substitutions and interstitials) and get_vacancy_charge_states() (for vacancies) functions.

These probability functions were found to give optimal performance in terms of efficiency and completeness when tested against other approaches (see the doped JOSS paper: https://doi.org/10.21105/joss.06433) but of course may not be perfect in all cases, so make sure to critically consider the estimated charge states for your system!

Parameters:
  • defect (Defect) – doped Defect object.

  • probability_threshold (float) – Probability threshold for including defect charge states (for substitutions and interstitials). Default is 0.0075.

  • padding (int) – Padding for vacancy charge states, such that the vacancy charge states are set to range(vacancy oxi state, padding) if vacancy oxidation state is negative, or to range(-padding, vacancy oxi state) if positive. Default is 1.

  • return_log (bool) – If true, returns a tuple of the defect charge states and a list of dictionaries of input & computed values used to determine charge state probability. Default is False.

Returns:

List of defect charge states (int) or a tuple of the defect charge states (list) and a list of dictionaries of input & computed values used to determine charge state probability.

doped.generation.name_defect_entries(defect_entries: list[DefectEntry | Defect], element_list: list[str] | None = None)[source]

Create a dictionary of {name: DefectEntry} from a list of DefectEntry objects, where the names are set according to the default doped algorithm; which is to use the pymatgen defect name (e.g. v_Cd, Cd_Te etc.) for vacancies/antisites/substitutions, unless there are multiple inequivalent sites for the defect, in which case the point group of the defect site is appended (e.g. v_Cd_Td, Cd_Te_Td etc.), and if this is still not unique, then element identity and distance to the nearest neighbour of the defect site is appended (e.g. v_Cd_Td_Te2.83, Cd_Te_Td_Cd2.83 etc.). Names do not yet have charge states included.

For interstitials, the same naming scheme is used, but the point group is always appended to the pymatgen defect name.

If still not unique after the 3rd nearest neighbour info, then “a, b, c” etc. is appended to the name of different defects to distinguish.

Parameters:
  • defect_entries (list) – List of DefectEntry or Defect objects to name.

  • element_list (list) – Sorted list of elements in the host structure, so that closest_site_info() returns deterministic results (in case two different elements located at the same distance from defect site). Default is None.

Returns:

Dictionary of {name: DefectEntry} objects.

Return type:

dict

doped.generation.sort_defect_entries(defect_entries: dict | list, element_list: list | None = None) dict | list[source]

Sort defect entries for deterministic behaviour; for output and when reloading DefectsGenerator objects, with DefectThermodynamics entries (particularly for deterministic plotting behaviour), and with DefectsParser objects.

Sorts defect entries by defect type (vacancies, substitutions, interstitials), then by order of appearance of elements in the host composition, then by periodic group (main groups 1, 2, 13-18 first, then TMs), then by atomic number, then (for defect entries of the same type) sort by name and charge state (from positive to negative).

Parameters:
  • defect_entries (dict or list) – Dictionary (in the format: {defect_entry_name: defect_entry}) or list of defect entries to sort.

  • element_list (list, optional) – List of elements present, used to determine preferential ordering. If None, determined by _get_element_list(), which orders by appearance of elements in the host composition, then by periodic group (main groups 1, 2, 13-18 first, then TMs), then by atomic number.

Returns:

Sorted dictionary or list of defect entries.