sumo.electronic_structure package¶
Module contents¶
Package containing functions for manipulating electron structure data.
Submodules¶
sumo.electronic_structure.bandstructure module¶
Module containing helper functions for dealing with band structures.
Todo
Extend get projections to allow specifying lm orbitals and atomic sites.
- sumo.electronic_structure.bandstructure.force_branches(bandstructure, return_mapping=False)[source]¶
Force a linemode band structure to contain branches.
Branches give a specific portion of the path from one high-symmetry point to another. Branches are required for the plotting methods to function correctly. Unfortunately, due to the pymatgen BandStructure implementation they require duplicate k-points in the band structure path. To avoid this unnecessary computational expense, this function can reconstruct branches in band structures without the duplicate k-points.
- Parameters:
bandstructure – A band structure object.
return_mapping – If True, return a mapping of the new k-points (with branches) to the original k-points.
- Returns:
A band structure with branches. If return_forced_branch_kpt_map is True, then a tuple is returned containing the band structure and the mapping from the new k-points to the original k-points.
- sumo.electronic_structure.bandstructure.get_projections(bs, selection, normalise=None)[source]¶
Returns orbital projections from a band structure.
- Parameters:
bs (
BandStructureSymmLine) – The band structure.selection (list) –
A list of
tupleorstringidentifying which projections to return. Projections can be specified by both element and orbital, for example:[('Bi', 's'), ('Bi', 'p'), ('S', 'p')]
If just the element is specified then all the orbitals of that element are combined. For example, the following will combine all the S orbitals into a single projection:
[('Bi', 's'), ('Bi', 'p'), 'S']
Particular orbitals can also be combined, for example:
[('Bi', 's'), ('Bi', 'p'), ('S', ('s', 'p', 'd'))]
normalise (
str, optional) –Normalisation the projections. Options are:
'all': Projections normalised against the sum of allother projections.
'select': Projections normalised against the sum of theselected projections.
None: No normalisation performed.
Defaults to
None.
- Returns:
A
listof orbital projections, in the same order as specified inselection, with the format:[ {spin: projections}, {spin: projections} ... ]
Where spin is a
pymatgen.electronic_structure.core.Spinobject and projections is anumpy.arrayof:projections[band_index][kpoint_index]
If there are no projections in the band structure, then an array of zeros is returned for each spin.
- Return type:
- sumo.electronic_structure.bandstructure.get_projections_by_branches(bs, selection, normalise=None)[source]¶
Returns orbital projections for each branch in a band structure.
- Parameters:
bs (
BandStructureSymmLine) – The band structure.selection (list) –
A list of
tupleorstringidentifying which projections to return. Projections can be specified by both element and orbital, for example:[('Sn', 's'), ('Bi', 'p'), ('S', 'p')]
If just the element is specified then all the orbitals of that element are combined. For example, the following will combine all the S orbitals into a single projection:
[('Bi', 's'), ('Bi', 'p'), 'S']
Particular orbitals can also be combined, for example:
[('Bi', 's'), ('Bi', 'p'), ('S', ('s', 'p', 'd'))]
normalise (
str, optional) –Normalisation the projections. Options are:
'all': Projections normalised against the sum of allother projections.
'select': Projections normalised against the sum of theselected projections.
None: No normalisation performed.
Defaults to
None.
- Returns:
A
listof orbital projections for each branch of the band structure, in the same order as specified inselection, with the format:[ [ {spin: projections} ], [ {spin: projections} ], ... ]
Where spin is a
pymatgen.electronic_structure.core.Spinobject and projections is anumpy.arrayof:projections[band_index][kpoint_index]
If there are no projections in the band structure, then an array of zeros is returned for each spin.
- Return type:
- sumo.electronic_structure.bandstructure.get_reconstructed_band_structure(list_bs, efermi=None, force_kpath_branches=True, return_forced_branch_kpt_map=False)[source]¶
Combine a list of band structures into a single band structure.
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.
This method will also ensure that any BandStructure objects will contain branches.
- Parameters:
list_bs (
listofBandStructureorBandStructureSymmLine) – The band structures.efermi (
float, optional) – The Fermi energy of the reconstructed band structure. If None, an average of all the Fermi energies across all band structures is used.force_kpath_branches (bool) – Force a linemode band structure to contain branches by adding repeated high-symmetry k-points in the path.
return_forced_branch_kpt_map (bool) – If True, return a mapping of the the new k-points to the original k-points.
- Returns:
A band structure object. The type depends on the type of the band structures in
list_bs. If return_forced_branch_kpt_map is True, then a tuple is returned containing the band structure and the mapping from the new k-points to the original k-points.- Return type:
pymatgen.electronic_structure.bandstructure.BandStructureorpymatgen.electronic_structure.bandstructureBandStructureSymmLine
sumo.electronic_structure.dos module¶
Module containing helper functions for dealing with
Dos and
CompleteDos objects.
- sumo.electronic_structure.dos.get_element_pdos(dos, element, sites, lm_orbitals=None, orbitals=None)[source]¶
Get the projected density of states for an element.
- Parameters:
dos (
CompleteDos) – The density of states.element (str) – Element symbol. E.g. ‘Zn’.
sites (tuple) – The atomic indices over which to sum the density of states, as a
tuple. Indices are zero based for each element. For example,(0, 1, 2)will sum the density of states for the 1st, 2nd and 3rd sites of the element specified.lm_orbitals (
tuple, optional) – The orbitals to decompose into their lm contributions (e.g. p -> px, py, pz). Should be provided as atupleofstr. For example,('p'), will extract the projected density of states for the px, py, and pz orbitals. Defaults toNone.orbitals (
tuple, optional) – The orbitals to extract from the projected density of states. Should be provided as atupleofstr. For example,('s', 'px', 'dx2')will extract the s, px, and dx2 orbitals, only. IfNone, all orbitals will be extracted. Defaults toNone.
- Returns:
The projected density of states. Formatted as a
dictmapping the orbitals toDosobjects. For example:{ 's': Dos, 'p': Dos }
- Return type:
- sumo.electronic_structure.dos.get_pdos(dos, lm_orbitals=None, atoms=None, elements=None)[source]¶
Extract the projected density of states from a CompleteDos object.
- Parameters:
dos (
CompleteDos) – The density of states.elements (
dict, optional) –The elements and orbitals to extract from the projected density of states. Should be provided as a
dictwith the keys as the element names and corresponding values as atupleof orbitals. For example, the following would extract the Bi s, px, py and d orbitals:{'Bi': ('s', 'px', 'py', 'd')}
If an element is included with an empty
tuple, all orbitals for that species will be extracted. Ifelementsis not set or set toNone, all elements for all species will be extracted.lm_orbitals (
dict, optional) –The orbitals to decompose into their lm contributions (e.g. p -> px, py, pz). Should be provided as a
dict, with the elements names as keys and atupleof orbitals as the corresponding values. For example, the following would be used to decompose the oxygen p and d orbitals:{'O': ('p', 'd')}
atoms (
dict, optional) –Which atomic sites to use when calculating the projected density of states. Should be provided as a
dict, with the element names as keys and atupleofintspecifying the atomic indices as the corresponding values. The elemental projected density of states will be summed only over the atom indices specified. If an element is included with an emptytuple, then all sites for that element will be included. The indices are 0 based for each element specified in the POSCAR. For example, the following will calculate the density of states for the first 4 Sn atoms and all O atoms in the structure:{'Sn': (1, 2, 3, 4), 'O': (, )}
If
atomsis not set or set toNonethen all atomic sites for all elements will be considered.
- Returns:
The projected density of states. Formatted as a
dictofdictmapping the elements and their orbitals toDosobjects. For example:{ 'Bi': {'s': Dos, 'p': Dos ... }, 'S': {'s': Dos} }
- Return type:
- sumo.electronic_structure.dos.load_dos(vasprun, elements=None, lm_orbitals=None, atoms=None, gaussian=None, total_only=False, log=False, adjust_fermi=True, scissor=None)[source]¶
Load a vasprun and extract the total and projected density of states.
- Parameters:
vasprun (str) – Path to a vasprun.xml or vasprun.xml.gz file or a
pymatgen.io.vasp.outputs.Vasprunobject.elements (
dict, optional) –The elements and orbitals to extract from the projected density of states. Should be provided as a
dictwith the keys as the element names and corresponding values as atupleof orbitals. For example, the following would extract the Bi s, px, py and d orbitals:{'Bi': ('s', 'px', 'py', 'd')}
If an element is included with an empty
tuple, all orbitals for that species will be extracted. Ifelementsis not set or set toNone, all elements for all species will be extracted.lm_orbitals (
dict, optional) –The orbitals to decompose into their lm contributions (e.g. p -> px, py, pz). Should be provided as a
dict, with the elements names as keys and atupleof orbitals as the corresponding values. For example, the following would be used to decompose the oxygen p and d orbitals:{'O': ('p', 'd')}
atoms (
dict, optional) –Which atomic sites to use when calculating the projected density of states. Should be provided as a
dict, with the element names as keys and atupleofintspecifying the atomic indices as the corresponding values. The elemental projected density of states will be summed only over the atom indices specified. If an element is included with an emptytuple, then all sites for that element will be included. The indices are 0 based for each element specified in the POSCAR. For example, the following will calculate the density of states for the first 4 Sn atoms and all O atoms in the structure:{'Sn': (1, 2, 3, 4), 'O': (, )}
If
atomsis not set or set toNonethen all atomic sites for all elements will be considered.gaussian (
float, optional) – Broaden the density of states using convolution with a gaussian function. This parameter controls the sigma or standard deviation of the gaussian distribution.total_only (
bool, optional) – Only extract the total density of states. Defaults toFalse.log (
bool) – Print logging messages. Defaults toFalse.adjust_fermi (
bool, optional) – Shift the Fermi level to sit at the valence band maximum (does not affect metals).scissor (
float, optional) – Apply a scissor operator to shift the band gap to equal scissor (rigid shift of the densities above the mid gap). Not compatible with metallic systems.
- Returns:
The total and projected density of states. Formatted as a
tupleof(dos, pdos), wheredosis aDosobject containing the total density of states andpdosis adictofdictmapping the elements and their orbitals toDosobjects. For example:{ 'Bi': {'s': Dos, 'p': Dos ... }, 'S': {'s': Dos} }
- Return type:
- sumo.electronic_structure.dos.sort_orbitals(element_pdos)[source]¶
Sort the orbitals of an element’s projected density of states.
Sorts the orbitals based on a standard format. E.g. s < p < d. Will also sort lm decomposed orbitals. This is useful for plotting/saving.
- sumo.electronic_structure.dos.write_files(dos, pdos, prefix=None, directory=None, zero_to_efermi=True)[source]¶
Write the density of states data to disk.
- Parameters:
dos (
DosorCompleteDos) – The total density of states.pdos (dict) –
The projected density of states. Formatted as a
dictofdictmapping the elements and their orbitals toDosobjects. For example:{ 'Bi': {'s': Dos, 'p': Dos}, 'S': {'s': Dos} }
prefix (
str, optional) – A prefix for file names.directory (
str, optional) – The directory in which to save files.zero_to_efermi (
bool, optional) – Normalise the energy such that the Fermi level is set as 0 eV.
sumo.electronic_structure.effective_mass module¶
Module containing helper functions for calculating band effective masses.
- sumo.electronic_structure.effective_mass.fit_effective_mass(distances, energies, parabolic=True)[source]¶
Fit the effective masses using either a parabolic or nonparabolic fit.
- Parameters:
distances (
numpy.ndarray) – The x-distances between k-points in reciprocal Angstroms, normalised to the band extrema.energies (
numpy.ndarray) – The band eigenvalues normalised to the eigenvalue of the band extrema.parabolic (
bool, optional) – Use a parabolic fit of the band edges. IfFalsethen nonparabolic fitting will be attempted. Defaults toTrue.
- Returns:
The effective mass in units of electron rest mass, \(m_0\).
- Return type:
- sumo.electronic_structure.effective_mass.get_fitting_data(bs, spin, band_id, kpoint_id, num_sample_points=3)[source]¶
Extract fitting data for band extrema based on spin, kpoint and band.
Searches forward and backward from the extrema point, but will only sample there data if there are enough points in that direction.
- Parameters:
bs (
BandStructureSymmLine) – The band structure.spin (
Spin) – Which spin channel to sample.band_id (int) – Index of the band to sample.
kpoint_id (int) – Index of the kpoint to sample.
- Returns:
The data necessary to calculate the effective mass, along with some metadata. Formatted as a
listofdict, each with the keys:- ’energies’ (
numpy.ndarray) Band eigenvalues in eV.
- ’distances’ (
numpy.ndarray) Distances of the k-points in reciprocal space.
- ’band_id’ (
int) The index of the band,
- ’spin’ (
Spin) The spin channel
- ’start_kpoint’ (
int) The index of the k-point at which the band extrema occurs
- ’end_kpoint’ (
int) The k-point towards which the data has been sampled.
- ’energies’ (
- Return type:
sumo.electronic_structure.optics module¶
Module containing functions to process dielectric and optical absorption data.
Todo
Remove magic values
- sumo.electronic_structure.optics.broaden_eps(dielectric, sigma)[source]¶
Apply gaussian broadening to the dielectric response.
- Parameters:
dielectric (tuple) –
The high-frequency dielectric data, following the same format as
pymatgen.io.vasp.outputs.Vasprun.dielectric. This is atuplecontaining the energy, the real part of the dielectric tensor, and the imaginary part of the tensor, as alistoffloats. E.g.:( [energies], [[real_xx, real_yy, real_zz, real_xy, real_yz, real_xz]], [[imag_xx, imag_yy, imag_zz, imag_xy, imag_yz, imag_xz]] )
sigma (float) – Standard deviation for gaussian broadening.
- Returns:
The broadened dielectric response. Returned as a tuple containing the energy, the real part of the dielectric tensor, and the imaginary part of the tensor. E.g.:
( [energies], [[real_xx, real_yy, real_zz, real_xy, real_yz, real_xz]], [[imag_xx, imag_yy, imag_zz, imag_xy, imag_yz, imag_xz]] )
- Return type:
- sumo.electronic_structure.optics.calculate_dielectric_properties(dielectric, properties, mode='average')[source]¶
Calculate optical properties from the dielectric function
Supported properties:
Absorption
The unit of alpha is \(\mathrm{cm}^{-1}\).
Refractive index \(n\) has real and imaginary parts:
\[n = [(e^\prime + ie^{\prime\prime} / e_0]^{1/2} = n^\prime + in^{\prime\prime}\]Relationship between \(a\) and imaginary \(n^{\prime\prime}\):
\[a = 4 \pi n^{\prime\prime} / \lambda\]Where:
\[\lambda = hc/E\]- Parameters:
dielectric (tuple) –
The high-frequency dielectric data, following the same format as
pymatgen.io.vasp.Vasprun.dielectric. This is atuplecontaining the energy, the real part of the dielectric tensor, and the imaginary part of the tensor, as alistoffloats. E.g.:( [energies], [[real_xx, real_yy, real_zz, real_xy, real_yz, real_xz]], [[imag_xx, imag_yy, imag_zz, imag_xy, imag_yz, imag_xz]] )
properties (set) – The set of properties to return. Intermediate properties will be calculated as needed. Accepted values: ‘eps_real’, ‘eps_imag’, ‘absorption’, ‘loss’, ‘n_real’, ‘n_imag’
mode (str) –
The format of the results: Options are:
”average”: Average the dielectric response for all directions.
”trace”: The trace of the dielectric tensor (i.e., along xx, yy, zz).
”full”: The full dielectric tensor (i.e., xx, xy, xz, yx, yy, yz, zx, zy, zz).
”eigs”: Calculate the eigenvalues of the dielectric response. The eigenvalues are sorted from smallest to largest.
- Returns:
The optical absorption is given in \(\mathrm{cm}^{-1}\). If
modeis"average", the property_value will be returned as:[property]
If
modeis"trace", the data will be returned as:[property_xx, property_yy, property_zz]
If
modeis"full", the data will be returned as:[xx, xy, xz, yx, yy, yz, zx, zy, zz]
If
modeis"eigs", the data will be returned as:[eig_1, eig_2, eig_3]
In all cases these are collected in a results dictionary with keys corresponding to the selected properties, e.g.:
{'absorption': [absorption], 'eps_real': [eps_real]}
- Return type:
tuple of
(energies, {property_name: property_value})
- sumo.electronic_structure.optics.kkr(de, eps_imag, cshift=1e-06)[source]¶
Kramers Kronig transformation of imaginary dielectric function
- Parameters:
de (
float) – Energy difference between evenly-spaced energy values corresponding to dielectric dataeps_imag (
listornp.array) – Evenly-spaced sequence of frequency-dependent 3x3 dielectric matrices (imaginary component only)cshift (
float, optional) – imaginary finite shift used in integration; this should be small (and results should not be very sensitive)
- Returns:
(
numpy.array) Real part of frequency-dependent dielectric function corresponding to eps_imag. Array shape (NEDOS, 3, 3)
- sumo.electronic_structure.optics.write_files(abs_data, basename='absorption', prefix=None, directory=None)[source]¶
Write the absorption or loss spectra to a file.
Note that this function expects to receive an iterable series of spectra.
- Parameters:
abs_data (tuple) –
Series (either
listortuple) of optical absorption or loss spectra. Each spectrum should be formatted as atupleoflistoffloat. If the data has been averaged, each spectrum should be:([energies], [alpha])
Else, if the data has not been averaged, each spectrum should be:
([energies], [alpha_xx, alpha_yy, alpha_zz]).
prefix (
str, optional) – Prefix for file names.directory (
str, optional) – The directory in which to save files.