pymatgen.electronic_structure.plotter module

class BSDOSPlotter(bs_projection='elements', dos_projection='elements', vb_energy_range=4, cb_energy_range=4, fixed_cb_energy=False, egrid_interval=1, font='Times New Roman', axis_fontsize=20, tick_fontsize=15, legend_fontsize=14, bs_legend='best', dos_legend='best', rgb_legend=True, fig_size=(11, 8.5))[source]

Bases: object

A joint, aligned band structure and density of states plot. Contributions from Jan Pohls as well as the online example from Germain Salvato-Vallverdu: http://gvallver.perso.univ-pau.fr/?p=587

Instantiate plotter settings.

Parameters:
  • bs_projection (str) – “elements” or None
  • dos_projection (str) – “elements”, “orbitals”, or None
  • vb_energy_range (float) – energy in eV to show of valence bands
  • cb_energy_range (float) – energy in eV to show of conduction bands
  • fixed_cb_energy (bool) – If true, the cb_energy_range will be interpreted as constant (i.e., no gap correction for cb energy)
  • egrid_interval (float) – interval for grid marks
  • font (str) – font family
  • axis_fontsize (float) – font size for axis
  • tick_fontsize (float) – font size for axis tick labels
  • legend_fontsize (float) – font size for legends
  • bs_legend (str) – matplotlib string location for legend or None
  • dos_legend (str) – matplotlib string location for legend or None
  • rgb_legend (bool) – (T/F) whether to draw RGB triangle/bar for element proj.
  • fig_size (tuple) – dimensions of figure size (width, height)
get_plot(bs, dos)[source]

Get a matplotlib plot object. :param bs: the bandstructure to plot. Projection

data must exist for projected plots.
Parameters:dos (Dos) – the Dos to plot. Projection data must exist (i.e., CompleteDos) for projected plots.
Returns:matplotlib.pyplot object on which you can call commands like show() and savefig()
class BSPlotter(bs)[source]

Bases: object

Class to plot or get data to facilitate the plot of band structure objects.

Parameters:bs – A BandStructureSymmLine object.
bs_plot_data(zero_to_efermi=True)[source]

Get the data nicely formatted for a plot

Parameters:zero_to_efermi – Automatically subtract off the Fermi energy from the eigenvalues and plot.
Returns: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 branch) of energy for each kpoint. The data is stored by branch to facilitate the plotting 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).
Return type:A dict of the following format
get_plot(zero_to_efermi=True, ylim=None, smooth=False, vbm_cbm_marker=False, smooth_tol=None)[source]

Get a matplotlib object for the bandstructure plot. Blue lines are up spin, red lines are down spin.

Parameters:
  • 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
  • smooth_tol (float) – tolerance for fitting spline to band data. Default is None such that no tolerance will be used.
get_ticks()[source]

Get all ticks and labels for a band structure plot.

Returns:a list of distance at which ticks should be set and ‘label’: a list of label for each of those ticks.
Return type:A dict with ‘distance’
plot_brillouin()[source]

plot the Brillouin zone

plot_compare(other_plotter, legend=True)[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

Parameters:band structure object defined along the same symmetry lines (another) –
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.

Parameters:
  • 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, smooth_tol=None)[source]

Show the plot using matplotlib.

Parameters:
  • 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
  • smooth_tol (float) – tolerance for fitting spline to band data. Default is None such that no tolerance will be used.
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.

Parameters:bs – A BandStructureSymmLine object with projections.
get_elt_projected_plots(zero_to_efermi=True, ylim=None, vbm_cbm_marker=False)[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)[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

spin up and spin down are differientiated by a ‘-‘ and a ‘–’ line

Parameters: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, vbm_cbm_marker=False)[source]

Method returning a plot composed of subplots along different elements and orbitals.

Parameters: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.
get_projected_plots_dots_patom_pmorb(dictio, dictpa, sum_atoms=None, sum_morbs=None, zero_to_efermi=True, ylim=None, vbm_cbm_marker=False, selected_branches=None, w_h_size=(12, 8), num_column=None)[source]

