Thermochemistry

ASE contains a ase.thermochemistry module that lets the user derive commonly desired thermodynamic quantities of molecules and crystalline solids from ASE output and some user-specified parameters. Three cases are currently handled by this module: the ideal-gas limit (in which translational and rotational degrees of freedom are taken into account), the harmonic limit (generally used for adsorbates, in which all degrees of freedom are treated harmonically), and a crystalline solid model (in which a lattice of N atoms is treated as a system of 3N independent harmonic oscillators). The first two cases rely on good vibrational energies being fed to the calculators, which can be calculated with the ase.vibrations module. Likewise, the crystalline solid model depends on an accurate phonon density of states; this is readily calculated using the ase.phonons module.

Ideal-gas limit

The thermodynamic quantities of ideal gases are calculated by assuming that all spatial degrees of freedom are independent and separable into translational, rotational, and vibrational degrees of freedom. The IdealGasThermo class supports calculation of enthalpy (\(H\)), entropy (\(S\)), and Gibbs free energy (\(G\)), and has the interface listed below.

class ase.thermochemistry.IdealGasThermo(vib_energies, geometry, electronicenergy=None, atoms=None, symmetrynumber=None, spin=None, natoms=None)[source]

Class for calculating thermodynamic properties of a molecule based on statistical mechanical treatments in the ideal gas approximation.

Inputs for enthalpy calculations:

vib_energies : list
a list of the vibrational energies of the molecule (e.g., from ase.vibrations.Vibrations.get_energies). The number of vibrations used is automatically calculated by the geometry and the number of atoms. If more are specified than are needed, then the lowest numbered vibrations are neglected. If either atoms or natoms is unspecified, then uses the entire list. Units are eV.
geometry : ‘monatomic’, ‘linear’, or ‘nonlinear’
geometry of the molecule
electronicenergy : float
the electronic energy in eV (if electronicenergy is unspecified, then the methods of this class can be interpreted as the enthalpy and free energy corrections)
natoms : integer
the number of atoms, used along with ‘geometry’ to determine how many vibrations to use. (Not needed if an atoms object is supplied in ‘atoms’ or if the user desires the entire list of vibrations to be used.)

Extra inputs needed for for entropy / free energy calculations:

atoms : an ASE atoms object
used to calculate rotational moments of inertia and molecular mass
symmetrynumber : integer
symmetry number of the molecule. See, for example, Table 10.1 and Appendix B of C. Cramer “Essentials of Computational Chemistry”, 2nd Ed.
spin : float
the total electronic spin. (0 for molecules in which all electrons are paired, 0.5 for a free radical with a single unpaired electron, 1.0 for a triplet with two unpaired electrons, such as O_2.)
get_enthalpy(temperature, verbose=True)[source]

Returns the enthalpy, in eV, in the ideal gas approximation at a specified temperature (K).

get_entropy(temperature, pressure, verbose=True)[source]

Returns the entropy, in eV/K, in the ideal gas approximation at a specified temperature (K) and pressure (Pa).

get_gibbs_energy(temperature, pressure, verbose=True)[source]

Returns the Gibbs free energy, in eV, in the ideal gas approximation at a specified temperature (K) and pressure (Pa).

Example

The IdealGasThermo class would generally be called after an energy optimization and a vibrational analysis. The user needs to supply certain parameters if the entropy or free energy are desired, such as the geometry and symmetry number. An example on the nitrogen molecule is:

from ase.structure import molecule
from ase.calculators.emt import EMT
from ase.optimize import QuasiNewton
from ase.vibrations import Vibrations
from ase.thermochemistry import IdealGasThermo

atoms = molecule('N2')
atoms.set_calculator(EMT())
dyn = QuasiNewton(atoms)
dyn.run(fmax=0.01)
electronicenergy = atoms.get_potential_energy()

vib = Vibrations(atoms)
vib.run()
vib_energies = vib.get_energies()

thermo = IdealGasThermo(vib_energies=vib_energies,
                        electronicenergy=electronicenergy,
                        atoms=atoms,
                        geometry='linear',
                        symmetrynumber=2, spin=0)
G = thermo.get_gibbs_energy(temperature=298.15, pressure=101325.)

This will give the thermodynamic summary output:

Enthalpy components at T = 298.15 K:
===============================
E_elec                 0.263 eV
E_ZPE                  0.076 eV
Cv_trans (0->T)        0.039 eV
Cv_rot (0->T)          0.026 eV
Cv_vib (0->T)          0.000 eV
(C_v -> C_p)           0.026 eV
-------------------------------
H                      0.429 eV
===============================

