easyunfold.unfold#

The main module for unfolding workflow and algorithm

Module Contents#

Classes#

UnfoldKSet

High-level interface for unfolding and serialization.

Unfold

Low level interface for performing unfolding related calculations. obtain the effective band structure (EBS).

Functions#

get_symmetry_dataset

Get the symmetry dataset using spglib

find_K_from_k

Get the K vector of the supercell onto which the k vector of the primitive cell folds. The unfolding vector G, which satisfy the following equation, is also returned.

rotate_kpt

Apply rotation to a kpoint based on the rotations of the crystals (in the real space)

expand_K_by_symmetry

Expand the sampling of the PC kpoints due to symmetry breaking of the supercell cell.

LorentzSmearing

Simulate the Delta function by a Lorentzian shape function $\( \Delta(x) = \lim_{\sigma\to 0} Lorentzian \)$

GaussianSmearing

Simulate the Delta function by a Lorentzian shape function $\( \Delta(x) = \lim_{\sigma\to 0} Gaussian \)$

make_kpath

Return a list of kpoints defining the path between the given kpoints.

clean_latex_string

Clean up latex labels and convert if necessary

spectral_function_from_weight_sets

Generate the spectral function

spectral_weight_multiple_source

Calculate the spectral weight for a list of kpoints in the primitive BZ from a list of wave function files.

concatenate_scf_kpoints

Concatenate SCF kpoints (from IBZKPT) with zero-weighted kpoints

create_white_colormap

Create colormap from white to certain colour.

create_white_colormap_from_existing

Create a white-based color map from an existing one.

parse_atoms_idx

Expanding syntax like 1-2 (inclusive)

process_projection_options

Process commandline-style specifications for project

read_poscar_contcar_if_present

Return an ase Atoms() object of the POSCAR or CONTCAR file if present in the current directory.

parse_atoms

Parse the specified atoms (and orbitals if set) from a comma-separated string (e.g. “Na,Bi”) into a list of strings (e.g. [“Na”, “Bi”]), as well as a list of the corresponding atom indices in the structure and the parse orbital projections.

API#

easyunfold.unfold.get_symmetry_dataset(atoms: ase.Atoms, symprec: float = 1e-05)#

Get the symmetry dataset using spglib

easyunfold.unfold.find_K_from_k(k: numpy.ndarray, M: numpy.ndarray)#

Get the K vector of the supercell onto which the k vector of the primitive cell folds. The unfolding vector G, which satisfy the following equation, is also returned.

\[ \vec{k} = \vec{K} + \vec{G} \]

where G is a reciprocal space vector of supercell. NOTE: M is the transformation matrix for the cell matrix consistent of row vectors!

easyunfold.unfold.rotate_kpt(k: numpy.ndarray, opt: numpy.ndarray)#

Apply rotation to a kpoint based on the rotations of the crystals (in the real space)

NOTE: The rotation matrix should be the one that act on fractional coordinates, e.g. from spglib.

easyunfold.unfold.expand_K_by_symmetry(kpt: Union[list, numpy.ndarray], opts_pc: numpy.ndarray, opts_sc: numpy.ndarray, time_reversal: bool = True)#

Expand the sampling of the PC kpoints due to symmetry breaking of the supercell cell.

Returns:

Expanded kpoints and corresponding weights for each primitive cell kpoint.

class easyunfold.unfold.UnfoldKSet(M: numpy.ndarray, kpts_pc: list, pc_latt: numpy.ndarray, pc_opts: numpy.ndarray, sc_opts: numpy.ndarray, time_reversal: bool = True, expand: bool = True, metadata: Union[None, dict] = None, expansion_results: Union[None, dict] = None, calculated_quantities: Union[None, dict] = None, kpoint_labels: Union[None, list] = None, dft_code='vasp')#

Bases: monty.json.MSONable

High-level interface for unfolding and serialization.

Initialization

Instantiate an UnfoldKSet object.

Parameters:
  • M – The supercell transformation matrix

  • kpts_pc – A list of kpoints in the PC

  • pc_latt – A 3x3 matrix of row lattice vectors of the primitive cell

  • pc_opts – Symmetry operations of the primitive cell

  • sc_opts – Symmetry operations of the supercell

  • time_reversal – Whether to assume time-reversal symmetry or not

  • expand – Whether to expand the kpoint to take account of broken symmetry or not

  • expansion_results – Using existing results of symmetry expansion

  • calculated_quantities – Existing calculated quantities

  • kpoint_labels – Labels of the kpoints as passed in pc_latt as a list of tuples (<idx>, <label>)

  • dft_code – Name of the DFT code to be used.

_VERSION#

‘0.1.0’

check_version()#

Check the version of the program

property is_calculated: bool#

Show the status of the work

property has_averaged_spectral_weights: bool#

Return True if the spectral weights stored is averaged

classmethod from_atoms(M: numpy.ndarray, kpts_pc: list, pc: ase.Atoms, sc: ase.Atoms, time_reversal: bool = True, expand=True, symprec: float = 1e-05, dft_code='vasp')#

Initialise from primitive cell and supercell atoms

