1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 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
|
"""Quantum ESPRESSO Calculator
Run pw.x jobs.
"""
import os
import warnings
from ase.calculators.genericfileio import (
BaseProfile,
CalculatorTemplate,
GenericFileIOCalculator,
read_stdout,
)
from ase.io import read, write
from ase.io.espresso import Namelist
compatibility_msg = (
'Espresso calculator is being restructured. Please use e.g. '
"Espresso(profile=EspressoProfile(argv=['mpiexec', 'pw.x'])) "
'to customize command-line arguments.'
)
# XXX We should find a way to display this warning.
# warn_template = 'Property "%s" is None. Typically, this is because the ' \
# 'required information has not been printed by Quantum ' \
# 'Espresso at a "low" verbosity level (the default). ' \
# 'Please try running Quantum Espresso with "high" verbosity.'
class EspressoProfile(BaseProfile):
configvars = {'pseudo_dir'}
def __init__(self, command, pseudo_dir, **kwargs):
super().__init__(command, **kwargs)
# not Path object to avoid problems in remote calculations from Windows
self.pseudo_dir = str(pseudo_dir)
@staticmethod
def parse_version(stdout):
import re
match = re.match(r'\s*Program PWSCF\s*v\.(\S+)', stdout, re.M)
assert match is not None
return match.group(1)
def version(self):
stdout = read_stdout(self._split_command)
return self.parse_version(stdout)
def get_calculator_command(self, inputfile):
return ['-in', inputfile]
class EspressoTemplate(CalculatorTemplate):
_label = 'espresso'
def __init__(self):
super().__init__(
'espresso',
['energy', 'free_energy', 'forces', 'stress', 'magmoms', 'dipole'],
)
self.inputname = f'{self._label}.pwi'
self.outputname = f'{self._label}.pwo'
self.errorname = f"{self._label}.err"
def write_input(self, profile, directory, atoms, parameters, properties):
dst = directory / self.inputname
input_data = Namelist(parameters.pop("input_data", None))
input_data.to_nested("pw")
input_data["control"].setdefault("pseudo_dir", str(profile.pseudo_dir))
parameters["input_data"] = input_data
write(
dst,
atoms,
format='espresso-in',
properties=properties,
**parameters,
)
def execute(self, directory, profile):
profile.run(directory, self.inputname, self.outputname,
errorfile=self.errorname)
def read_results(self, directory):
path = directory / self.outputname
atoms = read(path, format='espresso-out')
return dict(atoms.calc.properties())
def load_profile(self, cfg, **kwargs):
return EspressoProfile.from_config(cfg, self.name, **kwargs)
def socketio_parameters(self, unixsocket, port):
return {}
def socketio_argv(self, profile, unixsocket, port):
if unixsocket:
ipi_arg = f'{unixsocket}:UNIX'
else:
ipi_arg = f'localhost:{port:d}' # XXX should take host, too
return profile.get_calculator_command(self.inputname) + [
'--ipi',
ipi_arg,
]
class Espresso(GenericFileIOCalculator):
def __init__(
self,
*,
profile=None,
command=GenericFileIOCalculator._deprecated,
label=GenericFileIOCalculator._deprecated,
directory='.',
**kwargs,
):
"""
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).
"""
if command is not self._deprecated:
raise RuntimeError(compatibility_msg)
if label is not self._deprecated:
warnings.warn(
'Ignoring label, please use directory instead', FutureWarning
)
if 'ASE_ESPRESSO_COMMAND' in os.environ and profile is None:
warnings.warn(compatibility_msg, FutureWarning)
template = EspressoTemplate()
super().__init__(
profile=profile,
template=template,
directory=directory,
parameters=kwargs,
)
|