EPWpy.EPWpy_run#

Classes

PyRun(procs, env, code)

The PyRun class provides an interface to execute external electronic structure and many-body codes such as Quantum ESPRESSO (QE) and BerkeleyGW (BGW) within the EPWpy framework.

class EPWpy.EPWpy_run.PyRun(procs, env, code)[source]#

Bases: object

The PyRun class provides an interface to execute external electronic structure and many-body codes such as Quantum ESPRESSO (QE) and BerkeleyGW (BGW) within the EPWpy framework.

It can be used both as part of a higher-level EPWpy workflow or as a standalone runner to execute specific computational tasks.

Parameters:
  • procs (int) -- Number of processors to be used for the execution.

  • env (str) -- Type of execution environment. Typical values include: - "mpirun" for MPI execution, - "srun" for SLURM-based parallel runs, - "ibrun" for TACC systems, - "local" for serial runs.

  • code (str) -- Path to the executable code or command to run. For example: "/path/to/pw.x" or "/path/to/bgw.x".

env#

Environment command used to launch the executable.

Type:

str

procs#

Number of processors assigned to the job.

Type:

int

code#

Full path or command to the simulation code being executed.

Type:

str

Examples

>>> from epwpy.run import PyRun
>>> run_qe = PyRun(procs=8, env="mpirun", code="/path/to/pw.x")
>>> run_qe
<PyRun object for QE execution>

Notes

  • The class is primarily used internally by EPWpy to manage code execution for different simulation stages (DFT, DFPT, EPW, etc.).

  • In HPC environments, env should correspond to the launcher provided by the scheduler (e.g., ibrun on TACC or srun on SLURM systems).

run_abinit(folder='./abinit', name='abi')[source]#

Runner for abinit

run_absorption(folder=' ', name=None, flavor='cplx')[source]#

Run the absorption calculation using the BGW (BerkeleyGW) absorption module.

This method executes the absorption.x program in BerkeleyGW, which computes the optical absorption spectrum of the system. It uses the previously computed electron-hole interaction kernel (from kernel.x) and quasiparticle energies (from sigma.x) to calculate excitonic effects and optical properties.

Parameters:
  • folder (str, optional) -- Directory where the absorption calculation will be executed. Default is ' ', indicating the current directory.

  • name (str, optional) -- Base name for the input/output files. If None, defaults to 'absorption'.

  • flavor (str, optional) -- Flavor of the executable, typically 'cplx' for complex or 'rpa' for real calculations. Default is 'cplx'.

  • Behavior

  • --------

:param - Sets the working directory (self.dir).: :param - Determines the input (self.filein) and output (self.fileout): filenames based on name. :param - Constructs the execution command (self.util) for absorption.flavor.x.: :param - Calls runner() to execute or script-write the absorption calculation.:

Returns:

The process handle returned by runner() if executed, or None if the job is written to a script instead.

Return type:

subprocess.Popen or None

Notes

  • Requires prior kernel and sigma calculations.

  • Produces optical absorption spectra, including excitonic effects.

  • Supports different flavors of the executable (complex or real).

Examples

>>> run = PyRun(procs=8, env='mpirun', code='/path/to/bgw')
>>> run.run_absorption(folder='./absorption', name='graphene_abs', flavor='cplx')
# Executes equivalent to:
# mpirun -np 8 /path/to/bgw/absorption.cplx.x -in graphene_abs.in > graphene_abs.out
run_bands(folder='./bs', name='bands')[source]#

Runner for bands

run_bs(folder='./bs', name='bs')[source]#

Run a band structure (BS) calculation using Quantum ESPRESSO.

This method executes the band structure step of a typical DFT workflow. It reads the converged charge density from a previous SCF calculation and computes eigenvalues along a predefined high-symmetry k-path.

Parameters:
  • folder (str, optional) -- Directory where the band structure calculation will be executed. Default is './bs'.

  • name (str, optional) -- Base name for the input and output files. Default is 'bs', corresponding to 'bs.in' and 'bs.out'.

  • Behavior

  • --------

