Skip to content

pwxml

Classes for reading/manipulating PWscf xml files.

PWxml(filename, ionic_step_skip=1, ionic_step_offset=0, parse_dos=False, fildos=None, parse_pdos=False, filpdos=None, parse_projected_eigen=False, filproj=None, occu_tol=1e-08, separate_spins=False, **_)

Bases: Vasprun

Parser for PWscf xml files. Almost all Vasprun features are implemented. This class exposes a public API that is practically identical to the Vasprun class, with values always returned in the same coordinate and unit systems used by VASP. Can be used interchangeably with Vasprun for most purposes with no changes to the calling code. This class will pull data from the output of projwfc.x and dos.x when necessary, attempting to guess the filenames if they are not provided.

Some attributes are not particularly meaningful for QE, and a few are not implemented.

Missing Vasprun features:

  • Projected magnetization. Not sure we can get this from QE (maybe pp.x?). bands.x supports computing \(\langle \sigma_\alpha \rangle\), but not projected.
  • Phononic properties: Need to write parser for ph.x output for dynamical matrix, interatomic force constants, phonon frequencies, and eigenvectors.
  • Dielectric properties: None of these are implemented. Some can be obtained from ph.x or QE's old epsilon.x. Some are RPA and thus not available in QE.

Not all attributes are listed below since many are inherited directly from Vasprun.

ATTRIBUTE DESCRIPTION
kpoints_frac

kpoints in fractional coordinates, (nk, 3) array

TYPE: ndarray

kpoints_cart

kpoints in cartesian coords (1/angst), (nk, 3) array

TYPE: ndarray

actual_kpoints

Same as kpoints_frac, maintained for compatibility with Vasprun.

TYPE: ndarray

actual_kpoints_weights

List of kpoint weights, normalized.

TYPE: ndarray

atomic_symbols

List of atomic symbols, e.g., ["Li", "Fe", "Fe", "P", "P", "P"].

TYPE: list[str]

efermi

Fermi energy, eV.

TYPE: float

vbm

Valence band maximum, as computed by QE, in eV.

TYPE: float

cbm

Conduction band minimum, as computed by QE, in eV.

TYPE: float

eigenvalues

Final eigenvalues (in eV) as a dict with keys Spin.up/Spin.down, and the values are numpy arrays with dimensions (nk, nb, 2), where nk is the number of kpoints, nb is the number of bands, and the last axis is the eigenvalue and occupation. Identical representation to Vasprun.

TYPE: dict[Spin, ndarray]

projected_eigenvalues

Final projected eigenvalues as a dict of {spin: nd-array}. To access a particular value, you need to index as [spin][kpoint index][band index][atom index][orbital_index]. This representation is identical to Vasprun.xml.

TYPE: dict[Spin, ndarray]

tdos

Total dos from fildos or filpdos. filpdos overrides fildos.

TYPE: Dos

idos

Integrated dos from fildos.

TYPE: Dos

pdos

List of Dos objects, indexed as [atom][orbital][spin].

TYPE: [List[dict[Orbital, Dict[spin, ndarray]]]]

ionic_steps

All ionic steps, list of

{ "structure": structure at end of run, "total_energy": {"ehart", "etxc", ...} "stress": pressure in eV/A3, "forces": forces in eV/A }.

TYPE: list[dict]

nionic_steps

The total number of ionic steps. This number is always equal to the total number of steps in the actual run even if ionic_step_skip is used. nionic_steps here has a slightly different meaning from the Vasprun class. VASP will first do an SCF calculation with the input structure, then perform geometry optimization until you hit EDIFFG or NSW, then it's done. QE does the same thing for relax, but it will also do a final SCF calculation with the optimized structure and a new basis set for vc-relax. In reality, converged QE vc-relax calculations take nionic_steps-1 to converge.

TYPE: int

projected_magnetisation

Not implemented

TYPE: dict[Spin, ndarray]

