Statistiques
| Révision :

root / ase / thermochemistry.py @ 1

Historique | Voir | Annoter | Télécharger (14,07 ko)

1 1 tkerber
#!/usr/bin/env python
2 1 tkerber
"""Modules for calculating thermochemical information from computational
3 1 tkerber
outputs."""
4 1 tkerber
5 1 tkerber
import os
6 1 tkerber
import sys
7 1 tkerber
import numpy as np
8 1 tkerber
9 1 tkerber
from ase import units
10 1 tkerber
11 1 tkerber
def rotationalinertia(atoms):
12 1 tkerber
    """Calculates the three principle moments of inertia for an ASE atoms
13 1 tkerber
    object. This uses the atomic masses from ASE, which (if not explicitly
14 1 tkerber
    specified by the user) gives an inexact approximation of an isotopically
15 1 tkerber
    averaged result. Units are in amu*angstroms**2."""
16 1 tkerber
17 1 tkerber
    # Calculate the center of mass.
18 1 tkerber
    xcm, ycm, zcm = atoms.get_center_of_mass()
19 1 tkerber
    masses = atoms.get_masses()
20 1 tkerber
21 1 tkerber
    # Calculate moments of inertia in the current frame of reference.
22 1 tkerber
    Ixx = 0.
23 1 tkerber
    Iyy = 0.
24 1 tkerber
    Izz = 0.
25 1 tkerber
    Ixy = 0.
26 1 tkerber
    Ixz = 0.
27 1 tkerber
    Iyz = 0.
28 1 tkerber
    for index, atom in enumerate(atoms):
29 1 tkerber
        m = masses[index]
30 1 tkerber
        x = atom.get_x() - xcm
31 1 tkerber
        y = atom.get_y() - ycm
32 1 tkerber
        z = atom.get_z() - zcm
33 1 tkerber
        Ixx += m * (y**2. + z**2.)
34 1 tkerber
        Iyy += m * (x**2. + z**2.)
35 1 tkerber
        Izz += m * (x**2. + y**2.)
36 1 tkerber
        Ixy += m * x * y
37 1 tkerber
        Ixz += m * x * z
38 1 tkerber
        Iyz += m * y * z
39 1 tkerber
    # Create the inertia tensor in the current frame of reference.
40 1 tkerber
    I_ = np.matrix([[ Ixx, -Ixy, -Ixz],
41 1 tkerber
                    [-Ixy,  Iyy, -Iyz],
42 1 tkerber
                    [-Ixz, -Iyz,  Izz]])
43 1 tkerber
    # Find the eigenvalues, which are the principle moments of inertia.
44 1 tkerber
    I = np.linalg.eigvals(I_)
45 1 tkerber
    return I
46 1 tkerber
47 1 tkerber
class ThermoChem:
48 1 tkerber
    """Base class containing common methods used in thermochemistry
49 1 tkerber
    calculations."""
50 1 tkerber
51 1 tkerber
    def get_ZPE_correction(self):
52 1 tkerber
        """Returns the zero-point vibrational energy correction in eV."""
53 1 tkerber
        zpe = 0.
54 1 tkerber
        for energy in self.vib_energies:
55 1 tkerber
            zpe += 0.5 * energy
56 1 tkerber
        return zpe
57 1 tkerber
58 1 tkerber
    def _vibrational_energy_contribution(self, temperature):
59 1 tkerber
        """Calculates the change in internal energy due to vibrations from
60 1 tkerber
        0K to the specified temperature for a set of vibrations given in
61 1 tkerber
        inverse centimeters and a temperature given in Kelvin. Returns the
62 1 tkerber
        energy change in eV."""
63 1 tkerber
        kT = units.kB * temperature
64 1 tkerber
        dU = 0.
65 1 tkerber
        for energy in self.vib_energies:
66 1 tkerber
            dU += energy / (np.exp(energy/kT) - 1.)
67 1 tkerber
        return dU
68 1 tkerber
69 1 tkerber
    def _vibrational_entropy_contribution(self, temperature):
70 1 tkerber
        """Calculates the entropy due to vibrations for a set of vibrations
71 1 tkerber
        given in inverse centimeters and a temperature given in Kelvin.
72 1 tkerber
        Returns the entropy in eV/K."""
73 1 tkerber
        kT = units.kB * temperature
74 1 tkerber
        S_v = 0.
75 1 tkerber
        for energy in self.vib_energies:
76 1 tkerber
            x = energy / kT
77 1 tkerber
            S_v += x / (np.exp(x)-1.) - np.log(1-np.exp(-x))
78 1 tkerber
        S_v *= units.kB
79 1 tkerber
        return S_v