:param - If self.proc_set is None: the number of k-point and thread pools using set_processors(). :param automatically determines: the number of k-point and thread pools using set_processors(). :param - Constructs the command line for Quantum ESPRESSO’s pw.x executable.: :param - Sets the working directory (self.dir): and output file (self.fileout). :param input file (self.filein): and output file (self.fileout). :param : and output file (self.fileout). :param - Executes the calculation via runner().:

Returns:

The process handle returned by runner() if executed, or None if the job is written to a script instead.

Return type:

subprocess.Popen or None

Notes

  • Requires a prior SCF calculation with a converged charge density file (*.save).

  • Typically used for generating electronic band structures along high-symmetry directions specified in the NSCF k-path.

  • Uses pw.x as the default executable unless overridden by self.proc_set.

Examples

>>> run = PyRun(procs=8, env='srun', code='/path/to/qe')
>>> run.run_bs(folder='./bs', name='graphene_bs')
# Executes equivalent to:
# srun -n 8 /path/to/qe/pw.x -nk 4 -nt 2 -in graphene_bs.in > graphene_bs.out
run_custom(folder='./custom', name='custom', util=None)[source]#

Run a custom executable within the EPWpy environment.

This method provides a flexible interface for executing any external code (e.g., a user-defined tool, a Quantum ESPRESSO utility, or a post-processing script) using the same EPWpy execution infrastructure. It automatically handles directory changes, input/output file naming, and command construction before delegating to runner().

Parameters:
  • folder (str, optional) -- Directory in which the calculation should be executed. Default is './custom'.

  • name (str, optional) -- Base name of the input/output files. Default is 'custom', corresponding to 'custom.in' and 'custom.out'.

  • util (str or None, optional) -- Name of the executable or utility to run. If None, attempts to use the value stored in self.util.

  • Behavior

  • --------

:param - Sets the working directory self.dir to folder.: :param - Defines self.filein and self.fileout based on name.: :param - Constructs self.util using either:

  • The provided util string, if specified.

  • The existing self.util if util is None.

:param - Invokes runner() to execute or script-write the command.:

Returns:

The process handle from runner(), or None if the command is written to a script instead of executed.

Return type:

subprocess.Popen or None

Raises:

ValueError -- If both util and self.util are None, meaning no executable is defined.

Examples

>>> run = PyRun(procs=4, env='srun', code='/path/to/qe')
>>> run.run_custom(folder='./pp', name='dos', util='dos.x')
# Executes equivalent to:
# srun -n 4 /path/to/qe/dos.x dos.in > dos.out

Notes

  • This method is useful for integrating third-party or auxiliary codes into automated EPWpy workflows.

  • If the environment variables (self.env, self.procs, etc.) are properly configured, run_custom() behaves identically to EPW-specific runners.

run_dos(folder='./nscf', name='dos')[source]#

Run DOS calculation using QE

run_dvscf_q2r(folder='./ph', name='dvscf_q2r')[source]#

Runner for dvscf_q2r

run_dynmat(folder='./ph', name='dynmat')[source]#

Run a dynmat calculation using Quantum ESPRESSO's dynmat.x module.

This method executes the dynmat step, which computes the full dynamical matrices for the system, typically used for phonon analysis, lattice dynamics, and further post-processing. It is usually performed after SCF and NSCF calculations, and optionally after phonon (ph.x) runs.

Parameters:
  • folder (str, optional) -- Directory where the dynmat calculation will be executed. Default is './ph'.

  • name (str, optional) -- Base name for the input/output files. Default is 'dynmat', corresponding to 'dynmat.in' and 'dynmat.out'.

  • Behavior

  • --------

:param - Sets the working directory (self.dir): and output file (self.fileout). :param input file (self.filein): and output file (self.fileout). :param : and output file (self.fileout). :param - Constructs the command string (self.util) as 'dynmat.x -in'.: :param - Calls runner() to execute or script-write the dynmat calculation.:

Returns:

The process handle returned by runner() if executed, or None if the job is written to a script instead.

Return type:

subprocess.Popen or None

Notes

  • Requires a prior SCF/NSCF calculation and, if applicable, phonon inputs.

  • Produces dynamical matrices used for phonon dispersion, EPW calculations, or other lattice dynamics analyses.

Examples