Method returns a plot composed of subplots for different atoms and orbitals (subshell orbitals such as ‘s’, ‘p’, ‘d’ and ‘f’ defined by azimuthal quantum numbers l = 0, 1, 2 and 3, respectively or individual orbitals like ‘px’, ‘py’ and ‘pz’ defined by magnetic quantum numbers m = -1, 1 and 0, respectively). This is an extension of “get_projected_plots_dots” method.

Parameters:
  • dictio – The elements and the orbitals you need to project on. The format is {Element:[Orbitals]}, for instance: {‘Cu’:[‘dxy’,’s’,’px’],’O’:[‘px’,’py’,’pz’]} will give projections for Cu on orbitals dxy, s, px and for O on orbitals px, py, pz. If you want to sum over all individual orbitals of subshell orbitals, for example, ‘px’, ‘py’ and ‘pz’ of O, just simply set {‘Cu’:[‘dxy’,’s’,’px’],’O’:[‘p’]} and set sum_morbs (see explainations below) as {‘O’:[p],...}. Otherwise, you will get an error.
  • dictpa – The elements and their sites (defined by site numbers) you need to project on. The format is {Element: [Site numbers]}, for instance: {‘Cu’:[1,5],’O’:[3,4]} will give projections for Cu on site-1 and on site-5, O on site-3 and on site-4 in the cell. Attention: The correct site numbers of atoms are consistent with themselves in the structure computed. Normally, the structure should be totally similar with POSCAR file, however, sometimes VASP can rotate or translate the cell. Thus, it would be safe if using Vasprun class to get the final_structure and as a result, correct index numbers of atoms.
  • sum_atoms

    Sum projection of the similar atoms together (e.g.: Cu on site-1 and Cu on site-5). The format is {Element: [Site numbers]}, for instance:

    {‘Cu’: [1,5], ‘O’: [3,4]} means summing projections over Cu on site-1 and Cu on site-5 and O on site-3 and on site-4. If you do not want to use this functional, just turn it off by setting sum_atoms = None.
  • sum_morbs – Sum projections of individual orbitals of similar atoms together (e.g.: ‘dxy’ and ‘dxz’). The format is {Element: [individual orbitals]}, for instance: {‘Cu’: [‘dxy’, ‘dxz’], ‘O’: [‘px’, ‘py’]} means summing projections over ‘dxy’ and ‘dxz’ of Cu and ‘px’ and ‘py’ of O. If you do not want to use this functional, just turn it off by setting sum_morbs = None.
  • selected_branches – The index of symmetry lines you chose for plotting. This can be useful when the number of symmetry lines (in KPOINTS file) are manny while you only want to show for certain ones. The format is [index of line], for instance: [1, 3, 4] means you just need to do projection along lines number 1, 3 and 4 while neglecting lines number 2 and so on. By default, this is None type and all symmetry lines will be plotted.
  • w_h_size – This variable help you to control the width and height of figure. By default, width = 12 and height = 8 (inches). The width/height ratio is kept the same for subfigures and the size of each depends on how many number of subfigures are plotted.
  • num_column – This variable help you to manage how the subfigures are arranged in the figure by setting up the number of columns of subfigures. The value should be an int number. For example, num_column = 3 means you want to plot subfigures in 3 columns. By default, num_column = None and subfigures are aligned in 2 columns.
Returns:

A pylab object with different subfigures for different projections. The blue and red colors lines are bands for spin up and spin down. The green and cyan dots are projections for spin up and spin down. The bigger the green or cyan dots in the projected band structures, the higher character for the corresponding elements and orbitals. List of individual orbitals and their numbers (set up by VASP and no special meaning): s = 0; py = 1 pz = 2 px = 3; dxy = 4 dyz = 5 dz2 = 6 dxz = 7 dx2 = 8; f_3 = 9 f_2 = 10 f_1 = 11 f0 = 12 f1 = 13 f2 = 14 f3 = 15

