2 * This file is part of the GROMACS molecular simulation package.
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.
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 * Tests for entropy code
39 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
43 #include "gromacs/gmxana/thermochemistry.h"
44 #include "gromacs/utility/arrayref.h"
46 #include "testutils/refdata.h"
47 #include "testutils/testasserts.h"
55 class Entropy
: public ::testing::Test
58 test::TestReferenceData refData_
;
59 test::TestReferenceChecker checker_
;
60 Entropy() : checker_(refData_
.rootChecker()) {}
63 void runSchlitter(real temperature
, gmx_bool bLinear
)
65 std::vector
<double> ev
= { 0.00197861, 0.000389439, 0.000316043, 0.000150392, 0.000110254,
66 8.99659e-05, 8.06572e-05, 5.14339e-05, 4.34268e-05, 2.16063e-05,
67 1.65182e-05, 1.3965e-05, 1.01937e-05, 5.83113e-06, 4.59067e-06,
68 3.39577e-06, 2.14837e-06, 1.20468e-06, 9.60817e-18, 1.48941e-19,
69 1.13939e-19, 5.02332e-20, -6.90708e-21, -1.91354e-18 };
70 std::vector
<real
> eigenvalue
;
71 std::copy(ev
.begin(), ev
.end(), std::back_inserter(eigenvalue
));
73 real S
= calcSchlitterEntropy(eigenvalue
, temperature
, bLinear
);
74 checker_
.setDefaultTolerance(test::relativeToleranceAsFloatingPoint(1, 1e-7));
75 checker_
.checkReal(S
, "entropy");
78 void runQuasiHarmonic(real temperature
, gmx_bool bLinear
)
80 std::vector
<double> ev
= { -31.403, -7.73169, -3.80315, -2.15659e-06, -1.70991e-07, 236,
81 4609.83, 6718.07, 6966.27, 8587.85, 10736.3, 13543.7,
82 17721.3, 22868, 35517.8, 44118.1, 75827.9, 106277,
83 115132, 120782, 445118, 451149, 481058, 484576 };
84 std::vector
<real
> eigenvalue
;
85 std::copy(ev
.begin(), ev
.end(), std::back_inserter(eigenvalue
));
87 real S
= calcQuasiHarmonicEntropy(eigenvalue
, temperature
, bLinear
, 1.0);
89 checker_
.setDefaultTolerance(test::relativeToleranceAsFloatingPoint(1, 1e-7));
90 checker_
.checkReal(S
, "entropy");
94 TEST_F(Entropy
, Schlitter_300_NoLinear
)
96 runSchlitter(300.0, FALSE
);
99 TEST_F(Entropy
, Schlitter_300_Linear
)
101 runSchlitter(300.0, TRUE
);
104 TEST_F(Entropy
, QuasiHarmonic_300_NoLinear
)
106 runQuasiHarmonic(300.0, FALSE
);
109 TEST_F(Entropy
, QuasiHarmonic_200_NoLinear
)
111 runQuasiHarmonic(200.0, FALSE
);
114 TEST_F(Entropy
, QuasiHarmonic_200_Linear
)
116 runQuasiHarmonic(200.0, TRUE
);