>>> run = PyRun(procs=4, env='mpirun', code='/path/to/qe')
>>> run.run_dynmat(folder='./ph', name='graphene_dynmat')
# Executes equivalent to:
# mpirun -np 4 /path/to/qe/dynmat.x -in graphene_dynmat.in > graphene_dynmat.out
run_eps(folder='./eps', name='eps')[source]#

Runner for epsilon

run_epsilon(folder=' ', name=None, flavor='cplx')[source]#

Run the epsilon calculation using the BGW (BerkeleyGW) epsilon module.

This method executes the dielectric function calculation (epsilon.x) in BerkeleyGW. It computes the screened Coulomb interaction, which is essential for many-body perturbation theory calculations such as GW and excitonic properties.

Parameters:
  • folder (str, optional) -- Directory where the epsilon calculation will be executed. Default is a blank string ' ', which indicates the current directory.

  • name (str, optional) -- Base name for the input/output files. If None, defaults to 'epsilon'.

  • flavor (str, optional) -- Flavor of the executable, typically 'cplx' for complex or 'rpa' for real calculations. Default is 'cplx'.

  • Behavior

  • --------

:param - Sets the working directory (self.dir).: :param - Determines the input (self.filein) and output (self.fileout): filenames based on name. :param - Constructs the execution command (self.util) for epsilon.flavor.x.: :param - Calls runner() to execute or script-write the epsilon calculation.:

Returns:

The process handle returned by runner() if executed, or None if the job is written to a script instead.

Return type:

subprocess.Popen or None

Notes

  • Requires prior GW-related input preparation (e.g., wavefunctions and screening files).

  • This runner supports different flavors of the executable (complex or real).

  • Typically used as the first step in a BerkeleyGW GW workflow.

Examples

>>> run = PyRun(procs=8, env='mpirun', code='/path/to/bgw')
>>> run.run_epsilon(folder='./epsilon', name='epsilon.inp', flavor='cplx')
# Executes equivalent to:
# mpirun -np 8 /path/to/bgw/epsilon.cplx.x -in epsilon.inp
run_epw1(folder='./epw', name='epw1')[source]#

Execute the first stage of an EPW (Electron–Phonon Wannier) calculation.

This method configures and runs the EPW step 1 task, which typically computes the electron–phonon coupling matrix elements on coarse Brillouin-zone grids. It sets up the appropriate working directory, input/output filenames, and executable command, then delegates execution to the runner() method.

Parameters:
  • folder (str, optional) -- Directory in which the EPW 1 calculation will be executed. Default is './epw'.

  • name (str, optional) -- Base name of the EPW input/output files. Default is 'epw1', corresponding to input file 'epw1.in' and output file 'epw1.out'.

  • Behavior

  • --------

:param - Sets the working directory attribute self.dir to folder.: :param - Defines self.filein and self.fileout using name.: :param - Constructs self.util: whether self.proc_set is defined:

  • If self.proc_set is None: builds epw.x -nk <nprocs> -in.

  • Otherwise: uses the explicit process configuration stored in self.proc_set.

Parameters:
  • epw.x (the command to run) --

    whether self.proc_set is defined:

    • If self.proc_set is None: builds epw.x -nk <nprocs> -in.

    • Otherwise: uses the explicit process configuration stored in self.proc_set.

  • on (depending) --

    whether self.proc_set is defined:

    • If self.proc_set is None: builds epw.x -nk <nprocs> -in.

    • Otherwise: uses the explicit process configuration stored in self.proc_set.

:param - Calls runner() to execute (or script-write) the command.:

Returns:

The process handle from runner(), or None if the run command is only written to a script.

Return type:

subprocess.Popen or None

Notes

  • Relies on pre-defined attributes such as self.procs, self.env, and self.code, which should be set before invoking this method.

  • Designed to integrate into EPWpy’s automated workflow; users normally call higher-level workflow functions rather than this routine directly.

Examples

>>> run = PyRun(procs=8, env='mpirun', code='/path/to/qe')
>>> run.proc_set = None
>>> run.run_epw1(folder='./epw', name='epw1')
# Executes equivalent to:
# mpirun -np 8 /path/to/qe/epw.x -nk 8 -in epw1.in > epw1.out
run_epw2(folder='./epw', name='epw2')[source]#

Execute the second stage of an EPW (Electron–Phonon Wannier) calculation.

