Espresso

../../_images/espresso.png

Quantum ESPRESSO (QE) is an integrated suite of Open-Source computer codes for electronic-structure calculations and materials modeling at the nanoscale. It is based on density-functional theory, plane waves, and pseudopotentials. The ASE calculator is an interface to the pw.x executable, however, input/output operations can be managed for other executables (see the Input/Output section below). Users interested to run calculation with other executables, can have a look in the ase ecosystem pages for other packages.

Calculator

Basic usage

Calculations using the Espresso calculator can be run by defining an EspressoProfile and setting up an Espresso calculator with it:

Any calculation will need pseudopotentials for the elements involved. The directory for the pseudopotential files can be set with the pseudo_dir parameter in the EspressoProfile. Pseudopotential files can be passed to the calculator as a dictionary, where the keys are the element and values are the pseudopotential file names. The pseudopotential file names must be in the directory specified by the pseudo_dir parameter.

pseudopotentials = {"Na": "na_pbe_v1.5.uspp.F.UPF", "Cl": "cl_pbe_v1.4.uspp.F.UPF"}

A recommended list of pseudopotentials is the SSSP curated set, which can be found on Materials Cloud. With the all of the above in mind, a simple calculation can be set up like so:

from ase.build import bulk
from ase.calculators.espresso import Espresso, EspressoProfile
from ase.optimize import LBFGS

rocksalt = bulk("NaCl", crystalstructure="rocksalt", a=6.0)

# Pseudopotentials from SSSP Efficiency v1.3.0
pseudopotentials = {"Na": "na_pbe_v1.5.uspp.F.UPF", "Cl": "cl_pbe_v1.4.uspp.F.UPF"}

profile = EspressoProfile(
    binary="/path/to/pw.x", pseudo_dir="/path/to/pseudopotentials"
)

calc = Espresso(profile=profile, pseudopotentials=pseudopotentials)

rocksalt.calc = calc

rocksalt.get_potential_energy()  # This will run a single point calculation

opt = LBFGS(rocksalt)

opt.run(fmax=0.005)  # This will run a geometry optimization using ASE's LBFGS algorithm

# Print lattice constant...
print((8 * rocksalt.get_volume() / len(rocksalt)) ** (1.0 / 3.0))

The above example will run a single point calculation and a geometry optimization using default QE settings. In the next sections, we will describe how to customize the QE settings.

Parameters

The calculator will interpret most documented options for pw.x found at the PW input description.

Parameters must be passed in the input_data parameter, which is a dictionary , previously accepted kwargs keywords are now deprecated. The input_data dictionary can be used both in flat or nested format. In the flat format, keywords will be put into the correct section of the input file. In the nested format, the input_data dictionary will be used as is. Currently the nested format can be used to specify exotic input sections manually such as &FCP or &RISM.

All parameters must be given in QE units, usually Ry or atomic units in line with the documentation. ASE does not add any defaults over the defaults of QE.

input_data = {
    "system": {"ecutwfc": 60, "ecutrho": 480},
    "disk_io": "low",  # Automatically put into the 'control' section
}

calc = Espresso(
    profile=profile
    pseudopotentials=pseudopotentials,
    tstress=True,  # deprecated, put in input_data
    tprnfor=True,  # deprecated, put in input_data
    input_data=input_data,
)

