pymatgen.analysis.defects.corrections module

class BandEdgeShiftingCorrection[source]

Bases: pymatgen.analysis.defects.core.DefectCorrection

A class for BandEdgeShiftingCorrection class. Largely adapted from PyCDT code

Initializes the BandEdgeShiftingCorrection class

get_correction(entry)[source]

Gets the BandEdge correction for a defect entry :param entry: defect entry to compute BandFilling correction on.

Requires some parameters in the DefectEntry to properly function:
hybrid_cbm (float)

CBM of HYBRID bulk calculation one wishes to shift to

hybrid_vbm (float)

VBM of HYBRID bulk calculation one wishes to shift to

cbm (float)

CBM of bulk calculation (or band structure calculation of bulk); calculated on same level of theory as the defect (ex. GGA defects -> requires GGA cbm)

vbm (float)

VBM of bulk calculation (or band structure calculation of bulk); calculated on same level of theory as the defect (ex. GGA defects -> requires GGA vbm)

Returns

BandfillingCorrection value as a dictionary

class BandFillingCorrection(resolution=0.01)[source]

Bases: pymatgen.analysis.defects.core.DefectCorrection

A class for BandFillingCorrection class. Largely adapted from PyCDT code

Initializes the Bandfilling correction

Parameters

resolution (float) – energy resolution to maintain for gap states

get_correction(entry)[source]

Gets the BandFilling correction for a defect entry :param entry: defect entry to compute BandFilling correction on.

Requires following parameters in the DefectEntry to exist:
eigenvalues

dictionary of defect eigenvalues, as stored in a Vasprun object

kpoint_weights (list of floats)

kpoint weights corresponding to the dictionary of eigenvalues, as stored in a Vasprun object

potalign (float)

potential alignment for the defect calculation Only applies to non-zero charge, When using potential alignment correction (freysoldt or kumagai), need to divide by -q

cbm (float)

CBM of bulk calculation (or band structure calculation of bulk); calculated on same level of theory as the defect (ex. GGA defects -> requires GGA cbm)

vbm (float)

VBM of bulk calculation (or band structure calculation of bulk); calculated on same level of theory as the defect (ex. GGA defects -> requires GGA vbm)

Returns

Bandfilling Correction value as a dictionary

perform_bandfill_corr(eigenvalues, kpoint_weights, potalign, vbm, cbm)[source]

This calculates the band filling correction based on excess of electrons/holes in CB/VB…

Note that the total free holes and electrons may also be used for a “shallow donor/acceptor”
correction with specified band shifts:

+num_elec_cbm * Delta E_CBM (or -num_hole_vbm * Delta E_VBM)

class FreysoldtCorrection(dielectric_const, q_model=None, energy_cutoff=520, madetol=0.0001, axis=None)[source]

Bases: pymatgen.analysis.defects.core.DefectCorrection

A class for FreysoldtCorrection class. Largely adapated from PyCDT code

If this correction is used, please reference Freysoldt’s original paper. doi: 10.1103/PhysRevLett.102.016402

Initializes the FreysoldtCorrection class :param dielectric_const: Dielectric constant for the structure :type dielectric_const: float or 3x3 matrix :param q_model: instantiated QModel object or None.

Uses default parameters to instantiate QModel if None supplied

Parameters
  • energy_cutoff (int) – Maximum energy in eV in reciprocal space to perform integration for potential correction.

  • madeltol (float) – Convergence criteria for the Madelung energy for potential correction

  • axis (int) – Axis to calculate correction. If axis is None, then averages over all three axes is performed.

get_correction(entry)[source]

Gets the Freysoldt correction for a defect entry :param entry: defect entry to compute Freysoldt correction on.

Requires following keys to exist in DefectEntry.parameters dict:

axis_grid (3 x NGX where NGX is the length of the NGX grid in the x,y and z axis directions. Same length as planar average lists):

A list of 3 numpy arrays which contain the cartesian axis values (in angstroms) that correspond to each planar avg potential supplied.

bulk_planar_averages (3 x NGX where NGX is the length of the NGX grid in the x,y and z axis directions.):

A list of 3 numpy arrays which contain the planar averaged electrostatic potential for the bulk supercell.

defect_planar_averages (3 x NGX where NGX is the length of the NGX grid in the x,y and z axis directions.):

A list of 3 numpy arrays which contain the planar averaged electrostatic potential for the defective supercell.

initial_defect_structure (Structure) structure corresponding to

initial defect supercell structure (uses Lattice for charge correction)

defect_frac_sc_coords (3 x 1 array) Fractional co-ordinates of

defect location in supercell structure

Returns

FreysoldtCorrection values as a dictionary

perform_es_corr(lattice, q, step=0.0001)[source]

Peform Electrostatic Freysoldt Correction :param lattice: Pymatgen lattice object :param q: Charge of defect :type q: int :param step: step size for numerical integration :type step: float

Returns

Electrostatic Point Charge contribution to Freysoldt Correction (float)

perform_pot_corr(axis_grid, pureavg, defavg, lattice, q, defect_frac_position, axis, widthsample=1.0)[source]

For performing planar averaging potential alignment :param axis_grid (1 x NGX where NGX is the length of the NGX grid:

in the axis direction. Same length as pureavg list):

