2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2017,2018,2019, 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 zero point energy from an array of eigenvalues.
51 * This routine first converts the eigenvalues from a normal mode
52 * analysis to frequencies and then computes the zero point energy.
54 * \param[in] eigval The eigenvalues
55 * \param[in] scale_factor Factor to scale frequencies by before computing cv
56 * \return The zero point energy (kJ/mol)
58 double calcZeroPointEnergy(gmx::ArrayRef<const real> eigval, real scale_factor);
60 /*! \brief Compute heat capacity due to vibrational motion
62 * \param[in] eigval The eigenvalues
63 * \param[in] temperature Temperature (K)
64 * \param[in] linear TRUE if this is a linear molecule
65 * \param[in] scale_factor Factor to scale frequencies by before computing cv
66 * \return The heat capacity at constant volume (J/mol K)
68 double calcVibrationalHeatCapacity(gmx::ArrayRef<const real> eigval,
73 /*! \brief Compute entropy due to translational motion
75 * Following the equations in J. W. Ochterski,
76 * Thermochemistry in Gaussian, Gaussian, Inc., 2000
79 * \param[in] mass Molecular mass (Dalton)
80 * \param[in] temperature Temperature (K)
81 * \param[in] pressure Pressure (bar) at which to compute
82 * \returns The translational entropy (J/mol K)
84 double calcTranslationalEntropy(real mass, real temperature, real pressure);
86 /*! \brief Compute entropy due to rotational motion
88 * Following the equations in J. W. Ochterski,
89 * Thermochemistry in Gaussian, Gaussian, Inc., 2000
92 * \param[in] temperature Temperature (K)
93 * \param[in] natom Number of atoms
94 * \param[in] linear TRUE if this is a linear molecule
95 * \param[in] theta The principal moments of inertia (unit of Energy)
96 * \param[in] sigma_r Symmetry factor, should be >= 1
97 * \returns The rotational entropy (J/mol K)
99 double calcRotationalEntropy(real temperature, int natom, gmx_bool linear, const rvec theta, real sigma_r);
101 /*! \brief Compute internal energy due to vibrational motion
103 * \param[in] eigval The eigenvalues
104 * \param[in] temperature Temperature (K)
105 * \param[in] linear TRUE if this is a linear molecule
106 * \param[in] scale_factor Factor to scale frequencies by before computing E
107 * \return The internal energy (J/mol K)
109 double calcVibrationalInternalEnergy(gmx::ArrayRef<const real> eigval,
114 /*! \brief Compute entropy using Schlitter formula
116 * Computes entropy for a molecule / molecular system using the
117 * algorithm due to Schlitter (Chem. Phys. Lett. 215 (1993)
119 * The input should be eigenvalues from a covariance analysis,
120 * the units of the eigenvalues are those of energy.
122 * \param[in] eigval The eigenvalues
123 * \param[in] temperature Temperature (K)
124 * \param[in] linear True if this is a linear molecule (typically a diatomic molecule).
125 * \return the entropy (J/mol K)
127 double calcSchlitterEntropy(gmx::ArrayRef<const real> eigval, real temperature, gmx_bool linear);
129 /*! \brief Compute entropy using Quasi-Harmonic formula
131 * Computes entropy for a molecule / molecular system using the
132 * Quasi-harmonic algorithm (Macromolecules 1984, 17, 1370).
133 * The input should be eigenvalues from a normal mode analysis.
134 * In both cases the units of the eigenvalues are those of energy.
136 * \param[in] eigval The eigenvalues
137 * \param[in] temperature Temperature (K)
138 * \param[in] linear True if this is a linear molecule (typically a diatomic molecule).
139 * \param[in] scale_factor Factor to scale frequencies by before computing S0
140 * \return the entropy (J/mol K)
142 double calcQuasiHarmonicEntropy(gmx::ArrayRef<const real> eigval, real temperature, gmx_bool linear, real scale_factor);