This method launches EPW step 2, which typically performs fine-grid interpolation of electron–phonon coupling data, transport property calculations, or other post-processing tasks following EPW step 1. It configures the appropriate run parameters and calls the runner() method to execute the job.

Parameters:
  • folder (str, optional) -- Directory where the EPW step 2 calculation will be executed. Default is './epw'.

  • name (str, optional) -- Base name of the EPW input/output files. Default is 'epw2', corresponding to input file 'epw2.in' and output file 'epw2.out'.

  • Behavior

  • --------

:param - Sets the working directory attribute self.dir to folder.: :param - Defines self.filein and self.fileout using name.: :param - Constructs self.util: whether self.proc_set is defined:

  • If self.proc_set is None: builds epw.x -nk <nprocs> -in.

  • Otherwise: uses the explicit process setup stored in self.proc_set.

Parameters:
  • string (the EPW command) --

    whether self.proc_set is defined:

    • If self.proc_set is None: builds epw.x -nk <nprocs> -in.

    • Otherwise: uses the explicit process setup stored in self.proc_set.

  • on (depending) --

    whether self.proc_set is defined:

    • If self.proc_set is None: builds epw.x -nk <nprocs> -in.

    • Otherwise: uses the explicit process setup stored in self.proc_set.

:param - Calls runner() to execute or write the command script.:

Returns:

The process handle from runner(), or None if the run command is only written to a job script.

Return type:

subprocess.Popen or None

Notes

  • Assumes that self.procs, self.env, self.code, and self.proc_set (if applicable) have been initialized.

  • Typically called automatically by higher-level EPWpy workflow managers after successful completion of run_epw1().

Examples

>>> run = PyRun(procs=8, env='mpirun', code='/path/to/qe')
>>> run.proc_set = None
>>> run.run_epw2(folder='./epw', name='epw2')
# Executes equivalent to:
# mpirun -np 8 /path/to/qe/epw.x -nk 8 -in epw2.in > epw2.out
run_epw3(folder='./epw', name='epw3')[source]#

Runner for epw3

run_epw_outerbands(folder='./epw_outerbands', name='epw_outerbands')[source]#
run_fbw(folder='./fbw', name='fbw')[source]#
run_fbw_mu(folder='./fbw_mu', name='fbw_mu')[source]#
run_kernel(folder=' ', name=None, flavor='cplx')[source]#

Run the kernel calculation using the BGW (BerkeleyGW) kernel module.

This method executes the kernel.x program in BerkeleyGW, which computes the electron-hole interaction kernel required for solving the Bethe-Salpeter Equation (BSE). The kernel is essential for excitonic and optical properties calculations following GW and epsilon computations.

Parameters:
  • folder (str, optional) -- Directory where the kernel calculation will be executed. Default is ' ', indicating the current directory.

  • name (str, optional) -- Base name for the input/output files. If None, defaults to 'kernel'.

  • flavor (str, optional) -- Flavor of the executable, typically 'cplx' for complex or 'rpa' for real calculations. Default is 'cplx'.

  • Behavior

  • --------

:param - Sets the working directory (self.dir).: :param - Determines the input (self.filein) and output (self.fileout): filenames based on name. :param - Constructs the execution command (self.util) for kernel.flavor.x.: :param - Calls runner() to execute or script-write the kernel calculation.:

Returns:

The process handle returned by runner() if executed, or None if the job is written to a script instead.

Return type:

subprocess.Popen or None

Notes

  • Requires prior epsilon and sigma calculations.

  • Produces the electron-hole kernel necessary for BSE excitonic calculations.

  • Supports different flavors of the executable (complex or real).

Examples

>>> run = PyRun(procs=8, env='mpirun', code='/path/to/bgw')
>>> run.run_kernel(folder='./kernel', name='graphene_kernel', flavor='cplx')
# Executes equivalent to:
# mpirun -np 8 /path/to/bgw/kernel.cplx.x -in graphene_kernel.in > graphene_kernel.out
run_matdyn(folder='./ph', name='matdyn')[source]#

Run a matdyn calculation using Quantum ESPRESSO's matdyn.x module.

This method executes the matdyn step, which calculates phonon frequencies and eigenvectors along a specified path in the Brillouin zone, producing phonon band structures and related properties. It is typically used after q2r to obtain interpolated phonon bands.

