ifermi package¶
ifermi.analysis module¶
Isosurface and isoline analysis functions.
-
ifermi.analysis.
average_properties
(vertices, faces, properties, norm=False)[source]¶ Average property across an isosurface.
- Parameters
vertices (
ndarray
) – A (n, 3) float array of the vertices in the isosurface.faces (
ndarray
) – A (m, 3) int array of the faces of the isosurface.properties (
ndarray
) – A (m, …) array of the face properties as scalars or vectors.norm (
bool
) – Whether to average the norm of the properties (vector properties only).
- Return type
- Returns
The averaged property.
-
ifermi.analysis.
connected_images
(fractional_vertices, tol=1e-05)[source]¶ Find the images a set of vertices is connected to.
Note: This function expects the vertices to only belong to a single connected mesh, and the vertices should cover a 3x3x3 supercell, where the coordinates from -0.5 to 0.5 indicate the center image.
- Parameters
fractional_vertices (
ndarray
) – A (n, 3) float array of the vertices in the isosurface in fractional coordinates.tol – A tolerance for evaluating whether two vertices are equivalent.
- Return type
- Returns
A (n, 3) int array of the images across which the mesh is connected periodically.
-
ifermi.analysis.
connected_subsurfaces
(vertices, faces)[source]¶ Find connected sub-surfaces (those that share edges).
-
ifermi.analysis.
equivalent_surfaces
(surfaces_vertices, tol=1e-05)[source]¶ Note: This function expects the vertices of each surface to only belong to a single connected mesh, and the vertices should cover a 3x3x3 supercell, where the coordinates from -0.5 to 0.5 indicate the center image.
- Parameters
- Return type
- Returns
A (m, ) int array that maps each each surface in the original surface array to its equivalent.
-
ifermi.analysis.
equivalent_vertices
(vertices, tol=1e-05)[source]¶ Find vertices that are equivalent (closer than a tolerance).
Note that the algorithm used is effectively recursive. If vertex a is within the tolerance of b, and b is within the tolerance of c, even if a and c and not within the tolerance, a, b, and c will be considered equivalent.
-
ifermi.analysis.
isosurface_dimensionality
(fractional_vertices, faces)[source]¶ Calculate isosurface properties a single isosurface (fully connected).
The vertices must cover a 3x3x3 supercell and must not have been trimmed to fit inside the reciprocal lattice.
Note: This function expects the vertices to only belong to a single connected mesh, and the vertices should cover a 3x3x3 supercell, where the coordinates from -0.5 to 0.5 indicate the center image.
- Parameters
- Return type
- Returns
The dimensionality and (n, 3) int array orientation of the isosurface.
-
ifermi.analysis.
line_orientation
(images)[source]¶ Get the orientation (direction vector) from a list of rank 1 connected images.
-
ifermi.analysis.
longest_simple_paths
(vertices, edges)[source]¶ Find the shortest paths that go through all vertices.
The lines are broken up into the connected sublines. Note this function is only designed to work with connected sublines that are either simple cycles or chains with no branches.
-
ifermi.analysis.
plane_orientation
(images)[source]¶ Get the orientation (surface normal) from a list of rank 2 connected images.
-
ifermi.analysis.
sample_line_uniform
(segments, spacing)[source]¶ Sample line segments to a consistent density.
Note: the segments must be ordered so that they are adjacent.
-
ifermi.analysis.
sample_surface_uniform
(vertices, faces, grid_size)[source]¶ Sample isosurface faces uniformly.
The algorithm works by:
Splitting the mesh into a uniform grid with block sizes determined by
grid_size
.For each cell in the grid, finds whether the center of any faces falls within the cell.
If multiple face centers fall within the cell, it picks the closest one to the center of the cell. If no face centers fall within the cell then the cell is ignored.
Returns the indices of all the faces that have been selected.
This algorithm is not well optimised for small grid sizes.
ifermi.brillouin_zone module¶
Brillouin zone and slice geometries.
-
class
ifermi.brillouin_zone.
ReciprocalCell
(reciprocal_lattice, vertices, faces, centers, normals)[source]¶ Bases:
monty.json.MSONable
A parallelepiped reciprocal lattice cell.
-
reciprocal_lattice
¶ A (3, 3) float array of the reciprocal lattice vectors.
-
vertices
¶ A (n, 3) float array of the vertices of the Brillouin zone edges.
-
faces
¶ A (m, 3) int array of the faces of the Brillouin zone.
-
centers
¶ A (m, 3) float array of the face centers.
-
normals
¶ A (m, 3) float array of the face normal.
-
classmethod
from_structure
(structure)[source]¶ Initialise the reciprocal cell from a structure.
- Parameters
structure (
Structure
) – A structure.- Return type
- Returns
The reciprocal cell.
-
-
class
ifermi.brillouin_zone.
ReciprocalSlice
(reciprocal_space, vertices, transformation)[source]¶ Bases:
monty.json.MSONable
A slice along a pane in reciprocal space.
-
reciprocal_space
¶ The reciprocal space that the slice belongs to.
-
vertices
¶ A (n, 2) float array of the vertices defining the intersection of the slice with the Brillouin zone.
-
transformation
¶ (n, 4, 4) float array of the transformation that maps points in the 3D Brillouin zone to points on the 2D slice.
-
-
class
ifermi.brillouin_zone.
WignerSeitzCell
(reciprocal_lattice, vertices, faces, centers, normals)[source]¶ Bases:
ifermi.brillouin_zone.ReciprocalCell
WignerSeitz cell of the reciprocal lattice.
-
reciprocal_lattice
¶ A (3, 3) float array of the reciprocal lattice vectors.
-
vertices
¶ A (n, 3) float array of the vertices of the Brillouin zone edges.
-
faces
¶ A (m, 3) int array of the faces of the Brillouin zone.
-
centers
¶ A (m, 3) float array of the face centers.
-
normals
¶ A (m, 3) float array of the face normal.
-
ifermi.interpolate module¶
Tools for Fourier and linear interpolation.
-
class
ifermi.interpolate.
FourierInterpolator
(band_structure, magmom=None, mommat=None)[source]¶ Bases:
object
Class to perform Fourier interpolation of electronic band structures.
Interpolation is performed using BoltzTraP2.
- Parameters
band_structure (
BandStructure
) – The Bandstructure object to be interpolated.mommat (
Optional
[ndarray
]) – Momentum matrix, as supported by BoltzTraP2.
-
interpolate_bands
(interpolation_factor=5, return_velocities=False, nworkers=- 1)[source]¶ Get an interpolated pymatgen band structure.
Note, the interpolation mesh is determined using by
interpolate_factor
option in theFourierInterpolator
constructor.The degree of parallelization is controlled by the
nworkers
option.- Parameters
interpolation_factor (
float
) – The factor by which the band structure will be interpolated.return_velocities (
bool
) – Whether to return the group velocities.nworkers (
int
) – The number of processors used to perform the interpolation. If set to-1
, the number of workers will be set to the number of CPU cores.
- Returns
The interpolated electronic structure. If
return_velocities
is True, the group velocities will also be returned as a dict of{Spin: velocities}
where velocities is a numpy array with the shape (nbands, nkpoints, 3) and has units of m/s.
-
class
ifermi.interpolate.
LinearInterpolator
(kpoints, data)[source]¶ Bases:
object
Class to perform linear interpolation of periodic properties.
- Parameters
kpoints (
ndarray
) – The k-points in fractional coordinates as a numpy array. with the shape (nkpoints, 3). Note, the k-points must cover the full Brillouin zone, not just the irreducible part.data (
Dict
[Spin
,ndarray
]) – The data to interpolate. Should be given for spin up and spin down bands. If the system is not spin polarized then only spin up should be set. The data for each spin channel should be a numpy array with the shape (nbands, nkpoints, …). The values to interpolate can be scalar or multidimensional.
-
ifermi.interpolate.
trim_bandstructure
(energy_cutoff, band_structure)[source]¶ Trim the number of bands in a band structure object based on a cutoff.
- Parameters
energy_cutoff (
float
) – An energy cutoff within which to keep the bands. If the system is metallic then the bands to keep will fall within +/- the cutoff around the Fermi level. If the system has a band gap, the bands from the VBM - energy_cutoff to CBM + energy_cutoff will be kept.band_structure (
BandStructure
) – A band structure.
- Return type
- Returns
A trimmed band structure.
ifermi.kpoints module¶
k-point manipulation functions.
-
ifermi.kpoints.
kpoints_from_bandstructure
(bandstructure, cartesian=False)[source]¶ Extract the k-points from a band structure.
- Parameters
bandstructure (
BandStructure
) – A band structure object.cartesian (
bool
) – Whether to return the k-points in cartesian coordinates.
- Return type
- Returns
A (n, 3) float array of the k-points.
-
ifermi.kpoints.
kpoints_to_first_bz
(kpoints, tol=1e-05)[source]¶ Translate fractional k-points to the first Brillouin zone.
- I.e. all k-points will have fractional coordinates:
-0.5 <= fractional coordinates < 0.5
ifermi.plot module¶
Tools to plot FermiSurface and FermiSlice objects.
-
class
ifermi.plot.
FermiSlicePlotter
(fermi_slice, symprec=0.001)[source]¶ Bases:
object
Class to plot 2D isolines through a FermiSurface.
- Parameters
fermi_slice (
FermiSlice
) – A slice through a Fermi surface.symprec (
float
) – The symmetry precision in Angstrom for determining the high symmetry k-point labels.
-
get_plot
(ax=None, spin=None, colors=None, color_properties=True, vector_properties=False, projection_axis=None, scale_linewidth=False, vector_spacing=0.2, cmin=None, cmax=None, vnorm=None, hide_slice=False, hide_labels=False, hide_cell=False, arrow_pivot='tail', slice_kwargs=None, cbar_kwargs=None, quiver_kwargs=None, bz_kwargs=None, sym_pt_kwargs=None, sym_label_kwargs=None)[source]¶ Plot the Fermi slice.
- Parameters
ax (
Optional
[Any
]) – Matplotlib axes object on which to plot.spin (
Optional
[Spin
]) – Which spin channel to plot. By default plot both spin channels if available.colors (
Union
[str
,dict
,list
,None
]) –The color specification for the iso-surfaces. Valid options are:
A single color to use for all Fermi isolines, specified as a tuple of rgb values from 0 to 1. E.g., red would be
(1, 0, 0)
.A list of colors, specified as above.
A dictionary of
{Spin.up: color1, Spin.down: color2}
, where the colors are specified as above.A string specifying which matplotlib colormap to use. See https://matplotlib.org/tutorials/colors/colormaps.html for more information.
None
, in which case the default colors will be used.
color_properties (
Union
[str
,bool
]) – Whether to use the properties to color the Fermi isolines. If the properties is a vector then the norm of the properties will be used. Note, this will only take effect if the Fermi slice has properties. If set to True, the viridis colormap will be used. Alternative colormaps can be selected by settingcolor_properties
to a matplotlib colormap name. This setting will override thecolors
option. For vector properties, the arrows are colored according to the norm of the properties by default. If used in combination with theprojection_axis
option, the color will be determined by the dot product of the properties with the projection axis.vector_properties (
Union
[str
,bool
]) – Whether to plot arrows for vector properties. Note, this will only take effect if the Fermi slice has vector properties. If set to True, the viridis colormap will be used. Alternative colormaps can be selected by settingvector_properties
to a matplotlib colormap name. By default, the arrows are colored according to the norm of the properties. If used in combination with theprojection_axis
option, the color will be determined by the dot product of the properties with the projection axis.projection_axis (
Optional
[Tuple
[int
,int
,int
]]) – Projection axis that can be used to calculate the color of vector properties. If None, the norm of the properties will be used, otherwise the color will be determined by the dot product of the properties with the projection axis. Only has an effect when used with thevector_properties
option.scale_linewidth (
Union
[bool
,float
]) – Scale the linewidth by the absolute value of the segment properties. Can be true, false or a number. If a number, then this will be used as the max linewidth for scaling.vector_spacing (
float
) – The rough spacing between arrows. Uses a custom algorithm for resampling the Fermi surface to ensure that arrows are not too close together. Only has an effect when used with thevector_properties
option.cmin (
Optional
[float
]) – Minimum intensity for normalising properties colors (including vector colors). Only has an effect when used withcolor_properties
orvector_properties
options.cmax (
Optional
[float
]) – Maximum intensity for normalising properties colors (including vector colors). Only has an effect when used withcolor_properties
orvector_properties
options.vnorm (
Optional
[float
]) – The value by which to normalize the vector lengths. For example, spin properties should typically have a norm of 1 whereas group velocity properties can have larger or smaller norms depending on the structure. By changing this number, the size of the vectors will be scaled. Note that the properties of two materials can only be compared quantitatively if a fixed values is used for both plots. Only has an effect when used with thevector_properties
option.hide_slice (
bool
) – Whether to hide the Fermi surface. Only recommended in combination with thevector_properties
option.hide_labels (
bool
) – Whether to show the high-symmetry k-point labels.hide_cell (
bool
) – Whether to show the reciprocal cell boundary.arrow_pivot (
str
) – The part of the arrow that is anchored to the X, Y grid. The arrow rotates about this point, options are: tail, middle, tip.slice_kwargs (
Optional
[Dict
[str
,Any
]]) – Optional arguments that are passed toLineCollection
and are used to style the iso slice.cbar_kwargs (
Optional
[Dict
[str
,Any
]]) – Optional arguments that are passed tofig.colorbar
.quiver_kwargs (
Optional
[Dict
[str
,Any
]]) – Optional arguments that are passed toax.quiver
and are used to style the arrows.bz_kwargs (
Optional
[Dict
[str
,Any
]]) – Optional arguments that passed toLineCollection
and used to style the Brillouin zone boundary.sym_pt_kwargs (
Optional
[Dict
[str
,Any
]]) – Optional arguments that are passed toax.scatter
and are used to style the high-symmetry k-point symbols.sym_label_kwargs (
Optional
[Dict
[str
,Any
]]) – Optional arguments that are passed toax.text
and are used to style the high-symmetry k-point labels.
- Returns
matplotlib pyplot object.
-
class
ifermi.plot.
FermiSurfacePlotter
(fermi_surface, symprec=0.001)[source]¶ Bases:
object
Class to plot a FermiSurface.
- Parameters
fermi_surface (
FermiSurface
) – A FermiSurface object.symprec (
float
) – The symmetry precision in Angstrom for determining the high symmetry k-point labels.
-
get_plot
(plot_type='plotly', spin=None, colors=None, azimuth=45.0, elevation=35.0, color_properties=True, vector_properties=False, projection_axis=None, vector_spacing=0.2, cmin=None, cmax=None, vnorm=None, hide_surface=False, hide_labels=False, hide_cell=False, **plot_kwargs)[source]¶ Plot the Fermi surface.
- Parameters
plot_type (
str
) – Method used for plotting. Valid options are: “matplotlib”, “plotly”, “mayavi”, “crystal_toolkit”.spin (
Optional
[Spin
]) – Which spin channel to plot. By default plot both spin channels if available.azimuth (
float
) – The azimuth of the viewpoint in degrees. i.e. the angle subtended by the position vector on a sphere projected on to the x-y plane.elevation (
float
) – The zenith angle of the viewpoint in degrees, i.e. the angle subtended by the position vector and the z-axis.colors (
Union
[str
,dict
,list
,None
]) –The color specification for the iso-surfaces. Valid options are:
A single color to use for all Fermi surfaces, specified as a tuple of rgb values from 0 to 1. E.g., red would be
(1, 0, 0)
.A list of colors, specified as above.
A dictionary of
{Spin.up: color1, Spin.down: color2}
, where the colors are specified as above.A string specifying which matplotlib colormap to use. See https://matplotlib.org/tutorials/colors/colormaps.html for more information.
None
, in which case the default colors will be used.
color_properties (
Union
[str
,bool
]) – Whether to use the properties to color the Fermi surface. If the properties is a vector then the norm of the properties will be used. Note, this will only take effect if the Fermi surface has properties. If set to True, the viridis colormap will be used. Alternative colormaps can be selected by settingcolor_properties
to a matplotlib colormap name. This setting will override thecolors
option. For vector properties, the arrows are colored according to the norm of the properties by default. If used in combination with theprojection_axis
option, the color will be determined by the dot product of the properties with the projection axis.vector_properties (
Union
[str
,bool
]) – Whether to plot arrows for vector properties. Note, this will only take effect if the Fermi surface has vector properties. If set to True, the viridis colormap will be used. Alternative colormaps can be selected by settingvector_properties
to a matplotlib colormap name. By default, the arrows are colored according to the norm of the properties. If used in combination with theprojection_axis
option, the color will be determined by the dot product of the properties with the projection axis.projection_axis (
Optional
[Tuple
[int
,int
,int
]]) – Projection axis that can be used to calculate the color of vector properties. If None, the norm of the properties will be used, otherwise the color will be determined by the dot product of the properties with the projection axis. Only has an effect when used with thevector_properties
option.vector_spacing (
float
) – The rough spacing between arrows. Uses a custom algorithm for resampling the Fermi surface to ensure that arrows are not too close together. Only has an effect when used with thevector_properties
option.cmin (
Optional
[float
]) – Minimum intensity for normalising properties colors (including vector colors). Only has an effect when used withcolor_properties
orvector_properties
options.cmax (
Optional
[float
]) – Maximum intensity for normalising properties colors (including vector colors). Only has an effect when used withcolor_properties
orvector_properties
options.vnorm (
Optional
[float
]) – The value by which to normalize the vector lengths. For example, spin properties should typically have a norm of 1 whereas group velocity properties can have larger or smaller norms depending on the structure. By changing this number, the size of the vectors will be scaled. Note that the properties of two materials can only be compared quantitatively if a fixed values is used for both plots. Only has an effect when used with thevector_properties
option.hide_surface (
bool
) – Whether to hide the Fermi surface. Only recommended in combination with thevector_properties
option.hide_labels (
bool
) – Whether to show the high-symmetry k-point labels.hide_cell (
bool
) – Whether to show the reciprocal cell boundary.**plot_kwargs – Other keyword arguments supported by the individual plotting methods.
-
ifermi.plot.
cmap_to_plotly
(colormap)[source]¶ Convert a matplotlib colormap to plotly colorscale format.
-
ifermi.plot.
get_face_arrows
(fermi_surface, spins, vector_spacing, vnorm, projection_axis)[source]¶ Get face arrows from vector properties.
- Parameters
fermi_surface (
FermiSurface
) – The fermi surface containing the isosurfaces and properties.spins (
List
[Spin
]) – Spin channels from which to extract arrows.vector_spacing (
float
) – The rough spacing between arrows. Uses a custom algorithm for resampling the Fermi surface to ensure that arrows are not too close together.vnorm (
Optional
[float
]) – The value by which to normalize the vector lengths. For example, spin properties should typically have a norm of 1 whereas group velocity properties can have larger or smaller norms depending on the structure. By changing this number, the size of the vectors will be scaled. Note that the properties of two materials can only be compared quantitatively if a fixed values is used for both plots.projection_axis (
Optional
[Tuple
[int
,int
,int
]]) – Projection axis that can be used to calculate the color of vector projections. If None, the norm of the properties will be used, otherwise the color will be determined by the dot product of the properties with the projection axis.
- Return type
- Returns
The arrows, as a list of (starts, stops, intensities) for each face. The starts and stops are numpy arrays with the shape (narrows, 3) and intensities is a numpy array with the shape (narrows, ). The intensities are used to color the arrows during plotting.
-
ifermi.plot.
get_isosurface_colors
(colors, fermi_object, spins)[source]¶ Get colors for each Fermi surface.
- Parameters
colors (
Union
[str
,dict
,list
,None
]) –The color specification. Valid options are:
A single color to use for all Fermi surfaces, specified as a tuple of rgb values from 0 to 1. E.g., red would be
(1, 0, 0)
.A list of colors, specified as above.
A dictionary of
{Spin.up: color1, Spin.down: color2}
, where the colors are specified as above.A string specifying which matplotlib colormap to use. See https://matplotlib.org/tutorials/colors/colormaps.html for more information.
None
, in which case the default colors will be used.
fermi_object (
Union
[FermiSurface
,FermiSlice
]) – A Fermi surface or Fermi slice object.spins (
List
[Spin
]) – A list of spins for which colors will be generated.
- Return type
- Returns
The colors as a list of tuples, where each color is specified as the rgb values from 0 to 1. E.g., red would be
(1, 0, 0)
.
-
ifermi.plot.
get_segment_arrows
(fermi_slice, spins, vector_spacing, vnorm, projection_axis)[source]¶ Get segment arrows from vector properties.
- Parameters
fermi_slice (
FermiSlice
) – The Fermi slice containing the isolines and properties.spins (
Collection
[Spin
]) – Spin channels from which to extract arrows.vector_spacing (
float
) – The rough spacing between arrows. Uses a custom algorithm for resampling the Fermi slic to ensure that arrows are not too close together.vnorm (
Optional
[float
]) – The value by which to normalize the vector lengths. For example, spin properties should typically have a norm of 1 whereas group velocity properties can have larger or smaller norms depending on the structure. By changing this number, the size of the vectors will be scaled. Note that the properties of two materials can only be compared quantitatively if a fixed values is used for both plots.projection_axis (
Optional
[Tuple
[int
,int
,int
]]) – Projection axis that can be used to calculate the color of vector projects. If None, the norm of the properties will be used, otherwise the color will be determined by the dot product of the properties with the properties axis.
- Return type
- Returns
The arrows, as a list of (starts, stops, intensities) for each face. The starts and stops are numpy arrays with the shape (narrows, 3) and intensities is a numpy array with the shape (narrows, ). The intensities are used to color the arrows during plotting.
-
ifermi.plot.
plotly_arrow
(start, stop, color, line_kwargs=None, cone_kwargs=None)[source]¶ Create an arrow object.
- Parameters
start (
ndarray
) – The starting coordinates.stop (
ndarray
) – The ending coordinates.color (
Tuple
[float
,float
,float
]) – The arrow color in rgb format as a tuple of floats from 0 to 1.line_kwargs (
Optional
[Dict
[str
,Any
]]) – Additional keyword arguments used to style the arrow shaft and that are passed toScatter3d
.cone_kwargs (
Optional
[Dict
[str
,Any
]]) – Additional keyword arguments used to style the arrow cone and that are passed toCone
.
- Return type
- Returns
The arrow, formed by a line and cone.
-
ifermi.plot.
save_plot
(plot, filename, scale=4)[source]¶ Save a plot to file.
ifermi.slice module¶
Tools to generate Isolines and Fermi slices.
-
class
ifermi.slice.
FermiSlice
(isolines, reciprocal_slice, structure)[source]¶ Bases:
monty.json.MSONable
A FermiSlice object is a 2D slice through a Fermi surface.
-
isolines
¶ A dict containing a list of isolines for each spin channel.
-
reciprocal_space
¶ A reciprocal slice defining the intersection of the slice with the Brillouin zone edges.
-
structure
¶ The structure.
-
all_properties
(spins=None, projection_axis=None, norm=False)[source]¶ Get the properties for all isolines.
- Parameters
spins (
Union
[Spin
,Collection
[Spin
],None
]) – One or more spin channels to select. Default is all spins available.projection_axis (
Optional
[Tuple
[int
,int
,int
]]) – A (3, ) in array of the axis to project the properties onto (vector properties only).norm (
bool
) – Calculate the norm of the properties (vector properties only). Ignored ifprojection_axis
is set.
- Return type
- Returns
A list of properties arrays for each isosurface.
-
classmethod
from_fermi_surface
(fermi_surface, plane_normal, distance=0)[source]¶ Get a slice through the Fermi surface.
The slice is defined by the intersection of a plane with the Fermi surface.
- Parameters
- Return type
- Returns
The Fermi slice.
-
property
n_lines_per_band
¶ Get number of lines for each band index for each spin channel.
Returned as a dict of
{spin: {band_idx: count}}
.
-
property
n_lines_per_spin
¶ Get number of lines per spin channel.
Returned as a dict of
{spin: count}
.
-
-
class
ifermi.slice.
Isoline
(segments, band_idx, properties=None)[source]¶ Bases:
monty.json.MSONable
An isoline object contains line segments mesh and line properties.
-
segments
¶ A (n, 2, 2) float array of the line segments..
-
band_idx
¶ The band index to which the slice belongs.
-
properties
¶ An optional (n, …) float array containing segment properties as scalars or vectors.
-
-
ifermi.slice.
interpolate_segments
(segments, properties, max_spacing)[source]¶ Resample a series of line segments to a consistent density.
Note: the segments must be ordered so that they are adjacent.
- Parameters
- Returns
The interpolated segments and properties.
- Return type
(segments, properties)
-
ifermi.slice.
process_lines
(segments, face_idxs)[source]¶ Process segments and face_idxs from mesh_multiplane.
The key issue is that the segments from mesh_multiplane do not correspond to individual lines, nor are they sorted in a continuous order. Instead they are just a list of randomly ordered segments. This causes trouble later on when trying to interpolate the lines or add equally spaced arrows.
The goal of this function is to identify the separate paths in the segments (a path is a collection of segments that are connected together), and return these segments in the order that they are connected. By looping through each segment, you will be proceeding along the line in a certain direction.
Because the original segments may contain multiple paths, a list of segments and their corresponding face indices are returned.
Lastly, note that there is no guarantee that all of the original segments will be used in the processed segments. This is because some of the original segments are very small and will be filtered out.
ifermi.surface module¶
Tools to generate isosurfaces and Fermi surfaces.
-
class
ifermi.surface.
FermiSurface
(isosurfaces, reciprocal_space, structure)[source]¶ Bases:
monty.json.MSONable
A FermiSurface object contains isosurfaces and the reciprocal lattice definition.
-
isosurfaces
¶ A dict containing a list of isosurfaces for each spin channel.
-
reciprocal_space
¶ A reciprocal space defining periodic boundary conditions.
-
structure
¶ The structure.
-
all_properties
(spins=None, projection_axis=None, norm=False)[source]¶ Get the properties for all isosurfaces.
- Parameters
spins (
Union
[Spin
,Collection
[Spin
],None
]) – One or more spin channels to select. Default is all spins available.projection_axis (
Optional
[Tuple
[int
,int
,int
]]) – A (3, ) in array of the axis to project the properties onto (vector properties only).norm (
bool
) – Calculate the norm of the properties (vector properties only). Ignored ifprojection_axis
is set.
- Return type
- Returns
A list of properties arrays for each isosurface.
-
property
area_surfaces
¶ Area of each isosurface in the Fermi surface in Å-2.
-
average_properties
(norm=False, projection_axis=None)[source]¶ Average property across the full Fermi surface.
-
average_properties_surfaces
(norm=False, projection_axis=None)[source]¶ Average property for each isosurface in the Fermi surface.
- Parameters
- Return type
- Returns
The averaged property for each surface in each spin channel.
-
classmethod
from_band_structure
(band_structure, mu=0.0, wigner_seitz=False, decimate_factor=None, decimate_method='quadric', smooth=False, property_data=None, property_kpoints=None, calculate_dimensionality=False)[source]¶ Create a FermiSurface from a pymatgen band structure object.
- Parameters
band_structure (
BandStructure
) – A band structure. The k-points must cover the full Brillouin zone (i.e., not just be the irreducible mesh). Use theifermi.interpolator.FourierInterpolator
class to expand the k-points to the full Brillouin zone if required.mu (
float
) – Energy offset from the Fermi energy at which the isosurface is calculated.wigner_seitz (
bool
) – Controls whether the cell is the Wigner-Seitz cell or the reciprocal unit cell parallelepiped.decimate_factor (
Optional
[float
]) – If method is “quadric”, and factor is a floating point value then factor is the scaling factor by which to reduce the number of faces. I.e., final # faces = initial # faces * factor. If method is “quadric” but factor is an integer then factor is the target number of final faces. If method is “cluster”, factor is the voxel size in which to cluster points. Default is None (no decimation).decimate_method (
str
) – Algorithm to use for decimation. Options are “quadric” or “cluster”.smooth (
bool
) – If True, will smooth resulting isosurface. Requires PyMCubes. See compute_isosurfaces for more information.property_data (
Optional
[Dict
[Spin
,ndarray
]]) – A property to project onto the Fermi surface. It should be given as a dict of{spin: properties}
, where properties is numpy array with shape (nbands, nkpoints, …). The number of bands should equal the number of bands in the band structure but the k-point mesh can be different. Must be used in combination withkpoints
.property_kpoints (
Optional
[ndarray
]) – The k-points on which the data is generated. Must be used in combination withdata
.calculate_dimensionality (
bool
) – Whether to calculate isosurface dimensionalities.
- Return type
- Returns
A Fermi surface.
-
get_fermi_slice
(plane_normal, distance=0)[source]¶ Get a slice through the Fermi surface.
Slice defined by the intersection of a plane with the Fermi surface.
-
property
n_surfaces_per_band
¶ Get number of surfaces for each band index for each spin channel.
Returned as a dict of
{spin: {band_idx: count}}
.
-
property
n_surfaces_per_spin
¶ Get number of surfaces per spin channel.
Returned as a dict of
{spin: count}
.
-
property
spins
¶ The spin channels in the Fermi surface.
-
-
class
ifermi.surface.
Isosurface
(vertices, faces, band_idx, properties=None, dimensionality=None, orientation=None)[source]¶ Bases:
monty.json.MSONable
An isosurface object contains a triangular mesh and surface properties.
-
vertices
¶ A (n, 3) float array of the vertices in the isosurface.
-
faces
¶ A (m, 3) int array of the faces of the isosurface.
-
band_idx
¶ The band index to which the surface belongs.
-
properties
¶ An optional (m, …) float array containing face properties as scalars or vectors.
-
dimensionality
¶ The dimensionality of the surface.
-
orientation
¶ The orientation of the surface (for 1D and 2D surfaces only).
-
-
ifermi.surface.
compute_isosurfaces
(bands, kpoints, fermi_level, reciprocal_space, decimate_factor=None, decimate_method='quadric', smooth=False, calculate_dimensionality=False, property_interpolator=None)[source]¶ Compute the isosurfaces at a particular energy level.
- Parameters
bands (
Dict
[Spin
,ndarray
]) – The band energies, given as a dictionary of{spin: energies}
, where energies has the shape (nbands, nkpoints).kpoints (
ndarray
) – The k-points in fractional coordinates.fermi_level (
float
) – The energy at which to calculate the Fermi surface.reciprocal_space (
ReciprocalCell
) – The reciprocal space representation.decimate_factor (
Optional
[float
]) – If method is “quadric”, and factor is a floating point value then factor is the scaling factor by which to reduce the number of faces. I.e., final # faces = initial # faces * factor. If method is “quadric” but factor is an integer then factor is the target number of final faces. If method is “cluster”, factor is the voxel size in which to cluster points. Default is None (no decimation).decimate_method (
str
) – Algorithm to use for decimation. Options are “quadric” or “cluster”.smooth (
bool
) – If True, will smooth resulting isosurface. Requires PyMCubes. Smoothing algorithm will use constrained smoothing algorithm to preserve fine details if input dimension is lower than (500, 500, 500), otherwise will apply a Gaussian filter.calculate_dimensionality (
bool
) – Whether to calculate isosurface dimensionality.property_interpolator (
Optional
[LinearInterpolator
]) – An interpolator class for interpolating properties onto the surface. IfNone
, no properties will be calculated.
- Return type
Dict
[Spin
,List
[Isosurface
]]- Returns
A dictionary containing a list of isosurfaces for each spin channel.
-
ifermi.surface.
decimate_mesh
(vertices, faces, factor, method='quadric')[source]¶ Decimate mesh to reduce the number of triangles and vertices.
The open3d package is required for decimation.
- Parameters
vertices (
ndarray
) – A (n, 3) float array of the vertices in the isosurface.faces (
ndarray
) – A (m, 3) int array of the faces of the isosurface.factor (
Union
[int
,float
]) – If method is “quadric”, and factor is a floating point value then factor is the scaling factor by which to reduce the number of faces. I.e., final # faces = initial # faces * factor. If method is “quadric” but factor is an integer then factor is the target number of final faces. If method is “cluster”, factor is the voxel size in which to cluster points.method – Algorithm to use for decimation. Options are “quadric” or “cluster”.
- Returns
(vertices, faces) of the decimated mesh.
-
ifermi.surface.
expand_bands
(bands, fractional_kpoints, supercell_dim=3, 3, 3, center=0, 0, 0)[source]¶ Expand the band energies and k-points with periodic boundary conditions.
- Parameters
bands (
Dict
[Spin
,ndarray
]) – The band energies, given as a dictionary of{spin: energies}
, where energies has the shape (nbands, nkpoints).fractional_kpoints (
ndarray
) – A (n, 3) float array of the fractional k-point coordinates.supercell_dim (
Tuple
[int
,int
,int
]) – The supercell mesh dimensions.center (
Tuple
[int
,int
,int
]) – The cell on which the supercell is centered.
- Return type
- Returns
(energies, kpoints) The expanded band energies and k-points.
-
ifermi.surface.
face_properties
(interpolator, vertices, faces, band_idx, spin, reciprocal_lattice)[source]¶ Interpolate properties data onto the Fermi surfaces.
- Parameters
interpolator (
LinearInterpolator
) – Periodic interpolator for property data.vertices (
ndarray
) – A (n, 3) float array of the vertices in the isosurface.faces (
ndarray
) – A (m, 3) int array of the faces of the isosurface.band_idx (
int
) – The band that the isosurface belongs to.spin (
Spin
) – The spin channel the isosurface belongs to.reciprocal_lattice (
ndarray
) – Reciprocal lattice matrix
- Return type
- Returns
A (m, …) float array of the interpolated properties at the center of each face in the isosurface.
-
ifermi.surface.
trim_surface
(reciprocal_cell, vertices, faces)[source]¶ Trim the surface to remove parts outside the cell boundaries.
Will add new triangles at the boundary edges as necessary to produce a smooth surface.
- Parameters
reciprocal_cell (
ReciprocalCell
) – The reciprocal space object.vertices (
ndarray
) – A (n, 3) float array of the vertices in the isosurface.faces (
ndarray
) – A (m, 3) int array of the faces of the isosurface.
- Return type
- Returns
(vertices, faces) of the trimmed surface.