Statistiques
| Révision :

root / ase / calculators / turbomole.py @ 1

Historique | Voir | Annoter | Télécharger (6,95 ko)

1
"""This module defines an ASE interface to Turbomole
2

3
http://www.turbomole.com/
4
"""
5

    
6

    
7
import os, sys, string
8

    
9
from ase.data import chemical_symbols
10
from ase.units import Hartree, Bohr
11
from ase.io.turbomole import write_turbomole
12

    
13
import numpy as np
14

    
15
class Turbomole:
16
    """Class for doing Turbomole calculations.
17
    """
18
    def __init__(self, label='turbomole'):
19
        """Construct TURBOMOLE-calculator object.
20

21
        Parameters
22
        ==========
23
        label: str
24
            Prefix to use for filenames (label.txt, ...).
25
            Default is 'turbomole'.
26

27
        Examples
28
        ========
29
        This is poor man's version of ASE-Turbomole:
30

31
        First you do a normal turbomole preparation using turbomole's
32
        define-program for the initial and final states.
33

34
        (for instance in subdirectories Initial and Final)
35

36
        Then relax the initial and final coordinates with desired constraints
37
        using standard turbomole.
38

39
        Copy the relaxed initial and final coordinates 
40
        (initial.coord and final.coord)
41
        and the turbomole related files (from the subdirectory Initial) 
42
        control, coord, alpha, beta, mos, basis, auxbasis etc to the directory
43
        you do the diffusion run.
44
        
45
        For instance:
46
        cd $My_Turbomole_diffusion_directory
47
        cd Initial; cp control alpha beta mos basis auxbasis ../;
48
        cp coord ../initial.coord;
49
        cd ../;
50
        cp Final/coord ./final.coord;
51
        mkdir Gz ; cp * Gz ; gzip -r Gz
52
        
53
        from ase.io import read
54
        from ase.calculators.turbomole import Turbomole
55
        a = read('coord', index=-1, format='turbomole')
56
        calc = Turbomole()
57
        a.set_calculator(calc)
58
        e = a.get_potential_energy()
59
        
60
        """
61
        
62
        # get names of executables for turbomole energy and forces
63
        # get the path
64

    
65
        os.system('rm -f sysname.file; sysname > sysname.file')
66
        f = open('sysname.file')
67
        architechture = f.readline()[:-1]
68
        f.close()
69
        tmpath = os.environ['TURBODIR']
70
        pre = tmpath + '/bin/' + architechture + '/'
71

    
72

    
73
        if os.path.isfile('control'):
74
            f = open('control')
75
        else:
76
            print 'File control is missing'
77
            raise RuntimeError, \
78
                'Please run Turbomole define and come thereafter back'
79
        lines = f.readlines()
80
        f.close()
81
        self.tm_program_energy=pre+'dscf'
82
        self.tm_program_forces=pre+'grad'        
83
        for line in lines:
84
            if line.startswith('$ricore'):
85
                self.tm_program_energy=pre+'ridft'
86
                self.tm_program_forces=pre+'rdgrad'
87

    
88
        self.label = label
89
        self.converged = False
90
        #clean up turbomole energy file
91
        os.system('rm -f energy; touch energy')
92

    
93
        #turbomole has no stress
94
        self.stress = np.empty((3, 3))
95
        
96

    
97
    def update(self, atoms):
98
        """Energy and forces are calculated when atoms have moved
99
        by calling self.calculate
100
        """
101
        if (not self.converged or
102
            len(self.numbers) != len(atoms) or
103
            (self.numbers != atoms.get_atomic_numbers()).any()):
104
            self.initialize(atoms)
105
            self.calculate(atoms)
106
        elif ((self.positions != atoms.get_positions()).any() or
107
              (self.pbc != atoms.get_pbc()).any() or
108
              (self.cell != atoms.get_cell()).any()):
109
            self.calculate(atoms)