other_dielectric

Not implemented

TYPE: dict

force_constants

Not implemented.

TYPE: ndarray

normalmode_eigenvals

Not implemented.

TYPE: ndarray

normalmode_eigenvecs

Not implemented.

TYPE: ndarray

md_data

Molecular dynamics data, not implemented for QE.

TYPE: dict

parameters

Parameters of the PWscf run from the XML. Naturally very different from Vasprun.

TYPE: dict[str, Any]

incar

Empty Incar object, for compatibility.

TYPE: Incar

kpoints

Empty Kpoints object, for compatibility.

TYPE: Kpoints

pseudo_filenames

List of pseudopotential filenames, e.g., ["Li.pbe-spn-kjpaw_psl.0.1.UPF", "Fe.pbe-n-kjpaw_psl.0.2.1.UPF", ...].

TYPE: list[str]

potcar_symbols

Maintained for compatibility with Vasprun, has Psuedo + element names, e.g., ['Si.pbesol-n-rrkjus_psl.1.0.0.UPF Si'].

TYPE: list[str]

Author: Omar A. Ashour

PARAMETER DESCRIPTION
filename

Path to XML file to parse

TYPE: str | PathLike

ionic_step_skip

If ionic_step_skip is a number > 1, only every ionic_step_skip ionic steps will be read for structure and energies. Unlike Vasprun, the final energy will always be the total energy of the scf calculation performed after ionic convergence. This is not very useful since PWscf xml files aren't as huge as Vasprun files. Mainly kept for consistency with the Vasprun class.

TYPE: int DEFAULT: 1

ionic_step_offset

Used together with ionic_step_skip. If set, the first ionic step read will be offset by the amount of ionic_step_offset. For example, if you want to start reading every 10th structure but only from the 3rd structure onwards, set ionic_step_skip to 10 and ionic_step_offset to 3. Main use case is when doing statistical structure analysis with extremely long time scale PWscf calculations of varying numbers of steps, and kept for consistency with the Vasprun class.

TYPE: int DEFAULT: 0

parse_dos

Whether to parse the dos. If True, PWxml will use fildos if provided, or attempt to guess the filename.

TYPE: bool DEFAULT: False

fildos

If provided, forces parse_dos to be True and uses the provided string as the path to the dos file. This is the same as in dos.x input, and shouldn't be a full path to a file. For example, filpdos="path/to/filpdos" will look for "path/to/filpdos.pdos_*"

TYPE: str DEFAULT: None

filpdos

If provided, forces parse_dos to be True and uses the provided string as 'filproj', same as in projwfc.x input. It shouldn't include the rest of the filename. For example, filpdos="path/to/filpdos" will look for "path/to/filpdos.pdos_*"

TYPE: str DEFAULT: None

parse_projected_eigen

Whether to parse the projected eigenvalues and (magnetisation, not implemented). Defaults to False. If True, PWxml will look for a "filproj" from projwfc.x and parse it. It will look for files with the same name as the XML (or same QE prefix) but with a .projwfc_up extension, or will use the filproj argument

TYPE: bool DEFAULT: False

filproj

If provided, forces parse_projected_eigen to be True and uses the provided string as the filepath to the .projwfc_up file. Note that this is the same as in projwfc.x input, so it shouldn't include the .projwfc_up/down extension. It can also include a directory, e.g., "path/to/filproj" will look for "path/to/filproj.projwfc_up"

TYPE: str DEFAULT: None

occu_tol

Sets the minimum tol for the determination of the vbm and cbm. Usually the default of 1e-8 works well enough, but there may be pathological cases. Note that, unlike VASP, QE actually reports the VBM and CBM (accessible via the vbm and cbm attributes), so this is only used to recompute them and check against the reported values.

TYPE: float DEFAULT: 1e-08

separate_spins

Whether the band gap, CBM, and VBM should be reported for each individual spin channel. Defaults to False, which computes the eigenvalue band properties independent of the spin orientation. If True, the calculation must be spin-polarized.