80 1 tkerber
81 1 tkerber
    def _vprint(self, text):
82 1 tkerber
        """Print output if verbose flag True."""
83 1 tkerber
        if self.verbose:
84 1 tkerber
            sys.stdout.write(text + os.linesep)
85 1 tkerber
86 1 tkerber
class HarmonicThermo(ThermoChem):
87 1 tkerber
    """Class for calculating thermodynamic properties in the approximation
88 1 tkerber
    that all degrees of freedom are treated harmonically. Often used for
89 1 tkerber
    adsorbates.
90 1 tkerber

91 1 tkerber
    Inputs:
92 1 tkerber

93 1 tkerber
    vib_energies : list
94 1 tkerber
        a list of the harmonic energies of the adsorbate (e.g., from
95 1 tkerber
        ase.vibrations.Vibrations.get_energies). The number of
96 1 tkerber
        energies should match the number of degrees of freedom of the
97 1 tkerber
        adsorbate; i.e., 3*n, where n is the number of atoms. Note that
98 1 tkerber
        this class does not check that the user has supplied the correct
99 1 tkerber
        number of energies. Units of energies are eV.
100 1 tkerber
    electronicenergy : float
101 1 tkerber
        the electronic energy in eV
102 1 tkerber
        (If the electronicenergy is unspecified, then the methods of this
103 1 tkerber
        class can be interpreted as the energy corrections.)
104 1 tkerber
    """
105 1 tkerber
106 1 tkerber
    def __init__(self, vib_energies, electronicenergy=None):
107 1 tkerber
        self.vib_energies = vib_energies
108 1 tkerber
        # Check for imaginary frequencies.
109 1 tkerber
        if sum(np.iscomplex(self.vib_energies)):
110 1 tkerber
            raise ValueError('Imaginary vibrational energies are present.')
111 1 tkerber
        else:
112 1 tkerber
            self.vib_energies = np.real(self.vib_energies) # clear +0.j
113 1 tkerber
114 1 tkerber
        if electronicenergy:
115 1 tkerber
            self.electronicenergy = electronicenergy
116 1 tkerber
        else:
117 1 tkerber
            self.electronicenergy = 0.
118 1 tkerber
119 1 tkerber
    def get_internal_energy(self, temperature, verbose=True):
120 1 tkerber
        """Returns the internal energy, in eV, in the harmonic approximation
121 1 tkerber
        at a specified temperature (K)."""
122 1 tkerber
123 1 tkerber
        self.verbose = verbose
124 1 tkerber
        write = self._vprint
125 1 tkerber
        fmt = '%-15s%13.3f eV'
126 1 tkerber
        write('Internal energy components at T = %.2f K:' % temperature)
127 1 tkerber
        write('='*31)
128 1 tkerber
129 1 tkerber
        U = 0.
130 1 tkerber
131 1 tkerber
        write(fmt % ('E_elec', self.electronicenergy))
132 1 tkerber
        U += self.electronicenergy
133 1 tkerber
134 1 tkerber
        zpe = self.get_ZPE_correction()
135 1 tkerber
        write(fmt % ('E_ZPE', zpe))
136 1 tkerber
        U += zpe
137 1 tkerber
138 1 tkerber
        dU_v = self._vibrational_energy_contribution(temperature)
139 1 tkerber
        write(fmt % ('Cv_harm (0->T)', dU_v))
140 1 tkerber
        U += dU_v
141 1 tkerber
142 1 tkerber
        write('-'*31)
143 1 tkerber
        write(fmt % ('U', U))
144 1 tkerber
        write('='*31)
145 1 tkerber
        return U
146 1 tkerber
147 1 tkerber
    def get_entropy(self, temperature, verbose=True):
148 1 tkerber
        """Returns the entropy, in eV/K, in the harmonic approximation
149 1 tkerber
        at a specified temperature (K)."""
150 1 tkerber
151 1 tkerber
        self.verbose = verbose
152 1 tkerber
        write = self._vprint
153 1 tkerber
        fmt = '%-15s%13.7f eV/K%13.3f eV'
154 1 tkerber
        write('Entropy components at T = %.2f K:' % temperature)
155 1 tkerber
        write('='*49)
156 1 tkerber
        write('%15s%13s     %13s'%('', 'S', 'T*S'))
157 1 tkerber
158 1 tkerber
        S = 0.
159 1 tkerber
160 1 tkerber
        S_v = self._vibrational_entropy_contribution(temperature)
161 1 tkerber
        write(fmt%('S_harm', S_v, S_v*temperature))
162 1 tkerber
        S += S_v
163 1 tkerber
164 1 tkerber
        write('-'*49)
