doped.core module
Core functions and classes for defects in doped.
- class doped.core.Defect(structure: Structure, site: PeriodicSite, multiplicity: int | None = None, oxi_state: float | int | str | None = None, equivalent_sites: list[PeriodicSite] | None = None, symprec: float = 0.01, angle_tolerance: float = 5, user_charges: list[int] | None = None, **doped_kwargs)[source]
Bases:
Defect
doped
Defect
object.Subclass of
pymatgen.analysis.defects.core.Defect
with additional attributes and methods used bydoped
.- Parameters:
structure (Structure) – The structure in which to create the defect. Typically the primitive structure of the host crystal for defect generation, and/or the calculation supercell for defect parsing.
site (PeriodicSite) – The defect site in the structure.
multiplicity (int) – The multiplicity of the defect in the structure.
oxi_state (float, int or str) – The oxidation state of the defect. If not specified, this will be determined automatically.
equivalent_sites (list[PeriodicSite]) – A list of equivalent sites for the defect in the structure.
symprec (float) – Tolerance for symmetry finding.
angle_tolerance (float) – Angle tolerance for symmetry finding.
user_charges (list[int]) – User specified charge states. If specified,
get_charge_states
will return this list. IfNone
or empty list the charge states will be determined automatically.**doped_kwargs – Additional keyword arguments to define doped-specific attributes (listed below), in the form
doped_attribute_name=value
. (e.g.wyckoff = "4a"
).
- as_dict()[source]
JSON-serializable dict representation of
Defect
.Needs to be redefined because attributes not explicitly specified in subclasses, which is required for monty functions.
- property defect_site: PeriodicSite
The defect site in the structure.
Re-written from
pymatgen-analysis-defects
version to be far more efficient, when used in loops (e.g. for calculating defect concentrations as functions of chemical potentials, temperature etc.).
- classmethod from_json(filename: str)[source]
Load a
Defect
object from a json file.- Parameters:
filename (str) – Filename of json file to load
Defect
from.- Returns:
Defect
object
- get_charge_states(padding: int = 1) list[int] [source]
Refactored version of
pymatgen-analysis-defects
’sget_charge_states
to not break whenoxi_state
is not set.
- get_supercell_structure(sc_mat: ndarray | None = None, target_frac_coords: ndarray | None = None, return_sites: bool = False, min_image_distance: float = 10.0, min_atoms: int = 50, force_cubic: bool = False, force_diagonal: bool = False, ideal_threshold: float = 0.1, min_length: float | None = None, dummy_species: str | None = None) Structure [source]
Generate the simulation supercell for a defect.
Redefined from the parent class to allow the use of
target_frac_coords
to place the defect at the closest equivalent site to the target fractional coordinates in the supercell, while keeping the supercell fixed (to avoid any issues with defect parsing). Also returns information about equivalent defect sites in the supercell.If
sc_mat
is None, then the supercell is generated automatically using thedoped
algorithm described in theget_ideal_supercell_matrix
function docstring indoped.generation
.- Parameters:
sc_mat (3x3 matrix) – Transformation matrix of
self.structure
to create the supercell. If None, then automatically computed usingget_ideal_supercell_matrix
fromdoped.generation
.target_frac_coords (3x1 matrix) – If set, the defect will be placed at the closest equivalent site to these fractional coordinates (using self.equivalent_sites).
return_sites (bool) – If True, returns a tuple of the defect supercell, defect supercell site and list of equivalent supercell sites.
dummy_species (str) – Dummy species to highlight the defect position (for visualizing vacancies).
min_image_distance (float) – Minimum image distance in Å of the generated supercell (i.e. minimum distance between periodic images of atoms/sites in the lattice), if
sc_mat
is None. (Default = 10.0)min_atoms (int) – Minimum number of atoms allowed in the generated supercell, if
sc_mat
is None. (Default = 50)force_cubic (bool) – Enforce usage of
CubicSupercellTransformation
frompymatgen
for supercell generation (ifsc_mat
is None). (Default = False)force_diagonal (bool) – If True, return a transformation with a diagonal transformation matrix (if
sc_mat
is None). (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 to1 + perfect_cell_threshold
times larger (rounded up) are trialled, and will instead be returned if they yield an ideal transformation matrix (ifsc_mat
is None). (Default = 0.1; i.e. 10% larger than the minimum size)min_length (float) – Same as
min_image_distance
(kept for compatibility).
- Returns:
The defect supercell structure. If
return_sites
is True, also returns the defect supercell site and list of equivalent supercell sites.
- to_json(filename: str | None = None)[source]
Save the
Defect
object to a json file, which can be reloaded with the `` Defect``.from_json()`` class method.- Parameters:
filename (str) – Filename to save json file as. If None, the filename will be set as “{Defect.name}.json”.
- property volume: float
The volume (in ų) of the structure in which the defect is created (i.e.
Defect.structure
).Ensures volume is only computed once when calculating defect concentrations in loops (e.g. for calculating defect concentrations as functions of chemical potentials, temperature etc.).
- class doped.core.DefectEntry(defect: ~doped.core.Defect, charge_state: int, sc_entry: ~pymatgen.entries.computed_entries.ComputedStructureEntry, corrections: dict[str, float] = <factory>, corrections_metadata: dict[str, ~typing.Any] = <factory>, sc_defect_frac_coords: tuple[float, float, float] | None = None, bulk_entry: ~pymatgen.entries.computed_entries.ComputedEntry | None = None, entry_id: str | None = None, name: str = '', calculation_metadata: dict = <factory>, degeneracy_factors: dict = <factory>, conventional_structure: ~pymatgen.core.structure.Structure | None = None, conv_cell_frac_coords: ~numpy.ndarray | None = None, equiv_conv_cell_frac_coords: list[~numpy.ndarray] = <factory>, _BilbaoCS_conv_cell_vector_mapping: list[int] = <factory>, wyckoff: str | None = None, charge_state_guessing_log: dict = <factory>, defect_supercell: ~pymatgen.core.structure.Structure | None = None, defect_supercell_site: ~pymatgen.core.sites.PeriodicSite | None = None, equivalent_supercell_sites: list[~pymatgen.core.sites.PeriodicSite] = <factory>, bulk_supercell: ~pymatgen.core.structure.Structure | None = None, _bulk_entry_energy: float | None = None, _bulk_entry_hash: int | None = None, _sc_entry_energy: float | None = None, _sc_entry_hash: int | None = None)[source]
Bases:
DefectEntry
Subclass of
pymatgen.analysis.defects.thermo.DefectEntry
with additional attributes used bydoped
.- Core Attributes:
- defect:
doped
/pymatgen
defect object corresponding to the defect in the entry.- charge_state:
Charge state of the defect.
- sc_entry:
pymatgen
ComputedStructureEntry
for the defect supercell.- sc_defect_frac_coords:
The fractional coordinates of the defect in the supercell.
- bulk_entry:
pymatgen
ComputedEntry
for the bulk supercell reference. Required for calculating the defect formation energy.- corrections:
A dictionary of energy corrections which are summed and added to the defect formation energy.
- corrections_metadata:
A dictionary that acts as a generic container for storing information about how the corrections were calculated. Only used for debugging and plotting purposes.
- Parsing Attributes:
- calculation_metadata:
A dictionary of calculation parameters and data, used to perform charge corrections and compute formation energies.
- degeneracy_factors:
A dictionary of degeneracy factors contributing to the total degeneracy of the defect species (such as spin and configurational degeneracy etc). This is an important factor in the defect concentration equation (see discussion in https://doi.org/10.1039/D2FD00043A and https://doi.org/10.1039/D3CS00432E), and so affects the output of the defect concentration / Fermi level functions. This can be edited by the user if the
doped
defaults are not appropriate (e.g.doped
assumes singlet (S=0) state for even-electron defects and doublet (S=1/2) state for odd-electron defects, which is typically the case but can have triplets (S=1) or other multiplets for e.g. bipolarons, quantum / d-orbital / magnetic defects etc.); see https://doped.readthedocs.io/en/latest/Tips.html#spin-polarisation for discussion.
- Generation Attributes:
- name:
The
doped
-generated name of the defect entry. See docstrings ofDefectsGenerator
for the doped naming algorithm.- conventional_structure:
Conventional cell structure of the host according to the Bilbao Crystallographic Server (BCS) definition, used to determine defect site Wyckoff labels and multiplicities.
- conv_cell_frac_coords:
Fractional coordinates of the defect in the conventional cell.
- equiv_conv_cell_frac_coords:
Symmetry-equivalent defect positions in fractional coordinates of the conventional cell.
- _BilbaoCS_conv_cell_vector_mapping:
A vector mapping the lattice vectors of the
spglib
-defined conventional cell to that of the Bilbao Crystallographic Server definition (for most space groups the definitions are the same).- wyckoff:
Wyckoff label of the defect site.
- charge_state_guessing_log:
A log of the input & computed values used to determine charge state probabilities.
- defect_supercell:
pymatgen
Structure
object of the defect supercell.- defect_supercell_site:
pymatgen
PeriodicSite
object of the defect in the defect supercell.- equivalent_supercell_sites:
List of
pymatgen
PeriodicSite
objects of symmetry-equivalent defect sites in the defect supercell.- bulk_supercell:
pymatgen
Structure
object of the bulk (pristine, defect-free) supercell.
- bulk_entry: ComputedEntry | None = None
- property bulk_entry_energy
Get the raw energy of the bulk supercell calculation (i.e.
bulk_entry.energy
).Refactored from
pymatgen-analysis-defects
to be more efficient, reducing compute times when looping over formation energy / concentration functions.
- property bulk_site_concentration
Return the site concentration (in cm^-3) of the corresponding atomic site of the defect in the pristine bulk material (e.g. if the defect is V_O in SrTiO3, returns the site concentration of (symmetry-equivalent) oxygen atoms in SrTiO3).
- bulk_supercell: Structure | None = None
- calculation_metadata: dict
- charge_state: int
- charge_state_guessing_log: dict
- conv_cell_frac_coords: ndarray | None = None
- conventional_structure: Structure | None = None
- property corrected_energy: float
Get the energy of the defect supercell calculation, with all corrections (in
DefectEntry.corrections
) applied.Refactored from
pymatgen-analysis-defects
to be more efficient, reducing compute times when looping over formation energy / concentration functions.
- corrections: dict[str, float]
- corrections_metadata: dict[str, Any]
- defect_supercell: Structure | None = None
- defect_supercell_site: PeriodicSite | None = None
- degeneracy_factors: dict
- entry_id: str | None = None
- equilibrium_concentration(chempots: dict | None = None, limit: str | None = None, el_refs: dict | None = None, temperature: float = 300, fermi_level: float = 0, vbm: float | None = None, per_site: bool = False, symprec: float | None = None, formation_energy: float | None = None) float [source]
Compute the equilibrium concentration (in cm^-3) for the
DefectEntry
at a given chemical potential limit, fermi_level and temperature, assuming the dilute limit approximation.Note that these are the equilibrium defect concentrations!
DefectThermodynamics.get_quenched_fermi_level_and_concentrations()
can instead be used to calculate the Fermi level and defect concentrations for a material grown/annealed at higher temperatures and then cooled (quenched) to room/operating temperature (where defect concentrations are assumed to remain fixed) - this is known as the frozen defect approach and is typically the most valid approximation (see its docstring for more information, and discussion in 10.1039/D3CS00432E).The degeneracy/multiplicity factor “g” is an important parameter in the defect concentration equation (see discussion in https://doi.org/10.1039/D2FD00043A and https://doi.org/10.1039/D3CS00432E), affecting the final concentration by up to 2 orders of magnitude. This factor is taken from the product of the
defect_entry.defect.multiplicity
anddefect_entry.degeneracy_factors
attributes.- Parameters:
chempots (dict) –
Dictionary of chemical potentials to use for calculating the defect formation energy (and thus concentration). This can have the form of:
{"limits": [{'limit': [chempot_dict]}]}
(the format generated bydoped
's chemical potential parsing functions (see tutorials)) and specific limits (chemical potential limits) can then be chosen usinglimit
.Alternatively this can be a dictionary of chemical potentials for a single limit (limit), in the format:
{element symbol: chemical potential}
. If manually specifying chemical potentials this way, you can set theel_refs
option with the DFT reference energies of the elemental phases, in which case it is the formal chemical potentials (i.e. relative to the elemental references) that should be given here, otherwise the absolute (DFT) chemical potentials should be given.If None (default), sets all chemical potentials to 0 (inaccurate formation energies and concentrations!). (Default: None)
limit (str) –
The chemical potential limit for which to calculate the formation energy and thus concentration. Can be either:
None (default), if
chempots
corresponds to a single chemical potential limit - otherwise will use the first chemical potential limit in thechempots
dict.”X-rich”/”X-poor” where X is an element in the system, in which case the most X-rich/poor limit will be used (e.g. “Li-rich”).
A key in the
(self.)chempots["limits"]
dictionary.
The latter two options can only be used if
chempots
is in thedoped
format (see chemical potentials tutorial). (Default: None)el_refs (dict) – Dictionary of elemental reference energies for the chemical potentials in the format:
{element symbol: reference energy}
(to determine the formal chemical potentials, whenchempots
has been manually specified as{element symbol: chemical potential}
). Unnecessary ifchempots
is provided/present in format generated bydoped
(see tutorials). (Default: None)temperature (float) – Temperature in Kelvin at which to calculate the equilibrium concentration.
vbm (float) – VBM eigenvalue to use as Fermi level reference point for calculating the formation energy. If
None
(default), will use"vbm"
from thecalculation_metadata
dict attribute if present – which corresponds to the VBM of the bulk supercell calculation by default, unlessbulk_band_gap_vr
is set during defect parsing).fermi_level (float) – Value corresponding to the electron chemical potential, referenced to the VBM. Default is 0 (i.e. the VBM).
per_site (bool) – Whether to return the concentration as fractional concentration per site, rather than the default of per cm^3. Multiply by 100 for concentration in percent. (default: False)
symprec (float) – Symmetry tolerance for
spglib
to use when determining relaxed defect point symmetries and thus orientational degeneracies. Default is0.1
which matches that used by theMaterials Project
and is larger than thepymatgen
default of0.01
(which is used bydoped
for unrelaxed/bulk structures) to account for residual structural noise in relaxed defect supercells. You may want to adjust for your system (e.g. if there are very slight octahedral distortions etc.). Ifsymprec
is set, then the point symmetries and corresponding orientational degeneracy will be re-parsed/computed even if already present in theDefectEntry
objectcalculation_metadata
.formation_energy (float) – Pre-calculated formation energy to use for the defect concentration calculation, in order to reduce compute times (e.g. when looping over many chemical potential / temperature / etc ranges). Only really intended for internal
doped
usage. IfNone
(default), will calculate the formation energy using the inputchempots
,limit
,el_refs
,vbm
andfermi_level
arguments. (Default: None)
- Returns:
Concentration in cm^-3 (or as fractional per site, if per_site = True) (float)
- equiv_conv_cell_frac_coords: list[ndarray]
- equivalent_supercell_sites: list[PeriodicSite]
- formation_energy(chempots: dict | None = None, limit: str | None = None, el_refs: dict | None = None, vbm: float | None = None, fermi_level: float = 0) float [source]
Compute the formation energy for the
DefectEntry
at a given chemical potential limit and fermi_level.- Parameters:
chempots (dict) –
Dictionary of chemical potentials to use for calculating the defect formation energy. This can have the form of:
{"limits": [{'limit': [chempot_dict]}]}
(the format generated bydoped
's chemical potential parsing functions (see tutorials)) and specific limits (chemical potential limits) can then be chosen usinglimit
.Alternatively this can be a dictionary of chemical potentials for a single limit (limit), in the format:
{element symbol: chemical potential}
. If manually specifying chemical potentials this way, you can set theel_refs
option with the DFT reference energies of the elemental phases, in which case it is the formal chemical potentials (i.e. relative to the elemental references) that should be given here, otherwise the absolute (DFT) chemical potentials should be given.If None (default), sets all chemical potentials to zero. (Default: None)
limit (str) –
The chemical potential limit for which to calculate the formation energy. Can be either:
None (default), if
chempots
corresponds to a single chemical potential limit - otherwise will use the first chemical potential limit in thechempots
dict.”X-rich”/”X-poor” where X is an element in the system, in which case the most X-rich/poor limit will be used (e.g. “Li-rich”).
A key in the
(self.)chempots["limits"]
dictionary.
The latter two options can only be used if
chempots
is in thedoped
format (see chemical potentials tutorial). (Default: None)el_refs (dict) – Dictionary of elemental reference energies for the chemical potentials in the format:
{element symbol: reference energy}
(to determine the formal chemical potentials, whenchempots
has been manually specified as{element symbol: chemical potential}
). Unnecessary ifchempots
is provided/present in format generated bydoped
(see tutorials). (Default: None)vbm (float) – VBM eigenvalue to use as Fermi level reference point for calculating formation energy. If
None
(default), will use"vbm"
from thecalculation_metadata
dict attribute if present – which corresponds to the VBM of the bulk supercell calculation by default, unlessbulk_band_gap_vr
is set during defect parsing).fermi_level (float) – Value corresponding to the electron chemical potential, referenced to the VBM eigenvalue. Default is 0 (i.e. the VBM).
- Returns:
Formation energy value (float)
- classmethod from_dict(d: dict)[source]
Class method to create a
DefectEntry
object from a dictionary.Defined to avoid unnecessary
vise
/pydefect
INFO messages.- Parameters:
d (dict) – Dictionary representation of the
DefectEntry
object.- Returns:
DefectEntry
object
- classmethod from_json(filename: str)[source]
Load a
DefectEntry
object from a json file.- Parameters:
filename (str) – Filename of json file to load
DefectEntry
from.- Returns:
DefectEntry
object
- get_ediff() float [source]
Get the energy difference between the defect and bulk supercell calculations, including finite-size corrections.
Refactored from
pymatgen-analysis-defects
to be more efficient, reducing compute times when looping over formation energy / concentration functions.
- get_eigenvalue_analysis(plot: bool = True, filename: str | None = None, bulk_vr: str | Path | Vasprun | None = None, bulk_procar: str | Path | EasyunfoldProcar | Procar | None = None, defect_vr: str | Path | Vasprun | None = None, defect_procar: str | Path | EasyunfoldProcar | Procar | None = None, force_reparse: bool = False, **kwargs)[source]
Returns information about the band edge and in-gap electronic states and their orbital character / localisation degree for the defect entry, as well as a plot of the single-particle electronic eigenvalues and their occupation (if
plot=True
).Can be used to determine if a defect is adopting a perturbed host state (PHS / shallow state), see https://doped.readthedocs.io/en/latest/Tips.html#perturbed-host-states.
If eigenvalue data has not already been parsed for
DefectEntry
(default indoped
is to parse this data withDefectsParser
/DefectParser
, as controlled by theparse_projected_eigen
flag), then this function will attempt to load the eigenvalue data from either the inputVasprun
/PROCAR
objects or files, or from thebulk/defect_path
s indefect_entry.calculation_metadata
. If so, will initially try to load orbital projections fromvasprun.xml(.gz)
files (slightly slower but more accurate), or failing that fromPROCAR(.gz)
files if present.This function uses code from
pydefect
: Citation: https://doi.org/10.1103/PhysRevMaterials.5.123803.- Parameters:
plot (bool) – Whether to plot the single-particle eigenvalues. (Default: True)
filename (str) – Filename to save the eigenvalue plot to (if
plot = True
). IfNone
(default), plots are not saved.bulk_vr (str, Path, Vasprun) – Not required if eigenvalue data has already been parsed for
DefectEntry
(default behaviour when parsing withdoped
, data indefect_entry.calculation_metadata["eigenvalue_data"]
). Either a path to theVASP
vasprun.xml(.gz)
output file or apymatgen
Vasprun
object, for the reference bulk supercell calculation. IfNone
(default), tries to load theVasprun
object fromself.calculation_metadata["run_metadata"]["bulk_vasprun_dict"]
, or, failing that, from avasprun.xml(.gz)
file atself.calculation_metadata["bulk_path"]
.bulk_procar (str, Path, EasyunfoldProcar, Procar) – Not required if eigenvalue data has already been parsed for
DefectEntry
(default behaviour when parsing withdoped
, data indefect_entry.calculation_metadata["eigenvalue_data"]
), or ifbulk_vr
was parsed withparse_projected_eigen = True
. Either a path to theVASP
PROCAR
output file (withLORBIT > 10
in theINCAR
) or aneasyunfold
/pymatgen
Procar
object, for the reference bulk supercell calculation. IfNone
(default), tries to load from aPROCAR(.gz)
file atself.calculation_metadata["bulk_path"]
.defect_vr (str, Path, Vasprun) – Not required if eigenvalue data has already been parsed for
DefectEntry
(default behaviour when parsing withdoped
, data indefect_entry.calculation_metadata["eigenvalue_data"]
). Either a path to theVASP
vasprun.xml(.gz)
output file or apymatgen
Vasprun
object, for the defect supercell calculation. IfNone
(default), tries to load theVasprun
object fromself.calculation_metadata["run_metadata"]["defect_vasprun_dict"]
, or, failing that, from avasprun.xml(.gz)
file atself.calculation_metadata["defect_path"]
.defect_procar (str, Path, EasyunfoldProcar, Procar) – Not required if eigenvalue data has already been parsed for
DefectEntry
(default behaviour when parsing withdoped
, data indefect_entry.calculation_metadata["eigenvalue_data"]
), or ifdefect_vr
was parsed withparse_projected_eigen = True
. Either a path to theVASP
PROCAR
output file (withLORBIT > 10
in theINCAR
) or aneasyunfold
/pymatgen
Procar
object, for the defect supercell calculation. IfNone
(default), tries to load from aPROCAR(.gz)
file atself.calculation_metadata["defect_path"]
.force_reparse (bool) – Whether to force re-parsing of the eigenvalue data, even if already present in the
calculation_metadata
.**kwargs – Additional kwargs to pass to
doped.utils.eigenvalues.get_eigenvalue_analysis
, such asstyle_file
,ks_levels
,ylims
,legend_kwargs
,similar_orb_criterion
,similar_energy_criterion
.
- Returns:
pydefect
BandEdgeStates
object andmatplotlib
Figure
object (ifplot=True
).
- get_freysoldt_correction(dielectric: float | int | ndarray | list | None = None, defect_locpot: str | Locpot | dict | None = None, bulk_locpot: str | Locpot | dict | None = None, plot: bool = False, filename: str | None = None, axis=None, return_correction_error: bool = False, error_tolerance: float = 0.05, style_file: str | None = None, **kwargs) CorrectionResult [source]
Compute the isotropic Freysoldt (FNV) correction for the defect_entry.
The correction is added to the
defect_entry.corrections
dictionary (to be used in following formation energy calculations). If this correction is used, please cite Freysoldt’s original paper; 10.1103/PhysRevLett.102.016402.The charge correction error is estimated by computing the average standard deviation of the planar-averaged potential difference in the sampling region, and multiplying by the defect charge. This is expected to be a lower bound estimate of the true charge correction error.
- Parameters:
dielectric (float or int or 3x1 matrix or 3x3 matrix) – Total dielectric constant of the host compound (including both ionic and (high-frequency) electronic contributions), in the same xyz Cartesian basis as the supercell calculations. If None, then the dielectric constant is taken from the
defect_entry
calculation_metadata
if available.defect_locpot – Path to the output VASP
LOCPOT
file from the defect supercell calculation, or the correspondingpymatgen
Locpot
object, or a dictionary of the planar-averaged potential in the form:{i: Locpot.get_average_along_axis(i) for i in [0,1,2]}
. If None, will try to usedefect_locpot
from thedefect_entry
calculation_metadata
if available.bulk_locpot – Path to the output VASP
LOCPOT
file from the bulk supercell calculation, or the correspondingpymatgen
Locpot
object, or a dictionary of the planar-averaged potential in the form:{i: Locpot.get_average_along_axis(i) for i in [0,1,2]}
. If None, will try to usebulk_locpot
from thedefect_entry
calculation_metadata
if available.plot (bool) – Whether to plot the FNV electrostatic potential plots (for manually checking the behaviour of the charge correction here).
filename (str) – Filename to save the FNV electrostatic potential plots to. If None, plots are not saved.
axis (int or None) – If int, then the FNV electrostatic potential plot along the specified axis (0, 1, 2 for a, b, c) will be plotted. Note that the output charge correction is still that for all axes. If None, then all three axes are plotted.
return_correction_error (bool) – If True, also returns the average standard deviation of the planar-averaged potential difference times the defect charge (which gives an estimate of the error range of the correction energy). Default is False.
error_tolerance (float) – If the estimated error in the charge correction, based on the variance of the potential in the sampling region, is greater than this value (in eV), then a warning is raised. (default: 0.05 eV)
style_file (str) – Path to a
.mplstyle
file to use for the plot. IfNone
(default), uses the default doped style (fromdoped/utils/doped.mplstyle
).**kwargs – Additional kwargs to pass to
pymatgen.analysis.defects.corrections.freysoldt.get_freysoldt_correction
(e.g. energy_cutoff, mad_tol, q_model, step).
- Returns:
CorrectionResults (summary of the corrections applied and metadata), and the matplotlib figure object (or axis object if axis specified) if
plot
is True, and the estimated charge correction error ifreturn_correction_error
is True.
- get_kumagai_correction(dielectric: float | int | ndarray | list | None = None, defect_region_radius: float | None = None, excluded_indices: list[int] | None = None, defect_outcar: str | Outcar | None = None, bulk_outcar: str | Outcar | None = None, plot: bool = False, filename: str | None = None, return_correction_error: bool = False, error_tolerance: float = 0.05, style_file: str | None = None, **kwargs)[source]
Compute the Kumagai (eFNV) finite-size charge correction for the defect_entry. Compatible with both isotropic/cubic and anisotropic systems.
The correction is added to the
defect_entry.corrections
dictionary (to be used in following formation energy calculations). If this correction is used, please cite the Kumagai & Oba paper: 10.1103/PhysRevB.89.195205Typically for reasonably well-converged supercell sizes, the default
defect_region_radius
works perfectly well. However, for certain materials at small/intermediate supercell sizes, you may want to adjust this (and/orexcluded_indices
) to ensure the best sampling of the plateau region away from the defect position -doped
should throw a warning in these cases (about the correction error being above the default tolerance (50 meV)). For example, with layered materials, the defect charge is often localised to one layer, so we may want to adjustdefect_region_radius
and/orexcluded_indices
to ensure that only sites in other layers are used for the sampling region (plateau) - see example on doped docsTips
page.The correction error is estimated by computing the standard error of the mean of the sampled site potential differences, multiplied by the defect charge. This is expected to be a lower bound estimate of the true charge correction error.
- Parameters:
dielectric (float or int or 3x1 matrix or 3x3 matrix) – Total dielectric constant of the host compound (including both ionic and (high-frequency) electronic contributions), in the same xyz Cartesian basis as the supercell calculations. If None, then the dielectric constant is taken from the
defect_entry
calculation_metadata
if available.defect_region_radius (float) – Radius of the defect region (in Å). Sites outside the defect region are used for sampling the electrostatic potential far from the defect (to obtain the potential alignment). If None (default), uses the Wigner-Seitz radius of the supercell.
excluded_indices (list) – List of site indices (in the defect supercell) to exclude from the site potential sampling in the correction calculation/plot. If None (default), no sites are excluded.
defect_outcar (str or Outcar) – Path to the output VASP OUTCAR file from the defect supercell calculation, or the corresponding
pymatgen
Outcar object. If None, will try to use thedefect_supercell_site_potentials
from thedefect_entry
calculation_metadata
if available.bulk_outcar (str or Outcar) – Path to the output VASP OUTCAR file from the bulk supercell calculation, or the corresponding
pymatgen
Outcar object. If None, will try to use thebulk_supercell_site_potentials
from thedefect_entry
calculation_metadata
if available.plot (bool) – Whether to plot the Kumagai site potential plots (for manually checking the behaviour of the charge correction here).
filename (str) – Filename to save the Kumagai site potential plots to. If None, plots are not saved.
return_correction_error (bool) – If True, also returns the standard error of the mean of the sampled site potential differences times the defect charge (which gives an estimate of the error range of the correction energy). Default is False.
error_tolerance (float) – If the estimated error in the charge correction, based on the variance of the potential in the sampling region, is greater than this value (in eV), then a warning is raised. (default: 0.05 eV)
style_file (str) – Path to a
.mplstyle
file to use for the plot. IfNone
(default), uses the default doped style (fromdoped/utils/doped.mplstyle
).**kwargs – Additional kwargs to pass to pydefect.corrections.efnv_correction.ExtendedFnvCorrection (e.g. charge, defect_region_radius, defect_coords).
- Returns:
CorrectionResults (summary of the corrections applied and metadata), and the matplotlib figure object if
plot
is True, and the estimated charge correction error ifreturn_correction_error
is True.
- name: str = ''
- plot_site_displacements(separated_by_direction: bool | None = False, relative_to_defect: bool | None = False, vector_to_project_on: list | None = None, use_plotly: bool | None = False, mpl_style_file: str | None = '')[source]
Plot the site displacements as a function of distance from the defect site.
- Parameters:
separated_by_direction (bool) – Whether to plot the site displacements separated by the x, y and z directions (True) or all together (False). Defaults to False.
relative_to_defect (bool) – Whether to plot the signed displacements along the line from the defect site to that atom. Negative values indicate the atom moves towards the defect (compressive strain), positive values indicate the atom moves away from the defect (tensile strain). Uses the relaxed defect position as reference.
vector_to_project_on – Direction to project the site displacements along (e.g. [0, 0, 1]). Defaults to None (e.g. the displacements are calculated in the cartesian basis x, y, z).
use_plotly (bool) – Whether to use plotly (True) or matplotlib (False).
mpl_style_file (str) – Path to a matplotlib style file to use for the plot. If None, uses the default doped style file.
- sc_defect_frac_coords: tuple[float, float, float] | None = None
- sc_entry: ComputedStructureEntry
- property sc_entry_energy
Get the raw energy of the defect supercell calculation (i.e.
sc_entry.energy
).Refactored from
pymatgen-analysis-defects
to be more efficient, reducing compute times when looping over formation energy / concentration functions.
- to_json(filename: str | None = None)[source]
Save the
DefectEntry
object to a json file, which can be reloaded with theDefectEntry.from_json()
class method.- Parameters:
filename (str) – Filename to save json file as. If None, the filename will be set as
"{DefectEntry.name}.json"
.
- wyckoff: str | None = None
- class doped.core.Interstitial(*args, **kwargs)[source]
Bases:
Defect
,Interstitial
Subclass of
pymatgen.analysis.defects.core.Interstitial
with additional attributes and methods used bydoped
.
- class doped.core.Substitution(*args, **kwargs)[source]
Bases:
Defect
,Substitution
Subclass of
pymatgen.analysis.defects.core.Substitution
with additional attributes and methods used bydoped
.
- class doped.core.Vacancy(*args, **kwargs)[source]
Bases:
Defect
,Vacancy
Subclass of
pymatgen.analysis.defects.core.Vacancy
with additional attributes and methods used bydoped
.
- doped.core.doped_defect_from_pmg_defect(defect: Defect, bulk_oxi_states=False, **doped_kwargs)[source]
Create the corresponding
doped
Defect
(Vacancy
,Interstitial
,Substitution
) from an inputpymatgen
Defect
object.- Parameters:
defect –
pymatgen
Defect
object.bulk_oxi_states – Either a dict of bulk oxidation states to use, or a boolean. If True, re-guesses the oxidation state of the defect (ignoring the
pymatgen
Defect
oxi_state
attribute), otherwise uses the already-setoxi_state
(default = 0). Used in doped defect generation to make defect setup more robust and efficient (particularly for odd input structures, such as defect supercells etc).**doped_kwargs – Additional keyword arguments to define doped-specific attributes (see class docstring).
- doped.core.guess_and_set_oxi_states_with_timeout(structure, timeout_1=10, timeout_2=15, break_early_if_expensive=False) bool [source]
Tries to guess (and set) the oxidation states of the input structure, with a timeout catch for cases where the structure is complex and oxi state guessing will take a very very long time.
Tries first without using the
pymatgen
BVAnalyzer
class, and if this fails, tries using the ICSD oxidation state probabilities (with timeouts) to guess.- Parameters:
structure (Structure) – The structure for which to guess the oxidation states.
timeout_1 (float) – Timeout in seconds for the second attempt to guess the oxidation states, using ICSD oxidation state probabilities (with
max_sites=-1
). Default is 10 seconds.timeout_2 (float) – Timeout in seconds for the third attempt to guess the oxidation states, using ICSD oxidation state probabilities (without
max_sites=-1
). Default is 15 seconds.break_early_if_expensive (bool) – Whether to stop the function if the first oxi state guessing attempt (with
BVAnalyzer
) fails and the cost estimate for the ICSD probability guessing is high (expected to take a long time; > 10 seconds). Default isFalse
.
- Returns:
The structure with oxidation states guessed and set, or
False
if oxidation states could not be guessed.- Return type:
Structure
- doped.core.guess_and_set_struct_oxi_states(structure, try_without_max_sites=False)[source]
Tries to guess (and set) the oxidation states of the input structure, first using the
pymatgen
BVAnalyzer
class, and if that fails, using the ICSD oxidation state probabilities to guess.- Parameters:
structure (Structure) – The structure for which to guess the oxidation states.
try_without_max_sites (bool) – Whether to try to guess the oxidation states without using the
max_sites=-1
argument (True
)(which attempts to use the reduced composition for guessing oxi states) or not (False
), when using the ICSD oxidation state probability guessing.
- Returns:
The structure with oxidation states guessed and set, or
False
if oxidation states could not be guessed.- Return type:
Structure