TYPE: bool DEFAULT: False

Source code in pymatgen/io/espresso/outputs/pwxml.py
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
def __init__(
    self,
    filename: str | os.PathLike,
    ionic_step_skip: int = 1,
    ionic_step_offset: int = 0,
    parse_dos: bool = False,
    fildos: str | os.PathLike | None = None,
    parse_pdos: bool = False,
    filpdos: str | os.PathLike | None = None,
    parse_projected_eigen: bool = False,
    filproj: str | os.PathLike | None = None,
    occu_tol: float = 1e-8,
    separate_spins: bool = False,
    **_,  # Ignored arguments for compatibility with Vasprun
):
    """
    Constructor for the PWxml class.

    Args:
        filename (str | os.PathLike): Path to XML file to parse
        ionic_step_skip (int): If ionic_step_skip is a number > 1,
            only every ionic_step_skip ionic steps will be read for
            structure and energies. Unlike Vasprun, the final energy
            will always be the total energy of the scf calculation
            performed after ionic convergence. This is not very useful
            since PWscf xml files aren't as huge as Vasprun files.
            Mainly kept for consistency with the Vasprun class.
        ionic_step_offset (int): Used together with ionic_step_skip. If set,
            the first ionic step read will be offset by the amount of
            ionic_step_offset. For example, if you want to start reading
            every 10th structure but only from the 3rd structure onwards,
            set ionic_step_skip to 10 and ionic_step_offset to 3. Main use
            case is when doing statistical structure analysis with
            extremely long time scale PWscf calculations of
            varying numbers of steps, and kept for consistency with the
            Vasprun class.
        parse_dos (bool): Whether to parse the dos. If True, PWxml will use
            fildos if provided, or attempt to guess the filename.
        fildos (str): If provided, forces parse_dos to be True and uses the
            provided string as the path to the dos file. This is
            the same as in dos.x input, and shouldn't be a full path to a file.
            For example, filpdos="path/to/filpdos" will look for
            "path/to/filpdos.pdos_*"
        filpdos (str): If provided, forces parse_dos to be True and uses
             the provided string as 'filproj', same as in projwfc.x input.
             It shouldn't include the rest of the filename. For example,
             filpdos="path/to/filpdos" will look for "path/to/filpdos.pdos_*"
        parse_projected_eigen (bool): Whether to parse the projected
            eigenvalues and (magnetisation, not implemented). Defaults to False.
            If True, PWxml will look for a "filproj" from projwfc.x and parse it.
            It will look for files with the same name as the XML (or same QE prefix)
            but with a .projwfc_up extension, or will use the filproj argument
        filproj (str): If provided, forces parse_projected_eigen to be True and
            uses the provided string as the filepath to the .projwfc_up file.
            Note that this is the same as in projwfc.x input, so it shouldn't
            include the .projwfc_up/down extension. It can also include
            a directory, e.g., "path/to/filproj" will look for
            "path/to/filproj.projwfc_up"
        occu_tol (float): Sets the minimum tol for the determination of the
            vbm and cbm. Usually the default of 1e-8 works well enough,
            but there may be pathological cases. Note that, unlike VASP, QE
            actually reports the VBM and CBM (accessible via the vbm and cbm
            attributes), so this is only used to recompute them
            and check against the reported values.
        separate_spins (bool): Whether the band gap, CBM, and VBM should be
            reported for each individual spin channel. Defaults to False,
            which computes the eigenvalue band properties independent of
            the spin orientation. If True, the calculation must be spin-polarized.
    """
    self._filename = filename
    self.ionic_step_skip = ionic_step_skip
    self.ionic_step_offset = ionic_step_offset
    self.occu_tol = occu_tol
    self.separate_spins = separate_spins

    if filproj:
        parse_projected_eigen = True
    if fildos:
        parse_dos = True
    if filpdos:
        parse_pdos = True

    with zopen(filename, "rt") as f:
        self._parse(
            f,
            parse_dos=parse_dos,
            fildos=fildos,
            parse_pdos=parse_pdos,
            filpdos=filpdos,
            parse_projected_eigen=parse_projected_eigen,
            filproj=filproj,
            ionic_step_skip=ionic_step_skip,
            ionic_step_offset=ionic_step_offset,
        )

    if not self.converged:
        warnings.warn(
            (
                f"{filename} is an unconverged PWscf run.\n"
                f"Electronic convergence reached: {self.converged_electronic}.\n"
                f"Ionic convergence reached: {self.converged_ionic}."
            ),
            UnconvergedPWxmlWarning,
        )