Parameters:
  • M – The supercell transformation matrix

  • kpts_pc – A list of kpoints in the PC

  • pc – Primitive cell structure

  • sc – Supercell structure

  • time_reversal – Whether to assume time-reversal symmetry or not

  • expand – Whether to expand the kpoint to take account of broken symmetry or not

  • symprec – Symmetry detection precision

classmethod from_file(fname: str)#

Load from a file

expand_pc_kpoints() None#

Compute the pc kpoints to be unfolded into

__repr__() str#
property nsymm_orig: int#

Number of symmetry operation in the original cell

property nsymm_expand: int#

Number of symmetry operation in the original cell

property nkpts_orig: int#

Total number of unexpanded kpoints

property nkpts_expand: int#

Total number of expanded kpoints

generate_sc_kpoints() None#

Generate the supercell kpoints to be calculated. Results are stored into self.expansion_results.

write_sc_kpoints(file: str, nk_per_split: Union[None, list] = None, scf_kpoints_and_weights: Union[None, list] = None, use_separate_folders=False, **kwargs)#

Write the supercell kpoints to a file.

Parameters:
  • file – Name of the file

  • nk_per_split – Number of kpoints per split along the path

  • scf_kpoints_and_weights – SCF kpoint and their weights needed for split-path calculations

write_pc_kpoints(file: str, expanded: bool = False, **kwargs)#

Write the primitive cell kpoints

_read_weights(wavefunction: Union[str, List[str]], gamma: bool, ncl: bool, gamma_half: str)#

Read the weights from the wave function files for all kpoints

Returns the averaged weights and the original weights per set of kpoints

load_procars(procars: Union[str, List[str]])#

Read in PROCAR for band-based projection

_construct_procar_kmap() list#

Construct mapping from the set of kpoints to that in the PROCAR

property procar: Union[None, easyunfold.procar.Procar]#

Loaded PROCARS

property procar_kmaps#

Loaded PROCAR kpoint mapping

_get_spectral_weights(wavefunction, npoints=2000, sigma=0.01, emin=None, emax=None, gamma=False, ncl=False, gamma_half='x', also_spectral_function=False, atoms_idx=None, orbitals=None, symm_average=True)#

Fetch spectral weights from a wavecar and compute spectral function is requested

Args:

atomic_projects (tuple): A tuple of atoms and orbitals whose projected weights should be used.
get_band_weight_sets(atoms_idx: List[int], orbitals: List[Union[List[str], str]], procars: Union[None, List[str], str] = None) list#

Get weights array sets for bands

Construct the weights of each band in same format of the kpoint set. Each item is a numpy array of (nspins, nbands), containing the summed weights over the passed atom indices and orbitals.

Parameters:
  • atoms_idx – Indices of the atoms to be selected

  • orbitals – Orbitals to be selected for each atom

  • procars – Names of the PROCAR files to be loaded

Returns:

A list of weights for each band at each expanded kpoint

get_spectral_function(wavefunction: Union[None, List[str]] = None, npoints: int = 2000, sigma: float = 0.1, gamma: bool = False, ncl: bool = False, gamma_half: str = 'x', symm_average: bool = True, atoms_idx: Union[None, List[int]] = None, orbitals: Union[None, str, List[str]] = None)#

Compute and return the spectral function

Parameters:
  • wavecar – The wavefunction files to be used

  • npoints – Number of points for the energy axis

  • sigma – Smearing width for the Gaussian smearing

  • gamma – Need to be set to True for Gamma-only calculation

  • ncl – Need to be set to True for non-collinear magnetism calculation

  • gamma_half – Flag used for reading WAVECAR

  • symm_average – Whether to perform symmetry averaging

  • atoms_idx – Indices for the atoms for projection

  • orbitals – Orbitals of the atoms for projection

Returns:

A tuple of the energies and the spectral functioin.

get_spectral_weights(wavefunction=None, gamma: bool = False, ncl: bool = False, gamma_half: str = 'x', symm_average: bool = True)#

Return the spectral weights calculated

Parameters:
  • wavefunction – The wavefunction file(s) for calculation. Not need if the weights has been calculated.

  • gamma – Whether the wavefunction files are from \(\Gamma\)-only calculation.

Returns:

An array storing the spectral weights.

as_dict() dict#

Convert the object into a dictionary representation

get_kpoint_distances()#

Distances between the kpoints along the path in the reciprocal space. This does not take account of the breaking of the path.

Note

The reciprocal lattice vectors includes the \(2\pi\) factor, e.g. np.linalg.inv(L).T * 2 * np.pi.

easyunfold.unfold.LorentzSmearing(x, x0, sigma=0.02)#

Simulate the Delta function by a Lorentzian shape function $\( \Delta(x) = \lim_{\sigma\to 0} Lorentzian \)$

easyunfold.unfold.GaussianSmearing(x, x0, sigma=0.02)#

Simulate the Delta function by a Lorentzian shape function $\( \Delta(x) = \lim_{\sigma\to 0} Gaussian \)$

easyunfold.unfold.make_kpath(kbound: List[float], nseg=40)#

Return a list of kpoints defining the path between the given kpoints.

