Welcome to DFT Tools’s documentation!¶
DFT Tools is a python library for parsing, post-processing and presenting numerical data generated by simulation codes in physics and chemistry. The abbreviation ‘DFT’ comes from the denisty function theory being implemented in these codes.
With DFT Tools you will be able to:
- parse numerical data from textual output such as band structures, data on the grid, etc.;
- manipulate the data: build supercells, calculate density of states, etc.;
- visualise the data;
Contents:
Introduction to DFT Tools¶
A number of codes implementing DFT and it’s flavors is available in the web, see Wikipedia for example. The developers of these codes are usually scientists who never aimed to develop a user-friendly application mainly because they are not get paid for that. Thus, to be able to use such codes one has to master several tools, among which is data post-processing and presentation.
An average DFT code produces a set of text and binary data during the run. Typically, the data cannot be plotted directly and one needs a program to collect this data and present it. Here is an example of a Quantum Espresso band structure:
k = 0.0000 0.0000 0.0000 band energies (ev):
-5.8099 6.2549 6.2549 6.2549 8.8221 8.8221 8.8221 9.7232
k = 0.0000 0.0000 0.1000 band energies (ev):
-5.7668 5.9810 6.0722 6.0722 8.7104 9.0571 9.0571 9.9838
k = 0.0000 0.0000 0.2000 band energies (ev):
-5.6337 5.3339 5.6601 5.6601 8.4238 9.6301 9.6301 10.5192
With DFT Tools it can be plotted as easy as is following script:
from dfttools.simple import parse
from dfttools import presentation
from matplotlib import pyplot
with open("plot.py.data",'r') as f:
# Read bands data
bands = parse(f, "band-structure")
# Plot bands
presentation.matplotlib_bands(bands,pyplot.gca())
pyplot.show()
(Source code, png, hires.png, pdf)

Not only the band structure can be plotted, but atomic structure, data on the grid, etc., see examples.
Getting started¶
DFT Tools package is written in python. To be able to use it you have to download and install it locally.
Installing¶
The easiest way to install DFT Tools is to use pip:
$ pip install dfttools
For a local user it can be done with a --user
option:
$ pip install dfttools --user
You may also download the package and use the bundled setup.py:
$ python setup.py install
$ python setup.py install --user
The package explicitly requires numpy and numericalunits which will be automatically installed if not yet present in your system. Also, it is recommended to install matplotlib and svgwrite for data visualisation and scipy to be able to use some other functions. All packages are available through pip:
$ pip install matplotlib
$ pip install svgwrite
$ pip install scipy
Using¶
Once installed you may start using it by importing the package in your python script:
import dfttools
or just using one of the pre-set scripts:
$ dft-plot-bands my_dft_output_file
Examples¶
To run the examples you have to install all recommended packages, see corresponding section.
Atomic structure¶
With DFT Tools you can manipulate crystal structures easily: only very few lines of code required.
Example: Si unit cell¶
from dfttools.types import Basis, UnitCell
from dfttools.presentation import svgwrite_unit_cell
from numericalunits import angstrom as a
si_basis = Basis((3.9*a/2, 3.9*a/2, 3.9*a/2, .5,.5,.5), kind = 'triclinic')
si_cell = UnitCell(si_basis, (.5,.5,.5), 'Si')
svgwrite_unit_cell(si_cell, 'output.svg', size = (440,360), show_cell = True)
One can obtain a supercell by repeating the unit cell:
mult_cell = si_cell.repeated(3,3,3)
svgwrite_unit_cell(mult_cell, 'output2.svg', size = (440,360), show_cell = True)
Arbitrary supercell is available via the corresponding function:
cubic_cell = si_cell.supercell(
(1,-1,1),
(1,1,-1),
(-1,1,1),
)
svgwrite_unit_cell(cubic_cell, 'output3.svg', size = (440,360), show_cell = True, camera = (1,1,1))
A slab is prepared easily:
slab_cell = cubic_cell.repeated(5,5,3).isolated(0,0,10*a)
svgwrite_unit_cell(slab_cell, 'output4.svg', size = (440,360), camera = (1,1,1))
Example: Monolayer MoS2 line defect¶
A more complex example: monolayer MoS2 with a line defect:
from dfttools.types import Basis, UnitCell
from dfttools.presentation import svgwrite_unit_cell
from numericalunits import angstrom as a
mos2_basis = Basis(
(3.19*a, 3.19*a, 20*a, 0,0,.5),
kind = 'triclinic'
)
d = 1.57722483162840/20
# Unit cell with 3 atoms
mos2_cell = UnitCell(mos2_basis, (
(1./3,1./3,.5),
(2./3,2./3,0.5+d),
(2./3,2./3,0.5-d),
), ('Mo','S','S'))
# Rectangular supercell with 6 atoms
mos2_rectangular = mos2_cell.supercell(
(1,0,0),
(-1,2,0),
(0,0,1)
)
# Rectangular sheet with a defect
mos2_defect = mos2_rectangular.normalized()
mos2_defect.discard((mos2_defect.values == "S") * (mos2_defect.coordinates[:,1] < .5) * (mos2_defect.coordinates[:,2] < .5))
# Prepare a sheet
mos2_sheet = UnitCell.stack(*((mos2_rectangular,)*3 + (mos2_defect,) + (mos2_rectangular,)*3), vector = 'y')
# Draw
svgwrite_unit_cell(mos2_sheet.repeated(10,1,1), 'output.svg', size = (440,360), camera = (1,1,0.3), camera_top = (0,0,1))
Example: parsing structure data¶
It is also possible to obtain atomic structure from the supported files. In this particular case the file source and format can be determined automatically (OpenMX input file).
from dfttools.presentation import svgwrite_unit_cell
from dfttools.simple import parse
# Parse
with open("plot.py.data", "r") as f:
cell = parse(f, "unit-cell")
# Draw
svgwrite_unit_cell(cell, 'output.svg', size = (440,360), camera = (1,0,0))
Example: Moire pattern¶
The Moire pattern is obtained using UnitCell.supercell
.
from dfttools.types import Basis, UnitCell
from dfttools.presentation import svgwrite_unit_cell
from numericalunits import angstrom as a
graphene_basis = Basis(
(2.46*a, 2.46*a, 6.7*a, 0,0,.5),
kind = 'triclinic'
)
# Unit cell
graphene_cell = UnitCell(graphene_basis, (
(1./3,1./3,.5),
(2./3,2./3,.5),
), ('C','C'))
# Moire matching vectors
moire = [1, 26, 6, 23]
# A top layer
l1 = graphene_cell.supercell(
(moire[0],moire[1],0),
(-moire[1],moire[0]+moire[1],0),
(0,0,1)
)
# A bottom layer
l2 = graphene_cell.supercell(
(moire[2],moire[3],0),
(-moire[3],moire[2]+moire[3],0),
(0,0,1)
)
# Make the basis fit
l2.vectors[:2] = l1.vectors[:2]
# Draw
svgwrite_unit_cell(l1.stack(l2, vector='z'), 'output.svg', size = (440,360), camera = (0,0,-1), camera_top = (0,1,0), show_atoms = False)
Band structure¶
The band structures can be easily plotted directly from the output files.
Example: OpenMX¶
In this case to retrieve the band structure we import parser
dfttools.parsers.openmx.bands
explicitly.
from dfttools.parsers.openmx import bands
from dfttools import presentation
from matplotlib import pyplot
with open("plot.py.data",'r') as f:
# Read bands data
b = bands(f.read()).bands()
# Plot bands
presentation.matplotlib_bands(b,pyplot.gca())
pyplot.show()
(Source code, png, hires.png, pdf)