is_spin: bool property

RETURNS DESCRIPTION
bool

True if the calculation is spin-polarized.

converged_electronic: bool property

Returns true if electronic convergence was reached in the last ionic step. Note that, even though QE internally considers NSCF calculations unconverged, this property is set to True for NSCF (and bands) calcs to maintain consistency with the Vasprun class.

RETURNS DESCRIPTION
bool

True if electronic step convergence has been reached in the final

bool

ionic step

converged_ionic: bool property

RETURNS DESCRIPTION
bool

True if ionic step convergence has been reached

To maintain consistency with the Vasprun class, we return True if the calculation didn't involve geometric optimization (scf, nscf, ...)

final_energy: float property

Final energy from the PWscf run, in eV.

hubbards: dict[str, float] property

Hubbard U values used if a vasprun is a GGA+U run. {} otherwise.

RETURNS DESCRIPTION
dict

{element: U_value}

TYPE: dict[str, float]

run_type property

Returns the run type.

Should be able to detect functional, Hubbard U terms and vdW corrections.

eigenvalue_band_properties: tuple[float, float, float, bool] | tuple[tuple[float, float], tuple[float, float], tuple[float, float], tuple[bool, bool]] property

Returns the band gap, CBM, VBM, and whether the gap is direct. Directly uses the Vasprun implementation, with some extra check against the CBM and VBM values reported by QE.

RETURNS DESCRIPTION
tuple[float, float, float, bool] | tuple[tuple[float, float], tuple[float, float], tuple[float, float], tuple[bool, bool]]

If separate_spins is False, returns a tuple of the band gap, CBM, VBM, and whether the gap is direct. If separate_spins is True, returns a tuple of tuples of the band gap, CBM, VBM, and whether the gap is direct for each spin channel.

projected_eigenvalues: dict[Spin, np.ndarray] | None property

Returns the projected eigenvalues in the same format Vasprun uses (i.e., the VASP convention)

In the case of SOC, QE uses the |LJJz> basis while VASP uses the |LLz> basis. This function will sum all p states into where py should be (index 1) and all d states into where dxy should be (index 4). The rest will be 0.

TODO: cleanup and rewrite

projected_magnetisation property

Returns the projected magnetisation for each atom in a format compatible with the Vasprun class. Not currently implemented for QE.

md_data property

Molecular dynamics data. Not currently implemented for QE.

epsilon_static: list[float] property

The static part of the dielectric constant. Not implemented for QE.

epsilon_static_wolfe: list[float] property

The static part of the dielectric constant without any local field effects. Not implemented for QE.

epsilon_ionic: list[float] property

The ionic part of the static dielectric constant. Not implemented for QE.

dielectric: tuple[list, list, list] property

The real and imaginary part of the dielectric constant (e.g., computed by RPA) in function of the energy (frequency). Optical properties (e.g. absorption coefficient) can be obtained through this.

Not implemented for QE.

optical_absorption_coeff: list[float] | None property

The optical absorption coefficient from the dielectric constants. Not implemented for QE.

get_band_structure(kpoints_filename=None, efermi='smart', line_mode=False, force_hybrid_mode=False)

Get the band structure as a BandStructure object.

