Plotting Fleur DOS/bandstructures

This section discusses how to obtain plots of data in the banddos.hdf for density of states and bandstructure calculations.

In the following bandstructure and DOS plots are explained. Each section leads with the names of the recipes from the recipes module that can be used with the explained visualization function.

All Fleur specific plotting routines are found in fleur have implementations for both the matplotlib and bokeh plotting libraries and can be customized heavily. For an explanation on customizing plots refer to General Plotting routines.

Reading banddos.hdf files

The process here is divided in two parts. First we extract and transform the data in a way to make it easy to plot via the reader.HDF5Reader. For a detailed explanation of the capabilities of this tool refer to General HDF5 file reader. Here we show the basic usage:

from masci_tools.io.parsers.hdf5 import HDF5Reader
from masci_tools.io.parsers.hdf5.recipes import FleurBands

with HDF5Reader('files/banddos_bands.hdf') as h5reader:
   data, attributes = h5reader.read(recipe=FleurBands)

print(f"The following datasets were read: {list(data.keys())}")
print(f"The following attributes were read: {list(attributes.keys())}")
The following datasets were read: ['eigenvalues_up', 'kpath', 'INT_up', 'MT:1d_up', 'MT:1f_up', 'MT:1p_up', 'MT:1s_up', 'Sym_up', 'Total_up', 'MT:1_up']
The following attributes were read: ['group_name', 'kpoints', 'nkpts', 'nbands', 'atoms_elements', 'n_types', 'atoms_position', 'atoms_groups', 'reciprocal_cell', 'bravais_matrix', 'special_kpoint_indices', 'special_kpoint_labels', 'fermi_energy', 'spins']

Bandstructures

Compatible Recipes for the reader.HDF5Reader:

  • FleurBands: Default recipe reading in the kpoints, eignevalues and weights for atom and orbital contributions

  • FleurSimpleBands: Reads in only the kpoints and eigenvalues and now weights

  • FleurOrbcompBands: In addition to the eigenvalues the weights from an orbital decomposition calculation are read in

  • FleurjDOSBands: In addition to the eigenvalues the weights from a jDOS calculation are read in

  • FleurMCDBands: In addition to the eigenvalues the weights from a MCD calculation are read in

  • recipes.get_fleur_bands_specific_weights(): Function to generate a recipe for reading in the eigenvalues+a provided list of weights

The bandstructure visualization plot_fleur_bands() can be used to plot

  1. Non-spinpolarized/spinpolarized bandstructures

  2. Bandstructures with emphasized weights on all eigenvalues (Also non-spinpolarized and spinpolarized)

Standard bandstructure

To plot a simple bandstructure without any weighting we just have to pass the data, that the HDF5Reader provided to the plot_fleur_bands().

The two examples below show the resulting plots for a non-psinpolarized system (bulk Si) and a spin-polarized system (Fe fcc). For both systems the necessary code is exactly the same and is shown above the plots. The shown plots are the ones for the matplotlib plotting backend:

from masci_tools.io.parsers.hdf5 import HDF5Reader
from masci_tools.io.parsers.hdf5.recipes import FleurBands
from masci_tools.vis.fleur import plot_fleur_bands

#Read in data
with HDF5Reader('files/banddos_bands.hdf') as h5reader:
   data, attributes = h5reader.read(recipe=FleurBands)

#Plot the data
#Notice that you get the axis object of this plot is returned
#if you want to make any special additions
ax = plot_fleur_bands(data,
                      attributes,
                      limits={'y': (-13, 5)},
                      markersize=10)
../_images/08d8fabce4a0049d76c0bf525015d8310b80da62b43d25e768ee1c25ee45fad3.png

Non spinpolarized bandstructure for a bulk Si structure

../_images/bb7bcaa300470b15d2db4428cbfd055e34f1fdc9e0ade8d46042cf4c803c7faf.png

Spinpolarized bandstructure for a bulk Fe fcc structure

Bandstructure with weights

To plot a simple bandstructure with weighting we do the same procedure as above, but we pass in the entry we want to use for weights. These correspond to the entries in the banddos.hdf file (for example the weight for the s-orbital on the first atom type is called MT:1s). The weights will be used to change the size and color (according to a colormap) to indicate regions of high weight.

The two examples below show the resulting plots for a non-psinpolarized system (bulk Si) weighted for the s-orbital on the first atom and a spin-polarized system (Fe fcc) with weights for the d-orbital on the first atom type. For both systems the necessary code is exactly the same and is shown above the plots. The shown plots are the ones for the matplotlib plotting backend:

from masci_tools.io.parsers.hdf5 import HDF5Reader
from masci_tools.io.parsers.hdf5.recipes import FleurBands
from masci_tools.vis.fleur import plot_fleur_bands

