Introduce SimulatorBuilder
[gromacs.git] / src / gromacs / gmxana / thermochemistry.h
blob5d8010563c020ea3ca9f78967a6fbbae76c6ca54
1 /*
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.
35 /*! \internal \file
36 * \brief
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,
59 real scale_factor);
61 /*! \brief Compute heat capacity due to vibrational motion
63 * \param[in] eigval The eigenvalues
64 * \param[in] temperature Temperature (K)
65 * \param[in] linear TRUE if this is a linear molecule
66 * \param[in] scale_factor Factor to scale frequencies by before computing cv
67 * \return The heat capacity at constant volume (J/mol K)
69 double calcVibrationalHeatCapacity(gmx::ArrayRef<const real> eigval,
70 real temperature,
71 gmx_bool linear,
72 real scale_factor);
74 /*! \brief Compute entropy due to translational motion
76 * Following the equations in J. W. Ochterski,
77 * Thermochemistry in Gaussian, Gaussian, Inc., 2000
78 * Pitssburg PA
80 * \param[in] mass Molecular mass (Dalton)
81 * \param[in] temperature Temperature (K)
82 * \param[in] pressure Pressure (bar) at which to compute
83 * \returns The translational entropy (J/mol K)
85 double calcTranslationalEntropy(real mass,
86 real temperature,
87 real pressure);
89 /*! \brief Compute entropy due to rotational motion
91 * Following the equations in J. W. Ochterski,
92 * Thermochemistry in Gaussian, Gaussian, Inc., 2000
93 * Pitssburg PA
95 * \param[in] temperature Temperature (K)
96 * \param[in] natom Number of atoms
97 * \param[in] linear TRUE if this is a linear molecule
98 * \param[in] theta The principal moments of inertia (unit of Energy)
99 * \param[in] sigma_r Symmetry factor, should be >= 1
100 * \returns The rotational entropy (J/mol K)
102 double calcRotationalEntropy(real temperature,
103 int natom,
104 gmx_bool linear,
105 const rvec theta,
106 real sigma_r);
108 /*! \brief Compute internal energy due to vibrational motion
110 * \param[in] eigval The eigenvalues
111 * \param[in] temperature Temperature (K)
112 * \param[in] linear TRUE if this is a linear molecule
113 * \param[in] scale_factor Factor to scale frequencies by before computing E
114 * \return The internal energy (J/mol K)
116 double calcVibrationalInternalEnergy(gmx::ArrayRef<const real> eigval,
117 real temperature,
118 gmx_bool linear,
119 real scale_factor);
121 /*! \brief Compute entropy using Schlitter formula
123 * Computes entropy for a molecule / molecular system using the
124 * algorithm due to Schlitter (Chem. Phys. Lett. 215 (1993)
125 * 617-621).
126 * The input should be eigenvalues from a covariance analysis,
127 * the units of the eigenvalues are those of energy.
129 * \param[in] eigval The eigenvalues
130 * \param[in] temperature Temperature (K)
131 * \param[in] linear True if this is a linear molecule (typically a diatomic molecule).
132 * \return the entropy (J/mol K)
134 double calcSchlitterEntropy(gmx::ArrayRef<const real> eigval,
135 real temperature,
136 gmx_bool linear);
138 /*! \brief Compute entropy using Quasi-Harmonic formula
140 * Computes entropy for a molecule / molecular system using the
141 * Quasi-harmonic algorithm (Macromolecules 1984, 17, 1370).
142 * The input should be eigenvalues from a normal mode analysis.
143 * In both cases the units of the eigenvalues are those of energy.
145 * \param[in] eigval The eigenvalues
146 * \param[in] temperature Temperature (K)
147 * \param[in] linear True if this is a linear molecule (typically a diatomic molecule).
148 * \param[in] scale_factor Factor to scale frequencies by before computing S0
149 * \return the entropy (J/mol K)
151 double calcQuasiHarmonicEntropy(gmx::ArrayRef<const real> eigval,
152 real temperature,
153 gmx_bool linear,
154 real scale_factor);
156 #endif