Entropy components at T = 298.15 K and P = 101325.0 Pa:
=================================================
                           S               T*S
S_trans (1 atm)    0.0015579 eV/K        0.464 eV
S_rot              0.0004101 eV/K        0.122 eV
S_elec             0.0000000 eV/K        0.000 eV
S_vib              0.0000016 eV/K        0.000 eV
S (1 atm -> P)    -0.0000000 eV/K       -0.000 eV
-------------------------------------------------
S                  0.0019695 eV/K        0.587 eV
=================================================

Free energy components at T = 298.15 K and P = 101325.0 Pa:
=======================
    H          0.429 eV
 -T*S         -0.587 eV
-----------------------
    G         -0.158 eV
=======================

Harmonic limit

In the harmonic limit, all degrees of freedom are treated harmonically. The HarmonicThermo class supports the calculation of internal energy, entropy, and Gibbs free energy. This class uses all of the energies given to it in the vib_energies list; this is a list as can be generated with the .get_energies() method of ase.vibrations.Vibrations, but the user should take care that all of these energies are real (non-imaginary). The class HarmonicThermo has the interface described below.

class ase.thermochemistry.HarmonicThermo(vib_energies, electronicenergy=None)[source]

Class for calculating thermodynamic properties in the approximation that all degrees of freedom are treated harmonically. Often used for adsorbates.

Inputs:

vib_energies : list
a list of the harmonic energies of the adsorbate (e.g., from ase.vibrations.Vibrations.get_energies). The number of energies should match the number of degrees of freedom of the adsorbate; i.e., 3*n, where n is the number of atoms. Note that this class does not check that the user has supplied the correct number of energies. Units of energies are eV.
electronicenergy : float
the electronic energy in eV (if electronicenergy is unspecified, then the methods of this class can be interpreted as the energy corrections)
get_entropy(temperature, verbose=True)[source]

Returns the entropy, in eV/K, in the harmonic approximation at a specified temperature (K).

get_gibbs_energy(temperature, verbose=True)[source]

Returns the Gibbs free energy, in eV, in the harmonic approximation at a specified temperature (K).

get_internal_energy(temperature, verbose=True)[source]

Returns the internal energy, in eV, in the harmonic approximation at a specified temperature (K).

Crystals

In this model a crystalline solid is treated as a periodic system of independent harmonic oscillators. The CrystalThermo class supports the calculation of internal energy (\(U\)), entropy (\(S\)) and Helmholtz free energy (\(F\)), and has the interface listed below.

class ase.thermochemistry.CrystalThermo(phonon_DOS, phonon_energies, formula_units=None, electronicenergy=None)[source]

Class for calculating thermodynamic properties of a crystalline solid in the approximation that a lattice of N atoms behaves as a system of 3N independent harmonic oscillators.

Inputs:

phonon_DOS : list
a list of the phonon density of states, where each value represents the phonon DOS at the vibrational energy value of the corresponding index in phonon_energies.
phonon_energies : list
a list of the range of vibrational energies (hbar*omega) over which the phonon density of states has been evaluated. This list should be the same length as phonon_DOS and integrating phonon_DOS over phonon_energies should yield approximately 3N, where N is the number of atoms per unit cell. If the first element of this list is zero-valued it will be deleted along with the first element of phonon_DOS. Units of vibrational energies are eV.
electronicenergy : float
the electronic energy in eV (if electronicenergy is unspecified, then the methods of this class can be interpreted as the phonon energy corrections.)
formula_units : int
the number of formula units per unit cell. If unspecified, the thermodynamic quantities calculated will be listed on a per-unit-cell basis.
get_entropy(temperature, verbose=True)[source]

Returns the entropy, in eV/K, of crystalline solid at a specified temperature (K).

get_helmholtz_energy(temperature, verbose=True)[source]

Returns the Helmholtz free energy, in eV, of crystalline solid at a specified temperature (K).

get_internal_energy(temperature, verbose=True)[source]

Returns the internal energy, in eV, of crystalline solid at a specified temperature (K).

Example

The CrystalThermo class will generally be called after an energy optimization and a phonon vibrational analysis of the crystal. An example for bulk gold is:

from ase.lattice.spacegroup import crystal
from ase.calculators.emt import EMT
from ase.optimize import QuasiNewton
from ase.phonons import Phonons
from ase.thermochemistry import CrystalThermo

# Set up gold bulk and attach EMT calculator
a = 4.078
atoms = crystal('Au', (0.,0.,0.),
                spacegroup = 225,
                cellpar = [a, a, a, 90, 90, 90],
                pbc = (1, 1, 1))