Example: Quantum Espresso¶
The Quantum Espresso files can be identified automatically via
dfttools.simple.parse
routine.
from dfttools.simple import parse
from dfttools import presentation
from matplotlib import pyplot
with open("plot.py.data",'r') as f:
# Read bands data
bands = parse(f, "band-structure")
# Plot bands
presentation.matplotlib_bands(bands,pyplot.gca())
pyplot.show()
(Source code, png, hires.png, pdf)

The density of states can be plotted directly from the band structure. However, one has to note that the density calculated from a k-point path is usually not the relevant one.
from dfttools.simple import parse
from dfttools import presentation
from matplotlib import pyplot
with open("plot.py.data",'r') as f:
# Read bands data
bands = parse(f, "band-structure")
# Prepare axes
ax_left = pyplot.subplot2grid((1,3), (0, 0), colspan=2)
ax_right = pyplot.subplot2grid((1,3), (0, 2))
# Plot bands
presentation.matplotlib_bands(bands,ax_left)
presentation.matplotlib_bands_density(bands, ax_right, 100, orientation = 'portrait')
ax_right.set_ylabel('')
pyplot.show()
(Source code, png, hires.png, pdf)

Example: Density of states¶
To plot an accurate density of states (DoS) a large enough grid is required. Following is an example of a density of states of graphene.
from dfttools.types import Basis, Grid
from dfttools import presentation
from matplotlib import pyplot
from numericalunits import eV
import numpy
# A reciprocal basis
basis = Basis((1,1,1,0,0,-0.5), kind = 'triclinic', meta = {"Fermi": 0})
# Grid shape
shape = (50,50,1)
# A dummy grid with correct grid coordinates and empty grid values
grid = Grid(
basis,
tuple(numpy.linspace(0,1,x, endpoint = False)+.5/x for x in shape),
numpy.zeros(shape+(2,), dtype = numpy.float64),
)
# Calculate graphene band
k = grid.cartesian()*numpy.pi/3.**.5*2
e = (1+4*numpy.cos(k[...,1])**2 + 4*numpy.cos(k[...,1])*numpy.cos(k[...,0]*3.**.5))**.5*eV
# Set the band values
grid.values[...,0] = -e
grid.values[...,1] = e
presentation.matplotlib_bands_density(grid, pyplot.gca(), 200, energy_range = (-1, 1))
pyplot.show()
(Source code, png, hires.png, pdf)

Example: K-point grids: density of states and interpolation¶
They key point of presenting the density of states from a file is
converting the band structure to grid via UnitCell.as_grid
. This
only works if you indeed calculated band energies on a grid. Note
that while both Grid
and UnitCell
can be used for DoS, the
former one is considerably more accurate.
from dfttools.simple import parse
from dfttools import presentation
from matplotlib import pyplot
with open("plot.py.data",'r') as f:
# Retrieve the last band structure from the file
bands = parse(f, "band-structure")
# Convert to a grid
grid = bands.as_grid()
# Plot both
presentation.matplotlib_bands_density(bands, pyplot.gca(), 200, energy_range = (-2, 2), label = "bands")
presentation.matplotlib_bands_density(grid, pyplot.gca(), 200, energy_range = (-2, 2), label = "grid")
pyplot.legend()
pyplot.show()
(Source code, png, hires.png, pdf)

One can also plot the bands by interpolating data on the grid. The quality of the figure depends on the grid size and interpolation methods.
from dfttools.simple import parse
from dfttools import presentation
import numpy
from matplotlib import pyplot
with open("plot.py.data",'r') as f:
# Retrieve the last band structure from the file
bands = parse(f, "band-structure")
# Convert to a grid
grid = bands.as_grid()
# Interpolate
kp_path = numpy.linspace(0,1)[:,numpy.newaxis] * ((1./3,2./3,0),)
bands = grid.interpolate_to_cell(kp_path)
# Plot
presentation.matplotlib_bands(bands, pyplot.gca())
pyplot.show()
(Source code, png, hires.png, pdf)

Example: Band structure with weights¶
The band structure with weights is plotted using weights
keyword
argument. The weights array is just numbers assigned to each k-point and
each band.
from dfttools.types import Basis, UnitCell
from dfttools import presentation
from matplotlib import pyplot
from numericalunits import eV
import numpy
# A reciprocal basis
basis = Basis((1,1,1,0,0,-0.5), kind = 'triclinic', meta = {"Fermi": 0})
# G-K path
kp = numpy.linspace(0,1,100)[:,numpy.newaxis] * numpy.array(((1./3,2./3,0),))
# A dummy grid UnitCell with correct kp-path
bands = UnitCell(
basis,
kp,
numpy.zeros((100,2), dtype = numpy.float64),
)
# Calculate graphene band
k = bands.cartesian()*numpy.pi/3.**.5*2
e = (1+4*numpy.cos(k[...,1])**2 + 4*numpy.cos(k[...,1])*numpy.cos(k[...,0]*3.**.5))**.5*eV
# Set the band values
bands.values[...,0] = -e
bands.values[...,1] = e
# Assign some weights
weights = bands.values.copy()
weights -= weights.min()
weights /= weights.max()
# Prepare axes
ax_left = pyplot.subplot2grid((1,3), (0, 0), colspan=2)
ax_right = pyplot.subplot2grid((1,3), (0, 2))
# Plot bands
p = presentation.matplotlib_bands(bands,ax_left,weights = weights)
presentation.matplotlib_bands_density(bands, ax_right, 100, orientation = 'portrait')
presentation.matplotlib_bands_density(bands, ax_right, 100, orientation = 'portrait', weights = weights, use_fill = True, color = "#AAAAFF")
ax_right.set_ylabel('')
pyplot.colorbar(p)
pyplot.show()
(Source code, png, hires.png, pdf)

Data on the grid¶
Plotting of data (charge, potential, density, etc.) on a 3D grid is very straightforward.
from dfttools.types import Basis, Grid
from dfttools import presentation
from numericalunits import angstrom
from matplotlib import pyplot
import numpy
grid = Grid(
Basis((1*angstrom,1*angstrom,1*angstrom,0,0,-0.5), kind = 'triclinic'),
(
numpy.linspace(0,1,30,endpoint = False),
numpy.linspace(0,1,30,endpoint = False),
numpy.linspace(0,1,30,endpoint = False),
),
numpy.zeros((30,30,30)),
)
grid.values = numpy.prod(numpy.sin(grid.explicit_coordinates()*2*numpy.pi), axis = -1)
presentation.matplotlib_scalar(grid, pyplot.gca(), (0.1,0.1,0.1), 'z', show_cell = True)
pyplot.show()
(Source code, png, hires.png, pdf)