Any FixAtoms or FixCartesian constraints are converted to Espresso constraints. Some parameters are used by ASE and have additional meaning:

  • kpts, is used to specify the k-point sampling * If kpts is a tuple (or list) of 3 integers kpts=(int, int, int), it is interpreted as the dimensions of a Monkhorst-Pack grid. * If kpts is set to None, only the Γ-point will be included and QE will use routines optimized for Γ-point-only calculations. Compared to Γ-point-only calculations without this optimization (i.e. with kpts=(1, 1, 1)), the memory and CPU requirements are typically reduced by half. * If kpts is a dict, it will either be interpreted as a path in the Brillouin zone (see ase.dft.kpoints.bandpath) if it contains the ‘path’ keyword, otherwise it is converted to a Monkhorst-Pack grid (see ase.calculators.calculator.kpts2sizeandoffsets)

  • koffset=(0, 0, 0) set to 0 or 1 to displace the kpoint grid by a half cell in that direction. Also accepts True and False.

  • kspacing=0.1 sets the minimum distance between kpoints in reciprocal space.

  • nspin does not typically need to be specified by the user. If any atom has a magnetic moment, it is turned on automatically.

Parallelism

The above calc object will run the pw.x executable on 1 core. parallelism can be set by providing a parallel_info dictionary to the EspressoProfile. Currently there is no way to specify the parallelization keywords for QE, such as -nk. This is a current limitation of the ASE calculator interface.

parallel_info = {"binary": "mpirun", "np": 16, "-v": True}

profile = EspressoProfile(
    binary="/path/to/pw.x",
    pseudo_dir="/path/to/pseudopotentials",
    parallel_info=parallel_info,
)

# This will run "mpirun -np 16 -v /path/to/pw.x"

Input/Output

ASE can read and write input files for the espresso executables. The procedure depends on the executable. Each subsection will describe the procedure for each executable.

pw.x

The pw.x executable is the main executable of the Quantum Espresso suite. ASE can fully read and write input files for this executable. This must be done by using the main ase.io.write() and ase.io.read() functions. The

These two main functions are wrappers around the ase.io.espresso.write_espresso_in() and ase.io.espresso.read_espresso_out() functions. The former will write an input file for the pw.x executable, while the latter will read an input file and return Atoms object(s) with available results.

ase.io.espresso.write_espresso_in(fd, atoms, input_data=None, pseudopotentials=None, kspacing=None, kpts=None, koffset=(0, 0, 0), crystal_coordinates=False, additional_cards=None, **kwargs)[source]

Create an input file for pw.x.

Use set_initial_magnetic_moments to turn on spin, if nspin is set to 2 with no magnetic moments, they will all be set to 0.0. Magnetic moments will be converted to the QE units (fraction of valence electrons) using any pseudopotential files found, or a best guess for the number of valence electrons.

Units are not converted for any other input data, so use Quantum ESPRESSO units (Usually Ry or atomic units).

Keys with a dimension (e.g. Hubbard_U(1)) will be incorporated as-is so the \(i\) should be made to match the output.

Implemented features:

  • Conversion of ase.constraints.FixAtoms and ase.constraints.FixCartesian.

  • starting_magnetization derived from the magmoms and pseudopotentials (searches default paths for pseudo files.)

  • Automatic assignment of options to their correct sections.

Not implemented:

  • Non-zero values of ibrav

  • Lists of k-points

  • Other constraints

  • Hubbard parameters

  • Validation of the argument types for input

  • Validation of required options

Parameters:
  • fd (file | str) – A file to which the input is written.

  • atoms (Atoms) – A single atomistic configuration to write to fd.

  • input_data (dict) – A flat or nested dictionary with input parameters for pw.x

  • pseudopotentials (dict) – A filename for each atomic species, e.g. {‘O’: ‘O.pbe-rrkjus.UPF’, ‘H’: ‘H.pbe-rrkjus.UPF’}. A dummy name will be used if none are given.

  • kspacing (float) – Generate a grid of k-points with this as the minimum distance, in A^-1 between them in reciprocal space. If set to None, kpts will be used instead.

  • kpts ((int, int, int) or dict) – If kpts is a tuple (or list) of 3 integers, it is interpreted as the dimensions of a Monkhorst-Pack grid. If kpts is set to None, only the Γ-point will be included and QE will use routines optimized for Γ-point-only calculations. Compared to Γ-point-only calculations without this optimization (i.e. with kpts=(1, 1, 1)), the memory and CPU requirements are typically reduced by half. If kpts is a dict, it will either be interpreted as a path in the Brillouin zone (*) if it contains the ‘path’ keyword, otherwise it is converted to a Monkhorst-Pack grid (**). (*) see ase.dft.kpoints.bandpath (**) see ase.calculators.calculator.kpts2sizeandoffsets

  • koffset ((int, int, int)) – Offset of kpoints in each direction. Must be 0 (no offset) or 1 (half grid offset). Setting to True is equivalent to (1, 1, 1).

  • crystal_coordinates (bool) – Whether the atomic positions should be written to the QE input file in absolute (False, default) or relative (crystal) coordinates (True).

