Remove unnecessary includes of arrayref.h
[alexxy/gromacs.git] / src / gromacs / gmxana / thermochemistry.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2017,2018,2019,2020, 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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 /*! \internal \file
36  * \brief
37  * Code for computing entropy and heat capacity from eigenvalues
38  *
39  * \author David van der Spoel <david.vanderspoel@icm.uu.se>
40  */
41 #ifndef GMXANA_THERMOCHEMISTRY_H
42 #define GMXANA_THERMOCHEMISTRY_H
43
44 #include "gromacs/math/units.h"
45 #include "gromacs/math/vec.h"
46 #include "gromacs/utility/basedefinitions.h"
47
48 namespace gmx
49 {
50 template<typename>
51 class ArrayRef;
52 }
53
54 /*! \brief Compute zero point energy from an array of eigenvalues.
55  *
56  * This routine first converts the eigenvalues from a normal mode
57  * analysis to frequencies and then computes the zero point energy.
58  *
59  * \param[in] eigval       The eigenvalues
60  * \param[in] scale_factor Factor to scale frequencies by before computing cv
61  * \return The zero point energy (kJ/mol)
62  */
63 double calcZeroPointEnergy(gmx::ArrayRef<const real> eigval, real scale_factor);
64
65 /*! \brief Compute heat capacity due to vibrational motion
66  *
67  * \param[in] eigval       The eigenvalues
68  * \param[in] temperature  Temperature (K)
69  * \param[in] linear       TRUE if this is a linear molecule
70  * \param[in] scale_factor Factor to scale frequencies by before computing cv
71  * \return The heat capacity at constant volume (J/mol K)
72  */
73 double calcVibrationalHeatCapacity(gmx::ArrayRef<const real> eigval,
74                                    real                      temperature,
75                                    gmx_bool                  linear,
76                                    real                      scale_factor);
77
78 /*! \brief Compute entropy due to translational motion
79  *
80  * Following the equations in J. W. Ochterski,
81  * Thermochemistry in Gaussian, Gaussian, Inc., 2000
82  * Pitssburg PA
83  *
84  * \param[in] mass         Molecular mass (Dalton)
85  * \param[in] temperature  Temperature (K)
86  * \param[in] pressure     Pressure (bar) at which to compute
87  * \returns The translational entropy (J/mol K)
88  */
89 double calcTranslationalEntropy(real mass, real temperature, real pressure);
90
91 /*! \brief Compute entropy due to rotational motion
92  *
93  * Following the equations in J. W. Ochterski,
94  * Thermochemistry in Gaussian, Gaussian, Inc., 2000
95  * Pitssburg PA
96  *
97  * \param[in] temperature  Temperature (K)
98  * \param[in] natom        Number of atoms
99  * \param[in] linear       TRUE if this is a linear molecule
100  * \param[in] theta        The principal moments of inertia (unit of Energy)
101  * \param[in] sigma_r      Symmetry factor, should be >= 1
102  * \returns The rotational entropy (J/mol K)
103  */
104 double calcRotationalEntropy(real temperature, int natom, gmx_bool linear, const rvec theta, real sigma_r);
105
106 /*! \brief Compute internal energy due to vibrational motion
107  *
108  * \param[in] eigval       The eigenvalues
109  * \param[in] temperature  Temperature (K)
110  * \param[in] linear       TRUE if this is a linear molecule
111  * \param[in] scale_factor Factor to scale frequencies by before computing E
112  * \return The internal energy (J/mol K)
113  */
114 double calcVibrationalInternalEnergy(gmx::ArrayRef<const real> eigval,
115                                      real                      temperature,
116                                      gmx_bool                  linear,
117                                      real                      scale_factor);
118
119 /*! \brief Compute entropy using Schlitter formula
120  *
121  * Computes entropy for a molecule / molecular system using the
122  * algorithm due to Schlitter (Chem. Phys. Lett. 215 (1993)
123  * 617-621).
124  * The input should be eigenvalues from a covariance analysis,
125  * the units of the eigenvalues are those of energy.
126  *
127  * \param[in] eigval       The eigenvalues
128  * \param[in] temperature  Temperature (K)
129  * \param[in] linear       True if this is a linear molecule (typically a diatomic molecule).
130  * \return the entropy (J/mol K)
131  */
132 double calcSchlitterEntropy(gmx::ArrayRef<const real> eigval, real temperature, gmx_bool linear);
133
134 /*! \brief Compute entropy using Quasi-Harmonic formula
135  *
136  * Computes entropy for a molecule / molecular system using the
137  * Quasi-harmonic algorithm (Macromolecules 1984, 17, 1370).
138  * The input should be eigenvalues from a normal mode analysis.
139  * In both cases the units of the eigenvalues are those of energy.
140  *
141  * \param[in] eigval       The eigenvalues
142  * \param[in] temperature  Temperature (K)
143  * \param[in] linear       True if this is a linear molecule (typically a diatomic molecule).
144  * \param[in] scale_factor Factor to scale frequencies by before computing S0
145  * \return the entropy (J/mol K)
146  */
147 double calcQuasiHarmonicEntropy(gmx::ArrayRef<const real> eigval, real temperature, gmx_bool linear, real scale_factor);
148
149 #endif