#Read in data
with HDF5Reader('files/banddos_bands.hdf') as h5reader:
   data, attributes = h5reader.read(recipe=FleurBands)

#Plot the data
#Notice that you get the axis object of this plot is returned
#if you want to make any special additions
ax = plot_fleur_bands(data,
                      attributes,
                      weight='MT:1s',
                      limits={'y': (-13,5)})
../_images/27b30e71e1d278ea2396f3d7887f2600da292710f22c08ee2408143c0ada150a.png

Non spinpolarized bandstructure for a bulk Si structure. The s-like character inside the Muffin-tin sphere is highlighted

../_images/9c7bab2053ff5a850234c296aeff73959e61c60a02c736a6b46168bc1392c400.png

Spinpolarized bandstructure for a bulk Fe fcc structure. The d-like character inside the Muffin-tin sphere is highlighted

Plotting options for bandstructure plots

The plot_fleur_bands() function has a couple of options to modify, what is being displayed from the banddos.hdf file. Below we show a few examples of ways to use these options, together with examples of resulting plots.

Plotting bandstructure without spinpolarization

Providing spinpol=False will display the bandstructure as non spinpolarized, even if there are two spins in the data. Works for both non-weighted and weighted bandstructures.

from masci_tools.io.parsers.hdf5 import HDF5Reader
from masci_tools.io.parsers.hdf5.recipes import FleurBands
from masci_tools.vis.fleur import plot_fleur_bands

#Read in data
with HDF5Reader('files/banddos_spinpol_bands.hdf') as h5reader:
   data, attributes = h5reader.read(recipe=FleurBands)

#Plot the data
#Notice that you get the axis object of this plot is returned
#if you want to make any special additions
ax = plot_fleur_bands(data,
                      attributes,
                      limits={'y': (-9,4)},
                      markersize=10,
                      spinpol=False)
../_images/89f09f131c2744c078c1a56df980c62c46ed70be716752a782ea1602eca96dbb.png

Non spinpolarized bandstructure for a bulk Fe fcc structure.

../_images/d39032799b23b483eef679a54459663a85ba3da9329a6bd225fea1fea9935b22.png

Non spinpolarized bandstructure for a bulk Fe fcc structure. The d-like character inside the Muffin-tin sphere is highlighted

Selecting a specific spin channel

Providing only_spin='up' or 'down' will plot only the given spin channel

from masci_tools.io.parsers.hdf5 import HDF5Reader
from masci_tools.io.parsers.hdf5.recipes import FleurBands
from masci_tools.vis.fleur import plot_fleur_bands

#Read in data
with HDF5Reader('files/banddos_spinpol_bands.hdf') as h5reader:
   data, attributes = h5reader.read(recipe=FleurBands)

#Plot the data
#Notice that you get the axis object of this plot is returned
#if you want to make any special additions
ax = plot_fleur_bands(data,
                      attributes,
                      limits={'y': (-9,4)},
                      only_spin='up',
                      markersize=10,
                      color='darkblue')
../_images/72490d2a43fa874fbaa3c97e1d12ab78cd68eb4170d56f12852bbd0d37c8682f.png

Bandstructure for a bulk Fe fcc structure (only spin up).

../_images/ec90dd6dbecb69406f88491521f8c114a16a06aa720c88fcb09e8589fd494ab7.png

Bandstructure for a bulk Fe fcc structure (only spin down). The color is changed manually

Density of States

Compatible Recipes for the HDF5Reader:

  • FleurDOS: Default recipe reading in the total, interstitial, vacuum, atom and l-channel resolved DOS

  • FleurORBCOMP: Read in the DOS from an orbital decomposition calculation

  • FleurJDOS: Read in the DOS from a jDOS calculation

  • FleurMCD: Read in the DOS from a MCD calculation

The dos visualization plot_fleur_dos() can be used to plot non spinpolarized and spinpolarized DOS, with selection of which components to plot.

Standard density of states plot

from masci_tools.io.parsers.hdf5 import HDF5Reader
from masci_tools.io.parsers.hdf5.recipes import FleurDOS
from masci_tools.vis.fleur import plot_fleur_dos

#Read in data
with HDF5Reader('files/banddos_dos.hdf') as h5reader:
   data, attributes = h5reader.read(recipe=FleurDOS)

#Plot the data
#Notice that you get the axis object of this plot is returned
#if you want to make any special additions
ax = plot_fleur_dos(data,
                    attributes,
                    limits={'energy': (-13,5)})
../_images/d1dc8a7d52483c2f134b07eaf48af4896b821d25761188fdb72f8ac24019d689.png