Example:

from ase.io.espresso import write_espresso_in

input_data = {
    "calculation": "relax",
    "restart_mode": "from_scratch",
    "tprnfor": True,
    "etot_conv_thr": 1e-5,
    "forc_conv_thr": 1e-4,
    "ecutwfc": 60,
    "ecutrho": 480,
    "input_dft": "rpbe",
    "vdw_corr": "dft-d3",
    "occupations": "smearing",
    "degauss": 0.01,
    "smearing": "cold",
    "conv_thr": 1e-8,
    "mixing_mode": "local-TF",
    "mixing_beta": 0.35,
    "diagonalization": "david",
    "ion_dynamics": "bfgs",
    "bfgs_ndim": 6,
    "startingwfc": "random",
}  # This flat dictionary will be converted to a nested dictionary where, for example, "calculation" will be put into the "control" section

write("pw.in", atoms, input_data=input_data, format="espresso-in")

# Hubbard is not implemented in the write_espresso_in function, we can add it manually
additional_cards = ["HUBBARD (ortho-atomic)", "U Mn-3d 5.0", "U Ni-3d 6.0"]

write(
    "pw_hubbard.in",
    atoms,
    input_data=input_data,
    additional_cards=additional_cards,
    format="espresso-in",
)
ase.io.espresso.read_espresso_out(fileobj, index=slice(None, None, None), results_required=True)[source]

Reads Quantum ESPRESSO output files.

The atomistic configurations as well as results (energy, force, stress, magnetic moments) of the calculation are read for all configurations within the output file.

Will probably raise errors for broken or incomplete files.

Parameters:
  • fileobj (file|str) – A file like object or filename

  • index (slice) – The index of configurations to extract.

  • results_required (bool) – If True, atomistic configurations that do not have any associated results will not be included. This prevents double printed configurations and incomplete calculations from being returned as the final configuration with no results data.

Yields:

structure (Atoms) – The next structure from the index slice. The Atoms has a SinglePointCalculator attached with any results parsed from the file.

ph.x

As of January 2024, ph.x has custom read and write functions in ASE. The function write_espresso_ph() can be used to write a Fortran namelist file for ph.x.

ase.io.espresso.write_espresso_ph(fd, input_data=None, qpts=None, nat_todo_indices=None, **kwargs) None[source]

Function that write the input file for a ph.x calculation. Normal namelist cards are passed in the input_data dictionary. Which can be either nested or flat, ASE style. The q-points are passed in the qpts list. If qplot is set to True then qpts is expected to be a list of list|tuple of length 4. Where the first three elements are the coordinates of the q-point in units of 2pi/alat and the last element is the weight of the q-point. if qplot is set to False then qpts is expected to be a simple list of length 4 (single q-point). Finally if ldisp is set to True, the above is discarded and the q-points are read from the nq1, nq2, nq3 cards in the input_data dictionary.

Additionally, a nat_todo_indices kwargs (list[int]) can be specified in the kwargs. It will be used if nat_todo is set to True in the input_data dictionary.