calc = EMT()
atoms.set_calculator(calc)
qn = QuasiNewton(atoms)
qn.run(fmax = 0.05)
electronicenergy = atoms.get_potential_energy()

# Phonon analysis
N = 5
ph = Phonons(atoms, calc, supercell=(N, N, N), delta=0.05)
ph.run()
ph.read(acoustic=True)
phonon_energies, phonon_DOS = ph.dos(kpts=(40, 40, 40), npts=3000,
                                     delta=5e-4)

# Calculate the Helmholtz free energy
thermo = CrystalThermo(phonon_energies=phonon_energies,
                       phonon_DOS = phonon_DOS,
                       electronicenergy = electronicenergy,
                       formula_units = 4)
F = thermo.get_helmholtz_energy(temperature=298.15)

This will give the thermodynamic summary output:

Internal energy components at T = 298.15 K,
on a per-formula-unit basis:
===============================
E_elec                0.0022 eV
E_ZPE                 0.0135 eV
E_phonon              0.0629 eV
-------------------------------
U                     0.0786 eV
===============================

Entropy components at T = 298.15 K,
on a per-formula-unit basis:
=================================================
                           S               T*S
-------------------------------------------------
S                  0.0005316 eV/K       0.1585 eV
=================================================

Helmholtz free energy components at T = 298.15 K,
on a per-formula-unit basis:
=======================
    U         0.0786 eV
 -T*S        -0.1585 eV
-----------------------
    F        -0.0799 eV
=======================

Background

Ideal gas. The conversion of electronic structure calculations to thermodynamic properties in the ideal-gas limit is well documented; see, for example, Chapter 10 of Cramer, 2004. The key equations used in the IdealGasThermo class are summarized here.

C.J. Cramer. Essentials of Computational Chemistry, Second Edition. Wiley, 2004.

The ideal-gas enthalpy is calculated from extrapolation of the energy at 0 K to the relevant temperature (for an ideal gas, the enthalpy is not a function of pressure):

\[H(T) = E_\text{elec} + E_\text{ZPE} + \int_0^\text{T} C_P \, \text{d}T\]

where the first two terms are the electronic energy and the zero-point energy, and the integral is over the constant-pressure heat capacity. The heat capacity is separable into translational, rotational, vibrational, and electronic parts (plus a term of \(k_\text{B}\) to switch from constant-volume to constant-pressure):

\[C_P = k_\text{B} + C_{V\text{,trans}} + C_{V\text{,rot}} + C_{V\text{,vib}} + C_{V\text{,elec}}\]

The translational heat capacity is 3/2 \(k_\text{B}\) for a 3-dimensional gas. The rotational heat capacity is 0 for a monatomic species, \(k_\text{B}\) for a linear molecule, and 3/2 \(k_\text{B}\) for a nonlinear molecule. In this module, the electronic component of the heat capacity is assumed to be 0. The vibrational heat capacity contains \(3N-6\) degrees of freedom for nonlinear molecules and \(3N-5\) degrees of freedom for linear molecules (where \(N\) is the number of atoms). The integrated form of the vibrational heat capacity is:

\[\int_0^T C_{V,\text{vib}} \text{d}T = \sum_i^\text{vib DOF} \frac{\epsilon_i}{e^{\epsilon_i / k_\text{B} T} - 1 }\]

where \(\epsilon_i\) are the energies associated with the vibrational frequencies, \(\epsilon_i = h \omega_i\).

The ideal gas entropy can be calculated as a function of temperature and pressure as:

\[\begin{split}S(T,P) &= S(T,P^\circ) - k_\text{B} \ln \frac{P}{P^\circ} \\ &= S_\text{trans} + S_\text{rot} + S_\text{elec} + S_\text{vib} - k_\text{B} \ln \frac{P}{P^\circ}\end{split}\]

where the translational, rotational, electronic, and vibrational components are calculated as below. (Note that the translational component also includes components from the Stirling approximation, and that the vibrational degrees of freedom are enumerated the same as in the above.)