Non spinpolarized DOS for a bulk Si structure

../_images/a823e46ed019a756d9a2cc96cae20c364597c60eafc9311ff82c544d4829d357.png

Spinpolarized DOS for a bulk Fe fcc structure

Plotting options for DOS plots

The plot_fleur_dos() function has a couple of options to modify, what is being displayed from the banddos.hdf file. Below we show a few examples of ways to use these options, together with examples of resulting plots.

Selecting specific DOS components

The DOS is made up of a lot of contributions that can be displayed separately.

Here we list the options that are available and show example plots for only selecting the atom projected compinents of the density of states

  • plot_keys: Can be used to provide a explicit list of keys you want to display (Same format as in the banddos.hdf)

  • show_total: Control, whether to show the total density of states (default True)

  • show_interstitial: Control, whether to show the interstitial contribution of the density of states (default True)

  • show_atoms: Control, which total atom projected DOS to show. Can be either the string all (All components are shown), the value None (no components are shown) or a list of the integer indices of the atom types that should be displayed (default all)

  • show_lresolved: Control, on which atoms to show the orbital projected DOS. Can be either the string all (All components are shown), the value None (no components are shown) or a list of the integer indices of the atom types for which to display the orbital components (default None)

Below an example of only displaying the atom projected DOS together with their orbital contributions is shown.

from masci_tools.io.parsers.hdf5 import HDF5Reader
from masci_tools.io.parsers.hdf5.recipes import FleurDOS
from masci_tools.vis.fleur import plot_fleur_dos

#Read in data
with HDF5Reader('files/banddos_dos.hdf') as h5reader:
   data, attributes = h5reader.read(recipe=FleurDOS)

ax = plot_fleur_dos(data, attributes,
                    show_total=False,
                    show_interstitial=False,
                    show_lresolved='all',
                    limits={'energy': (-13,5)})
../_images/0429be05bf52c2b1b688c9f6113dce19c2bbc2aa895832de2c4a8f9282f45032.png

Non spinpolarized DOS for a bulk Si structure. Only the atom and l-channel projected DOS is shown

../_images/e0c25d92bd1e4c0c966fcfca9d77b1e7de5e76e6ab3ba48de6c0b325fe132c23.png

Non spinpolarized DOS for a bulk Fe fcc structure. Only the atom and l-channel projected DOS is shown

Plotting DOS without spinpolarization

Providing spinpol=False will display the DOS as non spinpolarized, even if there are two spins in the data.

from masci_tools.io.parsers.hdf5 import HDF5Reader
from masci_tools.io.parsers.hdf5.recipes import FleurDOS
from masci_tools.vis.fleur import plot_fleur_dos

#Read in data
with HDF5Reader('files/banddos_spinpol_dos.hdf') as h5reader:
   data, attributes = h5reader.read(recipe=FleurDOS)

#Plot the data
#Notice that you get the axis object of this plot is returned
#if you want to make any special additions
ax = plot_fleur_dos(data,
                    attributes,
                    limits={'energy': (-9,4)},
                    spinpol=False)
/home/docs/checkouts/readthedocs.org/user_builds/masci-tools/envs/latest/lib/python3.10/site-packages/masci_tools/vis/fleur.py:384: FutureWarning: DataFrame.groupby with axis=1 is deprecated. Do `frame.T.groupby(...)` without axis instead.
  new_dosdata = complete_spin.groupby(complete_spin.columns, axis=1).sum()
../_images/fb7d22e5284bf69fb3ca471d87a3085011f97c1f0767608a6d8127192424938e.png

Non spinpolarized DOS for a bulk Fe fcc structure.

Selecting a specific spin channel

Providing only_spin='up' or 'down' will plot only the given spin channel

from masci_tools.io.parsers.hdf5 import HDF5Reader
from masci_tools.io.parsers.hdf5.recipes import FleurDOS
from masci_tools.vis.fleur import plot_fleur_dos

#Read in data
with HDF5Reader('files/banddos_spinpol_dos.hdf') as h5reader:
   data, attributes = h5reader.read(recipe=FleurDOS)

#Plot the data
#Notice that you get the axis object of this plot is returned
#if you want to make any special additions
ax = plot_fleur_dos(data,
                    attributes,
                    limits={'energy': (-9,4)},
                    only_spin='up')
../_images/c0e680298f3e31ce609270e696818df9804c6b6e69e34d8db048f9c7764b5a3d.png

DOS for a bulk Fe fcc structure (only spin up).

../_images/6bd724512fb1198a70ebffabe4885dd2e9b88ab2d45f0beff23e9b057039a796.png

DOS for a bulk Fe fcc structure (only spin down).