root / ase / calculators / turbomole.py @ 3
Historique | Voir | Annoter | Télécharger (6,27 ko)
1 |
"""This module defines an ASE interface to Turbomole
|
---|---|
2 |
|
3 |
http://www.turbomole.com/
|
4 |
"""
|
5 |
import os, sys, string |
6 |
|
7 |
from ase.units import Hartree, Bohr |
8 |
from ase.io.turbomole import write_turbomole |
9 |
from ase.calculators.general import Calculator |
10 |
|
11 |
import numpy as np |
12 |
|
13 |
class Turbomole(Calculator): |
14 |
"""Class for doing Turbomole calculations.
|
15 |
"""
|
16 |
def __init__(self, label='turbomole', calculate_energy="dscf", calculate_forces="grad", post_HF = False): |
17 |
"""Construct TURBOMOLE-calculator object.
|
18 |
|
19 |
Parameters
|
20 |
==========
|
21 |
label: str
|
22 |
Prefix to use for filenames (label.txt, ...).
|
23 |
Default is 'turbomole'.
|
24 |
|
25 |
Examples
|
26 |
========
|
27 |
This is poor man's version of ASE-Turbomole:
|
28 |
|
29 |
First you do a normal turbomole preparation using turbomole's
|
30 |
define-program for the initial and final states.
|
31 |
|
32 |
(for instance in subdirectories Initial and Final)
|
33 |
|
34 |
Then relax the initial and final coordinates with desired constraints
|
35 |
using standard turbomole.
|
36 |
|
37 |
Copy the relaxed initial and final coordinates
|
38 |
(initial.coord and final.coord)
|
39 |
and the turbomole related files (from the subdirectory Initial)
|
40 |
control, coord, alpha, beta, mos, basis, auxbasis etc to the directory
|
41 |
you do the diffusion run.
|
42 |
|
43 |
For instance:
|
44 |
cd $My_Turbomole_diffusion_directory
|
45 |
cd Initial; cp control alpha beta mos basis auxbasis ../;
|
46 |
cp coord ../initial.coord;
|
47 |
cd ../;
|
48 |
cp Final/coord ./final.coord;
|
49 |
mkdir Gz ; cp * Gz ; gzip -r Gz
|
50 |
|
51 |
from ase.io import read
|
52 |
from ase.calculators.turbomole import Turbomole
|
53 |
a = read('coord', index=-1, format='turbomole')
|
54 |
calc = Turbomole()
|
55 |
a.set_calculator(calc)
|
56 |
e = a.get_potential_energy()
|
57 |
|
58 |
"""
|
59 |
|
60 |
self.label = label
|
61 |
self.converged = False |
62 |
|
63 |
# set calculators for energy and forces
|
64 |
self.calculate_energy = calculate_energy
|
65 |
self.calculate_forces = calculate_forces
|
66 |
|
67 |
# turbomole has no stress
|
68 |
self.stress = np.empty((3, 3)) |
69 |
|
70 |
# storage for energy and forces
|
71 |
self.e_total = None |
72 |
self.forces = None |
73 |
self.updated = False |
74 |
|
75 |
# atoms must be set
|
76 |
self.atoms = None |
77 |
|
78 |
# POST-HF method
|
79 |
self.post_HF = post_HF
|
80 |
|
81 |
def initialize(self, atoms): |
82 |
self.numbers = atoms.get_atomic_numbers().copy()
|
83 |
self.species = []
|
84 |
for a, Z in enumerate(self.numbers): |
85 |
self.species.append(Z)
|
86 |
self.converged = False |
87 |
|
88 |
def execute(self, command): |
89 |
from subprocess import Popen, PIPE |
90 |
try:
|
91 |
proc = Popen([command], shell=True, stderr=PIPE)
|
92 |
error = proc.communicate()[1]
|
93 |
if "abnormally" in error: |
94 |
raise OSError(error) |
95 |
print "TM command: ", command, "successfully executed" |
96 |
except OSError, e: |
97 |
print >>sys.stderr, "Execution failed:", e |
98 |
sys.exit(1)
|
99 |
|
100 |
def get_potential_energy(self, atoms): |
101 |
# update atoms
|
102 |
self.set_atoms(atoms)
|
103 |
# if update of energy is neccessary
|
104 |
if self.update_energy: |
105 |
# calculate energy
|
106 |
self.execute(self.calculate_energy + " > ASE.TM.energy.out") |
107 |
# check for convergence of dscf cycle
|
108 |
if os.path.isfile('dscf_problem'): |
109 |
print 'Turbomole scf energy calculation did not converge' |
110 |
raise RuntimeError, \ |
111 |
'Please run Turbomole define and come thereafter back'
|
112 |
# read energy
|
113 |
self.read_energy()
|
114 |
else :
|
115 |
print "taking old values (E)" |
116 |
self.update_energy = False |
117 |
return self.e_total |
118 |
|
119 |
def get_forces(self, atoms): |
120 |
# update atoms
|
121 |
self.set_atoms(atoms)
|
122 |
# complete energy calculations
|
123 |
if self.update_energy: |
124 |
self.get_potential_energy(atoms)
|
125 |
# if update of forces is neccessary
|
126 |
if self.update_forces: |
127 |
# calculate forces
|
128 |
self.execute(self.calculate_forces + " > ASE.TM.forces.out") |
129 |
# read forces
|
130 |
self.read_forces()
|
131 |
else :
|
132 |
print "taking old values (F)" |
133 |
self.update_forces = False |
134 |
return self.forces.copy() |
135 |
|
136 |
def get_stress(self, atoms): |
137 |
return self.stress.copy() |
138 |
|
139 |
def set_atoms(self, atoms): |
140 |
if self.atoms == atoms: |
141 |
return
|
142 |
# performs an update of the atoms
|
143 |
super(Turbomole, self).set_atoms(atoms) |
144 |
write_turbomole('coord', atoms)
|
145 |
# energy and forces must be re-calculated
|
146 |
self.update_energy = True |
147 |
self.update_forces = True |
148 |
|
149 |
def read_energy(self): |
150 |
"""Read Energy from Turbomole energy file."""
|
151 |
text = open('energy', 'r').read().lower() |
152 |
lines = iter(text.split('\n')) |
153 |
|
154 |
# Energy:
|
155 |
for line in lines: |
156 |
if line.startswith('$end'): |
157 |
break
|
158 |
elif line.startswith('$'): |
159 |
pass
|
160 |
else:
|
161 |
# print 'LINE',line
|
162 |
energy_tmp = float(line.split()[1]) |
163 |
if self.post_HF: |
164 |
energy_tmp += float(line.split()[4]) |
165 |
self.e_total = energy_tmp * Hartree
|
166 |
|
167 |
|
168 |
def read_forces(self): |
169 |
"""Read Forces from Turbomole gradient file."""
|
170 |
file = open('gradient','r') |
171 |
lines=file.readlines()
|
172 |
file.close()
|
173 |
|
174 |
forces = np.array([[0, 0, 0]]) |
175 |
|
176 |
nline = len(lines)
|
177 |
iline = -1
|
178 |
|
179 |
for i in xrange(nline): |
180 |
if 'cycle' in lines[i]: |
181 |
iline = i |
182 |
|
183 |
if iline < 0: |
184 |
print 'Gradients not found' |
185 |
raise RuntimeError("Please check TURBOMOLE gradients") |
186 |
|
187 |
iline += len(self.atoms) + 1 |
188 |
nline -= 1
|
189 |
for i in xrange(iline, nline): |
190 |
line = string.replace(lines[i],'D','E') |
191 |
#tmp=np.append(forces_tmp,np.array\
|
192 |
# ([[float(f) for f in line.split()[0:3]]]))
|
193 |
tmp=np.array([[float(f) for f in line.split()[0:3]]]) |
194 |
forces=np.concatenate((forces,tmp)) |
195 |
#note the '-' sign for turbomole, to get forces
|
196 |
self.forces = (-np.delete(forces, np.s_[0:1], axis=0))*Hartree/Bohr |