pymatgen.analysis.defects.point_defects module¶
-
class
Defect[source]¶ Bases:
objectAbstract class for point defects
-
get_coordinated_elements(n)[source]¶ Elements of sites in structure surrounding the defect site.
Parameters: n – Index of defect list
-
get_coordinated_sites(n)[source]¶ The sites in structure surrounding the defect site.
Parameters: n – Index of defects list
-
get_defectsite_coordination_number(n)[source]¶ Coordination number of interstitial site.
Parameters: n – Index of interstitial list
-
get_defectsite_multiplicity(n)[source]¶ Returns the symmtric multiplicity of the defect site at the index.
-
make_supercells_with_defects(scaling_matrix)[source]¶ Generate the supercell with input multipliers and create the defect. First supercell has no defects. To create unit cell with defect pass unit matrix.
-
struct_radii¶ Radii of elements in the structure
-
struct_valences¶ Valence of elements in the structure
-
structure¶ Returns the structure without any defects Useful for Mott-Littleton calculations.
-
-
class
Interstitial(structure, valences, radii, site_type='voronoi_vertex', accuracy='Normal', symmetry_flag=True, oxi_state=False)[source]¶ Bases:
pymatgen.analysis.defects.point_defects.DefectSubclass of Defect to generate interstitial sites
Given a structure, generate symmetrically distinct interstitial sites. For a non-ionic structure, use oxi_state=True and give atomic radii.
Parameters: - structure – pymatgen.core.structure.Structure
- valences – Dictionary of oxidation states of elements in {el:valence} form
- radii – Radii of elemnts in the structure
- site_type – “voronoi_vertex” uses voronoi nodes “voronoi_edgecenter” uses voronoi polyhedra edge centers “voronoi_facecenter” uses voronoi polyhedra face centers “all” combines vertices, edgecenters and facecenters. Default is “voronoi_vertex”
- accuracy – Flag denoting whether to use high accuracy version of Zeo++. Options are “Normal” and “High”. Default is normal.
- symmetry_flag – If True, only returns symmetrically distinct sites
- oxi_state – If False, input structure is considered devoid of oxidation-state decoration. And oxi-state for each site is determined. Use True, if input structure is oxi-state decorated. This option is useful when the structure is not electro-neutral after deleting/adding sites. In that case oxi-decorate the structure before deleting/adding the sites.
-
append_defectsite(site)[source]¶ Append a site to list of possible interstitials
Parameters: site – pymatgen.core.sites.Site
-
delete_defectsite(n)[source]¶ Remove a symmetrically distinct interstitial site
Parameters: n – Index of interstitial site
-
enumerate_defectsites()[source]¶ Enumerate all the symmetrically distinct interstitial sites. The defect site has “X” as occupied specie.
-
get_coordsites_charge_sum(n)[source]¶ Total charge of the interstitial coordinated sites.
Parameters: n – Index of interstitial list
-
get_coordsites_min_max_charge(n)[source]¶ Minimum and maximum charge of sites surrounding the interstitial site.
Parameters: n – Index of symmetrical distinct interstitial site
-
get_radius(n)[source]¶ Volume of the nth interstitial
Parameters: n – Index of symmetrically distinct vacancies list Returns: floating number representing radius of interstitial sphere
-
make_supercells_with_defects(scaling_matrix, element)[source]¶ Returns sequence of supercells in pymatgen.core.structure.Structure format, with each supercell containing an interstitial. First supercell has no defects.
-
prune_defectsites(el='C', oxi_state=4, dlta=0.1)[source]¶ Prune all the defect sites which can’t acoomodate the input elment with the input oxidation state.
-
class
InterstitialAnalyzer(inter, el, oxi_state, scd=2)[source]¶ Bases:
objectUse GULP to compute the interstitial formation energy, relaxed structures. Works only for metal oxides due to the use of Buckingham Potentials.
Parameters: - inter – pymatgen.defects.point_defects.Interstitial
- el – Element name in short hand notation (“El”)
- oxi_state – Oxidtation state
- scd – Super cell dimension as number. The scaling is equal along xyz.
-
get_percentage_bond_distance_change(n)[source]¶ Bond distance change after the introduction of interstitial
Parameters: n – Symmetrically distinct interstitial index
-
get_percentage_lattice_parameter_change(n)[source]¶ Lattice parameter change after the introduction of interstitial
Parameters: n – Symmetrically distinct interstitial index
-
get_percentage_volume_change(n)[source]¶ Volume change after the introduction of interstitial
Parameters: n – Symmetrically distinct interstitial index
-
class
InterstitialStructureRelaxer(interstitial, el, oxi_state, supercell_dim=2)[source]¶ Bases:
objectPerforms structural relaxation for each interstitial supercell.
Parameters: - interstitial – Unrelaxed interstitial
- el – Species string in short notation
- oxi_state – Oxidation state of the element
- supercell_dim – Dimensions of super cell
-
get_relaxed_energy(n)[source]¶ Get the relaxed structure of nth symmetrically distinct interstitial.
Parameters: n – Symmetrically distinct interstitial index Note
0 corresponds to relaxed bulk structure
-
get_relaxed_interstitial()[source]¶ Get the relaxed structure of nth symmetrically distinct interstitial.
Parameters: n – Symmetrically distinct interstitial index
-
get_relaxed_structure(n)[source]¶ Get the relaxed structure of nth symmetrically distinct interstitial.
Parameters: n – Symmetrically distinct interstitial index Note
0 corresponds to relaxed bulk structure
-
class
RelaxedInterstitial(struct_list, energy_list, valence_dict)[source]¶ Bases:
objectStores the relaxed supercell structures for each interstitial Used to compute formation energies, displacement of atoms near the the interstitial.
Parameters: - struct_list – List of structures(supercells). The first structure should represent relaxed bulk structure and the subsequent ones interstitial structures (with the extra interstitial site appended at the end).
- energy_list – List of energies for relaxed interstitial structures. The first energy should correspond to bulk structure
- valence_dict – Valences of elements in dictionary form
-
formation_energy(n, chem_pot=0)[source]¶ Compute the interstitial formation energy
Parameters: - n – Index of interstitials
- chem_pot – Chemical potential of interstitial site element. If not given, assumed as zero. The user is strongly urged to supply the chemical potential value
-
get_charge_coordination_number(n)[source]¶ Charge coordination number for nth interstitial.
Parameters: n – Index of interstitials
-
get_coordinated_bulk_sites(n)[source]¶ Bulk sites corresponding to the coordinated sites for nth interstitial.
Parameters: n – Index of interstitials
-
get_coordinated_site_displacement(n)[source]¶ Compute the total displacement of coordinated sites from the interstitial sites during the relaxation
Parameters: n – Index of defect site
-
get_coordinated_sites(n)[source]¶ Coordinated sites for nth interstitial.
Parameters: n – Index of interstitials
-
get_coordination_number(n)[source]¶ Coordination number for nth interstitial.
Parameters: n – Index of interstitials
-
get_defectsite(n)[source]¶ Returns the defect site of nth interstitial.
Parameters: n – Index of interstitial
-
get_percentage_bond_distance_change(n)[source]¶ Bond distance change after the introduction of interstitial.
Parameters: n – index of interstitials
-
class
StructureMotifInterstitial(struct, inter_elem, motif_types=('tet', 'oct'), op_threshs=(0.3, 0.5), dl=0.2, doverlap=1.0, facmaxdl=1.01, verbose=False)[source]¶ Bases:
pymatgen.analysis.defects.point_defects.DefectSubclass of Defect to generate interstitial sites at positions where the interstitialcy is coordinated by nearest neighbors in a way that resembles basic structure motifs (e.g., tetrahedra, octahedra). The algorithm will be formally introducted in an upcoming publication by Nils E. R. Zimmermann, Anubhav Jain, and Maciej Haranczyk, and it is already used by the Python Charged Defect Toolkit (PyCDT, https://arxiv.org/abs/1611.07481).
Generates symmetrically distinct interstitial sites at positions where the interstitial is coordinated by nearest neighbors in a pattern that resembles a supported structure motif (e.g., tetrahedra, octahedra).
Parameters: - struct (Structure) – input structure for which symmetrically distinct interstitial sites are to be found.
- inter_elem (string) – element symbol of desired interstitial.
- motif_types ([string]) – list of structure motif types that are to be considered. Permissible types are: tet (tetrahedron), oct (octahedron).
- op_threshs ([float]) – threshold values for the underlying order parameters to still recognize a given structural motif (i.e., for an OP value >= threshold the coordination pattern match is positive, for OP < threshold the match is negative.
- dl (float) – grid fineness in Angstrom. The input structure is divided into a grid of dimension a/dl x b/dl x c/dl along the three crystallographic directions, with a, b, and c being the lengths of the three lattice vectors of the input unit cell.
- doverlap (float) – distance that is considered to flag an overlap between any trial interstitial site and a host atom.
- facmaxdl (float) – factor to be multiplied with the maximum grid width that is then used as a cutoff distance for the clustering prune step.
- verbose (bool) – flag indicating whether (True) or not (False; default) to print additional information to screen.
-
enumerate_defectsites()[source]¶ Get all defect sites.
Returns: - list of periodic sites
- representing the interstitials.
Return type: defect_sites ([PeriodicSite])
-
get_coordinating_elements_cns(i)[source]¶ Get element-specific coordination numbers of defect with index i.
Returns: - dictionary storing the coordination numbers (int)
- with string representation of elements as keys. (i.e., {elem1 (string): cn1 (int), …}).
Return type: elem_cn (dict)
-
get_motif_type(i)[source]¶ Get the motif type of defect with index i (e.g., “tet”).
Returns: motif type. Return type: motif (string)
-
get_op_value(i)[source]¶ Get order-parameter value of defect with index i.
Returns: OP value. Return type: opval (float)
-
make_supercells_with_defects(scaling_matrix)[source]¶ Generate a sequence of supercells in which each supercell contains a single interstitial, except for the first supercell in the sequence which is a copy of the defect-free input structure.
Parameters: scaling_matrix (3x3 integer array) – scaling matrix to transform the lattice vectors. Returns: sequence of supercells. Return type: scs ([Structure])
-
class
Vacancy(structure, valences, radii)[source]¶ Bases:
pymatgen.analysis.defects.point_defects.DefectSubclass of Defect to generate vacancies and their analysis.
Parameters: - structure – pymatgen.core.structure.Structure
- valences – valences of elements as a dictionary
- radii – Radii of elements as a dictionary
-
get_coordsites_min_max_charge(n)[source]¶ Minimum and maximum charge of sites surrounding the vacancy site.
Parameters: n – Index of vacancy list
-
get_defectsite_effective_charge(n)[source]¶ Effective charge (In Kroger-Vink notation, cation vacancy has effectively -ve charge and anion vacancy has +ve charge.)
Parameters: n – Index of vacancy list Returns: Effective charnge of defect site
-
get_defectsite_structure_index(n)[source]¶ index of the vacacy site in the structure.sites list
Parameters: n – Index of vacancy list
-
get_surface_area(n)[source]¶ Surface area of the nth vacancy
Parameters: n – Index of symmetrically distinct vacancies list Returns: floating number representing volume of vacancy
-
get_volume(n)[source]¶ Volume of the nth vacancy
Parameters: n – Index of symmetrically distinct vacancies list Returns: floating number representing volume of vacancy
-
make_supercells_with_defects(scaling_matrix, species=None, limit_return_structures=False)[source]¶ Generate sequence of supercells in pymatgen.core.structure.Structure format, with each supercell containing one vacancy.
Parameters: - scaling_matrix – super cell scale parameters in matrix forms
- species – Species in list format only for which vacancy supercells are required. If not specified all the species are considered.
- limit_return_structures – Boolean or positive number If number, only that many structures are returned.
Returns: Supercells with vacancies. First supercell has no defects.
-
class
VacancyFormationEnergy(vacancy)[source]¶ Bases:
objectUsing GULP compute the vacancy formation energy. Works only for binary metal oxides due to the use of Buckingham Potentials
-
class
ValenceIonicRadiusEvaluator(structure)[source]¶ Bases:
objectComputes site valences and ionic radii for a structure using bond valence analyzer
Parameters: structure – pymatgen.core.structure.Structure -
radii¶ List of ionic radii of elements in the order of sites.
-
structure¶ Returns oxidation state decorated structure.
-
valences¶ List of oxidation states of elements in the order of sites.
-
-
symmetry_reduced_voronoi_nodes(structure, rad_dict, high_accuracy_flag=False, symm_flag=True)[source]¶ Obtain symmetry reduced voronoi nodes using Zeo++ and pymatgen.symmetry.finder.SpacegroupAnalyzer
Parameters: - strucutre – pymatgen Structure object
- rad_dict – Dictionary containing radii of spcies in the structure
- high_accuracy_flag – Flag denotting whether to use high accuracy version of Zeo++
- symm_flag – Flag denoting whether to return symmetrically distinct sites only
Returns: Symmetrically distinct voronoi nodes as pymatgen Strucutre