Hot-keys on this page
r m x p toggle line displays
j k next/prev highlighted chunk
0 (zero) top of page
1 (one) first highlighted chunk
1"""Quantum ESPRESSO Calculator
3export ASE_ESPRESSO_COMMAND="/path/to/pw.x -in PREFIX.pwi > PREFIX.pwo"
5Run pw.x jobs.
6"""
9import warnings
10from ase import io
11from ase.calculators.calculator import FileIOCalculator, PropertyNotPresent
14error_template = 'Property "%s" not available. Please try running Quantum\n' \
15 'Espresso first by calling Atoms.get_potential_energy().'
17warn_template = 'Property "%s" is None. Typically, this is because the ' \
18 'required information has not been printed by Quantum ' \
19 'Espresso at a "low" verbosity level (the default). ' \
20 'Please try running Quantum Espresso with "high" verbosity.'
23class Espresso(FileIOCalculator):
24 """
25 """
26 implemented_properties = ['energy', 'forces', 'stress', 'magmoms']
27 command = 'pw.x -in PREFIX.pwi > PREFIX.pwo'
28 discard_results_on_any_change = True
30 def __init__(self, restart=None,
31 ignore_bad_restart_file=FileIOCalculator._deprecated,
32 label='espresso', atoms=None, **kwargs):
33 """
34 All options for pw.x are copied verbatim to the input file, and put
35 into the correct section. Use ``input_data`` for parameters that are
36 already in a dict, all other ``kwargs`` are passed as parameters.
38 Accepts all the options for pw.x as given in the QE docs, plus some
39 additional options:
41 input_data: dict
42 A flat or nested dictionary with input parameters for pw.x
43 pseudopotentials: dict
44 A filename for each atomic species, e.g.
45 ``{'O': 'O.pbe-rrkjus.UPF', 'H': 'H.pbe-rrkjus.UPF'}``.
46 A dummy name will be used if none are given.
47 kspacing: float
48 Generate a grid of k-points with this as the minimum distance,
49 in A^-1 between them in reciprocal space. If set to None, kpts
50 will be used instead.
51 kpts: (int, int, int), dict, or BandPath
52 If kpts is a tuple (or list) of 3 integers, it is interpreted
53 as the dimensions of a Monkhorst-Pack grid.
54 If ``kpts`` is set to ``None``, only the Γ-point will be included
55 and QE will use routines optimized for Γ-point-only calculations.
56 Compared to Γ-point-only calculations without this optimization
57 (i.e. with ``kpts=(1, 1, 1)``), the memory and CPU requirements
58 are typically reduced by half.
59 If kpts is a dict, it will either be interpreted as a path
60 in the Brillouin zone (*) if it contains the 'path' keyword,
61 otherwise it is converted to a Monkhorst-Pack grid (**).
62 (*) see ase.dft.kpoints.bandpath
63 (**) see ase.calculators.calculator.kpts2sizeandoffsets
64 koffset: (int, int, int)
65 Offset of kpoints in each direction. Must be 0 (no offset) or
66 1 (half grid offset). Setting to True is equivalent to (1, 1, 1).
69 .. note::
70 Set ``tprnfor=True`` and ``tstress=True`` to calculate forces and
71 stresses.
73 .. note::
74 Band structure plots can be made as follows:
77 1. Perform a regular self-consistent calculation,
78 saving the wave functions at the end, as well as
79 getting the Fermi energy:
81 >>> input_data = {<your input data>}
82 >>> calc = Espresso(input_data=input_data, ...)
83 >>> atoms.calc = calc
84 >>> atoms.get_potential_energy()
85 >>> fermi_level = calc.get_fermi_level()
87 2. Perform a non-self-consistent 'band structure' run
88 after updating your input_data and kpts keywords:
90 >>> input_data['control'].update({'calculation':'bands',
91 >>> 'restart_mode':'restart',
92 >>> 'verbosity':'high'})
93 >>> calc.set(kpts={<your Brillouin zone path>},
94 >>> input_data=input_data)
95 >>> calc.calculate(atoms)
97 3. Make the plot using the BandStructure functionality,
98 after setting the Fermi level to that of the prior
99 self-consistent calculation:
101 >>> bs = calc.band_structure()
102 >>> bs.reference = fermi_energy
103 >>> bs.plot()
105 """
106 FileIOCalculator.__init__(self, restart, ignore_bad_restart_file,
107 label, atoms, **kwargs)
108 self.calc = None
110 def write_input(self, atoms, properties=None, system_changes=None):
111 FileIOCalculator.write_input(self, atoms, properties, system_changes)
112 io.write(self.label + '.pwi', atoms, **self.parameters)
114 def read_results(self):
115 output = io.read(self.label + '.pwo')
116 self.calc = output.calc
117 self.results = output.calc.results
119 def get_fermi_level(self):
120 if self.calc is None:
121 raise PropertyNotPresent(error_template % 'Fermi level')
122 return self.calc.get_fermi_level()
124 def get_ibz_k_points(self):
125 if self.calc is None:
126 raise PropertyNotPresent(error_template % 'IBZ k-points')
127 ibzkpts = self.calc.get_ibz_k_points()
128 if ibzkpts is None:
129 warnings.warn(warn_template % 'IBZ k-points')
130 return ibzkpts
132 def get_k_point_weights(self):
133 if self.calc is None:
134 raise PropertyNotPresent(error_template % 'K-point weights')
135 k_point_weights = self.calc.get_k_point_weights()
136 if k_point_weights is None:
137 warnings.warn(warn_template % 'K-point weights')
138 return k_point_weights
140 def get_eigenvalues(self, **kwargs):
141 if self.calc is None:
142 raise PropertyNotPresent(error_template % 'Eigenvalues')
143 eigenvalues = self.calc.get_eigenvalues(**kwargs)
144 if eigenvalues is None:
145 warnings.warn(warn_template % 'Eigenvalues')
146 return eigenvalues
148 def get_number_of_spins(self):
149 if self.calc is None:
150 raise PropertyNotPresent(error_template % 'Number of spins')
151 nspins = self.calc.get_number_of_spins()
152 if nspins is None:
153 warnings.warn(warn_template % 'Number of spins')
154 return nspins