import numpy as np
from EPWpy.utilities import Band as bnd
from EPWpy.plotting import plot_bands as plb
from EPWpy.QE import PW_util as pwu
from EPWpy.utilities.k_gen import k_mesh
from EPWpy.utilities.read_QE_xml import QE_XML
[docs]
class BandUtil:
"""
Utility class for handling band structure data from Quantum ESPRESSO (QE) outputs.
This class provides convenient access to file and directory information
associated with band structure, phonon, and EPW calculations, allowing
downstream tools to locate and process relevant outputs.
Attributes
----------
system : str
Name of the system or material (used as the main working directory).
prefix : str
Common prefix used in Quantum ESPRESSO input/output files.
scf_file : str
Name of the self-consistent field (SCF) input/output file (default: 'scf').
scf_fold : str
Folder containing SCF calculation results (default: 'scf').
nscf_file : str
Name of the non-self-consistent field (NSCF) input/output file (default: 'nscf').
nscf_fold : str
Folder containing NSCF calculation results (default: 'nscf').
ph_file : str
Name of the phonon input/output file (default: 'ph').
ph_fold : str
Folder containing phonon calculation results (default: 'ph').
bs_file : str
Name of the band structure input/output file (default: 'bs').
bs_fold : str
Folder containing band structure calculation results (default: 'bs').
epw_file : str
Name of the EPW input/output file (default: 'epw').
epw_fold : str
Folder containing EPW calculation results (default: 'epw').
Examples
--------
>>> bands = BandUtil(system='si', prefix='si')
>>> print(bands.bs_fold)
'bs'
>>> print(bands.scf_file)
'scf'
"""
def __init__(self, system: str, prefix: str,
scf_file: str = 'scf',
scf_fold: str = 'scf',
nscf_file: str = 'nscf',
nscf_fold: str = 'nscf',
ph_file: str = 'ph',
ph_fold: str = 'ph',
bs_fold: str = 'bs',
bs_file: str = 'bs',
epw_fold: str = 'epw',
epw_file: str = 'epw'):
"""
Initialize the BandUtil class with the system and directory structure.
Parameters
----------
system : str
Name of the material system.
prefix : str
File prefix used in Quantum ESPRESSO outputs.
scf_file, scf_fold, nscf_file, nscf_fold, ph_file, ph_fold, bs_file, bs_fold, epw_file, epw_fold : str, optional
File and directory names for SCF, phonon, band structure, and EPW runs.
"""
self.system = system
self.prefix = prefix
self.scf_file = scf_file
self.scf_fold = scf_fold
self.nscf_file = nscf_file
self.nscf_fold = nscf_fold
self.ph_file = ph_file
self.ph_fold = ph_fold
self.bs_file = bs_file
self.bs_fold = bs_fold
self.epw_file = epw_file
self.epw_fold = epw_fold
[docs]
def get_band_QE(self, file: str = None):
"""
Retrieve the band structure data from a Quantum ESPRESSO output file.
This method reads the QE band structure output (typically from ``bs.out``)
and extracts eigenvalues and k-points using the ``bnd.extract_band_scf``
utility function.
Parameters
----------
file : str, optional
Path to the QE bandstructure output file.
If not provided, it defaults to ``{system}/{bs_fold}/{bs_file}.out``.
Returns
-------
Band_QE : object
Extracted band structure data with shape ``(nk, nband)``.
Examples
--------
>>> bu = BandUtil(system='si', prefix='si')
>>> bands = bu.get_band_QE()
>>> print(bands.eigenvalues.shape)
(50, 10)
"""
if file is None:
file = f'{self.system}/{self.bs_fold}/{self.bs_file}.out'
Band_QE = bnd.extract_band_scf(f'{file}')
return(Band_QE)
[docs]
def get_band_EPW(self, file: str = None):
"""
Retrieve the interpolated band structure from an EPW calculation.
This method reads the EPW-generated ``band.eig`` file and extracts
eigenvalues and k-points using the ``bnd.extract_band_eig`` utility function.
Parameters
----------
file : str, optional
Path to the EPW band structure eigenvalue file.
Defaults to ``{system}/{epw_fold}/band.eig``.
Returns
-------
Band_EPW : object
Extracted EPW band structure data with shape ``(nk, nband)``.
Examples
--------
>>> bu = BandUtil(system='si', prefix='si')
>>> bands_epw = bu.get_band_EPW()
>>> print(bands_epw.eigenvalues[:5])
[ -5.78, -5.77, -5.75, -5.70, -5.68 ]
"""
if file is None:
file = f'{self.system}/{self.epw_fold}/band.eig'
Band_EPW = bnd.extract_band_eig(f'{file}')
return(Band_EPW)
[docs]
def get_phonon_QE(self, file: str = None):
"""
Retrieve the phonon dispersion data from Quantum ESPRESSO.
This method reads the phonon frequency data from the ``matdyn.freq`` file
(typically produced by ``matdyn.x``) and converts frequencies from meV
to THz using a conversion factor of 0.124.
Parameters
----------
file : str, optional
Path to the phonon frequency file.
Defaults to ``{system}/{ph_fold}/{prefix}.freq``.
Returns
-------
ph_band : ndarray
Extracted phonon band frequencies with shape ``(nq, nmode)``,
expressed in meV.
Notes
-----
The conversion factor ``0.124`` converts from cm^-1 → meV
Examples
--------
>>> bu = BandUtil(system='si', prefix='si')
>>> ph = bu.get_phonon_QE()
>>> print(ph.shape)
(120, 6)
"""
if file is None:
file = f'{self.system}/{self.ph_fold}/{self.prefix}.freq'
ph_band = bnd.extract_band_freq(f'{file}')
return(ph_band*0.124)
[docs]
def get_phonon_EPW(self, file: str = None):
"""
Retrieve the phonon dispersion data interpolated by EPW.
This method reads the phonon frequency file ``phband.freq`` generated by EPW,
which contains interpolated phonon frequencies along a specified q-path.
Parameters
----------
file : str, optional
Path to the EPW phonon frequency file.
Defaults to ``{system}/{epw_fold}/phband.freq``.
Returns
-------
epw_ph_band : ndarray
Extracted EPW phonon band frequencies with shape ``(nq, nmode)`` in THz.
Examples
--------
>>> bu = BandUtil(system='si', prefix='si')
>>> epw_ph = bu.get_phonon_EPW()
>>> print(epw_ph.shape)
(200, 6)
"""
if file is None:
file = f'{self.system}/{self.epw_fold}/phband.freq'
epw_ph_band = bnd.extract_band_eig(f'{file}')
return(epw_ph_band)
[docs]
def plot_band_QE(self, file: str = None, **kwargs):
"""
Plot the electronic band structure obtained from Quantum ESPRESSO.
This method reads the Quantum ESPRESSO bandstructure output file
(typically ``bs.out``), extracts the band data, aligns it with the Fermi
level, and plots the band structure using the built-in plotting utilities.
Parameters
----------
file : str, optional
Path to the Quantum ESPRESSO bandstructure output file.
Defaults to ``{system}/{bs_fold}/{bs_file}.out``.
**kwargs : dict, optional
Additional keyword arguments passed to the plotting function
``plb.plot_band_prod`` (e.g., color, linewidth, figsize).
Supported keys include:
- ``xlabel`` : str, label for x-axis (default: "Wavevector")
- ``ylabel`` : str, label for y-axis (default: "Energy (eV)")
Notes
-----
- The Fermi level is automatically extracted from the SCF output file
(``{system}/{scf_fold}/{scf_file}.out``).
- The energy axis is shifted such that the Fermi level is set to 0 eV.
Example
-------
>>> bu = BandUtil(system='si', prefix='si')
>>> bu.plot_band_QE(color_c='blue', linewidth=1.2)
"""
if 'xlabel' not in kwargs.keys():
xlabel = 'Wavevector'
if 'ylabel' not in kwargs.keys():
ylabel = 'Energy (eV)'
Band_QE = self.get_band_QE(file = file)
ef = pwu.find_fermi(f'{self.system}/{self.scf_fold}/{self.scf_file}.out')
plb.plot_band_prod(Band_QE,
ef0=ef,
**kwargs)
[docs]
def plot_band_EPW(self, file: str = None, **kwargs):
"""
Plot the interpolated electronic band structure obtained from EPW.
This method reads the EPW-interpolated bandstructure data
(typically from ``band.eig``), aligns the bands to the Fermi level,
and visualizes the interpolated dispersion using the internal plotting
utility.
Parameters
----------
file : str, optional
Path to the EPW bandstructure file.
Defaults to ``{system}/{epw_fold}/band.eig``.
**kwargs : dict, optional
Additional keyword arguments passed to the plotting function
``plb.plot_band_prod`` (e.g., color, linewidth, figsize).
Supported keys include:
- ``xlabel`` : str, label for x-axis (default: "Wavevector")
- ``ylabel`` : str, label for y-axis (default: "Energy (eV)")
Notes
-----
- The band structure is read from EPW’s interpolation output file
(``band.eig``), typically generated after Wannier interpolation.
- The Fermi level is automatically extracted from the SCF output file
(``{system}/{scf_fold}/{scf_file}.out``).
- The energy axis is shifted such that the Fermi level corresponds to 0 eV.
Example
-------
>>> bu = BandUtil(system='graphene', prefix='gr')
>>> bu.plot_band_EPW(color_c='red', linewidth=1.0, xlabel='k-path', ylabel='Energy (eV)')
"""
if 'xlabel' not in kwargs.keys():
kwargs['xlabel'] = 'Wavevector'
if 'ylabel' not in kwargs.keys():
kwargs['ylabel'] = 'Energy (eV)'
Band_EPW = self.get_band_EPW(file = file)
ef = pwu.find_fermi(f'{self.system}/{self.scf_fold}/{self.scf_file}.out')
plb.plot_band_prod(Band_EPW,
ef0=ef,
**kwargs)
[docs]
def plot_phonon_QE(self, file: str = None, **kwargs):
"""
Plot the phonon dispersion obtained from Quantum ESPRESSO.
This method reads the phonon frequencies calculated by QE
(typically from ``matdyn.freq`` or ``{prefix}.freq``),
converts them to meV, and visualizes the phonon bandstructure.
Parameters
----------
file : str, optional
Path to the QE phonon frequency file.
Defaults to ``{system}/{ph_fold}/{prefix}.freq``.
**kwargs : dict, optional
Additional keyword arguments passed to the plotting function
``plb.plot_band_prod`` (e.g., color, linewidth, figsize).
Supported keys include:
- ``xlabel`` : str, label for x-axis (default: "Wavevector")
- ``ylabel`` : str, label for y-axis (default: "Energy (meV)")
Notes
-----
- The phonon frequencies are extracted from QE’s ``matdyn.freq`` output.
- Frequencies are converted to meV using the factor 0.124.
- The zero energy reference is set at 0 meV.
Example
-------
>>> bu = BandUtil(system='graphene', prefix='gr')
>>> bu.plot_phonon_QE(color='blue', linewidth=1.2)
"""
if 'xlabel' not in kwargs.keys():
kwargs['xlabel'] = 'Wavevector'
if 'ylabel' not in kwargs.keys():
kwargs['ylabel'] = 'Energy (meV)'
ph_QE = self.get_phonon_QE(file = file)
plb.plot_band_prod(ph_QE,
ef0=0,
**kwargs)
[docs]
def plot_phonon_EPW(self, file: str = None, **kwargs):
"""
Plot the phonon dispersion interpolated using EPW.
This method reads the EPW-generated phonon dispersion file
(typically ``phband.freq``) and visualizes the interpolated
phonon bandstructure.
Parameters
----------
file : str, optional
Path to the EPW phonon frequency file.
Defaults to ``{system}/{epw_fold}/phband.freq``.
**kwargs : dict, optional
Additional keyword arguments passed to the plotting function
``plb.plot_band_prod`` (e.g., color, linewidth, figsize).
Supported keys include:
- ``xlabel`` : str, label for x-axis (default: "Wavevector")
- ``ylabel`` : str, label for y-axis (default: "Energy (meV)")
Notes
-----
- The phonon frequencies are obtained from EPW’s ``phband.freq`` output.
- Frequencies are assumed to be in meV.
- The zero energy reference is set at 0 meV.
Example
-------
>>> bu = BandUtil(system='graphene', prefix='gr')
>>> bu.plot_phonon_EPW(color='red', linewidth=1.2)
"""
if 'xlabel' not in kwargs.keys():
kwargs['xlabel'] = 'Wavevector'
if 'ylabel' not in kwargs.keys():
kwargs['ylabel'] = 'Energy (meV)'
ph_EPW = self.get_phonon_EPW(file = file)
plb.plot_band_prod(ph_EPW,
ef0=0,
**kwargs)
[docs]
def get_eff_mass(self, grid, file: str = None, nmax: int = 5):
"""
Get the effective mass tensor using FT method
nmax = nmax
nk1 = grid[0]
nk2 = grid[1]
nk3 = grid[2]
if file is None:
file = f'{self.system}/{self.nscf_fold}/{self.nscf_file}.out'
band = self.get_band_QE(f'{file}')
kx=[1,0,0]
ky=[0,1,0]
kz=[0,0,0]
h=np.linspace(0.0,1,nk1,endpoint=False)
k=np.linspace(0.0,1,nk2,endpoint=False)
l=np.linspace(0.0,1,nk3,endpoint=False)
kpoints = np.array(k_mesh(h,k,l,kx,ky,kz))
print(np.shape(band),np.shape(kpoints))
interp_func, coeffs, Gints = fourier_fit(band,kpoints[:,:3], nmax = nmax)
qe = QE_XML("si/bs/si.xml") # or your file name
qe.summary()
kpts = qe.get_kpoints()
eigs = qe.get_eigenvalues()
cell = qe.get_lattice()
kpoints_crys = kpoints_2pi_alat_to_fractional(kpts, cell, 10.262000)
band_dft = extract_band_scf('si/bs/bs.out')
band = []
for kpt in kpoints_crys:
band.append(interp_func(kpt))
band = np.array(band)
for ib in range(len(band[0,0,:])):
plt.plot(band[:,0,ib],color = 'k')
plt.plot(band_dft[:,ib],color = 'crimson')
bvecs = recoproc_lat_cart(cell)
min_arg = np.argmin(band[:,0,4])
print(effective_mass_tensor(kpoints_crys[min_arg,:], coeffs, Gints, cell*bohr2m, band = 3))
print(np.shape(band))
"""
pass