A numpy array which contain the cartesian axis values (in angstroms) that correspond to each planar avg potential supplied.

Parameters
  • (1 x NGX where NGX is the length of the NGX grid in (defavg) –

    the axis direction.):

    A numpy array for the planar averaged electrostatic potential of the bulk supercell.

  • (1 x NGX where NGX is the length of the NGX grid in

    the axis direction.):

    A numpy array for the planar averaged electrostatic potential of the defect supercell.

  • lattice – Pymatgen Lattice object of the defect supercell

  • q (float or int) – charge of the defect

  • defect_frac_position – Fracitional Coordinates of the defect in the supercell

  • axis (int) – axis for performing the freysoldt correction on

  • widthsample (float) – width (in Angstroms) of the region in between defects where the potential alignment correction is averaged. Default is 1 Angstrom.

Returns

Potential Alignment contribution to Freysoldt Correction (float)

plot(axis, title=None, saved=False)[source]

Plots the planar average electrostatic potential against the Long range and short range models from Freysoldt. Must run perform_pot_corr or get_correction (to load metadata) before this can be used. :param axis: axis to plot :type axis: int :param title: Title to be given to plot. Default is no title. :type title: str :param saved: Whether to save file or not. If False then returns plot

object. If True then saves plot as str(title) + “FreyplnravgPlot.pdf”

class KumagaiCorrection(dielectric_tensor, sampling_radius=None, gamma=None)[source]

Bases: pymatgen.analysis.defects.core.DefectCorrection

A class for KumagaiCorrection class. Largely adapated from PyCDT code

If this correction is used, please reference Kumagai and Oba’s original paper (doi: 10.1103/PhysRevB.89.195205) as well as Freysoldt’s original paper (doi: 10.1103/PhysRevLett.102.016402)

NOTE that equations 8 and 9 from Kumagai et al. reference are divided by (4 pi) to get SI units

Initializes the Kumagai Correction :param dielectric_tensor: Dielectric constant for the structure :type dielectric_tensor: float or 3x3 matrix :param optional data that can be tuned:

sampling_radius (float): radius (in Angstrom) which sites must be outside

of to be included in the correction. Publication by Kumagai advises to use Wigner-Seitz radius of defect supercell, so this is default value.

gamma (float): convergence parameter for gamma function.

Code will automatically determine this if set to None.

get_correction(entry)[source]

Gets the Kumagai correction for a defect entry :param entry: defect entry to compute Kumagai correction on.

Requires following parameters in the DefectEntry to exist:

bulk_atomic_site_averages (list): list of bulk structure”s atomic site averaged ESPs * charge,

in same order as indices of bulk structure note this is list given by VASP’s OUTCAR (so it is multiplied by a test charge of -1)

defect_atomic_site_averages (list): list of defect structure”s atomic site averaged ESPs * charge,

in same order as indices of defect structure note this is list given by VASP’s OUTCAR (so it is multiplied by a test charge of -1)

site_matching_indices (list): list of corresponding site index values for

bulk and defect site structures EXCLUDING the defect site itself (ex. [[bulk structure site index, defect structure”s corresponding site index], … ]

initial_defect_structure (Structure): Pymatgen Structure object representing un-relaxed defect

structure

defect_frac_sc_coords (array): Defect Position in fractional coordinates of the supercell

given in bulk_structure

Returns

KumagaiCorrection values as a dictionary

get_potential_shift(gamma, volume)[source]
get_real_summation(gamma, real_vectors)[source]

Get real summation term from list of real-space vectors

get_recip_summation(gamma, recip_vectors, volume, r=[0.0, 0.0, 0.0])[source]

Get Reciprocal summation term from list of reciprocal-space vectors

get_self_interaction(gamma)[source]
perform_es_corr(gamma, prec, lattice, charge)[source]

Peform Electrostatic Kumagai Correction :param gamma: Ewald parameter :type gamma: float :param prec: Precision parameter for reciprical/real lattice vector generation :type prec: int :param lattice: Pymatgen Lattice object corresponding to defect supercell :param charge: Defect charge :type charge: int

Returns

Electrostatic Point Charge contribution to Kumagai Correction (float)

perform_pot_corr(defect_structure, defect_frac_coords, site_list, sampling_radius, q, r_vecs, g_vecs, gamma)[source]

For performing potential alignment in manner described by Kumagai et al. :param defect_structure: Pymatgen Structure object corrsponding to the defect supercell :param defect_frac_coords: Defect Position in fractional coordinates of the supercell

given in bulk_structure

Parameters
  • site_list – list of corresponding site index values for bulk and defect site structures EXCLUDING the defect site itself (ex. [[bulk structure site index, defect structure”s corresponding site index], … ]

  • sampling_radius (float) – radius (in Angstrom) which sites must be outside of to be included in the correction. Publication by Kumagai advises to use Wigner-Seitz radius of defect supercell, so this is default value.

  • q (int) – Defect charge

  • r_vecs – List of real lattice vectors to use in summation

  • g_vecs – List of reciprocal lattice vectors to use in summation

  • gamma (float) – Ewald parameter

Returns

Potential alignment contribution to Kumagai Correction (float)

plot(title=None, saved=False)[source]

Plots the AtomicSite electrostatic potential against the Long range and short range models from Kumagai and Oba (doi: 10.1103/PhysRevB.89.195205)