Parameters:
  • folder (str, optional) -- Directory where the matdyn calculation will be executed. Default is './ph'.

  • name (str, optional) -- Base name for the input/output files. Default is 'matdyn', corresponding to 'matdyn.in' and 'matdyn.out'.

  • Behavior

  • --------

:param - Sets the working directory (self.dir): and output file (self.fileout). :param input file (self.filein): and output file (self.fileout). :param : and output file (self.fileout). :param - Constructs the command string (self.util) as 'matdyn.x -in'.: :param - Calls runner() to execute or script-write the matdyn calculation.:

Returns:

The process handle returned by runner() if executed, or None if the job is written to a script instead.

Return type:

subprocess.Popen or None

Notes

  • Requires prior q2r calculations to generate real-space force constants.

  • Produces interpolated phonon band structures, which can be used for plotting or electron-phonon calculations in EPW.

Examples

>>> run = PyRun(procs=4, env='mpirun', code='/path/to/qe')
>>> run.run_matdyn(folder='./ph', name='graphene_matdyn')
# Executes equivalent to:
# mpirun -np 4 /path/to/qe/matdyn.x -in graphene_matdyn.in > graphene_matdyn.out
run_nesting(folder='./nesting', name='nesting')[source]#
run_nscf(folder='./nscf', name='nscf')[source]#

Run a non-self-consistent field (NSCF) calculation using Quantum ESPRESSO.

This method executes the NSCF stage, which computes the electronic structure on a denser k-point grid using the converged charge density from the SCF run. It constructs the appropriate pw.x command and runs it via runner().

Parameters:
  • folder (str, optional) -- Directory where the NSCF calculation will be executed. Default is './nscf'.

  • name (str, optional) -- Base name for the input and output files. Default is 'nscf', corresponding to 'nscf.in' and 'nscf.out'.

  • Behavior

  • --------

  • using (- Determines the number of thread and k-point pools) -- set_processors() if self.proc_set is None.

  • executable. (- Constructs the NSCF command using Quantum ESPRESSO’s pw.x)

:param - Sets the execution directory (self.dir): and output (self.fileout) file paths. :param input (self.filein): and output (self.fileout) file paths. :param : and output (self.fileout) file paths. :param - Calls runner() to execute or write the NSCF command.:

Returns:

The process handle returned by runner() if executed, or None if the job is written to a script instead.

Return type:

subprocess.Popen or None

Notes

  • The NSCF run should follow a completed SCF calculation in the same system folder.

  • This step is essential for generating the wavefunctions required for phonon or EPW calculations.

  • Uses pw.x by default, unless overridden via self.proc_set.

Examples

>>> run = PyRun(procs=8, env='srun', code='/path/to/qe')
>>> run.run_nscf(folder='./nscf', name='graphene_nscf')
# Executes equivalent to:
# srun -n 8 /path/to/qe/pw.x -nk 4 -nt 2 -in graphene_nscf.in > graphene_nscf.out
run_nscf2supercond(folder='./nscf', name='nscf2supercond')[source]#

Runner for nscf2supercond

run_nscf_tetra(folder='./nscf_tetra', name='nscf')[source]#
run_pdos(folder='./nscf', name='pdos')[source]#

Run PDOS calculation

run_ph(folder='./ph', name='ph')[source]#

Run a phonon calculation using Quantum ESPRESSO's ph.x module.

This method executes the phonon calculation stage, which computes phonon frequencies, dynamical matrices, and related properties required for electron-phonon coupling analysis. It sets up the appropriate working directory, input/output files, and constructs the parallel execution command before calling runner().

Parameters:
  • folder (str, optional) -- Directory where the phonon calculation will be executed. Default is './ph'.

  • name (str, optional) -- Base name of the input/output files. Default is 'ph', corresponding to 'ph.in' and 'ph.out'.

  • Behavior

  • --------

:param - If self.proc_set is None: thread and k-point pools using set_processors(). :param determines the number of: thread and k-point pools using set_processors(). :param - Constructs the execution command using ph.x: -pd .true. flag for parallel diagonalization. :param including the: -pd .true. flag for parallel diagonalization. :param - Sets the working directory (self.dir): (self.filein), and output file (self.fileout). :param input file: (self.filein), and output file (self.fileout). :param - Calls runner() to execute or script-write the phonon job.:

Returns:

The process handle returned by runner() if executed, or None if the job is written to a script instead.

Return type:

subprocess.Popen or None

Notes

  • Requires a prior SCF/NSCF calculation for the converged charge density.

  • The -pd .true. flag enables parallel phonon calculations.

  • This step is essential for subsequent EPW calculations involving electron-phonon interactions.

Examples

>>> run = PyRun(procs=8, env='srun', code='/path/to/qe')
>>> run.run_ph(folder='./ph', name='graphene_ph')
# Executes equivalent to:
# srun -n 8 /path/to/qe/ph.x -pd .true. -nk 4 -nt 2 -in graphene_ph.in > graphene_ph.out
run_phdos(folder='./ph', name='phdos')[source]#
run_phselfen(folder='./phselfen', name='phselfen')[source]#
run_postahc(folder='./ph', name='postahc')[source]#

Runner for postahc

run_pp(folder='./pp', name='pp')[source]#

Run a post-processing (PP) calculation using Quantum ESPRESSO's pp.x module.

This method executes a generic post-processing step, which can be used to compute charge densities, potential maps, densities of states, or other derived quantities from previously calculated wavefunctions or densities. It sets up the appropriate working directory, input/output files, and constructs the execution command before delegating to runner().

Parameters:
  • folder (str, optional) -- Directory where the post-processing calculation will be executed. Default is './pp'.

  • name (str, optional) -- Base name of the input/output files. Default is 'pp', corresponding to 'pp.in' and 'pp.out'.

  • Behavior

  • --------

:param - If self.proc_set is None: basic parallel pp.x command. :param builds a default serial or: basic parallel pp.x command. :param - Otherwise: :param constructs the command using the specified processor setup.: :param - Sets the working directory (self.dir): (self.filein), and output file (self.fileout). :param input file: (self.filein), and output file (self.fileout). :param - Calls runner() to execute or script-write the post-processing job.:

Returns:

The process handle returned by runner() if executed, or None if the job is written to a script instead.

Return type:

subprocess.Popen or None

Notes

  • Requires a prior SCF or NSCF calculation providing the necessary wavefunctions or densities.

  • The pp.x module is flexible; this runner assumes the user has prepared the input file (pp.in) appropriately for the desired analysis.

Examples

>>> run = PyRun(procs=4, env='mpirun', code='/path/to/qe')
>>> run.run_pp(folder='./pp', name='graphene_pp')
# Executes equivalent to:
# mpirun -np 4 /path/to/qe/pp.x -in graphene_pp.in > graphene_pp.out
run_pw2bgw(folder=' ', name=None)[source]#

Runner for pw2BGW

run_q2r(folder='./ph', name='q2r')[source]#

Run a q2r calculation using Quantum ESPRESSO's q2r.x module.

This method executes the q2r step, which converts dynamical matrices (computed by ph.x) from reciprocal space to real-space interatomic force constants. The output is typically used for phonon interpolation or subsequent electron-phonon calculations.

Parameters:
  • folder (str, optional) -- Directory where the q2r calculation will be executed. Default is './ph'.

  • name (str, optional) -- Base name for the input/output files. Default is 'q2r', corresponding to 'q2r.in' and 'q2r.out'.

  • Behavior

  • --------

:param - Sets the working directory (self.dir): and output file (self.fileout). :param input file (self.filein): and output file (self.fileout). :param : and output file (self.fileout). :param - Constructs the command string (self.util) as 'q2r.x -in '.: :param - Delegates execution or script-writing to runner().:

Returns:

The process handle returned by runner() if executed, or None if the job is written to a script instead.

Return type:

subprocess.Popen or None

Notes

  • Requires prior phonon calculations using ph.x to generate the necessary dynamical matrices.

  • The generated real-space force constants are essential for EPW or post-processing phonon calculations.

Examples

>>> run = PyRun(procs=4, env='mpirun', code='/path/to/qe')
>>> run.run_q2r(folder='./ph', name='graphene_q2r')
# Executes equivalent to:
# mpirun -np 4 /path/to/qe/q2r.x -in graphene_q2r.in > graphene_q2r.out
run_scf(folder='./scf', name='scf')[source]#

