root / ase / calculators / aims.py @ 16
Historique | Voir | Annoter | Télécharger (17,02 ko)
1 | 1 | tkerber | """This module defines an ASE interface to FHI-aims.
|
---|---|---|---|
2 | 1 | tkerber |
|
3 | 1 | tkerber | Felix Hanke hanke@liverpool.ac.uk
|
4 | 1 | tkerber | Jonas Bjork j.bjork@liverpool.ac.uk
|
5 | 1 | tkerber | """
|
6 | 1 | tkerber | |
7 | 1 | tkerber | from general import Calculator |
8 | 1 | tkerber | import os |
9 | 1 | tkerber | import sys |
10 | 1 | tkerber | from os.path import join, isfile, islink |
11 | 1 | tkerber | |
12 | 1 | tkerber | import numpy as np |
13 | 1 | tkerber | |
14 | 1 | tkerber | import ase |
15 | 1 | tkerber | |
16 | 1 | tkerber | float_keys = [ |
17 | 1 | tkerber | 'charge',
|
18 | 1 | tkerber | 'charge_mix_param',
|
19 | 1 | tkerber | 'hartree_convergence_parameter',
|
20 | 1 | tkerber | 'ini_linear_mix_param',
|
21 | 1 | tkerber | 'ini_spin_mix_parma',
|
22 | 1 | tkerber | 'initial_moment',
|
23 | 1 | tkerber | 'MD_MB_init',
|
24 | 1 | tkerber | 'MD_time_step',
|
25 | 1 | tkerber | 'prec_mix_param',
|
26 | 1 | tkerber | 'set_vacuum_level',
|
27 | 1 | tkerber | 'spin_mix_param',
|
28 | 1 | tkerber | ] |
29 | 1 | tkerber | |
30 | 1 | tkerber | exp_keys = [ |
31 | 1 | tkerber | 'sc_accuracy_eev',
|
32 | 1 | tkerber | 'sc_accuracy_etot',
|
33 | 1 | tkerber | 'sc_accuracy_forces',
|
34 | 1 | tkerber | 'sc_accuracy_rho',
|
35 | 1 | tkerber | ] |
36 | 1 | tkerber | |
37 | 1 | tkerber | string_keys = [ |
38 | 1 | tkerber | 'communication_type',
|
39 | 1 | tkerber | 'KS_method',
|
40 | 1 | tkerber | 'mixer',
|
41 | 1 | tkerber | 'output_level',
|
42 | 1 | tkerber | 'packed_matrix_format',
|
43 | 1 | tkerber | 'restart',
|
44 | 1 | tkerber | 'restart_read_only',
|
45 | 1 | tkerber | 'restart_write_only',
|
46 | 1 | tkerber | 'spin',
|
47 | 1 | tkerber | 'total_energy_method',
|
48 | 1 | tkerber | 'qpe_calc',
|
49 | 1 | tkerber | 'xc',
|
50 | 1 | tkerber | ] |
51 | 1 | tkerber | |
52 | 1 | tkerber | int_keys = [ |
53 | 1 | tkerber | 'empty_states',
|
54 | 1 | tkerber | 'ini_linear_mixing',
|
55 | 1 | tkerber | 'max_relaxation_steps',
|
56 | 1 | tkerber | 'n_max_pulay',
|
57 | 1 | tkerber | 'sc_iter_limit',
|
58 | 1 | tkerber | 'walltime',
|
59 | 1 | tkerber | ] |
60 | 1 | tkerber | |
61 | 1 | tkerber | bool_keys = [ |
62 | 1 | tkerber | 'collect_eigenvectors',
|
63 | 1 | tkerber | 'compute_forces',
|
64 | 1 | tkerber | 'compute_kinetic',
|
65 | 1 | tkerber | 'distributed_spline_storage',
|
66 | 1 | tkerber | 'evaluate_work_function',
|
67 | 1 | tkerber | 'final_forces_cleaned',
|
68 | 1 | tkerber | 'hessian_to_restart_geometry',
|
69 | 1 | tkerber | 'MD_clean_rotations',
|
70 | 1 | tkerber | 'restart_relaxations',
|
71 | 1 | tkerber | 'squeeze_memory',
|
72 | 1 | tkerber | 'use_density_matrix',
|
73 | 1 | tkerber | 'use_dipole_correction',
|
74 | 1 | tkerber | 'use_local_index',
|
75 | 1 | tkerber | 'vdw_correction_hirshfeld',
|
76 | 1 | tkerber | ] |
77 | 1 | tkerber | |
78 | 1 | tkerber | list_keys = [ |
79 | 1 | tkerber | 'k_grid',
|
80 | 1 | tkerber | 'k_offset',
|
81 | 1 | tkerber | 'MD_run',
|
82 | 1 | tkerber | 'MD_schedule',
|
83 | 1 | tkerber | 'MD_segment',
|
84 | 1 | tkerber | 'mixer_threshold',
|
85 | 1 | tkerber | 'occupation_type',
|
86 | 1 | tkerber | 'output',
|
87 | 1 | tkerber | 'preconditioner',
|
88 | 1 | tkerber | 'relativistic',
|
89 | 1 | tkerber | 'relax_geometry',
|
90 | 1 | tkerber | ] |
91 | 1 | tkerber | |
92 | 1 | tkerber | input_keys = [ |
93 | 1 | tkerber | 'run_command',
|
94 | 1 | tkerber | 'run_dir',
|
95 | 1 | tkerber | 'species_dir',
|
96 | 1 | tkerber | 'cubes',
|
97 | 1 | tkerber | 'output_template',
|
98 | 1 | tkerber | 'track_output',
|
99 | 1 | tkerber | ] |
100 | 1 | tkerber | |
101 | 1 | tkerber | input_parameters_default = {'run_command':None, |
102 | 1 | tkerber | 'run_dir':None, |
103 | 1 | tkerber | 'species_dir':None, |
104 | 1 | tkerber | 'cubes':None, |
105 | 1 | tkerber | 'output_template':'aims', |
106 | 1 | tkerber | 'track_output':False} |
107 | 1 | tkerber | |
108 | 1 | tkerber | class Aims(Calculator): |
109 | 1 | tkerber | def __init__(self, **kwargs): |
110 | 1 | tkerber | self.name = 'Aims' |
111 | 1 | tkerber | self.float_params = {}
|
112 | 1 | tkerber | self.exp_params = {}
|
113 | 1 | tkerber | self.string_params = {}
|
114 | 1 | tkerber | self.int_params = {}
|
115 | 1 | tkerber | self.bool_params = {}
|
116 | 1 | tkerber | self.list_params = {}
|
117 | 1 | tkerber | self.input_parameters = {}
|
118 | 1 | tkerber | for key in float_keys: |
119 | 1 | tkerber | self.float_params[key] = None |
120 | 1 | tkerber | for key in exp_keys: |
121 | 1 | tkerber | self.exp_params[key] = None |
122 | 1 | tkerber | for key in string_keys: |
123 | 1 | tkerber | self.string_params[key] = None |
124 | 1 | tkerber | for key in int_keys: |
125 | 1 | tkerber | self.int_params[key] = None |
126 | 1 | tkerber | for key in bool_keys: |
127 | 1 | tkerber | self.bool_params[key] = None |
128 | 1 | tkerber | for key in list_keys: |
129 | 1 | tkerber | self.list_params[key] = None |
130 | 1 | tkerber | for key in input_keys: |
131 | 1 | tkerber | self.input_parameters[key] = input_parameters_default[key]
|
132 | 1 | tkerber | if os.environ.has_key('AIMS_SPECIES_DIR'): |
133 | 1 | tkerber | self.input_parameters['species_dir'] = os.environ['AIMS_SPECIES_DIR'] |
134 | 1 | tkerber | if os.environ.has_key('AIMS_COMMAND'): |
135 | 1 | tkerber | self.input_parameters['run_command'] = os.environ['AIMS_COMMAND'] |
136 | 1 | tkerber | |
137 | 1 | tkerber | self.positions = None |
138 | 1 | tkerber | self.atoms = None |
139 | 1 | tkerber | self.run_counts = 0 |
140 | 1 | tkerber | self.set(**kwargs)
|
141 | 1 | tkerber | |
142 | 1 | tkerber | def set(self, **kwargs): |
143 | 1 | tkerber | if 'control' in kwargs: |
144 | 1 | tkerber | fname = kwargs['control']
|
145 | 1 | tkerber | from ase.io.aims import read_aims_calculator |
146 | 1 | tkerber | file = open(fname, 'r') |
147 | 1 | tkerber | calc_temp = None
|
148 | 1 | tkerber | while True: |
149 | 1 | tkerber | line = file.readline()
|
150 | 1 | tkerber | if "List of parameters used to initialize the calculator:" in line: |
151 | 1 | tkerber | file.readline()
|
152 | 1 | tkerber | calc_temp = read_aims_calculator(file)
|
153 | 1 | tkerber | break
|
154 | 1 | tkerber | if calc_temp is not None: |
155 | 1 | tkerber | self.float_params = calc_temp.float_params
|
156 | 1 | tkerber | self.exp_params = calc_temp.exp_params
|
157 | 1 | tkerber | self.string_params = calc_temp.string_params
|
158 | 1 | tkerber | self.int_params = calc_temp.int_params
|
159 | 1 | tkerber | self.bool_params = calc_temp.bool_params
|
160 | 1 | tkerber | self.list_params = calc_temp.list_params
|
161 | 1 | tkerber | self.input_parameters = calc_temp.input_parameters
|
162 | 1 | tkerber | else:
|
163 | 1 | tkerber | raise TypeError("Control file to be imported can not be read by ASE = " + fname) |
164 | 1 | tkerber | for key in kwargs: |
165 | 1 | tkerber | if self.float_params.has_key(key): |
166 | 1 | tkerber | self.float_params[key] = kwargs[key]
|
167 | 1 | tkerber | elif self.exp_params.has_key(key): |
168 | 1 | tkerber | self.exp_params[key] = kwargs[key]
|
169 | 1 | tkerber | elif self.string_params.has_key(key): |
170 | 1 | tkerber | self.string_params[key] = kwargs[key]
|
171 | 1 | tkerber | elif self.int_params.has_key(key): |
172 | 1 | tkerber | self.int_params[key] = kwargs[key]
|
173 | 1 | tkerber | elif self.bool_params.has_key(key): |
174 | 1 | tkerber | self.bool_params[key] = kwargs[key]
|
175 | 1 | tkerber | elif self.list_params.has_key(key): |
176 | 1 | tkerber | self.list_params[key] = kwargs[key]
|
177 | 1 | tkerber | elif self.input_parameters.has_key(key): |
178 | 1 | tkerber | self.input_parameters[key] = kwargs[key]
|
179 | 1 | tkerber | elif key is not 'control': |
180 | 1 | tkerber | raise TypeError('Parameter not defined: ' + key) |
181 | 1 | tkerber | |
182 | 1 | tkerber | def update(self, atoms): |
183 | 1 | tkerber | if self.calculation_required(atoms,[]): |
184 | 1 | tkerber | self.calculate(atoms)
|
185 | 1 | tkerber | |
186 | 1 | tkerber | def calculation_required(self, atoms,quantities): |
187 | 1 | tkerber | if (self.positions is None or |
188 | 1 | tkerber | (self.atoms != atoms) or |
189 | 1 | tkerber | (self.atoms != self.old_atoms) or |
190 | 1 | tkerber | (self.float_params != self.old_float_params) or |
191 | 1 | tkerber | (self.exp_params != self.old_exp_params) or |
192 | 1 | tkerber | (self.string_params != self.old_string_params) or |
193 | 1 | tkerber | (self.int_params != self.old_int_params) or |
194 | 1 | tkerber | (self.bool_params != self.old_bool_params) or |
195 | 1 | tkerber | (self.list_params != self.old_list_params) or |
196 | 1 | tkerber | (self.input_parameters != self.old_input_parameters)): |
197 | 1 | tkerber | return True |
198 | 1 | tkerber | else:
|
199 | 1 | tkerber | return False |
200 | 1 | tkerber | |
201 | 1 | tkerber | def calculate(self, atoms): |
202 | 1 | tkerber | """Generate necessary files in the working directory.
|
203 | 1 | tkerber |
|
204 | 1 | tkerber | If the directory does not exist it will be created.
|
205 | 1 | tkerber |
|
206 | 1 | tkerber | """
|
207 | 1 | tkerber | positions = atoms.get_positions() |
208 | 1 | tkerber | have_lattice_vectors = atoms.get_pbc().any() |
209 | 1 | tkerber | have_k_grid = self.list_params['k_grid'] |
210 | 1 | tkerber | if have_lattice_vectors and not have_k_grid: |
211 | 1 | tkerber | raise RuntimeError("Found lattice vectors but no k-grid!") |
212 | 1 | tkerber | if not have_lattice_vectors and have_k_grid: |
213 | 1 | tkerber | raise RuntimeError("Found k-grid but no lattice vectors!") |
214 | 1 | tkerber | from ase.io.aims import write_aims |
215 | 1 | tkerber | write_aims('geometry.in', atoms)
|
216 | 1 | tkerber | self.write_control()
|
217 | 1 | tkerber | self.write_species()
|
218 | 1 | tkerber | self.run()
|
219 | 1 | tkerber | self.converged = self.read_convergence() |
220 | 1 | tkerber | if not self.converged: |
221 | 1 | tkerber | os.system("tail -20 "+self.out) |
222 | 1 | tkerber | raise RuntimeError("FHI-aims did not converge!\n"+ |
223 | 1 | tkerber | "The last lines of output are printed above "+
|
224 | 1 | tkerber | "and should give an indication why.")
|
225 | 1 | tkerber | |
226 | 1 | tkerber | self.set_results(atoms)
|
227 | 1 | tkerber | |
228 | 1 | tkerber | def set_results(self,atoms): |
229 | 1 | tkerber | self.read(atoms)
|
230 | 1 | tkerber | self.old_float_params = self.float_params.copy() |
231 | 1 | tkerber | self.old_exp_params = self.exp_params.copy() |
232 | 1 | tkerber | self.old_string_params = self.string_params.copy() |
233 | 1 | tkerber | self.old_int_params = self.int_params.copy() |
234 | 1 | tkerber | self.old_input_parameters = self.input_parameters.copy() |
235 | 1 | tkerber | self.old_bool_params = self.bool_params.copy() |
236 | 1 | tkerber | self.old_list_params = self.list_params.copy() |
237 | 1 | tkerber | self.old_atoms = self.atoms.copy() |
238 | 1 | tkerber | |
239 | 1 | tkerber | def run(self): |
240 | 1 | tkerber | if self.input_parameters['track_output']: |
241 | 1 | tkerber | self.out = self.input_parameters['output_template']+str(self.run_counts)+'.out' |
242 | 1 | tkerber | self.run_counts += 1 |
243 | 1 | tkerber | else:
|
244 | 1 | tkerber | self.out = self.input_parameters['output_template']+'.out' |
245 | 1 | tkerber | if self.input_parameters['run_command']: |
246 | 1 | tkerber | aims_command = self.input_parameters['run_command'] |
247 | 1 | tkerber | elif os.environ.has_key('AIMS_COMMAND'): |
248 | 1 | tkerber | aims_command = os.environ['AIMS_COMMAND']
|
249 | 1 | tkerber | else:
|
250 | 1 | tkerber | raise RuntimeError("No specification for running FHI-aims. Aborting!") |
251 | 1 | tkerber | aims_command = aims_command + ' >> '
|
252 | 1 | tkerber | if self.input_parameters['run_dir']: |
253 | 1 | tkerber | aims_command = aims_command + self.input_parameters['run_dir'] + '/' |
254 | 1 | tkerber | aims_command = aims_command + self.out
|
255 | 1 | tkerber | self.write_parameters('#',self.out) |
256 | 1 | tkerber | exitcode = os.system(aims_command) |
257 | 1 | tkerber | if exitcode != 0: |
258 | 1 | tkerber | raise RuntimeError('FHI-aims exited with exit code: %d. ' % exitcode) |
259 | 1 | tkerber | if self.input_parameters['cubes'] and self.input_parameters['track_output']: |
260 | 1 | tkerber | self.input_parameters['cubes'].move_to_base_name(self.input_parameters['output_template']+str(self.run_counts-1)) |
261 | 1 | tkerber | |
262 | 1 | tkerber | def write_parameters(self,prefix,filename): |
263 | 1 | tkerber | output = open(filename,'w') |
264 | 1 | tkerber | output.write(prefix+'=======================================================\n')
|
265 | 1 | tkerber | output.write(prefix+'FHI-aims file: '+filename+'\n') |
266 | 1 | tkerber | output.write(prefix+'Created using the Atomic Simulation Environment (ASE)\n'+prefix+'\n') |
267 | 1 | tkerber | output.write(prefix+'List of parameters used to initialize the calculator:\n')
|
268 | 1 | tkerber | output.write(prefix+'=======================================================\n')
|
269 | 1 | tkerber | for key, val in self.float_params.items(): |
270 | 1 | tkerber | if val is not None: |
271 | 1 | tkerber | output.write('%-35s%5.6f\n' % (key, val))
|
272 | 1 | tkerber | for key, val in self.exp_params.items(): |
273 | 1 | tkerber | if val is not None: |
274 | 1 | tkerber | output.write('%-35s%5.2e\n' % (key, val))
|
275 | 1 | tkerber | for key, val in self.string_params.items(): |
276 | 1 | tkerber | if val is not None: |
277 | 1 | tkerber | output.write('%-35s%s\n' % (key, val))
|
278 | 1 | tkerber | for key, val in self.int_params.items(): |
279 | 1 | tkerber | if val is not None: |
280 | 1 | tkerber | output.write('%-35s%d\n' % (key, val))
|
281 | 1 | tkerber | for key, val in self.bool_params.items(): |
282 | 1 | tkerber | if val is not None: |
283 | 1 | tkerber | if key == 'vdw_correction_hirshfeld' and val: |
284 | 1 | tkerber | output.write('%-35s\n' % (key))
|
285 | 1 | tkerber | elif val:
|
286 | 1 | tkerber | output.write('%-35s.true.\n' % (key))
|
287 | 1 | tkerber | elif key != 'vdw_correction_hirshfeld': |
288 | 1 | tkerber | output.write('%-35s.false.\n' % (key))
|
289 | 1 | tkerber | for key, val in self.list_params.items(): |
290 | 1 | tkerber | if val is not None: |
291 | 1 | tkerber | if key == 'output': |
292 | 1 | tkerber | if not isinstance(val,(list,tuple)): |
293 | 1 | tkerber | val = [val] |
294 | 1 | tkerber | for output_type in val: |
295 | 1 | tkerber | output.write('%-35s%s\n' % (key,str(output_type))) |
296 | 1 | tkerber | else:
|
297 | 1 | tkerber | output.write('%-35s' % (key))
|
298 | 1 | tkerber | if isinstance(val,str): |
299 | 1 | tkerber | output.write(val) |
300 | 1 | tkerber | else:
|
301 | 1 | tkerber | for sub_value in val: |
302 | 1 | tkerber | output.write(str(sub_value)+' ') |
303 | 1 | tkerber | output.write('\n')
|
304 | 1 | tkerber | for key, val in self.input_parameters.items(): |
305 | 1 | tkerber | if key is 'cubes': |
306 | 1 | tkerber | if val:
|
307 | 1 | tkerber | val.write(output) |
308 | 1 | tkerber | elif val and val != input_parameters_default[key]: |
309 | 1 | tkerber | output.write(prefix+'%-34s%s\n' % (key,val))
|
310 | 1 | tkerber | output.write(prefix+'=======================================================\n\n')
|
311 | 1 | tkerber | output.close() |
312 | 1 | tkerber | |
313 | 1 | tkerber | def write_control(self, file = 'control.in'): |
314 | 1 | tkerber | """Writes the control.in file."""
|
315 | 1 | tkerber | self.write_parameters('#',file) |
316 | 1 | tkerber | |
317 | 1 | tkerber | def write_species(self, file = 'control.in'): |
318 | 1 | tkerber | from ase.data import atomic_numbers |
319 | 1 | tkerber | |
320 | 1 | tkerber | if not self.input_parameters['species_dir']: |
321 | 1 | tkerber | raise RuntimeError('Missing species directory, THIS MUST BE SPECIFIED!') |
322 | 1 | tkerber | |
323 | 1 | tkerber | control = open(file, 'a') |
324 | 1 | tkerber | species_path = self.input_parameters['species_dir'] |
325 | 1 | tkerber | symbols = self.atoms.get_chemical_symbols()
|
326 | 1 | tkerber | symbols2 = [] |
327 | 1 | tkerber | for n, symbol in enumerate(symbols): |
328 | 1 | tkerber | if symbol not in symbols2: |
329 | 1 | tkerber | symbols2.append(symbol) |
330 | 1 | tkerber | for symbol in symbols2: |
331 | 1 | tkerber | fd = join(species_path, '%02i_%s_default' % (atomic_numbers[symbol], symbol))
|
332 | 1 | tkerber | for line in open(fd, 'r'): |
333 | 1 | tkerber | control.write(line) |
334 | 1 | tkerber | control.close() |
335 | 1 | tkerber | |
336 | 1 | tkerber | def get_dipole_moment(self, atoms): |
337 | 1 | tkerber | if self.list_params['output'] is None or 'dipole' not in self.list_params['output']: |
338 | 1 | tkerber | raise RuntimeError('output=[\'dipole\'] has to be set.') |
339 | 1 | tkerber | elif atoms.get_pbc().any():
|
340 | 1 | tkerber | raise RuntimeError('FHI-aims does not allow this for systems with periodic boundary conditions.') |
341 | 1 | tkerber | self.update(atoms)
|
342 | 1 | tkerber | return self.dipole |
343 | 1 | tkerber | |
344 | 1 | tkerber | def read_dipole(self): |
345 | 1 | tkerber | """Method that reads the electric dipole moment from the output file."""
|
346 | 1 | tkerber | |
347 | 1 | tkerber | dipolemoment=np.zeros([1,3]) |
348 | 1 | tkerber | for line in open(self.out, 'r'): |
349 | 1 | tkerber | if line.rfind('Total dipole moment [eAng]') > -1: |
350 | 1 | tkerber | dipolemoment=np.array([float(f) for f in line.split()[6:10]]) |
351 | 1 | tkerber | return dipolemoment
|
352 | 1 | tkerber | |
353 | 1 | tkerber | def read_energy(self, all=None): |
354 | 1 | tkerber | for line in open(self.out, 'r'): |
355 | 1 | tkerber | if line.rfind('Total energy corrected') > -1: |
356 | 1 | tkerber | E0 = float(line.split()[-2]) |
357 | 1 | tkerber | elif line.rfind('Total energy uncorrected') > -1: |
358 | 1 | tkerber | F = float(line.split()[-2]) |
359 | 1 | tkerber | energy_free, energy_zero = F, E0 |
360 | 1 | tkerber | return [energy_free, energy_zero]
|
361 | 1 | tkerber | |
362 | 1 | tkerber | def read_forces(self, atoms, all=False): |
363 | 1 | tkerber | """Method that reads forces from the output file.
|
364 | 1 | tkerber |
|
365 | 1 | tkerber | If 'all' is switched on, the forces for all ionic steps
|
366 | 1 | tkerber | in the output file will be returned, in other case only the
|
367 | 1 | tkerber | forces for the last ionic configuration are returned."""
|
368 | 1 | tkerber | lines = open(self.out, 'r').readlines() |
369 | 1 | tkerber | forces = np.zeros([len(atoms), 3]) |
370 | 1 | tkerber | for n, line in enumerate(lines): |
371 | 1 | tkerber | if line.rfind('Total atomic forces') > -1: |
372 | 1 | tkerber | for iatom in range(len(atoms)): |
373 | 1 | tkerber | data = lines[n+iatom+1].split()
|
374 | 1 | tkerber | for iforce in range(3): |
375 | 1 | tkerber | forces[iatom, iforce] = float(data[2+iforce]) |
376 | 1 | tkerber | return forces
|
377 | 1 | tkerber | |
378 | 1 | tkerber | def get_stress(self, atoms): |
379 | 1 | tkerber | raise NotImplementedError('Stresses are not currently available in FHI-aims, sorry. ') |
380 | 1 | tkerber | |
381 | 1 | tkerber | # methods that should be quickly implemented some time, haven't had time yet:
|
382 | 1 | tkerber | def read_fermi(self): |
383 | 1 | tkerber | """Method that reads Fermi energy from output file"""
|
384 | 1 | tkerber | return
|
385 | 1 | tkerber | |
386 | 1 | tkerber | def read_magnetic_moment(self): |
387 | 1 | tkerber | return
|
388 | 1 | tkerber | |
389 | 1 | tkerber | def read_convergence(self): |
390 | 1 | tkerber | converged = False
|
391 | 1 | tkerber | lines = open(self.out, 'r').readlines() |
392 | 1 | tkerber | for n, line in enumerate(lines): |
393 | 1 | tkerber | if line.rfind('Have a nice day') > -1: |
394 | 1 | tkerber | converged = True
|
395 | 1 | tkerber | return converged
|
396 | 1 | tkerber | |
397 | 1 | tkerber | def read_eigenvalues(self, kpt=0, spin=0): |
398 | 1 | tkerber | return
|
399 | 1 | tkerber | |
400 | 1 | tkerber | class AimsCube: |
401 | 1 | tkerber | """ object to ensure the output of cube files, can be attached to Aims object"""
|
402 | 1 | tkerber | def __init__(self,origin=(0,0,0), |
403 | 1 | tkerber | edges=[(0.1,0.0,0.0),(0.0,0.1,0.0),(0.0,0.0,0.1)], |
404 | 1 | tkerber | points=(50,50,50),plots=None): |
405 | 1 | tkerber | """ parameters:
|
406 | 1 | tkerber | origin, edges, points = same as in the FHI-aims output
|
407 | 1 | tkerber | plots: what to print, same names as in FHI-aims """
|
408 | 1 | tkerber | |
409 | 1 | tkerber | self.name = 'AimsCube' |
410 | 1 | tkerber | self.origin = origin
|
411 | 1 | tkerber | self.edges = edges
|
412 | 1 | tkerber | self.points = points
|
413 | 1 | tkerber | self.plots = plots
|
414 | 1 | tkerber | |
415 | 1 | tkerber | def ncubes(self): |
416 | 1 | tkerber | """returns the number of cube files to output """
|
417 | 1 | tkerber | if self.plots: |
418 | 1 | tkerber | number = len(self.plots) |
419 | 1 | tkerber | else:
|
420 | 1 | tkerber | number = 0
|
421 | 1 | tkerber | return number
|
422 | 1 | tkerber | |
423 | 1 | tkerber | def set(self,**kwargs): |
424 | 1 | tkerber | """ set any of the parameters ... """
|
425 | 1 | tkerber | # NOT IMPLEMENTED AT THE MOMENT!
|
426 | 1 | tkerber | |
427 | 1 | tkerber | def move_to_base_name(self,basename): |
428 | 1 | tkerber | """ when output tracking is on or the base namem is not standard,
|
429 | 1 | tkerber | this routine will rename add the base to the cube file output for
|
430 | 1 | tkerber | easier tracking """
|
431 | 1 | tkerber | for plot in self.plots: |
432 | 1 | tkerber | found = False
|
433 | 1 | tkerber | cube = plot.split() |
434 | 1 | tkerber | if cube[0] == 'total_density' or cube[0] == 'spin_density' or cube[0] == 'delta_density': |
435 | 1 | tkerber | found = True
|
436 | 1 | tkerber | old_name = cube[0]+'.cube' |
437 | 1 | tkerber | new_name = basename+'.'+old_name
|
438 | 1 | tkerber | if cube[0] == 'eigenstate' or cube[0] == 'eigenstate_density': |
439 | 1 | tkerber | found = True
|
440 | 1 | tkerber | state = int(cube[1]) |
441 | 1 | tkerber | s_state = cube[1]
|
442 | 1 | tkerber | for i in [10,100,1000,10000]: |
443 | 1 | tkerber | if state < i:
|
444 | 1 | tkerber | s_state = '0'+s_state
|
445 | 1 | tkerber | old_name = cube[0]+'_'+s_state+'_spin_1.cube' |
446 | 1 | tkerber | new_name = basename+'.'+old_name
|
447 | 1 | tkerber | if found:
|
448 | 1 | tkerber | os.system("mv "+old_name+" "+new_name) |
449 | 1 | tkerber | |
450 | 1 | tkerber | def add_plot(self,name): |
451 | 1 | tkerber | """ in case you forgot one ... """
|
452 | 1 | tkerber | plots += [name] |
453 | 1 | tkerber | |
454 | 1 | tkerber | def write(self,file): |
455 | 1 | tkerber | """ write the necessary output to the already opened control.in """
|
456 | 1 | tkerber | file.write('output cube '+self.plots[0]+'\n') |
457 | 1 | tkerber | file.write(' cube origin ') |
458 | 1 | tkerber | for ival in self.origin: |
459 | 1 | tkerber | file.write(str(ival)+' ') |
460 | 1 | tkerber | file.write('\n') |
461 | 1 | tkerber | for i in range(3): |
462 | 1 | tkerber | file.write(' cube edge '+str(self.points[i])+' ') |
463 | 1 | tkerber | for ival in self.edges[i]: |
464 | 1 | tkerber | file.write(str(ival)+' ') |
465 | 1 | tkerber | file.write('\n') |
466 | 1 | tkerber | if self.ncubes() > 1: |
467 | 1 | tkerber | for i in range(self.ncubes()-1): |
468 | 1 | tkerber | file.write('output cube '+self.plots[i+1]+'\n') |
469 | 1 | tkerber | |
470 | 1 | tkerber |