Coverage for /builds/ase/ase/ase/calculators/onetep.py : 10.46%

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"""This module defines an interface to ONETEP for use by the ASE.
3Authors:
4 Edward Tait, ewt23@cam.ac.uk
5 Nicholas Hine, n.d.m.hine@warwick.ac.uk (current maintainer)
7 Based on castep.py by:
8 Max Hoffmann, max.hoffmann@ch.tum.de
9 Jörg Meyer, joerg.meyer@ch.tum.de
10"""
12from copy import deepcopy
13from os.path import isfile
14from warnings import warn
16from numpy import array
18from ase import Atoms
19from ase.calculators.calculator import FileIOCalculator, ReadError
20from ase.parallel import paropen
21from ase.units import Bohr, Hartree
24__all__ = ['Onetep']
27class Onetep(FileIOCalculator):
28 """Implements the calculator for the onetep linear
29 scaling DFT code. Recommended ASE_ONETEP_COMMAND format
30 is "onetep_executable_name PREFIX.dat > PREFIX.out 2> PREFIX.err" """
32 implemented_properties = ['energy', 'forces', 'dipole', 'magmom']
34 # Used to indicate 'parameters' which shouldn't be written to
35 # the onetep input file in the standard <key> : <value> format
36 # for example the NGWF radius is used in the species block and isn't
37 # written elsewhere in the input file
38 _dummy_parameters = ['ngwf_radius', 'xc', 'species_ngwf_radius',
39 'species_ngwf_number', 'species_solver',
40 'ngwf_radius_cond', 'pseudo_suffix',
41 'species_pseudo', 'species_core_wf',
42 'species_solver_cond', 'species_ngwf_number_cond',
43 'species_ngwf_radius_cond']
45 # Used to indicate which parameters are a kpoint path and should be
46 # written as such
47 _path_parameters = ['bsunfld_kpoint_path', 'bs_kpoint_path']
49 # Used to indicate which parameters are a block listing atom
50 # groupings for a variety of purposes
51 _group_parameters = ['species_bsunfld_groups', 'species_ldos_groups',
52 'species_locdipole_groups',
53 'species_bsunfld_projatoms',
54 'species_pdos_groups', 'species_tddft_ct',
55 'species_tddft_kernel', 'nbo_write_species',
56 'species_ngwf_plot']
58 # Used to indicate which parameters are a block of any other sort
59 # other than those above (the contents of the parameter is reproduced
60 # verbatim within the block)
61 _block_parameters = _path_parameters + _group_parameters + [
62 'species_constraints', 'nbo_species_ngwflabel',
63 'ddec_rmse_vdw', 'vdw_params', 'sol_ions', 'swri']
65 default_parameters = {'cutoff_energy': '1000 eV',
66 'kernel_cutoff': '1000 bohr',
67 'ngwf_radius': 12.0,
68 'ngwf_radius_cond': -1.0}
70 name = 'onetep'
72 def __init__(self, restart=None,
73 ignore_bad_restart_file=FileIOCalculator._deprecated,
74 label=None, command=None, atoms=None, **kwargs):
75 FileIOCalculator.__init__(self, restart, ignore_bad_restart_file,
76 label, atoms, command, **kwargs)
78 self.species = []
79 self.species_cond = []
80 self.pseudos = []
81 self.core_wfs = []
82 self.solvers = []
83 self.solvers_cond = []
84 self.restart = False
85 self.prefix = label
86 self.directory = '.'
88 def read(self, label):
89 """Read a onetep .out file into the current instance."""
91 FileIOCalculator.read(self, label)
93 onetep_file = self.label + '.out'
95 warnings = []
97 try:
98 out = paropen(onetep_file, 'r')
99 except IOError:
100 raise ReadError('Could not open output file "%s"' % onetep_file)
102 # keep track of what we've read in
103 read_lattice = False
104 read_species = False
105 read_positions = False
107 line = out.readline()
109 if self.atoms is None:
110 self.atoms = Atoms()
111 self.atoms.calc = self
113 while line:
114 clean_line = line.strip().lower()
115 if '%block lattice_cart' in clean_line:
116 self._read_lattice(out)
117 read_lattice = True
118 elif '%block species_pot' in clean_line:
119 self._read_species_pot(out)
120 elif clean_line.endswith('%block species_atomic_set'):
121 self._read_species_solver(out)
122 elif clean_line.endswith('%block species'):
123 self._read_species(out)
124 read_species = True
125 elif '%block positions_abs' in clean_line:
126 self._read_positions(out)
127 read_positions = True
128 elif '%block species_cond' in clean_line:
129 self._read_species(out, cond=True)
130 elif '%block species_atomic_set_cond' in clean_line:
131 self._read_species_solver(out, cond=True)
132 elif 'warn' in line.lower():
133 warnings.append(line)
134 line = out.readline()
135 out.close()
137 if warnings:
138 warn('WARNING: %s contains warnings' % onetep_file)
139 for warning in warnings:
140 warn(warning)
142 if not (read_lattice and read_species and read_positions):
143 raise ReadError('Failed to read in essential calculation'
144 ' data from output file "%s"' % onetep_file)
146 self.read_results(label)
148 def read_results(self, label=None):
149 FileIOCalculator.read_results(self)
151 if label is None:
152 onetep_file = self.label + '.out'
153 else:
154 onetep_file = label + '.out'
156 warnings = []
158 try:
159 out = paropen(onetep_file, 'r')
160 except IOError:
161 raise ReadError('Could not open output file "%s"' % onetep_file)
163 line = out.readline()
164 while line:
165 if '| Total' in line:
166 self.results['energy'] = Hartree * float(line.split()[-2])
167 elif ('Element Atom Cartesian components (Eh/a)'
168 in line):
169 self._read_forces(out)
170 elif ('Final Configuration' in line):
171 self._read_geom_output(out)
172 elif ('Integrated spin density' in line):
173 self.results['magmom'] = self._read_magmom(line)
174 elif ('|Excitation| Energy (in Ha) | Oscillator Str'
175 in line):
176 self._read_excitations(out)
177 elif ('Dipole Moment Calculation' in line):
178 self.results['dipole'] = self._read_dipole(out)
179 elif 'warn' in line.lower():
180 warnings.append(line)
181 line = out.readline()
183 if warnings:
184 warn('WARNING: %s contains warnings' % onetep_file)
185 for warning in warnings:
186 warn(warning)
188 def _read_lattice(self, out):
189 """ read the lattice parameters out of a onetep .out formatted file
190 stream"""
192 axes = []
194 l = out.readline()
195 # onetep assumes lengths are in atomic units by default
196 conv_fac = Bohr
197 if 'ang' in l:
198 l = out.readline()
199 conv_fac = 1.0
200 elif 'bohr' in l:
201 l = out.readline()
203 for _ in range(0, 3):
204 l = l.strip()
205 p = l.split()
206 if len(p) != 3:
207 raise ReadError('Malformed Lattice block line "%s"' % l)
208 try:
209 axes.append([conv_fac * float(comp) for comp in p[0:3]])
210 except ValueError:
211 raise ReadError("Can't parse line \"%s\" in axes block" % l)
212 l = out.readline()
213 self.atoms.set_cell(axes)
215 def _read_positions(self, out):
216 """Read the contents of a positions_abs block into the calculator's
217 atoms object, setting both species and positions. Tries to strip out
218 comment lines and is aware of angstom vs. bohr"""
220 line = out.readline()
221 # onetep assumes lengths are in atomic units by default
222 conv_fac = Bohr
223 if 'ang' in line:
224 line = out.readline()
225 conv_fac = 1.0
226 elif 'bohr' in line:
227 line = out.readline()
228 symbols = []
229 positions = []
230 while '%endblock' not in line.lower():
231 line = line.strip()
232 if line[0] != '#':
233 atom, suffix = line.split(None, 1)
234 pos = suffix.split(None, 3)[0:3]
235 try:
236 pos = [conv_fac * float(p) for p in pos]
237 except ValueError:
238 raise ReadError('Malformed position line "%s"', line)
239 symbols.append(atom)
240 positions.append(pos)
241 line = out.readline()
242 tags = deepcopy(symbols)
243 for j in range(len(symbols)):
244 symbols[j] = ''.join(i for i in symbols[j] if not i.isdigit())
245 for j in range(len(tags)):
246 tags[j] = ''.join(i for i in tags[j] if not i.isalpha())
247 if tags[j] == '':
248 tags[j] = '0'
249 tags[j] = int(tags[j])
250 if len(self.atoms) != len(symbols):
251 self.atoms = Atoms(symbols=symbols, positions=positions)
252 self.atoms.set_chemical_symbols(symbols)
253 self.atoms.set_tags(tags)
254 self.atoms.set_positions(positions)
256 def _read_dipole(self, out):
257 """Reads total dipole moment from ONETEP output file"""
259 # Find start of total dipole moment block
260 line = ()
261 while 'Total dipole moment' not in line:
262 line = out.readline()
264 # Read total dipole moment
265 dipolemoment = []
266 for label, pos in sorted({'dx': 6, 'dy': 2, 'dz': 2}.items()):
267 assert label in line.split()
268 value = float(line.split()[pos]) * Bohr
269 dipolemoment.append(value)
270 line = out.readline()
272 return array(dipolemoment)
274 def _read_magmom(self, line):
275 """Reads magnetic moment from Integrated Spin line"""
276 return float(line.split()[4])
278 def _read_geom_output(self, out):
279 """Reads geometry optimisation output from ONETEP output file"""
280 conv_fac = Bohr
282 # Find start of atom positions
283 while 'x-----' not in out.readline():
284 pass
285 symbols = []
286 positions = []
287 # Read atom positions
288 line = out.readline()
289 while 'xxxxxx' not in line:
290 line = line.strip()
291 pos = line.split()[3:6]
292 pos = [conv_fac * float(p) for p in pos]
293 atom = line.split()[1]
294 positions.append(pos)
295 symbols.append(atom)
296 line = out.readline()
297 if len(positions) != len(self.atoms):
298 raise ReadError('Wrong number of atoms found in output geometry'
299 'block')
300 if len(symbols) != len(self.atoms):
301 raise ReadError('Wrong number of atoms found in output geometry'
302 'block')
304 # Update atoms object with new positions (and symbols)
305 self.atoms.set_positions(positions)
306 self.atoms.set_chemical_symbols(symbols)
308 def _read_species(self, out, cond=False):
309 """ Read in species block from a onetep output file"""
310 line = out.readline().strip()
311 species = []
312 while '%endblock' not in line.lower():
313 atom, element, z, nngwf, ngwf_radius = line.split(None, 5)
314 z = int(z)
315 nngwf = int(nngwf)
316 ngwf_radius = float(ngwf_radius)
317 species.append((atom, element, z, nngwf, ngwf_radius,))
318 line = out.readline().strip()
319 if not cond:
320 self.set_species(species)
321 else:
322 self.set_species_cond(species)
324 def _read_species_pot(self, out):
325 """ Read in pseudopotential information from a onetep output file"""
326 line = out.readline().strip()
327 pots = []
328 while '%endblock' not in line.lower() and len(line) > 0:
329 atom, suffix = line.split(None, 1)
330 filename = suffix.split('#', 1)[0].strip()
331 filename = filename.replace('"', '') # take out quotes
332 filename = filename.replace("'", '')
333 pots.append((atom, filename,))
334 line = out.readline().strip()
335 if len(line) == 0:
336 raise ReadError('End of file while reading potential block')
337 self.set_pseudos(pots)
339 def _read_species_solver(self, out, cond=False):
340 """ Read in pseudopotential information from a onetep output file"""
341 line = out.readline().strip()
342 solvers = []
343 while '%endblock' not in line.lower() and len(line) > 0:
344 atom, suffix = line.split(None, 1)
345 solver_str = suffix.split('#', 1)[0].strip()
346 solvers.append((atom, solver_str))
347 line = out.readline().strip()
348 if len(line) == 0:
349 raise ReadError('End of file while reading solver block')
350 if not cond:
351 self.set_solvers(solvers)
352 else:
353 self.set_solvers_cond(solvers)
355 def _read_forces(self, out):
356 """ Extract the computed forces from a onetep output file"""
357 forces = []
358 atomic2ang = Hartree / Bohr
359 while True:
360 line = out.readline()
361 fields = line.split()
362 if len(fields) > 6:
363 break
364 while len(fields) == 7:
365 force = [float(fcomp) * atomic2ang for fcomp in fields[-4:-1]]
366 forces.append(force)
367 line = out.readline()
368 fields = line.split()
369 self.results['forces'] = array(forces)
371 def _read_excitations(self, out):
372 """ Extract the computed electronic excitations from a onetep output
373 file."""
374 excitations = []
375 line = out.readline()
376 while line:
377 words = line.split()
378 if len(words) == 0:
379 break
380 excitations.append([float(words[0]), float(words[1]) * Hartree,
381 float(words[2])])
382 line = out.readline()
383 self.results['excitations'] = array(excitations)
385 def _generate_species_block(self, cond=False):
386 """Create a default onetep species block, use -1 for the NGWF number
387 to trigger automatic NGWF number assigment using onetep's internal
388 routines."""
390 # check if we need to do anything.
391 if len(self.species) == len(self.atoms.get_chemical_symbols()):
392 return
394 parameters = self.parameters
396 atoms = self.atoms
397 if not cond:
398 self.species = []
399 default_ngwf_radius = self.parameters['ngwf_radius']
400 species_ngwf_rad_var = 'species_ngwf_radius'
401 species_ngwf_num_var = 'species_ngwf_number'
402 else:
403 self.species_cond = []
404 default_ngwf_radius = self.parameters['ngwf_radius_cond']
405 species_ngwf_rad_var = 'species_ngwf_radius_cond'
406 species_ngwf_num_var = 'species_ngwf_number_cond'
407 for sp in set(zip(atoms.get_atomic_numbers(),
408 atoms.get_chemical_symbols(),
409 ["" if i == 0 else str(i) for i in
410 atoms.get_tags()])):
411 try:
412 ngrad = parameters[species_ngwf_rad_var][sp[1]]
413 except KeyError:
414 ngrad = default_ngwf_radius
415 try:
416 ngnum = parameters[species_ngwf_num_var][sp[1]]
417 except KeyError:
418 ngnum = -1
419 if not cond:
420 self.species.append((sp[1] + sp[2], sp[1], sp[0], ngnum, ngrad))
421 else:
422 self.species_cond.append(
423 (sp[1] + sp[2], sp[1], sp[0], ngnum, ngrad))
425 def _generate_pseudo_block(self):
426 """Create a default onetep pseudopotentials block, using the
427 element name with the variable pseudo_suffix appended to it by default,
428 unless the user has set overrides for specific species by setting
429 specific entries in species_pseudo"""
431 for sp in self.species:
432 try:
433 pseudo_string = self.parameters['species_pseudo'][sp[0]]
434 except KeyError:
435 try:
436 pseudo_string = sp[1] + self.parameters['pseudo_suffix']
437 except KeyError:
438 # bare elem name if pseudo suffix empty
439 pseudo_string = sp[1]
440 self.pseudos.append((sp[0], pseudo_string))
442 def _generate_solver_block(self, cond=False):
443 """Create a default onetep pseudoatomic solvers block, using 'SOLVE'
444 unless the user has set overrides for specific species by setting
445 specific entries in species_solver (_cond)"""
447 if not cond:
448 solver_var = 'species_solver'
449 else:
450 solver_var = 'species_solver_cond'
451 for sp in self.species:
452 try:
453 atomic_string = self.parameters[solver_var][sp[0]]
454 except KeyError:
455 atomic_string = 'SOLVE'
456 if not cond:
457 self.solvers.append((sp[0], atomic_string))
458 else:
459 self.solvers_cond.append((sp[0], atomic_string))
461 def _generate_core_wf_block(self):
462 """Create a default onetep core wavefunctions block, using 'NONE'
463 unless the user has set overrides for specific species by setting
464 specific entries in species_core_wf. If all are NONE, no block
465 will be printed"""
467 any_core_wfs = False
468 for sp in self.species:
469 try:
470 core_wf_string = self.parameters['species_core_wf'][sp[0]]
471 any_core_wfs = True
472 except KeyError:
473 core_wf_string = 'NONE'
475 self.core_wfs.append((sp[0], core_wf_string))
477 # if no species core wavefunction definitions were set to anything
478 # other than 'NONE', delete the block entirely
479 if not any_core_wfs:
480 self.core_wfs = []
482 def set_pseudos(self, pots):
483 """ Sets the pseudopotential files used in this dat file """
484 self.pseudos = deepcopy(pots)
486 def set_solvers(self, solvers):
487 """ Sets the solver strings used in this dat file """
488 self.solvers = deepcopy(solvers)
490 def set_solvers_cond(self, solvers):
491 """ Sets the solver strings used in this dat file """
492 self.solvers_cond = deepcopy(solvers)
494 def set_atoms(self, atoms):
495 self.atoms = atoms
497 def set_species(self, sp):
498 """ Sets the species in the current dat instance,
499 in onetep this includes both atomic number information
500 as well as NGWF parameters like number and cut off radius"""
501 self.species = deepcopy(sp)
503 def set_species_cond(self, spc):
504 """ Sets the conduction species in the current dat instance,
505 in onetep this includes both atomic number information
506 as well as NGWF parameters like number and cut off radius"""
507 self.species_cond = deepcopy(spc)
509 def write_input(self, atoms, properties=None, system_changes=None):
510 """Only writes the input .dat file and return
511 This can be useful if one quickly needs to prepare input files
512 for a cluster where no python or ASE is available. One can than
513 upload the file manually and read out the results using
514 Onetep().read().
515 """
517 if atoms is None:
518 atoms = self.atoms
520 if self.restart:
521 self.parameters['read_tightbox_ngwfs'] = True
522 self.parameters['read_denskern'] = True
524 self._generate_species_block()
526 if len(self.pseudos) < len(self.species):
527 if 'pseudo_suffix' in self.parameters:
528 self._generate_pseudo_block()
530 if len(self.solvers) < len(self.species):
531 self._generate_solver_block()
533 if 'ngwf_radius_cond' in self.parameters:
534 if len(self.species_cond) < len(self.species):
535 self._generate_species_block(cond=True)
536 if len(self.solvers_cond) < len(self.species):
537 self._generate_solver_block(cond=True)
539 if len(self.core_wfs) < len(self.species):
540 self._generate_core_wf_block()
542 self._write_dat()
544 def get_dipole_moment(self, atoms=None):
545 self.parameters['polarisation_calculate'] = True
546 self.parameters['do_properties'] = True
547 return FileIOCalculator.get_dipole_moment(self, atoms)
549 def get_forces(self, atoms=None):
550 self.parameters['write_forces'] = True
551 return FileIOCalculator.get_forces(self, atoms)
553 def _write_dat(self, force_write=True):
554 """This export function write minimal information to
555 a .dat file. If the atoms object is a trajectory, it will
556 take the last image.
557 """
558 filename = self.label + '.dat'
560 if self.atoms is None:
561 raise Exception('No associated atoms object.')
563 atoms = self.atoms
564 parameters = self.parameters
566 if isfile(filename) and not force_write:
567 raise Exception('Target input file already exists.')
569 if 'xc' in parameters and 'xc_functional' in parameters \
570 and parameters['xc'] != parameters['xc_functional']:
571 raise Exception('Conflicting functionals defined! %s vs. %s' %
572 (parameters['xc'], parameters['xc_functional']))
574 fd = open(filename, 'w')
575 fd.write('######################################################\n')
576 fd.write('#ONETEP .dat file: %s\n' % filename)
577 fd.write('#Created using the Atomic Simulation Environment (ASE)\n')
578 fd.write('######################################################\n\n')
580 fd.write('%BLOCK LATTICE_CART\n')
581 fd.write('ang\n')
582 for line in atoms.get_cell():
583 fd.write(' %.10f %.10f %.10f\n' % tuple(line))
584 fd.write('%ENDBLOCK LATTICE_CART\n\n\n')
586 keyword = 'POSITIONS_ABS'
587 positions = atoms.get_positions()
588 tags = ["" if i == 0 else str(i) for i in atoms.get_tags()]
589 pos_block = [('%s %8.6f %8.6f %8.6f' %
590 (x + z, y[0], y[1], y[2])) for (x, y, z)
591 in zip(atoms.get_chemical_symbols(), positions, tags)]
593 fd.write('%%BLOCK %s\n' % keyword)
594 fd.write('ang\n')
595 for line in pos_block:
596 fd.write(' %s\n' % line)
597 fd.write('%%ENDBLOCK %s\n\n' % keyword)
599 keyword = 'SPECIES'
600 sp_block = [('%s %s %d %d %8.6f' % sp) for sp in self.species]
601 fd.write('%%BLOCK %s\n' % keyword)
602 for line in sorted(sp_block):
603 fd.write(' %s\n' % line)
604 fd.write('%%ENDBLOCK %s\n\n' % keyword)
606 if ((self.parameters['ngwf_radius_cond'] > 0) or
607 len(self.species_cond) == len(self.species)):
608 keyword = 'SPECIES_COND'
609 sp_block = [('%s %s %d %d %8.6f' % sp) for sp in self.species_cond]
610 fd.write('%%BLOCK %s\n' % keyword)
611 for line in sorted(sp_block):
612 fd.write(' %s\n' % line)
613 fd.write('%%ENDBLOCK %s\n\n' % keyword)
615 keyword = 'SPECIES_POT'
616 fd.write('%%BLOCK %s\n' % keyword)
617 for sp in sorted(self.pseudos):
618 fd.write(' %s "%s"\n' % (sp[0], sp[1]))
619 fd.write('%%ENDBLOCK %s\n\n' % keyword)
621 keyword = 'SPECIES_ATOMIC_SET'
622 fd.write('%%BLOCK %s\n' % keyword)
623 for sp in sorted(self.solvers):
624 fd.write(' %s "%s"\n' % (sp[0], sp[1]))
625 fd.write('%%ENDBLOCK %s\n\n' % keyword)
627 if ((self.parameters['ngwf_radius_cond'] > 0) or
628 len(self.solvers_cond) == len(self.species)):
629 keyword = 'SPECIES_ATOMIC_SET_COND'
630 fd.write('%%BLOCK %s\n' % keyword)
631 for sp in sorted(self.solvers_cond):
632 fd.write(' %s "%s"\n' % (sp[0], sp[1]))
633 fd.write('%%ENDBLOCK %s\n\n' % keyword)
635 if self.core_wfs:
636 keyword = 'SPECIES_CORE_WF'
637 fd.write('%%BLOCK %s\n' % keyword)
638 for sp in sorted(self.core_wfs):
639 fd.write(' %s "%s"\n' % (sp[0], sp[1]))
640 fd.write('%%ENDBLOCK %s\n\n' % keyword)
642 if 'bsunfld_calculate' in self.parameters:
643 if 'species_bsunfld_groups' not in self.parameters:
644 self.parameters['species_bsunfld_groups'] = \
645 self.atoms.get_chemical_symbols()
647 # Loop over parameters entries in alphabetal order, outputting
648 # them as keywords or blocks as appropriate
649 for p, param in sorted(parameters.items()):
650 if param is not None and \
651 p.lower() not in self._dummy_parameters:
652 if p.lower() in self._block_parameters:
653 keyword = p.upper()
654 fd.write('\n%%BLOCK %s\n' % keyword)
655 if p.lower() in self._path_parameters:
656 self.write_kpt_path(fd, param)
657 elif p.lower() in self._group_parameters:
658 self.write_groups(fd, param)
659 else:
660 fd.write('%s\n' % str(param))
661 fd.write('%%ENDBLOCK %s\n\n' % keyword)
662 else:
663 fd.write('%s : %s\n' % (p, param))
664 if p.upper() == 'XC':
665 # Onetep calls XC something else...
666 fd.write('xc_functional : %s\n' % param)
667 fd.close()
669 def write_kpt_path(self, fd, path):
670 """Writes a k-point path to a ONETEP input file"""
671 for kpt in array(path):
672 fd.write(' %8.6f %8.6f %8.6f\n' % (kpt[0], kpt[1], kpt[2]))
674 def write_groups(self, fd, groups):
675 """Writes multiple groups of atom labels to a ONETEP input file"""
676 for grp in groups:
677 fd.write(" ".join(map(str, grp)))
678 fd.write('\n')
680 def __repr__(self):
681 """Returns generic, fast to capture representation of
682 ONETEP settings along with atoms object.
683 """
684 expr = ''
685 expr += '-----------------Atoms--------------------\n'
686 if self.atoms is not None:
687 expr += str('%20s\n' % self.atoms)
688 else:
689 expr += 'None\n'
691 expr += '\n-----------------Species---------------------\n'
692 expr += str(self.species)
693 expr += '\n-----------------Pseudos---------------------\n'
694 expr += str(self.pseudos)
695 expr += '\n-----------------Options------------\n'
696 for key in self.parameters:
697 expr += '%20s : %s\n' % (key, self.parameters[key])
699 return expr
701 def set_label(self, label):
702 """The label is part of each seed, which in turn is a prefix
703 in each ONETEP related file.
704 """
705 self.label = label
706 self.prefix = label