PARAMETER DESCRIPTION
kpoints_filename

Path of the PWscf input file from which the band structure is generated. If none is provided, the code will try to guess the appropriate file, see the docstring of FileGuesser.

TYPE: str | None DEFAULT: None

efermi

The Fermi energy associated with the bandstructure, in eV. By default (None), uses the value reported by PWscf in the xml. To manually set the Fermi energy, pass a float. Pass 'smart' to use the calculate_efermi() method, which is identical for metals but more accurate for insulators (mid-gap).

TYPE: float | Literal['smart'] | None DEFAULT: 'smart'

line_mode

Force the band structure to be considered as a run along symmetry lines. (Default: False)

TYPE: bool DEFAULT: False

force_hybrid_mode

Not Yet Implemented (Default: False)

TYPE: bool DEFAULT: False

RETURNS DESCRIPTION
BandStructureSymmLine | BandStructure

a BandStructure object (or more specifically a BandStructureSymmLine object if the run is detected to be a run along symmetry lines) NSCF (calc='nscf' or 'bands') calculations are accepted for Line-Mode with explicit PWscf input file, and 'crystal', 'crystal_b', 'tpiba' or 'tpiba_b' K_POINTS card. The k-points needs to have data on the kpoint label as a comment.

Source code in pymatgen/io/espresso/outputs/pwxml.py
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
def get_band_structure(
    self,
    kpoints_filename: str | None = None,
    efermi: float | Literal["smart"] | None = "smart",
    line_mode: bool = False,
    force_hybrid_mode: bool = False,
) -> BandStructureSymmLine | BandStructure:
    """Get the band structure as a BandStructure object.

    Args:
        kpoints_filename: Path of the PWscf input file from which the band
            structure is generated. If none is provided, the code will try to
            guess the appropriate file, see the docstring of `FileGuesser`.
        efermi: The Fermi energy associated with the bandstructure, in eV. By
            default (None), uses the value reported by PWscf in the xml. To
            manually set the Fermi energy, pass a float. Pass 'smart' to use the
            `calculate_efermi()` method, which is identical for metals but more
            accurate for insulators (mid-gap).
        line_mode: Force the band structure to be considered as
            a run along symmetry lines. (Default: False)
        force_hybrid_mode: Not Yet Implemented (Default: False)

    Returns:
        a BandStructure object (or more specifically a
            BandStructureSymmLine object if the run is detected to be a run
            along symmetry lines) NSCF (calc='nscf' or 'bands') calculations are accepted for Line-Mode with explicit PWscf input file, and 'crystal', 'crystal_b', 'tpiba' or 'tpiba_b' K_POINTS card. The k-points needs to have data on the kpoint label as a comment.
    """
    factor = 1 if self.noncolin else 2
    if self.nbands <= self.nelec / factor:
        warnings.warn(
            f"Number of bands ({self.nbands}) <= number of electrons/{factor} "
            f"({self.nelec / factor:.4f}). BSPlotter may not work properly.",
            DifferentFromVASPWarning,
        )

    if not kpoints_filename:
        kpoints_filename = FileGuesser("pwin", self._filename, self.prefix).guess()

    if kpoints_filename and not os.path.exists(kpoints_filename) and line_mode:
        raise PWxmlParserError(
            "PW input file needed to obtain band structure along symmetry lines."
        )

    if efermi == "smart":
        e_fermi = self.calculate_efermi()
    elif efermi is None:
        e_fermi = self.efermi
    else:
        e_fermi = efermi

    lattice_new = Lattice(self.final_structure.lattice.reciprocal_lattice.matrix)

    p_eigenvals: defaultdict[Spin, list] = defaultdict(list)
    eigenvals: defaultdict[Spin, list] = defaultdict(list)

    for spin, v in self.eigenvalues.items():
        v = np.swapaxes(v, 0, 1)
        eigenvals[spin] = v[:, :, 0]

        if self.projected_eigenvalues:
            peigen = self.projected_eigenvalues[spin]
            # Original axes for self.projected_eigenvalues are kpoints,
            # band, ion, orb.
            # For BS input, we need band, kpoints, orb, ion.
            peigen = np.swapaxes(peigen, 0, 1)  # Swap kpoint and band axes
            peigen = np.swapaxes(peigen, 2, 3)  # Swap ion and orb axes

            p_eigenvals[spin] = peigen

    k_card = None
    if kpoints_filename and os.path.exists(kpoints_filename):
        k_card = PWin.from_file(kpoints_filename).k_points
    coords_are_cartesian = False if k_card is None else k_card.coords_are_cartesian
    if coords_are_cartesian:
        kpoints = [np.array(kpt) for kpt in self.kpoints_cart]
    else:
        kpoints = [np.array(kpt) for kpt in self.kpoints_frac]

    if k_card.line_mode or line_mode:
        labels_dict = {}
        # TODO: check how hybrid band structs work in QE
        hybrid_band = False
        if hybrid_band or force_hybrid_mode:
            raise NotImplementedError(
                "Hybrid band structures not yet supported in line mode."
            )
        kpoints, eigenvals, p_eigenvals, labels_dict = self._vaspify_kpts_bands(
            kpoints, eigenvals, p_eigenvals, k_card, self.alat
        )
        return BandStructureSymmLine(
            kpoints,
            eigenvals,
            lattice_new,
            e_fermi,
            labels_dict,
            structure=self.final_structure,
            projections=p_eigenvals,
            coords_are_cartesian=coords_are_cartesian,
        )
    return BandStructure(
        kpoints,
        eigenvals,
        lattice_new,
        e_fermi,
        structure=self.final_structure,
        projections=p_eigenvals,
        coords_are_cartesian=coords_are_cartesian,
    )