class BoltztrapPlotter(bz)[source]

Bases: object

class containing methods to plot the data from Boltztrap.

Parameters:bz – a BoltztrapAnalyzer object
plot_carriers(temp=300)[source]

Plot the carrier concentration in function of Fermi level

Parameters:temp – the temperature
Returns:a matplotlib object
plot_conductivity_dop(temps='all', output='average', relaxation_time=1e-14)[source]

Plot the conductivity in function of doping levels for different temperatures.

Parameters:
  • temps – the default ‘all’ plots all the temperatures in the analyzer. Specify a list of temperatures if you want to plot only some.
  • output – with ‘average’ you get an average of the three directions with ‘eigs’ you get all the three directions.
  • relaxation_time – specify a constant relaxation time value
Returns:

a matplotlib object

plot_conductivity_mu(temp=600, output='eig', relaxation_time=1e-14, xlim=None)[source]

Plot the conductivity in function of Fermi level. Semi-log plot

Parameters:
  • temp – the temperature
  • xlim – a list of min and max fermi energy by default (0, and band gap)
  • tau – A relaxation time in s. By default none and the plot is by units of relaxation time
Returns:

a matplotlib object

plot_conductivity_temp(doping='all', output='average', relaxation_time=1e-14)[source]

Plot the conductivity in function of temperature for different doping levels.

Parameters:
  • dopings – the default ‘all’ plots all the doping levels in the analyzer. Specify a list of doping levels if you want to plot only some.
  • output – with ‘average’ you get an average of the three directions with ‘eigs’ you get all the three directions.
  • relaxation_time – specify a constant relaxation time value
Returns:

a matplotlib object

plot_dos(sigma=0.05)[source]

plot dos

Parameters:sigma – a smearing
Returns:a matplotlib object
plot_eff_mass_dop(temps='all', output='average')[source]

Plot the average effective mass in function of doping levels for different temperatures.

Parameters:
  • temps – the default ‘all’ plots all the temperatures in the analyzer. Specify a list of temperatures if you want to plot only some.
  • output – with ‘average’ you get an average of the three directions with ‘eigs’ you get all the three directions.
  • relaxation_time – specify a constant relaxation time value
Returns:

a matplotlib object

plot_eff_mass_temp(doping='all', output='average')[source]

Plot the average effective mass in function of temperature for different doping levels.

Parameters:
  • dopings – the default ‘all’ plots all the doping levels in the analyzer. Specify a list of doping levels if you want to plot only some.
  • output – with ‘average’ you get an average of the three directions with ‘eigs’ you get all the three directions.
Returns:

a matplotlib object

plot_hall_carriers(temp=300)[source]

Plot the Hall carrier concentration in function of Fermi level

Parameters:temp – the temperature
Returns:a matplotlib object
plot_power_factor_dop(temps='all', output='average', relaxation_time=1e-14)[source]

Plot the Power Factor in function of doping levels for different temperatures.

Parameters:
  • temps – the default ‘all’ plots all the temperatures in the analyzer. Specify a list of temperatures if you want to plot only some.
  • output – with ‘average’ you get an average of the three directions with ‘eigs’ you get all the three directions.
  • relaxation_time – specify a constant relaxation time value
Returns:

a matplotlib object

plot_power_factor_mu(temp=600, output='eig', relaxation_time=1e-14, xlim=None)[source]

Plot the power factor in function of Fermi level. Semi-log plot

Parameters:
  • temp – the temperature
  • xlim – a list of min and max fermi energy by default (0, and band gap)
  • tau – A relaxation time in s. By default none and the plot is by units of relaxation time
Returns:

a matplotlib object

plot_power_factor_temp(doping='all', output='average', relaxation_time=1e-14)[source]

Plot the Power Factor in function of temperature for different doping levels.

Parameters:
  • dopings – the default ‘all’ plots all the doping levels in the analyzer. Specify a list of doping levels if you want to plot only some.
  • output – with ‘average’ you get an average of the three directions with ‘eigs’ you get all the three directions.
  • relaxation_time – specify a constant relaxation time value