110

    
111
    def initialize(self, atoms):
112
        self.numbers = atoms.get_atomic_numbers().copy()
113
        self.species = []
114
        for a, Z in enumerate(self.numbers):
115
            self.species.append(Z)
116
        self.converged = False
117
        
118
    def get_potential_energy(self, atoms):
119
        self.update(atoms)
120
        return self.etotal
121

    
122
    def get_forces(self, atoms):
123
        self.update(atoms)
124
        return self.forces.copy()
125
    
126
    def get_stress(self, atoms):
127
        self.update(atoms)
128
        return self.stress.copy()
129

    
130
    def calculate(self, atoms):
131
        """Total Turbomole energy is calculated (to file 'energy'
132
        also forces are calculated (to file 'gradient')
133
        """
134
        self.positions = atoms.get_positions().copy()
135
        self.cell = atoms.get_cell().copy()
136
        self.pbc = atoms.get_pbc().copy()
137

    
138
        #write current coordinates to file 'coord' for Turbomole
139
        write_turbomole('coord', atoms)
140

    
141

    
142
        #Turbomole energy calculation
143
        os.system('rm -f output.energy.dummy; \
144
                      '+ self.tm_program_energy +'> output.energy.dummy')
145

    
146
        #check that the energy run converged
147
        if os.path.isfile('dscf_problem'):
148
            print 'Turbomole scf energy calculation did not converge'
149
            print 'issue command t2x -c > last.xyz'
150
            print 'and check geometry last.xyz and job.xxx or statistics'
151
            raise RuntimeError, \
152
                'Please run Turbomole define and come thereafter back'
153

    
154

    
155
        self.read_energy()
156

    
157
        #Turbomole atomic forces calculation
158
        #killing the previous gradient file because 
159
        #turbomole gradients are affected by the previous values
160
        os.system('rm -f gradient; rm -f output.forces.dummy; \
161
                      '+ self.tm_program_forces +'> output.forces.dummy')
162

    
163
        self.read_forces(atoms)
164

    
165
        self.converged = True
166

    
167
        
168
    def read_energy(self):
169
        """Read Energy from Turbomole energy file."""
170
        text = open('energy', 'r').read().lower()
171
        lines = iter(text.split('\n'))
172

    
173
        # Energy:
174
        for line in lines:
175
            if line.startswith('$end'):
176
                break
177
            elif line.startswith('$'):
178
                pass
179
            else:
180
                #print 'LINE',line
181
                energy_tmp = float(line.split()[1])
182
        #print 'energy_tmp',energy_tmp
183
        self.etotal = energy_tmp * Hartree
184

    
185

    
186
    def read_forces(self,atoms):
187
        """Read Forces from Turbomole gradient file."""
188

    
189
        file = open('gradient','r')
190
        line=file.readline()
191
        line=file.readline()
192
        tmpforces = np.array([[0, 0, 0]])
193
        while line:
194
            if 'cycle' in line:
195
                for i, dummy in enumerate(atoms):
196
                            line=file.readline()
197
                forces_tmp=[]
198
                for i, dummy in enumerate(atoms):
199
                            line=file.readline()
200
                            line2=string.replace(line,'D','E')
201
                            #tmp=np.append(forces_tmp,np.array\
202
                            #      ([[float(f) for f in line2.split()[0:3]]]))
203
                            tmp=np.array\
204
                                ([[float(f) for f in line2.split()[0:3]]])
205
                            tmpforces=np.concatenate((tmpforces,tmp))  
206
            line=file.readline()
207
            
208

    
209
        #note the '-' sign for turbomole, to get forces
210
        self.forces = (-np.delete(tmpforces, np.s_[0:1], axis=0))*Hartree/Bohr
211

    
212
        #print 'forces', self.forces
213

    
214
    def read(self):
215
        """Dummy stress for turbomole"""
216
        self.stress = np.empty((3, 3))