calculate_efermi(**_)

Calculate the Fermi level

PWscf returns the Fermi level for all calculations and the cbm and vbm for all insulators. These are stored in PWxml.efermi, PWxml.cbm, and PWxml.vbm. However, for insulators, the Fermi level is often slightly off from the exact mid-gap value.

If vbm and cbm are both undefined (metallic system), return the Fermi level if vbm is defined and cbm isn't, it's usually a sign of an insulator with as many bands as electrons (often nbnd isn't set in input) Such calculations don't work with BSPlotter()

RETURNS DESCRIPTION
float

The Fermi level, in eV

Source code in pymatgen/io/espresso/outputs/pwxml.py
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
def calculate_efermi(self, **_) -> float:
    """
    Calculate the Fermi level

    PWscf returns the Fermi level for all calculations and the cbm and vbm for
    all insulators. These are stored in PWxml.efermi, PWxml.cbm, and PWxml.vbm.
    However, for insulators, the Fermi level is often slightly off from the
    exact mid-gap value.

    If vbm and cbm are both undefined (metallic system), return the Fermi level
    if vbm is defined and cbm isn't, it's usually a sign of an insulator
    with as many bands as electrons (often nbnd isn't set in input)
    Such calculations don't work with BSPlotter()

    Returns:
        The Fermi level, in eV
    """
    if self.vbm is None or self.cbm is None:
        return self.efermi
    return (self.vbm + self.cbm) / 2

UnconvergedPWxmlWarning

Bases: Warning

Warning for unconverged PWscf run from xml file

PWxmlParserError

Bases: Exception

Exception class for PWxml parsing.

InconsistentWithXMLError

Bases: Exception

Exception class for data from external files that is inconsistent with the XML file.

DifferentFromVASPWarning

Bases: Warning

Warning for differences between QE and VASP outputs

NotImplementedForQE

Bases: Warning

Warning for differences between QE and VASP outputs

ZeroTotalEnergyWarning

Bases: Warning

Warning for zero total energy. Happens with bands/NSCF calcs in QE