Returns:

a matplotlib object

plot_seebeck_dop(temps='all', output='average')[source]

Plot the Seebeck in function of doping levels for different temperatures.

Parameters:
  • temps – the default ‘all’ plots all the temperatures in the analyzer. Specify a list of temperatures if you want to plot only some.
  • output – with ‘average’ you get an average of the three directions with ‘eigs’ you get all the three directions.
Returns:

a matplotlib object

plot_seebeck_mu(temp=600, output='eig', xlim=None)[source]

Plot the seebeck coefficient in function of Fermi level

Parameters:
  • temp – the temperature
  • xlim – a list of min and max fermi energy by default (0, and band gap)
Returns:

a matplotlib object

plot_seebeck_temp(doping='all', output='average')[source]

Plot the Seebeck coefficient in function of temperature for different doping levels.

Parameters:
  • dopings – the default ‘all’ plots all the doping levels in the analyzer. Specify a list of doping levels if you want to plot only some.
  • output – with ‘average’ you get an average of the three directions with ‘eigs’ you get all the three directions.
Returns:

a matplotlib object

plot_zt_dop(temps='all', output='average', relaxation_time=1e-14)[source]

Plot the figure of merit zT in function of doping levels for different temperatures.

Parameters:
  • temps – the default ‘all’ plots all the temperatures in the analyzer. Specify a list of temperatures if you want to plot only some.
  • output – with ‘average’ you get an average of the three directions with ‘eigs’ you get all the three directions.
  • relaxation_time – specify a constant relaxation time value
Returns:

a matplotlib object

plot_zt_mu(temp=600, output='eig', relaxation_time=1e-14, xlim=None)[source]

Plot the ZT in function of Fermi level.

Parameters:
  • temp – the temperature
  • xlim – a list of min and max fermi energy by default (0, and band gap)
  • tau – A relaxation time in s. By default none and the plot is by units of relaxation time
Returns:

a matplotlib object

plot_zt_temp(doping='all', output='average', relaxation_time=1e-14)[source]

Plot the figure of merit zT in function of temperature for different doping levels.

Parameters:
  • dopings – the default ‘all’ plots all the doping levels in the analyzer. Specify a list of doping levels if you want to plot only some.
  • output – with ‘average’ you get an average of the three directions with ‘eigs’ you get all the three directions.
  • relaxation_time – specify a constant relaxation time value
Returns:

a matplotlib object

class DosPlotter(zero_at_efermi=True, stack=False, sigma=None)[source]

Bases: object

Class for plotting DOSs. Note that the interface is extremely flexible given that there are many different ways in which people want to view DOS. The typical usage is:

# Initializes plotter with some optional args. Defaults are usually
# fine,
plotter = DosPlotter()

# Adds a DOS with a label.
plotter.add_dos("Total DOS", dos)

# Alternatively, you can add a dict of DOSs. This is the typical
# form returned by CompleteDos.get_spd/element/others_dos().
plotter.add_dos_dict({"dos1": dos1, "dos2": dos2})
plotter.add_dos_dict(complete_dos.get_spd_dos())
Parameters:
  • 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.

Parameters:
  • 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.

Parameters:
  • 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:{‘energies’:.., ‘densities’: {‘up’:...}, ‘efermi’:efermi}}
