#
"""
QE class for interfacing with quantum espresso
"""
from __future__ import annotations
__author__= "Sabyasachi Tiwari"
__copyright__= "Copyright 2024, EPWpy project"
__version__= "1.0"
__maintainer__= "Sabyasachi Tiwari"
__maintainer_email__= "sabyasachi.tiwari@austin.utexas.edu"
__status__= "Production"
__date__= "May 03, 2024"
import numpy as np
import argparse
import os
import glob
from EPWpy.utilities.set_files import SetDir
from EPWpy.default.default import *
from EPWpy.QE.write_QE import WriteQEFiles
from EPWpy.structure.lattice import *
from EPWpy.utilities.k_gen import k_mesh, k_path
from EPWpy.utilities.epw_pp import run_pp
from EPWpy.utilities.save_state import *
from EPWpy.utilities.EPW_util import *
from EPWpy.QE.PW_util import *
from EPWpy.QE.set_QE import *
from EPWpy.QE.PH_util import *
from EPWpy.utilities.printing import *# printing
[docs]
class PyQE(SetDir, SetQE, WriteQEFiles):
"""
High-level interface for building Quantum ESPRESSO and EPW input files.
The `PyQE` class automates setup and input generation for Quantum ESPRESSO
and EPW simulations. It inherits functionality from three base classes:
- **SetDir** → manages directory structure for calculations
- **SetQE** → provides default Quantum ESPRESSO parameters and configuration
- **WriteQEFiles** → writes the actual QE and EPW input files
This class can be used **standalone**, independent of EPWpy, to prepare
Quantum ESPRESSO workflows (e.g., SCF, NSCF, Phonon, or EPW inputs).
Attributes
----------
system : str
Name of the material or simulation system (used to create a working directory).
code : str
Path to the compiled Quantum ESPRESSO executables (default: current directory).
home : str
The current working directory when the class is initialized.
env : str
Parallel execution environment (e.g., ``"mpirun"``, ``"srun"``, or none).
"""
def __init__(self, type_input: dict, system: str = 'si', code: str = '.', env: str = 'mpirun'):
"""
Initialize the `PyQE` class and prepare the environment.
This constructor initializes Quantum ESPRESSO input setup by creating
a system folder, loading default parameters, and applying user-defined
settings. It also reads any master configuration file if present.
Parameters
----------
type_input : dict
Dictionary specifying the type of calculation and related parameters.
Example:
``{"type": "scf", "ecutwfc": 50, "kpoints": [6,6,6]}``
system : str, optional
Name of the material or system to simulate (default: ``"si"``).
code : str, optional
Path to the Quantum ESPRESSO executables (default: ``"."``).
env : str, optional
Execution command for parallel runs, such as ``"mpirun"`` or ``"srun"`` (default: ``"mpirun"``).
Notes
-----
- Automatically creates a folder named after ``system`` if it doesn’t exist.
- Loads default Quantum ESPRESSO parameters via ``default_values()``.
- Applies calculation-specific settings from ``set_values()``.
- Reads master configuration data (if available) using ``_read_master()``.
Example
-------
>>> from epwpy import PyQE
>>> params = {"type": "scf", "ecutwfc": 60, "kpoints": [8,8,8]}
>>> qe = PyQE(params, system="graphene", env="mpirun")
>>> qe.write_input() # Generates SCF input file for Quantum ESPRESSO
"""
self.system = system
self.code = code
self.home = os.getcwd()
self.env = env
os.system('mkdir ' + self.system)
self.default_values()
self.set_values(type_input)
self._read_master()
[docs]
def reset(self):
"""
Reset the working environment for a new calculation.
This method restores the internal state of the `PyQE` instance
to prepare for a fresh calculation run. It performs the following steps:
1. Calls internal ``set_work()`` to configure the working directory.
2. Reloads master configuration data via ``_read_master()``.
3. Returns to the original home directory using ``set_home()``.
Notes
-----
Use this method before starting a new calculation if the
previous one modified the working state (e.g., after an SCF or NSCF run).
Example
-------
>>> qe.reset()
>>> qe.run_scf()
"""
print("Resetting environment: obtaining NSCF and PH attributes")
self.set_work()
self._read_master()
self.set_home()
[docs]
def custom(self, parameters: dict = None, name: str = 'custom'):
"""
Run a fully custom Quantum ESPRESSO calculation.
This method enables running a calculation that is **not predefined**
(e.g., not SCF, NSCF, PH, or EPW). The user can directly specify
input parameters in the form of a dictionary.
Parameters
----------
parameters : dict, optional
Custom Quantum ESPRESSO parameters to include in the input file.
If not provided, defaults to an empty dictionary.
name : str, optional
Name of the calculation type (used for naming input files).
Default is ``"custom"``.
Notes
-----
- Creates a custom file in working system directory.
- Applies the provided parameters via ``set_custom()``.
- Writes the custom input file with ``write_custom()``.
- Returns to the home directory using ``set_home()``.
Example
-------
>>> custom_input = {
... "calculation": "bands",
... "ecutwfc": 60,
... "kpoints": [10, 10, 10]
... }
>>> qe.custom(custom_input, name="bands")
"""
if parameters is None:
parameters = {}
self.set_work()
self.set_custom(parameters)
self.write_custom()
self.set_home()
[docs]
def scf(self,control=None,
system=None,
electrons=None,
ions=None,
cell=None,
kpoints=None,
cell_parameters=None,
hubbard=None,
fcp=None,
solvents=None,
occupations=None,
atomic_velocities=None,
constraints=None,
atomic_forces=None,
rism=None,
cards=None,
name='scf'):
"""
Prepare a **self-consistent field (SCF)** calculation for Quantum ESPRESSO.
This method assembles all necessary input sections for a standard SCF
calculation by merging user-supplied dictionaries with default parameters.
The input file is written automatically in the working directory.
Parameters
----------
control, system, electrons, ions, cell, kpoints, cell_parameters : dict, optional
Dictionaries corresponding to Quantum ESPRESSO namelists such as
``&CONTROL``, ``&SYSTEM``, ``&ELECTRONS``, etc.
hubbard, fcp, solvents, occupations, atomic_velocities, constraints,
atomic_forces, rism, cards : dict, optional
Optional dictionaries specifying less-common QE sections (e.g.,
``HUBBARD``, ``FCP``, ``SOLVENTS``). Each is optional and defaults to `{}`.
name : str, optional
Name of the calculation (used for input/output file naming).
Default is ``"scf"``.
Workflow
--------
1. Creates or moves to the working directory using ``set_work()``.
2. Calls ``pw()`` to populate QE input variables.
3. Writes the SCF input file using ``write_scf()``.
4. Stores the SCF filename in ``self.scf_file``.
5. Returns to the home directory using ``set_home()``.
Example
-------
>>> qe = PyQE({"type": "scf"}, system="graphene")
>>> qe.scf(
... control={"calculation": "scf", "prefix": "graphene"},
... system={"ecutwfc": 60, "occupations": "smearing"},
... kpoints={"type": "automatic", "grid": [12, 12, 1]}
... )
"""
if control is None:
control = {}
if system is None:
system = {}
if electrons is None:
electrons = {}
if ions is None:
ions = {}
if cell is None:
cell = {}
if kpoints is None:
kpoints = {}
if cell_parameters is None:
cell_parameters = {}
if hubbard is None:
hubbard = {}
if fcp is None:
fcp = {}
if solvents is None:
solvents = {}
if occupations is None:
occupations = {}
if atomic_velocities is None:
atomic_velocities = {}
if constraints is None:
constraints = {}
if atomic_forces is None:
atomic_forces = {}
if rism is None:
rism = {}
if cards is None:
cards = {}
self.set_work()
"""
set the pw values
"""
self.pw(control,
system,
electrons,
ions,
cell,
kpoints,
cell_parameters,
hubbard,
fcp,
solvents,
occupations,
atomic_velocities,
constraints,
atomic_forces,
rism,
cards)
self.write_scf(name=name)
self.scf_file = name
self.set_home()
[docs]
def ph(self,phonons=None,
qpoints=None,
name='ph'):
"""
Prepare a **phonon (PH)** calculation input for Quantum ESPRESSO.
This method sets up and writes the input files required for phonon
calculations using the `ph.x` module of Quantum ESPRESSO. It combines
user-provided phonon parameters and q-point settings with the default
configuration stored in the class.
Parameters
----------
phonons : dict, optional
Dictionary containing keywords for the phonon calculation, corresponding
to the ``ph.x`` input parameters. Example:
``{"tr2_ph": 1.0e-14, "nq1": 4, "nq2": 4, "nq3": 4}``
qpoints : dict, optional
Dictionary defining the q-point grid or specific q-points for the
phonon calculation. Example:
``{"q_grid": [4, 4, 4]}``
name : str, optional
Base name for the generated phonon input file. Default is ``"ph"``.
Workflow
--------
1. Prepares the working directory via ``set_work()``.
2. Applies phonon and q-point parameters using ``set_phonons()``.
3. Writes the phonon input file using ``write_ph()``.
4. Saves a phonon configuration snapshot with ``save_PH()``.
5. Returns to the home directory via ``set_home()``.
Example
-------
>>> qe = PyQE({"type": "ph"}, system="graphene")
>>> qe.ph(
... phonons={"tr2_ph": 1.0e-14, "ldisp": True, "nq1": 4, "nq2": 4, "nq3": 1},
... qpoints={"grid": [4, 4, 1]},
... name="graphene_phonon"
... )
>>> print(qe.ph_file)
graphene_phonon
"""
if phonons is None:
phonons = {}
if qpoints is None:
qpoints = {}
self.set_work()
self.set_phonons(phonons,qpoints)
self.write_ph(name=name)
self.ph_file = name
self.save_PH()
self.set_home()
[docs]
def nscf(self,control=None,
system=None,
electrons=None,
ions=None,
cell=None,
kpoints=None,
cell_parameters=None,
hubbard=None,
fcp=None,
solvents=None,
occupations=None,
atomic_velocities=None,
constraints=None,
atomic_forces=None,
rism=None,
cards=None,
name='nscf'):
"""
Prepare a **non-self-consistent field (NSCF)** calculation input for Quantum ESPRESSO.
The NSCF step is typically performed after an SCF calculation to compute
electronic properties (e.g., band structure or Fermi surface) using a
denser k-point grid, without updating the charge density.
Parameters
----------
control : dict, optional
Dictionary of `&CONTROL` parameters for Quantum ESPRESSO.
system : dict, optional
Dictionary of `&SYSTEM` parameters, such as `occupations`, `ecutwfc`, etc.
electrons : dict, optional
Dictionary of `&ELECTRONS` parameters for electronic minimization.
ions : dict, optional
Dictionary of `&IONS` parameters (only relevant for relax/dynamics runs).
cell : dict, optional
Dictionary of `&CELL` parameters (for variable-cell calculations).
kpoints : dict, optional
Dictionary specifying the k-point grid or explicit list of k-points.
cell_parameters : dict, optional
Dictionary defining lattice vectors for non-cubic systems.
hubbard : dict, optional
Dictionary specifying DFT+U corrections.
fcp : dict, optional
Parameters for fixed-charge potential calculations.
solvents : dict, optional
Parameters for implicit solvent models (e.g., ENVIRON, SCCS).
occupations : dict, optional
Occupation-related settings (e.g., smearing type, degauss).
atomic_velocities : dict, optional
Initial atomic velocities for molecular dynamics.
constraints : dict, optional
Constraints on atomic motion or structure optimization.
atomic_forces : dict, optional
External forces on atoms.
rism : dict, optional
Settings for solvation calculations using 3D-RISM.
cards : dict, optional
Additional QE cards (e.g., ATOMIC_POSITIONS, K_POINTS, CELL_PARAMETERS).
name : str, optional
Base name for the NSCF input file (default: ``"nscf"``).
Workflow
--------
1. Sets up the working directory using ``set_work()``.
2. Forces the calculation type to ``'nscf'`` in ``pw_control``.
3. Applies all input dictionaries using ``pw()``.
4. Writes the NSCF input file via ``write_scf()``.
5. Saves current PW configuration using ``save_PW()``.
6. Returns to the home directory using ``set_home()``.
Example
-------
>>> qe = PyQE({"type": "nscf"}, system="graphene")
>>> qe.nscf(
... control={"verbosity": "high"},
... system={"occupations": "tetrahedra"},
... kpoints={"grid": [12, 12, 1]},
... name="graphene_nscf"
... )
>>> print(qe.nscf_file)
graphene_nscf
"""
if control is None:
control = {}
if system is None:
system = {}
if electrons is None:
electrons = {}
if ions is None:
ions = {}
if cell is None:
cell = {}
if kpoints is None:
kpoints = {}
if cell_parameters is None:
cell_parameters = {}
if hubbard is None:
hubbard = {}
if fcp is None:
fcp = {}
if solvents is None:
solvents = {}
if occupations is None:
occupations = {}
if atomic_velocities is None:
atomic_velocities = {}
if constraints is None:
constraints = {}
if atomic_forces is None:
atomic_forces = {}
if rism is None:
rism = {}
if cards is None:
cards = {}
self.set_work()
self.pw_control['calculation']='\'nscf\''
self.pw(control,
system,
electrons,
ions,
cell,
kpoints,
cell_parameters,
hubbard,
fcp,
solvents,
occupations,
atomic_velocities,
constraints,
atomic_forces,
rism,
cards)
self.write_scf(name)
self.nscf_file = name
self.save_PW()
self.set_home()
[docs]
@decorated_warning
def epw(self, epwin=None, name='epw'):
"""
Prepare an **EPW (Electron–Phonon Wannier)** calculation input.
This method generates the input file required for EPW, a post-processing
tool of Quantum ESPRESSO used for calculating electron–phonon coupling,
transport, superconductivity, and polaron properties.
The routine automatically detects previously generated Quantum ESPRESSO
(PW and PH) files, sets up the EPW environment, and writes the input
based on user-defined or default parameters.
Parameters
----------
epwin : dict, optional
Dictionary containing EPW-specific input parameters.
Example:
``{"ep_coupling": True, "wannierize": True, "degaussw": 0.01}``
name : str, optional
Base name for the EPW input file (default: ``"epw"``).
Workflow
--------
1. Sets up the working directory using ``set_work()``.
2. Reads stored Quantum ESPRESSO (PW) and Phonon (PH) data if available.
3. Determines whether to restart from a previous EPW run or start fresh:
- If ``name`` is ``"epw1"``, ``"fbw"``, or ``"fbw_mu"``, restart is disabled.
- Otherwise, EPW restart is attempted by loading saved EPW data.
4. Applies EPW parameters using ``set_epw()``.
5. Writes the EPW input file using ``write_epw()``.
6. Saves configuration with ``save_EPW()`` and returns to home directory.
Notes
-----
- The method attempts to read PW and PH data through ``read_PW()`` and ``read_PH()``.
- Missing files trigger a warning but do not stop execution.
- EPW restart handling depends on the provided input name.
Example
-------
>>> qe = PyQE({"type": "epw"}, system="graphene")
>>> epw_params = {
... "elph": '.true.',
... "wannierize": '.true.',
... "degaussw": 0.02,
... "nkf1": 60, "nkf2": 60, "nkf3": 1
... }
>>> qe.epw(epwin=epw_params, name="epw1")
>>> print(qe.epw_file)
epw1
"""
if epwin is None:
epwin = {}
"""
EPW calculation preparation
"""
self.set_work()
self.epw_file = name
self.epw_refresh = None
try:
self.read_PW()
self.read_PH()
except FileNotFoundError:
print(f'PW PH files not found')
if((name =='epw1') or (name =='fbw') or (name =='fbw_mu')):
self.epw_restart = False
else:
self.epw_restart = True
try:
self.read_EPW()
except FileNotFoundError:
self.epw_refresh = None
self.set_epw(epwin)
self.write_epw(name)
self.save_EPW()
self.set_home()
[docs]
def QE_pp(self, inputpp=None, plot=None, name='pp'):
"""
Post-processing calculation of QE
"""
if inputpp is None:
inputpp = {}
if plot is None:
plot = {}
self.set_work()
self.set_QE_pp(inputpp,plot)
self.write_QE_pp(name = name)
self.QE_pp_file = name
self.set_home()
[docs]
def q2r(self, q2r=None, name='q2r'):
"""
Prepare a **Q2R (q-points to real-space)** calculation input.
The `q2r()` method generates input files for the Quantum ESPRESSO
``q2r.x`` utility, which transforms dynamical matrices from reciprocal
space (q-space) — computed via ``ph.x`` — into real-space force constants.
These real-space force constants are typically used as input for
``matdyn.x`` to compute phonon dispersions.
Parameters
----------
q2r : dict, optional
Dictionary containing user-defined parameters for ``q2r.x`` input.
Example:
``{"fildyn": "si.dyn", "flfrc": "si.fc"}``
name : str, optional
Base name for the generated Q2R input file (default: ``"q2r"``).
Workflow
--------
1. Initializes a working directory with ``set_work()``.
2. Applies provided ``q2r`` parameters using ``set_q2r()``.
3. Writes the Q2R input file via ``write_q2r()``.
4. Stores file metadata (``self.q2r_file``).
5. Returns to the home directory using ``set_home()``.
Notes
-----
- This method assumes that phonon calculations (``ph.x``) have already
been completed and that the dynamical matrix files are available.
- The resulting real-space force constants file is usually named
``*.fc`` and used by ``matdyn.x``.
Example
-------
>>> qe = PyQE({"type": "q2r"}, system="si")
>>> params = {"fildyn": "si.dyn", "flfrc": "si.fc"}
>>> qe.q2r(q2r=params, name="q2r_si")
>>> print(qe.q2r_file)
q2r_si
"""
if q2r is None:
q2r = {}
self.set_work()
self.set_q2r(q2r)
self.write_q2r(name=name)
self.q2r_file = name
self.set_home()
[docs]
def dvscf_q2r(self, dvscf_q2r=None, name='dvscf_q2r'):
if dvscf_q2r is None:
dvscf_q2r = {}
"""
dvscf_q2r calculation preparation
"""
self.set_work()
self.set_dvscf_q2r(dvscf_q2r)
self.write_dvscf_q2r(name=name)
self.dvscf_q2r_file = name
self.set_home()
[docs]
def postahc(self, postahc=None, name='postahc'):
if postahc is None:
postahc = {}
"""
postahc calculation preparation
"""
self.set_work()
self.set_postahc(postahc)
self.write_postahc(name=name)
self.postahc_file = name
self.set_home()
[docs]
def zg(self, zg=None, azg=None, name='zg'):
"""
Prepare a **ZG (zone-averaged or auxiliary calculation)** input.
The `zg()` method generates input files for the Quantum ESPRESSO `zg.x` utility
(or any custom ZG post-processing module). This is typically used for
zone-averaged properties, auxiliary phonon data, or other custom calculations
derived from prior PW/PH runs.
Parameters
----------
zg : dict, optional
Dictionary of primary ZG parameters.
Example:
``{"fildyn": "si.dyn", "output_format": "cube"}``
azg : dict, optional
Dictionary of auxiliary ZG parameters (e.g., projection settings, filters).
Example:
``{"projection": True, "energy_window": [0, 10]}``
name : str, optional
Base name for the generated ZG input file (default: ``"zg"``).
Workflow
--------
1. Sets up the working directory using ``set_work()``.
2. Applies ZG parameters via ``set_zg()``.
3. Writes the ZG input file using ``write_zg()``.
4. Stores file metadata in ``self.zg_file``.
5. Returns to the home directory with ``set_home()``.
Notes
-----
- Assumes required PW/PH data is available.
- Used primarily for zone-averaged or auxiliary post-processing.
Example
-------
>>> qe = PyQE({"type": "zg"}, system="si")
>>> zg_params = {"fildyn": "si.dyn", "output_format": "cube"}
>>> azg_params = {"projection": True, "energy_window": [0, 10]}
>>> qe.zg(zg=zg_params, azg=azg_params, name="zg_si")
>>> print(qe.zg_file)
zg_si
"""
if zg is None:
zg = {}
if azg is None:
azg = {}
self.set_work()
self.set_zg(zg,azg)
self.write_zg(name=name)
self.zg_file = name
self.set_home()
[docs]
def eps(self, inputpp=None, energy_grid=None, name='eps'):
if inputpp is None:
inputpp = {}
if energy_grid is None:
energy_grid = {}
"""
epsilon from qe preparation
"""
self.set_work()
self.set_eps(inputpp,energy_grid)
self.write_eps(name = name)
self.eps_file = name
self.set_home()
[docs]
def wannier(self, win=None, pw2wan=None, name='win'):
if win is None:
win = {}
if pw2wan is None:
pw2wan = {}
"""
Wannier calculation
"""
self.set_work()
self.set_wannier(win,pw2wan)
self.write_wann(name=name)
self.wannier_file = name
self.set_home()
[docs]
def EPW_pp(self):
"""
Run PP for EPW if everything fails
"""
self.set_work()
os.chdir('./epw')
run_pp(self.prefix)
self.set_home()
[docs]
def pw(self, control,
system,
electrons,
ions,
cell,
kpoints,
cell_parameters,
hubbard,
fcp,
solvents,
occupations,
atomic_velocities,
constraints,
atomic_forces,
rism,
cards):
"""
PW setup and setting with a changing state
"""
self.set_control(control)
self.set_system(system)
self.set_electrons(electrons)
self.set_ions(ions)
self.set_cell(cell)
self.set_cell_parameters(cell_parameters)
self.set_kpoints(kpoints)
#self.cards = cards
#for card in cards:
if (len(cards) != 0):
self.set_cards(cards)
[docs]
def bands(self, bands=None, name='bands'):
if bands is None:
bands = {}
self.set_work()
self.write_bands(bands=bands,name=name)
self.set_home()
[docs]
def nscf2supercond(self, nscf2supercond=None, name='nscf2supercond'):
if nscf2supercond is None:
nscf2supercond = {}
self.set_work()
self.write_nscf2supercond(nscf2supercond=nscf2supercond,name=name)
self.set_home()
[docs]
def matdyn(self, matdyn=None, kpoints=None, name='matdyn'):
"""
Prepare a **Matdyn (phonon dispersion)** calculation input.
The `matdyn()` method generates input files for the Quantum ESPRESSO
``matdyn.x`` utility, which computes phonon frequencies and eigenvectors
at arbitrary k-points using real-space interatomic force constants
obtained from ``q2r.x``.
Parameters
----------
matdyn : dict, optional
Dictionary containing user-defined parameters for the ``matdyn.x`` input.
Example:
``{"flfrc": "si.fc", "dos": False, "asr": "simple"}``
kpoints : dict, optional
Dictionary specifying phonon k-points or path for dispersion.
Example:
``{"path": "G X L G", "npoints": 100}``
name : str, optional
Base name for the generated Matdyn input file (default: ``"matdyn"``).
Workflow
--------
1. Initializes working directory via ``set_work()``.
2. Applies phonon and k-point parameters with ``set_matdyn()``.
3. Writes the input file using ``write_matdyn()``.
4. Returns to the home directory with ``set_home()``.
Notes
-----
- ``matdyn.x`` typically follows ``q2r.x`` in the phonon workflow.
- The output may include:
- ``.freq`` file containing phonon frequencies.
- ``.modes`` file containing phonon eigenvectors.
- Can be used to compute phonon dispersion relations or DOS.
Example
-------
>>> qe = PyQE({"type": "matdyn"}, system="si")
>>> matdyn_params = {"flfrc": "si.fc", "dos": False, "asr": "simple"}
>>> kpts = {["kpoints": ['0.5','0.50','0.50','20'],
... ['0.0','0.00','0.00','20'],
... ['0.5','0.25','0.75','20']],
... "kpoints_type" : " crystal_b "}
>>> qe.matdyn(matdyn=matdyn_params, kpoints=kpts, name="matdyn_si")
"""
if matdyn is None:
matdyn = {}
if kpoints is None:
kpoints = {}
self.set_work()
self.set_matdyn(matdyn,kpoints)
self.write_matdyn(matdyn=matdyn,name=name)
self.set_home()
[docs]
def dynmat(self, dynmat=None, name='dynmat'):
"""
Prepare a **Dynmat (dynamical matrix)** calculation input.
The `dynmat()` method generates input files for the Quantum ESPRESSO
``dynmat.x`` utility, which is used to compute phonon dynamical matrices
from the outputs of previous phonon calculations or force constants.
These matrices are then used to analyze phonon modes or as input for
further post-processing (e.g., ``matdyn.x``).
Parameters
----------
dynmat : dict, optional
Dictionary of user-defined parameters for ``dynmat.x``.
Example:
``{"fildyn": "si.dyn", "qpoints": [[0,0,0], [0.5,0,0]]}``
name : str, optional
Base name for the generated Dynmat input file (default: ``"dynmat"``).
Workflow
--------
1. Sets up the working directory using ``set_work()``.
2. Applies Dynmat parameters with ``set_dynmat()``.
3. Writes the Dynmat input file via ``write_dynmat()``.
4. Returns to the home directory using ``set_home()``.
Notes
-----
- Assumes prior phonon or force constant calculations are available.
- Output is typically used for phonon analysis or as input for ``matdyn.x``.
Example
-------
>>> qe = PyQE({"type": "dynmat"}, system="si")
>>> dynmat_params = {"fildyn": "si.dyn", "qpoints": [[0,0,0], [0.5,0,0]]}
>>> qe.dynmat(dynmat=dynmat_params, name="dynmat_si")
"""
if dynmat is None:
dynmat = {}
self.set_work()
self.set_dynmat(dynmat)
self.write_dynmat(name=name)
self.set_home()
[docs]
def phdos(self,phdos=None,name='phdos'):
if phdos is None:
phdos = {}
self.set_work()
self.write_phdos(phdos=phdos,name=name)
self.set_home()
[docs]
def dos(self,dos=None,name='dos'):
if dos is None:
dos = {}
self.set_work()
self.write_dos(dos=dos, name=name)
self.set_home()
[docs]
def pdos(self,pdos=None,name='pdos'):
if pdos is None:
pdos = {}
self.set_work()
self.write_pdos(pdos=pdos, name=name)
self.set_home()
[docs]
def save_PW(self):
"""
Saves PW state
"""
dict1={'system':self.pw_system,
'kpoints':self.pw_kpoints,
'control':self.pw_control,
'electrons':self.pw_electrons,
'atomic_species':self.pw_atomic_species,
'cell_parameters':self.pw_cell_parameters}
if (len(self.cards) !=0):
for card in self.cards.keys():
dict1.update({card:self.cards[card]})
self._save_json(self.nscf_file+'_save',dict1)#self.pw_system)
[docs]
def PH_utilities(self):
"""
Access **Phonon (PH)** post-processing utilities.
This method initializes and returns a `PHProperties` object for analyzing
phonon calculation outputs from Quantum ESPRESSO. It automatically sets
relevant attributes such as the atomic species, number of atoms, and the
dynamical matrix file location.
Returns
-------
PHProperties
An object that provides utilities to analyze phonon data such as:
- Phonon eigenvectors
- Frequencies
- Force constants
Notes
-----
- Assumes that a phonon calculation (`ph.x`) has been performed and the output
file exists at ``{prefix}/{ph_fold}/{ph_file}.out``.
- `dynmat_file` is automatically set to ``{prefix}/{ph_fold}/{prefix}.axsf`` for plotting.
Example
-------
>>> ph_util = qe.PH_utilities()
>>> ph_util.get_frequencies()
"""
prefix = self.prefix.replace('\'','')
self.ph_util = PHProperties(f'{prefix}/{self.ph_fold}/{self.ph_file}.out', version = self.QE_ver)
self.ph_util.atoms = self.pw_atomic_species['atomic_species']
self.ph_util.natoms = self.pw_system['nat']
self.ph_util.dynmat_file = f'{prefix}/{self.ph_fold}/{prefix}.axsf'
return self.ph_util
[docs]
def PW_utilities(self):
"""
Access **PW (Plane-Wave)** post-processing utilities.
This method initializes and returns a `PWProperties` object for analyzing
self-consistent field (SCF) calculation outputs from Quantum ESPRESSO.
Returns
-------
PWProperties
An object that provides utilities to analyze electronic structure data such as:
- Total energies
- Band structures
- Charge densities
Notes
-----
- Assumes that an SCF calculation has been performed and the output file exists
at ``{prefix}/{scf_fold}/{scf_file}.out``.
Example
-------
>>> pw_util = qe.PW_utilities()
>>> pw_util.get_total_energy()
"""
prefix = self.prefix.replace('\'','')
self.pw_util = PWProperties(f'{prefix}/{self.scf_fold}/{self.scf_file}.out')
return self.pw_util
[docs]
def EPW_utilities(self):
"""
Access **EPW (Electron-Phonon Wannier) properties** utility.
This method initializes and returns an `EPWProperties` object for analyzing
EPW calculation outputs from Quantum ESPRESSO. It allows users to access
polaron properties, electron-phonon matrix elements, and Wannier functions.
Returns
-------
EPWProperties
An object that provides utilities to analyze EPW data such as:
- Polaron displacement (`dtau`)
- Polaron wavefunction (`psir_plrn`)
- Electron-phonon matrix (`gkk`)
- Wannier cube files
- Transport properties
Notes
-----
- Assumes that an EPW calculation has been performed and the output file exists
at ``{self.epw_file}.out`` inside ``{self.system}/{self.epw_fold}/``.
- EPW parameters used in the calculation are passed via ``epw_params``.
Example
-------
>>> epw_util = qe.EPW_utilities()
>>> dtau = epw_util.dtau
>>> gkk_matrix = epw_util.gkk
"""
prefix=self.prefix.replace('\'','')
self.EPW_util = EPWProperties(file = f'{self.epw_file}.out',
system = f'{self.system}',
epw_fold = f'{self.epw_fold}',
epw_params = self.epw_params,
version = self.EPW_ver)
return self.EPW_util
[docs]
def save_EPW(self):
"""
Saves PW state
"""
dict1={'epw_data':self.epw_params,
'refresh': '.true.'}
self._save_json('epw_save',dict1)#self.pw_system)
[docs]
def save_PH(self):
"""
Saves PH state
"""
dict1={'ph_params':self.ph_params}
self._save_json(self.ph_file+'_save',dict1)#self.pw_system)
def _save_json(self,folder,data):
self.save=Save_state(' ',' ')
self.save.folder = folder
self.save.data = data
self.save.save_state()
def _read_master(self):
for file in glob.glob("*.json"):
#print(file)
self.save=Save_state(' ',' ')
dict1=self._read_json(file.split('.json')[0])
#print(dict1.keys())
if('control' in dict1.keys()):
# print(dict1['control']['calculation'])
if(dict1['control']['calculation'] == '\'nscf\''):
self.nscf_file = file.split('_')[0]
if(dict1['control']['calculation'] == '\'bands\''):
self.nscf_file = file.split('_')[0]
if('ph_params' in dict1.keys()):
self.ph_file = file.split('_')[0]
if('epw_params' in dict1.keys()):
self.epw_file = file.split('_')[0]
if('epwpy_params' in dict1.keys()):
self._get_folds(dict1)
self._get_files(dict1)
# print(self.nscf_file,self.ph_file)
[docs]
def read_PW(self):
"""
Reads PW state
"""
self.save=Save_state(' ',' ')
dict1=self._read_json(self.nscf_file+'_save')#self.pw_system)
self.pw_kpoints = dict1['kpoints']
[docs]
def read_PH(self):
"""
Reads PH state
"""
self.save=Save_state(' ',' ')
dict1=self._read_json(self.ph_file+'_save')#self.pw_system)
self.ph_params = dict1['ph_params']
[docs]
def read_EPW(self):
"""
Reads PH state
"""
self.save=Save_state(' ',' ')
dict1=self._read_json('epw_save')#self.pw_system)
self.epw_refresh = dict1['refresh']
def _read_json(self,folder):
self.save.folder = folder
return(self.save.read_state())
def _get_folds(self,dict1):
""" Gets the folders for various calculations """
self.state = dict1
self.scf_fold = dict1['scf_fold']
try:
self.nscf_fold = dict1['nscf_fold']
except KeyError:
pass
try:
self.ph_fold = dict1['ph_fold']
except KeyError:
pass
try:
self.epw_fold = dict1['epw_fold']
except KeyError:
pass
def _get_files(self,dict1):
""" Gets the folders for various calculations """
self.state = dict1
self.scf_file = dict1['scf_file']
try:
self.nscf_file = dict1['nscf_file']
except KeyError:
pass
try:
self.ph_file = dict1['ph_file']
except KeyError:
pass
try:
self.epw_file = dict1['epw_file']
except KeyError:
pass