2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2017,2018, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
37 * Code for computing entropy and heat capacity from eigenvalues
39 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
41 #ifndef GMXANA_THERMOCHEMISTRY_H
42 #define GMXANA_THERMOCHEMISTRY_H
44 #include "gromacs/math/units.h"
45 #include "gromacs/math/vec.h"
46 #include "gromacs/utility/arrayref.h"
47 #include "gromacs/utility/basedefinitions.h"
49 /*! \brief Compute heat capacity due to vibrational motion
51 * \param[in] eigval The eigenvalues
52 * \param[in] temperature Temperature (K)
53 * \param[in] linear TRUE if this is a linear molecule
54 * \param[in] scale_factor Factor to scale frequencies by before computing cv
55 * \return The heat capacity at constant volume
57 double calcVibrationalHeatCapacity(gmx::ArrayRef<const real> eigval,
62 /*! \brief Compute entropy due to translational motion
64 * Following the equations in J. W. Ochterski,
65 * Thermochemistry in Gaussian, Gaussian, Inc., 2000
68 * \param[in] mass Molecular mass (Dalton)
69 * \param[in] temperature Temperature (K)
70 * \param[in] pressure Pressure (bar) at which to compute
71 * \returns The translational entropy
73 double calcTranslationalEntropy(real mass,
77 /*! \brief Compute entropy due to rotational motion
79 * Following the equations in J. W. Ochterski,
80 * Thermochemistry in Gaussian, Gaussian, Inc., 2000
83 * \param[in] temperature Temperature (K)
84 * \param[in] natom Number of atoms
85 * \param[in] linear TRUE if this is a linear molecule
86 * \param[in] theta The principal moments of inertia (unit of Energy)
87 * \param[in] sigma_r Symmetry factor, should be >= 1
88 * \returns The rotational entropy
90 double calcRotationalEntropy(real temperature,
96 /*! \brief Compute internal energy due to vibrational motion
98 * \param[in] eigval The eigenvalues
99 * \param[in] temperature Temperature (K)
100 * \param[in] linear TRUE if this is a linear molecule
101 * \param[in] scale_factor Factor to scale frequencies by before computing E
102 * \return The internal energy
104 double calcVibrationalInternalEnergy(gmx::ArrayRef<const real> eigval,
109 /*! \brief Compute entropy using Schlitter formula
111 * Computes entropy for a molecule / molecular system using the
112 * algorithm due to Schlitter (Chem. Phys. Lett. 215 (1993)
114 * The input should be eigenvalues from a covariance analysis,
115 * the units of the eigenvalues are those of energy.
117 * \param[in] eigval The eigenvalues
118 * \param[in] temperature Temperature (K)
119 * \param[in] linear True if this is a linear molecule (typically a diatomic molecule).
120 * \return the entropy (J/mol K)
122 double calcSchlitterEntropy(gmx::ArrayRef<const real> eigval,
126 /*! \brief Compute entropy using Quasi-Harmonic formula
128 * Computes entropy for a molecule / molecular system using the
129 * Quasi-harmonic algorithm (Macromolecules 1984, 17, 1370).
130 * The input should be eigenvalues from a normal mode analysis.
131 * In both cases the units of the eigenvalues are those of energy.
133 * \param[in] eigval The eigenvalues
134 * \param[in] temperature Temperature (K)
135 * \param[in] linear True if this is a linear molecule (typically a diatomic molecule).
136 * \param[in] scale_factor Factor to scale frequencies by before computing S0
137 * \return the entropy (J/mol K)
139 double calcQuasiHarmonicEntropy(gmx::ArrayRef<const real> eigval,