1# type: ignore 2""" 3This module defines an ASE interface to Turbomole: http://www.turbomole.com/ 4 5QMMM functionality provided by Markus Kaukonen <markus.kaukonen@iki.fi>. 6 7Please read the license file (../../LICENSE) 8 9Contact: Ivan Kondov <ivan.kondov@kit.edu> 10""" 11import os 12import re 13import warnings 14from subprocess import Popen, PIPE 15from math import log10, floor 16import numpy as np 17from ase import Atoms 18from ase.units import Ha, Bohr 19from ase.io import read, write 20from ase.calculators.calculator import FileIOCalculator 21from ase.calculators.calculator import PropertyNotImplementedError, ReadError 22 23 24def read_output(regex): 25 """collects all matching strings from the output""" 26 hitlist = [] 27 checkfiles = [] 28 for filename in os.listdir('.'): 29 if filename.startswith('job.') or filename.endswith('.out'): 30 checkfiles.append(filename) 31 for filename in checkfiles: 32 with open(filename, 'rt') as fd: 33 lines = fd.readlines() 34 for line in lines: 35 match = re.search(regex, line) 36 if match: 37 hitlist.append(match.group(1)) 38 return hitlist 39 40 41def execute(args, input_str=None, error_test=True, 42 stdout_tofile=True): 43 """executes a turbomole executable and process the outputs""" 44 45 if isinstance(args, str): 46 args = args.split() 47 48 if stdout_tofile: 49 stdout_file = 'ASE.TM.' + args[0] + '.out' 50 stdout = open(stdout_file, 'w') 51 else: 52 stdout = PIPE 53 54 if input_str: 55 stdin = input_str.encode() 56 else: 57 stdin = None 58 59 message = 'TM command "' + args[0] + '" execution failed' 60 try: 61 proc = Popen(args, stdin=PIPE, stderr=PIPE, stdout=stdout) 62 res = proc.communicate(input=stdin) 63 if error_test: 64 error = res[1].decode() 65 if 'abnormally' in error or 'ended normally' not in error: 66 message += ' with error:\n' + error 67 message += '\nSee file ' + stdout_file + ' for details.\n' 68 raise RuntimeError(message) 69 except RuntimeError as err: 70 raise err 71 except OSError as err: 72 raise OSError(err.args[1] + '\n' + message) 73 else: 74 print('TM command: "' + args[0] + '" successfully executed') 75 76 if not stdout_tofile: 77 return res[0].decode() 78 79 80def add_data_group(data_group, string=None, raw=False): 81 """write a turbomole data group to control file""" 82 if raw: 83 data = data_group 84 else: 85 data = '$' + data_group 86 if string: 87 data += ' ' + string 88 data += '\n' 89 fd = open('control', 'r+') 90 lines = fd.readlines() 91 fd.seek(0) 92 fd.truncate() 93 lines.insert(2, data) 94 fd.write(''.join(lines)) 95 fd.close() 96 97 98def read_data_group(data_group): 99 """read a turbomole data group from control file""" 100 args = ['sdg', data_group] 101 dg = execute(args, error_test=False, stdout_tofile=False) 102 return dg.strip() 103 104 105def delete_data_group(data_group): 106 """delete a turbomole data group from control file""" 107 command = ['kdg', data_group] 108 execute(command, error_test=False, stdout_tofile=False) 109 110 111class TurbomoleOptimizer: 112 def __init__(self, atoms, calc): 113 self.atoms = atoms 114 self.calc = calc 115 self.atoms.calc = self.calc 116 117 def todict(self): 118 return {'type': 'optimization', 119 'optimizer': 'TurbomoleOptimizer'} 120 121 def run(self, fmax=None, steps=None): 122 if fmax is not None: 123 self.calc.parameters['force convergence'] = fmax 124 self.calc.verify_parameters() 125 if steps is not None: 126 self.calc.parameters['geometry optimization iterations'] = steps 127 self.calc.verify_parameters() 128 self.calc.calculate() 129 self.atoms.positions[:] = self.calc.atoms.positions 130 self.calc.parameters['task'] = 'energy' 131 132 133class Turbomole(FileIOCalculator): 134 135 """constants""" 136 name = 'Turbomole' 137 138 implemented_properties = ['energy', 'forces', 'dipole', 'free_energy', 139 'charges'] 140 141 available_functionals = [ 142 'slater-dirac-exchange', 's-vwn', 'vwn', 's-vwn_Gaussian', 'pwlda', 143 'becke-exchange', 'b-lyp', 'b-vwn', 'lyp', 'b-p', 'pbe', 'tpss', 144 'bh-lyp', 'b3-lyp', 'b3-lyp_Gaussian', 'pbe0', 'tpssh', 'lhf', 'oep', 145 'b97-d', 'b2-plyp' 146 ] 147 tm_files = [ 148 'control', 'coord', 'basis', 'auxbasis', 'energy', 'gradient', 'mos', 149 'alpha', 'beta', 'statistics', 'GEO_OPT_CONVERGED', 'GEO_OPT_FAILED', 150 'not.converged', 'nextstep', 'hessapprox', 'job.last', 'job.start', 151 'optinfo', 'statistics', 'converged', 'vibspectrum', 152 'vib_normal_modes', 'hessian', 'dipgrad', 'dscf_problem', 'pc.txt', 153 'pc_gradients.txt' 154 ] 155 tm_tmp_files = [ 156 'errvec', 'fock', 'oldfock', 'dens', 'ddens', 'diff_densmat', 157 'diff_dft_density', 'diff_dft_oper', 'diff_fockmat', 'diis_errvec', 158 'diis_oldfock' 159 ] 160 spec_names = { 161 'default': 'default_parameters', 162 'comment': 'parameter_comment', 163 'updateable': 'parameter_updateable', 164 'type': 'parameter_type', 165 'key': 'parameter_key', 166 'group': 'parameter_group', 167 'units': 'parameter_units', 168 'mapping': 'parameter_mapping', 169 'non-define': 'parameter_no_define' 170 } 171 172 # flat dictionaries with parameters attributes 173 default_parameters = {} 174 parameter_comment = {} 175 parameter_updateable = {} 176 parameter_type = {} 177 parameter_key = {} 178 parameter_group = {} 179 parameter_units = {} 180 parameter_mapping = {} 181 parameter_no_define = {} 182 183 # nested dictionary with parameters attributes 184 parameter_spec = { 185 'automatic orbital shift': { 186 'comment': None, 187 'default': 0.1, 188 'group': 'scforbitalshift', 189 'key': 'automatic', 190 'mapping': { 191 'to_control': lambda a: a / Ha, 192 'from_control': lambda a: a * Ha 193 }, 194 'type': float, 195 'units': 'eV', 196 'updateable': True 197 }, 198 'basis set definition': { 199 'comment': 'used only in restart', 200 'default': None, 201 'group': 'basis', 202 'key': None, 203 'type': dict, 204 'units': None, 205 'updateable': False 206 }, 207 'basis set name': { 208 'comment': 'current default from module "define"', 209 'default': 'def-SV(P)', 210 'group': 'basis', 211 'key': None, 212 'type': str, 213 'units': None, 214 'updateable': False 215 }, 216 'closed-shell orbital shift': { 217 'comment': 'does not work with automatic', 218 'default': None, 219 'group': 'scforbitalshift', 220 'key': 'closedshell', 221 'mapping': { 222 'to_control': lambda a: a / Ha, 223 'from_control': lambda a: a * Ha 224 }, 225 'type': float, 226 'units': 'eV', 227 'updateable': True 228 }, 229 'damping adjustment step': { 230 'comment': None, 231 'default': None, 232 'group': 'scfdamp', 233 'key': 'step', 234 'type': float, 235 'units': None, 236 'updateable': True 237 }, 238 'density convergence': { 239 'comment': None, 240 'default': None, 241 'group': 'denconv', 242 'key': 'denconv', 243 'mapping': { 244 'to_control': lambda a: int(-log10(a)), 245 'from_control': lambda a: 10**(-a) 246 }, 247 'non-define': True, 248 'type': float, 249 'units': None, 250 'updateable': True 251 }, 252 'density functional': { 253 'comment': None, 254 'default': 'b-p', 255 'group': 'dft', 256 'key': 'functional', 257 'type': str, 258 'units': None, 259 'updateable': True 260 }, 261 'energy convergence': { 262 'comment': 'jobex -energy <int>', 263 'default': None, 264 'group': None, 265 'key': None, 266 'mapping': { 267 'to_control': lambda a: a / Ha, 268 'from_control': lambda a: a * Ha 269 }, 270 'type': float, 271 'units': 'eV', 272 'updateable': True 273 }, 274 'fermi annealing factor': { 275 'comment': None, 276 'default': 0.95, 277 'group': 'fermi', 278 'key': 'tmfac', 279 'type': float, 280 'units': None, 281 'updateable': True 282 }, 283 'fermi final temperature': { 284 'comment': None, 285 'default': 300, 286 'group': 'fermi', 287 'key': 'tmend', 288 'type': float, 289 'units': 'Kelvin', 290 'updateable': True 291 }, 292 'fermi homo-lumo gap criterion': { 293 'comment': None, 294 'default': 0.1, 295 'group': 'fermi', 296 'key': 'hlcrt', 297 'mapping': { 298 'to_control': lambda a: a / Ha, 299 'from_control': lambda a: a * Ha 300 }, 301 'type': float, 302 'units': 'eV', 303 'updateable': True 304 }, 305 'fermi initial temperature': { 306 'comment': None, 307 'default': 300, 308 'group': 'fermi', 309 'key': 'tmstrt', 310 'type': float, 311 'units': 'Kelvin', 312 'updateable': True 313 }, 314 'fermi stopping criterion': { 315 'comment': None, 316 'default': 0.001, 317 'group': 'fermi', 318 'key': 'stop', 319 'mapping': { 320 'to_control': lambda a: a / Ha, 321 'from_control': lambda a: a * Ha 322 }, 323 'type': float, 324 'units': 'eV', 325 'updateable': True 326 }, 327 'force convergence': { 328 'comment': 'jobex -gcart <int>', 329 'default': None, 330 'group': None, 331 'key': None, 332 'mapping': { 333 'to_control': lambda a: a / Ha * Bohr, 334 'from_control': lambda a: a * Ha / Bohr 335 }, 336 'type': float, 337 'units': 'eV/Angstrom', 338 'updateable': True 339 }, 340 'geometry optimization iterations': { 341 'comment': 'jobex -c <int>', 342 'default': None, 343 'group': None, 344 'key': None, 345 'type': int, 346 'units': None, 347 'updateable': True 348 }, 349 'grid size': { 350 'comment': None, 351 'default': 'm3', 352 'group': 'dft', 353 'key': 'gridsize', 354 'type': str, 355 'units': None, 356 'updateable': True 357 }, 358 'ground state': { 359 'comment': 'only this is currently supported', 360 'default': True, 361 'group': None, 362 'key': None, 363 'type': bool, 364 'units': None, 365 'updateable': False 366 }, 367 'initial damping': { 368 'comment': None, 369 'default': None, 370 'group': 'scfdamp', 371 'key': 'start', 372 'type': float, 373 'units': None, 374 'updateable': True 375 }, 376 'initial guess': { 377 'comment': '"eht", "hcore" or {"use": "<path/to/control>"}', 378 'default': 'eht', 379 'group': None, 380 'key': None, 381 'type': None, 382 'units': None, 383 'updateable': False 384 }, 385 'minimal damping': { 386 'comment': None, 387 'default': None, 388 'group': 'scfdamp', 389 'key': 'min', 390 'type': float, 391 'units': None, 392 'updateable': True 393 }, 394 'multiplicity': { 395 'comment': None, 396 'default': None, 397 'group': None, 398 'key': None, 399 'type': int, 400 'units': None, 401 'updateable': False 402 }, 403 'non-automatic orbital shift': { 404 'comment': None, 405 'default': False, 406 'group': 'scforbitalshift', 407 'key': 'noautomatic', 408 'type': bool, 409 'units': None, 410 'updateable': True 411 }, 412 'point group': { 413 'comment': 'only c1 supported', 414 'default': 'c1', 415 'group': 'symmetry', 416 'key': 'symmetry', 417 'type': str, 418 'units': None, 419 'updateable': False 420 }, 421 'ri memory': { 422 'comment': None, 423 'default': 1000, 424 'group': 'ricore', 425 'key': 'ricore', 426 'type': int, 427 'units': 'Megabyte', 428 'updateable': True 429 }, 430 'rohf': { 431 'comment': 'used only in restart', 432 'default': None, 433 'group': None, 434 'key': None, 435 'type': bool, 436 'units': None, 437 'updateable': False 438 }, 439 'scf energy convergence': { 440 'comment': None, 441 'default': None, 442 'group': 'scfconv', 443 'key': 'scfconv', 444 'mapping': { 445 'to_control': lambda a: int(floor(-log10(a / Ha))), 446 'from_control': lambda a: 10**(-a) * Ha 447 }, 448 'type': float, 449 'units': 'eV', 450 'updateable': True 451 }, 452 'scf iterations': { 453 'comment': None, 454 'default': 60, 455 'group': 'scfiterlimit', 456 'key': 'scfiterlimit', 457 'type': int, 458 'units': None, 459 'updateable': True 460 }, 461 'task': { 462 'comment': '"energy calculation" = "energy", ' 463 '"gradient calculation" = "gradient", ' 464 '"geometry optimization" = "optimize", ' 465 '"normal mode analysis" = "frequencies"', 466 'default': 'energy', 467 'group': None, 468 'key': None, 469 'type': str, 470 'units': None, 471 'updateable': True 472 }, 473 'title': { 474 'comment': None, 475 'default': '', 476 'group': 'title', 477 'key': 'title', 478 'type': str, 479 'units': None, 480 'updateable': False 481 }, 482 'total charge': { 483 'comment': None, 484 'default': 0, 485 'group': None, 486 'key': None, 487 'type': int, 488 'units': None, 489 'updateable': False 490 }, 491 'uhf': { 492 'comment': None, 493 'default': None, 494 'group': 'uhf', 495 'key': 'uhf', 496 'type': bool, 497 'units': None, 498 'updateable': False 499 }, 500 'use basis set library': { 501 'comment': 'only true implemented', 502 'default': True, 503 'group': 'basis', 504 'key': None, 505 'type': bool, 506 'units': None, 507 'updateable': False 508 }, 509 'use dft': { 510 'comment': None, 511 'default': True, 512 'group': 'dft', 513 'key': 'dft', 514 'type': bool, 515 'units': None, 516 'updateable': False 517 }, 518 'use fermi smearing': { 519 'comment': None, 520 'default': False, 521 'group': 'fermi', 522 'key': 'fermi', 523 'type': bool, 524 'units': None, 525 'updateable': True 526 }, 527 'use redundant internals': { 528 'comment': None, 529 'default': False, 530 'group': 'redundant', 531 'key': None, 532 'type': bool, 533 'units': None, 534 'updateable': False 535 }, 536 'use resolution of identity': { 537 'comment': None, 538 'default': False, 539 'group': 'rij', 540 'key': 'rij', 541 'type': bool, 542 'units': None, 543 'updateable': False 544 }, 545 'numerical hessian': { 546 'comment': 'NumForce will be used if dictionary exists', 547 'default': None, 548 'group': None, 549 'key': None, 550 'type': dict, 551 'units': None, 552 'updateable': True 553 }, 554 'esp fit': { 555 'comment': 'ESP fit', 556 'default': None, 557 'group': 'esp_fit', 558 'key': 'esp_fit', 559 'type': str, 560 'units': None, 561 'updateable': True, 562 'non-define': True 563 } 564 } 565 566 # initialize attributes 567 parameters = {} 568 results = {} 569 initialized = False 570 pc_initialized = False 571 converged = False 572 updated = False 573 update_energy = None 574 update_forces = None 575 update_geometry = None 576 update_hessian = None 577 atoms = None 578 forces = None 579 e_total = None 580 dipole = None 581 charges = None 582 version = None 583 runtime = None 584 datetime = None 585 hostname = None 586 pcpot = None 587 588 def __init__(self, label=None, calculate_energy='dscf', 589 calculate_forces='grad', post_HF=False, atoms=None, 590 restart=False, define_str=None, control_kdg=None, 591 control_input=None, reset_tolerance=1e-2, **kwargs): 592 593 FileIOCalculator.__init__(self) 594 595 self.label = label 596 self.calculate_energy = calculate_energy 597 self.calculate_forces = calculate_forces 598 self.post_HF = post_HF 599 self.restart = restart 600 self.define_str = define_str 601 self.control_kdg = control_kdg 602 self.control_input = control_input 603 self.reset_tolerance = reset_tolerance 604 605 # construct flat dictionaries with parameter attributes 606 for p in self.parameter_spec: 607 for k in self.spec_names: 608 if k in list(self.parameter_spec[p].keys()): 609 subdict = getattr(self, self.spec_names[k]) 610 subdict.update({p: self.parameter_spec[p][k]}) 611 612 if self.restart: 613 self._set_restart(kwargs) 614 else: 615 self.set_parameters(kwargs) 616 self.verify_parameters() 617 self.reset() 618 619 if atoms is not None: 620 atoms.calc = self 621 self.set_atoms(atoms) 622 623 def __getitem__(self, item): 624 return getattr(self, item) 625 626 def _set_restart(self, params_update): 627 """constructs atoms, parameters and results from a previous 628 calculation""" 629 630 # read results, key parameters and non-key parameters 631 self.read_restart() 632 params_old = self.read_parameters() 633 634 # filter out non-updateable parameters 635 for p in list(params_update.keys()): 636 if not self.parameter_updateable[p]: 637 del params_update[p] 638 warnings.warn('"' + p + '"' + ' cannot be changed') 639 640 # update and verify parameters 641 params_new = params_old.copy() 642 params_new.update(params_update) 643 self.set_parameters(params_new) 644 self.verify_parameters() 645 646 # if a define string is specified then run define 647 if self.define_str: 648 execute('define', input_str=self.define_str) 649 650 # updates data groups in the control file 651 if params_update or self.control_kdg or self.control_input: 652 self._update_data_groups(params_old, params_update) 653 654 self.initialized = True 655 # more precise convergence tests are necessary to set these flags: 656 self.update_energy = True 657 self.update_forces = True 658 self.update_geometry = True 659 self.update_hessian = True 660 661 def _update_data_groups(self, params_old, params_update): 662 """updates data groups in the control file""" 663 # construct a list of data groups to update 664 grps = [] 665 for p in list(params_update.keys()): 666 if self.parameter_group[p] is not None: 667 grps.append(self.parameter_group[p]) 668 669 # construct a dictionary of data groups and update params 670 dgs = {} 671 for g in grps: 672 dgs[g] = {} 673 for p in self.parameter_key: 674 if g == self.parameter_group[p]: 675 if self.parameter_group[p] == self.parameter_key[p]: 676 if p in list(params_update.keys()): 677 val = params_update[p] 678 pmap = list(self.parameter_mapping.keys()) 679 if val is not None and p in pmap: 680 fun = self.parameter_mapping[p]['to_control'] 681 val = fun(params_update[p]) 682 dgs[g] = val 683 else: 684 if p in list(params_old.keys()): 685 val = params_old[p] 686 pmap = list(self.parameter_mapping.keys()) 687 if val is not None and p in pmap: 688 fun = self.parameter_mapping[p]['to_control'] 689 val = fun(params_old[p]) 690 dgs[g][self.parameter_key[p]] = val 691 if p in list(params_update.keys()): 692 val = params_update[p] 693 pmap = list(self.parameter_mapping.keys()) 694 if val is not None and p in pmap: 695 fun = self.parameter_mapping[p]['to_control'] 696 val = fun(params_update[p]) 697 dgs[g][self.parameter_key[p]] = val 698 699 # write dgs dictionary to a data group 700 for g in dgs: 701 delete_data_group(g) 702 if isinstance(dgs[g], dict): 703 string = '' 704 for key in list(dgs[g].keys()): 705 if dgs[g][key] is None: 706 continue 707 elif isinstance(dgs[g][key], bool): 708 if dgs[g][key]: 709 string += ' ' + key 710 else: 711 string += ' ' + key + '=' + str(dgs[g][key]) 712 add_data_group(g, string=string) 713 else: 714 if isinstance(dgs[g], bool): 715 if dgs[g]: 716 add_data_group(g, string='') 717 else: 718 add_data_group(g, string=str(dgs[g])) 719 720 self._set_post_define() 721 722 def _set_post_define(self): 723 """non-define keys, user-specified changes in the control file""" 724 # process key parameters that are not written with define 725 for p in list(self.parameters.keys()): 726 if p in list(self.parameter_no_define.keys()): 727 if self.parameter_no_define[p]: 728 if self.parameters[p]: 729 if p in list(self.parameter_mapping.keys()): 730 fun = self.parameter_mapping[p]['to_control'] 731 val = fun(self.parameters[p]) 732 else: 733 val = self.parameters[p] 734 delete_data_group(self.parameter_group[p]) 735 add_data_group(self.parameter_group[p], str(val)) 736 else: 737 delete_data_group(self.parameter_group[p]) 738 739 # delete user-specified data groups 740 if self.control_kdg: 741 for dg in self.control_kdg: 742 delete_data_group(dg) 743 744 # append user-defined input to control 745 if self.control_input: 746 for inp in self.control_input: 747 add_data_group(inp, raw=True) 748 749 # add point charges if pcpot defined: 750 if self.pcpot: 751 self.set_point_charges() 752 753 def set_parameters(self, params): 754 """loads the default parameters and updates with actual values""" 755 self.parameters = self.default_parameters.copy() 756 self.parameters.update(params) 757 if self.parameters['use resolution of identity']: 758 self.calculate_energy = 'ridft' 759 self.calculate_forces = 'rdgrad' 760 761 def verify_parameters(self): 762 """detect wrong or not implemented parameters""" 763 764 # kwargs parameters are ignored if user provides define_str 765 if self.define_str is not None: 766 assert isinstance(self.define_str, str) 767 assert len(self.define_str) != 0 768 return 769 770 for par in self.parameters: 771 assert par in self.parameter_spec, 'invalid parameter: ' + par 772 773 if self.parameters['use dft']: 774 func_list = [x.lower() for x in self.available_functionals] 775 func = self.parameters['density functional'] 776 assert func.lower() in func_list, ( 777 'density functional not available / not supported' 778 ) 779 780 assert self.parameters['multiplicity'], 'multiplicity not defined' 781 782 if self.parameters['rohf']: 783 raise NotImplementedError('ROHF not implemented') 784 if self.parameters['initial guess'] not in ['eht', 'hcore']: 785 if not (isinstance(self.parameters['initial guess'], dict) and 786 'use' in self.parameters['initial guess'].keys()): 787 raise ValueError('Wrong input for initial guess') 788 if not self.parameters['use basis set library']: 789 raise NotImplementedError('Explicit basis set definition') 790 if self.parameters['point group'] != 'c1': 791 raise NotImplementedError('Point group not impemeneted') 792 793 def reset(self): 794 """removes all turbomole input, output and scratch files, 795 and deletes results dict and the atoms object""" 796 self.atoms = None 797 self.results = {} 798 self.results['calculation parameters'] = {} 799 ase_files = [f for f in os.listdir('.') if f.startswith('ASE.TM.')] 800 for f in self.tm_files + self.tm_tmp_files + ase_files: 801 if os.path.exists(f): 802 os.remove(f) 803 self.initialized = False 804 self.pc_initialized = False 805 self.converged = False 806 807 def set_atoms(self, atoms): 808 """Create the self.atoms object and writes the coord file. If 809 self.atoms exists a check for changes and an update of the atoms 810 is performed. Note: Only positions changes are tracked in this 811 version. 812 """ 813 changes = self.check_state(atoms, tol=1e-13) 814 if self.atoms == atoms or 'positions' not in changes: 815 # print('two atoms obj are (almost) equal') 816 if self.updated and os.path.isfile('coord'): 817 self.updated = False 818 a = read('coord').get_positions() 819 if np.allclose(a, atoms.get_positions(), rtol=0, atol=1e-13): 820 return 821 else: 822 return 823 824 changes = self.check_state(atoms, tol=self.reset_tolerance) 825 if 'positions' in changes: 826 # print(two atoms obj are different') 827 self.reset() 828 else: 829 # print('two atoms obj are slightly different') 830 if self.parameters['use redundant internals']: 831 self.reset() 832 833 write('coord', atoms) 834 self.atoms = atoms.copy() 835 self.update_energy = True 836 self.update_forces = True 837 self.update_geometry = True 838 self.update_hessian = True 839 840 def get_define_str(self): 841 """construct a define string from the parameters dictionary""" 842 define_str_tpl = ( 843 '\n__title__\na coord\n__inter__\n' 844 'bb all __basis_set__\n*\neht\ny\n__charge_str____occ_str__' 845 '__single_atom_str____norb_str____dft_str____ri_str__' 846 '__scfiterlimit____fermi_str____damp_str__q\n' 847 ) 848 849 params = self.parameters 850 851 if params['use redundant internals']: 852 internals_str = 'ired\n*' 853 else: 854 internals_str = '*\nno' 855 charge_str = str(params['total charge']) + '\n' 856 857 if params['multiplicity'] == 1: 858 if params['uhf']: 859 occ_str = 'n\ns\n*\n' 860 else: 861 occ_str = 'y\n' 862 elif params['multiplicity'] == 2: 863 occ_str = 'y\n' 864 elif params['multiplicity'] == 3: 865 occ_str = 'n\nt\n*\n' 866 else: 867 unpaired = params['multiplicity'] - 1 868 if params['use fermi smearing']: 869 occ_str = 'n\nuf ' + str(unpaired) + '\n*\n' 870 else: 871 occ_str = 'n\nu ' + str(unpaired) + '\n*\n' 872 873 if len(self.atoms) != 1: 874 single_atom_str = '' 875 else: 876 single_atom_str = '\n' 877 878 if params['multiplicity'] == 1 and not params['uhf']: 879 norb_str = '' 880 else: 881 norb_str = 'n\n' 882 883 if params['use dft']: 884 dft_str = 'dft\non\n*\n' 885 else: 886 dft_str = '' 887 888 if params['density functional']: 889 dft_str += 'dft\nfunc ' + params['density functional'] + '\n*\n' 890 891 if params['grid size']: 892 dft_str += 'dft\ngrid ' + params['grid size'] + '\n*\n' 893 894 if params['use resolution of identity']: 895 ri_str = 'ri\non\nm ' + str(params['ri memory']) + '\n*\n' 896 else: 897 ri_str = '' 898 899 if params['scf iterations']: 900 scfmaxiter = params['scf iterations'] 901 scfiter_str = 'scf\niter\n' + str(scfmaxiter) + '\n\n' 902 else: 903 scfiter_str = '' 904 if params['scf energy convergence']: 905 conv = floor(-log10(params['scf energy convergence'] / Ha)) 906 scfiter_str += 'scf\nconv\n' + str(int(conv)) + '\n\n' 907 908 fermi_str = '' 909 if params['use fermi smearing']: 910 fermi_str = 'scf\nfermi\n' 911 if params['fermi initial temperature']: 912 par = str(params['fermi initial temperature']) 913 fermi_str += '1\n' + par + '\n' 914 if params['fermi final temperature']: 915 par = str(params['fermi final temperature']) 916 fermi_str += '2\n' + par + '\n' 917 if params['fermi annealing factor']: 918 par = str(params['fermi annealing factor']) 919 fermi_str += '3\n' + par + '\n' 920 if params['fermi homo-lumo gap criterion']: 921 par = str(params['fermi homo-lumo gap criterion']) 922 fermi_str += '4\n' + par + '\n' 923 if params['fermi stopping criterion']: 924 par = str(params['fermi stopping criterion']) 925 fermi_str += '5\n' + par + '\n' 926 fermi_str += '\n\n' 927 928 damp_str = '' 929 damp_keys = ('initial damping', 'damping adjustment step', 930 'minimal damping') 931 damp_pars = [params[k] for k in damp_keys] 932 if any(damp_pars): 933 damp_str = 'scf\ndamp\n' 934 for par in damp_pars: 935 par_str = str(par) if par else '' 936 damp_str += par_str + '\n' 937 damp_str += '\n' 938 939 define_str = define_str_tpl 940 define_str = re.sub('__title__', params['title'], define_str) 941 define_str = re.sub('__basis_set__', params['basis set name'], 942 define_str) 943 define_str = re.sub('__charge_str__', charge_str, define_str) 944 define_str = re.sub('__occ_str__', occ_str, define_str) 945 define_str = re.sub('__norb_str__', norb_str, define_str) 946 define_str = re.sub('__dft_str__', dft_str, define_str) 947 define_str = re.sub('__ri_str__', ri_str, define_str) 948 define_str = re.sub('__single_atom_str__', single_atom_str, 949 define_str) 950 define_str = re.sub('__inter__', internals_str, define_str) 951 define_str = re.sub('__scfiterlimit__', scfiter_str, define_str) 952 define_str = re.sub('__fermi_str__', fermi_str, define_str) 953 define_str = re.sub('__damp_str__', damp_str, define_str) 954 955 return define_str 956 957 def initialize(self): 958 """prepare turbomole control file by running module 'define'""" 959 if self.initialized: 960 return 961 self.verify_parameters() 962 if not self.atoms: 963 raise RuntimeError('atoms missing during initialization') 964 if not os.path.isfile('coord'): 965 raise IOError('file coord not found') 966 967 if self.define_str is not None: 968 define_str = self.define_str 969 else: 970 define_str = self.get_define_str() 971 972 # run define 973 execute('define', input_str=define_str) 974 975 # process non-default initial guess 976 iguess = self.parameters['initial guess'] 977 if isinstance(iguess, dict) and 'use' in iguess.keys(): 978 # "use" initial guess 979 if self.parameters['multiplicity'] != 1 or self.parameters['uhf']: 980 define_str = '\n\n\ny\nuse ' + iguess['use'] + '\nn\nn\nq\n' 981 else: 982 define_str = '\n\n\ny\nuse ' + iguess['use'] + '\nn\nq\n' 983 execute('define', input_str=define_str) 984 elif self.parameters['initial guess'] == 'hcore': 985 # "hcore" initial guess 986 if self.parameters['multiplicity'] != 1 or self.parameters['uhf']: 987 delete_data_group('uhfmo_alpha') 988 delete_data_group('uhfmo_beta') 989 add_data_group('uhfmo_alpha', 'none file=alpha') 990 add_data_group('uhfmo_beta', 'none file=beta') 991 else: 992 delete_data_group('scfmo') 993 add_data_group('scfmo', 'none file=mos') 994 995 self._set_post_define() 996 997 self.initialized = True 998 self.converged = False 999 1000 def calculation_required(self, atoms, properties): 1001 if self.atoms != atoms: 1002 return True 1003 for prop in properties: 1004 if prop == 'energy' and self.e_total is None: 1005 return True 1006 elif prop == 'forces' and self.forces is None: 1007 return True 1008 return False 1009 1010 def calculate(self, atoms=None): 1011 """execute the requested job""" 1012 if atoms is None: 1013 atoms = self.atoms 1014 if self.parameters['task'] in ['energy', 'energy calculation']: 1015 self.get_potential_energy(atoms) 1016 if self.parameters['task'] in ['gradient', 'gradient calculation']: 1017 self.get_forces(atoms) 1018 if self.parameters['task'] in ['optimize', 'geometry optimization']: 1019 self.relax_geometry(atoms) 1020 if self.parameters['task'] in ['frequencies', 'normal mode analysis']: 1021 self.normal_mode_analysis(atoms) 1022 self.read_results() 1023 1024 def relax_geometry(self, atoms=None): 1025 """execute geometry optimization with script jobex""" 1026 if atoms is None: 1027 atoms = self.atoms 1028 self.set_atoms(atoms) 1029 if self.converged and not self.update_geometry: 1030 return 1031 self.initialize() 1032 jobex_flags = '' 1033 if self.parameters['use resolution of identity']: 1034 jobex_flags += ' -ri' 1035 if self.parameters['force convergence']: 1036 par = self.parameters['force convergence'] 1037 conv = floor(-log10(par / Ha * Bohr)) 1038 jobex_flags += ' -gcart ' + str(int(conv)) 1039 if self.parameters['energy convergence']: 1040 par = self.parameters['energy convergence'] 1041 conv = floor(-log10(par / Ha)) 1042 jobex_flags += ' -energy ' + str(int(conv)) 1043 geom_iter = self.parameters['geometry optimization iterations'] 1044 if geom_iter is not None: 1045 assert isinstance(geom_iter, int) 1046 jobex_flags += ' -c ' + str(geom_iter) 1047 self.converged = False 1048 execute('jobex' + jobex_flags) 1049 # check convergence 1050 self.converged = self.read_convergence() 1051 if self.converged: 1052 self.update_energy = False 1053 self.update_forces = False 1054 self.update_geometry = False 1055 self.update_hessian = True 1056 # read results 1057 new_struct = read('coord') 1058 atoms.set_positions(new_struct.get_positions()) 1059 self.atoms = atoms.copy() 1060 self.read_energy() 1061 1062 def normal_mode_analysis(self, atoms=None): 1063 """execute normal mode analysis with modules aoforce or NumForce""" 1064 from ase.constraints import FixAtoms 1065 if atoms is None: 1066 atoms = self.atoms 1067 self.set_atoms(atoms) 1068 self.initialize() 1069 if self.update_energy: 1070 self.get_potential_energy(atoms) 1071 if self.update_hessian: 1072 fixatoms = [] 1073 for constr in atoms.constraints: 1074 if isinstance(constr, FixAtoms): 1075 ckwargs = constr.todict()['kwargs'] 1076 if 'indices' in ckwargs.keys(): 1077 fixatoms.extend(ckwargs['indices']) 1078 if self.parameters['numerical hessian'] is None: 1079 if len(fixatoms) > 0: 1080 define_str = '\n\ny\n' 1081 for index in fixatoms: 1082 define_str += 'm ' + str(index + 1) + ' 999.99999999\n' 1083 define_str += '*\n*\nn\nq\n' 1084 execute('define', input_str=define_str) 1085 dg = read_data_group('atoms') 1086 regex = r'(mass\s*=\s*)999.99999999' 1087 dg = re.sub(regex, r'\g<1>9999999999.9', dg) 1088 dg += '\n' 1089 delete_data_group('atoms') 1090 add_data_group(dg, raw=True) 1091 execute('aoforce') 1092 else: 1093 optstr = '' 1094 pdict = self.parameters['numerical hessian'] 1095 if self.parameters['use resolution of identity']: 1096 optstr += ' -ri' 1097 if len(fixatoms) > 0: 1098 optstr += ' -frznuclei -central -c' 1099 if 'central' in pdict.keys(): 1100 optstr += ' -central' 1101 if 'delta' in pdict.keys(): 1102 optstr += ' -d ' + str(pdict['delta'] / Bohr) 1103 execute('NumForce' + optstr) 1104 self.update_hessian = False 1105 1106 def read_restart(self): 1107 """read a previous calculation from control file""" 1108 self.atoms = read('coord') 1109 self.atoms.calc = self 1110 self.converged = self.read_convergence() 1111 read_methods = [ 1112 self.read_energy, 1113 self.read_gradient, 1114 self.read_forces, 1115 self.read_basis_set, 1116 self.read_ecps, 1117 self.read_mos, 1118 self.read_occupation_numbers, 1119 self.read_dipole_moment, 1120 self.read_ssquare, 1121 self.read_hessian, 1122 self.read_vibrational_reduced_masses, 1123 self.read_normal_modes, 1124 self.read_vibrational_spectrum, 1125 self.read_charges, 1126 self.read_point_charges, 1127 self.read_run_parameters 1128 ] 1129 for method in read_methods: 1130 try: 1131 method() 1132 except ReadError as err: 1133 warnings.warn(err.args[0]) 1134 1135 def read_parameters(self): 1136 """read parameters from control file""" 1137 1138 def parse_data_group(dg, dg_name): 1139 """parse a data group""" 1140 if len(dg) == 0: 1141 return None 1142 lsep = None 1143 ksep = None 1144 ndg = dg.replace('$' + dg_name, '').strip() 1145 if '\n' in ndg: 1146 lsep = '\n' 1147 if '=' in ndg: 1148 ksep = '=' 1149 if not lsep and not ksep: 1150 return ndg 1151 result = {} 1152 lines = ndg.split(lsep) 1153 for line in lines: 1154 fields = line.strip().split(ksep) 1155 if len(fields) == 2: 1156 result[fields[0]] = fields[1] 1157 elif len(fields) == 1: 1158 result[fields[0]] = True 1159 return result 1160 1161 params = {} 1162 pdgs = {} 1163 for p in self.parameter_group: 1164 if self.parameter_group[p] and self.parameter_key[p]: 1165 pdgs[p] = parse_data_group( 1166 read_data_group(self.parameter_group[p]), 1167 self.parameter_group[p] 1168 ) 1169 1170 for p in self.parameter_key: 1171 if self.parameter_key[p]: 1172 if self.parameter_key[p] == self.parameter_group[p]: 1173 if pdgs[p] is None: 1174 if self.parameter_type[p] is bool: 1175 params[p] = False 1176 else: 1177 params[p] = None 1178 else: 1179 if self.parameter_type[p] is bool: 1180 params[p] = True 1181 else: 1182 typ = self.parameter_type[p] 1183 val = typ(pdgs[p]) 1184 mapping = self.parameter_mapping 1185 if p in list(mapping.keys()): 1186 fun = mapping[p]['from_control'] 1187 val = fun(val) 1188 params[p] = val 1189 else: 1190 if pdgs[p] is None: 1191 params[p] = None 1192 elif isinstance(pdgs[p], str): 1193 if self.parameter_type[p] is bool: 1194 params[p] = (pdgs[p] == self.parameter_key[p]) 1195 else: 1196 if self.parameter_key[p] not in list(pdgs[p].keys()): 1197 if self.parameter_type[p] is bool: 1198 params[p] = False 1199 else: 1200 params[p] = None 1201 else: 1202 typ = self.parameter_type[p] 1203 val = typ(pdgs[p][self.parameter_key[p]]) 1204 mapping = self.parameter_mapping 1205 if p in list(mapping.keys()): 1206 fun = mapping[p]['from_control'] 1207 val = fun(val) 1208 params[p] = val 1209 1210 # non-group or non-key parameters 1211 1212 # per-element and per-atom basis sets not implemented in calculator 1213 basis_sets = set([bs['nickname'] for bs in self.results['basis set']]) 1214 assert len(basis_sets) == 1 1215 params['basis set name'] = list(basis_sets)[0] 1216 params['basis set definition'] = self.results['basis set'] 1217 1218 # rohf, multiplicity and total charge 1219 orbs = self.results['molecular orbitals'] 1220 params['rohf'] = (bool(len(read_data_group('rohf'))) or 1221 bool(len(read_data_group('roothaan')))) 1222 core_charge = 0 1223 if self.results['ecps']: 1224 for ecp in self.results['ecps']: 1225 for symbol in self.atoms.get_chemical_symbols(): 1226 if symbol.lower() == ecp['element'].lower(): 1227 core_charge -= ecp['number of core electrons'] 1228 if params['uhf']: 1229 alpha_occ = [o['occupancy'] for o in orbs if o['spin'] == 'alpha'] 1230 beta_occ = [o['occupancy'] for o in orbs if o['spin'] == 'beta'] 1231 spin = (np.sum(alpha_occ) - np.sum(beta_occ)) * 0.5 1232 params['multiplicity'] = int(2 * spin + 1) 1233 nuclear_charge = np.sum(self.atoms.numbers) 1234 electron_charge = -int(np.sum(alpha_occ) + np.sum(beta_occ)) 1235 electron_charge += core_charge 1236 params['total charge'] = nuclear_charge + electron_charge 1237 elif not params['rohf']: # restricted HF (closed shell) 1238 params['multiplicity'] = 1 1239 nuclear_charge = np.sum(self.atoms.numbers) 1240 electron_charge = -int(np.sum([o['occupancy'] for o in orbs])) 1241 electron_charge += core_charge 1242 params['total charge'] = nuclear_charge + electron_charge 1243 else: 1244 raise NotImplementedError('ROHF not implemented') 1245 1246 # task-related parameters 1247 if os.path.exists('job.start'): 1248 with open('job.start', 'r') as log: 1249 lines = log.readlines() 1250 for line in lines: 1251 if 'CRITERION FOR TOTAL SCF-ENERGY' in line: 1252 en = int(re.search(r'10\*{2}\(-(\d+)\)', line).group(1)) 1253 params['energy convergence'] = en 1254 if 'CRITERION FOR MAXIMUM NORM OF SCF-ENERGY GRADIENT' in line: 1255 gr = int(re.search(r'10\*{2}\(-(\d+)\)', line).group(1)) 1256 params['force convergence'] = gr 1257 if 'AN OPTIMIZATION WITH MAX' in line: 1258 cy = int(re.search(r'MAX. (\d+) CYCLES', line).group(1)) 1259 params['geometry optimization iterations'] = cy 1260 return params 1261 1262 def read_convergence(self): 1263 """perform convergence checks""" 1264 if self.restart: 1265 if bool(len(read_data_group('restart'))): 1266 return False 1267 if bool(len(read_data_group('actual'))): 1268 return False 1269 if not bool(len(read_data_group('energy'))): 1270 return False 1271 if (os.path.exists('job.start') and 1272 os.path.exists('GEO_OPT_FAILED')): 1273 return False 1274 return True 1275 1276 if self.parameters['task'] in ['optimize', 'geometry optimization']: 1277 if os.path.exists('GEO_OPT_CONVERGED'): 1278 return True 1279 elif os.path.exists('GEO_OPT_FAILED'): 1280 # check whether a failed scf convergence is the reason 1281 checkfiles = [] 1282 for filename in os.listdir('.'): 1283 if filename.startswith('job.'): 1284 checkfiles.append(filename) 1285 for filename in checkfiles: 1286 for line in open(filename): 1287 if 'SCF FAILED TO CONVERGE' in line: 1288 # scf did not converge in some jobex iteration 1289 if filename == 'job.last': 1290 raise RuntimeError('scf failed to converge') 1291 else: 1292 warnings.warn('scf failed to converge') 1293 warnings.warn('geometry optimization failed to converge') 1294 return False 1295 else: 1296 raise RuntimeError('error during geometry optimization') 1297 else: 1298 if os.path.isfile('dscf_problem'): 1299 raise RuntimeError('scf failed to converge') 1300 else: 1301 return True 1302 1303 def read_results(self): 1304 """read all results and load them in the results entity""" 1305 self.read_energy() 1306 self.read_mos() 1307 self.read_basis_set() 1308 self.read_occupation_numbers() 1309 self.read_dipole_moment() 1310 self.read_ssquare() 1311 self.read_run_parameters() 1312 if self.parameters['task'] in ['gradient', 'optimize', 1313 'gradient calculation', 1314 'geometry optimization']: 1315 self.read_gradient() 1316 self.read_forces() 1317 if self.parameters['task'] in ['frequencies', 'normal mode analysis']: 1318 self.read_hessian() 1319 self.read_vibrational_reduced_masses() 1320 self.read_normal_modes() 1321 self.read_vibrational_spectrum() 1322 self.read_charges() 1323 1324 def read_run_parameters(self): 1325 """read parameters set by define and not in self.parameters""" 1326 1327 if 'calculation parameters' not in self.results.keys(): 1328 self.results['calculation parameters'] = {} 1329 parameters = self.results['calculation parameters'] 1330 dg = read_data_group('symmetry') 1331 parameters['point group'] = str(dg.split()[1]) 1332 parameters['uhf'] = '$uhf' in read_data_group('uhf') 1333 # Gaussian function type 1334 gt = read_data_group('pople') 1335 if gt == '': 1336 parameters['gaussian type'] = 'spherical harmonic' 1337 else: 1338 gt = gt.split()[1] 1339 if gt == 'AO': 1340 parameters['gaussian type'] = 'spherical harmonic' 1341 elif gt == 'CAO': 1342 parameters['gaussian type'] = 'cartesian' 1343 else: 1344 parameters['gaussian type'] = None 1345 1346 nvibro = read_data_group('nvibro') 1347 if nvibro: 1348 parameters['nuclear degrees of freedom'] = int(nvibro.split()[1]) 1349 1350 def read_energy(self): 1351 """Read energy from Turbomole energy file.""" 1352 try: 1353 with open('energy', 'r') as enf: 1354 text = enf.read().lower() 1355 except IOError: 1356 raise ReadError('failed to read energy file') 1357 if text == '': 1358 raise ReadError('empty energy file') 1359 1360 lines = iter(text.split('\n')) 1361 1362 for line in lines: 1363 if line.startswith('$end'): 1364 break 1365 elif line.startswith('$'): 1366 pass 1367 else: 1368 energy_tmp = float(line.split()[1]) 1369 if self.post_HF: 1370 energy_tmp += float(line.split()[4]) 1371 # update energy units 1372 self.e_total = energy_tmp * Ha 1373 self.results['total energy'] = self.e_total 1374 1375 def read_forces(self): 1376 """Read Forces from Turbomole gradient file.""" 1377 dg = read_data_group('grad') 1378 if len(dg) == 0: 1379 return 1380 file = open('gradient', 'r') 1381 lines = file.readlines() 1382 file.close() 1383 1384 forces = np.array([[0, 0, 0]]) 1385 1386 nline = len(lines) 1387 iline = -1 1388 1389 for i in range(nline): 1390 if 'cycle' in lines[i]: 1391 iline = i 1392 1393 if iline < 0: 1394 raise RuntimeError('Please check TURBOMOLE gradients') 1395 1396 # next line 1397 iline += len(self.atoms) + 1 1398 # $end line 1399 nline -= 1 1400 # read gradients 1401 for i in range(iline, nline): 1402 line = lines[i].replace('D', 'E') 1403 tmp = np.array([[float(f) for f in line.split()[0:3]]]) 1404 forces = np.concatenate((forces, tmp)) 1405 # Note the '-' sign for turbomole, to get forces 1406 self.forces = -np.delete(forces, np.s_[0:1], axis=0) * Ha / Bohr 1407 self.results['energy gradient'] = (-self.forces).tolist() 1408 1409 def read_occupation_numbers(self): 1410 """read occupation numbers with module 'eiger' """ 1411 if 'molecular orbitals' not in self.results.keys(): 1412 return 1413 mos = self.results['molecular orbitals'] 1414 args = ['eiger', '--all', '--pview'] 1415 output = execute(args, error_test=False, stdout_tofile=False) 1416 lines = output.split('\n') 1417 for line in lines: 1418 regex = ( 1419 r'^\s+(\d+)\.*\s+(\w*)\s+(\d+)\s+(\S+)' 1420 r'\s+(\d*\.*\d*)\s+([-+]?\d+\.\d*)' 1421 ) 1422 match = re.search(regex, line) 1423 if match: 1424 orb_index = int(match.group(3)) 1425 if match.group(2) == 'a': 1426 spin = 'alpha' 1427 elif match.group(2) == 'b': 1428 spin = 'beta' 1429 else: 1430 spin = None 1431 ar_index = next( 1432 index for (index, molecular_orbital) in enumerate(mos) 1433 if (molecular_orbital['index'] == orb_index and 1434 molecular_orbital['spin'] == spin) 1435 ) 1436 mos[ar_index]['index by energy'] = int(match.group(1)) 1437 irrep = str(match.group(4)) 1438 mos[ar_index]['irreducible representation'] = irrep 1439 if match.group(5) != '': 1440 mos[ar_index]['occupancy'] = float(match.group(5)) 1441 else: 1442 mos[ar_index]['occupancy'] = float(0) 1443 1444 def read_mos(self): 1445 """read the molecular orbital coefficients and orbital energies 1446 from files mos, alpha and beta""" 1447 1448 self.results['molecular orbitals'] = [] 1449 mos = self.results['molecular orbitals'] 1450 keywords = ['scfmo', 'uhfmo_alpha', 'uhfmo_beta'] 1451 spin = [None, 'alpha', 'beta'] 1452 1453 for index, keyword in enumerate(keywords): 1454 flen = None 1455 mo = {} 1456 orbitals_coefficients_line = [] 1457 mo_string = read_data_group(keyword) 1458 if mo_string == '': 1459 continue 1460 mo_string += '\n$end' 1461 lines = mo_string.split('\n') 1462 for line in lines: 1463 if re.match(r'^\s*#', line): 1464 continue 1465 if 'eigenvalue' in line: 1466 if len(orbitals_coefficients_line) != 0: 1467 mo['eigenvector'] = orbitals_coefficients_line 1468 mos.append(mo) 1469 mo = {} 1470 orbitals_coefficients_line = [] 1471 regex = (r'^\s*(\d+)\s+(\S+)\s+' 1472 r'eigenvalue=([\+\-\d\.\w]+)\s') 1473 match = re.search(regex, line) 1474 mo['index'] = int(match.group(1)) 1475 mo['irreducible representation'] = str(match.group(2)) 1476 eig = float(re.sub('[dD]', 'E', match.group(3))) * Ha 1477 mo['eigenvalue'] = eig 1478 mo['spin'] = spin[index] 1479 mo['degeneracy'] = 1 1480 continue 1481 if keyword in line: 1482 # e.g. format(4d20.14) 1483 regex = r'format\(\d+[a-zA-Z](\d+)\.\d+\)' 1484 match = re.search(regex, line) 1485 if match: 1486 flen = int(match.group(1)) 1487 if ('scfdump' in line or 'expanded' in line or 1488 'scfconv' not in line): 1489 self.converged = False 1490 continue 1491 if '$end' in line: 1492 if len(orbitals_coefficients_line) != 0: 1493 mo['eigenvector'] = orbitals_coefficients_line 1494 mos.append(mo) 1495 break 1496 sfields = [line[i:i + flen] 1497 for i in range(0, len(line), flen)] 1498 ffields = [float(f.replace('D', 'E').replace('d', 'E')) 1499 for f in sfields] 1500 orbitals_coefficients_line += ffields 1501 1502 def read_basis_set(self): 1503 """read the basis set""" 1504 self.results['basis set'] = [] 1505 self.results['basis set formatted'] = {} 1506 bsf = read_data_group('basis') 1507 self.results['basis set formatted']['turbomole'] = bsf 1508 lines = bsf.split('\n') 1509 basis_set = {} 1510 functions = [] 1511 function = {} 1512 primitives = [] 1513 read_tag = False 1514 read_data = False 1515 for line in lines: 1516 if len(line.strip()) == 0: 1517 continue 1518 if '$basis' in line: 1519 continue 1520 if '$end' in line: 1521 break 1522 if re.match(r'^\s*#', line): 1523 continue 1524 if re.match(r'^\s*\*', line): 1525 if read_tag: 1526 read_tag = False 1527 read_data = True 1528 else: 1529 if read_data: 1530 # end primitives 1531 function['primitive functions'] = primitives 1532 function['number of primitives'] = len(primitives) 1533 primitives = [] 1534 functions.append(function) 1535 function = {} 1536 # end contracted 1537 basis_set['functions'] = functions 1538 functions = [] 1539 self.results['basis set'].append(basis_set) 1540 basis_set = {} 1541 read_data = False 1542 read_tag = True 1543 continue 1544 if read_tag: 1545 match = re.search(r'^\s*(\w+)\s+(.+)', line) 1546 if match: 1547 basis_set['element'] = match.group(1) 1548 basis_set['nickname'] = match.group(2) 1549 else: 1550 raise RuntimeError('error reading basis set') 1551 else: 1552 match = re.search(r'^\s+(\d+)\s+(\w+)', line) 1553 if match: 1554 if len(primitives) > 0: 1555 # end primitives 1556 function['primitive functions'] = primitives 1557 function['number of primitives'] = len(primitives) 1558 primitives = [] 1559 functions.append(function) 1560 function = {} 1561 # begin contracted 1562 function['shell type'] = str(match.group(2)) 1563 continue 1564 regex = ( 1565 r'^\s*([-+]?[0-9]*\.?[0-9]+([eE][-+]?[0-9]+)?)' 1566 r'\s+([-+]?[0-9]*\.?[0-9]+([eE][-+]?[0-9]+)?)' 1567 ) 1568 match = re.search(regex, line) 1569 if match: 1570 exponent = float(match.group(1)) 1571 coefficient = float(match.group(3)) 1572 primitives.append( 1573 {'exponent': exponent, 'coefficient': coefficient} 1574 ) 1575 1576 def read_ecps(self): 1577 """read the effective core potentials""" 1578 ecpf = read_data_group('ecp') 1579 if not bool(len(ecpf)): 1580 self.results['ecps'] = None 1581 self.results['ecps formatted'] = None 1582 return 1583 self.results['ecps'] = [] 1584 self.results['ecps formatted'] = {} 1585 self.results['ecps formatted']['turbomole'] = ecpf 1586 lines = ecpf.split('\n') 1587 ecp = {} 1588 groups = [] 1589 group = {} 1590 terms = [] 1591 read_tag = False 1592 read_data = False 1593 for line in lines: 1594 if len(line.strip()) == 0: 1595 continue 1596 if '$ecp' in line: 1597 continue 1598 if '$end' in line: 1599 break 1600 if re.match(r'^\s*#', line): 1601 continue 1602 if re.match(r'^\s*\*', line): 1603 if read_tag: 1604 read_tag = False 1605 read_data = True 1606 else: 1607 if read_data: 1608 # end terms 1609 group['terms'] = terms 1610 group['number of terms'] = len(terms) 1611 terms = [] 1612 groups.append(group) 1613 group = {} 1614 # end group 1615 ecp['groups'] = groups 1616 groups = [] 1617 self.results['ecps'].append(ecp) 1618 ecp = {} 1619 read_data = False 1620 read_tag = True 1621 continue 1622 if read_tag: 1623 match = re.search(r'^\s*(\w+)\s+(.+)', line) 1624 if match: 1625 ecp['element'] = match.group(1) 1626 ecp['nickname'] = match.group(2) 1627 else: 1628 raise RuntimeError('error reading ecp') 1629 else: 1630 regex = r'ncore\s*=\s*(\d+)\s+lmax\s*=\s*(\d+)' 1631 match = re.search(regex, line) 1632 if match: 1633 ecp['number of core electrons'] = int(match.group(1)) 1634 ecp['maximum angular momentum number'] = \ 1635 int(match.group(2)) 1636 continue 1637 match = re.search(r'^(\w(\-\w)?)', line) 1638 if match: 1639 if len(terms) > 0: 1640 # end terms 1641 group['terms'] = terms 1642 group['number of terms'] = len(terms) 1643 terms = [] 1644 groups.append(group) 1645 group = {} 1646 # begin group 1647 group['title'] = str(match.group(1)) 1648 continue 1649 regex = (r'^\s*([-+]?[0-9]*\.?[0-9]+([eE][-+]?[0-9]+)?)\s+' 1650 r'(\d)\s+([-+]?[0-9]*\.?[0-9]+([eE][-+]?[0-9]+)?)') 1651 match = re.search(regex, line) 1652 if match: 1653 terms.append( 1654 { 1655 'coefficient': float(match.group(1)), 1656 'power of r': float(match.group(3)), 1657 'exponent': float(match.group(4)) 1658 } 1659 ) 1660 1661 def read_gradient(self): 1662 """read all information in file 'gradient'""" 1663 from ase import Atom 1664 grad_string = read_data_group('grad') 1665 if len(grad_string) == 0: 1666 return 1667# try to reuse ase: 1668# structures = read('gradient', index=':') 1669 lines = grad_string.split('\n') 1670 history = [] 1671 image = {} 1672 gradient = [] 1673 atoms = Atoms() 1674 (cycle, energy, norm) = (None, None, None) 1675 for line in lines: 1676 # cycle lines 1677 regex = ( 1678 r'^\s*cycle =\s*(\d+)\s+' 1679 r'SCF energy =\s*([-+]?[0-9]*\.?[0-9]+([eE][-+]?[0-9]+)?)\s+' 1680 r'\|dE\/dxyz\| =\s*([-+]?[0-9]*\.?[0-9]+([eE][-+]?[0-9]+)?)' 1681 ) 1682 match = re.search(regex, line) 1683 if match: 1684 if len(atoms): 1685 image['optimization cycle'] = cycle 1686 image['total energy'] = energy 1687 image['gradient norm'] = norm 1688 image['energy gradient'] = gradient 1689 history.append(image) 1690 image = {} 1691 atoms = Atoms() 1692 gradient = [] 1693 cycle = int(match.group(1)) 1694 energy = float(match.group(2)) * Ha 1695 norm = float(match.group(4)) * Ha / Bohr 1696 continue 1697 # coordinate lines 1698 regex = ( 1699 r'^\s*([-+]?[0-9]*\.?[0-9]+([eEdD][-+]?[0-9]+)?)' 1700 r'\s+([-+]?[0-9]*\.?[0-9]+([eEdD][-+]?[0-9]+)?)' 1701 r'\s+([-+]?[0-9]*\.?[0-9]+([eEdD][-+]?[0-9]+)?)' 1702 r'\s+(\w+)' 1703 ) 1704 match = re.search(regex, line) 1705 if match: 1706 x = float(match.group(1)) * Bohr 1707 y = float(match.group(3)) * Bohr 1708 z = float(match.group(5)) * Bohr 1709 symbol = str(match.group(7)).capitalize() 1710 1711 if symbol == 'Q': 1712 symbol = 'X' 1713 atoms += Atom(symbol, (x, y, z)) 1714 1715 continue 1716 # gradient lines 1717 regex = ( 1718 r'^\s*([-+]?[0-9]*\.?[0-9]+([eEdD][-+]?[0-9]+)?)' 1719 r'\s+([-+]?[0-9]*\.?[0-9]+([eEdD][-+]?[0-9]+)?)' 1720 r'\s+([-+]?[0-9]*\.?[0-9]+([eEdD][-+]?[0-9]+)?)' 1721 ) 1722 match = re.search(regex, line) 1723 if match: 1724 gradx = float(match.group(1).replace('D', 'E')) * Ha / Bohr 1725 grady = float(match.group(3).replace('D', 'E')) * Ha / Bohr 1726 gradz = float(match.group(5).replace('D', 'E')) * Ha / Bohr 1727 gradient.append([gradx, grady, gradz]) 1728 1729 image['optimization cycle'] = cycle 1730 image['total energy'] = energy 1731 image['gradient norm'] = norm 1732 image['energy gradient'] = gradient 1733 history.append(image) 1734 self.results['geometry optimization history'] = history 1735 1736 def read_hessian(self, noproj=False): 1737 """Read in the hessian matrix""" 1738 self.results['hessian matrix'] = {} 1739 self.results['hessian matrix']['array'] = [] 1740 self.results['hessian matrix']['units'] = '?' 1741 self.results['hessian matrix']['projected'] = True 1742 self.results['hessian matrix']['mass weighted'] = True 1743 dg = read_data_group('nvibro') 1744 if len(dg) == 0: 1745 return 1746 nvibro = int(dg.split()[1]) 1747 self.results['hessian matrix']['dimension'] = nvibro 1748 row = [] 1749 key = 'hessian' 1750 if noproj: 1751 key = 'npr' + key 1752 self.results['hessian matrix']['projected'] = False 1753 lines = read_data_group(key).split('\n') 1754 for line in lines: 1755 if key in line: 1756 continue 1757 fields = line.split() 1758 row.extend(fields[2:len(fields)]) 1759 if len(row) == nvibro: 1760 # check whether it is mass-weighted 1761 float_row = [float(element) for element in row] 1762 self.results['hessian matrix']['array'].append(float_row) 1763 row = [] 1764 1765 def read_normal_modes(self, noproj=False): 1766 """Read in vibrational normal modes""" 1767 self.results['normal modes'] = {} 1768 self.results['normal modes']['array'] = [] 1769 self.results['normal modes']['projected'] = True 1770 self.results['normal modes']['mass weighted'] = True 1771 self.results['normal modes']['units'] = '?' 1772 dg = read_data_group('nvibro') 1773 if len(dg) == 0: 1774 return 1775 nvibro = int(dg.split()[1]) 1776 self.results['normal modes']['dimension'] = nvibro 1777 row = [] 1778 key = 'vibrational normal modes' 1779 if noproj: 1780 key = 'npr' + key 1781 self.results['normal modes']['projected'] = False 1782 lines = read_data_group(key).split('\n') 1783 for line in lines: 1784 if key in line: 1785 continue 1786 if '$end' in line: 1787 break 1788 fields = line.split() 1789 row.extend(fields[2:len(fields)]) 1790 if len(row) == nvibro: 1791 # check whether it is mass-weighted 1792 float_row = [float(element) for element in row] 1793 self.results['normal modes']['array'].append(float_row) 1794 row = [] 1795 1796 def read_vibrational_reduced_masses(self): 1797 """Read vibrational reduced masses""" 1798 self.results['vibrational reduced masses'] = [] 1799 dg = read_data_group('vibrational reduced masses') 1800 if len(dg) == 0: 1801 return 1802 lines = dg.split('\n') 1803 for line in lines: 1804 if '$vibrational' in line: 1805 continue 1806 if '$end' in line: 1807 break 1808 fields = [float(element) for element in line.split()] 1809 self.results['vibrational reduced masses'].extend(fields) 1810 1811 def read_vibrational_spectrum(self, noproj=False): 1812 """Read the vibrational spectrum""" 1813 self.results['vibrational spectrum'] = [] 1814 key = 'vibrational spectrum' 1815 if noproj: 1816 key = 'npr' + key 1817 lines = read_data_group(key).split('\n') 1818 for line in lines: 1819 dictionary = {} 1820 regex = ( 1821 r'^\s+(\d+)\s+(\S*)\s+([-+]?\d+\.\d*)' 1822 r'\s+(\d+\.\d*)\s+(\S+)\s+(\S+)' 1823 ) 1824 match = re.search(regex, line) 1825 if match: 1826 dictionary['mode number'] = int(match.group(1)) 1827 dictionary['irreducible representation'] = str(match.group(2)) 1828 dictionary['frequency'] = { 1829 'units': 'cm^-1', 1830 'value': float(match.group(3)) 1831 } 1832 dictionary['infrared intensity'] = { 1833 'units': 'km/mol', 1834 'value': float(match.group(4)) 1835 } 1836 1837 if match.group(5) == 'YES': 1838 dictionary['infrared active'] = True 1839 elif match.group(5) == 'NO': 1840 dictionary['infrared active'] = False 1841 else: 1842 dictionary['infrared active'] = None 1843 1844 if match.group(6) == 'YES': 1845 dictionary['Raman active'] = True 1846 elif match.group(6) == 'NO': 1847 dictionary['Raman active'] = False 1848 else: 1849 dictionary['Raman active'] = None 1850 1851 self.results['vibrational spectrum'].append(dictionary) 1852 1853 def read_ssquare(self): 1854 """Read the expectation value of S^2 operator""" 1855 s2_string = read_data_group('ssquare from dscf') 1856 if s2_string == '': 1857 return 1858 string = s2_string.split('\n')[1] 1859 ssquare = float(re.search(r'^\s*(\d+\.*\d*)', string).group(1)) 1860 self.results['ssquare from scf calculation'] = ssquare 1861 1862 def read_dipole_moment(self): 1863 """Read the dipole moment""" 1864 dip_string = read_data_group('dipole') 1865 if dip_string == '': 1866 return 1867 lines = dip_string.split('\n') 1868 for line in lines: 1869 regex = ( 1870 r'^\s+x\s+([-+]?\d+\.\d*)\s+y\s+([-+]?\d+\.\d*)' 1871 r'\s+z\s+([-+]?\d+\.\d*)\s+a\.u\.' 1872 ) 1873 match = re.search(regex, line) 1874 if match: 1875 dip_vec = [float(match.group(c)) for c in range(1, 4)] 1876 regex = r'^\s+\| dipole \| =\s+(\d+\.*\d*)\s+debye' 1877 match = re.search(regex, line) 1878 if match: 1879 dip_abs_val = float(match.group(1)) 1880 self.results['electric dipole moment'] = {} 1881 self.results['electric dipole moment']['vector'] = { 1882 'array': dip_vec, 1883 'units': 'a.u.' 1884 } 1885 self.results['electric dipole moment']['absolute value'] = { 1886 'value': dip_abs_val, 1887 'units': 'Debye' 1888 } 1889 self.dipole = np.array(dip_vec) * Bohr 1890 1891 def read_version(self): 1892 """read the version from the tm output if stored in a file""" 1893 versions = read_output(r'TURBOMOLE\s+V(\d+\.\d+)\s+') 1894 if len(set(versions)) > 1: 1895 warnings.warn('different turbomole versions detected') 1896 self.version = list(set(versions)) 1897 elif len(versions) == 0: 1898 warnings.warn('no turbomole version detected') 1899 self.version = None 1900 else: 1901 self.version = versions[0] 1902 1903 def read_datetime(self): 1904 """read the datetime of the most recent calculation 1905 from the tm output if stored in a file 1906 """ 1907 datetimes = read_output( 1908 r'(\d{4}-[01]\d-[0-3]\d([T\s][0-2]\d:[0-5]' 1909 r'\d:[0-5]\d\.\d+)?([+-][0-2]\d:[0-5]\d|Z)?)') 1910 if len(datetimes) == 0: 1911 warnings.warn('no turbomole datetime detected') 1912 self.datetime = None 1913 else: 1914 # take the most recent time stamp 1915 self.datetime = sorted(datetimes, reverse=True)[0] 1916 1917 def read_runtime(self): 1918 """read the total runtime of calculations""" 1919 hits = read_output(r'total wall-time\s+:\s+(\d+.\d+)\s+seconds') 1920 if len(hits) == 0: 1921 warnings.warn('no turbomole runtimes detected') 1922 self.runtime = None 1923 else: 1924 self.runtime = np.sum([float(a) for a in hits]) 1925 1926 def read_hostname(self): 1927 """read the hostname of the computer on which the calc has run""" 1928 hostnames = read_output(r'hostname is\s+(.+)') 1929 if len(set(hostnames)) > 1: 1930 warnings.warn('runs on different hosts detected') 1931 self.hostname = list(set(hostnames)) 1932 else: 1933 self.hostname = hostnames[0] 1934 1935 def get_optimizer(self, atoms, trajectory=None, logfile=None): 1936 """returns a TurbomoleOptimizer object""" 1937 self.parameters['task'] = 'optimize' 1938 self.verify_parameters() 1939 return TurbomoleOptimizer(atoms, self) 1940 1941 def get_results(self): 1942 """returns the results dictionary""" 1943 return self.results 1944 1945 def get_potential_energy(self, atoms, force_consistent=True): 1946 # update atoms 1947 self.updated = self.e_total is None 1948 self.set_atoms(atoms) 1949 self.initialize() 1950 # if update of energy is necessary 1951 if self.update_energy: 1952 # calculate energy 1953 execute(self.calculate_energy) 1954 # check convergence 1955 self.converged = self.read_convergence() 1956 if not self.converged: 1957 return None 1958 # read energy 1959 self.read_energy() 1960 1961 self.update_energy = False 1962 return self.e_total 1963 1964 def get_forces(self, atoms): 1965 # update atoms 1966 self.updated = self.forces is None 1967 self.set_atoms(atoms) 1968 # complete energy calculations 1969 if self.update_energy: 1970 self.get_potential_energy(atoms) 1971 # if update of forces is necessary 1972 if self.update_forces: 1973 # calculate forces 1974 execute(self.calculate_forces) 1975 # read forces 1976 self.read_forces() 1977 1978 self.update_forces = False 1979 return self.forces.copy() 1980 1981 def get_dipole_moment(self, atoms): 1982 self.get_potential_energy(atoms) 1983 self.read_dipole_moment() 1984 return self.dipole 1985 1986 def get_property(self, name, atoms=None, allow_calculation=True): 1987 """return the value of a property""" 1988 1989 if name not in self.implemented_properties: 1990 # an ugly work around; the caller should test the raised error 1991 # if name in ['magmom', 'magmoms', 'charges', 'stress']: 1992 # return None 1993 raise PropertyNotImplementedError(name) 1994 1995 if atoms is None: 1996 atoms = self.atoms.copy() 1997 1998 persist_property = { 1999 'energy': 'e_total', 2000 'forces': 'forces', 2001 'dipole': 'dipole', 2002 'free_energy': 'e_total', 2003 'charges': 'charges' 2004 } 2005 property_getter = { 2006 'energy': self.get_potential_energy, 2007 'forces': self.get_forces, 2008 'dipole': self.get_dipole_moment, 2009 'free_energy': self.get_potential_energy, 2010 'charges': self.get_charges 2011 } 2012 getter_args = { 2013 'energy': [atoms], 2014 'forces': [atoms], 2015 'dipole': [atoms], 2016 'free_energy': [atoms, True], 2017 'charges': [atoms] 2018 } 2019 2020 if allow_calculation: 2021 result = property_getter[name](*getter_args[name]) 2022 else: 2023 if hasattr(self, persist_property[name]): 2024 result = getattr(self, persist_property[name]) 2025 else: 2026 result = None 2027 2028 if isinstance(result, np.ndarray): 2029 result = result.copy() 2030 return result 2031 2032 def get_charges(self, atoms): 2033 """return partial charges on atoms from an ESP fit""" 2034 self.get_potential_energy(atoms) 2035 self.read_charges() 2036 return self.charges 2037 2038 def read_charges(self): 2039 """read partial charges on atoms from an ESP fit""" 2040 epsfit_defined = ('esp fit' in self.parameters and 2041 self.parameters['esp fit'] is not None) 2042 if epsfit_defined or len(read_data_group('esp_fit')) > 0: 2043 filename = 'ASE.TM.' + self.calculate_energy + '.out' 2044 with open(filename, 'r') as infile: 2045 lines = infile.readlines() 2046 oklines = None 2047 for n, line in enumerate(lines): 2048 if 'atom radius/au charge' in line: 2049 oklines = lines[n + 1:n + len(self.atoms) + 1] 2050 if oklines is not None: 2051 qm_charges = [float(line.split()[3]) for line in oklines] 2052 self.charges = np.array(qm_charges) 2053 2054 def get_forces_on_point_charges(self): 2055 """return forces acting on point charges""" 2056 self.get_forces(self.atoms) 2057 lines = read_data_group('point_charge_gradients').split('\n')[1:] 2058 forces = [] 2059 for line in lines: 2060 linef = line.strip().replace('D', 'E') 2061 forces.append([float(x) for x in linef.split()]) 2062 # Note the '-' sign for turbomole, to get forces 2063 return -np.array(forces) * Ha / Bohr 2064 2065 def set_point_charges(self, pcpot=None): 2066 """write external point charges to control""" 2067 if pcpot is not None and pcpot != self.pcpot: 2068 self.pcpot = pcpot 2069 if self.pcpot.mmcharges is None or self.pcpot.mmpositions is None: 2070 raise RuntimeError('external point charges not defined') 2071 2072 if not self.pc_initialized: 2073 if len(read_data_group('point_charges')) == 0: 2074 add_data_group('point_charges', 'file=pc.txt') 2075 if len(read_data_group('point_charge_gradients')) == 0: 2076 add_data_group( 2077 'point_charge_gradients', 2078 'file=pc_gradients.txt' 2079 ) 2080 drvopt = read_data_group('drvopt') 2081 if 'point charges' not in drvopt: 2082 drvopt += '\n point charges\n' 2083 delete_data_group('drvopt') 2084 add_data_group(drvopt, raw=True) 2085 self.pc_initialized = True 2086 2087 if self.pcpot.updated: 2088 with open('pc.txt', 'w') as pcfile: 2089 pcfile.write('$point_charges nocheck list\n') 2090 for (x, y, z), charge in zip( 2091 self.pcpot.mmpositions, self.pcpot.mmcharges): 2092 pcfile.write('%20.14f %20.14f %20.14f %20.14f\n' 2093 % (x / Bohr, y / Bohr, z / Bohr, charge)) 2094 pcfile.write('$end \n') 2095 self.pcpot.updated = False 2096 2097 def read_point_charges(self): 2098 """read point charges from previous calculation""" 2099 pcs = read_data_group('point_charges') 2100 if len(pcs) > 0: 2101 lines = pcs.split('\n')[1:] 2102 (charges, positions) = ([], []) 2103 for line in lines: 2104 columns = [float(col) for col in line.strip().split()] 2105 positions.append([col * Bohr for col in columns[0:3]]) 2106 charges.append(columns[3]) 2107 self.pcpot = PointChargePotential(charges, positions) 2108 2109 def embed(self, charges=None, positions=None): 2110 """embed atoms in an array of point-charges; function used in 2111 qmmm calculations.""" 2112 self.pcpot = PointChargePotential(charges, positions) 2113 return self.pcpot 2114 2115 2116class PointChargePotential: 2117 """Point-charge potential for Turbomole""" 2118 def __init__(self, mmcharges, mmpositions=None): 2119 self.mmcharges = mmcharges 2120 self.mmpositions = mmpositions 2121 self.mmforces = None 2122 self.updated = True 2123 2124 def set_positions(self, mmpositions): 2125 """set the positions of point charges""" 2126 self.mmpositions = mmpositions 2127 self.updated = True 2128 2129 def set_charges(self, mmcharges): 2130 """set the values of point charges""" 2131 self.mmcharges = mmcharges 2132 self.updated = True 2133 2134 def get_forces(self, calc): 2135 """forces acting on point charges""" 2136 self.mmforces = calc.get_forces_on_point_charges() 2137 return self.mmforces 2138