Skip to content

Commit

Permalink
Merge pull request #349 from paulneves77/patch-2
Browse files Browse the repository at this point in the history
Added optional arguments to plot fermi surface outside first Brillouin zone
  • Loading branch information
utf authored Oct 26, 2023
2 parents 8f38de7 + 0084bae commit b6a1a1a
Showing 1 changed file with 14 additions and 2 deletions.
16 changes: 14 additions & 2 deletions src/ifermi/surface.py
Original file line number Diff line number Diff line change
Expand Up @@ -421,6 +421,8 @@ def from_band_structure(
property_data: dict[Spin, np.ndarray] | None = None,
property_kpoints: np.ndarray | None = None,
calculate_dimensionality: bool = False,
supercell_dim: tuple[int, int, int] = (3, 3, 3),
trim_to_first_bz: bool = True,
) -> FermiSurface:
"""Create a FermiSurface from a pymatgen band structure object.
Expand Down Expand Up @@ -452,6 +454,9 @@ def from_band_structure(
property_kpoints: The k-points on which the data is generated.
Must be used in combination with ``data``.
calculate_dimensionality: Whether to calculate isosurface dimensionalities.
supercell_dim: The supercell mesh dimensions.
trim_to_first_bz: If true, only includes Fermi surface within one Brillouin
zone. If false, include Fermi surface within entire supercell.
Returns:
A Fermi surface.
Expand Down Expand Up @@ -486,7 +491,7 @@ def from_band_structure(
raise ValueError("Both data and kpoints must be specified.")

bands, kpoints = expand_bands(
bands, kpoints, supercell_dim=(3, 3, 3), center=(1, 1, 1)
bands, kpoints, supercell_dim=supercell_dim, center=(1, 1, 1)
)
if isinstance(decimate_factor, int):
# increase number of target faces to account for 3x3x3 supercell
Expand All @@ -502,6 +507,7 @@ def from_band_structure(
smooth=smooth,
calculate_dimensionality=calculate_dimensionality,
property_interpolator=interpolator,
trim_to_first_bz=trim_to_first_bz,
)

if sum(len(spin_iso) for spin_iso in isosurfaces.values()) == 0:
Expand Down Expand Up @@ -552,6 +558,7 @@ def compute_isosurfaces(
smooth: bool = False,
calculate_dimensionality: bool = False,
property_interpolator: LinearInterpolator | None = None,
trim_to_first_bz: bool = True,
) -> dict[Spin, list[Isosurface]]:
"""Compute the isosurfaces at a particular energy level.
Expand All @@ -576,6 +583,8 @@ def compute_isosurfaces(
calculate_dimensionality: Whether to calculate isosurface dimensionality.
property_interpolator: An interpolator class for interpolating properties
onto the surface. If ``None``, no properties will be calculated.
trim_to_first_bz: If true, only includes Fermi surface within one Brillouin
zone. If false, include Fermi surface within entire supercell.
Returns:
A dictionary containing a list of isosurfaces for each spin channel.
Expand Down Expand Up @@ -614,6 +623,7 @@ def compute_isosurfaces(
smooth,
calculate_dimensionality,
property_interpolator,
trim_to_first_bz,
)
spin_isosurface.extend(band_isosurfaces)

Expand All @@ -635,6 +645,7 @@ def _calculate_band_isosurfaces(
smooth: bool,
calculate_dimensionality: bool,
property_interpolator: LinearInterpolator | None,
trim_to_first_bz: bool = True,
):
"""Helper function to calculate the connected isosurfaces for a band."""
from skimage.measure import marching_cubes
Expand Down Expand Up @@ -698,7 +709,8 @@ def _calculate_band_isosurfaces(
for (subverts, subfaces), idx in zip(subsurfaces, mapping):
# convert vertices to cartesian coordinates
subverts = np.dot(subverts, rlat)
subverts, subfaces = trim_surface(reciprocal_space, subverts, subfaces)
if trim_to_first_bz:
subverts, subfaces = trim_surface(reciprocal_space, subverts, subfaces)

if len(subverts) == 0:
# skip surfaces that do not enter the reciprocal space boundaries
Expand Down

0 comments on commit b6a1a1a

Please sign in to comment.