Globally, this function follows the convention set in the ph.x documentation (https://www.quantum-espresso.org/Doc/INPUT_PH.html)

Parameters:
  • fd – The file descriptor of the input file.

  • kwargs

    kwargs dictionary possibly containing the following keys:

    • input_data: dict

    • qpts: list[list[float]] | list[tuple[float]] | list[float]

    • nat_todo_indices: list[int]

Return type:

None

Example:

from ase.io.espresso import write_espresso_ph, read_espresso_ph

input_data = {
  'tr2_ph': 1.0e-12,
  'prefix': 'pwscf',
  'verbosity': 'high',
  'ldisp': True,
  'qplot': True,
  'alpha_mix(1)': 0.1,
}

qpts = [(0.0, 0.0, 0.0), (0.5, 0.5, 0.5), (0.5, 0.5, 0.0), (0.0, 0.0, 0.5)]

write_espresso_ph("input_ph.in", input_data, qpts=qpts) # Will automatically be built to nested format

After running the calculation, the output can be read with the read_espresso_ph() method.

ase.io.espresso.read_espresso_ph(fileobj)[source]

Function that reads the output of a ph.x calculation. It returns a dictionary where each q-point number is a key and the value is a dictionary with the following keys if available:

  • qpoints: The q-point in cartesian coordinates.

  • kpoints: The k-points in cartesian coordinates.

  • dieltensor: The dielectric tensor.

  • borneffcharge: The effective Born charges.

  • borneffcharge_dfpt: The effective Born charges from DFPT.

  • polarizability: The polarizability tensor.

  • modes: The phonon modes.

  • eqpoints: The symmetrically equivalent q-points.

  • freqs: The phonon frequencies.

  • mode_symmetries: The symmetries of the modes.

  • atoms: The atoms object.

Some notes:

  • For some reason, the cell is not defined to high level of precision in ph.x outputs. Be careful when using the atoms object retrieved from this function.

  • This function can be called on incomplete calculations i.e. if the calculation couldn’t diagonalize the dynamical matrix for some q-points, the results for the other q-points will still be returned.

Parameters:

fileobj – The file descriptor of the output file.

Returns:

The results dictionnary as described above.

Return type:

dict

Example:

results = read_espresso_ph("ph.out")

print(results[0]['qpoints']) # Will probably be (0, 0, 0)
print(results[0]['freq']) # Eigenvalues of the dynamical matrix at the qpoint

Other binaries

Since January 2024, it is possible to use ASE to write input files for other QE executables.

The ase.io.espresso.write_fortran_namelist() method will write a fortran namelist file, which can be read by QE executables. The list of currently implemented executable is available in ase.io.espresso_namelist.keys

ase.io.espresso.write_fortran_namelist(fd, input_data=None, binary=None, additional_cards=None, **kwargs) None[source]

Function which writes input for simple espresso binaries. List of supported binaries are in the espresso_keys.py file. Non-exhaustive list (to complete)

Note: “EOF” is appended at the end of the file. (https://lists.quantum-espresso.org/pipermail/users/2020-November/046269.html)

Additional fields are written ‘as is’ in the input file. It is expected to be a string or a list of strings.

Parameters:
  • fd – The file descriptor of the input file.

  • input_data (dict) – A flat or nested dictionary with input parameters for the binary.x

  • binary (str) – Name of the binary

  • additional_cards (str | list[str]) – Additional fields to be written at the end of the input file, after the namelist. It is expected to be a string or a list of strings.

Return type:

None

Example:

from ase.io.espresso import Namelist
from ase.io.espresso import write_fortran_namelist

# input_data for pp.x, built automatically
input_data = {"outdir": "/path/to/output", "prefix": "prefix", "verbosity": "high"}

input_data = Namelist(input_data)
input_data.to_nested(binary="pp")

write_fortran_namelist("input_pp.in", input_data)

# Alternatively, the input_data can be written manually
# without the need to specify a binary

input_data = {
  'environ': {
    'environ_type': 'water',
  ...} # In this case, the ``input_data`` is a nested dictionary
  # since ASE cannot guess the input sections for exotic executables

additional_cards = [
      "EXTERNAL_CHARGES (bohr)",
      "-0.5 0. 0. 25.697 1.0 2 3",
      "-0.5 0. 0. 20.697 1.0 2 3"
] # This will be added at the end of the file

write_fortran_namelist("environ.in", input_data)

The read_fortran_namelist() method will read a fortran namelist file and return a tuple, where the first element is the input_data dictionary and the second element is a list of additional cards at the end of the file.

ase.io.espresso.read_fortran_namelist(fileobj)[source]

Takes a fortran-namelist formatted file and returns nested dictionaries of sections and key-value data, followed by a list of lines of text that do not fit the specifications. Behaviour is taken from Quantum ESPRESSO 5.3. Parses fairly convoluted files the same way that QE should, but may not get all the MANDATORY rules and edge cases for very non-standard files Ignores anything after ‘!’ in a namelist, split pairs on ‘,’ to include multiple key=values on a line, read values on section start and end lines, section terminating character, ‘/’, can appear anywhere on a line. All of these are ignored if the value is in ‘quotes’.

Parameters:

fileobj (file) – An open file-like object.

Returns:

  • data (dict[str, dict]) – Dictionary for each section in the namelist with key = value pairs of data.

  • additional_cards (list[str]) – Any lines not used to create the data, assumed to belong to ‘cards’ in the input file.

Example:

from ase.io.espresso import read_fortran_namelist

input_data, additional_cards = read_fortran_namelist("input_pp.in")

print(input_data)
print(additional_cards)

You can use both of these methods to jungle between your favorite QE executables.

Namelist Class

The Namelist class is a simple class to handle fortran namelist files. It derived from the dict class and is case case-insensitive.

class ase.io.espresso.Namelist(dict=None, /, **kwargs)[source]

A case-insensitive dictionary for storing Quantum Espresso namelists. This class is a subclass of UserDict, which is a wrapper around a regular dictionary. This allows us to define custom behavior for the dictionary methods, while still having access to the full dictionary API.

to_string() have been added to handle the conversion of the dictionary to a string for writing to a file or quick lookup using print().

to_nested() have been added to convert the dictionary to a nested dictionary with the correct structure for the specified binary.

Espresso Calculator Class

class ase.calculators.espresso.Espresso(*, profile=None, command=<object object>, label=<object object>, directory='.', parallel_info=None, parallel=True, **kwargs)[source]

All options for pw.x are copied verbatim to the input file, and put into the correct section. Use input_data for parameters that are already in a dict.

input_data: dict

A flat or nested dictionary with input parameters for pw.x

pseudopotentials: dict

A filename for each atomic species, e.g. {'O': 'O.pbe-rrkjus.UPF', 'H': 'H.pbe-rrkjus.UPF'}. A dummy name will be used if none are given.

kspacing: float

Generate a grid of k-points with this as the minimum distance, in A^-1 between them in reciprocal space. If set to None, kpts will be used instead.

kpts: (int, int, int), dict, or BandPath

If kpts is a tuple (or list) of 3 integers, it is interpreted as the dimensions of a Monkhorst-Pack grid. If kpts is set to None, only the Γ-point will be included and QE will use routines optimized for Γ-point-only calculations. Compared to Γ-point-only calculations without this optimization (i.e. with kpts=(1, 1, 1)), the memory and CPU requirements are typically reduced by half. If kpts is a dict, it will either be interpreted as a path in the Brillouin zone (*) if it contains the ‘path’ keyword, otherwise it is converted to a Monkhorst-Pack grid (**). (*) see ase.dft.kpoints.bandpath (**) see ase.calculators.calculator.kpts2sizeandoffsets

koffset: (int, int, int)

Offset of kpoints in each direction. Must be 0 (no offset) or 1 (half grid offset). Setting to True is equivalent to (1, 1, 1).