Run a self-consistent field (SCF) calculation using Quantum ESPRESSO.

This method constructs and executes the SCF calculation command, automatically setting the parallelization options (-nk and -nt) based on user input or preconfigured processor settings.

Parameters:
  • folder (str, optional) -- Directory where the SCF calculation will be executed. Default is './scf'.

  • name (str, optional) -- Base name for the input and output files. Default is 'scf', corresponding to 'scf.in' and 'scf.out'.

  • Behavior

  • --------

  • from (- Determines the number of k-point and thread pools) -- set_processors() if self.proc_set is None.

  • pw.x. (- Constructs the Quantum ESPRESSO SCF command using)

:param - Sets the execution directory (self.dir): and output (self.fileout) filenames. :param input (self.filein): and output (self.fileout) filenames. :param : and output (self.fileout) filenames. :param - Calls runner() to execute the command or write it to a script.:

Returns:

The process handle from runner(), or None if the job is written to a script instead of being executed directly.

Return type:

subprocess.Popen or None

Notes

  • Requires a properly configured Quantum ESPRESSO environment.

  • This is typically the first stage of an EPW workflow, preceding NSCF, PH, and EPW runs.

  • Uses pw.x as the default executable unless overridden by self.proc_set.

Examples

>>> run = PyRun(procs=8, env='srun', code='/path/to/qe')
>>> run.run_scf(folder='./scf', name='graphene_scf')
# Executes equivalent to:
# srun -n 8 /path/to/qe/pw.x -nk 4 -nt 2 -in graphene_scf.in > graphene_scf.out
run_sig2wan(folder=' ', name=None)[source]#

Run the sigma-to-Wannier conversion using BGW's sig2wan module.

This method executes the sig2wan.x program in BerkeleyGW, which converts the computed electron self-energy (sigma) into a format compatible with Wannier interpolation. This step is typically used to interface GW corrections with Wannier90 or EPW workflows.

Parameters:
  • folder (str, optional) -- Directory where the sig2wan calculation will be executed. Default is ' ', indicating the current directory.

  • name (str, optional) -- Base name for the input/output files. If None, defaults to 'sig2wan'.

  • Behavior

  • --------

:param - Sets the working directory (self.dir).: :param - Determines the input (self.filein) and output (self.fileout): filenames based on name. :param - Constructs the execution command (self.util) for sig2wan.x.: :param - Calls runner() to execute or script-write the sig2wan calculation.:

Returns:

The process handle returned by runner() if executed, or None if the job is written to a script instead.

Return type:

subprocess.Popen or None

Notes

  • Requires prior sigma calculations from sigma.x.

  • The output can be used in Wannier90 or EPW workflows for band structure interpolation including GW corrections.

Examples

>>> run = PyRun(procs=8, env='mpirun', code='/path/to/bgw')
>>> run.run_sig2wan(folder='./sig2wan', name='graphene_sig2wan')
# Executes equivalent to:
# mpirun -np 8 /path/to/bgw/sig2wan.x -in graphene_sig2wan.in > graphene_sig2wan.out
run_sigma(folder=' ', name=None, flavor='cplx')[source]#

Run the sigma calculation using the BGW (BerkeleyGW) sigma module.

This method executes the self-energy calculation (sigma.x) in BerkeleyGW. It computes the electron self-energy based on a prior dielectric screening calculation (epsilon) and wavefunctions, which is required for GW quasiparticle energy corrections.

Parameters:
  • folder (str, optional) -- Directory where the sigma calculation will be executed. Default is a blank string ' ', indicating the current directory.

  • name (str, optional) -- Base name for the input/output files. If None, defaults to 'sigma'.

  • flavor (str, optional) -- Flavor of the executable, typically 'cplx' for complex or 'rpa' for real calculations. Default is 'cplx'.

  • Behavior

  • --------

:param - Sets the working directory (self.dir).: :param - Determines the input (self.filein) and output (self.fileout): filenames based on name. :param - Constructs the execution command (self.util) for sigma.flavor.x.: :param - Calls runner() to execute or script-write the sigma calculation.:

Returns:

The process handle returned by runner() if executed, or None if the job is written to a script instead.

Return type:

subprocess.Popen or None

Notes

  • Requires prior epsilon and wavefunction calculations.

  • This calculation provides the GW quasiparticle self-energy.

  • Supports different flavors of the executable (complex or real).

Examples

>>> run = PyRun(procs=8, env='mpirun', code='/path/to/bgw')
>>> run.run_sigma(folder='./sigma', name='sigma.inp', flavor='cplx')
# Executes equivalent to:
# mpirun -np 8 /path/to/bgw/sigma.cplx.x -in sigma.inp
run_vasp(folder='./vasp', name=' ', flavor='std')[source]#

Runner for abinit

run_wannier(folder='./wannier', name='win')[source]#

Runner for Wannier

run_zg(folder='./zg', name='zg')[source]#

Run a ZG (Zeta-Gamma) calculation using the ZG.x executable.

This method executes the ZG step of the workflow, which is typically used for specialized post-processing or analysis tasks defined by the ZG.x module. It sets up the working directory, input/output files, and constructs the execution command before delegating to runner().

Parameters:
  • folder (str, optional) -- Directory where the ZG calculation will be executed. Default is './zg'.

  • name (str, optional) -- Base name for the input and output files. Default is 'zg', corresponding to 'zg.in' and 'zg.out'.

  • Behavior

  • --------

:param - Sets self.dir: :param self.filein: :param and self.fileout.: :param - Constructs the command string self.util as 'ZG.x -in '.: :param - Calls runner() to execute or script-write the ZG calculation.:

Returns:

The process handle returned by runner() if executed, or None if the job is written to a script instead.

Return type:

subprocess.Popen or None

Notes

  • Requires a properly prepared input file (zg.in) according to the ZG.x module specifications.

  • This runner does not currently support parallelization flags, but could be extended to include them if needed.

Examples

>>> run = PyRun(procs=4, env='mpirun', code='/path/to/zg')
>>> run.run_zg(folder='./zg', name='graphene_zg')
# Executes equivalent to:
# mpirun -np 4 /path/to/zg/ZG.x -in graphene_zg.in > graphene_zg.out
runner(changedir=True)[source]#

Main execution routine for running external simulation executables.

This method constructs and executes the appropriate command-line instruction to launch a target code (e.g., Quantum ESPRESSO, EPW, or BGW) based on the user-specified environment (env), number of processors, and runtime configuration. It supports parallel environments such as mpirun, srun, and ibrun, as well as serial execution.

If a run script writer is active (self.writer.script_params['write_script'] == True), the method writes the command to a script instead of executing it.

Parameters:
  • changedir (bool, optional) -- Whether to change to the working directory (self.dir) before execution. Default is True.

  • Behavior

  • --------

  • provided (- If no env is) -- the environment type (mpirun, srun, ibrun, etc.).

  • on (constructs a parallel execution command based) -- the environment type (mpirun, srun, ibrun, etc.).

  • provided

  • run. (assumes a direct serial)

  • self.fileout. (- Redirects standard output to)

  • True (- When self.writer.script_params['write_script'] is) -- writes the run command to a job script instead of executing.

:param : writes the run command to a job script instead of executing. :param - When self.serial is True: :param the process waits for completion (p1.wait()).:

Returns:

A Popen process object representing the executed command, unless the command is only written to a script (in which case no process is launched).

Return type:

subprocess.Popen

Raises:
  • FileNotFoundError -- If changedir=True but self.dir does not exist.

  • AttributeError -- If required attributes such as self.code, self.util, or self.filein are missing.

Examples

>>> run = PyRun(procs=8, env='mpirun', code='/path/to/qe')
>>> run.util = 'pw.x'
>>> run.filein = 'scf.in'
>>> run.fileout = 'scf.out'
>>> process = run.runner()
>>> process.wait()

Notes

  • The constructed command typically looks like:

    mpirun -np 8 /path/to/qe/pw.x scf.in > scf.out
    
  • The method can be adapted for other HPC launchers by extending the env conditionals.

  • For debugging, verbosity levels control printed output:
    • verbosity > 1 prints the command before execution

    • verbosity > 2 appends real-time log following (tail -f)

set_processors(procs)[source]#