165 1 tkerber
        write(fmt%('S', S, S*temperature))
166 1 tkerber
        write('='*49)
167 1 tkerber
        return S
168 1 tkerber
169 1 tkerber
    def get_free_energy(self, temperature, verbose=True):
170 1 tkerber
        """Returns the free energy, in eV, in the harmonic approximation
171 1 tkerber
        at a specified temperature (K)."""
172 1 tkerber
173 1 tkerber
        self.verbose = True
174 1 tkerber
        write = self._vprint
175 1 tkerber
176 1 tkerber
        U = self.get_internal_energy(temperature, verbose=verbose)
177 1 tkerber
        write('')
178 1 tkerber
        S = self.get_entropy(temperature, verbose=verbose)
179 1 tkerber
        G = U - temperature * S
180 1 tkerber
181 1 tkerber
        write('')
182 1 tkerber
        write('Free energy components at T = %.2f K:' % temperature)
183 1 tkerber
        write('='*23)
184 1 tkerber
        fmt = '%5s%15.3f eV'
185 1 tkerber
        write(fmt % ('U', U))
186 1 tkerber
        write(fmt % ('-T*S', -temperature * S))
187 1 tkerber
        write('-'*23)
188 1 tkerber
        write(fmt % ('G', G))
189 1 tkerber
        write('='*23)
190 1 tkerber
        return G
191 1 tkerber
192 1 tkerber
class IdealGasThermo(ThermoChem):
193 1 tkerber
    """Class for calculating thermodynamic properties of a molecule
194 1 tkerber
    based on statistical mechanical treatments in the ideal gas
195 1 tkerber
    approximation.
196 1 tkerber

197 1 tkerber
    Inputs for enthalpy calculations:
198 1 tkerber

199 1 tkerber
    vib_energies : list
200 1 tkerber
        a list of the vibrational energies of the molecule (e.g., from
201 1 tkerber
        ase.vibrations.Vibrations.get_energies). The number of vibrations
202 1 tkerber
        used is automatically calculated by the geometry and the number of
203 1 tkerber
        atoms. If more are specified than are needed, then the lowest
204 1 tkerber
        numbered vibrations are neglected. If either atoms or natoms is
205 1 tkerber
        unspecified, then uses the entire list. Units are eV.
206 1 tkerber
    geometry : 'monatomic', 'linear', or 'nonlinear'
207 1 tkerber
        geometry of the molecule
208 1 tkerber
    electronicenergy : float
209 1 tkerber
        the electronic energy in eV
210 1 tkerber
        (If electronicenergy is unspecified, then the methods of this
211 1 tkerber
        class can be interpreted as the enthalpy and free energy
212 1 tkerber
        corrections.)
213 1 tkerber
    natoms : integer
214 1 tkerber
        the number of atoms, used along with 'geometry' to determine how
215 1 tkerber
        many vibrations to use. (Not needed if an atoms object is supplied
216 1 tkerber
        in 'atoms' or if the user desires the entire list of vibrations
217 1 tkerber
        to be used.)
218 1 tkerber

219 1 tkerber
    Extra inputs needed for for entropy / free energy calculations:
220 1 tkerber

221 1 tkerber
    atoms : an ASE atoms object
222 1 tkerber
        used to calculate rotational moments of inertia and molecular mass
223 1 tkerber
    symmetrynumber : integer
224 1 tkerber
        symmetry number of the molecule. See, for example, Table 10.1 and
225 1 tkerber
        Appendix B of C. Cramer "Essentials of Computational Chemistry",
226 1 tkerber
        2nd Ed.
227 1 tkerber
    spin : float
228 1 tkerber
        the total electronic spin. (0 for molecules in which all electrons
229 1 tkerber
        are paired, 0.5 for a free radical with a single unpaired electron,
230 1 tkerber
        1.0 for a triplet with two unpaired electrons, such as O_2.)
231 1 tkerber
    """
232 1 tkerber
233 1 tkerber
    def __init__(self, vib_energies, geometry, electronicenergy=None,
234 1 tkerber
                 atoms=None, symmetrynumber=None, spin=None, natoms=None):
235 1 tkerber
        if electronicenergy == None:
236 1 tkerber
            self.electronicenergy = 0.
237 1 tkerber
        else:
238 1 tkerber
            self.electronicenergy = electronicenergy
239 1 tkerber
        self.geometry = geometry
240 1 tkerber
        self.atoms = atoms
241 1 tkerber
        self.sigma = symmetrynumber
242 1 tkerber
        self.spin = spin
243 1 tkerber
        if natoms == None:
244 1 tkerber
            if atoms:
