electronic_structure Package
electronic_structure Package
This package contains electronic structure related tools and analyses.
bandstructure Module
This module provides classes to define everything related to band structures.
-
class BandStructure(kpoints, eigenvals, lattice, efermi, labels_dict=None, coords_are_cartesian=False, structure=None, projections=None)[source]
Bases: object
This is the most generic band structure data possible
it’s defined by a list of kpoints + energies for each of them
- Args:
- kpoints:
- list of kpoint as numpy arrays, in frac_coords of the
given lattice by default
- eigenvals:
- dict of energies for spin up and spin down
{Spin.up:[][],Spin.down:[][]}, the first index of the array
[][] refers to the band and the second to the index of the
kpoint. The kpoints are ordered according to the order of the
kpoints array. If the band structure is not spin polarized, we
only store one data set under Spin.up
- lattice:
- The reciprocal lattice as a pymatgen Lattice object.
- label_dict:
- (dict) of {} this link a kpoint (in frac coords or cartesian
coordinates depending on the coords).
- coords_are_cartesian:
- Whether coordinates are cartesian.
- efermi:
- fermi energy
- labels_dict:
- (dict) of {} this links a kpoint (in frac coords or cartesian
coordinates depending on the coords) to a label.
- coords_are_cartesian:
- Whether coordinates are cartesian.
- structure:
- the crystal structure (as a pymatgen Structure object)
associated with the band structure. This is needed if we
provide projections to the band structure
- projections:
- dict of orbital projections for spin up and spin down
{Spin.up:[][{Orbital:[]}],Spin.down:[][{Orbital:[]}]. The
format follows the one from eigenvals: The first index of the
array refers to the band and the second to the index of the
kpoint. The kpoints are ordered according to the order of the
kpoints array. For each band and kpoint, we associate a
dictionary indicating projections on orbitals and on different
sites the keys of the dictionary are Orbital objects and the
values are the projections on each site ordered as in the
structure object. If the band structure is not spin polarized,
we only store one data set under Spin.up.
-
bands[source]
returns the eigenvalues for each kpoints as a dictionary
{Spin.up:[][],Spin.down:[][]}, the first index of the array
[][] refers to the band and the second to the index of the
kpoint. The kpoints are ordered according to the order of the
self.kpoints. If the band structure is not spin polarized, we
only store one data set under Spin.up
-
efermi[source]
the fermi energy
-
get_projection_on_elements()[source]
Method returning a dictionary of projections on elements.
- Returns:
- a dictionary in the {Spin.up:[][{Element:values}],
Spin.down:[][{Element:values}]} format
if there is no projections in the band structure
returns an empty dict
-
get_projections_on_elts_and_orbitals(dictio)[source]
Method returning a dictionary of projections on elements and specific
orbitals
- Args:
- dictio:
- A dictionary of Elements and Orbitals for which we want to have
projections on. It is given as: {Element:[orbitals]},
e.g., {‘Cu’:[‘d’,’s’]}
- Returns:
- A dictionary of projections on elements in the
{Spin.up:[][{Element:{orb:values}}],
Spin.down:[][{Element:{orb:values}}]} format
if there is no projections in the band structure returns an empty
dict.
-
is_spin_polarized[source]
True if the band structure is spin-polarized, False otherwise
-
kpoints[source]
the list of kpoints (as Kpoint objects) in the band structure
-
lattice[source]
the lattice of the band structure as a pymatgen Lattice object
-
nb_bands[source]
returns the number of bands in the band structure
-
class BandStructureSymmLine(kpoints, eigenvals, lattice, efermi, labels_dict, coords_are_cartesian=False, structure=None, projections=None)[source]
Bases: pymatgen.electronic_structure.bandstructure.BandStructure, pymatgen.serializers.json_coders.MSONable
This object stores band structures along selected (symmetry) lines in the
Brillouin zone. We call the different symmetry lines (ex: Gamma to Z)
“branches”.
- Args:
- kpoints:
- list of kpoint as numpy arrays, in frac_coords of the
given lattice by default
- eigenvals:
- dict of energies for spin up and spin down
{Spin.up:[][],Spin.down:[][]}, the first index of the array
[][] refers to the band and the second to the index of the
kpoint. The kpoints are ordered according to the order of the
kpoints array. If the band structure is not spin polarized, we
only store one data set under Spin.up.
- lattice:
- The reciprocal lattice.
- efermi:
- fermi energy
- label_dict:
- (dict) of {} this link a kpoint (in frac coords or cartesian
coordinates depending on the coords).
- coords_are_cartesian:
- Whether coordinates are cartesian.
- structure:
- the crystal structure (as a pymatgen Structure object)
associated with the band structure. This is needed if we
provide projections to the band structure.
- projections:
- dict of orbital projections for spin up and spin down
{Spin.up:[][{Orbital:[]}],Spin.down:[][{Orbital:[]}]. The
format follows the one from eigenvals: the first index of the
array refers to the band and the second to the index of the
kpoint. The kpoints are ordered according to the order of the
kpoints array. For each band and kpoint, we associate a
dictionary indicating projections on orbitals and on different
sites the keys of the dictionary are Orbital objects and the
values are the projections on each site ordered as in the
structure object. If the band structure is not spin polarized,
we only store one data set under Spin.up.
-
apply_scissor(new_band_gap)[source]
Apply a scissor operator (shift of the CBM) to fit the given band gap.
If it’s a metal. We look for the band crossing the fermi level
and shift this one up. This will not work all the time for metals!
- Args:
- new_band_gap:
- the band gap the scissor band structure need to have.
- Returns:
- a BandStructureSymmLine object with the applied scissor shift
-
classmethod from_dict(d)[source]
- Args:
- A dict with all data for a band structure symm line object.
- Returns:
- A BandStructureSymmLine object
-
get_band_gap()[source]
Returns band gap data.
- Returns:
- a dict {“energy”,”direct”,”transition”}:
- “energy”:
- the band gap energy
- “direct”:
- A boolean telling if the gap is direct (True)
or not (False)c
- “transition”:
- The kpoint labels of the transition (e.g., “Gamma-X”)
-
get_branch(index)[source]
Returns in what branch(es) is the kpoint. There can be several
branches.
- Args:
- index:
- the kpoint index
- Returns:
- A list of dictionaries [{“name”,”start_index”,”end_index”,”index”}]
indicating all branches in which the k_point is. It takes into
account the fact that one kpoint (e.g., Gamma) can be in several
branches
-
get_cbm()[source]
Returns data about the CBM.
- Returns:
- {“band_index”,”kpoint_index”,”kpoint”,”energy”}
- “band_index”:
- A dict with spin keys pointing to a list of the indices of
the band containing the CBM (please note that you can have
several bands sharing the CBM) {Spin.up:[],Spin.down:[]}
- “kpoint_index”:
- The list of indices in self._kpoints for the kpoint vbm.
Please note that there can be several kpoint_indices
relating to the same kpoint (e.g., Gamma can occur at
different spots in the band structure line plot)
- “kpoint”:
- The kpoint (as a kpoint object)
- “energy”:
- The energy of the CBM
- “projections”:
- The projections along sites and orbitals of the CBM
if any projection data is available (else it is an empty
dictionnary)
the format is similar to the projections field in
BandStructure: {spin: {‘Orbital’: [proj]}} where the array
[proj] is ordered according to the sites in structure
-
get_direct_band_gap()[source]
Returns the direct band gap.
- Returns:
- the value of the direct band gap
-
get_equivalent_kpoints(index)[source]
Returns the list of kpoint indices equivalent (meaning they are the
same frac coords) to the given one.
- Args:
- index:
- the kpoint index
- Returns:
- a list of equivalent indices
TODO: now it uses the label we might want to use coordinates instead
(in case there was a mislabel)
-
get_vbm()[source]
Returns data about the VBM.
- Returns:
- dict as {“band_index”,”kpoint_index”,”kpoint”,”energy”}
- “band_index”:
- A dict with spin keys pointing to a list of the indices of
the band containing the VBM (please note that you can have
several bands sharing the VBM) {Spin.up:[],Spin.down:[]}
- “kpoint_index”:
- The list of indices in self._kpoints for the kpoint vbm.
Please note that there can be several kpoint_indices
relating to the same kpoint (e.g., Gamma can occur at
different spots in the band structure line plot)
- “kpoint”:
- The kpoint (as a kpoint object)
- “energy”:
- The energy of the VBM
- “projections”:
- The projections along sites and orbitals of the VBM
if any projection data is available (else it is an empty
dictionnary)
the format is similar to the projections field in
BandStructure: {spin:{‘Orbital’: [proj]}} where the array
[proj] is ordered according to the sites in structure
-
is_metal()[source]
Check if the band structure indicates a metal by looking if the fermi
level crosses a band.
- Returns:
- True if a metal, False if not
-
to_dict[source]
Json-serializable dict representation of BandStructureSymmLine.
-
class Kpoint(coords, lattice, to_unit_cell=False, coords_are_cartesian=False, label=None)[source]
Bases: pymatgen.serializers.json_coders.MSONable
Class to store kpoint objects. A kpoint is defined with a lattice and frac
or cartesian coordinates syntax similar than the site object in
pymatgen.core.structure.
- Args:
- coords:
- coordinate of the kpoint as a numpy array
- lattice:
- a pymatgen.core.lattice.Lattice lattice object representing the
reciprocal lattice of the kpoint
- to_unit_cell:
- translates fractional coordinate to the basic unit cell, i.e.
all fractional coordinates satisfy 0 <= a < 1.
Defaults to False.
- coords_are_cartesian:
- Boolean indicating if the coordinates given are in cartesian or
fractional coordinates (by default fractional)
- label:
- the label of the kpoint if any (None by default)
-
a[source]
Fractional a coordinate of the kpoint
-
b[source]
Fractional b coordinate of the kpoint
-
c[source]
Fractional c coordinate of the kpoint
-
cart_coords[source]
The cartesian coordinates of the kpoint as a numpy array
-
frac_coords[source]
The fractional coordinates of the kpoint as a numpy array
-
label[source]
The label associated with the kpoint
-
lattice[source]
The lattice associated with the kpoint. It’s a
pymatgen.core.lattice.Lattice object
-
to_dict[source]
Json-serializable dict representation of a kpoint
-
get_reconstructed_band_structure(list_bs, efermi=None)[source]
This method takes a list of band structures
and reconstruct one band structure object from all of them
this is typically very useful when you split non self consistent
band structure runs in several independent jobs and want to merge back
the results
- Args:
- list_bs:
- A list of BandStructure
- efermi:
- The fermi energy of the reconstructed band structure. If none
is assigned an average of all the fermi energy in each object
in the list_bs is used.
- Returns:
- A BandStructure or BandStructureSymmLine object (depending on
the type of the list_bs objects)
core Module
This module provides core classes needed by all define electronic structure,
such as the Spin, Orbital, etc.
-
class Orbital[source]
Bases: object
Enum type for OrbitalType. Indices are basically the azimutal quantum
number, l.
-
all_orbitals = (s, py, pz, px, dxy, dyz, dz2, dxz, dx2, f_3, f_2, f_1, f0, f1, f2, f3)
-
dx2 = dx2
-
dxy = dxy
-
dxz = dxz
-
dyz = dyz
-
dz2 = dz2
-
f0 = f0
-
f1 = f1
-
f2 = f2
-
f3 = f3
-
f_1 = f_1
-
f_2 = f_2
-
f_3 = f_3
-
static from_string(orb_str)[source]
Returns an orbital from a string representation, e.g., “s”, “px”.
-
static from_vasp_index(i)[source]
Returns an orbital based on the index of the orbital in VASP runs.
-
px = px
-
py = py
-
pz = pz
-
s = s
-
class Spin[source]
Bases: object
Enum type for Spin. Only up and down.
Usage: Spin.up, Spin.down.
-
all_spins = (1, -1)
-
down = -1
-
static from_int(i)[source]
Provides the spin from an int. +1 == Spin.up, -1 == Spin.down.
- Args:
- i:
- integer representing direction of spin.
-
up = 1
dos Module
This module defines classes to represent the density of states, etc.
-
class CompleteDos(structure, total_dos, pdoss)[source]
Bases: pymatgen.electronic_structure.dos.Dos
This wrapper class defines a total dos, and also provides a list of PDos.
Mainly used by pymatgen.io.vaspio.Vasprun to create a complete Dos from
a vasprun.xml file. You are unlikely to try to generate this object
manually.
-
structure
Structure associated with the CompleteDos.
-
pdos
Dict of partial densities of the form {Site:{Orbital:{Spin:Densities}}}
- Args:
- structure:
- Structure associated with this particular DOS.
- total_dos:
- total Dos for structure
- pdoss:
- The pdoss are supplied as an {Site:{Orbital:{
Spin:Densities}}}
-
classmethod from_dict(d)[source]
Returns CompleteDos object from dict representation.
-
get_element_dos()[source]
Get element projected Dos.
- Returns:
- dict of {Element: Dos}
-
get_site_dos(site)[source]
Get the total Dos for a site (all orbitals).
- Args:
- site:
- Site in Structure associated with CompleteDos.
- Returns:
- Dos containing summed orbital densities for site.
-
get_site_orbital_dos(site, orbital)[source]
Get the Dos for a particular orbital of a particular site.
- Args:
- site:
- Site in Structure associated with CompleteDos.
- orbital:
- Orbital in the site.
- Returns:
- Dos containing densities for orbital of site.
-
get_site_t2g_eg_resolved_dos(site)[source]
Get the t2g, eg projected DOS for a particular site.
- Args:
- site:
- Site in Structure associated with CompleteDos.
- Returns:
- A dict {“e_g”: Dos, “t2g”: Dos} containing summed e_g and t2g DOS
for the site.
-
get_spd_dos()[source]
Get orbital projected Dos.
- Returns:
- dict of {orbital: Dos}, e.g. {“s”: Dos object, ...}
-
to_dict[source]
Json-serializable dict representation of CompleteDos.
-
class Dos(efermi, energies, densities)[source]
Bases: pymatgen.serializers.json_coders.MSONable
Basic DOS object. All other DOS objects are extended versions of this
object.
- Args:
- efermi:
- Fermi level energy
- energies:
- A sequences of energies
- densities:
- A dict of {Spin: np.array} representing the density of states
for each Spin.
-
classmethod from_dict(d)[source]
Returns Dos object from dict representation of Dos.
-
get_cbm_vbm(tol=0.001, abs_tol=False, spin=None)[source]
Expects a DOS object and finds the cbm and vbm.
- Args:
- tol:
- tolerance in occupations for determining the gap
- abs_tol:
- an absolute tolerance (True) and a relative one (False)
- spin:
- Possible values are:
- None - finds the gap in the summed densities
Up - finds the gap in the up spin channel
Down - finds the gap in teh down spin channel
- Returns:
- (cbm, vbm): float in eV corresponding to the gap
-
get_densities(spin=None)[source]
Returns the density of states for a particular spin.
- Args:
- spin:
- Spin
- Returns:
- Returns the density of states for a particular spin. If Spin is
None, the sum of all spins is returned.
-
get_gap(tol=0.001, abs_tol=False, spin=None)[source]
Expects a DOS object and finds the gap.
- Args:
- tol:
- tolerance in occupations for determining the gap
- abs_tol:
- an absolute tolerance (True) and a relative one (False)
- spin:
- Possible values are:
- None - finds the gap in the summed densities
Up - finds the gap in the up spin channel
Down - finds the gap in teh down spin channel
- Returns:
- gap in eV
-
get_interpolated_gap(tol=0.001, abs_tol=False, spin=None)[source]
Expects a DOS object and finds the gap
- Args:
- tol:
- tolerance in occupations for determining the gap
- abs_tol:
- Set to True for an absolute tolerance and False for a relative
one.
- spin:
- Possible values are:
- None - finds the ap in the summed densities
Up - finds the gap in the up spin channel
Down - finds the gap in teh down spin channel
- Returns:
- (gap, cbm, vbm):
- Tuple of floats in eV corresponding to the gap, cbm and vbm.
-
get_interpolated_value(energy)[source]
Returns interpolated density for a particular energy.
- Args:
- energy:
- Energy to return the density for.
-
get_smeared_densities(*args, **kwargs)[source]
Returns the Dict representation of the densities, {Spin: densities},
but with a Gaussian smearing of std dev sigma applied about the fermi
level.
- Args:
- sigma:
- Std dev of Gaussian smearing function.
- Returns:
- Dict of Gaussian-smeared densities.
-
to_dict[source]
Json-serializable dict representation of Dos.
-
add_densities(density1, density2)[source]
Method to sum two densities.
- Args:
- density1:
- First density.
- density2:
- Second density.
- Returns:
- Dict of {spin: density}.
plotter Module
This module implements plotter for DOS and band structure.
-
class BSPlotter(bs)[source]
Bases: object
Class to plot or get data to facilitate the plot of band structure objects.
- Args:
- bs:
- A BandStructureSymmLine object.
-
bs_plot_data(zero_to_efermi=True)[source]
Get the data nicely formatted for a plot
- Args:
- zero_to_efermi:
- Automatically subtract off the Fermi energy from the
eigenvalues and plot.
- Returns:
- A dict of the following format:
- ticks:
- A dict with the ‘distances’ at which there is a kpoint (the
x axis) and the labels (None if no label)
- energy:
- A dict storing bands for spin up and spin down data
{Spin:[band_index][k_point_index]} as a list (one element
for each band) of energy for each kpoint.
- vbm:
- A list of tuples (distance,energy) marking the vbms. The
energies are shifted with respect to the fermi level is the
option has been selected.
- cbm:
- A list of tuples (distance,energy) marking the cbms. The
energies are shifted with respect to the fermi level is the
option has been selected.
- lattice:
- The reciprocal lattice.
- zero_energy:
- This is the energy used as zero for the plot.
- band_gap:
- A string indicating the band gap and its nature (empty if
it’s a metal).
- is_metal:
- True if the band structure is metallic (i.e., there is at
least one band crossing the fermi level).
-
get_plot(zero_to_efermi=True, ylim=None, smooth=False)[source]
get a matplotlib object for the bandstructure plot.
Blue lines are up spin, red lines are down
spin.
- Args:
- zero_to_efermi:
- Automatically subtract off the Fermi energy from the
eigenvalues and plot (E-Ef).
- ylim:
- specify the y-axis (energy) limits;
by default None let the code choose.
It is vbm-4 and cbm+4 if insulator
efermi-10 and efermi+10 if metal
- smooth:
- interpolates the bands by a spline
cubic
-
get_ticks()[source]
Get all ticks and labels for a band structure plot.
- Returns:
- A dict with
- ‘distance’: a list of distance at which ticks should be set.
‘label’: a list of label for each of those ticks.
-
plot_brillouin()[source]
plot the Brillouin zone
-
plot_compare(other_plotter)[source]
plot two band structure for comparison. One is in red the other in blue
(no difference in spins). The two band structures need to be defined
on the same symmetry lines! and the distance between symmetry lines is
the one of the band structure used to build the BSPlotter
- Args:
- another band structure object defined along the same symmetry lines
- Returns:
- a matplotlib object with both band structures
-
save_plot(filename, img_format='eps', ylim=None, zero_to_efermi=True, smooth=False)[source]
Save matplotlib plot to a file.
- Args:
- filename:
- Filename to write to.
- img_format:
- Image format to use. Defaults to EPS.
- ylim:
- Specifies the y-axis limits.
-
show(zero_to_efermi=True, ylim=None, smooth=False)[source]
Show the plot using matplotlib.
- Args:
- zero_to_efermi:
- Automatically subtract off the Fermi
energy from the eigenvalues
and plot (E-Ef).
- ylim
- specify the y-axis (energy) limits; by default None
let the code choose.
It is vbm-4 and cbm+4 if insulator
efermi-10 and efermi+10 if metal
- smooth:
- interpolates the bands by a spline
cubic
-
class BSPlotterProjected(bs)[source]
Bases: pymatgen.electronic_structure.plotter.BSPlotter
Class to plot or get data to facilitate the plot of band structure objects
projected along orbitals, elements or sites
- Args:
- bs:
- A BandStructureSymmLine object with projections.
-
get_elt_projected_plots(zero_to_efermi=True, ylim=None)[source]
Method returning a plot composed of subplots along different elements
- Returns:
- a pylab object with different subfigures for each projection
The blue and red colors are for spin up and spin down
The bigger the red or blue dot in the band structure the higher
character for the corresponding element and orbital
-
get_elt_projected_plots_color(zero_to_efermi=True, elt_ordered=None, smooth=False)[source]
returns a pylab plot object with one plot where the band structure
line color depends on the character of the band (along different
elements). Each element is associated with red, green or blue
and the corresponding rgb color depending on the character of the band
is used. The method can only deal with binary and ternary compounds
the method does not make a difference for now between spin up
and spin down
- Args:
- elt_ordered:
- a list of Element ordered. The first one is red, second green,
last blue
- Returns:
- a pylab object
-
get_projected_plots_dots(dictio, zero_to_efermi=True, ylim=None)[source]
Method returning a plot composed of subplots along different elements
and orbitals.
- Args:
- dictio:
- the element and orbitals you want a projection on. The format
is {Element:[Orbitals]} for instance {‘Cu’:[‘d’,’s’],’O’:[‘p’]}
will give projections for Cu on d and s orbitals and on oxygen
p.
- Returns:
- a pylab object with different subfigures for each projection
The blue and red colors are for spin up and spin down.
The bigger the red or blue dot in the band structure the higher
character for the corresponding element and orbital.
-
class DosPlotter(zero_at_efermi=True, stack=False, sigma=None)[source]
Bases: object
Class for plotting DOSes.
- Args:
- zero_at_efermi:
- Whether to shift all Dos to have zero energy at the fermi
energy. Defaults to True.
- stack:
- Whether to plot the DOS as a stacked area graph
- key_sort_func:
- function used to sort the dos_dict keys.
- sigma:
- A float specifying a standard deviation for Gaussian smearing
the DOS for nicer looking plots. Defaults to None for no
smearing.
-
add_dos(label, dos)[source]
Adds a dos for plotting.
- Args:
- label:
- label for the DOS. Must be unique.
- dos:
- Dos object
-
add_dos_dict(dos_dict, key_sort_func=None)[source]
Add a dictionary of doses, with an optional sorting function for the
keys.
- Args:
- dos_dict:
- dict of {label: Dos}
- key_sort_func:
- function used to sort the dos_dict keys.
-
get_dos_dict()[source]
Returns the added doses as a json-serializable dict. Note that if you
have specified smearing for the DOS plot, the densities returned will
be the smeared densities, not the original densities.
- Returns:
- Dict of dos data. Generally of the form, {label: {‘energies’:..,
‘densities’: {‘up’:...}, ‘efermi’:efermi}}
-
get_plot(xlim=None, ylim=None)[source]
Get a matplotlib plot showing the DOS.
- Args:
- xlim:
- Specifies the x-axis limits. Set to None for automatic
determination.
- ylim:
- Specifies the y-axis limits.
-
save_plot(filename, img_format='eps', xlim=None, ylim=None)[source]
Save matplotlib plot to a file.
- Args:
- filename:
- Filename to write to.
- img_format:
- Image format to use. Defaults to EPS.
- xlim:
- Specifies the x-axis limits. Set to None for automatic
determination.
- ylim:
- Specifies the y-axis limits.
-
show(xlim=None, ylim=None)[source]
Show the plot using matplotlib.
- Args:
- xlim:
- Specifies the x-axis limits. Set to None for automatic
determination.
- ylim:
- Specifies the y-axis limits.
-
get_lines_voronoi(data)[source]
-
qvertex_target(data, index)[source]
Input data should be in the form of a list of a list of floats.
index is the index of the targeted point
Returns the vertices of the voronoi construction around this target point.