\[S_\text{trans} = k_\text{B} \left\{ \ln \left[ \left( \frac{2 \pi M k_\text{B} T}{h^2} \right)^{3/2} \frac{k_\text{B} T}{P^\circ} \right] + \frac{5}{2} \right\}\]
\[\begin{split}S_\text{rot} = \left\{ \begin{array}{ll} 0 & \text{, if monatomic} \\ k_\text{B} \left[ \ln \left( \frac{8\pi^2 I k_\text{B}T}{\sigma h^2}\right) + 1 \right] & \text{, if linear} \\ k_\text{B} \left\{ \ln \left[ \frac{\sqrt{\pi I_\text{A} I_\text{B} I_\text{C}}}{\sigma} \left(\frac{8\pi^2 k_\text{B} T}{h^2}\right)^{3/2}\right] + \frac{3}{2} \right\} & \text{, if nonlinear} \\ \end{array} \right.\end{split}\]
\[S_\text{vib} = k_\text{B} \sum_i^\text{vib DOF} \left[ \frac{\epsilon_i}{k_\text{B}T\left(e^{\epsilon_i/k_\text{B}T}-1\right)} - \ln \left( 1 - e^{-\epsilon_i/k_\text{B}T} \right)\right]\]
\[S_\text{elec} = k_\text{B} \ln \left[ 2 \times \left(\text{spin multiplicity}\right) + 1\right]\]

\(I_\text{A}\) through \(I_\text{C}\) are the three principle moments of inertia for a non-linear molecule. \(I\) is the degenerate moment of inertia for a linear molecule. \(\sigma\) is the symmetry number of the molecule.

The ideal-gas Gibbs free energy is then just calculated from the combination of the enthalpy and entropy:

\[G(T,P) = H(T) - T\, S(T,P)\]

Harmonic limit. The conversion of electronic structure calculation information into thermodynamic properties is less established for adsorbates. However, the simplest approach often taken is to treat all \(3N\) degrees of freedom of the adsorbate harmonically since the adsorbate often has no real translational or rotational degrees of freedom. This is the approach implemented in the HarmonicThermo class. Thus, the internal energy and entropy of the adsorbate are calculated as

\[U(T) = E_\text{elec} + E_\text{ZPE} + \sum_i^\text{harm DOF} \frac{\epsilon_i}{e^{\epsilon_i / k_\text{B} T} - 1 }\]
\[S = k_\text{B} \sum_i^\text{harm DOF} \left[ \frac{\epsilon_i}{k_\text{B}T\left(e^{\epsilon_i/k_\text{B}T}-1\right)} - \ln \left( 1 - e^{-\epsilon_i/k_\text{B}T} \right)\right]\]

and the Gibbs free energy is calculated as

\[G(T) = U(T) - T\, S(T)\]

In this case, the number of harmonic energies (\(\epsilon_i\)) used in the summation is generally \(3N\), where \(N\) is the number of atoms in the adsorbate.

Crystalline solid

The derivation of the partition function for a crystalline solid is fairly straight-forward and can be found, for example, in Chapter 11 of McQuarrie, 2000.

D.A. McQuarrie. Statistical Mechanics. University Science Books, 2000.

The treatment implemented in the CrystalThermo class depends on introducing normal coordinates to the entire crystal and treating each atom in the lattice as an independent harmonic oscillator. This yields the partition function

\[Z = \prod_{j=1}^\text{3N} \left( \frac{e^{-\frac{1}{2}\hbar\omega_j\beta}}{1 - e^{-\hbar\omega_j\beta}} \right) e^{-E_\text{elec} \beta}\]

where \(\omega_j\) are the \(3N\) vibrational frequencies, \(E_\text{elec}\) is the electronic energy of the crystalline solid, and \(\beta = \frac{1}{k_\text{B} T}\). Now, taking the logarithm of the partition function and replacing the resulting sum with an integral (assuming that the energy level spacing is essentially continuous) gives

\[-\ln Z = E_\text{elec}\beta + \int_0^\infty \left[ \ln \left( 1 - e^{-\hbar\omega\beta} \right) + \frac{\hbar\omega\beta}{2} \right]\sigma (\omega) \text{d}\omega\]

Here \(\sigma (\omega)\) represents the degeneracy or phonon density of states as a function of vibrational frequency. Once this function has been determined (i.e. using the ase.phonons module), it is a simple matter to calculate the canonical ensemble thermodynamic quantities; namely the internal energy, the entropy and the Helmholtz free energy.

\[\begin{split}U(T) &= -\left( \frac{\partial \ln Z}{\partial \beta} \right)_\text{N,V} \\ &= E_\text{elec} + \int_0^\infty \left[ \frac{\hbar \omega}{e^{\hbar \omega \beta} - 1} + \frac{\hbar \omega}{2} \right]\sigma (\omega) \text{d}\omega\end{split}\]
\[\begin{split}S(T) &= \frac{U}{T} + k_\text{B} \ln Z \\ &= \int_0^\infty \left[ \frac{\hbar \omega}{T} \frac{1}{e^{\hbar \omega \beta} - 1} - k_\text{B} \ln \left(1 - e^{-\hbar \omega \beta} \right) \right]\sigma (\omega) \text{d}\omega\end{split}\]
\[F(T) = U(T) - T\, S(T,P)\]