245 1 tkerber
                natoms = len(atoms)
246 1 tkerber
        # Cut the vibrations to those needed from the geometry.
247 1 tkerber
        if natoms:
248 1 tkerber
            if geometry == 'nonlinear':
249 1 tkerber
                self.vib_energies = vib_energies[-(3*natoms-6):]
250 1 tkerber
            elif geometry == 'linear':
251 1 tkerber
                self.vib_energies = vib_energies[-(3*natoms-5):]
252 1 tkerber
            elif geometry == 'monatomic':
253 1 tkerber
                self.vib_energies = []
254 1 tkerber
        else:
255 1 tkerber
            self.vib_energies = vib_energies
256 1 tkerber
        # Make sure no imaginary frequencies remain.
257 1 tkerber
        if sum(np.iscomplex(self.vib_energies)):
258 1 tkerber
            raise ValueError('Imaginary frequencies are present.')
259 1 tkerber
        else:
260 1 tkerber
            self.vib_energies = np.real(self.vib_energies) # clear +0.j
261 1 tkerber
        self.referencepressure = 101325. # Pa
262 1 tkerber
263 1 tkerber
    def get_enthalpy(self, temperature, verbose=True):
264 1 tkerber
        """Returns the enthalpy, in eV, in the ideal gas approximation
265 1 tkerber
        at a specified temperature (K)."""
266 1 tkerber
267 1 tkerber
        self.verbose = verbose
268 1 tkerber
        write = self._vprint
269 1 tkerber
        fmt = '%-15s%13.3f eV'
270 1 tkerber
        write('Enthalpy components at T = %.2f K:' % temperature)
271 1 tkerber
        write('='*31)
272 1 tkerber
273 1 tkerber
        H = 0.
274 1 tkerber
275 1 tkerber
        write(fmt % ('E_elec', self.electronicenergy))
276 1 tkerber
        H += self.electronicenergy
277 1 tkerber
278 1 tkerber
        zpe = self.get_ZPE_correction()
279 1 tkerber
        write(fmt % ('E_ZPE', zpe))
280 1 tkerber
        H += zpe
281 1 tkerber
282 1 tkerber
        Cv_t = 3./2. * units.kB # translational heat capacity (3-d gas)
283 1 tkerber
        write(fmt % ('Cv_trans (0->T)', Cv_t * temperature))
284 1 tkerber
        H += Cv_t * temperature
285 1 tkerber
286 1 tkerber
        if self.geometry == 'nonlinear': # rotational heat capacity
287 1 tkerber
            Cv_r = 3./2.*units.kB
288 1 tkerber
        elif self.geometry == 'linear':
289 1 tkerber
            Cv_r = units.kB
290 1 tkerber
        elif self.geometry == 'monatomic':
291 1 tkerber
            Cv_r = 0.
292 1 tkerber
        write(fmt % ('Cv_rot (0->T)', Cv_r * temperature))
293 1 tkerber
        H += Cv_r*temperature
294 1 tkerber
295 1 tkerber
        dH_v = self._vibrational_energy_contribution(temperature)
296 1 tkerber
        write(fmt % ('Cv_vib (0->T)', dH_v))
297 1 tkerber
        H += dH_v
298 1 tkerber
299 1 tkerber
        Cp_corr = units.kB * temperature
300 1 tkerber
        write(fmt % ('(C_v -> C_p)', Cp_corr))
301 1 tkerber
        H += Cp_corr
302 1 tkerber
303 1 tkerber
        write('-'*31)
304 1 tkerber
        write(fmt % ('H', H))
305 1 tkerber
        write('='*31)
306 1 tkerber
        return H
307 1 tkerber
308 1 tkerber
    def get_entropy(self, temperature, pressure, verbose=True):
309 1 tkerber
        """Returns the entropy, in eV/K, in the ideal gas approximation
310 1 tkerber
        at a specified temperature (K) and pressure (Pa)."""
311 1 tkerber
312 1 tkerber
        if self.atoms == None or self.sigma == None or self.spin == None:
313 1 tkerber
            raise RuntimeError('atoms, symmetrynumber, and spin must be '
314 1 tkerber
                               'specified for entropy and free energy '
315 1 tkerber
                               'calculations.')
316 1 tkerber
        self.verbose = verbose
317 1 tkerber
        write = self._vprint
318 1 tkerber
        fmt = '%-15s%13.7f eV/K%13.3f eV'
319 1 tkerber
        write('Entropy components at T = %.2f K and P = %.1f Pa:' %
320 1 tkerber
               (temperature, pressure))
321 1 tkerber
        write('='*49)
322 1 tkerber
        write('%15s%13s     %13s' % ('', 'S', 'T*S'))
