EPWpy.QE.EPWpy_QE#

QE class for interfacing with quantum espresso

Classes

PyQE(type_input[, system, code, env])

High-level interface for building Quantum ESPRESSO and EPW input files.

class EPWpy.QE.EPWpy_QE.PyQE(type_input, system='si', code='.', env='mpirun')[source]#

Bases: 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).

system#

Name of the material or simulation system (used to create a working directory).

Type:

str

code#

Path to the compiled Quantum ESPRESSO executables (default: current directory).

Type:

str

home#

The current working directory when the class is initialized.

Type:

str

env#

Parallel execution environment (e.g., "mpirun", "srun", or none).

Type:

str

EPW_pp()[source]#

Run PP for EPW if everything fails

EPW_utilities()[source]#

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:

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

Return type:

EPWProperties

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
PH_utilities()[source]#

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:

An object that provides utilities to analyze phonon data such as: - Phonon eigenvectors - Frequencies - Force constants

Return type:

PHProperties

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()
PW_utilities()[source]#

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:

An object that provides utilities to analyze electronic structure data such as: - Total energies - Band structures - Charge densities

Return type:

PWProperties

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()
QE_pp(inputpp=None, plot=None, name='pp')[source]#

Post-processing calculation of QE

bands(bands=None, name='bands')[source]#
custom(parameters=None, name='custom')[source]#

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")
dos(dos=None, name='dos')[source]#
dvscf_q2r(dvscf_q2r=None, name='dvscf_q2r')[source]#
dynmat(dynmat=None, name='dynmat')[source]#

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

  • --------

  • set_work(). (1. Sets up the working directory using)

  • set_dynmat(). (2. Applies Dynmat parameters with)

  • write_dynmat(). (3. Writes the Dynmat input file via)

  • set_home(). (4. Returns to the home directory using)

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")
eps(inputpp=None, energy_grid=None, name='eps')[source]#
epw(epwin=None, name='epw')[source]#

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

  • --------

  • set_work(). (1. Sets up the working directory using)

  • available. (2. Reads stored Quantum ESPRESSO (PW) and Phonon (PH) data if)

  • fresh (3. Determines whether to restart from a previous EPW run or start) --

    • If name is "epw1", "fbw", or "fbw_mu", restart is disabled.

    • Otherwise, EPW restart is attempted by loading saved EPW data.

  • set_epw(). (4. Applies EPW parameters using)

  • write_epw(). (5. Writes the EPW input file using)

  • directory. (6. Saves configuration with save_EPW() and returns to home)

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
matdyn(matdyn=None, kpoints=None, name='matdyn')[source]#

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

  • --------

  • set_work(). (1. Initializes working directory via)

  • set_matdyn(). (2. Applies phonon and k-point parameters with)

  • write_matdyn(). (3. Writes the input file using)

  • set_home(). (4. Returns to the home directory with)

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")
nscf(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')[source]#

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

  • --------

  • set_work(). (1. Sets up the working directory using)

  • pw_control. (2. Forces the calculation type to 'nscf' in)

  • pw(). (3. Applies all input dictionaries using)

  • write_scf(). (4. Writes the NSCF input file via)

  • save_PW(). (5. Saves current PW configuration using)

  • set_home(). (6. Returns to the home directory using)

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
nscf2supercond(nscf2supercond=None, name='nscf2supercond')[source]#
pdos(pdos=None, name='pdos')[source]#
ph(phonons=None, qpoints=None, name='ph')[source]#

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

  • --------

  • set_work(). (1. Prepares the working directory via)

  • set_phonons(). (2. Applies phonon and q-point parameters using)

  • write_ph(). (3. Writes the phonon input file using)

  • save_PH(). (4. Saves a phonon configuration snapshot with)

  • set_home(). (5. Returns to the home directory via)

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
phdos(phdos=None, name='phdos')[source]#
postahc(postahc=None, name='postahc')[source]#
pw(control, system, electrons, ions, cell, kpoints, cell_parameters, hubbard, fcp, solvents, occupations, atomic_velocities, constraints, atomic_forces, rism, cards)[source]#

PW setup and setting with a changing state

q2r(q2r=None, name='q2r')[source]#

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

  • --------

  • set_work(). (1. Initializes a working directory with)

  • set_q2r(). (2. Applies provided q2r parameters using)

  • write_q2r(). (3. Writes the Q2R input file via)

  • (self.q2r_file). (4. Stores file metadata)

  • set_home(). (5. Returns to the home directory using)

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
read_EPW()[source]#

Reads PH state

read_PH()[source]#

Reads PH state

read_PW()[source]#

Reads PW state

reset()[source]#

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()
save_EPW()[source]#

Saves PW state

save_PH()[source]#

Saves PH state

save_PW()[source]#

Saves PW state

scf(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')[source]#

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 (dict, optional) -- Dictionaries corresponding to Quantum ESPRESSO namelists such as &CONTROL, &SYSTEM, &ELECTRONS, etc.

  • system (dict, optional) -- Dictionaries corresponding to Quantum ESPRESSO namelists such as &CONTROL, &SYSTEM, &ELECTRONS, etc.

  • electrons (dict, optional) -- Dictionaries corresponding to Quantum ESPRESSO namelists such as &CONTROL, &SYSTEM, &ELECTRONS, etc.

  • ions (dict, optional) -- Dictionaries corresponding to Quantum ESPRESSO namelists such as &CONTROL, &SYSTEM, &ELECTRONS, etc.

  • cell (dict, optional) -- Dictionaries corresponding to Quantum ESPRESSO namelists such as &CONTROL, &SYSTEM, &ELECTRONS, etc.

  • kpoints (dict, optional) -- Dictionaries corresponding to Quantum ESPRESSO namelists such as &CONTROL, &SYSTEM, &ELECTRONS, etc.

  • cell_parameters (dict, optional) -- Dictionaries corresponding to Quantum ESPRESSO namelists such as &CONTROL, &SYSTEM, &ELECTRONS, etc.

  • hubbard

  • fcp

  • solvents

  • occupations

  • atomic_velocities

  • constraints

:param : :type atomic_forces: :param atomic_forces: Optional dictionaries specifying less-common QE sections (e.g.,

HUBBARD, FCP, SOLVENTS). Each is optional and defaults to {}.

Parameters:
  • rism (dict, optional) -- Optional dictionaries specifying less-common QE sections (e.g., HUBBARD, FCP, SOLVENTS). Each is optional and defaults to {}.

  • 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

  • --------

  • set_work(). (1. Creates or moves to the working directory using)

  • variables. (2. Calls pw() to populate QE input)

  • write_scf(). (3. Writes the SCF input file using)

  • self.scf_file. (4. Stores the SCF filename in)

  • set_home(). (5. Returns to the home directory using)

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]}
... )
wannier(win=None, pw2wan=None, name='win')[source]#
zg(zg=None, azg=None, name='zg')[source]#

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

  • --------

  • set_work(). (1. Sets up the working directory using)

  • set_zg(). (2. Applies ZG parameters via)

  • write_zg(). (3. Writes the ZG input file using)

  • self.zg_file. (4. Stores file metadata in)

  • set_home(). (5. Returns to the home directory with)

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