#include "gromacs/math/units.h"
#include "gromacs/math/vec.h"
#include "gromacs/utility/arrayref.h"
#include "gromacs/utility/basedefinitions.h"
Code for computing entropy and heat capacity from eigenvalues.
- Author
- David van der Spoel david.nosp@m..van.nosp@m.dersp.nosp@m.oel@.nosp@m.icm.u.nosp@m.u.se
|
double | calcZeroPointEnergy (gmx::ArrayRef< const real > eigval, real scale_factor) |
| Compute zero point energy from an array of eigenvalues. More...
|
|
double | calcVibrationalHeatCapacity (gmx::ArrayRef< const real > eigval, real temperature, gmx_bool linear, real scale_factor) |
| Compute heat capacity due to vibrational motion. More...
|
|
double | calcTranslationalEntropy (real mass, real temperature, real pressure) |
| Compute entropy due to translational motion. More...
|
|
double | calcRotationalEntropy (real temperature, int natom, gmx_bool linear, const rvec theta, real sigma_r) |
| Compute entropy due to rotational motion. More...
|
|
double | calcVibrationalInternalEnergy (gmx::ArrayRef< const real > eigval, real temperature, gmx_bool linear, real scale_factor) |
| Compute internal energy due to vibrational motion. More...
|
|
double | calcSchlitterEntropy (gmx::ArrayRef< const real > eigval, real temperature, gmx_bool linear) |
| Compute entropy using Schlitter formula. More...
|
|
double | calcQuasiHarmonicEntropy (gmx::ArrayRef< const real > eigval, real temperature, gmx_bool linear, real scale_factor) |
| Compute entropy using Quasi-Harmonic formula. More...
|
|
Compute entropy using Quasi-Harmonic formula.
Computes entropy for a molecule / molecular system using the Quasi-harmonic algorithm (Macromolecules 1984, 17, 1370). The input should be eigenvalues from a normal mode analysis. In both cases the units of the eigenvalues are those of energy.
- Parameters
-
[in] | eigval | The eigenvalues |
[in] | temperature | Temperature (K) |
[in] | linear | True if this is a linear molecule (typically a diatomic molecule). |
[in] | scale_factor | Factor to scale frequencies by before computing S0 |
- Returns
- the entropy (J/mol K)
double calcRotationalEntropy |
( |
real |
temperature, |
|
|
int |
natom, |
|
|
gmx_bool |
linear, |
|
|
const rvec |
theta, |
|
|
real |
sigma_r |
|
) |
| |
Compute entropy due to rotational motion.
Following the equations in J. W. Ochterski, Thermochemistry in Gaussian, Gaussian, Inc., 2000 Pitssburg PA
- Parameters
-
[in] | temperature | Temperature (K) |
[in] | natom | Number of atoms |
[in] | linear | TRUE if this is a linear molecule |
[in] | theta | The principal moments of inertia (unit of Energy) |
[in] | sigma_r | Symmetry factor, should be >= 1 |
- Returns
- The rotational entropy (J/mol K)
Compute entropy using Schlitter formula.
Computes entropy for a molecule / molecular system using the algorithm due to Schlitter (Chem. Phys. Lett. 215 (1993) 617-621). The input should be eigenvalues from a covariance analysis, the units of the eigenvalues are those of energy.
- Parameters
-
[in] | eigval | The eigenvalues |
[in] | temperature | Temperature (K) |
[in] | linear | True if this is a linear molecule (typically a diatomic molecule). |
- Returns
- the entropy (J/mol K)
double calcTranslationalEntropy |
( |
real |
mass, |
|
|
real |
temperature, |
|
|
real |
pressure |
|
) |
| |
Compute entropy due to translational motion.
Following the equations in J. W. Ochterski, Thermochemistry in Gaussian, Gaussian, Inc., 2000 Pitssburg PA
- Parameters
-
[in] | mass | Molecular mass (Dalton) |
[in] | temperature | Temperature (K) |
[in] | pressure | Pressure (bar) at which to compute |
- Returns
- The translational entropy (J/mol K)
Compute heat capacity due to vibrational motion.
- Parameters
-
[in] | eigval | The eigenvalues |
[in] | temperature | Temperature (K) |
[in] | linear | TRUE if this is a linear molecule |
[in] | scale_factor | Factor to scale frequencies by before computing cv |
- Returns
- The heat capacity at constant volume (J/mol K)
Compute internal energy due to vibrational motion.
- Parameters
-
[in] | eigval | The eigenvalues |
[in] | temperature | Temperature (K) |
[in] | linear | TRUE if this is a linear molecule |
[in] | scale_factor | Factor to scale frequencies by before computing E |
- Returns
- The internal energy (J/mol K)
Compute zero point energy from an array of eigenvalues.
This routine first converts the eigenvalues from a normal mode analysis to frequencies and then computes the zero point energy.
- Parameters
-
[in] | eigval | The eigenvalues |
[in] | scale_factor | Factor to scale frequencies by before computing cv |
- Returns
- The zero point energy (kJ/mol)