323 1 tkerber
324 1 tkerber
        S = 0.
325 1 tkerber
326 1 tkerber
        # Translational entropy (term inside the log is in SI units).
327 1 tkerber
        mass = sum(self.atoms.get_masses()) * units._amu # kg/molecule
328 1 tkerber
        S_t = (2*np.pi*mass*units._k*temperature/units._hplanck**2)**(3./2)
329 1 tkerber
        S_t *= units._k * temperature / self.referencepressure
330 1 tkerber
        S_t = units.kB * (np.log(S_t) + 5./2.)
331 1 tkerber
        write(fmt % ('S_trans (1 atm)', S_t, S_t * temperature))
332 1 tkerber
        S += S_t
333 1 tkerber
334 1 tkerber
        # Rotational entropy (term inside the log is in SI units).
335 1 tkerber
        if self.geometry == 'monatomic':
336 1 tkerber
            S_r = 0.
337 1 tkerber
        elif self.geometry == 'nonlinear':
338 1 tkerber
            inertias = (rotationalinertia(self.atoms) * units._amu /
339 1 tkerber
                        (10.**10)**2) # kg m^2
340 1 tkerber
            S_r = np.sqrt(np.pi*np.product(inertias)) / self.sigma
341 1 tkerber
            S_r *= (8. * np.pi**2 * units._k * temperature /
342 1 tkerber
                    units._hplanck**2)**(3./2.)
343 1 tkerber
            S_r = units.kB * (np.log(S_r) + 3./2.)
344 1 tkerber
        elif self.geometry == 'linear':
345 1 tkerber
            inertias = (rotationalinertia(self.atoms) * units._amu /
346 1 tkerber
                        (10.**10)**2) # kg m^2
347 1 tkerber
            inertia = max(inertias) # should be two identical and one zero
348 1 tkerber
            S_r = (8 * np.pi**2 * inertia * units._k * temperature /
349 1 tkerber
                   self.sigma / units._hplanck**2)
350 1 tkerber
            S_r = units.kB * ( np.log(S_r) + 1.)
351 1 tkerber
        write(fmt % ('S_rot', S_r, S_r*temperature))
352 1 tkerber
        S += S_r
353 1 tkerber
354 1 tkerber
        # Electronic entropy.
355 1 tkerber
        S_e = units.kB * np.log(2*self.spin+1)
356 1 tkerber
        write(fmt % ('S_elec', S_e, S_e * temperature))
357 1 tkerber
        S += S_e
358 1 tkerber
359 1 tkerber
        # Vibrational entropy.
360 1 tkerber
        S_v = self._vibrational_entropy_contribution(temperature)
361 1 tkerber
        write(fmt % ('S_vib', S_v, S_v* temperature))
362 1 tkerber
        S += S_v
363 1 tkerber
364 1 tkerber
        # Pressure correction to translational entropy.
365 1 tkerber
        S_p = - units.kB * np.log(pressure/self.referencepressure)
366 1 tkerber
        write(fmt % ('S (1 atm -> P)', S_p, S_p * temperature))
367 1 tkerber
        S += S_p
368 1 tkerber
369 1 tkerber
        write('-'*49)
370 1 tkerber
        write(fmt % ('S', S, S * temperature))
371 1 tkerber
        write('='*49)
372 1 tkerber
        return S
373 1 tkerber
374 1 tkerber
    def get_free_energy(self, temperature, pressure, verbose=True):
375 1 tkerber
        """Returns the free energy, in eV, in the ideal gas approximation
376 1 tkerber
        at a specified temperature (K) and pressure (Pa)."""
377 1 tkerber
378 1 tkerber
        self.verbose = verbose
379 1 tkerber
        write = self._vprint
380 1 tkerber
381 1 tkerber
        H = self.get_enthalpy(temperature, verbose=verbose)
382 1 tkerber
        write('')
383 1 tkerber
        S = self.get_entropy(temperature, pressure, verbose=verbose)
384 1 tkerber
        G = H - temperature * S
385 1 tkerber
386 1 tkerber
        write('')
387 1 tkerber
        write('Free energy components at T = %.2f K and P = %.1f Pa:' %
388 1 tkerber
               (temperature, pressure))
389 1 tkerber
        write('='*23)
390 1 tkerber
        fmt = '%5s%15.3f eV'
391 1 tkerber
        write(fmt % ('H', H))
392 1 tkerber
        write(fmt % ('-T*S', -temperature * S))
393 1 tkerber
        write('-'*23)
394 1 tkerber
        write(fmt % ('G', G))
395 1 tkerber
        write('='*23)
396 1 tkerber
        return G