easyunfold.unfold.clean_latex_string(label: str)#

Clean up latex labels and convert if necessary

Returns:

Cleaned tag string

easyunfold.unfold.spectral_function_from_weight_sets(spectral_weight_sets: numpy.ndarray, kweight_sets: list, nedos: int = 4000, sigma: float = 0.02, emin=None, emax=None, band_weight_sets=None)#

Generate the spectral function

\[ A(k_i, E) = \sum_m P_{Km}(k_i)\Delta(E - Em) \]

Where the \Delta function can be approximated by Lorentzian or Gaussian function.

Parameters:

band_weight_sets – Additional weighting for each band, used for generating projection onto atomic orbitals.

class easyunfold.unfold.Unfold(M=None, fname: str = 'WAVECAR', gamma: bool = False, lsorbit: bool = False, gamma_half: str = 'x', verbose=False, time_reversal=False, dft_code='vasp')#

Low level interface for performing unfolding related calculations. obtain the effective band structure (EBS).

Reference

“Extracting E versus k effective band structure from supercell calculations on alloys and impurities” Phys. Rev. B 85, 085201 (2012)

Initialization

Initialization.

M is the transformation matrix between supercell and primitive cell:

M = np.dot(A, np.linalg.inv(a))

In real space, the basis vectors of Supercell (A) and those of the primitive cell (a) satisfy:

A = np.dot(M, a);      a = np.dot(np.linalg.inv(M), A)

Whereas in reciprocal space

b = np.dot(M.T, B);    B = np.dot(np.linalg.inv(M).T, b)
get_vbm_cbm(thresh: float = 1e-08) Tuple[float, float]#

Locate the VBM from the wave function data

get_ovlap_G(ikpt: int = 1, epsilon: float = 1e-05) Tuple[numpy.ndarray, numpy.ndarray]#

Get subset of the reciprocal space vectors of the supercell, specifically the ones that match the reciprocal space vectors of the primitive cell.

find_K_index(K0: numpy.ndarray) int#

Find the (one-based) index of a point K0.

spectral_weight_k(k0, whichspin=1)#

Spectral weight for a given \(k\):

\[ P_{Km}(k) = \sum_n |<Km | kn>|^2 \]

which is equivalent to $\( P_{Km}(k) = \sum_{G} |C_{Km}(G + k - K)|^2 \)$

where \(G\) is a subset of the reciprocal space vectors of the supercell.

spectral_weight(kpoints: List)#

Calculate the spectral weight for a list of kpoints in the primitive BZ.

spectral_function(nedos: int = 4000, sigma: float = 0.02)#

Generate the spectral function

\[ A(k_i, E) = \sum_m P_{Km}(k_i)\Delta(E - Em) \]

Where the \(\Delta\) function can be approximated by Lorentzian or Gaussian function.

easyunfold.unfold.spectral_weight_multiple_source(kpoints: list, unfold_objs: List[easyunfold.unfold.Unfold], transform_matrix: numpy.ndarray)#

Calculate the spectral weight for a list of kpoints in the primitive BZ from a list of wave function files.

easyunfold.unfold.concatenate_scf_kpoints(scf_kpts: list, scf_weights: list, kpoints: list)#

Concatenate SCF kpoints (from IBZKPT) with zero-weighted kpoints

easyunfold.unfold.create_white_colormap(color: Union[str, tuple, list]) matplotlib.colors.ListedColormap#

Create colormap from white to certain colour.

Args: color (str, tuple, list): HEX color string or tuple/list of RGB color.

easyunfold.unfold.create_white_colormap_from_existing(name: str) matplotlib.colors.ListedColormap#

Create a white-based color map from an existing one.

easyunfold.unfold.parse_atoms_idx(atoms_idx: str) List[int]#

Expanding syntax like 1-2 (inclusive)

For example, 1,2,3,4-6 will be expanded as 1,2,3,4,5,6.

Parameters:

atoms_idx – A string encode atom indices

Returns:

A list of atom indices

easyunfold.unfold.process_projection_options(atoms_idx: str, orbitals: str) Tuple[list, list]#

Process commandline-style specifications for project

Parameters:
  • atoms_idx – A comma- or hyphen-separated string of atom projections

  • orbitals – A comma-separated string of orbital projections

Returns:

A tuple of atom indices and the orbitals selected for projection.

easyunfold.unfold.read_poscar_contcar_if_present(poscar: str = 'POSCAR')#

Return an ase Atoms() object of the POSCAR or CONTCAR file if present in the current directory.

Returns:

ASE Atoms() object

easyunfold.unfold.parse_atoms(atoms_to_project: str, orbitals: str, poscar: str)#

Parse the specified atoms (and orbitals if set) from a comma-separated string (e.g. “Na,Bi”) into a list of strings (e.g. [“Na”, “Bi”]), as well as a list of the corresponding atom indices in the structure and the parse orbital projections.

Parameters:
  • atoms_to_project – A comma-separted string of atom symbols to project

  • orbitals – A “|”-separated string of orbital projections

  • poscar – The POSCAR file to read atom indices from

Returns:

A tuple of lists of atoms, atom indices and the orbitals selected for projection.