Return type:Dict of dos data. Generally of the form, {label
get_plot(xlim=None, ylim=None)[source]

Get a matplotlib plot showing the DOS.

Parameters:
  • 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.

Parameters:
  • 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.

Parameters:
  • xlim – Specifies the x-axis limits. Set to None for automatic determination.
  • ylim – Specifies the y-axis limits.
fold_point(p, lattice, coords_are_cartesian=False)[source]

Folds a point with coordinates p inside the first Brillouin zone of the lattice.

Parameters:
  • p – coordinates of one point
  • lattice – Lattice object used to convert from reciprocal to cartesian coordinates
  • coords_are_cartesian – Set to True if you are providing coordinates in cartesian coordinates. Defaults to False.
Returns:

The cartesian coordinates folded inside the first Brillouin zone

plot_brillouin_zone(bz_lattice, lines=None, labels=None, kpoints=None, fold=False, coords_are_cartesian=False, ax=None, **kwargs)[source]

Plots a 3D representation of the Brillouin zone of the structure. Can add to the plot paths, labels and kpoints

Parameters:
  • bz_lattice – Lattice object of the Brillouin zone
  • lines – list of lists of coordinates. Each list represent a different path
  • labels – dict containing the label as a key and the coordinates as value.
  • kpoints – list of coordinates
  • fold – whether the points should be folded inside the first Brillouin Zone. Defaults to False. Requires lattice if True.
  • coords_are_cartesian – Set to True if you are providing coordinates in cartesian coordinates. Defaults to False.
  • ax – matplotlib Axes or None if a new figure should be created.
  • kwargs – provided by add_fig_kwargs decorator
Returns:

matplotlib figure

keyword arguments controlling the display of the figure:

kwargs Meaning
title Title of the plot (Default: None).
show True to show the figure (default: True).
savefig ‘abc.png’ or ‘abc.eps’ to save the figure to a file.
size_kwargs Dictionary with options passed to fig.set_size_inches example: size_kwargs=dict(w=3, h=4)
tight_layout True if to call fig.tight_layout (default: False)
plot_brillouin_zone_from_kpath(kpath, ax=None, **kwargs)[source]
Gives the plot (as a matplotlib object) of the symmetry line path in
the Brillouin Zone.
Parameters:
  • kpath (HighSymmKpath) – a HighSymmKPath object
  • ax – matplotlib Axes or None if a new figure should be created.
  • **kwargs – provided by add_fig_kwargs decorator
Returns:

matplotlib figure

keyword arguments controlling the display of the figure:

kwargs Meaning
title Title of the plot (Default: None).
show True to show the figure (default: True).
savefig ‘abc.png’ or ‘abc.eps’ to save the figure to a file.
size_kwargs Dictionary with options passed to fig.set_size_inches example: size_kwargs=dict(w=3, h=4)
tight_layout True if to call fig.tight_layout (default: False)
plot_ellipsoid(hessian, center, lattice=None, rescale=1.0, ax=None, coords_are_cartesian=False, arrows=False, **kwargs)[source]

Plots a 3D ellipsoid rappresenting the Hessian matrix in input. Useful to get a graphical visualization of the effective mass of a band in a single k-point.

Parameters:
  • hessian – the Hessian matrix
  • center – the center of the ellipsoid in reciprocal coords (Default)
  • lattice – Lattice object of the Brillouin zone
  • rescale – factor for size scaling of the ellipsoid
  • ax – matplotlib Axes or None if a new figure should be created.
  • coords_are_cartesian – Set to True if you are providing a center in cartesian coordinates. Defaults to False.
  • kwargs – kwargs passed to the matplotlib function ‘plot_wireframe’. Color defaults to blue, rstride and cstride default to 4, alpha defaults to 0.2.
Returns:

matplotlib figure and matplotlib ax

Example of use:
fig,ax=plot_wigner_seitz(struct.reciprocal_lattice) plot_ellipsoid(hessian,[0.0,0.0,0.0], struct.reciprocal_lattice,ax=ax)
plot_fermi_surface(data, structure, cbm, energy_levels=[], multiple_figure=True, mlab_figure=None, kpoints_dict={}, color=(0, 0, 1), transparency_factor=[], labels_scale_factor=0.05, points_scale_factor=0.02, interative=True)[source]

Plot the Fermi surface at specific energy value.

Parameters:
  • data – energy values in a 3D grid from a CUBE file via read_cube_file function, or from a BoltztrapAnalyzer.fermi_surface_data
  • structure – structure object of the material
  • energy_levels – list of energy value of the fermi surface. Default: max energy value + 0.01 eV
  • cbm – Boolean value to specify if the considered band is a conduction band or not
  • multiple_figure – if True a figure for each energy level will be shown. If False all the surfaces will be shown in the same figure. In this las case, tune the transparency factor.
  • mlab_figure – provide a previous figure to plot a new surface on it.
  • kpoints_dict – dictionary of kpoints to show in the plot. example: {“K”:[0.5,0.0,0.5]}, where the coords are fractional.
  • color – tuple (r,g,b) of integers to define the color of the surface.
  • transparency_factor – list of values in the range [0,1] to tune the opacity of the surfaces.
  • labels_scale_factor – factor to tune the size of the kpoint labels
  • points_scale_factor – factor to tune the size of the kpoint points
  • interative – if True an interactive figure will be shown. If False a non interactive figure will be shown, but it is possible to plot other surfaces on the same figure. To make it interactive, run mlab.show().
Returns:

a Mayavi figure and a mlab module to control the plot.

Note: Experimental.
Please, double check the surface shown by using some other software and report issues.
plot_labels(labels, lattice=None, coords_are_cartesian=False, ax=None, **kwargs)[source]

Adds labels to a matplotlib Axes

Parameters:
  • labels – dict containing the label as a key and the coordinates as value.
  • lattice – Lattice object used to convert from reciprocal to cartesian coordinates
  • coords_are_cartesian – Set to True if you are providing. coordinates in cartesian coordinates. Defaults to False. Requires lattice if False.
  • ax – matplotlib Axes or None if a new figure should be created.
  • kwargs – kwargs passed to the matplotlib function ‘text’. Color defaults to blue and size to 25.
Returns:

matplotlib figure and matplotlib ax

plot_lattice_vectors(lattice, ax=None, **kwargs)[source]

Adds the basis vectors of the lattice provided to a matplotlib Axes

Parameters:
  • lattice – Lattice object
  • ax – matplotlib Axes or None if a new figure should be created.
  • kwargs – kwargs passed to the matplotlib function ‘plot’. Color defaults to green and linewidth to 3.
Returns:

matplotlib figure and matplotlib ax

plot_path(line, lattice=None, coords_are_cartesian=False, ax=None, **kwargs)[source]

Adds a line passing through the coordinates listed in ‘line’ to a matplotlib Axes

Parameters:
  • line – list of coordinates.
  • lattice – Lattice object used to convert from reciprocal to cartesian coordinates
  • coords_are_cartesian – Set to True if you are providing coordinates in cartesian coordinates. Defaults to False. Requires lattice if False.
  • ax – matplotlib Axes or None if a new figure should be created.
  • kwargs – kwargs passed to the matplotlib function ‘plot’. Color defaults to red and linewidth to 3.
Returns:

matplotlib figure and matplotlib ax

plot_points(points, lattice=None, coords_are_cartesian=False, fold=False, ax=None, **kwargs)[source]

Adds Points to a matplotlib Axes

Parameters:
  • points – list of coordinates
  • lattice – Lattice object used to convert from reciprocal to cartesian coordinates
  • coords_are_cartesian – Set to True if you are providing coordinates in cartesian coordinates. Defaults to False. Requires lattice if False.
  • fold – whether the points should be folded inside the first Brillouin Zone. Defaults to False. Requires lattice if True.
  • ax – matplotlib Axes or None if a new figure should be created.
  • kwargs – kwargs passed to the matplotlib function ‘scatter’. Color defaults to blue
Returns:

matplotlib figure and matplotlib ax

plot_wigner_seitz(lattice, ax=None, **kwargs)[source]

Adds the skeleton of the Wigner-Seitz cell of the lattice to a matplotlib Axes

Parameters:
  • lattice – Lattice object
  • ax – matplotlib Axes or None if a new figure should be created.
  • kwargs – kwargs passed to the matplotlib function ‘plot’. Color defaults to black and linewidth to 1.
Returns:

matplotlib figure and matplotlib ax