Package contents¶
dfttools
¶
dfttools.formatters
¶
This submodule contains routines presenting data (unit cell) in various text formats.
-
dfttools.formatters.
openmx_input
(cell, populations, l=None, r=None, tolerance=1e-10, indent=4)[source]¶ Generates OpenMX minimal input file with atomic structure.
Args:
cell (UnitCell): input unit cell;
populations (dict): a dict with initial electronic populations data;
Kwargs:
l (UnitCell): left lead;
r (UnitCell): right lead;
tolerance (float): tolerance for checking whether left-center-right unit cells can be stacked;
indent (int): size of indent;
Returns:
String with OpenMX input file formatted data.
-
dfttools.formatters.
pyscf_cell
(cell, **kwargs)[source]¶ Constructs a unit cell object in pyscf.
Args:
cell (UnitCell): a unit cell object to convert from;Kwargs are passed to pyscf.pbc.gto.M.
Returns:
A Pyscf Cell object.
-
dfttools.formatters.
qe_input
(cell=None, relax_triggers=0, parameters={}, inline_parameters={}, pseudopotentials={}, indent=4)[source]¶ Generates Quantum Espresso input file.
Kwargs:
cell (UnitCell): a unit cell with atomic coordinates;
relax_triggers (array,int): array with triggers for relaxation written asadditional columns in the input file;
parameters (dict): parameters for the input file;
inline_parameters (dict): a dict of inline parameters such as
crystal_b
, etc;pseudopotentials (dict): a dict of pseudopotential file names;
indent (int): size of indent;
Returns:
String contating Quantum Espresso input file contents.
-
dfttools.formatters.
siesta_input
(cell, indent=4)[source]¶ Generates Siesta minimal input file with atomic structure.
Args:
cell (UnitCell): input unit cell;Kwargs:
indent (int): size of indent;Returns:
String with Siesta input file contents.
dfttools.presentation
¶
This submodule contains data visualization routines.
-
dfttools.presentation.
matplotlib2svgwrite
(fig, svg, insert, size, **kwargs)[source]¶ Saves a matplotlib image to an existing svgwrite object.
Args:
fig (matplotlib.figure.Figure): a figure to save;
svg (svgwrite.Drawing): an svg drawing to save to;
insert (tuple): a tuple of ints defining destination to insert a drawing;
size (tuple): size of the inserted image;
Kwargs:
The kwargs are passed tofig.savefig
used to print the plot.
-
dfttools.presentation.
matplotlib_bands
(cell, axes, show_fermi=True, energy_range=None, energy_units='eV', energy_units_name=None, coordinate_units=None, coordinate_units_name=None, threshold=0.01, weights=None, weights_color=None, weights_size=None, optimize_visible=False, edge_names=[], mark_points=None, project=None, **kwargs)[source]¶ Plots basic band structure using pyplot.
Args:
cell (UnitCell): cell with the band structure;
axes (matplotlib.axes.Axes): axes to plot on;
Kwargs:
show_fermi (bool): shows the Fermi level if specified;
energy_range (array): 2 floats defining plot energy range. The units of energy are defined by the
units
keyword;energy_units (str, float): either a field from
numericalunits
package or a float with energy units;energy_units_name (str): a string used for the units. Used only if the
energy_units
keyword is a float;coordinate_units (str, float): either a field from
numericalunits
package or a float with coordinate units or None;coordinate_units_name (str): a string used for the coordinate units. Used only if the
coordinate_units
keyword is a float;threshold (float): threshold for determining edges of k point path;
weights, weights_color (array): a 2D array with weights on the band structure which will be converted to color according to current colormap;
weights_size (array): a 2D array with weights on the band structure which will be converted to line thickness;
optimize_visible (bool): draw only visible lines;
edge_names (list): the edges names to appear on the band structure;
mark_points (list): marks specific points on the band structure, the first number in each list element is interpreted as k-point while the second number is band number;
project (array): projects k-points along specified direction instead of unfolding the entire bands path. If
coordinate_units
specified the direction is expressed in the unit cell vectors, otherwise cartesian basis is used;The rest of kwargs are passed to
matplotlib.collections.LineCollection
.Returns:
A plotted LineCollection.
-
dfttools.presentation.
matplotlib_bands_density
(cell, axes, energies, show_fermi=True, energy_range=None, units='eV', units_name=None, weights=None, on_top_of=None, use_fill=False, orientation='landscape', gaussian_spread=None, method='default', **kwargs)[source]¶ Plots density of bands (density of states).
The cell values are considered to be band energies.
Args:
cell (Grid,UnitCell): a unit cell with the band structure, possibly on the grid;
axes (matplotlib.axes.Axes): axes to plot on;
energies (int,array): energies to calculate density at. The integer value has the meaning of number of points to cover the range
energy_range
. Otherwise the units of energy are defined by theunits
keyword;Kwargs:
show_fermi (bool): shows the Fermi level if specified;
energy_range (array): 2 floats defining plot energy range. The units of energy are defined by the
units
keyword;units (str, float): either a field from
numericalunits
package or a float with energy units;units_name (str): a string used for the units. Used only if the
units
keyword is a float;weights (array): a 2D array with weights on the band structure;
on_top_of (array): a 2D array with weights on the band structure to plot on top of;
use_fill (bool): fill the area below plot;
orientation (str): either ‘portrait’ or ‘landscape’ - orientation of the plot;
gaussian_spread (float): the gaussian spread for the density of states. This value is used only if the provided
cell
is not a Grid;method (bool): method to calculate density: ‘default’, ‘gaussian’ or ‘optimal’;
The rest of kwargs are passed to pyplot plotting functions.
Returns:
A plotted Line2D or a PolyCollection, depending onuse_fill
.
-
dfttools.presentation.
matplotlib_scalar
(grid, axes, origin, plane, units='angstrom', units_name=None, show_cell=False, normalize=True, ppu=None, isolines=None, window=None, margins=0.1, scale_bar=None, scale_bar_location=1, **kwargs)[source]¶ Plots scalar values on the grid using imshow.
Args:
grid (Grid): a 3D grid to be plotted;
axes (matplotlib.axes.Axes): axes to plot on;
origin (array): origin of the 2D slice to be plotted in the units of
grid
;plane (str, int): the plotting plane: either ‘x’,’y’ or ‘z’ or a correspondint int.
Kwargs:
units (str, float): either a field from
numericalunits
package or a float with energy units;units_name (str): a string used for the units. Used only if the
units
keyword is a float;show_cell (bool): if True then projected unit cell boundaries are shown on the final image;
normalize (bool): normalize data before plotting such that the minimum is set at zero and the maximum is equal to one;
ppu (float): points per
unit
for the raster image;isolines (array): plot isolines at the specified levels;
window (array): 4 values representing a window to plot the data: minimum and maximum ‘x’ coordinate and minimum and maximum ‘y’ coordinate;
margins (float): adds margins to the grid where the data is interpolated;
scale_bar (int): adds a scal bar to the image at the specified location;
scale_bar_location (int): location of the scale bar;
The rest of kwargs are passed to
pyplot.imshow
orpyplot.contour
.Returns:
Amatplotlib.image.AxesImage
plotted.
-
dfttools.presentation.
svgwrite_unit_cell
(cell, svg, camera=None, camera_top=None, insert=(0, 0), size=(600, 600), circle_size=0.4, circle_opacity=None, margin=6, show_cell=False, show_atoms=True, show_bonds=True, show_legend=True, show_numbers=False, fadeout_strength=0.8, bg=(255, 255, 255), bond_ratio=1, hook_atomic_color=None, coordinates='right', invisible=None, title=None)[source]¶ Creates an svg drawing of a unit cell.
Args:
cell (UnitCell): the cell to be visualized;
svg (str, svgwrite.Drawing): either file name to save the drawing to or an
svgwrite.Drawing
object to draw with.Kwargs:
camera (str, array): the direction of a camera: either ‘x’,’y’ or ‘z’ or an arbitrary 3D vector;
camera_top (array): a vector pointing up;
insert (array): a top-left corner of the drawing;
size (array): size of the bounding box;
circle_size (float): size of the circles representing atoms, arbitrary units;
circle_opacity (float,array): opacity of circles;
margin (float): size of the margin in all directions;
show_cell (bool, str): if True draws the unit cell edges projected, if ‘invisible’ the unit cell is invisible;
show_atoms (bool): if True draws atoms;
show_bonds (bool): if True draws bonds;
show_legend (bool): if True draws legend;
show_numbers (bool): if True shows numbers corresponding to the atomic order in the unit cell;
fadeout_strength (float): amount of fadeout applied to more distant atoms;
bg (array): an integer array defining background color;
bond_ratio (float): scale factor to determine whether the bond is rendered;
coordinates (str): the coordinate system, either ‘left’ or ‘right’;
hook_atomic_color (function): a function accepting integer (atom ID) and a 3-element list (suggested RGB color) and returning a new color of the atom;
invisible (str,array): make specified atoms invisible. If ‘auto’ specified, creates a supercell and makes all cell replica invisible. The bonds of invisible atoms will still be present on the final image;
title (str): a title to the drawing presented in the top left corner;
Returns:
An`svgwrite.Drawing
object. The object is saved if it was created inside this method.
dfttools.simple
¶
This submodule contains commonly used shortcuts to parse the data.
-
dfttools.simple.
get_all_parsers
(*modules)[source]¶ Retrieves all parsers.
Kwargs:
modules (list): a list of names ofparsers
submodules to search at.Returns:
A list of parsing classes.
-
dfttools.simple.
guess_parser
(f)[source]¶ Guesses parsers for a given data.
Args:
f (file): a file to parse.Returns:
A list of parser candidates.
dfttools.types
¶
This submodule contains key types for handling coordinate-dependent data:
UnitCell
, Grid
and Basis
.
-
class
dfttools.types.
Basis
(vectors, kind='default', meta=None, units=None)[source]¶ A class describing a set of vectors representing a basis.
Args:
vectors (array): a 2D or a 1D array of floats representing vectors of the basis set.Kwargs:
kind (str): a shortcut keyword for several most common basis sets:
- ‘default’: expects
vectors
to be a 2D array with basis vectors in cartesian coordinates; - ‘orthorombic’: expects
vectors
to be a 1D array with dimensions of an orthorombic basis set; - ‘triclinic’: expects
vectors
to be a 1D array with 3 lengths of edges and 3 cosines of face angles.
units (str,float): optional units for the Basis. The units are stored in self.meta[‘units’] and are used only during save/load process. The string value has to correspond to one of the values in numericalunits package.
meta (dict): a metadata for this Basis.
-
edges
()[source]¶ Computes pairs of cartesian coordinates of all edges of the basis cell.
Returns:
A list of pairs with cartesian coordinates of vertices forming edges.
-
faces
()[source]¶ Computes faces and returns corresponding cartesian coordinates.
Returns:
A list of lists of coordinates defining face polygon coordinates.
-
static
from_json
(j)[source]¶ Restores a Basis from JSON data.
Args:
j (dict): JSON data.Returns:
A Basis object.
-
generate_path
(points, n, anchor=True)[source]¶ Generates a path given key points and the total number of points on the path.
Args:
points (array): key points of the path expressed in this basis;
n (int): the total number of points on the path.
Kwargs:
anchor (bool): force the specified points to be present in the final path. If True alters slightly the total point number;Returns:
Path coordinates expressed in this basis.
-
reciprocal
()[source]¶ Computes a reciprocal basis.
Returns:
A reciprocal basis.Note
The multiplier is not present.
-
reorder_vectors
(*args, **kwargs)[source]¶ Reorders vectors.
Args:
new (array): new mapping of vectors.Example:
>>> basis.reorder_vectors(0, 1, 2) # does nothing >>> basis.reorder_vectors(1, 0, 2) # swaps first and second vectors.
-
repeated
(*args, **kwargs)[source]¶ Produces a new basis from a given by repeating it in all directions.
Args:
times (array): array of ints specifying how much the basis should be repeated along each of the vectors.
-
rotated
(axis, angle, units='rad')[source]¶ Rotates this basis.
Args:
axis (array): axis to rotate around;
angle (float): angle to rotate;
Kwargs:
units (str): units of the angle: ‘rad’, ‘deg’ or ‘frac’.Returns:
A rotated copy of this basis.
-
stack
(*args, **kwargs)[source]¶ Stacks several basises along one of the vectors.
Args:
basises (list): basises to stack. Corresponding vectors of all basises being stacked should match.Kwargs:
vector (str,int): a vector along which to stack, either ‘x’, ‘y’, ‘z’ or an int specifying the vector;
tolerance (float): a largest possible error in input basises’ vectors;
restrict_collinear (bool): if True will raise an exception if the vectors to stack along are not collinear
Raises:
ArgumentError: in the case of vector mismatch.Returns:
A larger basis containing all argument cell stacked.
-
to_json
()[source]¶ Prepares a JSON-compatible object representing this Basis.
Returns:
A JSON-compatible dict.
-
transform_from
(basis, coordinates)[source]¶ Transforms coordinates from another basis set.
Args:
basis (Basis): a basis to transform from.
coordinates (array): an array of coordinates to be transformed.
Returns:
An array with transformed coordinates.
-
transform_from_cartesian
(coordinates)[source]¶ Transforms coordinates from cartesian.
Args:
coordinates (array): an array of coordinates to be transformed.Returns:
An array with transformed coordinates.
-
transform_to
(basis, coordinates)[source]¶ Transforms coordinates to another basis set.
Args:
basis (Basis): a new basis to transform to.
coordinates (array): an array of coordinates to be transformed.
Returns:
An array with transformed coordinates.
-
transform_to_cartesian
(coordinates)[source]¶ Transforms coordinates to cartesian.
Args:
coordinates (array): an array of coordinates to be transformed.Returns:
An array with transformed coordinates.
-
units_aware
()[source]¶ Checks if units for this Basis are defined.
Returns:
True if units were defined.
- ‘default’: expects
-
class
dfttools.types.
Grid
(basis, coordinates, values, units=None)[source]¶ A class describing a data on a grid in a periodic environment.
Args:
basis (Basis): a crystal basis.
coordinates (array): a list of arrays of coordinates specifying grid.
values (array): a multidimensional array with data on the grid.
-
add
(*args, **kwargs)[source]¶ Adds grid points from another grids to this one.
Args:
grids (arguments): grids to be merged with.Returns:
A new grid with merged data.
-
apply
(*args, **kwargs)[source]¶ Applies selection to this grid.
Inverse of
Grid.discard
.Args:
selection (array): seleted grid points.Example:
>>> selection = grid.select((0,0,0,0.5,1,1)) # Selects species in the 'left' part of the grid. >>> grid.apply(selection) # Applies selection. Species outside the 'left' part are discarded.
-
cartesian
()[source]¶ Computes cartesian coordinates.
Returns:
A multidimensional numpy array with cartesian coordinates at each grid point.
-
static
combine_arrays
(arrays)[source]¶ Transforms input 1D arrays of coordinates into (N+1)D mesh array where first N dimensions correspond to a particular grid point and the last dimension specifies all coordinates of this grid point.
Args:
arrays (list): a list of 1D arrays;Returns:
A meshgrid array with coordinates.
-
cut
(*args, **kwargs)[source]¶ Selects a piece of the grid and returns it as a smaller basis.
Kwargs:
piece (array): fraction of the grid to be selected. The order of coordinates inpiece
isx_from, y_from, ..., z_from, x_to, y_to, ..., z_to
.Returns:
A smaller grid selected.
-
discard
(*args, **kwargs)[source]¶ Removes specified points from this grid.
Inverse of
Grid.apply
.Args:
selection (array): points to remove.Example:
>>> selection = grid.select((0,0,0,0.5,1,1)) # Selects points in the 'left' part of the grid. >>> grid.discard(selection) # Discards selection. Points inside the 'left' part are removed.
-
explicit_coordinates
()[source]¶ Creates an (N+1)D array with explicit coordinates at each grid point.
Returns:
An (N+1)D array with coordinates.
-
interpolate_to_array
(points, driver=None, periodic=True, **kwargs)[source]¶ Interpolates values at specified points and returns an array of interpolated values. By default uses
scipy.interpolate.interpn
.Args:
points (array): points to interpolate at.Kwargs:
driver (func): interpolation driver;
periodic (bool): employs periodicity of a unit cell;
The rest of keyword arguments are passed to the driver.
Returns:
An array with values of corresponding shape.
-
interpolate_to_cell
(points, **kwargs)[source]¶ Interpolates values at specified points and returns a unit cell. By default uses
scipy.interpolate.interpn
.Args:
points (array): points to interpolate at.Kwargs are passed to
self.interpolate_to_array
.Returns:
A unit cell interpolated values.
-
interpolate_to_grid
(points, **kwargs)[source]¶ Interpolates values at specified points and returns a grid. By default uses
scipy.interpolate.interpn
.Args:
points (array): points to interpolate at.Kwargs are passed to
self.interpolate_to_array
.Returns:
A grid with interpolated values.
-
interpolate_to_path
(points, n, anchor=True, **kwargs)[source]¶ Interpolates values to a path: performs path generation using
Basis.generate_path
and interpolates values on it.Args:
points (array): key points of the path expressed in lattice coordinates;
n (int): the total number of points on the path.
Kwargs:
anchor (bool): force the specified points to be present in the final path. If True alters slightly the total point number;
kwargs: keywords to the ‘interpolate_to_cell’ routine.
Returns:
A unit cell interpolated values.
-
isolated
(*args, **kwargs)[source]¶ Generates an isolated representation of this grid.
Symmetrically adds vacuum along all basis vectors such that resulting grid basis vectors are parallel to the initial ones.
Args:
gaps (array): size of the vacuum layer in each direction either in cartesian or in crystal units.Kwargs:
units (str): units of the vacuum size, ‘cartesian’ or ‘crystal’Returns:
A new isolated grid.
-
normalized
()[source]¶ Moves all grid points respecting periodicity so that each coordinate becomes in the unit range 0<=x<1 in the cell basis. Sorts the data.
Returns:
A new grid with normalized data.
-
reorder_vectors
(*args, **kwargs)[source]¶ Reorders vectors. Does not change output of
Grid.cartesian()
.Args:
new (array): new mapping of vectors.Example:
>>> grid.reorder_vectors(0, 1, 2) # does nothing >>> grid.reorder_vectors(1, 0, 2) # swaps first and second vectors.
Note
A call to this method does not modify the output of
self.cartesian()
.
-
select
(*args, **kwargs)[source]¶ Selects a piece of this grid.
Args:
piece (array): fraction of the grid to be selected, see examples. The order of coordinates inpiece
isx_from, y_from, ..., z_from, x_to, y_to, ..., z_to
.Returns:
A list of numpy arrays with bools defining whether particular grid point is selected or not.Example:
>>> grid.select((0,0,0,1,1,1)) # select all grid points with coordinates within (0,1) range >>> grid.select(0,0,0,1,1,1) # a simplified version of above >>> grid.select(0,0,0,0.5,1,1) # select the 'left' part >>> grid.select(0.5,0,0,1,1,1) # select the 'right' part
-
size
()[source]¶ Retrieves the total size of points on the grid.
Returns:
Number of species in cell as an integer.
-
stack
(*args, **kwargs)[source]¶ Stacks several grids along one of the vectors.
Args:
grids (list): grids to stack. Corresponding vectors of all grids being stacked should match.Kwargs:
vector (str,int): a vector along which to stack, either ‘x’, ‘y’, ‘z’ or an int specifying the vector.The rest of kwargs are redirected to
Basis.stack
.Raises:
ArgumentError: in the case of vector mismatch.Returns:
A bigger grid containing all argument grids stacked.
-
tetrahedron_density
(points, resolved=False, weights=None)[source]¶ Convolves data to calculate density (of states). Uses the tetrahedron method from PRB 49, 16223 by E. Blochl et al. Works only in a 3D space.
Args:
points (array): values to calculate density at.Kwargs:
resolved (bool): if True returns a spacially and index resolved density. The dimensions of the returned array are
self.values.shape + points.shape
.weights (array): if specified and
resolved
is False convolves result with the specified weights.Returns:
A numpy array containing density: 1D ifresolved == False
or a corresponding Grid ifresolved == True
.
-
static
uniform
(size, endpoint=False)[source]¶ Transform positive integers size into a meshgrid array representing a grid where grid points span uniformly zero to one intervals.
Args:
size (array): an array with positive integers;Kwargs:
endpoint (bool): indicates whether to include x=1 into grids.Returns:
A meshgrid array with coordinates.
-
-
class
dfttools.types.
UnitCell
(basis, coordinates, values, c_basis=None, units=None)[source]¶ A class describing a crystal unit cell in a periodic environment.
Args:
basis (Basis): a crystal basis.
coordinates (array): a 2D array of coordinates of atoms (or any other instances)
values (array): an array of atoms (or any other instances) with the leading dimenstion being the same as the one of
coordinates
array.Kwargs:
c_basis (str,Basis): a Basis for input coordinates or ‘cartesian’ if coordinates are passed in the cartesian basis.-
add
(*args, **kwargs)[source]¶ Adds species from another unit cells to this one.
Args:
cells (arguments): unit cells to be merged with.Returns:
A new unit cell with merged data.
-
angles
(*args, **kwargs)[source]¶ Computes angles between cell specimens.
Args:
ids (array): a set of specimen IDs to compute angles between. Several shapes are accepted:
- nx3 array: computes n cosines of angles [n,0]-[n,1]-[n,2];
- 1D array of length n: computes n-2 cosines of angles [n-1]-[n]-[n+1];
Returns:
A numpy array containing cosines of angles specified.Example:
Following are the valid calls:
>>> cell.angles((0,1,2)) # angle between vectors connecting {second and first} and {second and third} species >>> cell.angles(0,1,2) # a simplified version of above >>> cell.angles(0,1,3,2) # two angles along path: 0-1-3 and 1-3-2 >>> cell.angles(tuple(0,1,3,2)) # same as above >>> cell.angles((0,1,3),(1,3,2)) # same as above
-
apply
(*args, **kwargs)[source]¶ Applies selection to this cell.
Inverse of
UnitCell.discard
.Args:
selection (array): seleted species.Example:
>>> selection = cell.select((0,0,0,0.5,1,1)) # Selects species in the 'left' part of the unit cell. >>> cell.apply(selection) # Applies selection. Species outside the 'left' part are discarded.
-
as_grid
(fill=nan)[source]¶ Converts this unit cell to grid.
Kwargs:
fill: default value to fill withReturns:
A new grid with data from initial cell.
-
cartesian
()[source]¶ Computes cartesian coordinates.
Returns:
A numpy array with cartesian coordinates
-
cut
(*args, **kwargs)[source]¶ Selects a piece of this unit cell and returns it as a smaller unit cell.
Kwargs:
piece (array): fraction of the cell to be selected. The order of coordinates inpiece
isx_from, y_from, ..., z_from, x_to, y_to, ..., z_to
.Returns:
A smaller unit cell selected.
-
discard
(*args, **kwargs)[source]¶ Removes specified species from cell.
Inverse of
Cell.apply
.Args:
selection (array): species to remove.Example:
>>> selection = cell.select((0,0,0,0.5,1,1)) # Selects species in the 'left' part of the unit cell. >>> cell.discard(selection) # Discards selection. Species inside the 'left' part are removed.
-
distances
(*args, **kwargs)[source]¶ Computes distances between species and specified points.
Args:
ids (array): a list of specimen IDs to compute distances between. Several shapes are accepted:
- empty: returns a 2D matrix of all possible distances
- nx2 array of ints: returns n distances between each pair of [n,0]-[n,1] species;
- 1D array of ints of length n: returns n-1 distances between each pair of [n-1]-[n] species;
Returns:
A numpy array containing list of distances.
-
static
from_json
(j)[source]¶ Restores a UnitCell from JSON data.
Args:
j (dict): JSON data.Returns:
A UnitCell object.
-
interpolate
(*args, **kwargs)[source]¶ Interpolates values at specified points. By default uses
scipy.interpolate.griddata
.Args:
points (array): points to interpolate at.Kwargs:
driver (func): interpolation driver.
periodic (bool): employs periodicity of a unit cell.
kwargs: keywords to the driver.
Returns:
A new unit cell with interpolated data.
-
isolated
(*args, **kwargs)[source]¶ Generates an isolated representation of this cell.
Symmetrically adds vacuum along all unit cell vectors such that resulting unit cell vectors are parallel to the initial ones.
Args:
gaps (array): size of the vacuum layer in each direction either in cartesian or in crystal units.Kwargs:
units (str): units of the vacuum size, ‘cartesian’ or ‘crystal’Returns:
A new unit cell with spacially isolated species.
-
isolated2
(gap)[source]¶ Generates an isolated representation of this cell.
The resulting cell is rectangular and contains space gaps of at least “gap” size.
Args:
gap (float): size of the space gap along allself.vectors
.Returns:
A new unit cell with spacially isolated species.
-
normalized
(sort=None)[source]¶ Moves all species respecting periodicity so that each coordinate becomes in the unit range 0<=x<1 in the cell basis. Sorts the data if
sort
provided.Kwargs:
sort: coordinates to sort with: either ‘x’, ‘y’, ‘z’ or 0,1,2.Returns:
A new grid with normalized data.
-
packed
()[source]¶ Moves all species as close to the origin as it is possible. Does not perform translation.
Returns:
A new unit cell with packed coordinates.
-
reorder_vectors
(*args, **kwargs)[source]¶ Reorders vectors.
Args:
new (array): new mapping of vectors.Example:
>>> cell.reorder_vectors(0, 1, 2) # does nothing >>> cell.reorder_vectors(1, 0, 2) # swaps first and second vectors.
Note
A call to this method does not modify the output of
self.cartesian()
.
-
select
(*args, **kwargs)[source]¶ Selects a piece of this cell.
Args:
piece (array): fraction of the cell to be selected, see examples. The order of coordinates inpiece
isx_from, y_from, ..., z_from, x_to, y_to, ..., z_to
.Returns:
A numpy array with bools defining whether particular specimen is selected or not.Example:
>>> cell.select((0,0,0,1,1,1)) # select all species with coordinates within (0,1) range >>> cell.select(0,0,0,1,1,1) # a simplified version of above >>> cell.select(0,0,0,0.5,1,1) # select the 'left' part >>> cell.select(0.5,0,0,1,1,1) # select the 'right' part
-
size
()[source]¶ Retrieves the number of points or species in this unit cell.
Returns:
Number of points or species in cell.
-
species
()[source]¶ Collects number of species of each kind in this cell.
Particularly useful for counting the number of atoms.
Returns:
A dictionary containing species as keys and number of atoms as values.
-
stack
(*args, **kwargs)[source]¶ Stacks several cells along one of the vectors.
Args:
cells (list): cells to stack. Corresponding vectors of all cells being stacked should match.Kwargs:
vector (str,int): a vector along which to stack, either ‘x’, ‘y’, ‘z’ or an int specifying the vector.The rest of kwargs are redirected to
Basis.stack
.Raises:
ArgumentError: in the case of vector mismatch.Returns:
A bigger unit cell containing all argument cell stacked.
-
dfttools.parsers
¶
dfttools.parsers.elk
¶
Parsing ELK files.
-
class
dfttools.parsers.elk.
Bands
(file)[source]¶ Class for parsing band structure from BAND.OUT file.
Args:
data (string): contents of Elk BAND.OUT file
-
class
dfttools.parsers.elk.
Input
(file)[source]¶ Class for parsing elk.in input file.
Args:
data (string): contents of elk.in file
-
class
dfttools.parsers.elk.
Output
(file)[source]¶ Class for parsing INFO.OUT of elk output.
Args:
data (string): contents of INFO.OUT file
-
class
dfttools.parsers.elk.
UnitCellsParser
(file)[source]¶ Class for parsing elk GEOMETRY_OPT.OUT.
Args:
data (string): contents of GEOMETRY_OPT.OUT file
-
dfttools.parsers.elk.
unitcells
¶ alias of
UnitCellsParser
dfttools.parsers.generic
¶
Contains helper routines to parse text.
-
class
dfttools.parsers.generic.
AbstractJSONParser
(data)[source]¶ A root class for JSON parsers.
Args:
data (str): text representation of JSON to parse.
-
class
dfttools.parsers.generic.
AbstractParser
(file)[source]¶ A root class for text parsers.
Args:
data (str): text to parse or a file to read.
-
class
dfttools.parsers.generic.
StringParser
(string)[source]¶ Simple parser for a string with position memory.
This class can be used to parse words, numbers, floats and arrays from a given string. Based on re, it provides the basic functionality for the rest of parsing libraries.
Args:
string (str): input string to be parsed.Note
The input string can be further accessed by
self.string
field. The contents of the string is not copied.-
closest
(exprs)[source]¶ Returns the closest match of a set of expressions.
Args:
exprs (list): a set of expressions being matched.Returns:
Index of the closest expression. The distance is measured to the beginnings of matches. Returns None if none of expressions matched.Example:
>>> sp = StringParser("This is a large string") >>> sp.closest(("a","string","this")) 2
-
distance
(expression, n=1, default=None)[source]¶ Calculates distance to nth occurrence of expression in characters.
Args:
expression (str,re.RegexObject): expression to match. If expression is str then the case is ignored.Kwargs:
n (int): consequetive number of expression to calculate distance to;
default: return value if StopIteration occurs. Ignored if None.
Returns:
Numbers of characters between caret position and nth occurrence of expression or default if too few occurrences found.Raises:
StopIteration: No occurrences left in the string.
-
floatAfter
(after, n=None)[source]¶ Reads floats from string after the next regular expression. Returns the caret to initial position. Particularly useful for getting value for parameter name.
Args:
after (re.RegexObject) - pattern to skip;Kwargs:
n (array,int,str,re.RegexObject): specifies either shape of the numpy array returned or the regular expression to stop matching before;Returns:
If n is specified returns a numpy array of a given shape filled with floats from string. Otherwise returns a single float.Raises:
StopIteration: Not enough floats left in the string.Example:
>>> sp = StringParser("apples = 3.4; bananas = 7") >>> sp.floatAfter("bananas") 7.0 >>> sp.floatAfter("apples") 3.4
-
goto
(expression, n=1)[source]¶ Goes to the beginning of nth occurrence of expression in the string.
Args:
expression (str,re.RegexObject): expression to match. If expression is str then the case is ignored.Kwargs:
n (int): number of occurrences to match.Raises:
StopIteration: No occurrences left in the string.
-
intAfter
(after, n=None)[source]¶ Reads integers from string after the next regular expression. Returns the caret to initial position. Particularly useful for getting value for parameter name.
Args:
after (re.RegexObject) - pattern to skip;Kwargs:
n (array,int,str,re.RegexObject): specifies either shape of the numpy array returned or the regular expression to stop matching before;Returns:
If n is specified returns a numpy array of a given shape filled with integers from string. Otherwise returns a single int.Raises:
StopIteration: Not enough integers left in the string.Example:
>>> sp = StringParser("cows = 3, rabbits = 5") >>> sp.intAfter("rabbits") 5 >>> sp.intAfter("cows") 3
-
matchAfter
(after, match, n=None)[source]¶ Matches pattern after another pattern and returns caret to initial position. Particularly useful for getting value for parameter name. Supports matching arrays via keyword parameter n.
Args:
after (re.RegexObject): pattern to skip;
match (re.RegexObject): pattern to match;
Kwargs:
n (array,int,str,re.RegexObject): specifies either shape of the numpy array returned or the regular expression to stop matching before;Returns:
If n is specified returns a numpy array of a given shape filled with matches from string. Otherwise returns a single match.Raises:
StopIteration: Not enough matches left in the string.The function is equal to
>>> sp = StringParser("Some string") >>> sp.save() >>> sp.skip(after) >>> result = sp.nextMatch(match, n = n) >>> sp.pop()
-
nextFloat
(n=None)[source]¶ Reads floats from string.
Kwargs:
n (array,int,str,re.RegexObject): specifies either shape of the numpy array returned or the regular expression to stop matching before;Returns:
If n is specified returns a numpy array of a given shape filled with floats from string. Otherwise returns a single float. The caret is put behind the last float read.Raises:
StopIteration: Not enough floats left in the string.Example:
>>> sp = StringParser("1.9 2.8 3.7 56.2E-2 abc") >>> sp.nextFloat(2) array([ 1.9, 2.8]) >>> sp.nextFloat("abc") array([ 3.7 , 0.562])
-
nextInt
(n=None)[source]¶ Reads integers from string.
Kwargs:
n (array,int,str,re.RegexObject): specifies either shape of the numpy array returned or the regular expression to stop matching before;Returns:
If n is specified returns a numpy array of a given shape filled with integers from string. Otherwise returns a single int. The caret is put behind the last integer read.Raises:
StopIteration: Not enough integers left in the string.Example:
>>> sp = StringParser("1 2 3 4 5 6 7 8 9 abc 10") >>> sp.nextInt((2,3)) array([[1, 2, 3], [4, 5, 6]]) >>> sp.nextInt("abc") array([ 7., 8., 9.])
-
nextLine
(n=None)[source]¶ Reads lines from string.
Kwargs:
n (array,int,str,re.RegexObject): specifies either shape of the numpy array returned or the regular expression to stop matching before;Returns:
If n is specified returns a numpy array of a given shape filled with lines from string. Otherwise returns a single line. The caret is put behind the last line read.Raises:
StopIteration: Not enough lines left in the string.
-
nextMatch
(match, n=None)[source]¶ Basic function for matching data.
Args:
match (re.RegexObject): object to match;Kwargs:
n (array,int,str,re.RegexObject): specifies either shape of the numpy array returned or the regular expression to stop matching before;Returns:
If n is specified returns a numpy array of a given shape filled with matches from string. Otherwise returns a single match. The caret is put behind the last match.Raises:
StopIteration: Not enough matches left in the string.
-
pop
()[source]¶ Returns to the previously saved position of the parser.
Raises:
IndexError: No saved positions left.
-
present
(expression)[source]¶ Test the string for the presence of expression.
Args:
expression (str,re.RegexObject): expression to match. If expression is str then the case is ignored.Returns:
True if expression is matched to the right of current position of the caret.
-
save
()[source]¶ Saves the current position of the parser.
Example:
sp = StringParser("A very important integer 123 describes something.") sp.skip("very") # The caret is set to the right of "very" sp.save() # The caret position is saved sp.skip("describes") # The caret is set to the right of "describes" # Now the call to StringParser.nextInt() will yield StopIteration. # To return the caret to the previously saved position # StringParser.pop() is used. sp.pop() # Now it is possible to read the integer sp.nextInt()
-
skip
(expression, n=1)[source]¶ Skips n occurrences of expression in the string.
Args:
expression (str,re.RegexObject): expression to match. If expression is str then the case is ignored.Kwargs:
n (int): number of occurrences to skip.Raises:
StopIteration: No occurrences left in the string.
-
-
dfttools.parsers.generic.
parse
¶ alias of
StringParser
dfttools.parsers.openmx
¶
Parsing OpenMX files.
-
class
dfttools.parsers.openmx.
Bands
(file)[source]¶ Class for parsing band structure from openmx.Band file.
Args:
data (string): contents of OpenMX Band file-
bands
()[source]¶ Retrieves bands.
Returns:
A UnitCell object with band energies.Note
This method can be shortcut
dfttools.simple.parse(file,"band-structure")
.
-
-
class
dfttools.parsers.openmx.
Input
(file)[source]¶ Class for parsing parameter values from OpenMX input files.
Args:
data (str): contents of OpenMX input file-
getFloat
(parameter)[source]¶ A shortcut to parser.StringParser.matchAfter designed to obtain float parameter values from textual configuration files.
Args:
parameter (str): parameter nameReturns:
parameter value
-
getInt
(parameter)[source]¶ A shortcut to parser.StringParser.matchAfter designed to obtain integer parameter values from textual configuration files.
Args:
parameter (str): parameter nameReturns:
parameter value
-
getNonSpaced
(parameter)[source]¶ A shortcut to parser.StringParser.matchAfter designed to obtain parameter values without spaces from textual configuration files.
Args:
parameter (str): parameter nameReturns:
parameter value
-
getWord
(parameter)[source]¶ A shortcut to parser.StringParser.matchAfter designed to obtain word-like parameter values from textual configuration files.
Args:
parameter (str): parameter nameReturns:
parameter value
-
unitCell
(l=None, r=None, tolerance=1e-12)[source]¶ Retrieves atomic position data.
Kwargs:
l,r (UnitCell): left lead and right lead cells. This information is required for parsing the cell from NEGF calculation input file;
tolerance (float): a tolerance for comparing atomic position data from the keywords and from the file itself in
aBohr
.Returns:
An input unit cell.Raises:
ValueError: left and right lead cells are not specified for NEGF input file.Note
This method can be shortcut
dfttools.simple.parse(file,"unit-cell")
.
-
-
class
dfttools.parsers.openmx.
JSON_DOS
(data)[source]¶ Parses JSON with OpenMX density of states.
Args:
data (str): contents of OpenMX JSON DoS file.
-
class
dfttools.parsers.openmx.
Output
(file)[source]¶ Class for parsing parameter values from OpenMX output files.
Args:
data (string): contents of OpenMX output file-
convergence
()[source]¶ Retrieves convergence error values.
Returns:
A numpy array of convergence errors.
-
neutral_charge
()[source]¶ Retrieves the number of valence electrons in the calculation for the charge neutral system.
Returns:
The number of electrons.
-
populations
()[source]¶ Retrieves Mulliken populations during scf process.
Returns:
A numpy array where the first index corresponds to iteration number and the second one is atomic ID. The populations are renormalized to reproduce the total charge.
-
total
()[source]¶ Retrieves total energy calculated.
Returns:
An array of floats with total energy per each SCF cycle.
-
unitCells
(startingCell, noraise=False)[source]¶ Retrieves atomic positions data for relax calculation.
Args:
startingCell (qetools.cell.Cell): a unit cell from the input file. It is required since no chemical captions are written in the output.Kwargs:
noraise (bool): retirieves as much structures as possible without raising exceptions.Returns:
A set of all unit cells found.
-
-
class
dfttools.parsers.openmx.
Transmission
(file)[source]¶ Class for parsing transmission from openmx.tran file.
Args:
data (string): contents of openmx.tran file
-
dfttools.parsers.openmx.
joint_populations
(files)[source]¶ Collects several files with Lowdin populations and parses them into a single object.
Args:
files (list of str): file contents in an array.Returns:
A dict with joint JSON data.
-
dfttools.parsers.openmx.
populations
(s)[source]¶ Parses JSON with Lowdin populations. Replaces corresponding arrays with numpy objects. Adds units where needed.
Args:
s (str): file contents.Returns:
A dict with JSON data.
-
dfttools.parsers.openmx.
transmission
¶ alias of
Transmission
dfttools.parsers.qe
¶
Parsing Quantum Espresso files.
-
class
dfttools.parsers.qe.
Bands
(file)[source]¶ Class for parsing output files created by bands.x binary of Quantum Espresso package.
Args:
data (str): string with the contents of the bands.x output file.-
bands
(basis)[source]¶ Retrieves the bands.
Args:
basis (types.Basis): the reciprocal unit cell of the band structure.Returns:
A unit cell containing band energies.
-
ne
()[source]¶ Retrieves number of bands from the output file header.
Returns:
Integer number of bands.
-
-
class
dfttools.parsers.qe.
Cond
(file)[source]¶ Class for parsing output files created by pwcond.x binary of Quantum Espresso package.
Args:
data (str): string with the contents of the output file.-
transmission
(kind='resolved')[source]¶ Retrives transmission data from pwcond output file.
Kwargs:
- kind (str): either “resolved”, “total”, “states_in” or
“states_out”.
- resolved: retrieves transmisson matrix elements btw pairs of states;
- total: retrieves total transmission for all incoming states;
- states_in, states_out: retrieves only incoming or outgoing states without forming pairs and obtaining transmissions.
Warning
The “resolved” mode essentially picks all transmission matrix elements available in the input. Therefore it will not record incident states without corresponding outgoing states. However these states will show up in “total” regime with zero transmission.
Returns:
A numpy array of records with states and transmission:
- energy (float): energy of the state in eV;
- kx,ky (float): x and y components of k-vectors;
- incoming (float): z component of k-vector of incoming state (only for kind == “resolved” or kind == “total” or kind == “states_in”);
- outgoing (float): z component of k-vector of outgoing state (only for kind == “resolved” or kind == “states_out”);
- transmission (float): corresponding transmission matrix element or total transmission (only for kind == “resolved” or kind == “total”).
The k vector projections are given in units of reciprocal lattice.
-
-
class
dfttools.parsers.qe.
Input
(file)[source]¶ Class for parsing input file for pw.x binary of a Quantum Espresso package.
Args:
data (str): string with the contents of the input file.
-
class
dfttools.parsers.qe.
Output
(file)[source]¶ Class for parsing output files created by pw.x binary of Quantum Espresso package.
Args:
data (str): string with the contents of the output file.-
bands
(index=-1, skipVCRelaxException=False)[source]¶ Retrieves bands.
Kwargs:
index (int or None): index of a band structure or
None
if all band structures need to be parsed. Supports negative indexing.skipVCRelaxException (bool): forces to skip variable cell relaxation exception. In this very special case no reciprocal lattice vectors are provided for the new cells in the output file.
Returns:
A set of Cell objects with bands data stored inCell.values
. Specifically,Cell.values
is a n by m array where n is a number of k points and m is a number of bands.Raises:
Exception: if a variable cell calculation data found.
-
fermi
()[source]¶ Retrieves Fermi energies.
Returns:
A numpy array containing Fermi energies for each MD step.
-
force
()[source]¶ Retrieves total force.
Returns:
A numpy array containing total forces for each self-consistent calculation.
-
routineError
()[source]¶ Checks “error in routine” entry in the file.
Returns:
String with textual information about the error. Returns None if no error recorded.
-
scf_accuracy
()[source]¶ Retrieves scf convergence history.
Returns:
A numpy array containing estimated errors after all scf steps during calculations. The energies are given in eV.
-
scf_failed
()[source]¶ Checks for “convergence NOT achieved” signature.
Returns:
True if the signature is present.
-
scf_steps
()[source]¶ Retrieved number of scf steps.
Returns:
A numpy array containing numbers of consequetive scf steps performed to reach convergences.
-
success
()[source]¶ Checks for success signature in the end of the file.
Returns:
True if the signature is present.
-
threads
()[source]¶ Retrieves the number of MPI threads.
Returns:
A number of MPI threads for this calculation as an integer.
-
time
()[source]¶ Retrieves cpu times.
Returns:
Time stamps measured by Quantum Espresso in a numpy array.
-
-
class
dfttools.parsers.qe.
Proj
(file)[source]¶ Class for parsing output files created by projwfc.x binary of Quantum Espresso package.
Args:
data (str): string with the contents of the output file.-
basis
()[source]¶ Retrieves the localized basis set.
Returns:
A numpy array of records:
- state (int): ID of the state as provided by Quantum Espresso;
- atom (int): ID of particular atom in the unit cell;
- atomName (str): chemical caption;
- wfc (int): particular wave function ID;
- m,l (float): quantum numbers for collinear calculation;
- j,m_j,l (float): quantum numbers for non-collinear calculation.
-
dfttools.parsers.structure
¶
Parsing various atomic structure files.
-
class
dfttools.parsers.structure.
GaussianCube
(file)[source]¶ Class for parsing Gaussian CUBE files.
Args:
data (str): string with the contents of the Gaussian CUBE file.
-
class
dfttools.parsers.structure.
XSF
(file)[source]¶ Class for parsing xcrysden files commonly used in solid state visualisations.
Args:
data (str): string with the contents of the xsf file.
-
class
dfttools.parsers.structure.
XYZ
(file)[source]¶ Class for parsing XYZ structure files.
Args:
data (str): string with the contents of the XYZ file.
-
dfttools.parsers.structure.
cube
¶ alias of
GaussianCube
dfttools.parsers.vasp
¶
Parsing VASP files.